Hall C ROOT/C++ Analyzer (hcana)
THcCherenkov.cxx
Go to the documentation of this file.
1 
8 #include "THcCherenkov.h"
9 #include "THcHodoscope.h"
10 #include "TClonesArray.h"
11 #include "THcSignalHit.h"
12 #include "THaEvData.h"
13 #include "THaDetMap.h"
14 #include "THcDetectorMap.h"
15 #include "THcGlobals.h"
16 #include "THaCutList.h"
17 #include "THcParmList.h"
18 #include "THcHitList.h"
19 #include "THaApparatus.h"
20 #include "VarDef.h"
21 #include "VarType.h"
22 #include "THaTrack.h"
23 #include "TClonesArray.h"
24 #include "TMath.h"
25 #include "THaTrackProj.h"
26 #include "THcHallCSpectrometer.h"
27 
28 #include <algorithm>
29 #include <cstring>
30 #include <cstdio>
31 #include <cstdlib>
32 #include <iostream>
33 #include <string>
34 #include <iomanip>
35 
36 
37 using namespace std;
38 
39 using std::cout;
40 using std::cin;
41 using std::endl;
42 using std::setw;
43 using std::setprecision;
44 
45 //_____________________________________________________________________________
46 THcCherenkov::THcCherenkov( const char* name, const char* description,
47  THaApparatus* apparatus ) :
48  THaNonTrackingDetector(name,description,apparatus)
49 {
50  // Normal constructor with name and description
51  frAdcPedRaw = new TClonesArray("THcSignalHit", MaxNumCerPmt*MaxNumAdcPulse);
52  frAdcPulseIntRaw = new TClonesArray("THcSignalHit", MaxNumCerPmt*MaxNumAdcPulse);
53  frAdcPulseAmpRaw = new TClonesArray("THcSignalHit", MaxNumCerPmt*MaxNumAdcPulse);
54  frAdcPulseTimeRaw = new TClonesArray("THcSignalHit", MaxNumCerPmt*MaxNumAdcPulse);
55  frAdcPed = new TClonesArray("THcSignalHit", MaxNumCerPmt*MaxNumAdcPulse);
56  frAdcPulseInt = new TClonesArray("THcSignalHit", MaxNumCerPmt*MaxNumAdcPulse);
57  frAdcPulseAmp = new TClonesArray("THcSignalHit", MaxNumCerPmt*MaxNumAdcPulse);
58  frAdcPulseTime = new TClonesArray("THcSignalHit", MaxNumCerPmt*MaxNumAdcPulse);
59  fAdcErrorFlag = new TClonesArray("THcSignalHit", MaxNumCerPmt*MaxNumAdcPulse);
60 
61  fNumAdcHits = vector<Int_t> (MaxNumCerPmt, 0.0);
62  fNumGoodAdcHits = vector<Int_t> (MaxNumCerPmt, 0.0);
63  fNumTracksMatched = vector<Int_t> (MaxNumCerPmt, 0.0);
64  fNumTracksFired = vector<Int_t> (MaxNumCerPmt, 0.0);
65  fNpe = vector<Double_t> (MaxNumCerPmt, 0.0);
66  fGoodAdcPed = vector<Double_t> (MaxNumCerPmt, 0.0);
67  fGoodAdcMult = vector<Double_t> (MaxNumCerPmt, 0.0);
68  fGoodAdcHitUsed = vector<Double_t> (MaxNumCerPmt, 0.0);
69  fGoodAdcPulseInt = vector<Double_t> (MaxNumCerPmt, 0.0);
70  fGoodAdcPulseIntRaw = vector<Double_t> (MaxNumCerPmt, 0.0);
71  fGoodAdcPulseAmp = vector<Double_t> (MaxNumCerPmt, 0.0);
72  fGoodAdcPulseTime = vector<Double_t> (MaxNumCerPmt, 0.0);
73  fGoodAdcTdcDiffTime = vector<Double_t> (MaxNumCerPmt, 0.0);
74 
75  InitArrays();
76 }
77 
78 //_____________________________________________________________________________
80  THaNonTrackingDetector()
81 {
82  // Constructor
83  frAdcPedRaw = NULL;
84  frAdcPulseIntRaw = NULL;
85  frAdcPulseAmpRaw = NULL;
86  frAdcPulseTimeRaw = NULL;
87  frAdcPed = NULL;
88  frAdcPulseInt = NULL;
89  frAdcPulseAmp = NULL;
90  frAdcPulseTime = NULL;
91  fAdcErrorFlag = NULL;
92 
93  InitArrays();
94 }
95 
96 //_____________________________________________________________________________
98 {
99  // Destructor
100  delete frAdcPedRaw; frAdcPedRaw = NULL;
101  delete frAdcPulseIntRaw; frAdcPulseIntRaw = NULL;
102  delete frAdcPulseAmpRaw; frAdcPulseAmpRaw = NULL;
103  delete frAdcPulseTimeRaw; frAdcPulseTimeRaw = NULL;
104  delete frAdcPed; frAdcPed = NULL;
105  delete frAdcPulseInt; frAdcPulseInt = NULL;
106  delete frAdcPulseAmp; frAdcPulseAmp = NULL;
107  delete frAdcPulseTime; frAdcPulseTime = NULL;
108  delete fAdcErrorFlag; fAdcErrorFlag = NULL;
109 
110  DeleteArrays();
111 }
112 
113 //_____________________________________________________________________________
115 {
116  fGain = NULL;
117  fPedSum = NULL;
118  fPedSum2 = NULL;
119  fPedLimit = NULL;
120  fPedMean = NULL;
121  fPedCount = NULL;
122  fPed = NULL;
123  fThresh = NULL;
124 }
125 //_____________________________________________________________________________
127 {
128  delete [] fGain; fGain = NULL;
129 
130  // 6 Gev variables
131  delete [] fPedSum; fPedSum = NULL;
132  delete [] fPedSum2; fPedSum2 = NULL;
133  delete [] fPedLimit; fPedLimit = NULL;
134  delete [] fPedMean; fPedMean = NULL;
135  delete [] fPedCount; fPedCount = NULL;
136  delete [] fPed; fPed = NULL;
137  delete [] fThresh; fThresh = NULL;
138 
139  delete [] fPedDefault; fPedDefault = 0;
140  delete [] fAdcTimeWindowMin; fAdcTimeWindowMin = 0;
141  delete [] fAdcTimeWindowMax; fAdcTimeWindowMax = 0;
142  delete [] fRegionValue; fRegionValue = 0;
143 }
144 
145 //_____________________________________________________________________________
146 THaAnalysisObject::EStatus THcCherenkov::Init( const TDatime& date )
147 {
148  // cout << "THcCherenkov::Init for: " << GetName() << endl;
149 
150  string EngineDID = string(GetApparatus()->GetName()).substr(0, 1) + GetName();
151  std::transform(EngineDID.begin(), EngineDID.end(), EngineDID.begin(), ::toupper);
152  if(gHcDetectorMap->FillMap(fDetMap, EngineDID.c_str()) < 0) {
153  static const char* const here = "Init()";
154  Error(Here(here), "Error filling detectormap for %s.", EngineDID.c_str());
155  return kInitError;
156  }
157 
158  EStatus status;
159  if((status = THaNonTrackingDetector::Init( date )))
160  return fStatus=status;
161 
162  // Should probably put this in ReadDatabase as we will know the
163  // maximum number of hits after setting up the detector map
164  InitHitList(fDetMap, "THcCherenkovHit", fDetMap->GetTotNumChan()+1,
165  0, fADC_RefTimeCut);
166 
167  THcHallCSpectrometer *app=dynamic_cast<THcHallCSpectrometer*>(GetApparatus());
168  if( !app ||
169  !(fglHod = dynamic_cast<THcHodoscope*>(app->GetDetector("hod"))) ) {
170  static const char* const here = "ReadDatabase()";
171  Warning(Here(here),"Hodoscope \"%s\" not found. ","hod");
172  }
173 
174  fPresentP = 0;
175  THaVar* vpresent = gHaVars->Find(Form("%s.present",GetApparatus()->GetName()));
176  if(vpresent) {
177  fPresentP = (Bool_t *) vpresent->GetValuePointer();
178  }
179  return fStatus = kOK;
180 }
181 
182 //_____________________________________________________________________________
184 {
185  // This function is called by THaDetectorBase::Init() once at the beginning
186  // of the analysis.
187 
188  // cout << "THcCherenkov::ReadDatabase for: " << GetName() << endl; // Ahmed
189 
190  string prefix = string(GetApparatus()->GetName()).substr(0, 1) + GetName();
191  std::transform(prefix.begin(), prefix.end(), prefix.begin(), ::tolower);
192 
193  CreateMissReportParms(prefix.c_str());
194 
195  fNRegions = 4; // Defualt if not set in paramter file
196 
197  DBRequest list_1[] = {
198  {"_tot_pmts", &fNelem, kInt},
199  {"_num_regions", &fNRegions, kInt},
200  {0}
201  };
202 
203  gHcParms->LoadParmValues(list_1, prefix.c_str());
204 
205  Bool_t optional = true;
206 
207  cout << "Created Cherenkov detector " << GetApparatus()->GetName() << "."
208  << GetName() << " with " << fNelem << " PMTs" << endl;
209 
210  // 6 GeV pedestal paramters
211  fPedLimit = new Int_t[fNelem];
212  fGain = new Double_t[fNelem];
213  fPedMean = new Double_t[fNelem];
214  fAdcTimeWindowMin = new Double_t[fNelem];
215  fAdcTimeWindowMax= new Double_t[fNelem];
216  fPedDefault= new Int_t[fNelem];
217  // Region parameters
218  fRegionsValueMax = fNRegions * 8;
220  fAdcGoodElem = new Int_t[fNelem];
221  fAdcPulseAmpTest = new Double_t[fNelem];
222 
223  DBRequest list[]={
224  {"_ped_limit", fPedLimit, kInt, (UInt_t) fNelem, optional},
225  {"_adc_to_npe", fGain, kDouble, (UInt_t) fNelem},
226  {"_red_chi2_min", &fRedChi2Min, kDouble},
227  {"_red_chi2_max", &fRedChi2Max, kDouble},
228  {"_beta_min", &fBetaMin, kDouble},
229  {"_beta_max", &fBetaMax, kDouble},
230  {"_enorm_min", &fENormMin, kDouble},
231  {"_enorm_max", &fENormMax, kDouble},
232  {"_dp_min", &fDpMin, kDouble},
233  {"_dp_max", &fDpMax, kDouble},
234  {"_mirror_zpos", &fMirrorZPos, kDouble},
235  {"_npe_thresh", &fNpeThresh, kDouble},
236  {"_debug_adc", &fDebugAdc, kInt, 0, 1},
237  {"_adcTimeWindowMin", fAdcTimeWindowMin, kDouble,(UInt_t) fNelem,1},
238  {"_adcTimeWindowMax", fAdcTimeWindowMax, kDouble, (UInt_t) fNelem,1},
239  {"_PedDefault", fPedDefault, kInt, (UInt_t) fNelem,1},
240  {"_adc_tdc_offset", &fAdcTdcOffset, kDouble, 0, 1},
241  {"_region", &fRegionValue[0], kDouble, (UInt_t) fRegionsValueMax},
242  {"_adcrefcut", &fADC_RefTimeCut, kInt, 0, 1},
243  {0}
244  };
245  for (Int_t i=0;i<fNelem;i++) {
246  fAdcTimeWindowMin[i]=-1000.;
247  fAdcTimeWindowMax[i]=1000.;
248  fPedDefault[i]=0;
249  }
250  fDebugAdc = 0; // Set ADC debug parameter to false unless set in parameter file
251  fAdcTdcOffset = 0.0;
252  fADC_RefTimeCut = 0;
253 
254  gHcParms->LoadParmValues((DBRequest*)&list, prefix.c_str());
255 
256  // if (fDebugAdc) cout << "Cherenkov ADC Debug Flag Set To TRUE" << endl;
257 
258  fIsInit = true;
259 
260  // cout << "Track Matching Parameters for: " << GetName() << endl;
261  // for (Int_t iregion = 0; iregion < fNRegions; iregion++) {
262  // cout << "Region = " << iregion + 1 << endl;
263  // for (Int_t ivalue = 0; ivalue < 8; ivalue++)
264  // cout << fRegionValue[GetIndex(iregion, ivalue)] << " ";
265  // cout << endl;
266  // }
267 
268  // Create arrays to hold pedestal results
270 
271  return kOK;
272 }
273 
274 //_____________________________________________________________________________
276 {
277  // Initialize global variables for histogramming and tree
278  // cout << "THcCherenkov::DefineVariables called for: " << GetName() << endl;
279 
280  if( mode == kDefine && fIsSetup ) return kOK;
281  fIsSetup = ( mode == kDefine );
282 
283  // Register variables in global list
284 
285  if (fDebugAdc) {
286  RVarDef vars[] = {
287  {"numAdcHits", "Number of ADC Hits Per PMT", "fNumAdcHits"}, // Cherenkov occupancy
288  {"totNumAdcHits", "Total Number of ADC Hits", "fTotNumAdcHits"}, // Cherenkov multiplicity
289  {"adcPedRaw", "Raw ADC pedestals", "frAdcPedRaw.THcSignalHit.GetData()"},
290  {"adcPulseIntRaw", "Raw ADC pulse integrals", "frAdcPulseIntRaw.THcSignalHit.GetData()"},
291  {"adcPulseAmpRaw", "Raw ADC pulse amplitudes", "frAdcPulseAmpRaw.THcSignalHit.GetData()"},
292  {"adcPulseTimeRaw", "Raw ADC pulse times", "frAdcPulseTimeRaw.THcSignalHit.GetData()"},
293  {"adcPed", "ADC pedestals", "frAdcPed.THcSignalHit.GetData()"},
294  {"adcPulseInt", "ADC pulse integrals", "frAdcPulseInt.THcSignalHit.GetData()"},
295  {"adcPulseAmp", "ADC pulse amplitudes", "frAdcPulseAmp.THcSignalHit.GetData()"},
296  {"adcPulseTime", "ADC pulse times", "frAdcPulseTime.THcSignalHit.GetData()"},
297  { 0 }
298  };
299  DefineVarsFromList( vars, mode);
300  } //end debug statement
301 
302  RVarDef vars[] = {
303  {"adcCounter", "ADC counter numbers", "frAdcPulseIntRaw.THcSignalHit.GetPaddleNumber()"},
304  {"adcErrorFlag", "Error Flag for When FPGA Fails", "fAdcErrorFlag.THcSignalHit.GetData()"},
305 
306  {"numGoodAdcHits", "Number of Good ADC Hits Per PMT", "fNumGoodAdcHits"}, // Cherenkov occupancy
307  {"totNumGoodAdcHits", "Total Number of Good ADC Hits", "fTotNumGoodAdcHits"}, // Cherenkov multiplicity
308 
309  {"numTracksMatched", "Number of Tracks Matched Per Region", "fNumTracksMatched"},
310  {"numTracksFired", "Number of Tracks that Fired Per Region", "fNumTracksFired"},
311  {"totNumTracksMatched", "Total Number of Tracks Matched Per Region", "fTotNumTracksMatched"},
312  {"totNumTracksFired", "Total Number of Tracks that Fired", "fTotNumTracksFired"},
313 
314  {"xAtCer", "Track X at Cherenkov mirror", "fXAtCer"},
315  {"yAtCer", "Track Y at Cherenkov mirror", "fYAtCer"},
316 
317  {"npe", "Number of PEs", "fNpe"},
318  {"npeSum", "Total Number of PEs", "fNpeSum"},
319 
320  {"goodAdcPed", "Good ADC pedestals", "fGoodAdcPed"},
321  {"goodAdcMult", "Good ADC Multiplicity", "fGoodAdcMult"},
322  {"goodAdcHitUsed", "Good ADC Hit Used", "fGoodAdcHitUsed"},
323  {"goodAdcPulseInt", "Good ADC pulse integrals", "fGoodAdcPulseInt"},
324  {"goodAdcPulseIntRaw", "Good ADC raw pulse integrals", "fGoodAdcPulseIntRaw"},
325  {"goodAdcPulseAmp", "Good ADC pulse amplitudes", "fGoodAdcPulseAmp"},
326  {"goodAdcPulseTime", "Good ADC pulse times", "fGoodAdcPulseTime"},
327  {"goodAdcTdcDiffTime", "Good Hodo Start - ADC pulse times", "fGoodAdcTdcDiffTime"},
328  {"RefTime", "Raw ADC RefTime (chan) ", "fRefTime"}, // Raw reference time
329  { 0 }
330  };
331 
332  return DefineVarsFromList(vars, mode);
333 }
334 //_____________________________________________________________________________
335 inline
337 {
338  // Clear the hit lists
339  fNhits = 0;
340  fTotNumAdcHits = 0;
341  fTotNumGoodAdcHits = 0;
343  fTotNumTracksFired = 0;
344 
345  fXAtCer = 0.0;
346  fYAtCer = 0.0;
347 
348  fNpeSum = 0.0;
349  fRefTime=kBig;
350 
351  frAdcPedRaw->Clear();
355 
356  frAdcPed->Clear();
357  frAdcPulseInt->Clear();
358  frAdcPulseAmp->Clear();
360  fAdcErrorFlag->Clear();
361 
362  for (UInt_t ielem = 0; ielem < fNumAdcHits.size(); ielem++)
363  fNumAdcHits.at(ielem) = 0;
364  for (UInt_t ielem = 0; ielem < fNumGoodAdcHits.size(); ielem++)
365  fNumGoodAdcHits.at(ielem) = 0;
366  for (UInt_t ielem = 0; ielem < fNumTracksMatched.size(); ielem++)
367  fNumTracksMatched.at(ielem) = 0;
368  for (UInt_t ielem = 0; ielem < fNumTracksFired.size(); ielem++)
369  fNumTracksFired.at(ielem) = 0;
370  for (UInt_t ielem = 0; ielem < fGoodAdcPed.size(); ielem++) {
371  fGoodAdcPed.at(ielem) = 0.0;
372  fGoodAdcMult.at(ielem) = 0.0;
373  fGoodAdcHitUsed.at(ielem) = 0.0;
374  fGoodAdcPulseInt.at(ielem) = 0.0;
375  fGoodAdcPulseIntRaw.at(ielem) = 0.0;
376  fGoodAdcPulseAmp.at(ielem) = 0.0;
377  fGoodAdcPulseTime.at(ielem) = kBig;
378  fGoodAdcTdcDiffTime.at(ielem) = kBig;
379  fNpe.at(ielem) = 0.0;
380  }
381 
382 }
383 
384 //_____________________________________________________________________________
386 {
387  // Get the Hall C style hitlist (fRawHitList) for this event
388  Bool_t present = kTRUE; // Suppress reference time warnings
389  if(fPresentP) { // if this spectrometer not part of trigger
390  present = *fPresentP;
391  }
392  fNhits = DecodeToHitList(evdata, !present);
393 
394  //THcHallCSpectrometer *app = dynamic_cast<THcHallCSpectrometer*>(GetApparatus());
395  // cout << "Cerenkov Event num = " << evdata.GetEvNum() << " spec = " << app->GetName() << endl;
396 
397  if(gHaCuts->Result("Pedestal_event")) {
399  fAnalyzePedestals = 1; // Analyze pedestals first normal events
400  return(0);
401  }
402 
403  if(fAnalyzePedestals) {
405  fAnalyzePedestals = 0; // Don't analyze pedestals next event
406  }
407 
408  Int_t ihit = 0;
409  UInt_t nrAdcHits = 0;
410 
411  while(ihit < fNhits) {
412 
414  Int_t npmt = hit->fCounter;
415  THcRawAdcHit& rawAdcHit = hit->GetRawAdcHitPos();
416  if (rawAdcHit.GetNPulses() >0 && rawAdcHit.HasRefTime()) {
417  fRefTime=rawAdcHit.GetRefTime() ;
418  }
419  //if (rawAdcHit.GetNPulses()>0) cout << "Cer npmt = " << " ped = " << rawAdcHit.GetPed() << endl;
420  for (UInt_t thit = 0; thit < rawAdcHit.GetNPulses(); thit++) {
421 
422  ((THcSignalHit*) frAdcPedRaw->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPedRaw());
423  ((THcSignalHit*) frAdcPed->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPed());
424 
425  ((THcSignalHit*) frAdcPulseIntRaw->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseIntRaw(thit));
426  ((THcSignalHit*) frAdcPulseInt->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseInt(thit));
427 
428  ((THcSignalHit*) frAdcPulseAmpRaw->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseAmpRaw(thit));
429  ((THcSignalHit*) frAdcPulseAmp->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseAmp(thit));
430 
431  ((THcSignalHit*) frAdcPulseTimeRaw->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseTimeRaw(thit));
432  ((THcSignalHit*) frAdcPulseTime->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseTime(thit)+fAdcTdcOffset);
433 
434  if (rawAdcHit.GetPulseAmpRaw(thit) > 0) ((THcSignalHit*) fAdcErrorFlag->ConstructedAt(nrAdcHits))->Set(npmt, 0);
435  if (rawAdcHit.GetPulseAmpRaw(thit) <= 0) ((THcSignalHit*) fAdcErrorFlag->ConstructedAt(nrAdcHits))->Set(npmt, 1);
436 
437  if (rawAdcHit.GetPulseAmpRaw(thit) <= 0) {
438  Double_t PeakPedRatio= rawAdcHit.GetF250_PeakPedestalRatio();
439  Int_t NPedSamples= rawAdcHit.GetF250_NPedestalSamples();
440  Double_t AdcToC = rawAdcHit.GetAdcTopC();
441  Double_t AdcToV = rawAdcHit.GetAdcTomV();
442  if (fPedDefault[npmt-1] !=0) {
443  Double_t tPulseInt = AdcToC*(rawAdcHit.GetPulseIntRaw(thit) - fPedDefault[npmt-1]*PeakPedRatio);
444  ((THcSignalHit*) frAdcPulseInt->ConstructedAt(nrAdcHits))->Set(npmt, tPulseInt);
445  ((THcSignalHit*) frAdcPedRaw->ConstructedAt(nrAdcHits))->Set(npmt, fPedDefault[npmt-1]);
446  ((THcSignalHit*) frAdcPed->ConstructedAt(nrAdcHits))->Set(npmt, float(fPedDefault[npmt-1])/float(NPedSamples)*AdcToV);
447 
448  }
449  ((THcSignalHit*) frAdcPulseAmp->ConstructedAt(nrAdcHits))->Set(npmt, 0.);
450 
451  }
452 
453 
454  ++nrAdcHits;
455  fTotNumAdcHits++;
456  fNumAdcHits.at(npmt-1) = npmt;
457  }
458  ihit++;
459  }
460  return ihit;
461 }
462 
463 //_____________________________________________________________________________
465 {
466  return(0);
467 }
468 
469 //_____________________________________________________________________________
471 {
472  Double_t StartTime = 0.0;
473  if( fglHod ) StartTime = fglHod->GetStartTime();
474  Double_t OffsetTime = 0.0;
475  if( fglHod ) OffsetTime = fglHod->GetOffsetTime();
476  for(Int_t ipmt = 0; ipmt < fNelem; ipmt++) {
477  fAdcPulseAmpTest[ipmt] = -1000.;
478  fAdcGoodElem[ipmt]=-1;
479  }
480  //
481  for(Int_t ielem = 0; ielem < frAdcPulseInt->GetEntries(); ielem++) {
482  Int_t npmt = ((THcSignalHit*) frAdcPulseInt->ConstructedAt(ielem))->GetPaddleNumber() - 1;
483  Double_t pulseTime = ((THcSignalHit*) frAdcPulseTime->ConstructedAt(ielem))->GetData();
484  Double_t pulseAmp = ((THcSignalHit*) frAdcPulseAmp->ConstructedAt(ielem))->GetData();
485  Bool_t errorFlag = ((THcSignalHit*) fAdcErrorFlag->ConstructedAt(ielem))->GetData();
486  Double_t adctdcdiffTime = StartTime-pulseTime+OffsetTime;
487  Bool_t pulseTimeCut = adctdcdiffTime > fAdcTimeWindowMin[npmt] && adctdcdiffTime < fAdcTimeWindowMax[npmt];
488  fGoodAdcMult.at(npmt) += 1;
489  if (!errorFlag) {
490  if (pulseTimeCut && pulseAmp > fAdcPulseAmpTest[npmt]) {
491  fAdcGoodElem[npmt]=ielem;
492  fAdcPulseAmpTest[npmt] = pulseAmp;
493  }
494  } else {
495  if (pulseTimeCut) fAdcGoodElem[npmt]=ielem;
496  }
497  }
498  // Loop over the npmt
499  for(Int_t npmt = 0; npmt < fNelem; npmt++) {
500  Int_t ielem = fAdcGoodElem[npmt];
501  if (ielem != -1) {
502  Double_t pulsePed = ((THcSignalHit*) frAdcPed->ConstructedAt(ielem))->GetData();
503  Double_t pulseInt = ((THcSignalHit*) frAdcPulseInt->ConstructedAt(ielem))->GetData();
504  Double_t pulseIntRaw = ((THcSignalHit*) frAdcPulseIntRaw->ConstructedAt(ielem))->GetData();
505  Double_t pulseAmp = ((THcSignalHit*) frAdcPulseAmp->ConstructedAt(ielem))->GetData();
506  Double_t pulseTime = ((THcSignalHit*) frAdcPulseTime->ConstructedAt(ielem))->GetData();
507  Double_t adctdcdiffTime = StartTime-pulseTime+OffsetTime;
508  fGoodAdcPed.at(npmt) = pulsePed;
509  fGoodAdcHitUsed.at(npmt) = ielem+1;
510  fGoodAdcPulseInt.at(npmt) = pulseInt;
511  fGoodAdcPulseIntRaw.at(npmt) = pulseIntRaw;
512  fGoodAdcPulseAmp.at(npmt) = pulseAmp;
513  fGoodAdcPulseTime.at(npmt) = pulseTime;
514  fGoodAdcTdcDiffTime.at(npmt) = adctdcdiffTime;
515 
516  fNpe.at(npmt) = fGain[npmt]*pulseInt;
517  fNpeSum += fNpe.at(npmt);
518 
520  fNumGoodAdcHits.at(npmt) = npmt + 1;
521  }
522  }
523  return 0;
524 }
525 
526 //_____________________________________________________________________________
528 {
529 
530  Int_t nTracks = tracks.GetLast() + 1;
531 
532  for (Int_t itrack = 0; itrack < nTracks; itrack++) {
533 
534  THaTrack* track = dynamic_cast<THaTrack*> (tracks[itrack]);
535  if (track->GetIndex() != 0) continue; // Select the best track
536 
537  Double_t trackChi2 = track->GetChi2();
538  Int_t trackNDoF = track->GetNDoF();
539  Double_t trackRedChi2 = trackChi2/trackNDoF;
540  Double_t trackBeta = track->GetBeta();
541  Double_t trackEnergy = track->GetEnergy();
542  Double_t trackMom = track->GetP();
543  Double_t trackDp = track->GetDp();
544  Double_t trackENorm = trackEnergy/trackMom;
545  Double_t trackXfp = track->GetX();
546  Double_t trackYfp = track->GetY();
547  Double_t trackTheta = track->GetTheta();
548  Double_t trackPhi = track->GetPhi();
549 
550  Bool_t trackRedChi2Cut = trackRedChi2 > fRedChi2Min && trackRedChi2 < fRedChi2Max;
551  Bool_t trackBetaCut = trackBeta > fBetaMin && trackBeta < fBetaMax;
552  Bool_t trackENormCut = trackENorm > fENormMin && trackENorm < fENormMax;
553  Bool_t trackDpCut = trackDp > fDpMin && trackDp < fDpMax;
554  fXAtCer = trackXfp + trackTheta * fMirrorZPos;
555  fYAtCer = trackYfp + trackPhi * fMirrorZPos;
556 
557  if (trackRedChi2Cut && trackBetaCut && trackENormCut && trackDpCut) {
558 
559  // Project the track to the Cherenkov mirror planes
560 
561  // cout << "Cherenkov Detector: " << GetName() << " has fNRegions = " << fNRegions << endl;
562  // cout << "nTracks = " << nTracks << "\t" << "trackChi2 = " << trackChi2
563  // << "\t" << "trackNDof = " << trackNDoF << "\t" << "trackRedChi2 = " << trackRedChi2 << endl;
564  // cout << "trackBeta = " << trackBeta << "\t" << "trackEnergy = " << trackEnergy << "\t"
565  // << "trackMom = " << trackMom << "\t" << "trackENorm = " << trackENorm << endl;
566  // cout << "trackXfp = " << trackXfp << "\t" << "trackYfp = " << trackYfp << "\t"
567  // << "trackTheta = " << trackTheta << "\t" << "trackPhi = " << trackPhi << endl;
568  // cout << "fMirrorZPos = " << fMirrorZPos << "\t" << "fXAtCer = " << fXAtCer << "\t" << "fYAtCer = " << fYAtCer << endl;
569  // cout << "=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:" << endl;
570 
571  for (Int_t iregion = 0; iregion < fNRegions; iregion++) {
572 
573  if ((TMath::Abs(fRegionValue[GetIndex(iregion, 0)] - fXAtCer) < fRegionValue[GetIndex(iregion, 4)]) &&
574  (TMath::Abs(fRegionValue[GetIndex(iregion, 1)] - fYAtCer) < fRegionValue[GetIndex(iregion, 5)]) &&
575  (TMath::Abs(fRegionValue[GetIndex(iregion, 2)] - trackTheta) < fRegionValue[GetIndex(iregion, 6)]) &&
576  (TMath::Abs(fRegionValue[GetIndex(iregion, 3)] - trackPhi) < fRegionValue[GetIndex(iregion, 7)])) {
577 
579  fNumTracksMatched.at(iregion) = iregion + 1;
580 
581  if (fNpeSum > fNpeThresh) {
583  fNumTracksFired.at(iregion) = iregion + 1;
584  } // NPE threshold cut
585  } // Regional cuts
586  } // Loop over regions
587  } // Tracking cuts
588  } // Track loop
589 
590  return 0;
591 }
592 
593 //_____________________________________________________________________________
595 {
596  fNPedestalEvents = 0;
597  fMinPeds = 500; // In engine, this is set in parameter file
598  fPedSum = new Int_t [fNelem];
599  fPedSum2 = new Int_t [fNelem];
600  fPedCount = new Int_t [fNelem];
601  fPed = new Double_t [fNelem];
602  fThresh = new Double_t [fNelem];
603  for(Int_t i = 0; i < fNelem; i++) {
604  fPedSum[i] = 0;
605  fPedSum2[i] = 0;
606  fPedCount[i] = 0;
607  }
608 }
609 
610 //_____________________________________________________________________________
612 {
613  // Extract data from the hit list, accumulating into arrays for
614  // calculating pedestals
615 
616  Int_t nrawhits = rawhits->GetLast()+1;
617 
618  Int_t ihit = 0;
619  while(ihit < nrawhits) {
620  THcCherenkovHit* hit = (THcCherenkovHit *) rawhits->At(ihit);
621 
622  Int_t element = hit->fCounter - 1;
623  Int_t nadc = hit->GetRawAdcHitPos().GetPulseIntRaw();
624  if(nadc <= fPedLimit[element]) {
625  fPedSum[element] += nadc;
626  fPedSum2[element] += nadc*nadc;
627  fPedCount[element]++;
628  if(fPedCount[element] == fMinPeds/5) {
629  fPedLimit[element] = 100 + fPedSum[element]/fPedCount[element];
630  }
631  }
632  ihit++;
633  }
634 
636 
637  return;
638 }
639 
640 //_____________________________________________________________________________
642 {
643  // Use the accumulated pedestal data to calculate pedestals
644  // Later add check to see if pedestals have drifted ("Danger Will Robinson!")
645  // cout << "Plane: " << fPlaneNum << endl;
646  for(Int_t i=0; i<fNelem;i++) {
647 
648  // PMT tubes
649  fPed[i] = ((Double_t) fPedSum[i]) / TMath::Max(1, fPedCount[i]);
650  fThresh[i] = fPed[i] + 15;
651 
652  // Just a copy for now, but allow the possibility that fXXXPedMean is set
653  // in a parameter file and only overwritten if there is a sufficient number of
654  // pedestal events. (So that pedestals are sensible even if the pedestal events were
655  // not acquired.)
656  if(fMinPeds > 0) {
657  if(fPedCount[i] > fMinPeds) {
658  fPedMean[i] = fPed[i];
659  }
660  }
661  }
662  // cout << " " << endl;
663 
664 }
665 
666 //_____________________________________________________________________________
668 {
669  return fNRegions * nValue + nRegion;
670 }
671 
672 
673 //_____________________________________________________________________________
674 void THcCherenkov::Print(const Option_t* opt) const
675 {
676  THaNonTrackingDetector::Print(opt);
677  // Print out the pedestals
678  cout << endl;
679  cout << "Cherenkov Pedestals" << endl;
680  // Ahmed
681  cout << "No. ADC" << endl;
682  for(Int_t i=0; i<fNelem; i++)
683  cout << " " << i << " " << fPed[i] << endl;
684  cout << endl;
685 }
686 
687 //_____________________________________________________________________________
689 {
690  return fNpeSum;
691 }
692 //_____________________________________________________________________________
693 Int_t THcCherenkov::End(THaRunBase* run)
694 {
695  MissReport(Form("%s.%s", GetApparatus()->GetName(), GetName()));
696  return 0;
697 }
Double_t GetF250_PeakPedestalRatio()
Definition: THcRawAdcHit.h:24
std::string GetName(const std::string &scope_name)
virtual void Clear(Option_t *opt="")
vector< Int_t > fNumTracksMatched
Definition: THcCherenkov.h:67
virtual Int_t FillMap(THaDetMap *detmap, const char *detectorname)
Double_t * fPed
Definition: THcCherenkov.h:107
Double_t fRedChi2Max
Definition: THcCherenkov.h:82
Class representing a Cherenkov PMT hit.
Double_t fENormMin
Definition: THcCherenkov.h:85
vector< Int_t > fNumAdcHits
Definition: THcCherenkov.h:65
Int_t GetLast() const
const char Option_t
Class for gas Cherenkov detectors.
Definition: THcCherenkov.h:16
Double_t fDpMin
Definition: THcCherenkov.h:87
Double_t fNpeSum
Definition: THcCherenkov.h:62
Int_t fRegionsValueMax
Definition: THcCherenkov.h:80
void MissReport(const char *name)
Definition: THcHitList.cxx:557
static const Int_t MaxNumAdcPulse
Definition: THcCherenkov.h:45
Int_t fAnalyzePedestals
Definition: THcCherenkov.h:50
Double_t GetOffsetTime() const
Definition: THcHodoscope.h:59
Int_t GetPulseTimeRaw(UInt_t iPulse=0) const
Gets raw pulse time. In subsamples.
int Int_t
void DeleteArrays()
Int_t fDebugAdc
Definition: THcCherenkov.h:51
bool Bool_t
vector< Int_t > fNumTracksFired
Definition: THcCherenkov.h:68
Double_t * fAdcPulseAmpTest
Definition: THcCherenkov.h:109
Int_t fTotNumAdcHits
Definition: THcCherenkov.h:57
STL namespace.
TClonesArray * frAdcPulseInt
Definition: THcCherenkov.h:118
TClonesArray * fRawHitList
Definition: THcHitList.h:51
Int_t fCounter
Definition: THcRawHit.h:55
Short_t Abs(Short_t d)
Double_t * fAdcTimeWindowMin
Definition: THcCherenkov.h:93
virtual Int_t ApplyCorrections(void)
Double_t * fAdcTimeWindowMax
Definition: THcCherenkov.h:94
Double_t fAdcTdcOffset
Definition: THcCherenkov.h:96
virtual Int_t ReadDatabase(const TDatime &date)
Bool_t * fPresentP
Definition: THcCherenkov.h:49
Double_t fNpeThresh
Definition: THcCherenkov.h:92
TObject * ConstructedAt(Int_t idx)
Double_t fRedChi2Min
Definition: THcCherenkov.h:81
TClonesArray * frAdcPulseTimeRaw
Definition: THcCherenkov.h:116
virtual Int_t Decode(const THaEvData &)
tuple app
Int_t * fPedSum2
Definition: THcCherenkov.h:103
Double_t fBetaMax
Definition: THcCherenkov.h:84
vector< Double_t > fGoodAdcPed
Definition: THcCherenkov.h:69
Double_t fRefTime
Definition: THcCherenkov.h:61
TClonesArray * frAdcPedRaw
Definition: THcCherenkov.h:113
virtual Int_t FineProcess(TClonesArray &tracks)
Double_t GetPed() const
Gets sample pedestal. In channels.
void Error(const char *location, const char *msgfmt,...)
TClonesArray * frAdcPulseAmpRaw
Definition: THcCherenkov.h:115
static const Int_t MaxNumCerPmt
Definition: THcCherenkov.h:44
vector< Double_t > fGoodAdcMult
Definition: THcCherenkov.h:70
virtual void Clear(Option_t *option="")
Double_t GetCerNPE()
void InitHitList(THaDetMap *detmap, const char *hitclass, Int_t maxhits, Int_t tdcref_cut=0, Int_t adcref_cut=0)
Save the electronics module to detector mapping and initialize a hit array of hits of class hitclass...
Definition: THcHitList.cxx:65
TClonesArray * frAdcPulseTime
Definition: THcCherenkov.h:120
vector< Double_t > fGoodAdcPulseIntRaw
Definition: THcCherenkov.h:73
THcRawAdcHit & GetRawAdcHitPos()
Double_t * fGain
Definition: THcCherenkov.h:63
Int_t GetPulseAmpRaw(UInt_t iPulse=0) const
Gets raw pulse amplitude. In channels.
Int_t fADC_RefTimeCut
Definition: THcCherenkov.h:54
vector< Double_t > fGoodAdcTdcDiffTime
Definition: THcCherenkov.h:76
Double_t GetPulseAmp(UInt_t iPulse=0) const
Gets pedestal subtracted pulse amplitude. In channels.
Double_t GetPulseInt(UInt_t iPulse=0) const
Gets pedestal subtracted pulse integral. In channels.
Int_t fTotNumTracksFired
Definition: THcCherenkov.h:60
virtual EStatus Init(const TDatime &run_time)
unsigned int UInt_t
char * Form(const char *fmt,...)
void CreateMissReportParms(const char *prefix)
Definition: THcHitList.cxx:544
TClonesArray * fAdcErrorFlag
Definition: THcCherenkov.h:121
Double_t GetStartTime() const
Definition: THcHodoscope.h:58
Int_t fNRegions
Definition: THcCherenkov.h:79
Int_t * fPedDefault
Definition: THcCherenkov.h:95
void Warning(const char *location, const char *msgfmt,...)
TClonesArray * frAdcPed
Definition: THcCherenkov.h:117
tuple list
Definition: SConscript.py:9
Double_t fBetaMin
Definition: THcCherenkov.h:83
vector< Double_t > fGoodAdcHitUsed
Definition: THcCherenkov.h:71
Int_t * fPedLimit
Definition: THcCherenkov.h:104
Int_t GetPulseIntRaw(UInt_t iPulse=0) const
Gets raw pulse integral. In channels.
Int_t * fPedSum
Definition: THcCherenkov.h:102
R__EXTERN class THcDetectorMap * gHcDetectorMap
Definition: THcGlobals.h:12
Int_t GetRefTime() const
Double_t fYAtCer
Definition: THcCherenkov.h:91
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
double Double_t
Int_t fNPedestalEvents
Definition: THcCherenkov.h:100
Bool_t HasRefTime() const
THcHodoscope * fglHod
Definition: THcCherenkov.h:125
Double_t fDpMax
Definition: THcCherenkov.h:88
UInt_t GetNPulses() const
Gets number of set pulses.
TClonesArray * frAdcPulseAmp
Definition: THcCherenkov.h:119
Int_t fTotNumTracksMatched
Definition: THcCherenkov.h:59
Double_t * fThresh
Definition: THcCherenkov.h:108
Int_t fTotNumGoodAdcHits
Definition: THcCherenkov.h:58
virtual void InitializePedestals()
virtual ~THcCherenkov()
Int_t End(THaRunBase *run)
Int_t GetEntries() const
Double_t fENormMax
Definition: THcCherenkov.h:86
ClassImp(THcDCLookupTTDConv) THcDCLookupTTDConv
virtual void AccumulatePedestals(TClonesArray *rawhits)
virtual void Print(const Option_t *opt) const
TClonesArray * frAdcPulseIntRaw
Definition: THcCherenkov.h:114
virtual Int_t CoarseProcess(TClonesArray &tracks)
virtual void CalculatePedestals()
Int_t GetPedRaw() const
Gets raw signal pedestal. In channels.
Short_t Max(Short_t a, Short_t b)
vector< Double_t > fGoodAdcPulseInt
Definition: THcCherenkov.h:72
Double_t fXAtCer
Definition: THcCherenkov.h:90
Double_t * fPedMean
Definition: THcCherenkov.h:106
Double_t GetPulseTime(UInt_t iPulse=0) const
vector< Double_t > fGoodAdcPulseAmp
Definition: THcCherenkov.h:74
Double_t * fRegionValue
Definition: THcCherenkov.h:97
Double_t fMirrorZPos
Definition: THcCherenkov.h:89
R__EXTERN class THcParmList * gHcParms
Definition: THcGlobals.h:11
Int_t * fAdcGoodElem
Definition: THcCherenkov.h:110
TObject * At(Int_t idx) const
Int_t GetF250_NPedestalSamples()
Definition: THcRawAdcHit.h:25
void GetData(std::string s, double *x, double *y, double *ey)
vector< Double_t > fGoodAdcPulseTime
Definition: THcCherenkov.h:75
static const Double_t kBig
Definition: THcFormula.cxx:31
virtual Int_t DefineVariables(EMode mode=kDefine)
const Bool_t kTRUE
Int_t * fPedCount
Definition: THcCherenkov.h:105
Double_t GetAdcTopC() const
Class representing a single raw ADC hit.
Definition: THcRawAdcHit.h:7
virtual Int_t DecodeToHitList(const THaEvData &evdata, Bool_t suppress=kFALSE)
Populate the hitlist from the raw event data.
Definition: THcHitList.cxx:191
vector< Double_t > fNpe
Definition: THcCherenkov.h:77
Double_t GetAdcTomV() const
Int_t GetIndex(Int_t nRegion, Int_t nValue)
A standard Hall C spectrometer apparatus.
vector< Int_t > fNumGoodAdcHits
Definition: THcCherenkov.h:66