Class Reference for E1039 Core & Analysis Software
GFTrack.cxx
Go to the documentation of this file.
1 #include "GFTrack.h"
2 
3 #include <iostream>
4 #include <algorithm>
5 
6 #include <TVectorD.h>
7 #include <TMatrixDSym.h>
8 
9 #include <GenFit/MeasuredStateOnPlane.h>
10 #include <GenFit/MeasurementOnPlane.h>
11 #include <GenFit/WireMeasurement.h>
12 #include <GenFit/PlanarMeasurement.h>
13 #include <GenFit/RKTrackRep.h>
14 #include <GenFit/AbsFitterInfo.h>
15 #include <GenFit/KalmanFitterInfo.h>
16 #include <GenFit/StateOnPlane.h>
17 
18 #include <phool/recoConsts.h>
19 
20 #include "SRawEvent.h"
21 
22 namespace
23 {
24  //static flag to indicate the initialized has been done
25  static bool inited = false;
26 
27  static double Z_TARGET;
28  static double Z_DUMP;
29  static double Z_UPSTREAM;
30 
31  static double X_BEAM;
32  static double Y_BEAM;
33  static double SIGX_BEAM;
34  static double SIGY_BEAM;
35 
36  //initialize global variables
37  void initGlobalVariables()
38  {
39  if(!inited)
40  {
41  inited = true;
42 
44  Z_TARGET = rc->get_DoubleFlag("Z_TARGET");
45  Z_DUMP = rc->get_DoubleFlag("Z_DUMP");
46  Z_UPSTREAM = rc->get_DoubleFlag("Z_UPSTREAM");
47 
48  X_BEAM = rc->get_DoubleFlag("X_BEAM");
49  Y_BEAM = rc->get_DoubleFlag("Y_BEAM");
50  SIGX_BEAM = rc->get_DoubleFlag("SIGX_BEAM");
51  SIGY_BEAM = rc->get_DoubleFlag("SIGY_BEAM");
52  }
53  }
54 };
55 
56 namespace SQGenFit
57 {
58 
59 GFTrack::GFTrack(): _track(nullptr), _trkrep(nullptr), _propState(nullptr), _virtMeas(nullptr), _trkcand(nullptr), _pdg(0)
60 {
61  initGlobalVariables();
62 }
63 
64 GFTrack::GFTrack(SRecTrack& recTrack): _track(nullptr), _trkrep(nullptr), _propState(nullptr), _virtMeas(nullptr), _trkcand(nullptr), _pdg(0)
65 {
66  initGlobalVariables();
67 
68  _pdg = recTrack.getCharge() > 0 ? -13 : 13;
69  _trkrep = new genfit::RKTrackRep(_pdg);
70 
71  TVector3 seed_mom = recTrack.getMomentumVecSt1();
72  TVector3 seed_pos = recTrack.getPositionVecSt1();
73 
74  TVectorD seed_state(6);
75  seed_state[0] = seed_pos.X();
76  seed_state[1] = seed_pos.Y();
77  seed_state[2] = seed_pos.Z();
78  seed_state[3] = seed_mom.Px();
79  seed_state[4] = seed_mom.Py();
80  seed_state[5] = seed_mom.Pz();
81 
82  TMatrixDSym seed_cov(6);
83  double uncertainty[6] = {10., 10., 10., 3., 3., 10.};
84  for(int i = 0; i < 6; i++)
85  {
86  for(int j = 0; j < 6; j++)
87  {
88  seed_cov[i][j] = uncertainty[i]*uncertainty[j];
89  }
90  }
91  _track = new genfit::Track(_trkrep, seed_state, seed_cov);
92 
93  _fitstates.clear();
94  int nHits = recTrack.getNHits();
95  for(int i = 0; i < nHits; ++i)
96  {
97  genfit::SharedPlanePtr detPlane(new genfit::DetPlane(recTrack.getGFPlaneO(i), recTrack.getGFPlaneU(i), recTrack.getGFPlaneV(i)));
98 
99  _fitstates.push_back(genfit::MeasuredStateOnPlane(recTrack.getGFState(i), recTrack.getGFCov(i), detPlane, _trkrep, recTrack.getGFAuxInfo(i)));
100  }
101 }
102 
103 GFTrack::GFTrack(Tracklet& tracklet): _track(nullptr), _trkrep(nullptr), _propState(nullptr), _virtMeas(nullptr), _trkcand(nullptr), _pdg(0)
104 {
105  initGlobalVariables();
106  setTracklet(tracklet, 590., false);
107 }
108 
110 {
111  if(_track != nullptr) delete _track;
112 }
113 
114 void GFTrack::setVerbosity(unsigned int v)
115 {
116  if(_trkrep != nullptr) _trkrep->setDebugLvl(v);
117 }
118 
119 void GFTrack::addMeasurements(std::vector<GFMeasurement*>& measurements)
120 {
121  for(auto iter = measurements.begin(); iter != measurements.end(); ++iter)
122  {
123  if(!(*iter)->isEnabled()) continue;
124  addMeasurement(*iter);
125  }
126 }
127 
129 {
130  measurement->setTrackPtr(this);
131 
132  genfit::TrackPoint* tp = new genfit::TrackPoint(measurement, _track);
133  tp->setSortingParameter(measurement->getZ());
134  _track->insertPoint(tp);
135 }
136 
138 {
139  genfit::AbsTrackRep* rep = _track->getCardinalRep();
140  if(rep)
141  {
142  genfit::FitStatus* fs = _track->getFitStatus(rep);
143  if(fs) return fs->getChi2();
144  }
145 
146  return -1.;
147 }
148 
150 {
151  genfit::AbsTrackRep* rep = _track->getCardinalRep();
152  if(rep)
153  {
154  genfit::FitStatus* fs = _track->getFitStatus(rep);
155  if(fs) return fs->getNdf();
156  }
157 
158  return -1.;
159 }
160 
162 {
163  double z = meas->getZ();
164  if(z < _measurements.front()->getZ())
165  {
166  return 0;
167  }
168  else if(z > _measurements.back()->getZ())
169  {
170  return _measurements.size() - 1;
171  }
172  else
173  {
174  auto iter = std::lower_bound(_measurements.begin(), _measurements.end(), meas, [](GFMeasurement* a, GFMeasurement* b) { return (a->getZ() < b->getZ()); });
175  int idx = std::distance(_measurements.begin(), iter);
176  return fabs(z - _measurements[idx]->getZ()) < fabs(z - _measurements[idx-1]->getZ()) ? idx : idx-1;
177  }
178 }
179 
180 double GFTrack::extrapolateToLine(TVector3& endPoint1, TVector3& endPoint2, const int startPtID)
181 {
182  TVector3 linePoint = (endPoint2 + endPoint1);
183  TVector3 lineDir = endPoint2 - endPoint1;
184 
185  if(!setInitialStateForExtrap(startPtID)) return -9999.;
186 
187  genfit::AbsTrackRep* rep = _track->getCardinalRep();
188  double len = rep->extrapolateToLine(*_propState, linePoint, lineDir);
189 
190  TVectorD hitcoord(7);
191  hitcoord[0] = endPoint1.X();
192  hitcoord[1] = endPoint1.Y();
193  hitcoord[2] = endPoint1.Z();
194  hitcoord[3] = endPoint2.X();
195  hitcoord[4] = endPoint2.Y();
196  hitcoord[5] = endPoint2.Z();
197  hitcoord[6] = 0.;
198  TMatrixDSym tempcov(7);
199  tempcov.Zero();
200  _virtMeas.reset(new genfit::WireMeasurement(hitcoord, tempcov, 999, 999, nullptr));
201 
202  return len;
203 }
204 
205 double GFTrack::extrapolateToPlane(TVector3& pO, TVector3& pU, TVector3& pV, const int startPtID)
206 {
207  genfit::SharedPlanePtr destPlane(new genfit::DetPlane(pO, pU, pV));
208 
209  if(!setInitialStateForExtrap(startPtID)) return -9999.;
210 
211  genfit::AbsTrackRep* rep = _track->getCardinalRep();
212  double len = rep->extrapolateToPlane(*_propState, destPlane);
213 
214  TVectorD hitcoord(2);
215  hitcoord[0] = pO.X();
216  hitcoord[1] = pO.Y();
217  TMatrixDSym tempcov(2);
218  tempcov.UnitMatrix();
219 
220  genfit::PlanarMeasurement* pMeas = new genfit::PlanarMeasurement(hitcoord, tempcov, 998, 998, nullptr);
221  pMeas->setPlane(destPlane);
222  _virtMeas.reset(pMeas);
223 
224  return len;
225 }
226 
227 double GFTrack::extrapolateToPoint(TVector3& point, bool update, const int startPtID)
228 {
229  if(!setInitialStateForExtrap(startPtID)) return -9999.;
230 
231  genfit::AbsTrackRep* rep = _track->getCardinalRep();
232  double len = rep->extrapolateToPoint(*_propState, point);
233 
234  //TODO: implement the virtual spacepoint measurement here
235 
236  return len;
237 }
238 
239 double GFTrack::updatePropState(const TVectorD& meas, const TMatrixDSym& V)
240 {
241  //get the H matrix from measurement
242  std::unique_ptr<const genfit::AbsHMatrix> H(_virtMeas->constructHMatrix(_propState->getRep()));
243 
244  //Handcrafted KF
245  TVectorD stateVec(_propState->getState());
246  TMatrixDSym cov(_propState->getCov());
247  double chi2 = 9999.;
248  double ndf = meas.GetNrows();
249 
250  TVectorD res(meas - H->Hv(stateVec));
251  TMatrixDSym covSumInv(cov);
252  H->HMHt(covSumInv);
253  covSumInv += V;
254  genfit::tools::invertMatrix(covSumInv);
255 
256  TMatrixD CHt(H->MHt(cov));
257  TVectorD update(TMatrixD(CHt, TMatrixD::kMult, covSumInv)*res);
258 
259  stateVec += update;
260  covSumInv.Similarity(CHt);
261  cov -= covSumInv;
262  _propState->setStateCov(stateVec, cov);
263 
264  //calc chi2
265  TVectorD resNew(meas - H->Hv(stateVec));
266  TMatrixDSym HCHt(cov);
267  H->HMHt(HCHt);
268  HCHt -= V;
269  HCHt *= -1;
270  genfit::tools::invertMatrix(HCHt);
271  chi2 = HCHt.Similarity(resNew);
272 
273  return chi2/ndf;
274 }
275 
276 bool GFTrack::setInitialStateForExtrap(const int startPtID)
277 {
278  genfit::AbsTrackRep* rep = _track->getCardinalRep();
279  if(_track->getNumPoints() > 0)
280  {
281  genfit::TrackPoint* tp = _track->getPointWithMeasurementAndFitterInfo(startPtID, rep);
282  if(tp == nullptr) return false;
283 
284  genfit::KalmanFitterInfo* info = static_cast<genfit::KalmanFitterInfo*>(tp->getFitterInfo(rep));
285  _propState.reset(new genfit::MeasuredStateOnPlane(info->getFittedState()));
286  }
287  else
288  {
289  _propState.reset(new genfit::MeasuredStateOnPlane(_fitstates[startPtID]));
290  }
291 
292  return true;
293 }
294 
295 double GFTrack::swimToVertex(double z, TVector3* pos, TVector3* mom, TMatrixDSym* cov)
296 {
297  //Basic constants
298  TVector3 pU(1., 0., 0.);
299  TVector3 pV(0., 1., 0.);
300  TVector3 pO(0., 0., z );
301 
302  TVectorD beamCenter(2);
303  beamCenter[0] = X_BEAM; beamCenter[1] = Y_BEAM;
304  TMatrixDSym beamCov(2);
305  beamCov.Zero();
306  beamCov(0, 0) = SIGX_BEAM*SIGX_BEAM; beamCov(1, 1) = SIGY_BEAM*SIGY_BEAM;
307 
308  double chi2 = -1.;
309  try
310  {
311  double len = extrapolateToPlane(pO, pU, pV);
312  if(fabs(len) > 6000.) throw len;
313 
314  chi2 = updatePropState(beamCenter, beamCov);
315  if(cov != nullptr)
316  {
317  getExtrapPosMomCov(*pos, *mom, *cov);
318  }
319  else if(mom != nullptr)
320  {
321  getExtrapPosMom(*pos, *mom);
322  }
323  }
324  catch(genfit::Exception& e)
325  {
326  std::cerr << __FILE__ << " " << __LINE__ << ": hypo test failed vertex @Z=" << z << ": " << e.what() << std::endl;
327  return -1.;
328  }
329  catch(double len)
330  {
331  std::cerr << __FILE__ << " " << __LINE__ << ": hypo test failed vertex @Z=" << z << ": " << len << std::endl;
332  return -1.;
333  }
334 
335  return chi2;
336 }
337 
338 void GFTrack::setTracklet(Tracklet& tracklet, double z_reference, bool wildseedcov)
339 {
340  _trkcand = &tracklet;
341  _pdg = tracklet.getCharge() > 0 ? -13 : 13;
342  _trkrep = new genfit::RKTrackRep(_pdg);
343 
344  TVectorD seed_state(6);
345 
346  TVector3 seed_mom = tracklet.getExpMomentum(z_reference);
347  TVector3 seed_pos(tracklet.getExpPositionX(z_reference), tracklet.getExpPositionY(z_reference), z_reference);
348  seed_state[0] = seed_pos.X();
349  seed_state[1] = seed_pos.Y();
350  seed_state[2] = seed_pos.Z();
351  seed_state[3] = seed_mom.Px();
352  seed_state[4] = seed_mom.Py();
353  seed_state[5] = seed_mom.Pz();
354 
355  TMatrixDSym seed_cov(6);
356  double uncertainty[6] = {10., 10., 10., 3., 3., 10.};
357  for(int i = 0; i < 6; i++)
358  {
359  for(int j = 0; j < 6; j++)
360  {
361  seed_cov[i][j] = uncertainty[i]*uncertainty[j];
362  }
363  }
364 
365  if(!wildseedcov)
366  {
367  //Implement a smaller cov based on the chi^2 fit?
368  }
369  _track = new genfit::Track(_trkrep, seed_state, seed_cov);
370 
371  _measurements.clear();
372  for(auto iter = tracklet.hits.begin(); iter != tracklet.hits.end(); ++iter)
373  {
374  if(iter->hit.index < 0) continue;
375 
376  GFMeasurement* meas = new GFMeasurement(*iter);
377  _measurements.push_back(meas);
378  }
379  std::sort(_measurements.begin(), _measurements.end(), [](GFMeasurement* a, GFMeasurement* b) { return (a->getZ() < b->getZ()); });
380 
381  addMeasurements(_measurements);
383 }
384 
385 void GFTrack::postFitUpdate(bool updateMeasurements)
386 {
387  if(!updateMeasurements) return;
388  for(auto iter = _measurements.begin(); iter != _measurements.end(); ++iter)
389  {
390  (*iter)->postFitUpdate();
391  }
392 }
393 
395 {
396  //postFitUpdate();
397  //The following steps are pretty hacky and should be considered as only a temporary solution
398  SRecTrack strack;
399  strack.setChisq(getChi2());
400  for(auto iter = _measurements.begin(); iter != _measurements.end(); ++iter)
401  {
402  strack.insertHitIndex((*iter)->getBeforeFitHit().hit.index);
403 
404  const genfit::MeasuredStateOnPlane& fitstate = (*iter)->getTrackPoint()->getKalmanFitterInfo()->getFittedState(true);
405  strack.insertGFState(fitstate);
406 
407  TVector3 pos, mom;
408  fitstate.getPosMom(pos, mom);
409  TMatrixD stateVec(5, 1);
410  stateVec[0][0] = getCharge()/mom.Mag();
411  stateVec[1][0] = mom.Px()/mom.Pz();
412  stateVec[2][0] = mom.Py()/mom.Pz();
413  stateVec[3][0] = pos.X();
414  stateVec[4][0] = pos.Y();
415  strack.insertStateVector(stateVec);
416  strack.insertZ(pos.Z());
417 
418  TMatrixD cov(fitstate.getCov());
419  strack.insertCovariance(cov);
420 
421  strack.insertChisq((*iter)->getTrackPoint()->getKalmanFitterInfo()->getSmoothedChi2());
422  }
423 
424  //Swim to various places and save info
425  strack.swimToVertex(nullptr, nullptr, false);
426 
427  //Hypothesis test should be implemented here
428  //test Z_UPSTREAM
429  strack.setChisqUpstream(swimToVertex(Z_UPSTREAM));
430 
431  //test Z_TARGET
432  strack.setChisqTarget(swimToVertex(Z_TARGET));
433 
434  //test Z_DUMP
435  strack.setChisqDump(swimToVertex(Z_DUMP));
436 
437  /*
438  //Find POCA to beamline -- it seems to be funky and mostly found some place way upstream or downstream
439  // most likely because the cross product of the track direction and beam line direction is way too small
440  // z axis to provide reasonable calculation of the POCA location. It's disabled for now.
441  TVector3 ep1(0., 0., -499.);
442  TVector3 ep2(0., 0., 0.);
443  try
444  {
445  extrapolateToLine(ep1, ep2);
446  TVectorD beamR(1); beamR(0) = 0.;
447  TMatrixDSym beamC(1); beamC(0, 0) = 1000.;
448  strack.setChisqVertex(updatePropState(beamR, beamC));
449  }
450  catch(genfit::Exception& e)
451  {
452  std::cerr << "Hypothesis test failed at beamline: " << e.what() << std::endl;
453  print(0);
454  }
455  */
456 
457  //test Z_VERTEX
458  TVector3 pos, mom;
459  strack.setChisqVertex(swimToVertex(strack.getVertexPos().Z(), &pos, &mom));
460  strack.setVertexPos(pos);
461  strack.setVertexMom(mom);
462 
463  strack.setKalmanStatus(1);
464  return strack;
465 }
466 
467 void GFTrack::print(unsigned int debugLvl)
468 {
469  std::cout << "=============== SGTrack ===============" << std::endl;
470  std::cout << "------------ Track candidate ----------" << std::endl;
471  _trkcand->print();
472 
473  std::cout << "------------ Fit status ---------------" << std::endl;
474  _track->getFitStatus(_trkrep)->Print();
475 
476  std::cout << "------------ Fit Result ---------------" << std::endl;
477  _track->getFittedState().Print(); //default to the first hit
478 
479  if(debugLvl < 1) return;
480  for(auto iter = _measurements.begin(); iter != _measurements.end(); ++iter)
481  {
482  (*iter)->print(debugLvl-1);
483  }
484 
485  if(debugLvl < 20) return;
486  std::cout << "------------ GenFit Track ---------------" << std::endl;
487  _track->Print();
488 }
489 
490 }
virtual double get_DoubleFlag(const std::string &name) const
Definition: PHFlag.cc:49
void setTrackPtr(GFTrack *trackPtr)
double getNDF()
Definition: GFTrack.cxx:149
void print(unsigned int debugLvl=0)
Definition: GFTrack.cxx:467
double extrapolateToPlane(TVector3 &pO, TVector3 &pU, TVector3 &pV, const int startPtID=0)
Definition: GFTrack.cxx:205
void getExtrapPosMomCov(TVector3 &pos, TVector3 &mom, TMatrixDSym &cov)
Definition: GFTrack.h:45
void addMeasurements(std::vector< GFMeasurement * > &measurements)
Definition: GFTrack.cxx:119
void setVerbosity(unsigned int v)
Definition: GFTrack.cxx:114
void postFitUpdate(bool updateMeasurements=true)
Definition: GFTrack.cxx:385
int getNearestMeasurementID(GFMeasurement *meas)
Definition: GFTrack.cxx:161
double swimToVertex(double z, TVector3 *pos=nullptr, TVector3 *mom=nullptr, TMatrixDSym *cov=nullptr)
Definition: GFTrack.cxx:295
void addMeasurement(GFMeasurement *measurement)
Definition: GFTrack.cxx:128
SRecTrack getSRecTrack()
Definition: GFTrack.cxx:394
double extrapolateToLine(TVector3 &endPoint1, TVector3 &endPoint2, const int startPtID=0)
Definition: GFTrack.cxx:180
void getExtrapPosMom(TVector3 &pos, TVector3 &mom)
Definition: GFTrack.h:46
double getChi2()
Definition: GFTrack.cxx:137
void setTracklet(Tracklet &tracklet, double z_reference=590., bool wildseedcov=false)
Definition: GFTrack.cxx:338
double extrapolateToPoint(TVector3 &point, bool update=false, const int startPtID=0)
Definition: GFTrack.cxx:227
double updatePropState(const TVectorD &meas, const TMatrixDSym &V)
Definition: GFTrack.cxx:239
int getCharge()
Definition: GFTrack.h:36
void checkConsistency()
Definition: GFTrack.h:50
Int_t getNHits() const
Definition: SRecEvent.h:102
Int_t getCharge() const
Gets.
Definition: SRecEvent.h:101
TVector3 getVertexPos()
Definition: SRecEvent.h:196
void setVertexPos(TVector3 pos)
Definition: SRecEvent.h:208
void setChisqDump(Double_t chisq)
Definition: SRecEvent.h:215
TVectorD getGFAuxInfo(Int_t i)
Definition: SRecEvent.h:117
void insertChisq(Double_t chisq)
Definition: SRecEvent.h:159
void setChisq(Double_t chisq)
Sets.
Definition: SRecEvent.h:154
void insertCovariance(TMatrixD covar)
Definition: SRecEvent.h:157
void insertGFState(const genfit::MeasuredStateOnPlane &msop)
Definition: SRecEvent.cxx:351
TVector3 getMomentumVecSt1() const
Definition: SRecEvent.h:129
void setChisqTarget(Double_t chisq)
Definition: SRecEvent.h:216
TVector3 getPositionVecSt1() const
Definition: SRecEvent.h:137
TVectorD getGFState(Int_t i)
Definition: SRecEvent.h:118
TVector3 getGFPlaneU(Int_t i)
Definition: SRecEvent.h:115
void setVertexMom(TVector3 mom)
Definition: SRecEvent.h:214
TVector3 getGFPlaneO(Int_t i)
Definition: SRecEvent.h:114
void insertZ(Double_t z)
Definition: SRecEvent.h:158
void setChisqUpstream(Double_t chisq)
Definition: SRecEvent.h:217
TMatrixDSym getGFCov(Int_t i)
Definition: SRecEvent.h:119
void swimToVertex(TVector3 *pos=nullptr, TVector3 *mom=nullptr, bool hyptest=true)
Simple swim to vertex.
Definition: SRecEvent.cxx:402
void insertHitIndex(Int_t index)
Definition: SRecEvent.h:155
void setKalmanStatus(Int_t status)
Definition: SRecEvent.h:148
TVector3 getGFPlaneV(Int_t i)
Definition: SRecEvent.h:116
void insertStateVector(TMatrixD state)
Definition: SRecEvent.h:156
void setChisqVertex(Double_t chisq)
Definition: SRecEvent.h:218
double getExpPositionY(double z) const
std::list< SignedHit > hits
Definition: FastTracklet.h:228
TVector3 getExpMomentum(double z) const
double getExpPositionX(double z) const
void print(std::ostream &os=std::cout)
int getCharge() const
Return the charge (+1 or -1) of this tracklet.
static recoConsts * instance()
Definition: recoConsts.cc:7