LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
CalWireAna_module.cc
Go to the documentation of this file.
1 //
3 // CalWireAna class designed to make histos.
4 //
5 //
6 //
8 
10 
11 // C++ includes
12 #include <algorithm>
13 #include <string>
14 #include <vector>
15 
16 // Framework includes
20 #include "fhiclcpp/ParameterSet.h"
25 
26 // LArSoft includes
31 #include "larevt/CalibrationDBI/Interface/ChannelStatusService.h"
32 #include "larevt/CalibrationDBI/Interface/ChannelStatusProvider.h"
34 
35 // ROOT includes
36 #include <TMath.h>
37 #include <TH1.h>
38 #include <TH2.h>
39 #include <TGraph.h>
40 
41 class TH2F;
42 class TH1F;
43 
44 namespace geo { class Geometry; }
45 
46 
47 namespace caldata {
48 
50  class CalWireAna : public art::EDAnalyzer {
51 
52  public:
53 
54  explicit CalWireAna(fhicl::ParameterSet const& pset);
55  virtual ~CalWireAna();
56 
58  void analyze (const art::Event& evt);
59  void beginJob();
60  void endJob();
61 
62  private:
63 
64  std::string fCalWireModuleLabel;
65  std::string fDetSimModuleLabel; //< name of module that produced the digits
66 
67  TH1F* fDiffsR;
68  TH1F* fDiffsW;
69  TH1F* fDiffsRW;
70  TH1F* fDiffsRWgph;
71  TH1F* fMin;
72  TH1F* fMax;
73  TH1F* fIR;
74  TH1F* fCR;
75  TH1F* fIW;
76  TH1F* fCW;
77  TH1F* fRawIndPeak; //Raw Induction Peak Values
78  TH1F* fRawColPeak; //Raw Collection Peak Values
79  TH1F* fCalIndPeak; //Calibrated Induction Peak Values
80  TH1F* fCalColPeak; //Calibrated Collection Peak Values
81  TH1F* fNoiseHist; //Noise Frequency Spectrum
82  TH1F* fNoiseRMS; //Noise RMS values
83  TH2F* fWireSig;
84  TH2F* fRawSig;
87  TH2F* fDiffsRWvsR;
89  TH2F* fWindow;
90 
91 
92  }; // class CalWireAna
93 
94 } // End caldata namespace.
95 
96 
97 namespace caldata{
98 
99  //-------------------------------------------------
100  CalWireAna::CalWireAna(fhicl::ParameterSet const& pset)
101  : EDAnalyzer (pset)
102  , fCalWireModuleLabel(pset.get< std::string >("CalWireModuleLabel"))
103  , fDetSimModuleLabel (pset.get< std::string >("DetSimModuleLabel"))
104  {
105 
106  }
107 
108  //-------------------------------------------------
110  {
111  }
112 
113  //-------------------------------------------------
115  {
116  // get access to the TFile service
118 
119  fDiffsR = tfs->make<TH1F>("One timestamp diffs in RawDigits", ";#Delta ADC;", 40, -19.5, 20.5);
120  fDiffsW = tfs->make<TH1F>("One timestamp diffs in Wires", ";#Delta ADC;", 20, -9.5, 10.5);
121  fDiffsRW = tfs->make<TH1F>("Same timestamp diffs in RD-Wires", ";#Delta ADC;", 20, -9.5, 10.5);
122  fDiffsRWvsR = tfs->make<TH2F>("Same timestamp diffs in RD-Wires vs R", ";#Delta ADC;", 481,-0.5,480.5, 20, -9.5, 10.5);
123  fDiffsRWgph = tfs->make<TH1F>("Same timestamp diffs in RD-Wires gph", ";#Delta ADC;", 20, -9.5, 10.5);
124  fDiffsRWvsRgph = tfs->make<TH2F>("Same timestamp diffs in RD-Wires vs R gph", ";#Delta ADC;", 481,-0.5,480.5, 20, -9.5, 10.5);
125  fRawSig = tfs->make<TH2F>("One event, one channel Raw", "timestamp", 481,-0.5,480.5, 21,-0.5, 20.5);
126  fWireSig = tfs->make<TH2F>("One event, one channel Wire", "timestamp", 481,-0.5,480.5, 21, -0.5, 20.5);
127 
128  fRD_WireMeanDiff2D = tfs->make<TH2F>("Mean (Raw-CALD)-over-Raw in Window gph","Wire number",481,-0.05,480.5, 40, -1., 1.);
129  fRD_WireRMSDiff2D = tfs->make<TH2F>("RMS (Raw-CALD)-over-Raw in Window gph","Wire number",481,-0.05,480.5, 10, 0., 2.);
130 
131  fWindow = tfs->make<TH2F>("tmax-tmin vs indMax", "ticks", 200, 0, 2000, 20, -2.5, 60.5);
132  fMin = tfs->make<TH1F>("Value of min", "ticks", 21, -20.5, 0.5);
133  fMax = tfs->make<TH1F>("Value of max", "ticks", 21, 0.5, 20.5);
134  fRawIndPeak = tfs->make<TH1F>("indPeakRaw", ";Induction Peaks Raw;",40,5,45);
135  fRawColPeak = tfs->make<TH1F>("colPeakRaw", ";Collection Peaks Raw;",40,5,45);
136  fCalIndPeak = tfs->make<TH1F>("indPeakCal", ";Induction Peaks Calibrated;",40,5,45);
137  fCalColPeak = tfs->make<TH1F>("colPeakCal", ";Collection Peaks Calibrated;",40,5,45);
138 
139  fIR = tfs->make<TH1F>("Raw Ind signal","time ticks",4096,0.0,4096.);
140  fCR = tfs->make<TH1F>("Raw Coll signal","time ticks",4096,0.0,4096.);
141  fIW = tfs->make<TH1F>("Wire Ind signal","time ticks",4096,0.0,4096.);
142  fCW = tfs->make<TH1F>("Wire Coll signal","time ticks",4096,0.0,4096.);
143  fNoiseHist = tfs->make<TH1F>("Noise Histogram","FFT Bins",2049,0,2049);
144  fNoiseRMS = tfs->make<TH1F>("Noise RMS","RMS",25,0,2.0);
145 
146  return;
147 
148  }
149 
150  //-------------------------------------------------
152  {
153  }
154 
155  //-------------------------------------------------
157  {
158 
159  // loop over the raw digits and get the adc vector for each, then compress it and uncompress it
160 
161  lariov::ChannelStatusProvider const& channelStatus
164  evt.getByLabel(fDetSimModuleLabel,rdHandle);
166  evt.getByLabel(fCalWireModuleLabel,wHandle);
167 
168  // get the raw::RawDigit associated by fCalWireModuleLabel to wires in wHandle;
169  // RawDigitsFromWire.at(index) will be a art::Ptr<raw::RawDigit>
170  art::FindOneP<raw::RawDigit> RawDigitsFromWire(wHandle, evt, fCalWireModuleLabel);
171 
173  for(unsigned int i = 0; i < wHandle->size(); ++i){
174  art::Ptr<recob::Wire> w(wHandle,i);
175  wvec.push_back(w);
176  }
178  for(unsigned int i = 0; i < rdHandle->size(); ++i){
179  art::Ptr<raw::RawDigit> r(rdHandle,i);
180  rdvec.push_back(r);
181  }
182 
185  double pedestal = rdvec[0]->GetPedestal();
186  double threshold = 9.0;
187  double signalSize = rdvec[0]->Samples();
188  uint32_t indChan0=64;
189  uint32_t indChan1=110;
190  uint32_t colChan0=312;
191  uint32_t colChan1=354;
192  std::vector<double> ir(fft->FFTSize()),iw(fft->FFTSize()),cr(fft->FFTSize()),cw(fft->FFTSize());
193  ir[signalSize]=iw[signalSize]=cr[signalSize]=cw[signalSize]=1.0;
195  for(unsigned int rd = 0; rd < rdvec.size(); ++rd){
196  // Find corresponding wire.
197  std::vector<double> signal(fft->FFTSize());
198  for(unsigned int wd = 0; wd < wvec.size(); ++wd){
199  // if (wvec[wd]->RawDigit() == rdvec[rd]){
200  if (RawDigitsFromWire.at(wd) == rdvec[rd]){
201  std::vector<float> wirSig = wvec[wd]->Signal();
202  if(wirSig.size() > signal.size()) {
203  LOG_DEBUG("CalWireAna")<<"Incompatible vector size "<<wirSig.size()
204  <<" "<<signal.size();
205  return;
206  }
207  for(unsigned int ii = 0; ii < wirSig.size(); ++ii) {
208  signal[ii] = wirSig[ii];
209  }
210  break;
211  }
212  if (wd == (wvec.size()-1) ){
213  LOG_DEBUG("CalWireAna") << "caldata::CalWireAna:Big problem! No matching Wire for RawDigit. Bailing." << rd;
214  return;
215  }
216  }
217 
218  std::vector<double> adc(fft->FFTSize());
219 
220  for(unsigned int t = 1; t < rdvec[rd]->Samples(); ++t){
221  fDiffsR->Fill(rdvec[rd]->ADC(t) - rdvec[rd]->ADC(t-1));
222  adc[t-1]=rdvec[rd]->ADC(t-1);
223  fRawSig->Fill(rd,rdvec[rd]->ADC(t));
224  }
225  //get the last one for the adc vector
226  adc[rdvec[rd]->Samples()-1] = rdvec[rd]->ADC(rdvec[rd]->Samples()-1);
227  if(!channelStatus.IsBad(rdvec[rd]->Channel()) &&
228  (*max_element(adc.begin(),adc.end()) < pedestal+threshold &&
229  *min_element(adc.begin(),adc.end()) >pedestal -threshold)) {
230  double sum=0;
231  for(int i = 0; i < signalSize; i++) sum+=pow(adc[i]-pedestal,2.0);
232  fNoiseRMS->Fill(TMath::Sqrt(sum/(double)signalSize));
233  std::vector<double> temp(fft->FFTSize());
234  std::vector<TComplex> fTemp(fft->FFTSize()/2+1);
235  for(int i = 0; i < signalSize; i++) temp[i]=(adc[i]-pedestal)*sin(TMath::Pi()*(double)i/signalSize);
236  fft->DoFFT(temp,fTemp);
237  for(int i = 0; i < fft->FFTSize()/2+1; i++) fNoiseHist->Fill(i,fTemp[i].Rho());
238  }
239  if(geom->SignalType(rdvec[rd]->Channel()) == geo::kInduction &&
240  rdvec[rd]->Channel() > indChan0 &&
241  rdvec[rd]->Channel() < indChan1){
242  fft->AlignedSum(ir,adc);
243  fft->AlignedSum(iw,signal);
244  }
245  if(geom->SignalType(rdvec[rd]->Channel()) == geo::kCollection &&
246  rdvec[rd]->Channel() > colChan0 &&
247  rdvec[rd]->Channel() < colChan1) {
248  fft->AlignedSum(cr,adc);
249  fft->AlignedSum(cw,signal);
250  }
251  if(geom->SignalType(rdvec[rd]->Channel()) == geo::kInduction) {
252  if(*max_element(adc.begin(),adc.end()) > pedestal+threshold)
253  fRawIndPeak->Fill(*max_element(adc.begin(),adc.end()));
254  if(*max_element(signal.begin(),signal.end()) > pedestal+threshold)
255  fCalIndPeak->Fill(*max_element(signal.begin(),signal.end()));
256  }
257  if(geom->SignalType(rdvec[rd]->Channel()) == geo::kCollection) {
258  if(*max_element(adc.begin(),adc.end()) > pedestal +threshold)
259  fRawColPeak->Fill(*max_element(adc.begin(),adc.end()));
260  if(*max_element(signal.begin(),signal.end()) > pedestal+threshold)
261  fCalColPeak->Fill(*max_element(signal.begin(),signal.end()));
262  }
263 
264  int window = 8;
265  static unsigned int pulseHeight = 5;
266  unsigned int tmin = 1;
267  unsigned int tmax = 1;
268  int indMax = TMath::LocMax(signalSize,&adc[0]);
269  double sigMin = 0.0;
270  double sigMax = TMath::MaxElement(signalSize,&adc[0]);
271  if(geom->SignalType(rdvec[rd]->Channel()) == geo::kInduction && sigMax>=pulseHeight) {
272  int indMin = TMath::LocMin(signalSize,&adc[0]);
273  sigMin = TMath::MinElement(signalSize,&adc[0]);
274  tmin = std::max(indMax-window,0);
275  tmax = std::min(indMin+window,(int)signalSize-1);
276  LOG_DEBUG("CalWireAna") << "Induction channel, indMin,tmin,tmax "
277  << rd<< " " << indMin<< " " << tmin << " " << tmax;
278  }
279  else if (sigMax>=pulseHeight){
280  tmin = std::max(indMax-window,0);
281  tmax = std::min(indMax+window,(int)signalSize-1);
282  LOG_DEBUG("CalWireAna") << "Collection channel, tmin,tmax "<< rd<< " " << tmin << " " << tmax;
283  }
284 
285  fMin->Fill(sigMin);
286  fMax->Fill(sigMax);
287  fWindow->Fill(indMax, tmax - tmin);
288 
289  std::vector<double> winDiffs;
290  int cnt=0;
291  // for(unsigned int t = tmin; t < tmax; ++t)
292  static unsigned int tRawLead = 0;
293  for(unsigned int t = 1; t < signalSize; ++t){
294  fDiffsW->Fill(signal[t]-signal[t-1]);
295  fWireSig->Fill(rd,signal[t]);
296 
297  if (t>=tmin && t<=tmax && tmax>=pulseHeight && (t+tRawLead)<signalSize){
298  cnt++;
299  winDiffs.push_back((adc[t+tRawLead]-signal[t])/adc[t+tRawLead]);
300  fDiffsRW->Fill(adc[t+tRawLead]-signal[t]);
301  fDiffsRWvsR->Fill(rd,adc[t+tRawLead]-signal[t]);
302  if (sigMax >= pulseHeight){
303  fDiffsRWgph->Fill(adc[t+tRawLead]-signal[t]);
304  fDiffsRWvsRgph->Fill(rd,adc[t+tRawLead]-signal[t]);
305  }
306  }
307  }
308 
309  LOG_DEBUG("CalWireAna") << "on channel " << rdvec[rd]->Channel();
310  // TMath::Mean with iterators doesn't work. EC,23-Sep-2010.
311  double tmp = TMath::Mean(winDiffs.size(),&winDiffs[0]);
312  double tmp2 = TMath::RMS(winDiffs.size(),&winDiffs[0]);
313  double tmp3=0;
314  for (int ii=0; ii<rdvec[rd]->Samples(); ii++) tmp3+=rdvec[rd]->ADC(ii);
315  for(int i = 0; i < fft->FFTSize(); i++) {
316  fIR->Fill(i,ir[i]);
317  fIW->Fill(i,iw[i]);
318  fCR->Fill(i,cr[i]);
319  fCW->Fill(i,cw[i]);
320  }
321  fRD_WireMeanDiff2D->Fill(rd,tmp);
322  fRD_WireRMSDiff2D->Fill(rd,tmp2);
323 
324  }//end loop over rawDigits
325 
326  return;
327  }//end analyze method
328 
330 
331 }//end namespace
TH2F * fRD_WireMeanDiff2D
histogram of difference between original tdc value and compressesed value vs original value ...
void analyze(const art::Event &evt)
read/write access to event
STL namespace.
Definition of basic raw digits.
Float_t tmp
Definition: plot.C:37
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:97
SigType_t SignalType(geo::PlaneID const &pid) const
Returns the type of signal on the channels of specified TPC plane.
creation of calibrated signals on wires
TH1F * fDiffsR
histogram of Raw tdc to tdc differences
int FFTSize() const
Definition: LArFFT.h:69
TH2F * fRD_WireRMSDiff2D
histogram of difference between original tdc value and compressesed value vs original value ...
Int_t max
Definition: plot.C:27
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
void beginJob()
Definition: Breakpoints.cc:14
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:441
Signal from induction planes.
Definition: geo_types.h:92
size_type size() const
Definition: PtrVector.h:308
std::string fCalWireModuleLabel
name of module that produced the wires
T * make(ARGS...args) const
Encapsulate the construction of a single detector plane.
void AlignedSum(std::vector< T > &input, std::vector< T > &output, bool add=true)
Definition: LArFFT.h:243
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
Int_t min
Definition: plot.C:26
std::string fDetSimModuleLabel
TH1F * fDiffsW
histogram of Wire (post-deconvoution) tdc to tdc differences
#define LOG_DEBUG(id)
Declaration of basic channel signal object.
Namespace collecting geometry-related classes utilities.
Float_t w
Definition: plot.C:23
Base class for creation of raw signals on wires.
Definition: fwd.h:25
art framework interface to geometry description
Signal from collection planes.
Definition: geo_types.h:93