Class Reference for E1039 Core & Analysis Software
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
VertexFit.cxx
Go to the documentation of this file.
1 /*
2 VertexFit.cxx
3 
4 Implementation of the VertexFit, in the beginning only a old-fashioned primary vertex finder
5 is implemented
6 
7 Author: Kun Liu, liuk@fnal.gov
8 Created: 2-8-2012
9 */
10 
11 #include "VertexFit.h"
12 #include "KalmanTrack.h"
13 #include "KalmanFilter.h"
14 #include "GenFitExtrapolator.h"
15 
17 #include <fun4all/PHTFileServer.h>
18 #include <phool/PHNodeIterator.h>
19 #include <phool/PHIODataNode.h>
20 #include <phool/getClass.h>
21 #include <phool/recoConsts.h>
22 #include <phfield/PHFieldConfig_v3.h>
23 #include <phfield/PHFieldUtility.h>
24 #include <phgeom/PHGeomUtility.h>
25 
26 #include <TMatrixD.h>
27 
28 #include <iostream>
29 #include <cmath>
30 #include <memory>
31 #include <fstream>
32 
33 using namespace std;
34 
35 namespace
36 {
37  //static flag to indicate the initialized has been done
38  static bool inited = false;
39 
40  //static flag of kmag strength
41  static double FMAGSTR;
42  static double KMAGSTR;
43 
44  //Beam position and shape
45  static double X_BEAM;
46  static double Y_BEAM;
47  static double SIGX_BEAM;
48  static double SIGY_BEAM;
49 
50  //Simple swimming settings
51  static int NSTEPS_TARGET;
52  static int NSTEPS_SHIELDING;
53  static int NSTEPS_FMAG;
54 
55  //Geometric constants
56  static double Z_TARGET;
57  static double Z_DUMP;
58  static double Z_UPSTREAM;
59  static double Z_DOWNSTREAM;
60 
61  //initialize global variables
62  void initGlobalVariables()
63  {
64  if(!inited)
65  {
66  inited = true;
67 
69  FMAGSTR = rc->get_DoubleFlag("FMAGSTR");
70  KMAGSTR = rc->get_DoubleFlag("KMAGSTR");
71 
72  X_BEAM = rc->get_DoubleFlag("X_BEAM");
73  Y_BEAM = rc->get_DoubleFlag("Y_BEAM");
74  SIGX_BEAM = rc->get_DoubleFlag("SIGX_BEAM");
75  SIGY_BEAM = rc->get_DoubleFlag("SIGY_BEAM");
76 
77  NSTEPS_TARGET = rc->get_IntFlag("NSTEPS_TARGET");
78  NSTEPS_SHIELDING = rc->get_IntFlag("NSTEPS_SHIELDING");
79  NSTEPS_FMAG = rc->get_IntFlag("NSTEPS_FMAG");
80 
81  Z_TARGET = rc->get_DoubleFlag("Z_TARGET");
82  Z_DUMP = rc->get_DoubleFlag("Z_DUMP");
83  Z_UPSTREAM = rc->get_DoubleFlag("Z_UPSTREAM");
84  Z_DOWNSTREAM = rc->get_DoubleFlag("Z_DOWNSTREAM");
85  }
86  }
87 }
88 
89 VertexFit::VertexFit(const std::string& name) :
90  SubsysReco(name)
91 {
92  initGlobalVariables();
93 
95  TMatrixD m(2, 1), cov(2, 2), proj(2, 5);
96  m[0][0] = X_BEAM;
97  m[1][0] = Y_BEAM;
98 
99  cov.Zero();
100  cov[0][0] = SIGX_BEAM*SIGX_BEAM;
101  cov[1][1] = SIGY_BEAM*SIGY_BEAM;
102 
103  proj.Zero();
104  proj[0][3] = 1.;
105  proj[1][4] = 1.;
106 
107  _node_vertex.getMeasurement().ResizeTo(2, 1);
108  _node_vertex.getMeasurementCov().ResizeTo(2, 2);
109  _node_vertex.getProjector().ResizeTo(2, 5);
110 
111  _node_vertex.getMeasurement() = m;
112  _node_vertex.getMeasurementCov() = cov;
113  _node_vertex.getProjector() = proj;
114 
115  _max_iteration = 200;
116  _tolerance = .05;
117 
119  optimize = false;
120  //disable targert center fit by default
121  fit_target_center = false;
122 
124  evalFileName = "vertex_eval.root";
125  evalFile = nullptr;
126  evalTree = nullptr;
127 
128 }
129 
131 {
132  if(evalFile != nullptr)
133  {
134  evalFile->cd();
135  evalTree->Write();
136  evalFile->Close();
137  }
138 }
139 
141 
142  this->bookEvaluation(evalFileName.c_str());
143 
145 }
146 
148 
151 
152  int ret = MakeNodes(topNode);
153  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
154 
155  ret = GetNodes(topNode);
156  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
157 
158  ret = InitField(topNode);
159  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
160 
161  ret = InitGeom(topNode);
162  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
163 
164  PHField* field = PHFieldUtility::GetFieldMapNode(nullptr, topNode);
165  assert(field);
166 
167  if(verbosity > 2) {
168  cout << "PHField check: " << "-------" << endl;
169  std::ofstream field_scan("field_scan.csv");
170  field->identify(field_scan);
171  field_scan.close();
172  }
173 
174  _kmfit = KalmanFilter::instance();
175  _kmfit->enableDumpCorrection();
176  _extrapolator.init(field, _t_geo_manager);
177 
179  _extrapolator.setPropCalc(false);
180  _extrapolator.setLengthCalc(false);
181 
183 }
184 
186 
187  if(verbosity > 0) _recEvent->identify();
188 
189  _recEvent->clearDimuons(); // In case the vertexing is re-done.
190 
191  int ret = -99;
192  try {
193  ret = setRecEvent(_recEvent);
194  _recEvent->setRecStatus(ret);
196  LogInfo("Vertexing Returned: " << ret);
197  }
198  } catch (...) {
199  LogInfo("VertexFitting failed");
200  }
201 
203 }
204 
207 }
208 
209 int VertexFit::MakeNodes(PHCompositeNode* topNode) {
211 }
212 
213 int VertexFit::GetNodes(PHCompositeNode* topNode) {
214 
215  _recEvent = findNode::getClass<SRecEvent>(topNode, "SRecEvent");
216  if (!_recEvent) {
217  if(Verbosity() > 2) LogInfo("!_recEvent");
219  }
220 
222 }
223 
224 
225 int VertexFit::InitField(PHCompositeNode *topNode)
226 {
227  if (verbosity > 1) cout << "VertexFit::InitField" << endl;
228  PHField * phfield = PHFieldUtility::GetFieldMapNode(nullptr, topNode);
229 
230  if(phfield && verbosity > 1)
231  cout << "VertexFit::InitField - use filed from NodeTree." << endl;
232 
233  if(!phfield) {
234  if (verbosity > 1) cout << "VertexFit::InitField - create magnetic field setup" << endl;
235  unique_ptr<PHFieldConfig> default_field_cfg(nullptr);
236  default_field_cfg.reset(new PHFieldConfig_v3(recoConsts::instance()->get_CharFlag("fMagFile"), recoConsts::instance()->get_CharFlag("kMagFile"), recoConsts::instance()->get_DoubleFlag("FMAGSTR"), recoConsts::instance()->get_DoubleFlag("KMAGSTR"), 5.));
237  phfield = PHFieldUtility::GetFieldMapNode(default_field_cfg.get(), topNode, 0);
238  }
239 
240  assert(phfield);
241 
243 }
244 
245 int VertexFit::InitGeom(PHCompositeNode *topNode)
246 {
247  if (verbosity > 1) cout << "VertexFit::InitGeom" << endl;
248 
249  _t_geo_manager = PHGeomUtility::GetTGeoManager(topNode);
250 
251  if(_t_geo_manager && verbosity > 1)
252  cout << "VertexFit::InitGeom - use geom from NodeTree." << endl;
253 
254  assert(_t_geo_manager);
255 
257 }
258 
259 int VertexFit::setRecEvent(SRecEvent* recEvent, int sign1, int sign2)
260 {
261  std::vector<int> idx_pos = recEvent->getChargedTrackIDs(sign1);
262  std::vector<int> idx_neg = recEvent->getChargedTrackIDs(sign2);
263 
264  nPos = idx_pos.size();
265  nNeg = idx_neg.size();
266  if(nPos*nNeg == 0) return VFEXIT_FAIL_DIMUONPAIR;
267 
268  //Prepare evaluation output
269  if(evalTree != nullptr)
270  {
271  runID = recEvent->getRunID();
272  eventID = recEvent->getEventID();
273  targetPos = recEvent->getTargetPos();
274  }
275 
276  //Loop over all possible combinations
277  for(int i = 0; i < nPos; ++i)
278  {
279  if(!recEvent->getTrack(idx_pos[i]).isValid()) continue;
281  LogInfo("pos track OK");
282  }
283 
284  int j = sign1 + sign2 == 0 ? 0 : i + 1; // this is to avoid using same track twice in like-sign mode
285  for(; j < nNeg; ++j)
286  {
287  //Only needed for like-sign muons
288  if(idx_pos[i] == idx_neg[j]) continue;
290  LogInfo("two track OK");
291  }
292 
293  SRecTrack track_neg = recEvent->getTrack(idx_neg[j]);
294  if(!track_neg.isValid()) continue;
296  LogInfo("neg track OK");
297  }
298  SRecTrack track_pos = recEvent->getTrack(idx_pos[i]);
299 
300  SRecDimuon dimuon;
301  dimuon.trackID_pos = idx_pos[i];
302  dimuon.trackID_neg = idx_neg[j];
303 
304  dimuon.p_pos_single = track_pos.getMomentumVertex();
305  dimuon.p_neg_single = track_neg.getMomentumVertex();
306  dimuon.vtx_pos = track_pos.getVertex();
307  dimuon.vtx_neg = track_neg.getVertex();
308  dimuon.chisq_single = track_pos.getChisqVertex() + track_neg.getChisqVertex();
309 
310  //Start prepare the vertex fit
311  init();
312  if(KMAGSTR > 0)
313  {
314  addTrack(0, track_pos);
315  addTrack(1, track_neg);
316  }
317  else
318  {
319  addTrack(0, track_neg);
320  addTrack(1, track_pos);
321  }
322  addHypothesis(0.5*(dimuon.vtx_pos[2] + dimuon.vtx_neg[2]), 50.);
323  addHypothesis(findDimuonVertexFast(track_pos, track_neg), 50.);
324  choice_eval = processOnePair();
325 
326  //Fill the dimuon info which are not related to track refitting first
327  dimuon.chisq_vx = getVXChisq();
328  dimuon.vtx.SetXYZ(_vtxpar_curr._r[0][0], _vtxpar_curr._r[1][0], _vtxpar_curr._r[2][0]);
329 
330  dimuon.proj_target_pos = track_pos.getTargetPos();
331  dimuon.proj_dump_pos = track_pos.getDumpPos();
332  dimuon.proj_target_neg = track_neg.getTargetPos();
333  dimuon.proj_dump_neg = track_neg.getDumpPos();
334 
335  //Retrieve the results
336  double z_vertex_opt = getVertexZ0();
337  if(optimize)
338  {
339  //if(z_vertex_opt < -80. && getKFChisq() < 10.) z_vertex_opt = Z_TARGET;
340  if(dimuon.proj_target_pos.Perp() < dimuon.proj_dump_pos.Perp() && dimuon.proj_target_neg.Perp() < dimuon.proj_dump_neg.Perp())
341  {
342  int nTry = 0;
343  double z_curr = 9999.;
344  while(fabs(z_curr - z_vertex_opt) > 0.5 && nTry < 100)
345  {
346  z_curr = z_vertex_opt;
347  ++nTry;
348 
349  track_pos.setZVertex(z_vertex_opt);
350  track_neg.setZVertex(z_vertex_opt);
351 
352  double m = (track_pos.getMomentumVertex() + track_neg.getMomentumVertex()).M();
353  //z_vertex_opt = -189.6 + 17.71*m - 1.159*m*m; //parameterization of r1.4.0
354  //z_vertex_opt = -305.465 + 104.731*m - 24.3589*m*m + 2.5564*m*m*m - 0.0978876*m*m*m*m; //for E906 geomtery
355 
356  z_vertex_opt = -310.0540 + 4.3539*m - 0.9518*m*m + 0.0803*m*m*m ;//use this parametrization for E1039 geometry for now, more correction is due (Abi)
357 
358  //std::cout << nTry << " " << z_curr << " " << z_vertex_opt << " " << (track_pos.getMomentumVertex() + track_neg.getMomentumVertex()).M() << std::endl;
359  }
360  }
361  }
362 
363  if(fit_target_center) z_vertex_opt = Z_TARGET; //fit in the target center
364 
365  track_pos.setZVertex(z_vertex_opt);
366  track_neg.setZVertex(z_vertex_opt);
367  dimuon.p_pos = track_pos.getMomentumVertex();
368  dimuon.p_neg = track_neg.getMomentumVertex();
369  dimuon.chisq_kf = track_pos.getChisqVertex() + track_neg.getChisqVertex();
370 
371  /*
372  //If we are running in the like-sign mode, reverse one sign of px
373  if(sign1 + sign2 != 0)
374  {
375  if(dimuon.p_pos.Px() < 0)
376  {
377  dimuon.p_pos.SetPx(-dimuon.p_pos.Px());
378  dimuon.p_pos_single.SetPx(-dimuon.p_pos_single.Px());
379  }
380  if(dimuon.p_neg.Px() > 0)
381  {
382  dimuon.p_neg.SetPx(-dimuon.p_neg.Px());
383  dimuon.p_neg_single.SetPx(-dimuon.p_neg_single.Px());
384  }
385  if(dimuon.p_pos.Py()*dimuon.p_neg.Py() > 0)
386  {
387  dimuon.p_pos.SetPy(-dimuon.p_pos.Py());
388  dimuon.p_pos_single.SetPy(-dimuon.p_pos_single.Py());
389  }
390  }
391  */
392  dimuon.calcVariables();
393 
394  //Test three fixed hypothesis
395  dimuon.chisq_dump = track_pos.getChisqDump() + track_neg.getChisqDump();
396  dimuon.chisq_target = track_pos.getChisqTarget() + track_neg.getChisqTarget();
397  dimuon.chisq_upstream = track_pos.getChisqUpstream() + track_neg.getChisqUpstream();
398 
399  //Fill the final data
400  recEvent->insertDimuon(dimuon);
401 
402  //Fill the evaluation data
403  p_idx_eval = dimuon.trackID_pos;
404  m_idx_eval = dimuon.trackID_neg;
405  fillEvaluation();
406  }
407  }
408 
409  if(recEvent->getNDimuons() > 0) return VFEXIT_SUCCESS;
410  return VFEXIT_FAIL_ITERATION;
411 }
412 
414 {
415  //Swim both tracks all the way down, and store the numbers
416  TVector3 pos1[NSTEPS_FMAG + NSTEPS_TARGET + 1];
417  TVector3 pos2[NSTEPS_FMAG + NSTEPS_TARGET + 1];
418  TVector3 mom1[NSTEPS_FMAG + NSTEPS_TARGET + 1];
419  TVector3 mom2[NSTEPS_FMAG + NSTEPS_TARGET + 1];
420  track1.swimToVertex(pos1, mom1);
421  track2.swimToVertex(pos2, mom2);
422 
423  int iStep_min = -1;
424  double dist_min = 1E6;
425  int charge1 = track1.getCharge();
426  int charge2 = track2.getCharge();
427  for(int iStep = 0; iStep < NSTEPS_FMAG + NSTEPS_TARGET + 1; ++iStep)
428  {
429  double dist = (pos1[iStep] - pos2[iStep]).Perp();
430  if(dist < dist_min && FMAGSTR*charge1*mom1[iStep].Px() > 0 && FMAGSTR*charge2*mom2[iStep].Px() > 0)
431  {
432  iStep_min = iStep;
433  dist_min = dist;
434  }
435  }
436 
437  if(iStep_min == -1) return Z_DUMP;
438  return pos1[iStep_min].Z();
439 }
440 
442 {
443  _trkpar_curr.clear();
444 
445  z_vertex.clear();
446  r_vertex.clear();
447  chisq_km.clear();
448  chisq_vx.clear();
449 
450  z_start.clear();
451  sig_z_start.clear();
452 
454  //addHypothesis(Z_DUMP, 50.);
455  //addHypothesis(Z_TARGET, 50.);
456 }
457 
458 void VertexFit::setStartingVertex(double z_start, double sigz_start)
459 {
461  _vtxpar_curr._r[0][0] = X_BEAM;
462  _vtxpar_curr._r[1][0] = Y_BEAM;
463  _vtxpar_curr._r[2][0] = z_start;
464 
465  _vtxpar_curr._cov.Zero();
466  _vtxpar_curr._cov[0][0] = SIGX_BEAM*SIGX_BEAM;
467  _vtxpar_curr._cov[1][1] = SIGY_BEAM*SIGY_BEAM;
468  _vtxpar_curr._cov[2][2] = sigz_start*sigz_start;
469 
470  _chisq_vertex = 0.;
471  _chisq_kalman = 0.;
472 }
473 
475 {
476  double chisq_min = 1E6;
477  int index_min = -1;
478 
479  nStart = z_start.size();
480  for(int i = 0; i < nStart; ++i)
481  {
482 #ifdef _DEBUG_ON
483  LogInfo("Testing starting point: " << z_start[i]);
484 #endif
485 
486  setStartingVertex(z_start[i], sig_z_start[i]);
487  nIter_eval[i] = findVertex();
488 
489  chisq_km.push_back(_chisq_kalman);
490  chisq_vx.push_back(_chisq_vertex);
491  z_vertex.push_back(_vtxpar_curr._r[2][0]);
492  r_vertex.push_back(sqrt(_vtxpar_curr._r[0][0]*_vtxpar_curr._r[0][0] + _vtxpar_curr._r[1][0]*_vtxpar_curr._r[1][0]));
493 
494  if(_chisq_kalman < chisq_min)
495  {
496  index_min = i;
497  chisq_min = _chisq_kalman;
498  }
499  }
500 
501  if(index_min < 0) return 0;
502 
503  _chisq_kalman = chisq_km[index_min];
504  _chisq_vertex = chisq_vx[index_min];
505  _vtxpar_curr._r[2][0] = z_vertex[index_min];
506  if(z_vertex[index_min] > 1E4)
507  {
508  index_min = z_start.size();
509  _vtxpar_curr._r[2][0] = z_start.back();
510  }
511 
512  return index_min;
513 }
514 
515 void VertexFit::addTrack(int index, KalmanTrack& _track)
516 {
517  if(_track.getNodeList().front().isSmoothDone())
518  {
519  _trkpar_curr.push_back(_track.getNodeList().front().getSmoothed());
520  }
521  else if(_track.getNodeList().front().isFilterDone())
522  {
523  _trkpar_curr.push_back(_track.getNodeList().front().getFiltered());
524  }
525  else
526  {
527  _trkpar_curr.push_back(_track.getNodeList().front().getPredicted());
528  }
529 }
530 
531 void VertexFit::addTrack(int index, SRecTrack& _track)
532 {
533  TrkPar _trkpar;
534  _trkpar._state_kf = _track.getStateVector(0);
535  _trkpar._covar_kf = _track.getCovariance(0);
536  _trkpar._z = _track.getZ(0);
537 
538  _trkpar_curr.push_back(_trkpar);
539 }
540 
541 void VertexFit::addTrack(int index, TrkPar& _trkpar)
542 {
543  _trkpar_curr.push_back(_trkpar);
544 }
545 
547 {
548  int nIter = 0;
549  double _chisq_kalman_prev = 1E6;
550  for(; nIter < _max_iteration; ++nIter)
551  {
552 #ifdef _DEBUG_ON
553  LogInfo("Iteration: " << nIter);
554 #endif
555 
556  _chisq_vertex = 0.;
557  _chisq_kalman = 0.;
558  for(unsigned int j = 0; j < _trkpar_curr.size(); j++)
559  {
560  _node_vertex.resetFlags();
561  //_node_vertex.getMeasurement() = _vtxpar_curr._r.GetSub(0, 1, 0, 0);
562  //_node_vertex.getMeasurementCov() = _vtxpar_curr._cov.GetSub(0, 1, 0, 1);
563  _node_vertex.setZ(_vtxpar_curr._r[2][0]);
564 
565  _kmfit->setCurrTrkpar(_trkpar_curr[j]);
566  if(!_kmfit->fit_node(_node_vertex))
567  {
568 #ifdef _DEBUG_ON
569  LogInfo("Vertex fit for this track failed!");
570 #endif
571  _chisq_kalman = 1E5;
572  break;
573  }
574  else
575  {
576  _chisq_kalman += _node_vertex.getChisq();
577  }
578 
579  updateVertex();
580  }
581 
583 #ifdef _DEBUG_ON
584  LogInfo("At this iteration: ");
585  LogInfo(_vtxpar_curr._r[2][0] << " === " << _node_vertex.getZ() << " : " << _chisq_kalman << " === " << _chisq_vertex << " : " << _vtxpar_curr._r[0][0] << " === " << _vtxpar_curr._r[1][0]);
586 #endif
587  if(_vtxpar_curr._r[2][0] < Z_UPSTREAM || _vtxpar_curr._r[2][0] > Z_DOWNSTREAM)
588  {
589  _chisq_kalman = 1E5;
590  _chisq_vertex = 1E5;
591  break;
592  }
593 
594  if(nIter > 0 && (fabs(_vtxpar_curr._r[2][0] - _node_vertex.getZ()) < _tolerance || _chisq_kalman > _chisq_kalman_prev)) break;
595  _chisq_kalman_prev = _chisq_kalman;
596  }
597 
598  return nIter+1;
599 }
600 
602 {
603  double p = fabs(1./_node_vertex.getFiltered()._state_kf[0][0]);
604  double px = _node_vertex.getFiltered()._state_kf[1][0];
605  double py = _node_vertex.getFiltered()._state_kf[2][0];
606  double pz = sqrt(p*p - px*px - py*py);
607 
609  TMatrixD H(2, 3);
610  H.Zero();
611 
612  H[0][0] = 1.;
613  H[1][1] = 1.;
614  H[0][2] = -px/pz;
615  H[1][2] = -py/pz;
616 
617  TMatrixD vertex_dummy(3, 1);
618  vertex_dummy.Zero();
619  vertex_dummy[2][0] = _vtxpar_curr._r[2][0];
620 
621  TMatrixD mxy = _node_vertex.getFiltered()._state_kf.GetSub(3, 4, 0, 0);
622  TMatrixD Vxy = _node_vertex.getFiltered()._covar_kf.GetSub(3, 4, 3, 4);
623  TMatrixD S = SMatrix::invertMatrix(Vxy + SMatrix::getABCt(H, _vtxpar_curr._cov, H));
624  TMatrixD K = SMatrix::getABtC(_vtxpar_curr._cov, H, S);
625  TMatrixD zeta = mxy - H*(_vtxpar_curr._r - vertex_dummy);
626  TMatrixD _r_filtered = _vtxpar_curr._r + K*zeta;
627  TMatrixD _cov_filtered = _vtxpar_curr._cov - K*H*(_vtxpar_curr._cov);
628 
629  _chisq_vertex += SMatrix::getAtBC(zeta, S, zeta)[0][0];
630 
631  _vtxpar_curr._r = _r_filtered;
632  _vtxpar_curr._cov = _cov_filtered;
633 }
634 
636 {
637  TrkPar _trkpar_start;
638  _trkpar_start._state_kf = _track.getStateVector(0);
639  _trkpar_start._covar_kf = _track.getCovariance(0);
640  _trkpar_start._z = _track.getZ(0);
641 
642  return findSingleMuonVertex(_trkpar_start);
643 }
644 
646 {
647  TrkPar _trkpar_start;
648  if(_node_start.isSmoothDone())
649  {
650  _trkpar_start = _node_start.getSmoothed();
651  }
652  else if(_node_start.isFilterDone())
653  {
654  _trkpar_start = _node_start.getFiltered();
655  }
656  else
657  {
658  _trkpar_start = _node_start.getPredicted();
659  }
660 
661  return findSingleMuonVertex(_trkpar_start);
662 }
663 
665 {
666  _extrapolator.setInitialStateWithCov(_trkpar_start._z, _trkpar_start._state_kf, _trkpar_start._covar_kf);
667  double z_vertex_single = _extrapolator.extrapolateToIP();
668 
669  return z_vertex_single;
670 }
671 
672 void VertexFit::bookEvaluation(std::string evalFileName)
673 {
674  evalFile = new TFile(evalFileName.c_str(), "recreate");
675  evalTree = new TTree("T", "VertexFit eval");
676 
677  evalTree->Branch("runID", &runID, "runID/I");
678  evalTree->Branch("eventID", &eventID, "eventID/I");
679  evalTree->Branch("targetPos", &targetPos, "targetPos/I");
680  evalTree->Branch("choice_eval", &choice_eval, "choice_eval/I");
681  evalTree->Branch("choice_by_kf_eval", &choice_by_kf_eval, "choice_by_kf_eval/I");
682  evalTree->Branch("choice_by_vx_eval", &choice_by_vx_eval, "choice_by_vx_eval/I");
683 
684  evalTree->Branch("nStart", &nStart, "nStart/I");
685  evalTree->Branch("nIter_eval", nIter_eval, "nIter_eval[nStart]/I");
686  evalTree->Branch("chisq_kf_eval", chisq_kf_eval, "chisq_kf_eval[nStart]/D");
687  evalTree->Branch("chisq_vx_eval", chisq_vx_eval, "chisq_vx_eval[nStart]/D");
688  evalTree->Branch("z_vertex_eval", z_vertex_eval, "z_vertex_eval[nStart]/D");
689  evalTree->Branch("r_vertex_eval", r_vertex_eval, "r_vertex_eval[nStart]/D");
690  evalTree->Branch("z_start_eval", z_start_eval, "z_start_eval[nStart]/D");
691 
692  evalTree->Branch("m_chisq_kf_eval", &m_chisq_kf_eval, "m_chisq_kf_eval/D");
693  evalTree->Branch("m_z_vertex_eval", &m_z_vertex_eval, "m_z_vertex_eval/D");
694  evalTree->Branch("s_chisq_kf_eval", &s_chisq_kf_eval, "s_chisq_kf_eval/D");
695  evalTree->Branch("s_z_vertex_eval", &s_z_vertex_eval, "s_z_vertex_eval/D");
696 }
697 
699 {
700  if(evalTree == nullptr) return;
701 
702  choice_by_kf_eval = -1;
703  choice_by_vx_eval = -1;
704  double chisq_kf_min = 1E6;
705  double chisq_vx_min = 1E6;
706 
707  m_chisq_kf_eval = 0.;
708  m_z_vertex_eval = 0.;
709  for(int i = 0; i < nStart; ++i)
710  {
711  chisq_kf_eval[i] = chisq_km[i];
712  chisq_vx_eval[i] = chisq_vx[i];
713  z_vertex_eval[i] = z_vertex[i];
714  r_vertex_eval[i] = r_vertex[i];
715  z_start_eval[i] = z_start[i];
716 
717  m_chisq_kf_eval += chisq_kf_eval[i];
718  m_z_vertex_eval += z_vertex_eval[i];
719 
720  if(chisq_kf_min > chisq_km[i])
721  {
722  choice_by_kf_eval = i;
723  chisq_kf_min = chisq_km[i];
724  }
725 
726  if(chisq_vx_min > chisq_vx[i])
727  {
728  choice_by_vx_eval = i;
729  chisq_vx_min = chisq_vx[i];
730  }
731  }
732  m_chisq_kf_eval /= nStart;
733  m_z_vertex_eval /= nStart;
734 
735  s_chisq_kf_eval = 0.;
736  s_z_vertex_eval = 0.;
737  if(nStart == 1)
738  {
739  evalTree->Fill();
740  return;
741  }
742 
743  for(int i = 0; i < nStart; ++i)
744  {
745  s_chisq_kf_eval += ((m_chisq_kf_eval - chisq_kf_eval[i])*(m_chisq_kf_eval - chisq_kf_eval[i]));
746  s_z_vertex_eval += ((m_z_vertex_eval - z_vertex_eval[i])*(m_z_vertex_eval - z_vertex_eval[i]));
747  }
748  s_chisq_kf_eval = sqrt(s_chisq_kf_eval/(nStart-1));
749  s_z_vertex_eval = sqrt(s_z_vertex_eval/(nStart-1));
750 
751  evalTree->Fill();
752 }
753 
755 {
756  using namespace std;
757 
758  for(unsigned int i = 0; i < z_start.size(); i++)
759  {
760  cout << "============= Hypothesis " << i << " ============" << endl;
761  cout << "z_start = " << z_start[i] << ", sigz_start = " << sig_z_start[i] << endl;
762  cout << "Found z_vertex = " << z_vertex[i] << endl;
763  cout << "With chisq_km = " << chisq_km[i] << " and chisq_vx = " << chisq_vx[i] << endl;
764  }
765 }
#define LogInfo(message)
Definition: DbSvc.cc:15
#define VFEXIT_FAIL_DIMUONPAIR
Definition: GlobalConsts.h:26
#define VFEXIT_SUCCESS
Definition: GlobalConsts.h:18
#define VFEXIT_FAIL_ITERATION
Definition: GlobalConsts.h:27
const double FMAGSTR
Definition: MacroCommon.h:2
const double KMAGSTR
Definition: MacroCommon.h:3
TFile clean handling.
int verbosity
The verbosity level. 0 means not verbose at all.
Definition: Fun4AllBase.h:75
@ VERBOSITY_A_LOT
Output a lot of messages.
Definition: Fun4AllBase.h:48
virtual int Verbosity() const
Gets the verbosity of this module.
Definition: Fun4AllBase.h:64
void setPropCalc(bool option)
External control of modes.
bool init(const PHField *field, const TGeoManager *geom)
Initialize geometry and physics.
double extrapolateToIP()
Extrapolate to the primary vertex.
void setLengthCalc(bool option)
void setInitialStateWithCov(double z_in, TMatrixD &state_in, TMatrixD &cov_in)
Set input initial state parameters.
static KalmanFilter * instance()
singlton instance
void setCurrTrkpar(Node &_node)
set the current track parameter using the current node
Definition: KalmanFilter.h:45
bool fit_node(Node &_node)
Fit one node.
void enableDumpCorrection()
Enable the dump mode: stop calc prop matrix, start calc travel length.
Definition: KalmanFilter.h:52
std::list< Node > & getNodeList()
Definition: KalmanTrack.h:84
bool isFilterDone()
Definition: KalmanUtil.h:132
TrkPar & getSmoothed()
Definition: KalmanUtil.h:108
void setZ(double z)
Definition: KalmanUtil.h:137
double getChisq()
Definition: KalmanUtil.h:127
TrkPar & getPredicted()
Gets.
Definition: KalmanUtil.h:106
TMatrixD & getMeasurement()
Definition: KalmanUtil.h:110
bool isSmoothDone()
Definition: KalmanUtil.h:133
void resetFlags()
Definition: KalmanUtil.cxx:381
double getZ()
Definition: KalmanUtil.h:126
TMatrixD & getProjector()
Definition: KalmanUtil.h:122
TrkPar & getFiltered()
Definition: KalmanUtil.h:107
TMatrixD & getMeasurementCov()
Definition: KalmanUtil.h:111
PHFieldConfig_v3 implements field configuration information for input a field map file.
static PHField * GetFieldMapNode(const PHFieldConfig *default_config=nullptr, PHCompositeNode *topNode=nullptr, const int verbosity=0)
Get transient PHField from DST nodes. If not found, make a new one based on default_config.
transient DST object for field storage and access
Definition: PHField.h:14
virtual void identify(std::ostream &os=std::cout) const
Definition: PHField.h:30
virtual double get_DoubleFlag(const std::string &name) const
Definition: PHFlag.cc:49
virtual int get_IntFlag(const std::string &name) const
Definition: PHFlag.cc:117
static TGeoManager * GetTGeoManager(PHCompositeNode *topNode)
Main user interface: DST node -> TGeoManager for downstream use.
static TMatrixD getABtC(const TMatrixD &A, const TMatrixD &B, const TMatrixD &C)
Definition: KalmanUtil.cxx:141
static TMatrixD getABCt(const TMatrixD &A, const TMatrixD &B, const TMatrixD &C)
Definition: KalmanUtil.cxx:125
static TMatrixD getAtBC(const TMatrixD &A, const TMatrixD &B, const TMatrixD &C)
Definition: KalmanUtil.cxx:133
static TMatrixD invertMatrix(const TMatrixD &m)
Definition: KalmanUtil.cxx:52
Double_t chisq_single
Definition: SRecEvent.h:398
TVector3 vtx
3-vector vertex position
Definition: SRecEvent.h:379
void calcVariables(int choice=0)
Calculate the kinematic vairables, 0 = vertex, 1 = target, 2 = dump.
Definition: SRecEvent.cxx:633
Double_t chisq_dump
Definition: SRecEvent.h:406
Double_t chisq_upstream
Definition: SRecEvent.h:407
TVector3 proj_target_neg
Definition: SRecEvent.h:386
Double_t chisq_vx
Definition: SRecEvent.h:402
Int_t trackID_pos
Index of muon track used in host SRecEvent.
Definition: SRecEvent.h:361
Double_t chisq_target
Chisq of three test position.
Definition: SRecEvent.h:405
TVector3 proj_dump_pos
Definition: SRecEvent.h:385
TLorentzVector p_neg_single
Definition: SRecEvent.h:376
TVector3 proj_dump_neg
Definition: SRecEvent.h:387
Int_t trackID_neg
Definition: SRecEvent.h:362
TLorentzVector p_neg
Definition: SRecEvent.h:366
TLorentzVector p_pos_single
4-momentum of the muon tracks before vertex fit
Definition: SRecEvent.h:375
Double_t chisq_kf
Vertex fit chisqs.
Definition: SRecEvent.h:401
TLorentzVector p_pos
4-momentum of the muon tracks after vertex fit
Definition: SRecEvent.h:365
TVector3 vtx_neg
Definition: SRecEvent.h:381
TVector3 vtx_pos
Definition: SRecEvent.h:380
TVector3 proj_target_pos
Track projections at different location.
Definition: SRecEvent.h:384
void setRecStatus(int status)
Definition: SRecEvent.h:432
Int_t getEventID()
Definition: SRecEvent.h:442
std::vector< Int_t > getChargedTrackIDs(Int_t charge)
Get track IDs.
Definition: SRecEvent.cxx:785
void identify(std::ostream &os=std::cout) const
PHObject virtual overloads.
Definition: SRecEvent.h:418
Int_t getRunID()
Definition: SRecEvent.h:440
Int_t getTargetPos()
Definition: SRecEvent.h:443
void insertDimuon(SRecDimuon dimuon)
Insert dimuon.
Definition: SRecEvent.h:470
void clearDimuons()
Clear the dimuon list.
Definition: SRecEvent.cxx:824
SRecTrack & getTrack(Int_t i)
Definition: SRecEvent.h:456
Int_t getNDimuons()
Get dimuons.
Definition: SRecEvent.h:462
TMatrixD getCovariance(Int_t i)
Definition: SRecEvent.h:110
TVector3 getTargetPos()
Definition: SRecEvent.h:188
Int_t getCharge() const
Gets.
Definition: SRecEvent.h:101
Double_t getChisqTarget()
Definition: SRecEvent.h:199
Double_t getChisqDump()
Definition: SRecEvent.h:198
TVector3 getVertex()
Definition: SRecEvent.h:181
TLorentzVector getMomentumVertex()
Get the vertex info.
Definition: SRecEvent.cxx:340
void setZVertex(Double_t z, bool update=true)
Definition: SRecEvent.cxx:168
TMatrixD getStateVector(Int_t i)
Definition: SRecEvent.h:109
Double_t getZ(Int_t i)
Definition: SRecEvent.h:111
Double_t getChisqUpstream()
Definition: SRecEvent.h:200
int isValid() const
isValid returns non zero if object contains vailid data
Definition: SRecEvent.cxx:371
void swimToVertex(TVector3 *pos=nullptr, TVector3 *mom=nullptr, bool hyptest=true)
Simple swim to vertex.
Definition: SRecEvent.cxx:402
TVector3 getDumpPos()
Get mom/pos at a given location.
Definition: SRecEvent.h:186
Double_t getChisqVertex()
Definition: SRecEvent.h:183
TMatrixD _state_kf
State vectors and its covariance.
Definition: KalmanUtil.h:85
TMatrixD _covar_kf
Definition: KalmanUtil.h:86
double _z
Definition: KalmanUtil.h:87
int process_event(PHCompositeNode *topNode)
Definition: VertexFit.cxx:185
int setRecEvent(SRecEvent *recEvent, int sign1=1, int sign2=-1)
Set the SRecEvent, main external call the use vertex fit.
Definition: VertexFit.cxx:259
void fillEvaluation()
Definition: VertexFit.cxx:698
void addHypothesis(double z, double sigz=50.)
Definition: VertexFit.h:90
double findDimuonVertexFast(SRecTrack &track1, SRecTrack &track2)
Definition: VertexFit.cxx:413
double findSingleMuonVertex(SRecTrack &_track)
Definition: VertexFit.cxx:635
int findVertex()
Find the primary vertex.
Definition: VertexFit.cxx:546
void print()
Debugging output.
Definition: VertexFit.cxx:754
void bookEvaluation(std::string evalFileName="vtx_eval.root")
Evaluation.
Definition: VertexFit.cxx:672
void addTrack(int index, SRecTrack &_track)
Add one track parameter set into the fit.
Definition: VertexFit.cxx:531
void updateVertex()
Core function, update the vertex prediction according to the track info.
Definition: VertexFit.cxx:601
double getVXChisq()
Definition: VertexFit.h:110
void init()
Initialize and reset.
Definition: VertexFit.cxx:441
int Init(PHCompositeNode *topNode)
Definition: VertexFit.cxx:140
int InitRun(PHCompositeNode *topNode)
Definition: VertexFit.cxx:147
int End(PHCompositeNode *topNode)
Called at the end of all processing.
Definition: VertexFit.cxx:205
void setStartingVertex(double z_start, double sigz_start)
Definition: VertexFit.cxx:458
VertexFit(const std::string &name="VertexFit")
Definition: VertexFit.cxx:89
double getVertexZ0()
Gets.
Definition: VertexFit.h:109
int processOnePair()
After setting both tracks and hypothesis, start the iteration.
Definition: VertexFit.cxx:474
TMatrixD _r
Definition: VertexFit.h:58
TMatrixD _cov
Definition: VertexFit.h:59
static recoConsts * instance()
Definition: recoConsts.cc:7