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

Classes

struct  recoHits
 
struct  removePairIndex
 
struct  sortPair
 
struct  tempHits
 
struct  tracksPair
 

Public Member Functions

 TwoCRTMatching (fhicl::ParameterSet const &p)
 
 TwoCRTMatching (TwoCRTMatching const &)=delete
 
 TwoCRTMatching (TwoCRTMatching &&)=delete
 
TwoCRTMatchingoperator= (TwoCRTMatching const &)=delete
 
TwoCRTMatchingoperator= (TwoCRTMatching &&)=delete
 
void analyze (art::Event const &e) override
 
bool moduleMatcher (int module1, int module2)
 
double signedPointToLineDistance (double firstPoint1, double firstPoint2, double secondPoint1, double secondPoint2, double trackPoint1, double trackPoint2)
 
double signed3dDistance (double firstPoint1, double firstPoint2, double firstPoint3, double secondPoint1, double secondPoint2, double secondPoint3, TVector3 trackPos)
 
void beginJob () override
 
void endJob () override
 
void createPNG (TH1D *histo)
 
double setAngle (double angle)
 
- 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
 
TTree * fCRTdQTree
 
TTree * fMCCMuon
 
TTree * fTrackInfo
 
int run
 
int subRun
 
bool fMCCSwitch
 
bool fCTBTriggerOnly
 
bool fSCECorrection
 
bool fModuleSwitch
 
int fADCThreshold
 
int fModuletoModuleTimingCut
 
int fFronttoBackTimingCut
 
double averageSignedDistanceYZ
 
double averageSignedDistanceXZ
 
double averageSignedDistance
 
double displAngleXZ
 
double displAngleYZ
 
double displAngleXY
 
double CRT_TOF
 
double deltaX_F
 
double deltaX_B
 
double deltaY_F
 
double deltaY_B
 
double dotCos
 
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
 
double trackX1
 
double trackX2
 
double trackY1
 
double trackY2
 
double trackZ1
 
double trackZ2
 
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 measuredXOffset
 
bool recoPandoraT0Check
 
int trackID
 
double truthEnergy
 
int mccTrackId
 
double mccT0
 
double SCECorrectX_F
 
double SCECorrectY_F
 
double SCECorrectZ_F
 
double SCECorrectX_B
 
double SCECorrectY_B
 
double SCECorrectZ_B
 
int endTPC
 
double trkPosX
 
double trkPosY
 
double trkPosZ
 
double displXZ
 
double displYZ
 
double mccCRTStartX
 
double mccCRTStartY
 
double mccCRTStartZ
 
double mccCRTEndX
 
double mccCRTEndY
 
double mccCRTEndZ
 
int candidateCRT
 
double trkhitx
 
double trkhity
 
double trkhitz
 
double trkhitt0
 
double crtt0
 
double trkdqdx
 
int trkhitIntegral
 
double sigmaHit
 
int TPCID
 
int WireID
 
int sumADC
 
int rangeTime
 
long long timeStamp
 
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 86 of file TwoCRTMatching_module.cc.

Constructor & Destructor Documentation

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

Definition at line 274 of file TwoCRTMatching_module.cc.

275  :
276  EDAnalyzer(p), fCRTLabel(p.get < art::InputTag > ("CRTLabel")), fCTBLabel(p.get<art::InputTag>("CTBLabel")) {
277  consumes < std::vector < CRT::Trigger >> (fCRTLabel);
278  consumes < std::vector < art::Assns < sim::AuxDetSimChannel, CRT::Trigger >>> (fCRTLabel);
279  fMCCSwitch=(p.get<bool>("MCC"));
280  fCTBTriggerOnly=(p.get<bool>("CTBOnly"));
281  fSCECorrection=(p.get<bool>("SCECorrection"));
282  }
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
p
Definition: test.py:223
CRT::TwoCRTMatching::TwoCRTMatching ( TwoCRTMatching const &  )
delete
CRT::TwoCRTMatching::TwoCRTMatching ( TwoCRTMatching &&  )
delete

Member Function Documentation

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

Implements art::EDAnalyzer.

Definition at line 346 of file TwoCRTMatching_module.cc.

349 {
350 
351  nEvents++;
352  run=event.run();
353  subRun=event.subRun();
354  mccTrackId=-1;
355  if (fMCCSwitch){
356  fModuleSwitch=1;
357  fADCThreshold=800;
360 
361 }
362  else {
363  fModuleSwitch=0;
364  fADCThreshold=20;
367  art::ValidHandle<std::vector<raw::RDTimeStamp>> timingHandle = event.getValidHandle<std::vector<raw::RDTimeStamp>>("timingrawdecoder:daq");
368  timeStamp=timingHandle->at(0).GetTimeStamp();
369 
370 }
371 
372 
373  auto const* SCE = lar::providerFrom<spacecharge::SpaceChargeService>();
374 
375 
376  if(!fMCCSwitch){
377  //const auto& pdspctbs = event.getValidHandle<std::vector<raw::ctb::pdspctb>>(fCTB_tag);
378  art::ValidHandle<std::vector<raw::RDTimeStamp>> timingHandle = event.getValidHandle<std::vector<raw::RDTimeStamp>>("timingrawdecoder:daq");
379 
380  const raw::RDTimeStamp& timeStamp = timingHandle->at(0);
381  if (fCTBTriggerOnly){
382  if(timeStamp.GetFlags()!= 13) return;}
383  }
384  int nHits = 0;
385 
386 
387 
388 
389  //Detector properties service
390  primaryHits_F.clear();
391  primaryHits_B.clear();
392  allTracksPair.clear();
393  tempHits_F.clear();
394  tempHits_B.clear(); // Arrays to compile hits and move them through
395 
396 
397  //Get triggers
398  cout << "Getting triggers" << endl;
399  const auto & triggers = event.getValidHandle < std::vector < CRT::Trigger >> (fCRTLabel);
400 
401  art::FindManyP < sim::AuxDetSimChannel > trigToSim(triggers, event, fCRTLabel);
402 
403  //Get a handle to the Geometry service to look up AuxDetGeos from module numbers
405 
406 
407 
408  //Mapping from channel to trigger
409  std::unordered_map < size_t, double > prevTimes;
410  int hitID = 0;
411  cout << "Looking for hits in Triggers" << endl;
412 
413  int trigID=0;
414  for (const auto & trigger: * triggers) {
415  const auto & hits = trigger.Hits();
416  for (const auto & hit: hits) { // Collect hits on all modules
417  //cout<<hits.size()<<','<<hit.ADC()<<endl;
418  if (hit.ADC() > fADCThreshold) { // Keep if they are above threshold
419 
420  tempHits tHits;
421  if (!fMCCSwitch){
422 
423  tHits.module = trigger.Channel(); // Values to add to array
424  tHits.channel=hit.Channel();
425  tHits.adc = hit.ADC();
426  tHits.triggerTime=trigger.Timestamp()-timeStamp;
427  }
428  else{
429  tHits.module = trigger.Channel(); // Values to add to array
430  tHits.channel=hit.Channel();
431  tHits.adc = hit.ADC();
432  tHits.triggerTime=trigger.Timestamp();
433  }
434  //cout<<trigger.Channel()<<','<<hit.Channel()<<','<<hit.ADC()<<endl;
435  nHits++;
436  tHits.triggerNumber=trigID;
437  const auto & trigGeo = geom -> AuxDet(trigger.Channel()); // Get geo
438  const auto & csens = trigGeo.SensitiveVolume(hit.Channel());
439  const auto center = csens.GetCenter();
440  if (center.Z() < 100) tempHits_F.push_back(tHits); // Sort F/B from Z
441  else tempHits_B.push_back(tHits);
442  hitID++;
443  }
444  }
445  trigID++;
446  }
447  nHitsPerEvent=nHits;
448  cout << "Hits compiled for event: " << nEvents << endl;
449  cout << "Number of Hits above Threshold: " << hitID << endl;
450 
451  for (unsigned int f = 0; f < tempHits_F.size(); f++) {
452  for (unsigned int f_test = 0; f_test < tempHits_F.size(); f_test++) {
453  if (fabs(tempHits_F[f_test].triggerTime-tempHits_F[f].triggerTime)>fModuletoModuleTimingCut) continue;
454  const auto & trigGeo = geom -> AuxDet(tempHits_F[f].module);
455  const auto & trigGeo2 = geom -> AuxDet(tempHits_F[f_test].module);
456 
457  const auto & hit1Geo = trigGeo.SensitiveVolume(tempHits_F[f].channel);
458  const auto hit1Center = hit1Geo.GetCenter();
459  // Create 2D hits from geo of the Y and X modules
460  const auto & hit2Geo = trigGeo2.SensitiveVolume(tempHits_F[f_test].channel);
461  const auto hit2Center = hit2Geo.GetCenter();
462  bool moduleMatched;
463  moduleMatched=moduleMatcher(tempHits_F[f_test].module, tempHits_F[f].module);
464  if (moduleMatched) {
465  // Get the center of the hits (CRT_Res=2.5 cm)
466  double hitX = hit1Center.X();
467  for (unsigned int a = 0; a < tempHits_F.size(); a++)
468  {
469  if(tempHits_F[a].module==tempHits_F[f].module && (tempHits_F[a].channel-1)==tempHits_F[f].channel) hitX=hit1Center.X()+1.25;
470  }
471  double hitYPrelim=hit2Center.Y();
472  for (unsigned int a = 0; a < tempHits_F.size(); a++)
473  {
474  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;
475  }
476 
477 
478 
479  double hitY=hitYPrelim;
480  double hitZ = (hit1Center.Z() + hit2Center.Z()) / 2.f;
481 
482  recoHits rHits;
483  rHits.hitPositionX = hitX;
484  rHits.hitPositionY = hitY;
485  rHits.hitPositionZ = hitZ;
486  rHits.geoX=tempHits_F[f].module;
487  rHits.geoY=tempHits_F[f_test].module;
488  rHits.stripX=tempHits_F[f].channel;
489  rHits.stripY=tempHits_F[f_test].channel;
490 
491  rHits.adcX=tempHits_F[f].adc;
492  rHits.adcY=tempHits_F[f_test].adc;
493  rHits.trigNumberX=tempHits_F[f].triggerNumber;
494  rHits.trigNumberY=tempHits_F[f_test].triggerNumber;
495  rHits.timeAvg = (tempHits_F[f_test].triggerTime+tempHits_F[f].triggerTime)/2.0;
496  primaryHits_F.push_back(rHits); // Add array
497  }
498  }
499  }
500  for (unsigned int f = 0; f < tempHits_B.size(); f++) {
501  for (unsigned int f_test = 0; f_test < tempHits_B.size(); f_test++) { // Same as above but for back CRT
502  if (fabs(tempHits_B[f_test].triggerTime-tempHits_B[f].triggerTime)>fModuletoModuleTimingCut) continue;
503 
504  const auto & trigGeo = geom -> AuxDet(tempHits_B[f].module);
505  const auto & trigGeo2 = geom -> AuxDet(tempHits_B[f_test].module);
506  const auto & hit1Geo = trigGeo.SensitiveVolume(tempHits_B[f].channel);
507  const auto hit1Center = hit1Geo.GetCenter();
508 
509  const auto & hit2Geo = trigGeo2.SensitiveVolume(tempHits_B[f_test].channel);
510  const auto hit2Center = hit2Geo.GetCenter();
511  bool moduleMatched;
512  moduleMatched=moduleMatcher(tempHits_B[f_test].module, tempHits_B[f].module);
513 
514  if (moduleMatched) {
515  double hitX = hit1Center.X();
516 
517 
518  for (unsigned int a = 0; a < tempHits_B.size(); a++)
519  {
520  if(tempHits_B[a].module==tempHits_B[f].module && (tempHits_B[a].channel-1)==tempHits_B[f].channel) hitX=hit1Center.X()+1.25;
521  }
522 
523  double hitYPrelim = hit2Center.Y();
524 
525  for (unsigned int a = 0; a < tempHits_B.size(); a++)
526  {
527  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;
528  }
529  double hitY=hitYPrelim;
530 
531 
532  double hitZ = (hit1Center.Z() + hit2Center.Z()) / 2.f;
533 
534  recoHits rHits;
535  rHits.hitPositionX = hitX;
536  rHits.hitPositionY = hitY;
537  rHits.hitPositionZ = hitZ;
538  rHits.geoX=tempHits_B[f].module;
539  rHits.geoY=tempHits_B[f_test].module;
540  rHits.adcX=tempHits_B[f].adc;
541  rHits.adcY=tempHits_B[f_test].adc;
542  rHits.stripX=tempHits_B[f].channel;
543  rHits.stripY=tempHits_B[f_test].channel;
544  rHits.trigNumberX=tempHits_B[f].triggerNumber;
545  rHits.trigNumberY=tempHits_B[f_test].triggerNumber;
546  rHits.timeAvg = (tempHits_B[f_test].triggerTime+tempHits_B[f].triggerTime)/2.0;
547  primaryHits_B.push_back(rHits);
548 
549  }
550  }
551  }
552 
553  std::cout<<primaryHits_F.size()<<','<<primaryHits_B.size()<<std::endl;
554  // Reconstruciton information
555 
556  vector < art::Ptr < recob::Track > > trackList;
557  auto trackListHandle = event.getHandle < vector < recob::Track > >(fTrackModuleLabel);
558  if (trackListHandle) {
559  art::fill_ptr_vector(trackList, trackListHandle);
560  }
561 
562  vector<art::Ptr<recob::PFParticle> > pfplist;
563  auto PFPListHandle = event.getHandle< std::vector<recob::PFParticle> >("pandora");
564  if (PFPListHandle) art::fill_ptr_vector(pfplist, PFPListHandle);
565 
566  if (pfplist.size()<1) return;
567 art::FindManyP<anab::T0> trk_t0_assn_v(PFPListHandle, event ,"pandora");
568  art::FindManyP<recob::PFParticle> pfp_trk_assn(trackListHandle,event,"pandoraTrack");
569 
570  art::FindManyP<recob::Hit> trackHits(trackListHandle, event, "pandoraTrack");
571  int nTracksReco = trackList.size();
572  //cout<<"Number of Potential CRT Reconstructed Through-Going Muon: "<<combTrackHits.size()<<endl;
573  std::vector<art::Ptr<recob::Hit>> hitlist; // to get information about the hits
574  auto hitListHandle = event.getHandle< std::vector<recob::Hit> >("hitpdune");
575  if (hitListHandle)
576  art::fill_ptr_vector(hitlist, hitListHandle);
577  art::FindManyP < recob::Hit > hitsFromTrack(trackListHandle, event, fTrackModuleLabel);
578 
579 
581 
582  art::FindManyP<anab::Calorimetry> fmcal(trackListHandle, event, "pandoracalo");
583  int tempId = 0;
584  allTracksPair.clear();
585  for (int iRecoTrack = 0; iRecoTrack < nTracksReco; ++iRecoTrack) {
586  if (pfplist.size()<1) break;
587  std::vector< art::Ptr<recob::Hit> > allHits = hitsFromTrack.at(iRecoTrack);
588 
589  art::Ptr<recob::Track> ptrack(trackListHandle, iRecoTrack);
590 
591 
592  std::vector<art::Ptr<recob::PFParticle>> pfps=pfp_trk_assn.at(iRecoTrack);
593  if(!pfps.size()) continue;
594  std::vector<art::Ptr<anab::T0>> t0s=trk_t0_assn_v.at(pfps[0].key());
595  /*if(t0s.size()){
596  auto t0=t0s.at(0);
597  int t_zero=t0->Time();
598  cout<<"Pandora T0: "<<t_zero<<endl;
599  }*/
600 
601 
602 
603 
604 
605 
606  double trackStartPositionZ_noSCE = trackList[iRecoTrack]->Vertex().Z();
607  double trackEndPositionZ_noSCE = trackList[iRecoTrack] -> End().Z();
608 
609  double trackStartPositionX_notCorrected = trackList[iRecoTrack]->Vertex().X();
610  double trackStartPositionY_noSCE = trackList[iRecoTrack]->Vertex().Y();
611 
612 
613  double trackEndPositionX_notCorrected = trackList[iRecoTrack] -> End().X();
614  double trackEndPositionY_noSCE = trackList[iRecoTrack] -> End().Y();
615 
616  int firstHit=0;
617  //int lastHit=allHits.size()-2;
618  if (trackStartPositionZ_noSCE>trackEndPositionZ_noSCE){
619  trackEndPositionZ_noSCE = trackList[iRecoTrack]->Vertex().Z();
620  trackStartPositionZ_noSCE = trackList[iRecoTrack] -> End().Z();
621  trackEndPositionX_notCorrected = trackList[iRecoTrack]->Vertex().X();
622  trackEndPositionY_noSCE = trackList[iRecoTrack]->Vertex().Y();
623 
624 
625  trackStartPositionX_notCorrected=trackList[iRecoTrack] -> End().X();
626  trackStartPositionY_noSCE = trackList[iRecoTrack] -> End().Y();
627 
628  }
629 
630 
631 
632 
633 
634 
635 
636 
637  if ((trackEndPositionZ_noSCE > 660 && trackStartPositionZ_noSCE < 50) || (trackStartPositionZ_noSCE > 660 && trackEndPositionZ_noSCE < 50)) {
638 
639 
640  double min_delta = DBL_MAX;
641  double best_XF = DBL_MAX;
642  double best_YF = DBL_MAX;
643  double best_ZF = DBL_MAX;
644  double best_XB = DBL_MAX;
645  double best_YB = DBL_MAX;
646  double best_ZB = DBL_MAX;
647  double best_dotProductCos = DBL_MAX;
648  double best_deltaXF = DBL_MAX;
649  double best_deltaYF = DBL_MAX;
650  double best_deltaXB = DBL_MAX;
651  double best_deltaYB = DBL_MAX;
652  double best_T=DBL_MAX;
653  int bestHitIndex_F=0;
654  int bestHitIndex_B=0;
655 
656  /*
657  int best_trigXF=0;
658  int best_trigYF=0;
659  int best_trigXB=0;
660  int best_trigYB=0;
661  */
662  for (unsigned int f = 0; f < primaryHits_F.size(); f++) {
663  for (unsigned int b = 0; b < primaryHits_B.size(); b++) {
664 
665  double X1 = primaryHits_F[f].hitPositionX;
666  double Y1 = primaryHits_F[f].hitPositionY;
667  double Z1 = primaryHits_F[f].hitPositionZ;
668  double X2 = primaryHits_B[b].hitPositionX;
669  double Y2 = primaryHits_B[b].hitPositionY;
670  double Z2= primaryHits_B[b].hitPositionZ;
671 
672 
673 
674  if (fabs(primaryHits_F[f].timeAvg-primaryHits_B[b].timeAvg)>fFronttoBackTimingCut) continue;
675 
676  //std::cout<<"FOUND A COMBO"<<std::endl;
677  double t0=(primaryHits_F[f].timeAvg+primaryHits_B[b].timeAvg)/2.f;
678 
679  int tempId = 0;
680  double xOffset=0;
681 
682 
683  double trackStartPositionX_noSCE=trackStartPositionX_notCorrected;
684  double trackEndPositionX_noSCE=trackEndPositionX_notCorrected;
685  if (t0s.empty()){
686 
687  int RDOffset=0;
688  if (!fMCCSwitch) RDOffset=111;
689  double ticksOffset=0;
690 
691  if (!fMCCSwitch) ticksOffset = (t0+RDOffset)/25.f+detProp.GetXTicksOffset(allHits[firstHit]->WireID().Plane, allHits[firstHit]->WireID().TPC, allHits[firstHit]->WireID().Cryostat);
692 
693  else if (fMCCSwitch) ticksOffset = (t0/500.f)+detProp.GetXTicksOffset(allHits[firstHit]->WireID().Plane, allHits[firstHit]->WireID().TPC, allHits[firstHit]->WireID().Cryostat);
694 
695  xOffset=detProp.ConvertTicksToX(ticksOffset,allHits[firstHit]->WireID().Plane, allHits[firstHit]->WireID().TPC, allHits[firstHit]->WireID().Cryostat);
696  //double xOffset=.08*ticksOffset
697 
698  trackStartPositionX_noSCE=trackStartPositionX_notCorrected-xOffset;
699  trackEndPositionX_noSCE=trackEndPositionX_notCorrected-xOffset;
700 
701 }
702 
703  double trackStartPositionX=trackStartPositionX_noSCE;
704  double trackStartPositionY=trackStartPositionY_noSCE;
705  double trackStartPositionZ=trackStartPositionZ_noSCE;
706 
707  double trackEndPositionX=trackEndPositionX_noSCE;
708  double trackEndPositionY=trackEndPositionY_noSCE;
709  double trackEndPositionZ=trackEndPositionZ_noSCE;
710  if (fSCECorrection && SCE->EnableCalSpatialSCE()){
711 if(geom->PositionToTPCID(geo::Point_t(trackEndPositionX, trackEndPositionY, trackEndPositionZ)).deepestIndex()<13 && geom->PositionToTPCID(geo::Point_t(trackStartPositionX, trackStartPositionY, trackStartPositionZ)).deepestIndex()<13){
712  auto const & posOffsets_F = SCE->GetCalPosOffsets(geo::Point_t(trackStartPositionX, trackStartPositionY, trackStartPositionZ), geom->PositionToTPCID(geo::Point_t(trackStartPositionX, trackStartPositionY, trackStartPositionZ)).deepestIndex());
713  trackStartPositionX -= posOffsets_F.X();
714  trackStartPositionY += posOffsets_F.Y();
715  trackStartPositionZ += posOffsets_F.Z();
716  auto const & posOffsets_B = SCE->GetCalPosOffsets(geo::Point_t(trackEndPositionX, trackEndPositionY, trackEndPositionZ), geom->PositionToTPCID(geo::Point_t(trackEndPositionX, trackEndPositionY, trackEndPositionZ)).deepestIndex());
717  trackEndPositionX -= posOffsets_B.X();
718  trackEndPositionY += posOffsets_B.Y();
719  trackEndPositionZ += posOffsets_B.Z();
720  }
721  }
722 
723 
724 
725 
726 
727  // Make metrics for a CRT pair to compare later
728  TVector3 trackStart(trackStartPositionX, trackStartPositionY, trackStartPositionZ);
729  TVector3 trackEnd(trackEndPositionX, trackEndPositionY, trackEndPositionZ);
730  TVector3 v1(X1,Y1,Z1);
731  TVector3 v2(X2, Y2, Z2);
732 
733  TVector3 v4(trackStartPositionX,
734  trackStartPositionY,
735  trackStartPositionZ);
736  TVector3 v5(trackEndPositionX,
737  trackEndPositionY,
738  trackEndPositionZ);
739  TVector3 trackVector = (v5-v4).Unit();
740  TVector3 hitVector=(v2-v1).Unit();
741 
742 
743 
744 
745  double predictedHitPositionY1 = (v1.Z()-v5.Z())/(v4.Z()-v5.Z())*(v4.Y()-v5.Y())+v5.Y();
746  double predictedHitPositionY2 = (v2.Z()-v5.Z())/(v4.Z()-v5.Z())*(v4.Y()-v5.Y())+v5.Y();
747 
748  double predictedHitPositionX1 = (v1.Z()-v5.Z())/(v4.Z()-v5.Z())*(v4.X()-v5.X())+v5.X();
749  double predictedHitPositionX2 = (v2.Z()-v5.Z())/(v4.Z()-v5.Z())*(v4.X()-v5.X())+v5.X();
750 
751  double dotProductCos=trackVector*hitVector;
752 
753  double deltaX1 = (predictedHitPositionX1-X1);
754  double deltaX2 = (predictedHitPositionX2-X2);
755 
756  double deltaY1 = (predictedHitPositionY1-Y1);
757 
758  double deltaY2 = (predictedHitPositionY2-Y2);
759 
760 
761 
762  tempId++;
763  //iRecoTrack
764  if (min_delta > std::abs(deltaX1)+std::abs(deltaX2) + std::abs(deltaY1)+std::abs(deltaY2) ){
765 
766  min_delta=std::abs(deltaX1)+std::abs(deltaX2) + std::abs(deltaY1)+std::abs(deltaY2);
767  best_XF = X1;
768  best_YF = Y1;
769  best_ZF = Z1;
770  best_XB = X2;
771  best_YB = Y2;
772  best_ZB = Z2;
773  best_dotProductCos = dotProductCos;
774 
775  best_deltaXF = deltaX1;
776  best_deltaYF = deltaY1;
777  best_deltaXB = deltaX2;
778  best_deltaYB = deltaY2;
779  /*
780  best_trigXF=primaryHits_F[f].trigNumberX;
781  best_trigYF=primaryHits_F[f].trigNumberY;
782  best_trigXB=primaryHits_B[b].trigNumberX;
783  best_trigYB=primaryHits_B[b].trigNumberY;
784  */
785  best_T = t0;
786  bestHitIndex_F=f;
787  bestHitIndex_B=b;
788  if (!fMCCSwitch) best_T=(111.f+t0)*20.f;
789  // Added 111 tick CRT-CTB offset
790  }
791  }
792  }
793  X_F=best_XF; Y_F=best_YF; Z_F=best_ZF; X_B=best_XB; Y_B=best_YB; Z_B=best_ZB; int f=bestHitIndex_F; int b=bestHitIndex_B;
794  double t0=best_T;
795  deltaX_F=best_deltaXF; deltaY_F=best_deltaYF; deltaX_B=best_deltaXB; deltaY_B=best_deltaYB;
796  std::cout<<deltaX_F<<','<<deltaX_B<<','<<deltaY_F<<','<<deltaY_B<<','<<t0<<std::endl;
797  if (deltaX_F==DBL_MAX && deltaY_F==DBL_MAX && deltaX_B==DBL_MAX && deltaY_B==DBL_MAX) continue;
798 
799 
800  double trackStartPositionX_noSCE=trackStartPositionX_notCorrected;
801  double trackEndPositionX_noSCE=trackEndPositionX_notCorrected;
802  double xOffset=0;
803  if (t0s.empty()){
804  double ticksOffset=0;
805 
806 
807  ticksOffset=t0/500.f+detProp.GetXTicksOffset(allHits[firstHit]->WireID().Plane, allHits[firstHit]->WireID().TPC, allHits[firstHit]->WireID().Cryostat);
808 
809  xOffset=detProp.ConvertTicksToX(ticksOffset,allHits[firstHit]->WireID().Plane, allHits[firstHit]->WireID().TPC, allHits[firstHit]->WireID().Cryostat);
810 
811  trackStartPositionX_noSCE=trackStartPositionX_notCorrected-xOffset;
812  trackEndPositionX_noSCE=trackEndPositionX_notCorrected-xOffset;
813  }
814 
815 
816  double trackStartPositionX=trackStartPositionX_noSCE;
817  double trackStartPositionY=trackStartPositionY_noSCE;
818  double trackStartPositionZ=trackStartPositionZ_noSCE;
819 
820  double trackEndPositionX=trackEndPositionX_noSCE;
821  double trackEndPositionY=trackEndPositionY_noSCE;
822  double trackEndPositionZ=trackEndPositionZ_noSCE;
823 
824 
825  TVector3 trackStart(trackStartPositionX, trackStartPositionY, trackStartPositionZ);
826  TVector3 trackEnd(trackEndPositionX, trackEndPositionY, trackEndPositionZ);
827 
828  tracksPair tPair;
829  tPair.tempId = tempId;
830  tPair.recoId = iRecoTrack;
831  tPair.deltaX_F = deltaX_F;
832 
833  tPair.deltaX_B = deltaX_B;
834  tPair.deltaY_F = deltaY_F;
835  tPair.deltaY_B = deltaY_B;
836 
837  tPair.dotProductCos=best_dotProductCos;
838 
839  tPair.moduleX1 = primaryHits_F[f].geoX;
840  tPair.moduleX2 = primaryHits_B[b].geoX;
841  tPair.moduleY1 = primaryHits_F[f].geoY;
842  tPair.moduleY2 = primaryHits_B[b].geoY;
843  tPair.timeDiff=primaryHits_F[f].timeAvg-primaryHits_B[b].timeAvg;
844  tPair.adcX2=primaryHits_B[b].adcX;
845  tPair.adcX1=primaryHits_F[f].adcX;
846  tPair.adcY2=primaryHits_B[b].adcY;
847  tPair.adcY1=primaryHits_F[f].adcY;
848  tPair.stripX1 = primaryHits_F[f].stripX;
849  tPair.stripX2 = primaryHits_B[b].stripX;
850  tPair.stripY1 = primaryHits_F[f].stripY;
851  tPair.stripY2 = primaryHits_B[b].stripY;
852  tPair.X1 = X_F;
853  tPair.Y1 = Y_F;
854  tPair.Z1 = Z_F;
855  tPair.X2 = X_B;
856  tPair.Y2 = Y_B;
857  tPair.Z2 = Z_B;
858 
859  tPair.endTPC=allHits[0]->WireID().TPC;
860 
861  tPair.xOffset=xOffset;
862  tPair.t0=t0;
863  tPair.trackStartPosition=trackStart;
864  tPair.trackEndPosition=trackEnd;
865 
866  if (t0s.empty()) tPair.pandoraT0Check=0;
867  else tPair.pandoraT0Check=1;
868  allTracksPair.push_back(tPair);
869 
870 
871 
872 
873  tempId++;
874  } //iRecoTrack
875  }
876 
877  //Sort pair by ascending order of absolute distance
878  sort(allTracksPair.begin(), allTracksPair.end(), sortPair());
879 
880  // Compare, sort, and eliminate CRT hits for just the best one
881 
882 
883  vector < tracksPair > allUniqueTracksPair;
884  while (allTracksPair.size()) {
885  allUniqueTracksPair.push_back(allTracksPair.front());
886  allTracksPair.erase(remove_if(allTracksPair.begin(), allTracksPair.end(), removePairIndex(allTracksPair.front())),
887  allTracksPair.end());
888  }
889 
890  cout<<"Number of reco and CRT pairs: "<<allUniqueTracksPair.size()<<endl;
891 // For the best one, add the validation metrics to a tree
892  if (allUniqueTracksPair.size() > 0) {
893  for (unsigned int u = 0; u < allUniqueTracksPair.size(); u++) {
894  trackID = allUniqueTracksPair[u].recoId;
895 
896 
897 
898 
899  deltaX_F=allUniqueTracksPair[u].deltaX_F;
900  deltaX_B=allUniqueTracksPair[u].deltaX_B;
901  deltaY_F=allUniqueTracksPair[u].deltaY_F;
902  deltaY_B=allUniqueTracksPair[u].deltaY_B;
903  dotCos=allUniqueTracksPair[u].dotProductCos;
904  trackX1=allUniqueTracksPair[u].trackStartPosition.X();
905  trackY1=allUniqueTracksPair[u].trackStartPosition.Y();
906  trackZ1=allUniqueTracksPair[u].trackStartPosition.Z();
907 
908  trackX2=allUniqueTracksPair[u].trackEndPosition.X();
909  trackY2=allUniqueTracksPair[u].trackEndPosition.Y();
910  trackZ2=allUniqueTracksPair[u].trackEndPosition.Z();
911 
912  moduleX_F=allUniqueTracksPair[u].moduleX1;
913  moduleX_B=allUniqueTracksPair[u].moduleX2;
914  moduleY_F=allUniqueTracksPair[u].moduleY1;
915  moduleY_B=allUniqueTracksPair[u].moduleY2;
916 
917  adcX_F=allUniqueTracksPair[u].adcX1;
918  adcY_F=allUniqueTracksPair[u].adcY1;
919  adcX_B=allUniqueTracksPair[u].adcX2;
920  adcY_B=allUniqueTracksPair[u].adcY2;
921 
922  stripX_F=allUniqueTracksPair[u].stripX1;
923  stripX_B=allUniqueTracksPair[u].stripX2;
924  stripY_F=allUniqueTracksPair[u].stripY1;
925  stripY_B=allUniqueTracksPair[u].stripY2;
926 
927  X_F=allUniqueTracksPair[u].X1;
928  Y_F=allUniqueTracksPair[u].Y1;
929  Z_F=allUniqueTracksPair[u].Z1;
930 
931  X_B=allUniqueTracksPair[u].X2;
932  Y_B=allUniqueTracksPair[u].Y2;
933  Z_B=allUniqueTracksPair[u].Z2;
934 
935  endTPC=allUniqueTracksPair[u].endTPC;
936 
937  recoPandoraT0Check=allUniqueTracksPair[u].pandoraT0Check;
938  mccT0=allUniqueTracksPair[u].mccT0;
939  measuredT0=allUniqueTracksPair[u].t0;
940  measuredXOffset=allUniqueTracksPair[u].xOffset;
941  CRT_TOF=allUniqueTracksPair[u].timeDiff;
942  if (fabs(allUniqueTracksPair[u].dotProductCos)>0.99 && fabs(deltaX_F)+fabs(deltaX_B)<40 && fabs(deltaY_F)+fabs(deltaY_B)<40) {
943  //cout<<allUniqueTracksPair[u].timeDiff<<endl;
944  //cout<<fabs(allUniqueTracksPair[u].dotProductCos)<<endl;
945 
946  cout<< "Delta X and Y: "<<deltaX_F<<','<<deltaY_F<<','<<deltaX_B<<','<<deltaY_B<<endl;
947  cout<< "Predicted X and Y: "<< deltaX_F+X_F<<','<<deltaY_F+Y_F<<','<<deltaX_B+X_B<<','<<deltaY_B+Y_B<<endl;
948  cout<< "Detected X and Y: "<< X_F<<','<<Y_F<<','<<X_B<<','<<Y_B<<endl;
949  cout<<moduleX_F<<','<<stripX_F<<endl;
950  cout<<moduleY_F<<','<<stripY_F<<endl;
951  cout<<moduleX_B<<','<<stripX_B<<endl;
952  cout<<moduleY_B<<','<<stripY_B<<endl;
953  cout<<"ADC Values: "<<adcX_F<<','<<adcY_F<<','<<adcX_B<<','<<adcY_B<<endl;
954 
955 
956  int iRecoTrack=trackID;
960  const size_t lastPoint=trackList[iRecoTrack]->NumberTrajectoryPoints();
961  for (size_t trackpoint = 0; trackpoint < lastPoint; ++trackpoint) {
962  double trackPosX=trackList[iRecoTrack] -> LocationAtPoint(trackpoint).X()-measuredXOffset;
963  double trackPosY=trackList[iRecoTrack] -> LocationAtPoint(trackpoint).Y();
964  double trackPosZ=trackList[iRecoTrack] -> LocationAtPoint(trackpoint).Z();
965 
966  if (trackPosY==-999) continue;
967  if (trackPosX==-999) continue;
968  TVector3 trackPos(trackPosX, trackPosY, trackPosZ);
969  double distanceYZ = signedPointToLineDistance( Y_F,Z_F, Y_B,Z_B, trackPos.Y(), trackPos.Z() ); //only the Y and Z of trackpos will be used
970  double distanceXZ = signedPointToLineDistance( X_F,Z_F, X_B,Z_B, trackPos.X(), trackPos.Z() );
971 
972 
973  double distance=signed3dDistance(X_F,Y_F,Z_F,X_B,Y_B,Z_B, trackPos);
974 
975 
976  trkPosX=trackPosX;
977  trkPosY=trackPosY;
978  trkPosZ=trackPosZ;
979 
980  displXZ=distanceXZ;
981 
982  displYZ =distanceYZ;
983 
984  averageSignedDistance += distance/(lastPoint+1);
985  averageSignedDistanceYZ += distanceYZ/(lastPoint+1);
986  averageSignedDistanceXZ += distanceXZ/(lastPoint+1);
987 //averageSignedDistanceXY += distanceXY/(lastPoint+1);
988  fTrackInfo->Fill();
989 
990  }
991 
992 
993  std::vector<art::Ptr<anab::Calorimetry>> calos=fmcal.at(iRecoTrack);
994  //std::vector< art::Ptr<recob::Hit> > allHits = hitsFromTrack.at(iRecoTrack);
995  for(size_t ical = 0; ical<calos.size(); ++ical){
996  if (calos[ical]->PlaneID().Plane!=2) continue;
997 
998  const size_t NHits=calos[ical]->dEdx().size();
999  for(size_t iHit = 0; iHit < NHits; ++iHit){
1000 
1001  const auto& TrkPos = (calos[ical] -> XYZ()[iHit]);
1002 
1003  trkdqdx=(calos[ical]->dQdx()[iHit]);
1004  size_t tpIndex=(calos[ical]->TpIndices()[iHit]);
1005  if (hitlist[tpIndex]->Multiplicity()>1) continue;
1006  trkhitt0=hitlist[tpIndex]->PeakTime()*500.f;
1007  trkhitIntegral=hitlist[tpIndex]->Integral();
1008  crtt0=allUniqueTracksPair[u].t0;
1009  trkhitx=TrkPos.X()-measuredXOffset;
1010  trkhity=TrkPos.Y();
1011  trkhitz=TrkPos.Z();
1012  WireID=hitlist[tpIndex]->WireID().Wire;
1013  TPCID=hitlist[tpIndex]->WireID().TPC;
1014  sumADC=hitlist[tpIndex]->SummedADC();
1015  sigmaHit=hitlist[tpIndex]->SigmaIntegral();
1016  rangeTime=hitlist[tpIndex]->EndTick()-hitlist[tpIndex]->StartTick();
1017  fCRTdQTree->Fill();
1018  }
1019 
1020  }
1021 
1022 
1023  fCRTTree->Fill();
1024 
1025 
1026  }
1027  }
1028  }
1029 
1030  }
code to link reconstructed objects back to the MC truth information
double signedPointToLineDistance(double firstPoint1, double firstPoint2, double secondPoint1, double secondPoint2, double trackPoint1, double trackPoint2)
bool moduleMatcher(int module1, int module2)
geo::TPCID PositionToTPCID(geo::Point_t const &point) const
Returns the ID of the TPC at specified location.
uint8_t channel
Definition: CRTFragment.hh:201
T abs(T value)
def key(type, name=None)
Definition: graph.py:13
const double a
double distance(double x1, double y1, double z1, double x2, double y2, double z2)
Detector simulation of raw signals on wires.
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:184
std::vector< tempHits > tempHits_F
std::vector< recoHits > primaryHits_B
void End(void)
Definition: gXSecComp.cxx:210
def center(depos, point)
Definition: depos.py:117
double signed3dDistance(double firstPoint1, double firstPoint2, double firstPoint3, double secondPoint1, double secondPoint2, double secondPoint3, TVector3 trackPos)
std::vector< tracksPair > allTracksPair
std::vector< tempHits > tempHits_B
static bool * b
Definition: config.cpp:1043
constexpr auto const & deepestIndex() const
Returns the value of the deepest ID available (TPC&#39;s).
Definition: geo_types.h:428
detail::Node< FrameID, bool > PlaneID
Definition: CRTID.h:125
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
recob::tracking::Plane Plane
Definition: TrackState.h:17
std::vector< recoHits > primaryHits_F
QTextStream & endl(QTextStream &s)
Event finding and building.
void CRT::TwoCRTMatching::beginJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 1034 of file TwoCRTMatching_module.cc.

1034  {
1035  art::ServiceHandle<art::TFileService> fileServiceHandle;
1036  fCRTTree = fileServiceHandle->make<TTree>("Displacement", "track by track info");
1037  fMCCMuon= fileServiceHandle->make<TTree>("MCCTruths", "event by event info");
1038 
1039  fTrackInfo= fileServiceHandle->make<TTree>("TrackInfo", "track by track info");
1040  fCRTdQTree=fileServiceHandle->make<TTree>("CRTdQ", "track by track info");
1041  fCRTTree->Branch("nEvents", &nEvents, "fnEvents/I");
1042  fCRTTree->Branch("hRun", &run, "run/I");
1043  fCRTTree->Branch("hSubRun", &subRun, "subRun/I");
1044  fCRTTree->Branch("htrackID", &trackID, "trackID/I");
1045  fCRTTree->Branch("hmccTrackId", &mccTrackId, "mccTrackId/I");
1046 
1047  fMCCMuon->Branch("nEvents", &nEvents, "nEvents/I");
1048  fMCCMuon->Branch("mccCRTStartX",&mccCRTStartX,"mccCRTStartX/D");
1049  fMCCMuon->Branch("mccCRTStartY",&mccCRTStartY,"mccCRTStartY/D");
1050  fMCCMuon->Branch("mccCRTStartZ",&mccCRTStartZ,"mccCRTStartZ/D");
1051  fMCCMuon->Branch("mccCRTEndX",&mccCRTEndX,"mccCRTEndX/D");
1052  fMCCMuon->Branch("mccCRTEndY",&mccCRTEndY,"mccCRTEndY/D");
1053  fMCCMuon->Branch("mccCRTEndZ",&mccCRTEndZ,"mccCRTEndZ/D");
1054  fMCCMuon->Branch("truthEnergy", &truthEnergy, "truthEnergy/D");
1055  fMCCMuon->Branch("hRun", &run, "run/I");
1056  fMCCMuon->Branch("hSubRun", &subRun, "subRun/I");
1057  fMCCMuon->Branch("hmccTrackId", &mccTrackId, "mccTrackId/I");
1058  fMCCMuon->Branch("hcandidateCRT", &candidateCRT, "candidateCRT/I");
1059  fMCCMuon->Branch("hmccT0", &mccT0, "mccT0/D");
1060 
1061 
1062  fCRTdQTree->Branch("trkhitx",&trkhitx,"trkhitx/D");
1063  fCRTdQTree->Branch("trkhity",&trkhity,"trkhity/D");
1064  fCRTdQTree->Branch("trkhitz",&trkhitz,"trkhitz/D");
1065  fCRTdQTree->Branch("trkdqdx",&trkdqdx,"trkdqdx/D");
1066  fCRTdQTree->Branch("trkhitt0",&trkhitt0,"trkhitt0/D");
1067  fCRTdQTree->Branch("trkhitIntegral",&trkhitIntegral,"trkhitIntegral/I");
1068  fCRTdQTree->Branch("crtt0",&crtt0,"crtt0/D");
1069  fCRTdQTree->Branch("hX_F", &X_F, "X_F/D");
1070  fCRTdQTree->Branch("hX_B", &X_B, "X_B/D");
1071  fCRTdQTree->Branch("hY_F", &Y_F, "Y_F/D");
1072  fCRTdQTree->Branch("hY_B", &Y_B, "Y_B/D");
1073  fCRTdQTree->Branch("hZ_F", &Z_F, "Z_F/D");
1074  fCRTdQTree->Branch("hZ_B", &Z_B, "Z_B/D");
1075 
1076  fCRTdQTree->Branch("hTPCID", &TPCID, "TPCID/I");
1077  fCRTdQTree->Branch("hWireID", &WireID, "WireID/I");
1078  fCRTdQTree->Branch("hsumADC",&sumADC,"sumADC/I");
1079  fCRTdQTree->Branch("hrangeTime",&rangeTime,"rangeTime/I");
1080  fCRTdQTree->Branch("hsigmaHit",&sigmaHit,"sigmaHit/D");
1081 
1082  fCRTTree->Branch("nHitsPerEvent", &nHitsPerEvent, "fnHitsPerEvent/I");
1083  fCRTTree->Branch("haverageSignedDistanceYZ", &averageSignedDistanceYZ, "averageSignedDistanceYZ/D");
1084  fCRTTree->Branch("haverageSignedDistanceXZ", &averageSignedDistanceXZ, "averageSignedDistanceXZ/D");
1085  fCRTTree->Branch("haverageSignedDistance", &averageSignedDistance, "averageSignedDistance/D");
1086  fCRTTree->Branch("hdisplAngleXY", &displAngleXY, "displAngleXY/D");
1087  fCRTTree->Branch("hdisplAngleYZ", &displAngleYZ, "displAngleYZ/D");
1088  fCRTTree->Branch("hdisplAngleXZ", &displAngleXZ, "displAngleXZ/D");
1089 
1090  fCRTTree->Branch("hdeltaX_F", &deltaX_F, "deltaX_F/D");
1091  fCRTTree->Branch("hdeltaX_B", &deltaX_B, "deltaX_B/D");
1092  fCRTTree->Branch("hdeltaY_F", &deltaY_F, "deltaY_F/D");
1093  fCRTTree->Branch("hdeltaY_B", &deltaY_B, "deltaY_B/D");
1094  fCRTTree->Branch("hdotProductCos", &dotCos, "dotCos/D");
1095 
1096 
1097  fCRTTree->Branch("hX_F", &X_F, "X_F/D");
1098  fCRTTree->Branch("hX_B", &X_B, "X_B/D");
1099  fCRTTree->Branch("hY_F", &Y_F, "Y_F/D");
1100  fCRTTree->Branch("hY_B", &Y_B, "Y_B/D");
1101  fCRTTree->Branch("hZ_F", &Z_F, "Z_F/D");
1102  fCRTTree->Branch("hZ_B", &Z_B, "Z_B/D");
1103 
1104  fCRTTree->Branch("hSCECorrectX_F", &SCECorrectX_F, "SCECorrectX_F/D");
1105  fCRTTree->Branch("hSCECorrectY_F", &SCECorrectY_F, "SCECorrectY_F/D");
1106  fCRTTree->Branch("hSCECorrectZ_F", &SCECorrectZ_F, "SCECorrectZ_F/D");
1107  fCRTTree->Branch("hSCECorrectX_B", &SCECorrectX_B, "SCECorrectX_B/D");
1108  fCRTTree->Branch("hSCECorrectY_B", &SCECorrectY_B, "SCECorrectY_B/D");
1109  fCRTTree->Branch("hSCECorrectZ_B", &SCECorrectZ_B, "SCECorrectZ_B/D");
1110  fCRTTree->Branch("hendTPC", &endTPC, "endTPC/I");
1111 
1112  fCRTTree->Branch("htrackStartX", &trackX1, "trackX1/D");
1113  fCRTTree->Branch("htrackStartY", &trackY1, "trackY1/D");
1114  fCRTTree->Branch("htrackStartZ", &trackZ1, "trackZ1/D");
1115 
1116  fCRTTree->Branch("htrackEndX", &trackX2, "trackX2/D");
1117  fCRTTree->Branch("htrackEndY", &trackY2, "trackY2/D");
1118  fCRTTree->Branch("htrackEndZ", &trackZ2, "trackZ2/D");
1119  fCRTTree->Branch("CRT_TOF", &CRT_TOF, "CRT_TOF/D");
1120 
1121  fCRTTree->Branch("hmoduleX_F", &moduleX_F, "moduleX_F/I");
1122  fCRTTree->Branch("hmoduleX_B", &moduleX_B, "moduleX_B/I");
1123  fCRTTree->Branch("hmoduleY_F", &moduleY_F, "moduleY_F/I");
1124  fCRTTree->Branch("hmoduleY_B", &moduleY_B, "moduleY_B/I");
1125 
1126  fCRTTree->Branch("hstripX_F", &stripX_F, "stripX_F/I");
1127  fCRTTree->Branch("hstripX_B", &stripX_B, "stripX_B/I");
1128  fCRTTree->Branch("hstripY_F", &stripY_F, "stripY_F/I");
1129  fCRTTree->Branch("hstripY_B", &stripY_B, "stripY_B/I");
1130 
1131  fCRTTree->Branch("hadcX_F", &adcX_F, "adcX_F/I");
1132  fCRTTree->Branch("hadcY_F", &adcY_F, "adcY_F/I");
1133  fCRTTree->Branch("hadcX_B", &adcX_B, "adcX_B/I");
1134  fCRTTree->Branch("hadcY_B", &adcY_B, "adcY_B/I");
1135 
1136  fCRTTree->Branch("hmeasuredT0", &measuredT0, "measuredT0/D");
1137  fCRTTree->Branch("hmeasuredXOffset", &measuredXOffset, "measuredXOffset/D");
1138  fCRTTree->Branch("hmccT0", &mccT0, "mccT0/D");
1139  fCRTTree->Branch("hrecoPandoraT0Check", &recoPandoraT0Check, "recoPandoraT0Check/B");
1140 
1141 
1142 
1143  fTrackInfo->Branch("nEvents", &nEvents, "nEvents/I");
1144  fTrackInfo->Branch("hdotCos", &dotCos, "dotCos/D");
1145 
1146 
1147  fTrackInfo->Branch("htrkPosX", &trkPosX, "trkPosX/D");
1148  fTrackInfo->Branch("htrkPosY", &trkPosY, "trkPosY/D");
1149  fTrackInfo->Branch("htrkPosZ", &trkPosZ, "trkPosZ/D");
1150  fTrackInfo->Branch("hdisplXZ", &displXZ, "displXZ/D");
1151  fTrackInfo->Branch("hdisplYZ", &displYZ, "displYZ/D");
1152  fTrackInfo->Branch("hmeasuredXOffset", &measuredXOffset, "measuredXOffset/D");
1153 
1154  fTrackInfo->Branch("hX_F", &X_F, "X_F/D");
1155  fTrackInfo->Branch("hX_B", &X_B, "X_B/D");
1156  fTrackInfo->Branch("hY_F", &Y_F, "Y_F/D");
1157  fTrackInfo->Branch("hY_B", &Y_B, "Y_B/D");
1158  fTrackInfo->Branch("hZ_F", &Z_F, "Z_F/D");
1159  fTrackInfo->Branch("hZ_B", &Z_B, "Z_B/D");
1160 
1161 
1162 
1163 
1164 }
void CRT::TwoCRTMatching::createPNG ( TH1D *  histo)

Definition at line 324 of file TwoCRTMatching_module.cc.

324  {
325  //Save important histograms as PNG
326  TCanvas * c = new TCanvas;
327  TH1D * h = histo;
328  h -> Draw();
329  TImage * img = TImage::Create();
330  img -> FromPad(c);
331  img -> WriteImage((std::string(histo -> GetName()) + ".png").c_str());
332  delete h;
333  delete c;
334  delete img;
335 
336 }
std::string string
Definition: nybbler.cc:12
void Draw(const char *plot, const char *title)
Definition: gXSecComp.cxx:580
void CRT::TwoCRTMatching::endJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 1166 of file TwoCRTMatching_module.cc.

1167 {
1168 
1169 
1170 
1171 }
bool CRT::TwoCRTMatching::moduleMatcher ( int  module1,
int  module2 
)

Definition at line 316 of file TwoCRTMatching_module.cc.

316  {
317  // Function checking if two hits could reasonably be matched into a 2D hit
318  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;
319  else return 0;
320 
321 }
TwoCRTMatching& CRT::TwoCRTMatching::operator= ( TwoCRTMatching const &  )
delete
TwoCRTMatching& CRT::TwoCRTMatching::operator= ( TwoCRTMatching &&  )
delete
double CRT::TwoCRTMatching::setAngle ( double  angle)

Definition at line 337 of file TwoCRTMatching_module.cc.

337  {
338  if (angle < 0) {
339  angle += 3.14159265359;
340  }
341  angle *= 180.00 / 3.14159265359;
342  return angle;
343 }
double CRT::TwoCRTMatching::signed3dDistance ( double  firstPoint1,
double  firstPoint2,
double  firstPoint3,
double  secondPoint1,
double  secondPoint2,
double  secondPoint3,
TVector3  trackPos 
)

Definition at line 293 of file TwoCRTMatching_module.cc.

293  {
294 
295 double denominator = sqrt( (secondPoint2-firstPoint2)*(secondPoint2-firstPoint2) + (secondPoint1-firstPoint1)*(secondPoint1-firstPoint1)+ (secondPoint3-firstPoint3)*(secondPoint3-firstPoint3));
296 
297 double X1=point.X()-firstPoint1;
298 double Y1=point.Y()-firstPoint2;
299 double Z1=point.Z()-firstPoint3;
300 
301 double X2=point.X()-secondPoint1;
302 double Y2=point.Y()-secondPoint2;
303 double Z2=point.Z()-secondPoint3;
304 
305 double numerator=(X1*Y2-Y1*X2)-(X1*Z2-X2*Z1)+(Y1*Z2-Z1*Y2);
306 
307 return numerator/denominator;
308 
309 
310 
311 
312 }
double CRT::TwoCRTMatching::signedPointToLineDistance ( double  firstPoint1,
double  firstPoint2,
double  secondPoint1,
double  secondPoint2,
double  trackPoint1,
double  trackPoint2 
)

Definition at line 286 of file TwoCRTMatching_module.cc.

286  {
287  double numerator = (secondPoint2-firstPoint2)*trackPoint1 - (secondPoint1-firstPoint1) * trackPoint2 + secondPoint1*firstPoint2 - firstPoint1*secondPoint2; //removed the absolute value here, so will give signed distance //the sign indicates a right-hand ruled sign: cross product of line-to-point vector and the direction of the vector (in that order) gives the sign of the result
288  double denominator = sqrt( (secondPoint2-firstPoint2)*(secondPoint2-firstPoint2) + (secondPoint1-firstPoint1)*(secondPoint1-firstPoint1) );
289 
290  return numerator/denominator;
291 
292 }

Member Data Documentation

int CRT::TwoCRTMatching::adcX_B
private

Definition at line 150 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::adcX_F
private

Definition at line 150 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::adcY_B
private

Definition at line 150 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::adcY_F
private

Definition at line 150 of file TwoCRTMatching_module.cc.

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

Definition at line 270 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::averageSignedDistance
private

Definition at line 140 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::averageSignedDistanceXZ
private

Definition at line 139 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::averageSignedDistanceYZ
private

Definition at line 138 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::candidateCRT
private

Definition at line 168 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::CRT_TOF
private

Definition at line 144 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::crtt0
private

Definition at line 169 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::deltaX_B
private

Definition at line 146 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::deltaX_F
private

Definition at line 145 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::deltaY_B
private

Definition at line 148 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::deltaY_F
private

Definition at line 147 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::displAngleXY
private

Definition at line 143 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::displAngleXZ
private

Definition at line 141 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::displAngleYZ
private

Definition at line 142 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::displXZ
private

Definition at line 165 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::displYZ
private

Definition at line 165 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::dotCos
private

Definition at line 149 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::endTPC
private

Definition at line 164 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::fADCThreshold
private

Definition at line 132 of file TwoCRTMatching_module.cc.

TTree* CRT::TwoCRTMatching::fCRTdQTree
private

Definition at line 124 of file TwoCRTMatching_module.cc.

art::InputTag CRT::TwoCRTMatching::fCRTLabel
private

Definition at line 120 of file TwoCRTMatching_module.cc.

TTree* CRT::TwoCRTMatching::fCRTTree
private

Definition at line 122 of file TwoCRTMatching_module.cc.

art::InputTag CRT::TwoCRTMatching::fCTBLabel
private

Definition at line 121 of file TwoCRTMatching_module.cc.

bool CRT::TwoCRTMatching::fCTBTriggerOnly
private

Definition at line 129 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::fFronttoBackTimingCut
private

Definition at line 134 of file TwoCRTMatching_module.cc.

TTree* CRT::TwoCRTMatching::fMCCMuon
private

Definition at line 125 of file TwoCRTMatching_module.cc.

bool CRT::TwoCRTMatching::fMCCSwitch
private

Definition at line 128 of file TwoCRTMatching_module.cc.

bool CRT::TwoCRTMatching::fModuleSwitch
private

Definition at line 131 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::fModuletoModuleTimingCut
private

Definition at line 133 of file TwoCRTMatching_module.cc.

bool CRT::TwoCRTMatching::fSCECorrection
private

Definition at line 130 of file TwoCRTMatching_module.cc.

TTree* CRT::TwoCRTMatching::fTrackInfo
private

Definition at line 126 of file TwoCRTMatching_module.cc.

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

Definition at line 90 of file TwoCRTMatching_module.cc.

ofstream CRT::TwoCRTMatching::logFile
private

Definition at line 117 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::mccCRTEndX
private

Definition at line 167 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::mccCRTEndY
private

Definition at line 167 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::mccCRTEndZ
private

Definition at line 167 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::mccCRTStartX
private

Definition at line 166 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::mccCRTStartY
private

Definition at line 166 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::mccCRTStartZ
private

Definition at line 166 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::mccT0
private

Definition at line 161 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::mccTrackId
private

Definition at line 160 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::measuredT0
private

Definition at line 155 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::measuredXOffset
private

Definition at line 156 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::moduleX_B
private

Definition at line 153 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::moduleX_F
private

Definition at line 153 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::moduleY_B
private

Definition at line 153 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::moduleY_F
private

Definition at line 153 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::nEvents = 0

Definition at line 113 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::nHaloMuons = 0

Definition at line 114 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::nHitsPerEvent =0

Definition at line 115 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::nMCCMuons =0

Definition at line 116 of file TwoCRTMatching_module.cc.

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

Definition at line 266 of file TwoCRTMatching_module.cc.

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

Definition at line 265 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::rangeTime
private

Definition at line 175 of file TwoCRTMatching_module.cc.

bool CRT::TwoCRTMatching::recoPandoraT0Check
private

Definition at line 157 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::run
private

Definition at line 127 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::SCECorrectX_B
private

Definition at line 163 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::SCECorrectX_F
private

Definition at line 162 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::SCECorrectY_B
private

Definition at line 163 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::SCECorrectY_F
private

Definition at line 162 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::SCECorrectZ_B
private

Definition at line 163 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::SCECorrectZ_F
private

Definition at line 162 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::sigmaHit
private

Definition at line 173 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::stripX_B
private

Definition at line 154 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::stripX_F
private

Definition at line 154 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::stripY_B
private

Definition at line 154 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::stripY_F
private

Definition at line 154 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::subRun
private

Definition at line 127 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::sumADC
private

Definition at line 175 of file TwoCRTMatching_module.cc.

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

Definition at line 269 of file TwoCRTMatching_module.cc.

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

Definition at line 268 of file TwoCRTMatching_module.cc.

long long CRT::TwoCRTMatching::timeStamp
private

Definition at line 176 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::TPCID
private

Definition at line 174 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::trackID
private

Definition at line 158 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trackX1
private

Definition at line 152 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trackX2
private

Definition at line 152 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trackY1
private

Definition at line 152 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trackY2
private

Definition at line 152 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trackZ1
private

Definition at line 152 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trackZ2
private

Definition at line 152 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trkdqdx
private

Definition at line 169 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::trkhitIntegral
private

Definition at line 170 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trkhitt0
private

Definition at line 169 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trkhitx
private

Definition at line 169 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trkhity
private

Definition at line 169 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trkhitz
private

Definition at line 169 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trkPosX
private

Definition at line 165 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trkPosY
private

Definition at line 165 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::trkPosZ
private

Definition at line 165 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::truthEnergy
private

Definition at line 159 of file TwoCRTMatching_module.cc.

int CRT::TwoCRTMatching::WireID
private

Definition at line 174 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::X_B
private

Definition at line 151 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::X_F
private

Definition at line 151 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::Y_B
private

Definition at line 151 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::Y_F
private

Definition at line 151 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::Z_B
private

Definition at line 151 of file TwoCRTMatching_module.cc.

double CRT::TwoCRTMatching::Z_F
private

Definition at line 151 of file TwoCRTMatching_module.cc.


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