Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DDetectorMatches_factory.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: DDetectorMatches_factory.cc
4 // Created: Tue Aug 9 14:29:24 EST 2011
5 // Creator: pmatt (on Linux ifarml6 2.6.18-128.el5 x86_64)
6 //
7 
9 
10 //------------------
11 // init
12 //------------------
14 {
15  return NOERROR;
16 }
17 
18 //------------------
19 // brun
20 //------------------
21 jerror_t DDetectorMatches_factory::brun(jana::JEventLoop *locEventLoop, int32_t runnumber)
22 {
23  //LEAVE THIS EMPTY!!! OR ELSE WON'T BE INITIALIZED PROPERLY WHEN "COMBO" FACTORY CALLS Create_DDetectorMatches ON REST DATA!!
24  return NOERROR;
25 }
26 
27 //------------------
28 // evnt
29 //------------------
30 jerror_t DDetectorMatches_factory::evnt(jana::JEventLoop* locEventLoop, uint64_t eventnumber)
31 {
32  vector<const DTrackTimeBased*> locTrackTimeBasedVector;
33  locEventLoop->Get(locTrackTimeBasedVector);
34 
35  DDetectorMatches* locDetectorMatches = Create_DDetectorMatches(locEventLoop, locTrackTimeBasedVector);
36  _data.push_back(locDetectorMatches);
37 
38  return NOERROR;
39 }
40 
41 DDetectorMatches* DDetectorMatches_factory::Create_DDetectorMatches(jana::JEventLoop* locEventLoop, vector<const DTrackTimeBased*>& locTrackTimeBasedVector)
42 {
43  const DParticleID* locParticleID = NULL;
44  locEventLoop->GetSingle(locParticleID);
45 
46  vector<const DSCHit*> locSCHits;
47  locEventLoop->Get(locSCHits);
48 
49  vector<const DTOFPoint*> locTOFPoints;
50  locEventLoop->Get(locTOFPoints);
51 
52  vector<const DFCALShower*> locFCALShowers;
53  locEventLoop->Get(locFCALShowers);
54 
55  vector<const DBCALShower*> locBCALShowers;
56  locEventLoop->Get(locBCALShowers);
57 
58  vector<const DDIRCPmtHit*> locDIRCHits;
59  locEventLoop->Get(locDIRCHits);
60 
61  // cheat and get truth info of track at bar
62  vector<const DDIRCTruthBarHit*> locDIRCBarHits;
63  locEventLoop->Get(locDIRCBarHits);
64 
65  DDetectorMatches* locDetectorMatches = new DDetectorMatches();
66 
67  //Match tracks to showers/hits
68  for(size_t loc_i = 0; loc_i < locTrackTimeBasedVector.size(); ++loc_i)
69  {
70  MatchToBCAL(locParticleID, locTrackTimeBasedVector[loc_i], locBCALShowers, locDetectorMatches);
71  MatchToTOF(locParticleID, locTrackTimeBasedVector[loc_i], locTOFPoints, locDetectorMatches);
72  MatchToFCAL(locParticleID, locTrackTimeBasedVector[loc_i], locFCALShowers, locDetectorMatches);
73  MatchToSC(locParticleID, locTrackTimeBasedVector[loc_i], locSCHits, locDetectorMatches);
74  MatchToDIRC(locParticleID, locTrackTimeBasedVector[loc_i], locDIRCHits, locDetectorMatches, locDIRCBarHits);
75  }
76 
77  //Find nearest tracks to showers
78  for(size_t loc_i = 0; loc_i < locBCALShowers.size(); ++loc_i)
79  MatchToTrack(locParticleID, locBCALShowers[loc_i], locTrackTimeBasedVector, locDetectorMatches);
80  for(size_t loc_i = 0; loc_i < locFCALShowers.size(); ++loc_i)
81  MatchToTrack(locParticleID, locFCALShowers[loc_i], locTrackTimeBasedVector, locDetectorMatches);
82 
83  //Set flight-time/p correlations
84  for(size_t loc_i = 0; loc_i < locTrackTimeBasedVector.size(); ++loc_i)
85  {
86  double locFlightTimePCorrelation = locParticleID->Calc_BCALFlightTimePCorrelation(locTrackTimeBasedVector[loc_i], locDetectorMatches);
87  if(isfinite(locFlightTimePCorrelation))
88  locDetectorMatches->Set_FlightTimePCorrelation(locTrackTimeBasedVector[loc_i], SYS_BCAL, locFlightTimePCorrelation);
89 
90  locFlightTimePCorrelation = locParticleID->Calc_TOFFlightTimePCorrelation(locTrackTimeBasedVector[loc_i], locDetectorMatches);
91  if(isfinite(locFlightTimePCorrelation))
92  locDetectorMatches->Set_FlightTimePCorrelation(locTrackTimeBasedVector[loc_i], SYS_TOF, locFlightTimePCorrelation);
93 
94  locFlightTimePCorrelation = locParticleID->Calc_FCALFlightTimePCorrelation(locTrackTimeBasedVector[loc_i], locDetectorMatches);
95  if(isfinite(locFlightTimePCorrelation))
96  locDetectorMatches->Set_FlightTimePCorrelation(locTrackTimeBasedVector[loc_i], SYS_FCAL, locFlightTimePCorrelation);
97 
98  locFlightTimePCorrelation = locParticleID->Calc_SCFlightTimePCorrelation(locTrackTimeBasedVector[loc_i], locDetectorMatches);
99  if(isfinite(locFlightTimePCorrelation))
100  locDetectorMatches->Set_FlightTimePCorrelation(locTrackTimeBasedVector[loc_i], SYS_START, locFlightTimePCorrelation);
101  }
102 
103  return locDetectorMatches;
104 }
105 
106 void DDetectorMatches_factory::MatchToBCAL(const DParticleID* locParticleID, const DTrackTimeBased* locTrackTimeBased, const vector<const DBCALShower*>& locBCALShowers, DDetectorMatches* locDetectorMatches) const
107 {
108  vector<DTrackFitter::Extrapolation_t>extrapolations=locTrackTimeBased->extrapolations.at(SYS_BCAL);
109  if (extrapolations.size()==0) return;
110 
111  double locInputStartTime = locTrackTimeBased->t0();
112  for(size_t loc_i = 0; loc_i < locBCALShowers.size(); ++loc_i)
113  {
114  shared_ptr<DBCALShowerMatchParams> locShowerMatchParams;
115  if(locParticleID->Cut_MatchDistance(extrapolations, locBCALShowers[loc_i], locInputStartTime, locShowerMatchParams))
116  locDetectorMatches->Add_Match(locTrackTimeBased, locBCALShowers[loc_i], locShowerMatchParams);
117  }
118 }
119 
120 void DDetectorMatches_factory::MatchToTOF(const DParticleID* locParticleID, const DTrackTimeBased* locTrackTimeBased, const vector<const DTOFPoint*>& locTOFPoints, DDetectorMatches* locDetectorMatches) const
121 {
122  vector<DTrackFitter::Extrapolation_t> extrapolations=locTrackTimeBased->extrapolations.at(SYS_TOF);
123  if (extrapolations.size()==0) return;
124 
125  double locInputStartTime = locTrackTimeBased->t0();
126  for(size_t loc_i = 0; loc_i < locTOFPoints.size(); ++loc_i)
127  {
128  shared_ptr<DTOFHitMatchParams> locTOFHitMatchParams;
129  if(locParticleID->Cut_MatchDistance(extrapolations, locTOFPoints[loc_i], locInputStartTime, locTOFHitMatchParams))
130  locDetectorMatches->Add_Match(locTrackTimeBased, locTOFPoints[loc_i], locTOFHitMatchParams);
131  }
132 }
133 
134 void DDetectorMatches_factory::MatchToFCAL(const DParticleID* locParticleID, const DTrackTimeBased* locTrackTimeBased, const vector<const DFCALShower*>& locFCALShowers, DDetectorMatches* locDetectorMatches) const
135 {
136  vector<DTrackFitter::Extrapolation_t> extrapolations=locTrackTimeBased->extrapolations.at(SYS_FCAL);
137  if (extrapolations.size()==0) return;
138 
139  double locInputStartTime = locTrackTimeBased->t0();
140  for(size_t loc_i = 0; loc_i < locFCALShowers.size(); ++loc_i)
141  {
142  shared_ptr<DFCALShowerMatchParams>locShowerMatchParams;
143  if(locParticleID->Cut_MatchDistance(extrapolations, locFCALShowers[loc_i], locInputStartTime, locShowerMatchParams))
144  locDetectorMatches->Add_Match(locTrackTimeBased, locFCALShowers[loc_i], locShowerMatchParams);
145  }
146 }
147 
148 void DDetectorMatches_factory::MatchToSC(const DParticleID* locParticleID, const DTrackTimeBased* locTrackTimeBased, const vector<const DSCHit*>& locSCHits, DDetectorMatches* locDetectorMatches) const
149 {
150  vector<DTrackFitter::Extrapolation_t> extrapolations=locTrackTimeBased->extrapolations.at(SYS_START);
151  if (extrapolations.size()==0) return;
152 
153  double locInputStartTime = locTrackTimeBased->t0();
154  for(size_t loc_i = 0; loc_i < locSCHits.size(); ++loc_i)
155  {
156  shared_ptr<DSCHitMatchParams>locSCHitMatchParams;
157  if(locParticleID->Cut_MatchDistance(extrapolations, locSCHits[loc_i], locInputStartTime, locSCHitMatchParams, true))
158  locDetectorMatches->Add_Match(locTrackTimeBased, locSCHits[loc_i], locSCHitMatchParams);
159  }
160 }
161 
162 void DDetectorMatches_factory::MatchToDIRC(const DParticleID* locParticleID, const DTrackTimeBased* locTrackTimeBased, const vector<const DDIRCPmtHit*>& locDIRCHits, DDetectorMatches* locDetectorMatches, const vector<const DDIRCTruthBarHit*>& locDIRCBarHits) const
163 {
164  vector<DTrackFitter::Extrapolation_t> extrapolations=locTrackTimeBased->extrapolations.at(SYS_DIRC);
165  if(extrapolations.size()==0) return;
166 
167  double locInputStartTime = locTrackTimeBased->t0();
168 
169  // objects to hold DIRC match parameters and links between tracks and DIRC hits
170  shared_ptr<DDIRCMatchParams> locDIRCMatchParams;
171  map<shared_ptr<const DDIRCMatchParams>, vector<const DDIRCPmtHit*> > locDIRCTrackMatchParams;
172  locDetectorMatches->Get_DIRCTrackMatchParamsMap(locDIRCTrackMatchParams);
173 
174  // run DIRC LUT algorithm and add detector match
175  if(locParticleID->Cut_MatchDIRC(extrapolations, locDIRCHits, locInputStartTime, locTrackTimeBased->PID(), locDIRCMatchParams, locDIRCBarHits, locDIRCTrackMatchParams))
176  locDetectorMatches->Add_Match(locTrackTimeBased, locDIRCMatchParams);
177 }
178 
179 void DDetectorMatches_factory::MatchToTrack(const DParticleID* locParticleID, const DBCALShower* locBCALShower, const vector<const DTrackTimeBased*>& locTrackTimeBasedVector, DDetectorMatches* locDetectorMatches) const
180 {
181  double locMinDistance = 999.0;
182  double locFinalDeltaPhi = 999.0, locFinalDeltaZ = 999.0;
183  for(size_t loc_i = 0; loc_i < locTrackTimeBasedVector.size(); ++loc_i)
184  {
185  shared_ptr<DBCALShowerMatchParams> locShowerMatchParams;
186  double locInputStartTime = locTrackTimeBasedVector[loc_i]->t0();
187 
188  map<DetectorSystem_t,vector<DTrackFitter::Extrapolation_t> >extrapolations=locTrackTimeBasedVector[loc_i]->extrapolations;
189  if (extrapolations.size()==0) continue;
190 
191  if(!locParticleID->Distance_ToTrack(extrapolations.at(SYS_BCAL), locBCALShower, locInputStartTime, locShowerMatchParams))
192  continue;
193 
194  double locRSq = locBCALShower->x*locBCALShower->x + locBCALShower->y*locBCALShower->y;
195  double locDeltaPhi = locShowerMatchParams->dDeltaPhiToShower;
196  double locDeltaZ = locShowerMatchParams->dDeltaZToShower;
197  double locDistance = sqrt(locDeltaZ*locDeltaZ + locDeltaPhi*locDeltaPhi*locRSq);
198  if(locDistance >= locMinDistance)
199  continue;
200 
201  locMinDistance = locDistance;
202  locFinalDeltaPhi = locDeltaPhi;
203  locFinalDeltaZ = locDeltaZ;
204  }
205  locDetectorMatches->Set_DistanceToNearestTrack(locBCALShower, locFinalDeltaPhi, locFinalDeltaZ);
206 }
207 
208 void DDetectorMatches_factory::MatchToTrack(const DParticleID* locParticleID, const DFCALShower* locFCALShower, const vector<const DTrackTimeBased*>& locTrackTimeBasedVector, DDetectorMatches* locDetectorMatches) const
209 {
210  double locMinDistance = 999.0;
211  for(size_t loc_i = 0; loc_i < locTrackTimeBasedVector.size(); ++loc_i)
212  {
213  map<DetectorSystem_t,vector<DTrackFitter::Extrapolation_t> >extrapolations=locTrackTimeBasedVector[loc_i]->extrapolations;
214  if (extrapolations.size()==0) return;
215 
216  shared_ptr<DFCALShowerMatchParams> locShowerMatchParams;
217  double locInputStartTime = locTrackTimeBasedVector[loc_i]->t0();
218  if(!locParticleID->Distance_ToTrack(extrapolations.at(SYS_FCAL), locFCALShower, locInputStartTime, locShowerMatchParams))
219  continue;
220  if(locShowerMatchParams->dDOCAToShower < locMinDistance)
221  locMinDistance = locShowerMatchParams->dDOCAToShower;
222  }
223  locDetectorMatches->Set_DistanceToNearestTrack(locFCALShower, locMinDistance);
224 }
double t0(void) const
Definition: DTrackingData.h:23
void MatchToTOF(const DParticleID *locParticleID, const DTrackTimeBased *locTrackTimeBased, const vector< const DTOFPoint * > &locTOFPoints, DDetectorMatches *locDetectorMatches) const
bool Get_DIRCTrackMatchParamsMap(map< shared_ptr< const DDIRCMatchParams >, vector< const DDIRCPmtHit * > > &locDIRCTrackMatchParamsMap)
bool Cut_MatchDIRC(const vector< DTrackFitter::Extrapolation_t > &extrapolations, const vector< const DDIRCPmtHit * > locDIRCHits, double locInputStartTime, Particle_t locPID, shared_ptr< DDIRCMatchParams > &locDIRCMatchParams, const vector< const DDIRCTruthBarHit * > locDIRCBarHits, map< shared_ptr< const DDIRCMatchParams >, vector< const DDIRCPmtHit * > > &locDIRCTrackMatchParams, DVector3 *locOutputProjPos=nullptr, DVector3 *locOutputProjMom=nullptr) const
double Calc_SCFlightTimePCorrelation(const DTrackingData *locTrack, const DDetectorMatches *locDetectorMatches) const
Definition: GlueX.h:19
void Add_Match(const DTrackingData *locTrack, const DBCALShower *locBCALShower, const shared_ptr< const DBCALShowerMatchParams > &locShowerMatchParams)
double Calc_TOFFlightTimePCorrelation(const DTrackingData *locTrack, DDetectorMatches *locDetectorMatches) const
double Calc_FCALFlightTimePCorrelation(const DTrackingData *locTrack, DDetectorMatches *locDetectorMatches) const
double Calc_BCALFlightTimePCorrelation(const DTrackingData *locTrack, DDetectorMatches *locDetectorMatches) const
Definition: GlueX.h:20
void Set_DistanceToNearestTrack(const DBCALShower *locBCALShower, double locDeltaPhi, double locDeltaZ)
Definition: GlueX.h:22
DDetectorMatches * Create_DDetectorMatches(jana::JEventLoop *locEventLoop, vector< const DTrackTimeBased * > &locTrackTimeBasedVector)
Definition: GlueX.h:26
jerror_t init(void)
Called once at program start.
void MatchToDIRC(const DParticleID *locParticleID, const DTrackTimeBased *locTrackTimeBased, const vector< const DDIRCPmtHit * > &locDIRCHits, DDetectorMatches *locDetectorMatches, const vector< const DDIRCTruthBarHit * > &locDIRCBarHits) const
bool Cut_MatchDistance(const DReferenceTrajectory *rt, const DBCALShower *locBCALShower, double locInputStartTime, shared_ptr< DBCALShowerMatchParams > &locShowerMatchParams, DVector3 *locOutputProjPos=nullptr, DVector3 *locOutputProjMom=nullptr) const
jerror_t evnt(jana::JEventLoop *locEventLoop, uint64_t eventnumber)
Called every event.
double sqrt(double)
void MatchToFCAL(const DParticleID *locParticleID, const DTrackTimeBased *locTrackTimeBased, const vector< const DFCALShower * > &locFCALShowers, DDetectorMatches *locDetectorMatches) const
void MatchToBCAL(const DParticleID *locParticleID, const DTrackTimeBased *locTrackTimeBased, const vector< const DBCALShower * > &locBCALShowers, DDetectorMatches *locDetectorMatches) const
void MatchToTrack(const DParticleID *locParticleID, const DBCALShower *locBCALShower, const vector< const DTrackTimeBased * > &locTrackTimeBasedVector, DDetectorMatches *locDetectorMatches) const
jerror_t brun(jana::JEventLoop *locEventLoop, int32_t runnumber)
Called everytime a new run number is detected.
map< DetectorSystem_t, vector< DTrackFitter::Extrapolation_t > > extrapolations
bool Distance_ToTrack(const DReferenceTrajectory *rt, const DFCALShower *locFCALShower, double locInputStartTime, shared_ptr< DFCALShowerMatchParams > &locShowerMatchParams, DVector3 *locOutputProjPos=nullptr, DVector3 *locOutputProjMom=nullptr) const
Definition: DParticleID.cc:737
void Set_FlightTimePCorrelation(const DTrackingData *locTrack, DetectorSystem_t locDetectorSystem, double locCorrelation)
Particle_t PID(void) const
void MatchToSC(const DParticleID *locParticleID, const DTrackTimeBased *locTrackTimeBased, const vector< const DSCHit * > &locSCHits, DDetectorMatches *locDetectorMatches) const