Class Reference for E1039 Core & Analysis Software
SRecEvent.cxx
Go to the documentation of this file.
1 /*
2 SRecEvent.cxx
3 
4 Implimentation of the calss SRecTrack and SRecEvent
5 
6 Author: Kun Liu, liuk@fnal.gov
7 Created: 01-21-2013
8 */
9 
10 #include <cmath>
11 
12 #include <TLorentzVector.h>
13 #include <TMatrixD.h>
14 
15 #include <GenFit/SharedPlanePtr.h>
16 #include <phool/recoConsts.h>
17 
18 #include "SRecEvent.h"
19 #include "KalmanUtil.h"
20 #include "KalmanFilter.h"
21 
25 
26 namespace
27 {
28  //static flag to indicate the initialized has been done
29  static bool inited = false;
30 
31  //static flag of kmag on/off
32  static bool KMAG_ON;
33 
34  //static flag of kmag strength
35  static double FMAGSTR;
36  static double KMAGSTR;
37 
38  static double PT_KICK_FMAG;
39  static double PT_KICK_KMAG;
40 
41  //Geometric positions
42  static double FMAG_LENGTH;
43  static double Z_TARGET;
44  static double Z_DUMP;
45  static double Z_UPSTREAM;
46  static double Z_DOWNSTREAM;
47  static double FMAG_HOLE_LENGTH;
48  static double FMAG_HOLE_RADIUS;
49 
50  //Beam position and shape
51  static double X_BEAM;
52  static double Y_BEAM;
53  static double SIGX_BEAM;
54  static double SIGY_BEAM;
55 
56  //Simple swimming settings
57  static int NSTEPS_TARGET;
58  static int NSTEPS_SHIELDING;
59  static int NSTEPS_FMAG;
60 
61  static double DEDX_UNIT_0;
62  static double DEDX_UNIT_1;
63  static double DEDX_UNIT_2;
64  static double DEDX_UNIT_3;
65  static double DEDX_UNIT_4;
66  static double PTKICK_UNIT;
67 
68  static double STEP_TARGET;
69  static double STEP_SHIELDING;
70  static double STEP_FMAG;
71 
72  //initialize global variables
73  void initGlobalVariables()
74  {
75  if(!inited)
76  {
77  inited = true;
78 
80  KMAG_ON = rc->get_BoolFlag("KMAG_ON");
81  FMAGSTR = rc->get_DoubleFlag("FMAGSTR");
82  KMAGSTR = rc->get_DoubleFlag("KMAGSTR");
83  PT_KICK_FMAG = rc->get_DoubleFlag("PT_KICK_FMAG")*FMAGSTR;
84  PT_KICK_KMAG = rc->get_DoubleFlag("PT_KICK_KMAG")*KMAGSTR;
85 
86  X_BEAM = rc->get_DoubleFlag("X_BEAM");
87  Y_BEAM = rc->get_DoubleFlag("Y_BEAM");
88  SIGX_BEAM = rc->get_DoubleFlag("SIGX_BEAM");
89  SIGY_BEAM = rc->get_DoubleFlag("SIGY_BEAM");
90 
91  FMAG_LENGTH = rc->get_DoubleFlag("FMAG_LENGTH");
92  Z_TARGET = rc->get_DoubleFlag("Z_TARGET");
93  Z_DUMP = rc->get_DoubleFlag("Z_DUMP");
94  Z_UPSTREAM = rc->get_DoubleFlag("Z_UPSTREAM");
95  Z_DOWNSTREAM = rc->get_DoubleFlag("Z_DOWNSTREAM");
96  FMAG_HOLE_RADIUS = rc->get_DoubleFlag("FMAG_HOLE_RADIUS");
97  FMAG_HOLE_LENGTH = rc->get_DoubleFlag("FMAG_HOLE_LENGTH");
98 
99  NSTEPS_TARGET = rc->get_IntFlag("NSTEPS_TARGET");
100  NSTEPS_SHIELDING = rc->get_IntFlag("NSTEPS_SHIELDING");
101  NSTEPS_FMAG = rc->get_IntFlag("NSTEPS_FMAG");
102 
103  DEDX_UNIT_0 = rc->get_DoubleFlag("DEDX_FE_P0")/FMAG_LENGTH;
104  DEDX_UNIT_1 = rc->get_DoubleFlag("DEDX_FE_P1")/FMAG_LENGTH;
105  DEDX_UNIT_2 = rc->get_DoubleFlag("DEDX_FE_P2")/FMAG_LENGTH;
106  DEDX_UNIT_3 = rc->get_DoubleFlag("DEDX_FE_P3")/FMAG_LENGTH;
107  DEDX_UNIT_4 = rc->get_DoubleFlag("DEDX_FE_P4")/FMAG_LENGTH;
108  PTKICK_UNIT = PT_KICK_FMAG/FMAG_LENGTH;
109 
110  STEP_TARGET = fabs(Z_UPSTREAM)/NSTEPS_TARGET;
111  STEP_SHIELDING = 0.;
112  STEP_FMAG = FMAG_LENGTH/NSTEPS_FMAG/2.;
113  }
114  }
115 }
116 
118 {
119  fChisq = -99.;
120 
121  fHitIndex.clear();
122  fState.clear();
123  fCovar.clear();
124  fZ.clear();
125  fChisqAtNode.clear();
126 
127  for(Int_t i = 0; i < 3; ++i) fGFDetPlaneVec[i].clear();
128  fGFAuxInfo.clear();
129  fGFStateVec.clear();
130  fGFCov.clear();
131 
132  fChisqVertex = -99.;
133  fVertexPos.SetXYZ(999., 999., 999.);
134  fStateVertex.ResizeTo(5, 1);
135  fCovarVertex.ResizeTo(5, 5);
136 
137  fChisqTarget = -99.;
138  fChisqDump = -99.;
139  fChisqUpstream = -99.;
140 
141  initGlobalVariables();
142 }
143 
144 bool SRecTrack::operator<(const SRecTrack& elem) const
145 {
146  return getNHits() == elem.getNHits() ? (fChisq < elem.fChisq) : (getProb() > elem.getProb());
147 }
148 
149 Int_t SRecTrack::getNHitsInStation(Int_t stationID)
150 {
151  if(stationID != 1 && stationID != 2 && stationID != 3) return 0;
152 
153  Double_t z_ref[4] = {0., 700., 1600., 2500.};
154  Int_t nHits = 0;
155  for(std::vector<Double_t>::iterator iter = fZ.begin(); iter != fZ.end(); ++iter)
156  {
157  if(*iter > z_ref[stationID-1] && *iter < z_ref[stationID]) ++nHits;
158  }
159 
160  return nHits;
161 }
162 
163 Double_t SRecTrack::getProb() const
164 {
165  return KMAG_ON == 1 ? TMath::Prob(fChisq, getNHits() - 5) : TMath::Prob(fChisq, getNHits() - 4);
166 }
167 
168 void SRecTrack::setZVertex(Double_t z, bool update)
169 {
170  Node _node_vertex;
171  _node_vertex.setZ(z);
172 
173  TMatrixD m(2, 1), cov(2, 2), proj(2, 5);
174  m[0][0] = X_BEAM;
175  m[1][0] = Y_BEAM;
176 
177  cov.Zero();
178  cov[0][0] = SIGX_BEAM*SIGX_BEAM;
179  cov[1][1] = SIGY_BEAM*SIGY_BEAM;
180 
181  proj.Zero();
182  proj[0][3] = 1.;
183  proj[1][4] = 1.;
184 
185  _node_vertex.getMeasurement().ResizeTo(2, 1);
186  _node_vertex.getMeasurementCov().ResizeTo(2, 2);
187  _node_vertex.getProjector().ResizeTo(2, 5);
188 
189  _node_vertex.getMeasurement() = m;
190  _node_vertex.getMeasurementCov() = cov;
191  _node_vertex.getProjector() = proj;
192 
193  TrkPar _trkpar_curr;
194  _trkpar_curr._state_kf = fState[0];
195  _trkpar_curr._covar_kf = fCovar[0];
196  _trkpar_curr._z = fZ[0];
197 
199  kmfit->enableDumpCorrection();
200  kmfit->setCurrTrkpar(_trkpar_curr);
201  kmfit->fit_node(_node_vertex);
202 
203  fChisqVertex = _node_vertex.getChisq();
204  if(!update) return;
205 
206  fVertexPos.SetXYZ(_node_vertex.getFiltered().get_x(), _node_vertex.getFiltered().get_y(), z);
207  fVertexMom = _node_vertex.getFiltered().get_mom_vec();
208 
209  fStateVertex = _node_vertex.getFiltered()._state_kf;
210  fCovarVertex = _node_vertex.getFiltered()._covar_kf;
211 }
212 
214 {
215  setZVertex(Z_TARGET, false);
216  setChisqTarget(fChisqVertex);
217 
218  setZVertex(Z_DUMP, false);
219  setChisqDump(fChisqVertex);
220 
221  setZVertex(Z_UPSTREAM+10., false);
222  setChisqUpstream(fChisqVertex);
223 
224  setZVertex(getZVertex(), true);
225 }
226 
227 void SRecTrack::setVertexFast(TVector3 mom, TVector3 pos)
228 {
229  fVertexPos = pos;
230  fVertexMom = mom;
231 
232  fStateVertex[0][0] = getCharge()/mom.Mag();
233  fStateVertex[1][0] = mom[0]/mom[2];
234  fStateVertex[2][0] = mom[1]/mom[2];
235  fStateVertex[3][0] = pos[0];
236  fStateVertex[4][0] = pos[1];
237 
238  fCovarVertex.UnitMatrix();
239 }
240 
242 {
243  if(fChisqVertex > 50.) return false;
244  if(fVertexPos.Z() < Z_UPSTREAM || fVertexPos.Z() > Z_DOWNSTREAM) return false;
245 
246  return true;
247 }
248 
249 Int_t SRecTrack::getNearestNode(Double_t z)
250 {
251  Int_t nNodes = getNHits();
252  Int_t index_min = 0;
253  Double_t dist_min = 1E6;
254  for(Int_t i = 0; i < nNodes; ++i)
255  {
256  Double_t dist = fabs(z - fZ[i]);
257  if(dist < dist_min)
258  {
259  dist_min = dist;
260  index_min = i;
261  }
262  }
263 
264  return index_min;
265 }
266 
267 void SRecTrack::getExpPositionFast(Double_t z, Double_t& x, Double_t& y, Int_t iNode)
268 {
269  if(iNode < 0 || iNode >= getNHits())
270  {
271  iNode = getNearestNode(z);
272  }
273 
274  TrkPar _trkpar;
275  _trkpar._state_kf = fState[iNode];
276  _trkpar._covar_kf = fCovar[iNode];
277 
278  Double_t z_ref = fZ[iNode];
279  Double_t x_ref = _trkpar.get_x();
280  Double_t y_ref = _trkpar.get_y();
281  Double_t axz = _trkpar.get_dxdz();
282  Double_t ayz = _trkpar.get_dydz();
283 
284  x = x_ref + axz*(z - z_ref);
285  y = y_ref + ayz*(z - z_ref);
286 }
287 
288 void SRecTrack::getExpPosErrorFast(Double_t z, Double_t& dx, Double_t& dy, Int_t iNode)
289 {
290  if(iNode < 0 || iNode >= getNHits())
291  {
292  iNode = getNearestNode(z);
293  }
294 
295  Double_t z_ref = fZ[iNode];
296  Double_t dx_ref = sqrt((fCovar[iNode])[3][3]);
297  Double_t dy_ref = sqrt((fCovar[iNode])[4][4]);
298  Double_t daxz = sqrt((fCovar[iNode])[1][1]);
299  Double_t dayz = sqrt((fCovar[iNode])[2][2]);
300 
301  dx = 2.*(dx_ref + fabs(daxz*(z - z_ref)));
302  dy = 2.*(dy_ref + fabs(dayz*(z - z_ref)));
303 }
304 
305 double SRecTrack::getExpMomentumFast(Double_t z, Int_t iNode)
306 {
307  Double_t px, py, pz;
308  return getExpMomentumFast(z, px, py, pz, iNode);
309 }
310 
311 Double_t SRecTrack::getExpMomentumFast(Double_t z, Double_t& px, Double_t& py, Double_t& pz, Int_t iNode)
312 {
313  if(iNode < 0 || iNode >= getNHits())
314  {
315  iNode = getNearestNode(z);
316  }
317 
318  return getMomentum(fState[iNode], px, py, pz);
319 
320 }
321 
322 Double_t SRecTrack::getMomentum(const TMatrixD& state, Double_t& px, Double_t& py, Double_t& pz) const
323 {
324  Double_t p = 1./fabs(state[0][0]);
325  pz = p/sqrt(1. + state[1][0]*state[1][0] + state[2][0]*state[2][0]);
326  px = pz*state[1][0];
327  py = pz*state[2][0];
328 
329  return p;
330 }
331 
332 Double_t SRecTrack::getPosition(const TMatrixD& state, Double_t& x, Double_t& y) const
333 {
334  x = state[3][0];
335  y = state[4][0];
336 
337  return sqrt(x*x + y*y);
338 }
339 
341 {
342  Double_t mmu = 0.10566;
343  Double_t px, py, pz, E;
344 
345  getMomentumVertex(px, py, pz);
346  E = sqrt(px*px + py*py + pz*pz + mmu*mmu);
347 
348  return TLorentzVector(px, py, pz, E);
349 }
350 
351 void SRecTrack::insertGFState(const genfit::MeasuredStateOnPlane& msop)
352 {
353  fGFStateVec.push_back(msop.getState());
354  fGFCov.push_back(msop.getCov());
355  fGFAuxInfo.push_back(msop.getAuxInfo());
356 
357  const genfit::SharedPlanePtr& plane = msop.getPlane();
358  fGFDetPlaneVec[0].push_back(plane->getO());
359  fGFDetPlaneVec[1].push_back(plane->getU());
360  fGFDetPlaneVec[2].push_back(plane->getV());
361 }
362 
363 void SRecTrack::adjustKMag(double kmagStr)
364 {
365  for(std::vector<TMatrixD>::iterator iter = fState.begin(); iter != fState.end(); ++iter)
366  {
367  (*iter)[0][0] = (*iter)[0][0]/kmagStr;
368  }
369 }
370 
372 {
373  //Vertex valid
374  if(!isVertexValid()) return false;
375 
376  //Number of hits cut
377  Int_t nHits = getNHits();
378  if(nHits < 14) return false;
379 
380  //Trigger road cut
381  //if(fTriggerID == 0) return false;
382 
383  //Total chisq, may change to cut on prob
384  if(getChisq()/(nHits - 5) > 15.) return false;
385 
386  //Check the px polarity
387  //if(FMAGSTR*getCharge()*fVertexMom.Px() < 0) return false;
388 
389  return true;
390 }
391 
393 {
394  return (fVertexPos.Z() > -300 && fVertexPos.Z() < 0. && fChisqDump - fChisqTarget > 10.);
395 }
396 
398 {
399  return (fVertexPos.Z() > 0. && fVertexPos.Z() < 150. && fChisqTarget - fChisqDump > 10.);
400 }
401 
402 void SRecTrack::swimToVertex(TVector3* pos, TVector3* mom, bool hyptest)
403 {
404  //Store the steps on each point (center of the interval)
405  bool cleanupPos = false;
406  bool cleanupMom = false;
407  if(pos == nullptr)
408  {
409  pos = new TVector3[NSTEPS_FMAG + NSTEPS_TARGET + 1];
410  cleanupPos = true;
411  }
412  if(mom == nullptr)
413  {
414  mom = new TVector3[NSTEPS_FMAG + NSTEPS_TARGET + 1];
415  cleanupMom = true;
416  }
417 
418  //track slope/location in upstream
419  double tx = fState.front()[1][0];
420  double ty = fState.front()[2][0];
421  double x0 = fState.front()[3][0];
422  double y0 = fState.front()[4][0];
423  double z0 = fZ.front();
424 
425  //Initial position should be on the downstream face of beam dump
426  pos[0].SetXYZ(x0 + tx*(FMAG_LENGTH - z0), y0 + ty*(FMAG_LENGTH - z0), FMAG_LENGTH);
427  mom[0] = getMomentumVecSt1();
428 
429  //Charge of the track
430  double charge = getCharge();
431 
432  //Now make the swim
433  int iStep = 1;
434  for(; iStep <= NSTEPS_FMAG; ++iStep)
435  {
436  //Make pT kick at the center of slice, add energy loss at both first and last half-slice
437  //Note that ty is the global class data member, which does not change during the entire swimming
438  double tx_i = mom[iStep-1].Px()/mom[iStep-1].Pz();
439  double tx_f = tx_i + 2.*charge*PTKICK_UNIT*STEP_FMAG/sqrt(mom[iStep-1].Px()*mom[iStep-1].Px() + mom[iStep-1].Pz()*mom[iStep-1].Pz());
440 
441  TVector3 trajVec1(tx_i*STEP_FMAG, ty*STEP_FMAG, STEP_FMAG);
442  TVector3 pos_b = pos[iStep-1] - trajVec1;
443 
444  double p_tot_i = mom[iStep-1].Mag();
445  double p_tot_b;
446  if(pos_b[2] > FMAG_HOLE_LENGTH || pos_b.Perp() > FMAG_HOLE_RADIUS)
447  {
448  p_tot_b = p_tot_i + (DEDX_UNIT_0 + p_tot_i*DEDX_UNIT_1 + p_tot_i*p_tot_i*DEDX_UNIT_2 + p_tot_i*p_tot_i*p_tot_i*DEDX_UNIT_3 + p_tot_i*p_tot_i*p_tot_i*p_tot_i*DEDX_UNIT_4)*trajVec1.Mag();
449  }
450  else
451  {
452  p_tot_b = p_tot_i;
453  }
454 
455  TVector3 trajVec2(tx_f*STEP_FMAG, ty*STEP_FMAG, STEP_FMAG);
456  pos[iStep] = pos_b - trajVec2;
457 
458  double p_tot_f;
459  if(pos[iStep][2] > FMAG_HOLE_LENGTH || pos[iStep].Perp() > FMAG_HOLE_RADIUS)
460  {
461  p_tot_f = p_tot_b + (DEDX_UNIT_0 + p_tot_b*DEDX_UNIT_1 + p_tot_b*p_tot_b*DEDX_UNIT_2 + p_tot_b*p_tot_b*p_tot_b*DEDX_UNIT_3 + p_tot_b*p_tot_b*p_tot_b*p_tot_b*DEDX_UNIT_4)*trajVec2.Mag();
462  }
463  else
464  {
465  p_tot_f = p_tot_b;
466  }
467 
468  //Now the final position and momentum in this step
469  double pz_f = p_tot_f/sqrt(1. + tx_f*tx_f + ty*ty);
470  mom[iStep].SetXYZ(pz_f*tx_f, pz_f*ty, pz_f);
471  pos[iStep] = pos[iStep-1] - trajVec1 - trajVec2;
472 
473  //Save the dump position when applicable
474  if(fabs(pos_b.Z() - Z_DUMP) < STEP_FMAG)
475  {
476  double dz = Z_DUMP - pos_b.Z();
477  if(dz < 0)
478  {
479  setDumpPos(pos_b + TVector3(tx_f*dz, ty*dz, dz));
480  setDumpMom(mom[iStep]);
481  }
482  else
483  {
484  setDumpPos(pos_b + TVector3(tx_i*dz, ty*dz, dz));
485  setDumpMom(mom[iStep-1]);
486  }
487  }
488 
489 #ifdef _DEBUG_ON_LEVEL_2
490  std::cout << "FMAG: " << iStep << ": " << pos[iStep-1][2] << " ==================>>> " << pos[iStep][2] << std::endl;
491  std::cout << mom[iStep-1][0]/mom[iStep-1][2] << " " << mom[iStep-1][1]/mom[iStep-1][2] << " " << mom[iStep-1][2] << " ";
492  std::cout << pos[iStep-1][0] << " " << pos[iStep-1][1] << " " << pos[iStep-1][2] << std::endl << std::endl;
493  std::cout << mom[iStep][0]/mom[iStep][2] << " " << mom[iStep][1]/mom[iStep][2] << " " << mom[iStep][2] << " ";
494  std::cout << pos[iStep][0] << " " << pos[iStep][1] << " " << pos[iStep][2] << std::endl << std::endl;
495 #endif
496  }
497 
498  for(; iStep < NSTEPS_FMAG+NSTEPS_TARGET+1; ++iStep)
499  {
500  //Simple straight line flight
501  double tx_i = mom[iStep-1].Px()/mom[iStep-1].Pz();
502  TVector3 trajVec(tx_i*STEP_TARGET, ty*STEP_TARGET, STEP_TARGET);
503 
504  mom[iStep] = mom[iStep-1];
505  pos[iStep] = pos[iStep-1] - trajVec;
506 
507 #ifdef _DEBUG_ON_LEVEL_2
508  std::cout << "TARGET: " << iStep << ": " << pos[iStep-1][2] << " ==================>>> " << pos[iStep][2] << std::endl;
509  std::cout << mom[iStep-1][0]/mom[iStep-1][2] << " " << mom[iStep-1][1]/mom[iStep-1][2] << " " << mom[iStep-1][2] << " ";
510  std::cout << pos[iStep-1][0] << " " << pos[iStep-1][1] << " " << pos[iStep-1][2] << std::endl << std::endl;
511  std::cout << mom[iStep][0]/mom[iStep][2] << " " << mom[iStep][1]/mom[iStep][2] << " " << mom[iStep][2] << " ";
512  std::cout << pos[iStep][0] << " " << pos[iStep][1] << " " << pos[iStep][2] << std::endl << std::endl;
513 #endif
514  }
515 
516  //Now the swimming is done, find the point with closest distance of approach, let iStep store the index of that step
517  double dca_min = 1E9;
518  double dca_xmin = 1E9;
519  double dca_ymin = 1E9;
520 
521  iStep = NSTEPS_FMAG+NSTEPS_TARGET; // set the default point to the most upstream
522  int iStep_x = iStep; // the point when track cross beam line in X and in Y
523  int iStep_y = iStep; // both intialized with the most upstream position
524  for(int i = 0; i < NSTEPS_FMAG+NSTEPS_TARGET+1; ++i)
525  {
526  if(FMAGSTR*charge*mom[i].Px() < 0.) continue; // this is the upstream accidental cross, ignore
527 
528  double dca = (pos[i] - TVector3(X_BEAM, Y_BEAM, pos[i].Z())).Perp();
529  if(dca < dca_min)
530  {
531  dca_min = dca;
532  iStep = i;
533  }
534 
535  double dca_x = fabs(pos[i].X() - X_BEAM);
536  if(dca_x < dca_xmin)
537  {
538  dca_xmin = dca_x;
539  iStep_x = i;
540  }
541 
542  double dca_y = fabs(pos[i].Y() - Y_BEAM);
543  if(dca_y < dca_ymin)
544  {
545  dca_ymin = dca_y;
546  iStep_y = i;
547  }
548  }
549 
550  setVertexFast(mom[iStep], pos[iStep]);
551  setDumpFacePos(pos[NSTEPS_FMAG]);
552  setDumpFaceMom(mom[NSTEPS_FMAG]);
553  setTargetPos(fDumpFacePos + TVector3(fDumpFaceMom.Px()/fDumpFaceMom.Pz()*Z_TARGET, fDumpFaceMom.Py()/fDumpFaceMom.Pz()*Z_TARGET, Z_TARGET));
554  setTargetMom(mom[NSTEPS_FMAG]);
555 
556  double dz_x = -pos[iStep_x].X()/mom[iStep_x].Px()*mom[iStep_x].Pz();
557  setXVertexPos(pos[iStep_x] + TVector3(mom[iStep_x].Px()/mom[iStep_x].Pz()*dz_x, mom[iStep_x].Py()/mom[iStep_x].Pz()*dz_x, dz_x));
558  setXVertexMom(mom[iStep_x]);
559 
560  double dz_y = -pos[iStep_y].Y()/mom[iStep_y].Py()*mom[iStep_y].Pz();
561  setYVertexPos(pos[iStep_y] + TVector3(mom[iStep_y].Px()/mom[iStep_y].Pz()*dz_y, mom[iStep_y].Py()/mom[iStep_y].Pz()*dz_y, dz_y));
562  setYVertexMom(mom[iStep_y]);
563 
564 #ifdef _DEBUG_ON_LEVEL_2
565  std::cout << "The one with minimum DCA is: " << iStep << ": " << std::endl;
566  std::cout << mom[iStep][0]/mom[iStep][2] << " " << mom[iStep][1]/mom[iStep][2] << " " << mom[iStep][2] << " ";
567  std::cout << pos[iStep][0] << " " << pos[iStep][1] << " " << pos[iStep][2] << std::endl << std::endl;
568 #endif
569 
570 #ifdef _ENABLE_KF
571  if(hyptest) updateVtxHypothesis();
572 #endif
573 
574  if(cleanupPos) delete[] pos;
575  if(cleanupMom) delete[] mom;
576 }
577 
578 void SRecTrack::print(std::ostream& os) const
579 {
580  os << "=============== Reconstructed track ==================" << std::endl;
581  os << "This candidate has " << fHitIndex.size() << " hits!" << std::endl;
582  os << "Most upstream momentum is: " << 1./fabs((fState[0])[0][0]) << std::endl;
583  os << "Chi square of the track is: " << fChisq << std::endl;
584 
585  os << "Current vertex position: " << std::endl;
586  for(Int_t i = 0; i < 3; i++) os << fVertexPos[i] << " ";
587  os << std::endl;
588 
589  os << "Momentum at vertex: " << 1./fabs(fStateVertex[0][0]) << std::endl;
590  os << "Chi square at vertex: " << fChisqVertex << std::endl;
591 }
592 
593 void SRecTrack::printGF(std::ostream& os) const
594 {
595  os << "SRecTrack::PrintGF():\n"
596  << " DetPlaneVecO";
597  for (unsigned int iv = 0; iv < fGFDetPlaneVec[0].size(); iv++) {
598  const TVector3* vec = &fGFDetPlaneVec[0][iv];
599  os << " " << iv << ":(" << vec->X() << " " << vec->Y() << " " << vec->Z() << ")" ;
600  }
601  os << "\n DetPlaneVecU";
602  for (unsigned int iv = 0; iv < fGFDetPlaneVec[1].size(); iv++) {
603  const TVector3* vec = &fGFDetPlaneVec[1][iv];
604  os << " " << iv << ":(" << vec->X() << " " << vec->Y() << " " << vec->Z() << ")" ;
605  }
606  os << "\n DetPlaneVecV";
607  for (unsigned int iv = 0; iv < fGFDetPlaneVec[2].size(); iv++) {
608  const TVector3* vec = &fGFDetPlaneVec[2][iv];
609  os << " " << iv << ":(" << vec->X() << " " << vec->Y() << " " << vec->Z() << ")" ;
610  }
611  os << "\n AuxInfo";
612  for (unsigned int iv = 0; iv < fGFAuxInfo.size(); iv++) {
613  os << " [" << iv << "]";
614  for (int ir = 0; ir < fGFAuxInfo[iv].GetNrows(); ir++) os << " " << fGFAuxInfo[iv](ir);
615  }
616  os << "\n StateVec";
617  for (unsigned int iv = 0; iv < fGFStateVec.size(); iv++) {
618  os << " [" << iv << "]";
619  for (int ir = 0; ir < fGFStateVec[iv].GetNrows(); ir++) os << " " << fGFStateVec[iv](ir);
620  }
621  os << "\n Cov";
622  for (unsigned int iv = 0; iv < fGFCov.size(); iv++) {
623  os << " [" << iv << "]";
624  for (int ir = 0; ir < fGFCov[iv].GetNrows(); ir++) {
625  os << " (";
626  for (int ic = 0; ic < fGFCov[iv].GetNcols(); ic++) os << " " << fGFCov[iv][ir][ic];
627  os << ")";
628  }
629  }
630  os << std::endl;
631 }
632 
634 {
635  TLorentzVector pos, neg;
636  if(choice == 0)
637  {
638  pos = p_pos;
639  neg = p_neg;
640  }
641  else if(choice == 1)
642  {
643  pos = p_pos_target;
644  neg = p_neg_target;
645  }
646  else if(choice == 2)
647  {
648  pos = p_pos_dump;
649  neg = p_neg_dump;
650  }
651  else
652  {
653  std::cout << "ERROR!!! In SRecDimuon::calcVariables(choice), only three dimuon vertex choice numbers are provided - ";
654  std::cout << " 0 (default) = fitted vertex position, 1 = target position, 2 = dump position" << std::endl;
655 
656  //Unphysical error values
657  mass = -1.;
658  xF = -10.;
659  x1 = -10.;
660  x2 = -10.;
661  pT = -1.;
662  costh = -10.;
663  phi = -10.;
664  mass_single = -1.;
665 
666  return;
667  }
668 
669  Double_t mp = 0.938;
670  Double_t ebeam = 120.;
671 
672  TLorentzVector p_beam(0., 0., sqrt(ebeam*ebeam - mp*mp), ebeam);
673  TLorentzVector p_target(0., 0., 0., mp);
674 
675  TLorentzVector p_cms = p_beam + p_target;
676  TLorentzVector p_sum = pos + neg;
677 
678  mass = p_sum.M();
679  pT = p_sum.Perp();
680 
681  x1 = (p_target*p_sum)/(p_target*p_cms);
682  x2 = (p_beam*p_sum)/(p_beam*p_cms);
683 
684  Double_t s = p_cms.M2();
685  TVector3 bv_cms = p_cms.BoostVector();
686  p_sum.Boost(-bv_cms);
687  xF = 2.*p_sum.Pz()/TMath::Sqrt(s)/(1. - mass*mass/s);
688 
689  costh = 2.*(p_neg.E()*pos.Pz() - pos.E()*neg.Pz())/mass/sqrt(mass*mass + pT*pT);
690  phi = atan2(2.*sqrt(mass*mass + pT*pT)*(neg.X()*pos.Y() - pos.X()*neg.Y()), mass*(pos.X()*pos.X() - neg.X()*neg.X() + pos.Y()*pos.Y() - neg.Y()*neg.Y()));
692 }
693 
695 {
696  //Chisq of vertex fit
697  if(chisq_kf > 15. || chisq_kf < 0.) return false;
698 
699  //Kinematic cuts
700  if(FMAGSTR*(p_pos.Px() - p_neg.Px()) < 0.) return false;
701  if(fabs(xF) > 1.) return false;
702  if(x1 < 0. || x1 > 1.) return false;
703  if(x2 < 0. || x2 > 1.) return false;
704  if(mass < 0. || mass > 10.) return false;
705  if(fabs(vtx.X()) > 2.) return false;
706  if(fabs(vtx.Y()) > 2.) return false;
707  if(fabs(p_pos.Px() + p_neg.Px()) > 3.) return false;
708  if(fabs(p_pos.Py() + p_neg.Py()) > 3.) return false;
709  if(vtx.Z() > 200. || vtx.Z() < -300.) return false;
710  if(p_pos.Pz() + p_neg.Pz() > 120. || p_pos.Pz() + p_neg.Pz() < 30.) return false;
711 
712  //Track separation cuts
713  if(fabs(vtx_pos.Z() - vtx_neg.Z()) > 250.) return false;
714 
715  //Everything is fine
716  return true;
717 }
718 
720 {
721  //single muon vertex
722  if(vtx_pos.Z() > -100. || vtx_pos.Z() < -500.) return false;
723  if(vtx_neg.Z() > -100. || vtx_neg.Z() < -500.) return false;
724 
725  //Track projection comparison
726  double pzp = p_pos.Pz();
727  if(proj_dump_pos.Perp() - proj_target_pos.Perp() < 9.44301-0.356141*pzp+0.00566071*pzp*pzp-3.05556e-05*pzp*pzp*pzp) return false;
728 
729  double pzm = p_neg.Pz();
730  if(proj_dump_neg.Perp() - proj_target_neg.Perp() < 9.44301-0.356141*pzm+0.00566071*pzm*pzm-3.05556e-05*pzm*pzm*pzm) return false;
731 
732  return true;
733 }
734 
736 {
737  //single muon vertex
738  if(vtx_pos.Z() > 150.) return false;
739  if(vtx_neg.Z() > 150.) return false;
740 
741  //Track projection comparison
742  double pzp = p_pos.Pz();
743  if(proj_target_pos.Perp() - proj_dump_pos.Perp() < 9.44301-0.356141*pzp+0.00566071*pzp*pzp-3.05556e-05*pzp*pzp*pzp) return false;
744 
745  double pzm = p_neg.Pz();
746  if(proj_target_neg.Perp() - proj_dump_neg.Perp() < 9.44301-0.356141*pzm+0.00566071*pzm*pzm-3.05556e-05*pzm*pzm*pzm) return false;
747 
748  return true;
749 }
750 
752 {
753  fRunID = -1;
754  fSpillID = -1;
755  fEventID = -1;
756  fRecStatus = 0;
757 
758  fSource1 = -1;
759  fSource2 = -1;
760 
761  clear();
762 }
763 
765 {
766  clear();
767  setEventInfo(rawEvent);
768 
769  for(Int_t i = 0; i < rawEvent->getNHitsAll(); i++)
770  {
771  fLocalID.insert(std::map<Int_t, Int_t>::value_type(rawEvent->getHit(i).index, i));
772  }
773 }
774 
776 {
777  fRunID = rawEvent->getRunID();
778  fSpillID = rawEvent->getSpillID();
779  fEventID = rawEvent->getEventID();
780 
781  fTriggerBits = rawEvent->getTriggerBits();
782  fTargetPos = rawEvent->getTargetPos();
783 }
784 
785 std::vector<Int_t> SRecEvent::getChargedTrackIDs(Int_t charge)
786 {
787  std::vector<Int_t> trkIDs;
788  trkIDs.clear();
789 
790  Int_t nTracks = getNTracks();
791  for(Int_t i = 0; i < nTracks; i++)
792  {
793  if(fAllTracks[i].getCharge() == charge)
794  {
795  trkIDs.push_back(i);
796  }
797  }
798 
799  return trkIDs;
800 }
801 
803 {
804  fAllTracks.clear();
805  fLocalID.clear();
806  fDimuons.clear();
807 
808  fRecStatus = 0;
809 }
810 
812 {
813  fAllTracks.clear();
814  fLocalID.clear();
815  fRecStatus = 0;
816 }
817 
819 
825 {
826  fDimuons.clear();
827 }
const double FMAGSTR
Definition: MacroCommon.h:2
const double KMAGSTR
Definition: MacroCommon.h:3
ClassImp(SRecTrack) ClassImp(SRecDimuon) ClassImp(SRecEvent) namespace
Definition: SRecEvent.cxx:22
Int_t index
Definition: SRawEvent.h:77
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
void setZ(double z)
Definition: KalmanUtil.h:137
double getChisq()
Definition: KalmanUtil.h:127
TMatrixD & getMeasurement()
Definition: KalmanUtil.h:110
TMatrixD & getProjector()
Definition: KalmanUtil.h:122
TrkPar & getFiltered()
Definition: KalmanUtil.h:107
TMatrixD & getMeasurementCov()
Definition: KalmanUtil.h:111
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
virtual bool get_BoolFlag(const std::string &name) const
Definition: PHFlag.cc:151
Int_t getSpillID()
Definition: SRawEvent.h:151
Int_t getEventID()
Definition: SRawEvent.h:150
Int_t getTriggerBits()
Set/get the trigger types.
Definition: SRawEvent.h:174
Int_t getNHitsAll()
Definition: SRawEvent.h:118
Hit getHit(Int_t index)
Definition: SRawEvent.h:144
Int_t getTargetPos()
Definition: SRawEvent.h:195
Int_t getRunID()
Definition: SRawEvent.h:149
Double_t costh
Definition: SRecEvent.h:395
Double_t mass_single
Definition: SRecEvent.h:397
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
TLorentzVector p_neg_target
Definition: SRecEvent.h:370
TVector3 proj_target_neg
Definition: SRecEvent.h:386
bool isDump()
Dump dimuon.
Definition: SRecEvent.cxx:735
int isValid() const
isValid returns non zero if object contains vailid data
Definition: SRecEvent.cxx:694
Double_t xF
Definition: SRecEvent.h:392
TVector3 proj_dump_pos
Definition: SRecEvent.h:385
TLorentzVector p_neg_single
Definition: SRecEvent.h:376
TVector3 proj_dump_neg
Definition: SRecEvent.h:387
TLorentzVector p_pos_target
Track momentum projections at different location.
Definition: SRecEvent.h:369
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
Double_t x2
Definition: SRecEvent.h:394
TLorentzVector p_pos
4-momentum of the muon tracks after vertex fit
Definition: SRecEvent.h:365
Double_t x1
Definition: SRecEvent.h:393
Double_t mass
Kinematic variables.
Definition: SRecEvent.h:390
TVector3 vtx_neg
Definition: SRecEvent.h:381
bool isTarget()
Target dimuon.
Definition: SRecEvent.cxx:719
Double_t pT
Definition: SRecEvent.h:391
TVector3 vtx_pos
Definition: SRecEvent.h:380
TVector3 proj_target_pos
Track projections at different location.
Definition: SRecEvent.h:384
TLorentzVector p_pos_dump
Definition: SRecEvent.h:371
TLorentzVector p_neg_dump
Definition: SRecEvent.h:372
Double_t phi
Definition: SRecEvent.h:396
std::vector< Int_t > getChargedTrackIDs(Int_t charge)
Get track IDs.
Definition: SRecEvent.cxx:785
void setRawEvent(SRawEvent *rawEvent)
directly setup everything by raw event
Definition: SRecEvent.cxx:764
void clearTracks()
Definition: SRecEvent.cxx:811
void setEventInfo(SRawEvent *rawEvent)
Set/Get event info.
Definition: SRecEvent.cxx:775
void clear()
Clear everything.
Definition: SRecEvent.cxx:802
Int_t getNTracks()
Get tracks.
Definition: SRecEvent.h:455
void clearDimuons()
Clear the dimuon list.
Definition: SRecEvent.cxx:824
void getExpPosErrorFast(Double_t z, Double_t &dx, Double_t &dy, Int_t iNode=-1)
Definition: SRecEvent.cxx:288
Double_t getChisq() const
Definition: SRecEvent.h:104
Int_t getNHits() const
Definition: SRecEvent.h:102
Int_t getCharge() const
Gets.
Definition: SRecEvent.h:101
Double_t getMomentum(const TMatrixD &state, Double_t &px, Double_t &py, Double_t &pz) const
Definition: SRecEvent.cxx:322
bool isDump()
Definition: SRecEvent.cxx:397
Double_t getPosition(const TMatrixD &state, Double_t &x, Double_t &y) const
Definition: SRecEvent.cxx:332
bool isVertexValid() const
Vertex stuff.
Definition: SRecEvent.cxx:241
void setChisqDump(Double_t chisq)
Definition: SRecEvent.h:215
TLorentzVector getMomentumVertex()
Get the vertex info.
Definition: SRecEvent.cxx:340
void setDumpPos(TVector3 pos)
Set mom/pos at a given location.
Definition: SRecEvent.h:203
void getExpPositionFast(Double_t z, Double_t &x, Double_t &y, Int_t iNode=-1)
Definition: SRecEvent.cxx:267
bool operator<(const SRecTrack &elem) const
Comparitor.
Definition: SRecEvent.cxx:144
void setTargetPos(TVector3 pos)
Definition: SRecEvent.h:205
void setZVertex(Double_t z, bool update=true)
Definition: SRecEvent.cxx:168
void print(std::ostream &os=std::cout) const
Debugging output.
Definition: SRecEvent.cxx:578
void insertGFState(const genfit::MeasuredStateOnPlane &msop)
Definition: SRecEvent.cxx:351
Double_t getProb() const
Definition: SRecEvent.cxx:163
void setYVertexPos(TVector3 pos)
Definition: SRecEvent.h:207
void adjustKMag(double kmagStr)
Fast-adjust of kmag.
Definition: SRecEvent.cxx:363
TVector3 getMomentumVecSt1() const
Definition: SRecEvent.h:129
void setChisqTarget(Double_t chisq)
Definition: SRecEvent.h:216
void setXVertexPos(TVector3 pos)
Definition: SRecEvent.h:206
void setYVertexMom(TVector3 mom)
Definition: SRecEvent.h:213
void printGF(std::ostream &os=std::cout) const
Definition: SRecEvent.cxx:593
Int_t getNHitsInStation(Int_t stationID)
Definition: SRecEvent.cxx:149
int isValid() const
isValid returns non zero if object contains vailid data
Definition: SRecEvent.cxx:371
bool isTarget()
Overall track quality cut.
Definition: SRecEvent.cxx:392
void setDumpFaceMom(TVector3 mom)
Definition: SRecEvent.h:210
Int_t getNearestNode(Double_t z)
Definition: SRecEvent.cxx:249
void setXVertexMom(TVector3 mom)
Definition: SRecEvent.h:212
Double_t getZVertex()
Definition: SRecEvent.h:179
void setTargetMom(TVector3 mom)
Definition: SRecEvent.h:211
void setChisqUpstream(Double_t chisq)
Definition: SRecEvent.h:217
Double_t getExpMomentumFast(Double_t z, Double_t &px, Double_t &py, Double_t &pz, Int_t iNode=-1)
Definition: SRecEvent.cxx:311
void setDumpMom(TVector3 mom)
Definition: SRecEvent.h:209
void swimToVertex(TVector3 *pos=nullptr, TVector3 *mom=nullptr, bool hyptest=true)
Simple swim to vertex.
Definition: SRecEvent.cxx:402
void updateVtxHypothesis()
Definition: SRecEvent.cxx:213
void setDumpFacePos(TVector3 pos)
Definition: SRecEvent.h:204
void setVertexFast(TVector3 mom, TVector3 pos)
Plain setting, no KF-related stuff.
Definition: SRecEvent.cxx:227
double get_y()
Definition: KalmanUtil.h:59
TVector3 get_mom_vec()
Definition: KalmanUtil.h:65
double get_dydz()
Definition: KalmanUtil.h:62
TMatrixD _state_kf
State vectors and its covariance.
Definition: KalmanUtil.h:85
TMatrixD _covar_kf
Definition: KalmanUtil.h:86
double get_dxdz()
Definition: KalmanUtil.h:61
double _z
Definition: KalmanUtil.h:87
double get_x()
Definition: KalmanUtil.h:58
static recoConsts * instance()
Definition: recoConsts.cc:7