Classes | Public Member Functions | Public Attributes | Private Attributes | List of all members
CRT::TwoCRTReco Class Reference
Inheritance diagram for CRT::TwoCRTReco:
art::EDAnalyzer art::detail::Analyzer art::detail::LegacyModule art::Observer art::ModuleBase

Classes

struct  recoHits
 
struct  tempHits
 
struct  tracksPair
 

Public Member Functions

 TwoCRTReco (fhicl::ParameterSet const &p)
 
 TwoCRTReco (TwoCRTReco const &)=delete
 
 TwoCRTReco (TwoCRTReco &&)=delete
 
TwoCRTRecooperator= (TwoCRTReco const &)=delete
 
TwoCRTRecooperator= (TwoCRTReco &&)=delete
 
void analyze (art::Event const &e) override
 
bool moduleMatcher (int module1, int module2)
 
void beginJob () override
 
void endJob () override
 
void createPNG (TH1D *histo)
 
double setAngle (double angle)
 
int moduletoCTB (int module2, int module1)
 
- Public Member Functions inherited from art::EDAnalyzer
 EDAnalyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDAnalyzer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Analyzer
virtual ~Analyzer () noexcept
 
 Analyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 Analyzer (Table< Config > const &config)
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Observer
 ~Observer () noexcept
 
 Observer (Observer const &)=delete
 
 Observer (Observer &&)=delete
 
Observeroperator= (Observer const &)=delete
 
Observeroperator= (Observer &&)=delete
 
void registerProducts (ProductDescriptions &, ModuleDescription const &)
 
void fillDescriptions (ModuleDescription const &)
 
fhicl::ParameterSetID selectorConfig () const
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Public Attributes

std::string fTrackModuleLabel = "pandoraTrack"
 
int nEvents = 0
 
int nHaloMuons = 0
 
int nHitsPerEvent =0
 
int nMCCMuons =0
 

Private Attributes

ofstream logFile
 
art::InputTag fCRTLabel
 
art::InputTag fCTBLabel
 
TTree * fCRTTree
 
int run
 
int subRun
 
bool fMCCSwitch
 
bool fCTBTriggerOnly
 
bool fSCECorrection
 
bool fModuleSwitch
 
int fADCThreshold
 
int fModuletoModuleTimingCut
 
int fFronttoBackTimingCut
 
int adcX_F
 
int adcX_B
 
int adcY_F
 
int adcY_B
 
double X_F
 
double X_B
 
double Y_F
 
double Y_B
 
double Z_F
 
double Z_B
 
int moduleX_F
 
int moduleX_B
 
int moduleY_F
 
int moduleY_B
 
int stripX_F
 
int stripX_B
 
int stripY_F
 
int stripY_B
 
double measuredT0
 
double CRT_TOF
 
long long timeStamp
 
int eventNum
 
std::vector< recoHitsprimaryHits_F
 
std::vector< recoHitsprimaryHits_B
 
std::vector< tempHitstempHits_F
 
std::vector< tempHitstempHits_B
 
std::vector< tracksPairallTracksPair
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
std::string const & processName () const
 
bool wantAllEvents () const noexcept
 
bool wantEvent (ScheduleID id, Event const &e) const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) const
 
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &select_paths, std::vector< std::string > const &reject_paths, fhicl::ParameterSet const &config)
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 87 of file TwoCRTReco_module.cc.

Constructor & Destructor Documentation

CRT::TwoCRTReco::TwoCRTReco ( fhicl::ParameterSet const &  p)
explicit

Definition at line 210 of file TwoCRTReco_module.cc.

211  :
212  EDAnalyzer(p), fCRTLabel(p.get < art::InputTag > ("CRTLabel")), fCTBLabel(p.get<art::InputTag>("CTBLabel")) {
213  consumes < std::vector < CRT::Trigger >> (fCRTLabel);
214  consumes < std::vector < art::Assns < sim::AuxDetSimChannel, CRT::Trigger >>> (fCRTLabel);
215  fMCCSwitch=(p.get<bool>("MCC"));
216  fCTBTriggerOnly=(p.get<bool>("CTBOnly"));
217  fSCECorrection=(p.get<bool>("SCECorrection"));
218  }
art::InputTag fCRTLabel
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
p
Definition: test.py:223
art::InputTag fCTBLabel
CRT::TwoCRTReco::TwoCRTReco ( TwoCRTReco const &  )
delete
CRT::TwoCRTReco::TwoCRTReco ( TwoCRTReco &&  )
delete

Member Function Documentation

void CRT::TwoCRTReco::analyze ( art::Event const &  e)
overridevirtual

Implements art::EDAnalyzer.

Definition at line 291 of file TwoCRTReco_module.cc.

294 {
295 
296  nEvents++;
297  eventNum=event.event();
298  run=event.run();
299  subRun=event.subRun();
300  if (fMCCSwitch){
301  fModuleSwitch=1;
302  fADCThreshold=800;
305 
306 }
307  else {
308  fModuleSwitch=0;
309  fADCThreshold=20;
312  art::ValidHandle<std::vector<raw::RDTimeStamp>> timingHandle = event.getValidHandle<std::vector<raw::RDTimeStamp>>("timingrawdecoder:daq");
313  timeStamp=timingHandle->at(0).GetTimeStamp();
314 
315 }
316 
317 
318 
319  if(!fMCCSwitch){
320  //const auto& pdspctbs = event.getValidHandle<std::vector<raw::ctb::pdspctb>>(fCTB_tag);
321  art::ValidHandle<std::vector<raw::RDTimeStamp>> timingHandle = event.getValidHandle<std::vector<raw::RDTimeStamp>>("timingrawdecoder:daq");
322 
323  const raw::RDTimeStamp& timeStamp = timingHandle->at(0);
324  if (fCTBTriggerOnly){
325  if(timeStamp.GetFlags()!= 13) return;}
326  }
327  int nHits = 0;
328 
329 
330 
331 
332  //Detector properties service
333  primaryHits_F.clear();
334  primaryHits_B.clear();
335  allTracksPair.clear();
336  tempHits_F.clear();
337  tempHits_B.clear(); // Arrays to compile hits and move them through
338 
339 
340  //Get triggers
341  cout << "Getting triggers" << endl;
342  const auto & triggers = event.getValidHandle < std::vector < CRT::Trigger >> (fCRTLabel);
343 
344  art::FindManyP < sim::AuxDetSimChannel > trigToSim(triggers, event, fCRTLabel);
345 
346  //Get a handle to the Geometry service to look up AuxDetGeos from module numbers
348 
349 
350 
351  //Mapping from channel to trigger
352  std::unordered_map < size_t, double > prevTimes;
353  int hitID = 0;
354  cout << "Looking for hits in Triggers" << endl;
355 
356  // Find CTB pixels
357  int pixel0 = -1;
358  int pixel1 = -1;
359  if (!fMCCSwitch)
360  {
361  const auto& pdspctbs = *event.getValidHandle<std::vector<raw::ctb::pdspctb>>(fCTBLabel);
362  std::vector<int> uS, dS;
363  const size_t npdspctbs = pdspctbs.size();
364  for(size_t j=0;j<npdspctbs;++j)
365  {
366  const std::vector<raw::ctb::Trigger> HLTriggers = pdspctbs[j].GetHLTriggers();
367  const std::vector<raw::ctb::ChStatus> chs = pdspctbs[j].GetChStatusAfterHLTs();
368 for (size_t k=0; k<HLTriggers.size(); ++k)
369  {
370  //cout<<chs[k].timestamp<<endl;
371  int num = chs[k].crt;
372  //cout<<num<<endl;
373 
374  const std::string binary = std::bitset<32>(num).to_string();
375  const auto crtmask=chs[k].crt;
376  pixel0 = -1;
377  pixel1 = -1;
378  //cout<<crtmask<<endl;
379  for (int i = 0; i<32; ++i){
380  if (crtmask & (1<<i)){
381  if (i<16){
382  pixel0 = i;
383  }
384  else {
385  pixel1 = i;
386  }
387  }
388  }
389  if (pixel0!=-1 && pixel1!=-1) {
390  //cout<<nEvents<<" TJYang Pixels: "<<pixel0<<","<<pixel1<<endl;
391  }
392  else if (fCTBTriggerOnly) return;
393  }
394  }
395  }
396 
397  int trigID=0;
398  for (const auto & trigger: * triggers) {
399  const auto & hits = trigger.Hits();
400  for (const auto & hit: hits) { // Collect hits on all modules
401  //cout<<hits.size()<<','<<hit.ADC()<<endl;
402  if (hit.ADC() > fADCThreshold) { // Keep if they are above threshold
403 
404  tempHits tHits;
405  if (!fMCCSwitch){
406 
407  tHits.module = trigger.Channel(); // Values to add to array
408  tHits.channel=hit.Channel();
409  tHits.adc = hit.ADC();
410  tHits.triggerTime=trigger.Timestamp()-timeStamp;
411  }
412  else{
413  tHits.module = trigger.Channel(); // Values to add to array
414  tHits.channel=hit.Channel();
415  tHits.adc = hit.ADC();
416  tHits.triggerTime=trigger.Timestamp();
417  }
418  //cout<<trigger.Channel()<<','<<hit.Channel()<<','<<hit.ADC()<<endl;
419  nHits++;
420  tHits.triggerNumber=trigID;
421  const auto & trigGeo = geom -> AuxDet(trigger.Channel()); // Get geo
422  const auto & csens = trigGeo.SensitiveVolume(hit.Channel());
423  const auto center = csens.GetCenter();
424  if (center.Z() < 100) tempHits_F.push_back(tHits); // Sort F/B from Z
425  else tempHits_B.push_back(tHits);
426  hitID++;
427  }
428  }
429  trigID++;
430  }
431  nHitsPerEvent=nHits;
432  cout << "Hits compiled for event: " << nEvents << endl;
433  cout << "Number of Hits above Threshold: " << hitID << endl;
434 
435  for (unsigned int f = 0; f < tempHits_F.size(); f++) {
436  for (unsigned int f_test = 0; f_test < tempHits_F.size(); f_test++) {
437  if (fabs(tempHits_F[f_test].triggerTime-tempHits_F[f].triggerTime)>fModuletoModuleTimingCut) continue;
438  const auto & trigGeo = geom -> AuxDet(tempHits_F[f].module);
439  const auto & trigGeo2 = geom -> AuxDet(tempHits_F[f_test].module);
440 
441  const auto & hit1Geo = trigGeo.SensitiveVolume(tempHits_F[f].channel);
442  const auto hit1Center = hit1Geo.GetCenter();
443  // Create 2D hits from geo of the Y and X modules
444  const auto & hit2Geo = trigGeo2.SensitiveVolume(tempHits_F[f_test].channel);
445  const auto hit2Center = hit2Geo.GetCenter();
446  bool moduleMatched;
447  moduleMatched=moduleMatcher(tempHits_F[f_test].module, tempHits_F[f].module);
448  if (moduleMatched) {
449  // Get the center of the hits (CRT_Res=2.5 cm)
450  double hitX = hit1Center.X();
451  for (unsigned int a = 0; a < tempHits_F.size(); a++)
452  {
453  if(tempHits_F[a].module==tempHits_F[f].module && (tempHits_F[a].channel-1)==tempHits_F[f].channel) hitX=hit1Center.X()+1.25;
454  }
455  double hitYPrelim=hit2Center.Y();
456  for (unsigned int a = 0; a < tempHits_F.size(); a++)
457  {
458  if(tempHits_F[a].module==tempHits_F[f_test].module && (tempHits_F[a].channel-1)==tempHits_F[f_test].channel) hitYPrelim=hit2Center.Y()+1.25;
459  }
460 
461 
462 
463  double hitY=hitYPrelim;
464  double hitZ = (hit1Center.Z() + hit2Center.Z()) / 2.f;
465 
466  recoHits rHits;
467  rHits.hitPositionX = hitX;
468  rHits.hitPositionY = hitY;
469  rHits.hitPositionZ = hitZ;
470  rHits.geoX=tempHits_F[f].module;
471  rHits.geoY=tempHits_F[f_test].module;
472  rHits.stripX=tempHits_F[f].channel;
473  rHits.stripY=tempHits_F[f_test].channel;
474 
475  rHits.adcX=tempHits_F[f].adc;
476  rHits.adcY=tempHits_F[f_test].adc;
477  rHits.trigNumberX=tempHits_F[f].triggerNumber;
478  rHits.trigNumberY=tempHits_F[f_test].triggerNumber;
479  rHits.timeAvg = (tempHits_F[f_test].triggerTime+tempHits_F[f].triggerTime)/2.0;
480  primaryHits_F.push_back(rHits); // Add array
481  }
482  }
483  }
484  for (unsigned int f = 0; f < tempHits_B.size(); f++) {
485  for (unsigned int f_test = 0; f_test < tempHits_B.size(); f_test++) { // Same as above but for back CRT
486  if (fabs(tempHits_B[f_test].triggerTime-tempHits_B[f].triggerTime)>fModuletoModuleTimingCut) continue;
487 
488  const auto & trigGeo = geom -> AuxDet(tempHits_B[f].module);
489  const auto & trigGeo2 = geom -> AuxDet(tempHits_B[f_test].module);
490  const auto & hit1Geo = trigGeo.SensitiveVolume(tempHits_B[f].channel);
491  const auto hit1Center = hit1Geo.GetCenter();
492 
493  const auto & hit2Geo = trigGeo2.SensitiveVolume(tempHits_B[f_test].channel);
494  const auto hit2Center = hit2Geo.GetCenter();
495  bool moduleMatched;
496  moduleMatched=moduleMatcher(tempHits_B[f_test].module, tempHits_B[f].module);
497 
498  if (moduleMatched) {
499  double hitX = hit1Center.X();
500 
501 
502  for (unsigned int a = 0; a < tempHits_B.size(); a++)
503  {
504  if(tempHits_B[a].module==tempHits_B[f].module && (tempHits_B[a].channel-1)==tempHits_B[f].channel) hitX=hit1Center.X()+1.25;
505  }
506 
507  double hitYPrelim = hit2Center.Y();
508 
509  for (unsigned int a = 0; a < tempHits_B.size(); a++)
510  {
511  if(tempHits_B[a].module==tempHits_B[f_test].module && (tempHits_B[a].channel-1)==tempHits_B[f_test].channel) hitYPrelim=hit2Center.Y()+1.25;
512  }
513  double hitY=hitYPrelim;
514 
515 
516  double hitZ = (hit1Center.Z() + hit2Center.Z()) / 2.f;
517 
518  recoHits rHits;
519  rHits.hitPositionX = hitX;
520  rHits.hitPositionY = hitY;
521  rHits.hitPositionZ = hitZ;
522  rHits.geoX=tempHits_B[f].module;
523  rHits.geoY=tempHits_B[f_test].module;
524  rHits.adcX=tempHits_B[f].adc;
525  rHits.adcY=tempHits_B[f_test].adc;
526  rHits.stripX=tempHits_B[f].channel;
527  rHits.stripY=tempHits_B[f_test].channel;
528  rHits.trigNumberX=tempHits_B[f].triggerNumber;
529  rHits.trigNumberY=tempHits_B[f_test].triggerNumber;
530  rHits.timeAvg = (tempHits_B[f_test].triggerTime+tempHits_B[f].triggerTime)/2.0;
531  primaryHits_B.push_back(rHits);
532 
533  }
534  }
535  }
536 
537  std::cout<<"Number of Hits: "<<primaryHits_F.size()<<','<<primaryHits_B.size()<<std::endl;
538 
539  std::cout<<"CTB Pixes: "<<pixel0<<','<<pixel1<<std::endl;
540  int tempId = 0;
541  adcX_F=0;
542  adcX_B=0;
543  adcY_F=0;
544  adcY_B=0;
545  int crtPixel0=-1;
546  int crtPixel1=-1;
547  for (unsigned int f = 0; f < primaryHits_F.size(); f++) {
548  if (pixel0==-1 || pixel1==-1) break;
549  for (unsigned int b = 0; b < primaryHits_B.size(); b++) {
550  if (fabs(primaryHits_F[f].timeAvg-primaryHits_B[b].timeAvg)>fFronttoBackTimingCut) continue;
551  crtPixel0=moduletoCTB(primaryHits_F[f].geoX,primaryHits_F[f].geoY);
552  crtPixel1=moduletoCTB(primaryHits_B[b].geoX, primaryHits_B[b].geoY);
553  if (crtPixel0!=pixel0 || crtPixel1!=pixel1) continue;
554  //std::cout<<"HEY"<<std::endl;
555  if (adcX_F<primaryHits_F[f].adcX && adcY_F<primaryHits_F[f].adcY && adcX_B<primaryHits_B[b].adcX && adcY_B<primaryHits_B[b].adcY){
556  X_F = primaryHits_F[f].hitPositionX;
557  Y_F = primaryHits_F[f].hitPositionY;
558  Z_F = primaryHits_F[f].hitPositionZ;
559  X_B = primaryHits_B[b].hitPositionX;
560  Y_B = primaryHits_B[b].hitPositionY;
561  Z_B= primaryHits_B[b].hitPositionZ;
562  adcX_F=primaryHits_F[f].adcX;
563  adcX_B=primaryHits_B[b].adcX;
564  adcY_F=primaryHits_F[f].adcY;
565  adcY_B=primaryHits_B[b].adcY;
566  moduleX_F=primaryHits_F[f].geoX;
567  moduleY_F=primaryHits_F[f].geoY;
568  moduleX_B=primaryHits_B[b].geoX;
569  moduleY_B=primaryHits_B[b].geoY;
570 
571  stripX_F=primaryHits_F[f].stripX;
572  stripY_F=primaryHits_F[f].stripY;
573  stripX_B=primaryHits_B[b].stripX;
574  stripY_B=primaryHits_B[b].stripY;
575  //std::cout<<"FOUND A COMBO"<<std::endl;
576  measuredT0=(primaryHits_F[f].timeAvg+primaryHits_B[b].timeAvg)/2.f;
577  CRT_TOF=(primaryHits_F[f].timeAvg-primaryHits_B[b].timeAvg);
578  tempId++;
579 
580  }
581  }
582  }
583  if (adcX_F>0) fCRTTree->Fill();
584  // Filter return if (adcX_F==0) return false;
585  // Filter return
586  // return true;
587  }
std::vector< recoHits > primaryHits_B
std::string string
Definition: nybbler.cc:12
art::InputTag fCRTLabel
std::vector< tempHits > tempHits_B
std::vector< recoHits > primaryHits_F
uint8_t channel
Definition: CRTFragment.hh:201
std::vector< tempHits > tempHits_F
int moduletoCTB(int module2, int module1)
const double a
std::vector< tracksPair > allTracksPair
Detector simulation of raw signals on wires.
bool moduleMatcher(int module1, int module2)
def center(depos, point)
Definition: depos.py:117
art::InputTag fCTBLabel
static bool * b
Definition: config.cpp:1043
std::string to_string(ModuleType const mt)
Definition: ModuleType.h:34
QTextStream & endl(QTextStream &s)
Event finding and building.
void CRT::TwoCRTReco::beginJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 591 of file TwoCRTReco_module.cc.

591  {
592  art::ServiceHandle<art::TFileService> fileServiceHandle;
593  fCRTTree = fileServiceHandle->make<TTree>("CRTCandidate", "track by track info");
594 
595  fCRTTree->Branch("nEvents", &nEvents, "fnEvents/I");
596  fCRTTree->Branch("hRun", &run, "run/I");
597  fCRTTree->Branch("hEvent",&eventNum,"eventNum/I");
598 
599 
600  fCRTTree->Branch("hX_F", &X_F, "X_F/D");
601  fCRTTree->Branch("hX_B", &X_B, "X_B/D");
602  fCRTTree->Branch("hY_F", &Y_F, "Y_F/D");
603  fCRTTree->Branch("hY_B", &Y_B, "Y_B/D");
604  fCRTTree->Branch("hZ_F", &Z_F, "Z_F/D");
605  fCRTTree->Branch("hZ_B", &Z_B, "Z_B/D");
606 
607 
608  fCRTTree->Branch("CRT_TOF", &CRT_TOF, "CRT_TOF/D");
609 
610  fCRTTree->Branch("hmoduleX_F", &moduleX_F, "moduleX_F/I");
611  fCRTTree->Branch("hmoduleX_B", &moduleX_B, "moduleX_B/I");
612  fCRTTree->Branch("hmoduleY_F", &moduleY_F, "moduleY_F/I");
613  fCRTTree->Branch("hmoduleY_B", &moduleY_B, "moduleY_B/I");
614 
615  fCRTTree->Branch("hstripX_F", &stripX_F, "stripX_F/I");
616  fCRTTree->Branch("hstripX_B", &stripX_B, "stripX_B/I");
617  fCRTTree->Branch("hstripY_F", &stripY_F, "stripY_F/I");
618  fCRTTree->Branch("hstripY_B", &stripY_B, "stripY_B/I");
619 
620  fCRTTree->Branch("hadcX_F", &adcX_F, "adcX_F/I");
621  fCRTTree->Branch("hadcY_F", &adcY_F, "adcY_F/I");
622  fCRTTree->Branch("hadcX_B", &adcX_B, "adcX_B/I");
623  fCRTTree->Branch("hadcY_B", &adcY_B, "adcY_B/I");
624 
625  fCRTTree->Branch("hmeasuredT0", &measuredT0, "measuredT0/D");
626 
627 
628 
629 
630 
631 }
void CRT::TwoCRTReco::createPNG ( TH1D *  histo)

Definition at line 269 of file TwoCRTReco_module.cc.

269  {
270  //Save important histograms as PNG
271  TCanvas * c = new TCanvas;
272  TH1D * h = histo;
273  h -> Draw();
274  TImage * img = TImage::Create();
275  img -> FromPad(c);
276  img -> WriteImage((std::string(histo -> GetName()) + ".png").c_str());
277  delete h;
278  delete c;
279  delete img;
280 
281 }
std::string string
Definition: nybbler.cc:12
void Draw(const char *plot, const char *title)
Definition: gXSecComp.cxx:580
void CRT::TwoCRTReco::endJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 633 of file TwoCRTReco_module.cc.

634 {
635 
636 
637 
638 }
bool CRT::TwoCRTReco::moduleMatcher ( int  module1,
int  module2 
)

Definition at line 224 of file TwoCRTReco_module.cc.

224  {
225  // Function checking if two hits could reasonably be matched into a 2D hit
226  if ((module1 == 6 && (module2 == 10 || module2 == 11)) || (module1 == 14 && (module2 == 10 || module2 == 11)) || (module1 == 19 && (module2 == 26 || module2 == 27)) || (module1 == 31 && (module2 == 26 || module2 == 27)) || (module1 == 7 && (module2 == 12 || module2 == 13)) || (module1 == 15 && (module2 == 12 || module2 == 13)) || (module1 == 18 && (module2 == 24 || module2 == 25)) || (module1 == 30 && (module2 == 24 || module2 == 25)) || (module1 == 1 && (module2 == 4 || module2 == 5)) || (module1 == 9 && (module2 == 4 || module2 == 5)) || (module1 == 16 && (module2 == 20 || module2 == 21)) || (module1 == 28 && (module2 == 20 || module2 == 21)) || (module1 == 0 && (module2 == 2 || module2 == 3)) || (module1 == 8 && (module2 == 2 || module2 == 3)) || (module1 == 17 && (module2 == 22 || module2 == 23)) || (module1 == 29 && (module2 == 22 || module2 == 23))) return 1;
227  else return 0;
228 
229 }
int CRT::TwoCRTReco::moduletoCTB ( int  module2,
int  module1 
)

Definition at line 232 of file TwoCRTReco_module.cc.

232  {
233  if (module1 == 15 && module2 == 12 ) return 15;
234  else if (module1 == 15 && module2 == 13) return 10;
235  else if (module1 == 7 && module2 == 12) return 8;
236  else if (module1 == 7 && module2 == 13) return 9;
237  else if (module1 == 19 && module2 == 26) return 4;
238  else if (module1 == 19 && module2 == 27) return 13;
239  else if (module1 == 31 && module2 == 26) return 3;
240  else if (module1 == 31 && module2 == 27) return 2;
241  else if (module1 == 30 && module2 == 24) return 1;
242  else if (module1 == 30 && module2 == 25) return 0;
243  else if (module1 == 18 && module2 == 24) return 12;
244  else if (module1 == 18 && module2 == 25) return 11;
245  else if (module1 == 6 && module2 == 11) return 7;
246  else if (module1 == 6 && module2 == 10) return 6;
247  else if (module1 == 14 && module2 == 11) return 14;
248  else if (module1 == 14 && module2 == 10) return 5;
249  else if (module1 == 0 && module2 == 3) return 25;
250  else if (module1 == 0 && module2 == 2) return 24;
251  else if (module1 == 8 && module2 == 3) return 26;
252  else if (module1 == 8 && module2 == 2) return 30;
253  else if (module1 == 17 && module2 == 23) return 27;
254  else if (module1 == 17 && module2 == 22) return 28;
255  else if (module1 == 29 && module2 == 23) return 16;
256  else if (module1 == 29 && module2 == 22) return 17;
257  else if (module1 == 28 && module2 == 21) return 18;
258  else if (module1 == 28 && module2 == 20) return 19;
259  else if (module1 == 16 && module2 == 21) return 29;
260  else if (module1 == 16 && module2 == 20) return 20;
261  else if (module1 == 9 && module2 == 5) return 31;
262  else if (module1 == 9 && module2 == 4) return 21;
263  else if (module1 == 1 && module2 == 5) return 23;
264  else if (module1 == 1 && module2 == 4) return 22;
265  else return -1;
266 }
TwoCRTReco& CRT::TwoCRTReco::operator= ( TwoCRTReco const &  )
delete
TwoCRTReco& CRT::TwoCRTReco::operator= ( TwoCRTReco &&  )
delete
double CRT::TwoCRTReco::setAngle ( double  angle)

Definition at line 282 of file TwoCRTReco_module.cc.

282  {
283  if (angle < 0) {
284  angle += 3.14159265359;
285  }
286  angle *= 180.00 / 3.14159265359;
287  return angle;
288 }

Member Data Documentation

int CRT::TwoCRTReco::adcX_B
private

Definition at line 134 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::adcX_F
private

Definition at line 134 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::adcY_B
private

Definition at line 134 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::adcY_F
private

Definition at line 134 of file TwoCRTReco_module.cc.

std::vector< tracksPair > CRT::TwoCRTReco::allTracksPair
private

Definition at line 206 of file TwoCRTReco_module.cc.

double CRT::TwoCRTReco::CRT_TOF
private

Definition at line 140 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::eventNum
private

Definition at line 142 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::fADCThreshold
private

Definition at line 127 of file TwoCRTReco_module.cc.

art::InputTag CRT::TwoCRTReco::fCRTLabel
private

Definition at line 118 of file TwoCRTReco_module.cc.

TTree* CRT::TwoCRTReco::fCRTTree
private

Definition at line 120 of file TwoCRTReco_module.cc.

art::InputTag CRT::TwoCRTReco::fCTBLabel
private

Definition at line 119 of file TwoCRTReco_module.cc.

bool CRT::TwoCRTReco::fCTBTriggerOnly
private

Definition at line 124 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::fFronttoBackTimingCut
private

Definition at line 129 of file TwoCRTReco_module.cc.

bool CRT::TwoCRTReco::fMCCSwitch
private

Definition at line 123 of file TwoCRTReco_module.cc.

bool CRT::TwoCRTReco::fModuleSwitch
private

Definition at line 126 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::fModuletoModuleTimingCut
private

Definition at line 128 of file TwoCRTReco_module.cc.

bool CRT::TwoCRTReco::fSCECorrection
private

Definition at line 125 of file TwoCRTReco_module.cc.

std::string CRT::TwoCRTReco::fTrackModuleLabel = "pandoraTrack"

Definition at line 91 of file TwoCRTReco_module.cc.

ofstream CRT::TwoCRTReco::logFile
private

Definition at line 115 of file TwoCRTReco_module.cc.

double CRT::TwoCRTReco::measuredT0
private

Definition at line 139 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::moduleX_B
private

Definition at line 137 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::moduleX_F
private

Definition at line 137 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::moduleY_B
private

Definition at line 137 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::moduleY_F
private

Definition at line 137 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::nEvents = 0

Definition at line 111 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::nHaloMuons = 0

Definition at line 112 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::nHitsPerEvent =0

Definition at line 113 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::nMCCMuons =0

Definition at line 114 of file TwoCRTReco_module.cc.

std::vector< recoHits > CRT::TwoCRTReco::primaryHits_B
private

Definition at line 202 of file TwoCRTReco_module.cc.

std::vector< recoHits > CRT::TwoCRTReco::primaryHits_F
private

Definition at line 201 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::run
private

Definition at line 122 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::stripX_B
private

Definition at line 138 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::stripX_F
private

Definition at line 138 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::stripY_B
private

Definition at line 138 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::stripY_F
private

Definition at line 138 of file TwoCRTReco_module.cc.

int CRT::TwoCRTReco::subRun
private

Definition at line 122 of file TwoCRTReco_module.cc.

std::vector< tempHits > CRT::TwoCRTReco::tempHits_B
private

Definition at line 205 of file TwoCRTReco_module.cc.

std::vector< tempHits > CRT::TwoCRTReco::tempHits_F
private

Definition at line 204 of file TwoCRTReco_module.cc.

long long CRT::TwoCRTReco::timeStamp
private

Definition at line 141 of file TwoCRTReco_module.cc.

double CRT::TwoCRTReco::X_B
private

Definition at line 135 of file TwoCRTReco_module.cc.

double CRT::TwoCRTReco::X_F
private

Definition at line 135 of file TwoCRTReco_module.cc.

double CRT::TwoCRTReco::Y_B
private

Definition at line 135 of file TwoCRTReco_module.cc.

double CRT::TwoCRTReco::Y_F
private

Definition at line 135 of file TwoCRTReco_module.cc.

double CRT::TwoCRTReco::Z_B
private

Definition at line 135 of file TwoCRTReco_module.cc.

double CRT::TwoCRTReco::Z_F
private

Definition at line 135 of file TwoCRTReco_module.cc.


The documentation for this class was generated from the following file: