Hall C ROOT/C++ Analyzer (hcana)
THcTrigDet.cxx
Go to the documentation of this file.
1 
96 //TODO: Check if fNumAdc < fMaxAdcChannels && fNumTdc < fMaxTdcChannels.
97 
98 #include "THcTrigDet.h"
99 
100 #include <algorithm>
101 #include <iostream>
102 #include <stdexcept>
103 
104 #include "TDatime.h"
105 #include "TString.h"
106 
107 #include "THaApparatus.h"
108 #include "THaEvData.h"
109 
110 #include "THcDetectorMap.h"
111 #include "THcGlobals.h"
112 #include "THcParmList.h"
113 #include "THcRawAdcHit.h"
114 #include "THcRawTdcHit.h"
115 #include "THcTrigApp.h"
116 #include "THcTrigRawHit.h"
117 
118 //_____________________________________________________________________________
120 
121 //_____________________________________________________________________________
123  const char* name, const char* description, THaApparatus* app
124 ) :
125  THaDetector(name, description, app), THcHitList(),
126  fKwPrefix(""),
127  fNumAdc(0), fNumTdc(0), fAdcNames(), fTdcNames(),
128  fTdcTimeRaw(), fTdcTime(),
129  fAdcPedRaw(), fAdcPulseIntRaw(), fAdcPulseAmpRaw(), fAdcPulseTimeRaw(),
130  fAdcPed(), fAdcPulseInt(), fAdcPulseAmp(), fAdcPulseTime(),
131  fTdcMultiplicity(), fAdcMultiplicity()
132 {
133  // Guess at spectrometer name that this trigger detector is associated with
134  // Can override with SetSpectName
135  fSpectName = name[0];
136 }
137 
138 //_____________________________________________________________________________
140  delete [] fAdcTimeWindowMin; fAdcTimeWindowMin = NULL;
141  delete [] fAdcTimeWindowMax; fAdcTimeWindowMax = NULL;
142  delete [] fTdcTimeWindowMin; fTdcTimeWindowMin = NULL;
143  delete [] fTdcTimeWindowMax; fTdcTimeWindowMax = NULL;
144 
145 }
146 
147 //_____________________________________________________________________________
148 THaAnalysisObject::EStatus THcTrigDet::Init(const TDatime& date) {
149  // Call `Setup` before everything else.
150  Setup(GetName(), GetTitle());
151 
152  // Initialize all variables.
153  for (int i=0; i<fMaxAdcChannels; ++i) {
154  fAdcPedRaw[i] = 0;
155  fAdcPulseIntRaw[i] = 0;
156  fAdcPulseAmpRaw[i] = 0;
157  fAdcPulseTimeRaw[i] = 0;
158  fAdcPulseTime[i] = kBig;
159  fAdcPed[i] = 0.0;
160  fAdcPulseInt[i] = 0.0;
161  fAdcPulseAmp[i] = 0.0;
162  fAdcMultiplicity[i] = 0;
163  };
164  for (int i=0; i<fMaxTdcChannels; ++i) {
165  fTdcTimeRaw[i] = 0;
166  fTdcTime[i] = 0.0;
167  fTdcMultiplicity[i] = 0;
168  };
169 
170  // Call initializer for base class.
171  // This also calls `ReadDatabase` and `DefineVariables`.
172  EStatus status = THaDetector::Init(date);
173  if (status) {
174  fStatus = status;
175  return fStatus;
176 
177 }
178  // Fill in detector map.
179  string EngineDID = string(GetApparatus()->GetName()).substr(0, 1) + GetName();
180  std::transform(EngineDID.begin(), EngineDID.end(), EngineDID.begin(), ::toupper);
181  if (gHcDetectorMap->FillMap(fDetMap, EngineDID.c_str()) < 0) {
182  static const char* const here = "Init()";
183  Error(Here(here), "Error filling detectormap for %s.", EngineDID.c_str());
184  return kInitError;
185  }
186  DBRequest listextra[]={
187  {"_trig_tdcrefcut", &fTDC_RefTimeCut, kInt, 0, 1},
188  {"_trig_adcrefcut", &fADC_RefTimeCut, kInt, 0, 1},
189  {0}
190  };
191  fTDC_RefTimeCut = -1000; // Minimum allowed reference times
192  fADC_RefTimeCut = -1000;
193  gHcParms->LoadParmValues((DBRequest*)&listextra,fKwPrefix.c_str());
194  // Initialize hitlist part of the class.
195  // printf(" Init trig det hitlist\n");
196  InitHitList(fDetMap, "THcTrigRawHit", 200,fTDC_RefTimeCut,fADC_RefTimeCut);
198 
199  fPresentP = 0;
200  THaVar* vpresent = gHaVars->Find(Form("%s.present",fSpectName.Data()));
201  if(vpresent) {
202  fPresentP = (Bool_t *) vpresent->GetValuePointer();
203  }
204 
205  fStatus = kOK;
206  return fStatus;
207 }
208 
209 //_____________________________________________________________________________
211  THaAnalysisObject::Clear(opt);
212 
213  // Reset all data.
214  fTdcRefTime = kBig;
215  for (int i=0; i<fNumAdc; ++i) {
216  fAdcPedRaw[i] = 0;
217  fAdcPulseIntRaw[i] = 0;
218  fAdcPulseAmpRaw[i] = 0;
219  fAdcPulseTimeRaw[i] = 0;
220  fAdcPulseTime[i] = kBig;
221  fAdcPed[i] = 0.0;
222  fAdcPulseInt[i] = 0.0;
223  fAdcPulseAmp[i] = 0.0;
224  fAdcMultiplicity[i] = 0;
225  };
226  for (int i=0; i<fNumTdc; ++i) {
227  fTdcTimeRaw[i] = 0;
228  fTdcTime[i] = 0.0;
229  fTdcMultiplicity[i] = 0;
230  };
231 }
232 
233 //_____________________________________________________________________________
235 
236  // Decode raw data for this event.
237  Bool_t present = kTRUE; // Don't suppress reference time warnings
238  if(HaveIgnoreList()) {
239  if(IsIgnoreType(evData.GetEvType())) {
240  present = kFALSE;
241  }
242  } else if(fPresentP) { // if this spectrometer not part of trigger
243  present = *fPresentP;
244  }
245  Int_t numHits = DecodeToHitList(evData, !present);
246 
247  // Process each hit and fill variables.
248  Int_t iHit = 0;
249  while (iHit < numHits) {
250  THcTrigRawHit* hit = dynamic_cast<THcTrigRawHit*>(fRawHitList->At(iHit));
251 
252  Int_t cnt = hit->fCounter-1;
253  if (hit->fPlane == 1) {
254  THcRawAdcHit& rawAdcHit = hit->GetRawAdcHit();
255  fAdcMultiplicity[cnt] = rawAdcHit.GetNPulses();
256  UInt_t good_hit=999;
257  for (UInt_t thit=0; thit<rawAdcHit.GetNPulses(); ++thit) {
258  Int_t TestTime=rawAdcHit.GetPulseTimeRaw(thit);
259  if (TestTime>=fAdcTimeWindowMin[cnt]&&TestTime<=fAdcTimeWindowMax[cnt]&&good_hit==999) {
260  good_hit=thit;
261  }
262  }
263  if (good_hit<rawAdcHit.GetNPulses()) {
264  fAdcPedRaw[cnt] = rawAdcHit.GetPedRaw();
265  fAdcPulseIntRaw[cnt] = rawAdcHit.GetPulseIntRaw(good_hit);
266  fAdcPulseAmpRaw[cnt] = rawAdcHit.GetPulseAmpRaw(good_hit);
267  fAdcPulseTimeRaw[cnt] = rawAdcHit.GetPulseTimeRaw(good_hit);
268  fAdcPulseTime[cnt] = rawAdcHit.GetPulseTime(good_hit)+fAdcTdcOffset;
269  fAdcPed[cnt] = rawAdcHit.GetPed();
270  fAdcPulseInt[cnt] = rawAdcHit.GetPulseInt(good_hit);
271  fAdcPulseAmp[cnt] = rawAdcHit.GetPulseAmp(good_hit);
272  }
273  }
274  else if (hit->fPlane == 2) {
275  THcRawTdcHit& rawTdcHit = hit->GetRawTdcHit();
276  if (rawTdcHit.GetNHits() >0 && rawTdcHit.HasRefTime() && fTdcRefTime == kBig) fTdcRefTime=rawTdcHit.GetRefTime() ;
277  UInt_t good_hit=999;
278  UInt_t closest_hit=999;
279  Int_t TimeDiff=1000000;
280  for (UInt_t thit=0; thit<rawTdcHit.GetNHits(); ++thit) {
281  Int_t TestTime= rawTdcHit.GetTimeRaw(thit);
282  if (abs(TestTime-fTdcTimeWindowMin[cnt]) < TimeDiff) {
283  closest_hit=thit;
284  TimeDiff=abs(TestTime-fTdcTimeWindowMin[cnt]);
285  }
286  if (TestTime>=fTdcTimeWindowMin[cnt]&&TestTime<=fTdcTimeWindowMax[cnt]&&good_hit==999) {
287  good_hit=thit;
288  }
289  }
290  if (good_hit == 999 and closest_hit != 999) good_hit=closest_hit;
291  if (good_hit<rawTdcHit.GetNHits()) {
292  fTdcTimeRaw[cnt] = rawTdcHit.GetTimeRaw(good_hit);
293  fTdcTime[cnt] = rawTdcHit.GetTime(good_hit)*fTdcChanperNS+fTdcOffset;
294  }
295  fTdcMultiplicity[cnt] = rawTdcHit.GetNHits();
296  }
297  else {
298  throw std::out_of_range(
299  "`THcTrigDet::Decode`: only planes `1` and `2` available!"
300  );
301  }
302 
303  ++iHit;
304  }
305 
306  return 0;
307 }
308 
309 //_____________________________________________________________________________
310 void THcTrigDet::Setup(const char* name, const char* description) {
311  // Prefix for parameters in `param` file.
312  string kwPrefix = string(GetApparatus()->GetName()) + "_" + name;
313  std::transform(kwPrefix.begin(), kwPrefix.end(), kwPrefix.begin(), ::tolower);
314  fKwPrefix = kwPrefix;
315 }
316 
317 //_____________________________________________________________________________
319  std::string adcNames, tdcNames;
320  std::string trigNames="pTRIG1_ROC1 pTRIG4_ROC1 pTRIG1_ROC2 pTRIG4_ROC2";
321  // SJDK 12/04/21 - Added new RF names for use in getter
322  std::string RFNames="pRF hRF";
323  DBRequest list[] = {
324  {"_numAdc", &fNumAdc, kInt}, // Number of ADC channels.
325  {"_numTdc", &fNumTdc, kInt}, // Number of TDC channels.
326  {"_adcNames", &adcNames, kString}, // Names of ADC channels.
327  {"_tdcNames", &tdcNames, kString}, // Names of TDC channels.
328  {"_trigNames", &trigNames, kString,0,1}, // Names of Triggers for coincidence time.
329  {"_RFNames", &RFNames, kString,0, 1}, // SJDK 12/04/21 - Names for RF time
330  {"_tdcoffset", &fTdcOffset, kDouble,0,1}, // Offset of tdc channels
331  {"_adc_tdc_offset", &fAdcTdcOffset, kDouble,0,1}, // Offset of Adc Pulse time (ns)
332  {"_tdcchanperns", &fTdcChanperNS, kDouble,0,1}, // Convert channesl to ns
333  {0}
334  };
335  fTdcChanperNS=0.09766;
336  fTdcOffset=300.;
337  fAdcTdcOffset=200.;
338  gHcParms->LoadParmValues(list, fKwPrefix.c_str());
339 
344  for(Int_t ip=0;ip<fNumAdc;ip++) { // Set a large default window
345  fAdcTimeWindowMin[ip]=0;
346  fAdcTimeWindowMax[ip]=100000;
347  }
348  for(Int_t ip=0;ip<fNumTdc;ip++) { // Set a large default window
349  fTdcTimeWindowMin[ip]=0;
350  fTdcTimeWindowMax[ip]=100000;
351  }
352  DBRequest list2[]={
353  {"_AdcTimeWindowMin", fAdcTimeWindowMin, kDouble, (UInt_t) fNumAdc, 1},
354  {"_AdcTimeWindowMax", fAdcTimeWindowMax, kDouble, (UInt_t) fNumAdc, 1},
355  {"_TdcTimeWindowMin", fTdcTimeWindowMin, kDouble, (UInt_t) fNumTdc, 1},
356  {"_TdcTimeWindowMax", fTdcTimeWindowMax, kDouble, (UInt_t) fNumTdc, 1},
357  {0}
358  };
359 
360  gHcParms->LoadParmValues(list2, fKwPrefix.c_str());
361  for(Int_t ip=0;ip<fNumTdc;ip++) {
362  // cout << ip << " " << fTdcNames.at(ip) << " " << fTdcTimeWindowMin[ip] << " " << fTdcTimeWindowMax[ip] << endl;
363  }
364  // Split the names to std::vector<std::string>.
365  fAdcNames = Podd::vsplit(adcNames);
366  fTdcNames = Podd::vsplit(tdcNames);
367  fTrigNames = Podd::vsplit(trigNames);
368  fRFNames = Podd::vsplit(RFNames); // SJDK 12/04/21 - For RF getter
369  //default index values
370 
371  //Assign an index to coincidence trigger times strings
372  for (UInt_t j = 0; j <fTrigNames.size(); j++) {
373  fTrigId[j]=-1;
374  }
375  for (int i = 0; i <fNumTdc; i++) {
376  for (UInt_t j = 0; j <fTrigNames.size(); j++) {
377  if(fTdcNames.at(i)==fTrigNames[j]) fTrigId[j]=i;
378  }
379  }
380 
381  cout << " Trig = " << fTrigNames.size() << endl;
382  for (UInt_t j = 0; j <fTrigNames.size(); j++) {
383  cout << fTrigNames[j] << " " << fTrigId[j] << endl;
384  }
385 
386  // SJDK - 12/04/21 - For RF getter
387  // Assign an index to RF times strings
388  for (UInt_t j = 0; j <fRFNames.size(); j++) {
389  fRFId[j]=-1;
390  }
391  for (int i = 0; i <fNumTdc; i++) {
392  for (UInt_t j = 0; j <fRFNames.size(); j++) {
393  if(fTdcNames.at(i)==fRFNames[j]) fRFId[j]=i;
394  }
395  }
396  for (UInt_t j = 0; j <fRFNames.size(); j++) {
397  cout << fRFNames[j] << " " << fRFId[j] << endl;
398  }
399 
400  return kOK;
401 }
402 
403 //_____________________________________________________________________________
404 Int_t THcTrigDet::DefineVariables(THaAnalysisObject::EMode mode) {
405 
406  if (mode == kDefine && fIsSetup) return kOK;
407  fIsSetup = (mode == kDefine);
408 
409  std::vector<RVarDef> vars;
410 
411  //Push the variable names for ADC channels.
412  std::vector<TString> adcPedRawTitle(fNumAdc), adcPedRawVar(fNumAdc);
413  std::vector<TString> adcPulseIntRawTitle(fNumAdc), adcPulseIntRawVar(fNumAdc);
414  std::vector<TString> adcPulseAmpRawTitle(fNumAdc), adcPulseAmpRawVar(fNumAdc);
415  std::vector<TString> adcPulseTimeRawTitle(fNumAdc), adcPulseTimeRawVar(fNumAdc);
416  std::vector<TString> adcPulseTimeTitle(fNumAdc), adcPulseTimeVar(fNumAdc);
417  std::vector<TString> adcPedTitle(fNumAdc), adcPedVar(fNumAdc);
418  std::vector<TString> adcPulseIntTitle(fNumAdc), adcPulseIntVar(fNumAdc);
419  std::vector<TString> adcPulseAmpTitle(fNumAdc), adcPulseAmpVar(fNumAdc);
420  std::vector<TString> adcMultiplicityTitle(fNumAdc), adcMultiplicityVar(fNumAdc);
421 
422  TString RefTimeTitle= "TdcRefTime";
423  TString RefTimeVar= "fTdcRefTime";
424  RVarDef entryRefTime {
425  RefTimeTitle.Data(),
426  RefTimeTitle.Data(),
427  RefTimeVar.Data()
428  };
429  vars.push_back(entryRefTime);
430 
431  for (int i=0; i<fNumAdc; ++i) {
432  adcPedRawTitle.at(i) = fAdcNames.at(i) + "_adcPedRaw";
433  adcPedRawVar.at(i) = TString::Format("fAdcPedRaw[%d]", i);
434  RVarDef entry1 {
435  adcPedRawTitle.at(i).Data(),
436  adcPedRawTitle.at(i).Data(),
437  adcPedRawVar.at(i).Data()
438  };
439  vars.push_back(entry1);
440 
441  adcPulseIntRawTitle.at(i) = fAdcNames.at(i) + "_adcPulseIntRaw";
442  adcPulseIntRawVar.at(i) = TString::Format("fAdcPulseIntRaw[%d]", i);
443  RVarDef entry2 {
444  adcPulseIntRawTitle.at(i).Data(),
445  adcPulseIntRawTitle.at(i).Data(),
446  adcPulseIntRawVar.at(i).Data()
447  };
448  vars.push_back(entry2);
449 
450  adcPulseAmpRawTitle.at(i) = fAdcNames.at(i) + "_adcPulseAmpRaw";
451  adcPulseAmpRawVar.at(i) = TString::Format("fAdcPulseAmpRaw[%d]", i);
452  RVarDef entry3 {
453  adcPulseAmpRawTitle.at(i).Data(),
454  adcPulseAmpRawTitle.at(i).Data(),
455  adcPulseAmpRawVar.at(i).Data()
456  };
457  vars.push_back(entry3);
458 
459  adcPulseTimeRawTitle.at(i) = fAdcNames.at(i) + "_adcPulseTimeRaw";
460  adcPulseTimeRawVar.at(i) = TString::Format("fAdcPulseTimeRaw[%d]", i);
461  RVarDef entry4 {
462  adcPulseTimeRawTitle.at(i).Data(),
463  adcPulseTimeRawTitle.at(i).Data(),
464  adcPulseTimeRawVar.at(i).Data()
465  };
466  vars.push_back(entry4);
467 
468  adcPedTitle.at(i) = fAdcNames.at(i) + "_adcPed";
469  adcPedVar.at(i) = TString::Format("fAdcPed[%d]", i);
470  RVarDef entry5 {
471  adcPedTitle.at(i).Data(),
472  adcPedTitle.at(i).Data(),
473  adcPedVar.at(i).Data()
474  };
475  vars.push_back(entry5);
476 
477  adcPulseIntTitle.at(i) = fAdcNames.at(i) + "_adcPulseInt";
478  adcPulseIntVar.at(i) = TString::Format("fAdcPulseInt[%d]", i);
479  RVarDef entry6 {
480  adcPulseIntTitle.at(i).Data(),
481  adcPulseIntTitle.at(i).Data(),
482  adcPulseIntVar.at(i).Data()
483  };
484  vars.push_back(entry6);
485 
486  adcPulseAmpTitle.at(i) = fAdcNames.at(i) + "_adcPulseAmp";
487  adcPulseAmpVar.at(i) = TString::Format("fAdcPulseAmp[%d]", i);
488  RVarDef entry7 {
489  adcPulseAmpTitle.at(i).Data(),
490  adcPulseAmpTitle.at(i).Data(),
491  adcPulseAmpVar.at(i).Data()
492  };
493  vars.push_back(entry7);
494 
495  adcMultiplicityTitle.at(i) = fAdcNames.at(i) + "_adcMultiplicity";
496  adcMultiplicityVar.at(i) = TString::Format("fAdcMultiplicity[%d]", i);
497  RVarDef entry8 {
498  adcMultiplicityTitle.at(i).Data(),
499  adcMultiplicityTitle.at(i).Data(),
500  adcMultiplicityVar.at(i).Data()
501  };
502  vars.push_back(entry8);
503 
504  adcPulseTimeTitle.at(i) = fAdcNames.at(i) + "_adcPulseTime";
505  adcPulseTimeVar.at(i) = TString::Format("fAdcPulseTime[%d]", i);
506  RVarDef entry9 {
507  adcPulseTimeTitle.at(i).Data(),
508  adcPulseTimeTitle.at(i).Data(),
509  adcPulseTimeVar.at(i).Data()
510  };
511  vars.push_back(entry9);
512  } // loop over fNumAdc
513  // Push the variable names for TDC channels.
514  std::vector<TString> tdcTimeRawTitle(fNumTdc), tdcTimeRawVar(fNumTdc);
515  std::vector<TString> tdcTimeTitle(fNumTdc), tdcTimeVar(fNumTdc);
516  std::vector<TString> tdcMultiplicityTitle(fNumTdc), tdcMultiplicityVar(fNumTdc);
517 
518  for (int i=0; i<fNumTdc; ++i) {
519  tdcTimeRawTitle.at(i) = fTdcNames.at(i) + "_tdcTimeRaw";
520  tdcTimeRawVar.at(i) = TString::Format("fTdcTimeRaw[%d]", i);
521 
522  RVarDef entry1 {
523  tdcTimeRawTitle.at(i).Data(),
524  tdcTimeRawTitle.at(i).Data(),
525  tdcTimeRawVar.at(i).Data()
526  };
527  vars.push_back(entry1);
528 
529  tdcTimeTitle.at(i) = fTdcNames.at(i) + "_tdcTime";
530  tdcTimeVar.at(i) = TString::Format("fTdcTime[%d]", i);
531  RVarDef entry2 {
532  tdcTimeTitle.at(i).Data(),
533  tdcTimeTitle.at(i).Data(),
534  tdcTimeVar.at(i).Data()
535  };
536  vars.push_back(entry2);
537 
538  tdcMultiplicityTitle.at(i) = fTdcNames.at(i) + "_tdcMultiplicity";
539  tdcMultiplicityVar.at(i) = TString::Format("fTdcMultiplicity[%d]", i);
540  RVarDef entry3 {
541  tdcMultiplicityTitle.at(i).Data(),
542  tdcMultiplicityTitle.at(i).Data(),
543  tdcMultiplicityVar.at(i).Data()
544  };
545  vars.push_back(entry3);
546  }
547 
548  RVarDef end {0};
549  vars.push_back(end);
550 
551  return DefineVarsFromList(vars.data(), mode);
552 }
553 void THcTrigDet::SetSpectName( const char* name)
554 {
555  fSpectName = name;
556 }
557 
558 void THcTrigDet::AddEvtType(int evtype) {
559  eventtypes.push_back(evtype);
560 }
561 
562 void THcTrigDet::SetEvtType(int evtype) {
563  eventtypes.clear();
564  AddEvtType(evtype);
565 }
566 
567 //_____________________________________________________________________________
569 {
570  for (UInt_t i=0; i < eventtypes.size(); i++) {
571  if (evtype == eventtypes[i]) return kTRUE;
572  }
573  return kFALSE;
574 }
575 
576 //_____________________________________________________________________________
578 {
579  return( (eventtypes.size()>0) ? kTRUE : kFALSE);
580 }
581 
582 //_____________________________________________________________________________
583 Int_t THcTrigDet::End(THaRunBase* run)
584 {
585  MissReport(Form("%s.%s", GetApparatus()->GetName(), GetName()));
586  return 0;
587 }
588 
Int_t GetTimeRaw(UInt_t iHit=0) const
Gets raw TDC time. In channels.
std::string GetName(const std::string &scope_name)
virtual Int_t ReadDatabase(const TDatime &date)
Definition: THcTrigDet.cxx:318
Double_t fTdcTime[fMaxTdcChannels]
Definition: THcTrigDet.h:66
Int_t fTrigId[4]
Definition: THcTrigDet.h:50
virtual Int_t FillMap(THaDetMap *detmap, const char *detectorname)
Int_t fAdcPulseTimeRaw[fMaxAdcChannels]
Definition: THcTrigDet.h:78
Bool_t * fPresentP
Definition: THcTrigDet.h:91
Int_t fAdcPedRaw[fMaxAdcChannels]
Definition: THcTrigDet.h:75
~THcTrigDet()
A destructor.
Definition: THcTrigDet.cxx:139
Int_t fAdcPulseAmpRaw[fMaxAdcChannels]
Definition: THcTrigDet.h:77
const char Option_t
static const int fMaxAdcChannels
Definition: THcTrigDet.h:62
THcRawTdcHit & GetRawTdcHit()
Gets reference to THcRawTdcHit.
void MissReport(const char *name)
Definition: THcHitList.cxx:557
virtual void AddEvtType(int evtype)
Definition: THcTrigDet.cxx:558
Builds a Hall C ENGINE style list of raw hits from raw data.
Definition: THcHitList.h:27
Double_t * fTdcTimeWindowMin
Definition: THcTrigDet.h:72
Int_t fAdcPulseIntRaw[fMaxAdcChannels]
Definition: THcTrigDet.h:76
std::vector< std::string > fTrigNames
Definition: THcTrigDet.h:59
Int_t GetPulseTimeRaw(UInt_t iPulse=0) const
Gets raw pulse time. In subsamples.
int Int_t
bool Bool_t
virtual EStatus Init(const TDatime &date)
Initializes the detector variables.
Definition: THcTrigDet.cxx:148
Double_t * fAdcTimeWindowMax
Definition: THcTrigDet.h:71
UInt_t GetNHits() const
Gets the number of set hits.
Class representing a single raw TDC hit.
Definition: THcRawTdcHit.h:7
TClonesArray * fRawHitList
Definition: THcHitList.h:51
Int_t fCounter
Definition: THcRawHit.h:55
std::vector< std::string > fRFNames
Definition: THcTrigDet.h:60
Double_t fAdcPulseInt[fMaxAdcChannels]
Definition: THcTrigDet.h:81
Int_t fRFId[2]
Definition: THcTrigDet.h:51
const char * Data() const
Int_t fNumAdc
Definition: THcTrigDet.h:48
Int_t fTdcTimeRaw[fMaxTdcChannels]
Definition: THcTrigDet.h:65
static TString Format(const char *fmt,...)
virtual Int_t DefineVariables(EMode mode=kDefine)
Definition: THcTrigDet.cxx:404
Double_t fTdcChanperNS
Definition: THcTrigDet.h:55
Int_t fTDC_RefTimeCut
Definition: THcTrigDet.h:67
virtual void SetEvtType(int evtype)
Definition: THcTrigDet.cxx:562
Double_t GetPed() const
Gets sample pedestal. In channels.
void Error(const char *location, const char *msgfmt,...)
A mock detector to hold trigger related data.
Definition: THcTrigDet.h:16
std::vector< std::string > fTdcNames
Definition: THcTrigDet.h:58
Double_t fAdcPulseAmp[fMaxAdcChannels]
Definition: THcTrigDet.h:82
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
virtual void SetSpectName(const char *name)
Definition: THcTrigDet.cxx:553
Int_t GetPulseAmpRaw(UInt_t iPulse=0) const
Gets raw pulse amplitude. In channels.
Int_t End(THaRunBase *run)
Definition: THcTrigDet.cxx:583
TString fSpectName
Definition: THcTrigDet.h:89
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.
Double_t fAdcPulseTime[fMaxAdcChannels]
Definition: THcTrigDet.h:83
Bool_t HasRefTime() const
Queries whether reference time has been set.
virtual Bool_t IsIgnoreType(Int_t evtype) const
Definition: THcTrigDet.cxx:568
unsigned int UInt_t
char * Form(const char *fmt,...)
void CreateMissReportParms(const char *prefix)
Definition: THcHitList.cxx:544
Int_t fADC_RefTimeCut
Definition: THcTrigDet.h:68
std::vector< Int_t > eventtypes
Definition: THcTrigDet.h:90
virtual Bool_t HaveIgnoreList() const
Definition: THcTrigDet.cxx:577
Class representing a single raw hit for the THcTrigDet.
Definition: THcTrigRawHit.h:10
Double_t * fTdcTimeWindowMax
Definition: THcTrigDet.h:73
tuple list
Definition: SConscript.py:9
const Bool_t kFALSE
Int_t GetPulseIntRaw(UInt_t iPulse=0) const
Gets raw pulse integral. In channels.
R__EXTERN class THcDetectorMap * gHcDetectorMap
Definition: THcGlobals.h:12
THcRawAdcHit & GetRawAdcHit()
Gets reference to THcRawAdcHit.
Int_t fAdcMultiplicity[fMaxAdcChannels]
Definition: THcTrigDet.h:86
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
double Double_t
Double_t fAdcPed[fMaxAdcChannels]
Definition: THcTrigDet.h:80
void TestTime(const Func &f)
Int_t GetTime(UInt_t iHit=0) const
Gets TDC time. In channels.
std::string fKwPrefix
Definition: THcTrigDet.h:46
Int_t Decode(const THaEvData &evData)
Decodes and processes events.
Definition: THcTrigDet.cxx:234
UInt_t GetNPulses() const
Gets number of set pulses.
Int_t fNumTdc
Definition: THcTrigDet.h:49
Double_t * fAdcTimeWindowMin
Definition: THcTrigDet.h:70
Int_t fPlane
Definition: THcRawHit.h:54
UInt_t GetEvType() const
ClassImp(THcDCLookupTTDConv) THcDCLookupTTDConv
Double_t fAdcTdcOffset
Definition: THcTrigDet.h:53
Int_t GetRefTime() const
Gets reference time. In channels.
Int_t GetPedRaw() const
Gets raw signal pedestal. In channels.
std::vector< std::string > fAdcNames
Definition: THcTrigDet.h:57
Int_t fTdcMultiplicity[fMaxTdcChannels]
Definition: THcTrigDet.h:85
Double_t GetPulseTime(UInt_t iPulse=0) const
R__EXTERN class THcParmList * gHcParms
Definition: THcGlobals.h:11
TObject * At(Int_t idx) const
Double_t fTdcOffset
Definition: THcTrigDet.h:54
Double_t fTdcRefTime
Definition: THcTrigDet.h:87
static const int fMaxTdcChannels
Definition: THcTrigDet.h:63
static const Double_t kBig
Definition: THcFormula.cxx:31
const Bool_t kTRUE
Class representing a single raw ADC hit.
Definition: THcRawAdcHit.h:7
void Setup(const char *name, const char *description)
Definition: THcTrigDet.cxx:310
virtual Int_t DecodeToHitList(const THaEvData &evdata, Bool_t suppress=kFALSE)
Populate the hitlist from the raw event data.
Definition: THcHitList.cxx:191
const char * cnt
virtual void Clear(Option_t *opt="")
Clears variables before next event.
Definition: THcTrigDet.cxx:210