CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFTrackTransformer.cc
Go to the documentation of this file.
1 //
2 // -*- C++ -*-
3 // Package: PFTracking
4 // Class: PFTrackTransformer
5 //
6 // Original Author: Michele Pioppi
7 // Other Author: Daniele Benedetti
8 
10 
13 
17 
20 
22 // Add by Daniele
27 
28 using namespace std;
29 using namespace reco;
30 using namespace edm;
31 
32 
33 
35  LogInfo("PFTrackTransformer")<<"PFTrackTransformer built";
36 
37  PFGeometry pfGeometry;
38  onlyprop_=false;
39 }
40 
42 
43 }
44 
45 
46 bool
48  const reco::Track& track,
49  const Trajectory& traj) const {
50 
51  LogDebug("PFTrackTransformer")<<"Trajectory propagation started";
52  using namespace reco;
53  using namespace std;
54 
55  float PT= track.pt();
56  float pfmass= (pftrack.algoType()==reco::PFRecTrack::KF_ELCAND) ? 0.0005 : 0.139;
57  float pfenergy=sqrt((pfmass*pfmass)+(track.p()*track.p()));
58  // closest approach
59  BaseParticlePropagator theParticle =
62  track.py(),
63  track.pz(),
64  pfenergy),
65  XYZTLorentzVector(track.vertex().x(),
66  track.vertex().y(),
67  track.vertex().z(),
68  0.)),
69  0.,0.,B_.z());
70 
71  theParticle.setCharge(track.charge());
72  float pfoutenergy=sqrt((pfmass*pfmass)+track.outerMomentum().Mag2());
73  BaseParticlePropagator theOutParticle =
76  track.outerMomentum().y(),
77  track.outerMomentum().z(),
78  pfoutenergy),
79  XYZTLorentzVector(track.outerPosition().x(),
80  track.outerPosition().y(),
81  track.outerPosition().z(),
82  0.)),
83  0.,0.,B_.z());
84  theOutParticle.setCharge(track.charge());
85 
86 
87  math::XYZTLorentzVector momClosest
88  = math::XYZTLorentzVector(track.px(), track.py(),
89  track.pz(), track.p());
90  math::XYZPoint posClosest = track.vertex();
91 
92  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ClosestApproach,
93  posClosest,momClosest));
94 
95 
96  //BEAMPIPE
97  theParticle.setPropagationConditions(PFGeometry::outerRadius(PFGeometry::BeamPipe),
99  theParticle.propagate();
100  if(theParticle.getSuccess()!=0)
101  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::BeamPipeOrEndVertex,
102  math::XYZPoint(theParticle.vertex()),
103  math::XYZTLorentzVector(theParticle.momentum())));
104  else {
105  PFTrajectoryPoint dummyMaxSh;
106  pftrack.addPoint(dummyMaxSh);
107  }
108 
109 
110 
111  //trajectory points
112 
113  if (!onlyprop_){
114  bool direction =(traj.direction() == alongMomentum);
115  vector<TrajectoryMeasurement> measurements =traj.measurements();
116  int iTrajFirst = (direction) ? 0 : measurements.size() - 1;
117  int increment = (direction) ? +1 : -1;
118  int iTrajLast = (direction) ? int(measurements.size()) : -1;
119 
120 
121  for (int iTraj = iTrajFirst; iTraj != iTrajLast; iTraj += increment) {
122  GlobalPoint v=measurements[iTraj].updatedState().globalPosition();
123  GlobalVector p=measurements[iTraj].updatedState().globalMomentum();
124  unsigned int iid=measurements[iTraj].recHit()->det()->geographicalId().rawId();
125  pftrack.addPoint(PFTrajectoryPoint(iid,-1,
126  math::XYZPoint(v.x(), v.y(), v.z()),
127  math::XYZTLorentzVector(p.x(),p.y(),p.z(),p.mag())));
128  }
129  }
130 
131  bool isBelowPS=false;
132  theOutParticle.propagateToPreshowerLayer1(false);
133  if(theOutParticle.getSuccess()!=0)
134  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS1,
135  math::XYZPoint(theOutParticle.vertex()),
136  math::XYZTLorentzVector(theOutParticle.momentum())));
137  else {
138  PFTrajectoryPoint dummyPS1;
139  pftrack.addPoint(dummyPS1);
140  }
141 
142 
143  theOutParticle.propagateToPreshowerLayer2(false);
144  if(theOutParticle.getSuccess()!=0){
145  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS2,
146  math::XYZPoint(theOutParticle.vertex()),
147  math::XYZTLorentzVector(theOutParticle.momentum())));
148  isBelowPS=true;
149  } else {
150  PFTrajectoryPoint dummyPS2;
151  pftrack.addPoint(dummyPS2);
152  }
153 
154  theOutParticle.propagateToEcalEntrance(false);
155 
156  if(theOutParticle.getSuccess()!=0){
157  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALEntrance,
158  math::XYZPoint(theOutParticle.vertex()),
159  math::XYZTLorentzVector(theOutParticle.momentum())));
160  double ecalShowerDepth
161  = PFCluster::getDepthCorrection(theOutParticle.momentum().E(),
162  isBelowPS,
163  false);
164 
165  math::XYZPoint meanShower=math::XYZPoint(theOutParticle.vertex())+
166  math::XYZTLorentzVector(theOutParticle.momentum()).Vect().Unit()*ecalShowerDepth;
167 
168  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALShowerMax,
169  meanShower,
170  math::XYZTLorentzVector(theOutParticle.momentum())));}
171  else {
172  if (PT>5.)
173  LogWarning("PFTrackTransformer")<<"KF TRACK "<<pftrack<< " PROPAGATION TO THE ECAL HAS FAILED";
174  PFTrajectoryPoint dummyECAL;
175  pftrack.addPoint(dummyECAL);
176  PFTrajectoryPoint dummyMaxSh;
177  pftrack.addPoint(dummyMaxSh);
178  }
179 
180 
181 
182  //HCAL entrance
183  theOutParticle.propagateToHcalEntrance(false);
184  if(theOutParticle.getSuccess()!=0)
185  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALEntrance,
186  math::XYZPoint(theOutParticle.vertex()),
187  math::XYZTLorentzVector(theOutParticle.momentum())));
188  else{
189  if (PT>5.)
190  LogWarning("PFTrackTransformer")<<"KF TRACK "<<pftrack<< " PROPAGATION TO THE HCAL ENTRANCE HAS FAILED";
191  PFTrajectoryPoint dummyHCALentrance;
192  pftrack.addPoint(dummyHCALentrance);
193  }
194 
195  //HCAL exit
196  theOutParticle.propagateToHcalExit(false);
197  if(theOutParticle.getSuccess()!=0)
198  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALExit,
199  math::XYZPoint(theOutParticle.vertex()),
200  math::XYZTLorentzVector(theOutParticle.momentum())));
201  else{
202  if (PT>5.)
203  LogWarning("PFTrackTransformer")<<"KF TRACK "<<pftrack<< " PROPAGATION TO THE HCAL EXIT HAS FAILED";
204  PFTrajectoryPoint dummyHCALexit;
205  pftrack.addPoint(dummyHCALexit);
206  }
207 
208  return true;
209 }
210 bool
212  const reco::Track& track,
213  const Trajectory& traj,
214  const bool& GetMode) const {
215 
216  float PT= track.pt();
217  // Trajectory for each trajectory point
218 
219  bool direction =(traj.direction() == alongMomentum);
220  vector<TrajectoryMeasurement> measurements =traj.measurements();
221  int iTrajFirst = (direction) ? 0 : measurements.size() - 1;
222  int increment = (direction) ? +1 : -1;
223  int iTrajLast = (direction) ? int(measurements.size()) : -1;
224 
225 
226  unsigned int iTrajPos = 0;
227  for (int iTraj = iTrajFirst; iTraj != iTrajLast; iTraj += increment) {
228 
229  GlobalPoint v=measurements[iTraj].updatedState().globalPosition();
230  PFGsfHelper* PFGsf = new PFGsfHelper(measurements[iTraj]);
231  //if (PFGsf->isValid()){
232  bool ComputeMODE = GetMode;
233  GlobalVector p = PFGsf->computeP(ComputeMODE);
234  double DP = PFGsf->fittedDP();
235  double SigmaDP = PFGsf->sigmafittedDP();
236  unsigned int iid=measurements[iTraj].recHit()->det()->geographicalId().rawId();
237  delete PFGsf;
238 
239  // -------------------------- Fill GSF Track -------------------------------------
240 
241 
242  float pfmass= (pftrack.algoType()==reco::PFRecTrack::KF_ELCAND) ? 0.0005 : 0.139;
243  float ptot = sqrt((p.x()*p.x())+(p.y()*p.y())+(p.z()*p.z()));
244  float pfenergy=sqrt((pfmass*pfmass)+(ptot *ptot));
245 
246  if (iTraj == iTrajFirst) {
247 
248  math::XYZTLorentzVector momClosest
249  = math::XYZTLorentzVector(p.x(), p.y(),
250  p.z(), ptot);
251  math::XYZPoint posClosest = track.vertex();
252  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ClosestApproach,
253  posClosest,momClosest));
254 
255  BaseParticlePropagator theInnerParticle =
258  p.y(),
259  p.z(),
260  pfenergy),
261  XYZTLorentzVector(track.vertex().x(),
262  track.vertex().y(),
263  track.vertex().z(),
264  0.)), //DANIELE Same thing v.x(),v.y(),v.()?
265  0.,0.,B_.z());
266  theInnerParticle.setCharge(track.charge());
267 
268  //BEAMPIPE
271  theInnerParticle.propagate();
272  if(theInnerParticle.getSuccess()!=0)
273  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::BeamPipeOrEndVertex,
274  math::XYZPoint(theInnerParticle.vertex()),
275  math::XYZTLorentzVector(theInnerParticle.momentum())));
276  else {
277  PFTrajectoryPoint dummyMaxSh;
278  pftrack.addPoint(dummyMaxSh);
279  }
280 
281  // First Point for the trajectory == Vertex ??
282  pftrack.addPoint(PFTrajectoryPoint(iid,-1,
283  math::XYZPoint(v.x(), v.y(), v.z()),
284  math::XYZTLorentzVector(p.x(),p.y(),p.z(),p.mag())));
285 
286 
287  }
288  if (iTraj != iTrajFirst && iTraj != (abs(iTrajLast)-1)) {
289  pftrack.addPoint(PFTrajectoryPoint(iid,-1,
290  math::XYZPoint(v.x(), v.y(), v.z()),
291  math::XYZTLorentzVector(p.x(),p.y(),p.z(),p.mag())));
292 
293 
294  }
295  if (iTraj == (abs(iTrajLast)-1)) {
296 
297  // Last Trajectory Meas
298  pftrack.addPoint(PFTrajectoryPoint(iid,-1,
299  math::XYZPoint(v.x(), v.y(), v.z()),
300  math::XYZTLorentzVector(p.x(),p.y(),p.z(),p.mag())));
301 
302 
303 
304 
305  BaseParticlePropagator theOutParticle =
308  p.y(),
309  p.z(),
310  pfenergy),
311  XYZTLorentzVector(v.x(),
312  v.y(),
313  v.z(),
314  0.)),
315  0.,0.,B_.z());
316  theOutParticle.setCharge(track.charge());
317  bool isBelowPS=false;
318  theOutParticle.propagateToPreshowerLayer1(false);
319  if(theOutParticle.getSuccess()!=0)
320  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS1,
321  math::XYZPoint(theOutParticle.vertex()),
322  math::XYZTLorentzVector(theOutParticle.momentum())));
323  else {
324  PFTrajectoryPoint dummyPS1;
325  pftrack.addPoint(dummyPS1);
326  }
327 
328 
329  theOutParticle.propagateToPreshowerLayer2(false);
330  if(theOutParticle.getSuccess()!=0){
331  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS2,
332  math::XYZPoint(theOutParticle.vertex()),
333  math::XYZTLorentzVector(theOutParticle.momentum())));
334  isBelowPS=true;
335  } else {
336  PFTrajectoryPoint dummyPS2;
337  pftrack.addPoint(dummyPS2);
338  }
339 
340  theOutParticle.propagateToEcalEntrance(false);
341 
342  if(theOutParticle.getSuccess()!=0){
343  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALEntrance,
344  math::XYZPoint(theOutParticle.vertex()),
345  math::XYZTLorentzVector(theOutParticle.momentum())));
346  double ecalShowerDepth
347  = PFCluster::getDepthCorrection(theOutParticle.momentum().E(),
348  isBelowPS,
349  false);
350 
351  math::XYZPoint meanShower=math::XYZPoint(theOutParticle.vertex())+
352  math::XYZTLorentzVector(theOutParticle.momentum()).Vect().Unit()*ecalShowerDepth;
353 
354  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALShowerMax,
355  meanShower,
356  math::XYZTLorentzVector(theOutParticle.momentum())));}
357  else {
358  if (PT>5.)
359  LogWarning("PFTrackTransformer")<<"GSF TRACK "<<pftrack<< " PROPAGATION TO THE ECAL HAS FAILED";
360  PFTrajectoryPoint dummyECAL;
361  pftrack.addPoint(dummyECAL);
362  PFTrajectoryPoint dummyMaxSh;
363  pftrack.addPoint(dummyMaxSh);
364  }
365 
366 
367 
368  //HCAL entrance
369  theOutParticle.propagateToHcalEntrance(false);
370  if(theOutParticle.getSuccess()!=0)
371  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALEntrance,
372  math::XYZPoint(theOutParticle.vertex()),
373  math::XYZTLorentzVector(theOutParticle.momentum())));
374  else{
375  if (PT>5.)
376  LogWarning("PFTrackTransformer")<<"GSF TRACK "<<pftrack<< " PROPAGATION TO THE HCAL ENTRANCE HAS FAILED";
377  PFTrajectoryPoint dummyHCALentrance;
378  pftrack.addPoint(dummyHCALentrance);
379  }
380 
381  //HCAL exit
382  theOutParticle.propagateToHcalExit(false);
383  if(theOutParticle.getSuccess()!=0)
384  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALExit,
385  math::XYZPoint(theOutParticle.vertex()),
386  math::XYZTLorentzVector(theOutParticle.momentum())));
387  else{
388  if (PT>5.)
389  LogWarning("PFTrackTransformer")<<"GSF TRACK "<<pftrack<< " PROPAGATION TO THE HCAL EXIT HAS FAILED";
390  PFTrajectoryPoint dummyHCALexit;
391  pftrack.addPoint(dummyHCALexit);
392  }
393 
394  }
395 
396  // -------------------------- END GSF Track -------------------------------------
397 
398  // -------------------------- Fill Brem "Track" ---------------------------------
399  // Fill the brem for each traj point
400 
401  //check that the vertex of the brem is in the tracker volume
402  if ((v.perp()>110) ||(fabs(v.z())>280)) continue;
403  unsigned int iTrajPoint = iTrajPos + 2;
404  if(iid%2 == 1) iTrajPoint = 99;
405 
406  PFBrem brem(DP,SigmaDP,iTrajPoint);
407 
408 
409  GlobalVector p_gamma= p*(fabs(DP)/p.mag()); // Direction from the electron (tangent), DP without any sign!;
410  float e_gamma = fabs(DP); // DP = pout-pin so could be negative
411  BaseParticlePropagator theBremParticle =
413  RawParticle(XYZTLorentzVector(p_gamma.x(),
414  p_gamma.y(),
415  p_gamma.z(),
416  e_gamma),
417  XYZTLorentzVector(v.x(),
418  v.y(),
419  v.z(),
420  0.)),
421  0.,0.,B_.z());
422  int gamma_charge = 0;
423  theBremParticle.setCharge(gamma_charge);
424 
425 
426  // add TrajectoryPoint for Brem, PS, ECAL, ECALShowMax, HCAL
427  // Brem Entrance PS Layer1
428 
429  PFTrajectoryPoint dummyClosest; // Added just to have the right number order in PFTrack.cc
430  brem.addPoint(dummyClosest);
431 
432 
433  PFTrajectoryPoint dummyBeamPipe; // Added just to have the right number order in PFTrack.cc
434  brem.addPoint(dummyBeamPipe);
435 
436 
437 
438  bool isBelowPS=false;
439  theBremParticle.propagateToPreshowerLayer1(false);
440  if(theBremParticle.getSuccess()!=0)
441  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS1,
442  math::XYZPoint(theBremParticle.vertex()),
443  math::XYZTLorentzVector(theBremParticle.momentum())));
444  else {
445  PFTrajectoryPoint dummyPS1;
446  brem.addPoint(dummyPS1);
447  }
448 
449  // Brem Entrance PS Layer 2
450 
451  theBremParticle.propagateToPreshowerLayer2(false);
452  if(theBremParticle.getSuccess()!=0){
453  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS2,
454  math::XYZPoint(theBremParticle.vertex()),
455  math::XYZTLorentzVector(theBremParticle.momentum())));
456  isBelowPS=true;
457  } else {
458  PFTrajectoryPoint dummyPS2;
459  brem.addPoint(dummyPS2);
460  }
461 
462  theBremParticle.propagateToEcalEntrance(false);
463 
464  if(theBremParticle.getSuccess()!=0){
465  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALEntrance,
466  math::XYZPoint(theBremParticle.vertex()),
467  math::XYZTLorentzVector(theBremParticle.momentum())));
468  double ecalShowerDepth
469  = PFCluster::getDepthCorrection(theBremParticle.momentum().E(),
470  isBelowPS,
471  false);
472 
473  math::XYZPoint meanShower=math::XYZPoint(theBremParticle.vertex())+
474  math::XYZTLorentzVector(theBremParticle.momentum()).Vect().Unit()*ecalShowerDepth;
475 
476  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALShowerMax,
477  meanShower,
478  math::XYZTLorentzVector(theBremParticle.momentum())));}
479  else {
480  if ((DP>5.) && ((DP/SigmaDP)>3))
481  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE ECAL HAS FAILED";
482  PFTrajectoryPoint dummyECAL;
483  brem.addPoint(dummyECAL);
484  PFTrajectoryPoint dummyMaxSh;
485  brem.addPoint(dummyMaxSh);
486  }
487 
488 
489 
490  //HCAL entrance
491  theBremParticle.propagateToHcalEntrance(false);
492  if(theBremParticle.getSuccess()!=0)
493  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALEntrance,
494  math::XYZPoint(theBremParticle.vertex()),
495  math::XYZTLorentzVector(theBremParticle.momentum())));
496  else{
497  if ((DP>5.) && ((DP/SigmaDP)>3))
498  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE HCAL ENTRANCE HAS FAILED";
499  PFTrajectoryPoint dummyHCALentrance;
500  brem.addPoint(dummyHCALentrance);
501  }
502 
503  //HCAL exit
504  theBremParticle.propagateToHcalExit(false);
505  if(theBremParticle.getSuccess()!=0)
506  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALExit,
507  math::XYZPoint(theBremParticle.vertex()),
508  math::XYZTLorentzVector(theBremParticle.momentum())));
509  else{
510  if ((DP>5.) && ((DP/SigmaDP)>3))
511  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE HCAL EXIT HAS FAILED";
512  PFTrajectoryPoint dummyHCALexit;
513  brem.addPoint(dummyHCALexit);
514  }
515 
516  brem.calculatePositionREP();
517  pftrack.addBrem(brem);
518  iTrajPos++;
519  }
520  return true;
521 }
522 
523 bool
525  const reco::GsfTrack& track,
526  const MultiTrajectoryStateTransform& mtjstate) const {
527 
528  float PT= track.pt();
529  unsigned int iTrajPos = 0;
530  unsigned int iid = 0; // not anymore saved
531 
532 
533  // ***************************** INNER State *************************************
534  TrajectoryStateOnSurface inTSOS = mtjstate.innerStateOnSurface((track));
535  TrajectoryStateOnSurface outTSOS = mtjstate.outerStateOnSurface((track));
536 
537  if(!inTSOS.isValid() || !outTSOS.isValid()) {
538  if(!inTSOS.isValid())
539  LogWarning("PFTrackTransformer")<<" INNER TSOS NOT VALID ";
540  if(!outTSOS.isValid())
541  LogWarning("PFTrackTransformer")<<" OUTER TSOS NOT VALID ";
542  return false;
543  }
544 
545  GlobalVector InMom;
546  GlobalPoint InPos;
547  if(inTSOS.isValid()) {
548  mtsMode_->momentumFromModeCartesian(inTSOS,InMom);
549  mtsMode_->positionFromModeCartesian(inTSOS,InPos);
550  }
551  else {
552  InMom = GlobalVector(track.pxMode(),track.pyMode(),track.pzMode());
553  InPos = GlobalPoint(0.,0.,0.);
554  }
555 
556  float pfmass= (pftrack.algoType()==reco::PFRecTrack::KF_ELCAND) ? 0.0005 : 0.139;
557  float ptot = sqrt((InMom.x()*InMom.x())+(InMom.y()*InMom.y())+(InMom.z()*InMom.z()));
558  float pfenergy=sqrt((pfmass*pfmass)+(ptot *ptot));
559 
560  math::XYZTLorentzVector momClosest
561  = math::XYZTLorentzVector(InMom.x(), InMom.y(),
562  InMom.z(), ptot);
563  math::XYZPoint posClosest = track.vertex();
564  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ClosestApproach,
565  posClosest,momClosest));
566 
567  BaseParticlePropagator theInnerParticle =
569  InMom.y(),
570  InMom.z(),
571  pfenergy),
572  XYZTLorentzVector(track.vertex().x(),
573  track.vertex().y(),
574  track.vertex().z(),
575  0.)), //DANIELE Same thing v.x(),v.y(),v.()?
576  0.,0.,B_.z());
577  theInnerParticle.setCharge(track.charge()); // Use the chargeMode ??
578  //BEAMPIPE
581  theInnerParticle.propagate();
582  if(theInnerParticle.getSuccess()!=0)
583  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::BeamPipeOrEndVertex,
584  math::XYZPoint(theInnerParticle.vertex()),
585  math::XYZTLorentzVector(theInnerParticle.momentum())));
586  else {
587  PFTrajectoryPoint dummyBeam;
588  pftrack.addPoint(dummyBeam);
589  }
590 
591 
592  // first tjpoint
593  pftrack.addPoint(PFTrajectoryPoint(iid,-1,
594  math::XYZPoint(InPos.x(),InPos.y(), InPos.z()),
595  math::XYZTLorentzVector(InMom.x(),InMom.y(),InMom.z(),InMom.mag())));
596 
597 
598  //######### Photon at INNER State ##########
599 
600 
601  unsigned int iTrajPoint = iTrajPos + 2;
602  double dp_tang = ptot;
603  double sdp_tang = track.ptModeError()*(track.pMode()/track.ptMode());
604  PFBrem brem(dp_tang,sdp_tang,iTrajPoint);
605  BaseParticlePropagator theBremParticle =
608  InMom.y(),
609  InMom.z(),
610  dp_tang),
611  XYZTLorentzVector(InPos.x(),
612  InPos.y(),
613  InPos.z(),
614  0.)),
615  0.,0.,B_.z());
616  int gamma_charge = 0;
617  theBremParticle.setCharge(gamma_charge);
618  // add TrajectoryPoint for Brem, PS, ECAL, ECALShowMax, HCAL
619  // Brem Entrance PS Layer1
620  PFTrajectoryPoint dummyClosest; // Added just to have the right number order in PFTrack.cc
621  brem.addPoint(dummyClosest);
622 
623 
624  PFTrajectoryPoint dummyBeamPipe; // Added just to have the right number order in PFTrack.cc
625  brem.addPoint(dummyBeamPipe);
626 
627 
628 
629  bool isBelowPS=false;
630  theBremParticle.propagateToPreshowerLayer1(false);
631  if(theBremParticle.getSuccess()!=0)
632  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS1,
633  math::XYZPoint(theBremParticle.vertex()),
634  math::XYZTLorentzVector(theBremParticle.momentum())));
635  else {
636  PFTrajectoryPoint dummyPS1;
637  brem.addPoint(dummyPS1);
638  }
639 
640  // Brem Entrance PS Layer 2
641 
642  theBremParticle.propagateToPreshowerLayer2(false);
643  if(theBremParticle.getSuccess()!=0){
644  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS2,
645  math::XYZPoint(theBremParticle.vertex()),
646  math::XYZTLorentzVector(theBremParticle.momentum())));
647  isBelowPS=true;
648  } else {
649  PFTrajectoryPoint dummyPS2;
650  brem.addPoint(dummyPS2);
651  }
652 
653  theBremParticle.propagateToEcalEntrance(false);
654 
655  if(theBremParticle.getSuccess()!=0){
656  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALEntrance,
657  math::XYZPoint(theBremParticle.vertex()),
658  math::XYZTLorentzVector(theBremParticle.momentum())));
659 
660  // for the first brem give a low default DP of 100 MeV.
661  double EDepthCorr = 0.01;
662  double ecalShowerDepth
663  = PFCluster::getDepthCorrection(EDepthCorr,
664  isBelowPS,
665  false);
666 
667  math::XYZPoint meanShower=math::XYZPoint(theBremParticle.vertex())+
668  math::XYZTLorentzVector(theBremParticle.momentum()).Vect().Unit()*ecalShowerDepth;
669 
670  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALShowerMax,
671  meanShower,
672  math::XYZTLorentzVector(theBremParticle.momentum())));}
673  else {
674  if ((dp_tang>5.) && ((dp_tang/sdp_tang)>3))
675  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE ECAL HAS FAILED";
676  PFTrajectoryPoint dummyECAL;
677  brem.addPoint(dummyECAL);
678  PFTrajectoryPoint dummyMaxSh;
679  brem.addPoint(dummyMaxSh);
680  }
681 
682 
683 
684  //HCAL entrance
685  theBremParticle.propagateToHcalEntrance(false);
686  if(theBremParticle.getSuccess()!=0)
687  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALEntrance,
688  math::XYZPoint(theBremParticle.vertex()),
689  math::XYZTLorentzVector(theBremParticle.momentum())));
690  else{
691  if ((dp_tang>5.) && ((dp_tang/sdp_tang)>3))
692  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE HCAL ENTRANCE HAS FAILED";
693  PFTrajectoryPoint dummyHCALentrance;
694  brem.addPoint(dummyHCALentrance);
695  }
696 
697  //HCAL exit
698  theBremParticle.propagateToHcalExit(false);
699  if(theBremParticle.getSuccess()!=0)
700  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALExit,
701  math::XYZPoint(theBremParticle.vertex()),
702  math::XYZTLorentzVector(theBremParticle.momentum())));
703  else{
704  if ((dp_tang>5.) && ((dp_tang/sdp_tang)>3))
705  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE HCAL EXIT HAS FAILED";
706  PFTrajectoryPoint dummyHCALexit;
707  brem.addPoint(dummyHCALexit);
708  }
709 
710  brem.calculatePositionREP();
711  pftrack.addBrem(brem);
712  iTrajPos++;
713 
714 
715 
716 
717  // ***************************** INTERMIDIATE State *************************************
718  //From the new Wolfgang code
719 
720  // To think if the cout should be removed.
721  if(track.gsfExtra()->tangentsSize() == 0)
722  LogError("PFTrackTransformer")
723  <<"BE CAREFUL: Gsf Tangents not stored in the event. You need to re-reco the particle-flow with RecoToDisplay_cfg.py and not RecoToDisplay_NoTracking_cfg.py ";
724 
725 
726  vector<GsfTangent> gsftang = track.gsfExtra()->tangents();
727  for(unsigned int iTang = 0; iTang < track.gsfExtra()->tangentsSize(); iTang++) {
728 
729  dp_tang = gsftang[iTang].deltaP().value();
730  sdp_tang = gsftang[iTang].deltaP().error();
731 
732  //check that the vertex of the brem is in the tracker volume
733  if ((sqrt(gsftang[iTang].position().x()*gsftang[iTang].position().x()
734  + gsftang[iTang].position().y()*gsftang[iTang].position().y())>110)
735  ||(fabs(gsftang[iTang].position().z())>280)) continue;
736 
737  iTrajPoint = iTrajPos + 2;
738  PFBrem brem(dp_tang,sdp_tang,iTrajPoint);
739 
740 
741 
742  GlobalVector p_tang= GlobalVector(gsftang[iTang].momentum().x(),
743  gsftang[iTang].momentum().y(),
744  gsftang[iTang].momentum().z());
745 
746 
747  // ###### track tj points
748  pftrack.addPoint(PFTrajectoryPoint(iid,-1,
749  math::XYZPoint(gsftang[iTang].position().x(),gsftang[iTang].position().y(),gsftang[iTang].position().z()),
750  math::XYZTLorentzVector(p_tang.x(),p_tang.y(),p_tang.z(),p_tang.mag())));
751 
752 
753  //rescale
754  GlobalVector p_gamma = p_tang *(fabs(dp_tang)/p_tang.mag());
755 
756  // GlobalVector
757 
758 
759  double e_gamma = fabs(dp_tang); // DP = pout-pin so could be negative
760  theBremParticle = BaseParticlePropagator(
761  RawParticle(XYZTLorentzVector(p_gamma.x(),
762  p_gamma.y(),
763  p_gamma.z(),
764  e_gamma),
765  XYZTLorentzVector(gsftang[iTang].position().x(),
766  gsftang[iTang].position().y(),
767  gsftang[iTang].position().z(),
768  0.)),
769  0.,0.,B_.z());
770 
771  theBremParticle.setCharge(gamma_charge);
772 
773 
774  PFTrajectoryPoint dummyClosest; // Added just to have the right number order in PFTrack.cc
775  brem.addPoint(dummyClosest);
776 
777 
778  PFTrajectoryPoint dummyBeamPipe; // Added just to have the right number order in PFTrack.cc
779  brem.addPoint(dummyBeamPipe);
780 
781 
782 
783  isBelowPS=false;
784  theBremParticle.propagateToPreshowerLayer1(false);
785  if(theBremParticle.getSuccess()!=0)
786  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS1,
787  math::XYZPoint(theBremParticle.vertex()),
788  math::XYZTLorentzVector(theBremParticle.momentum())));
789  else {
790  PFTrajectoryPoint dummyPS1;
791  brem.addPoint(dummyPS1);
792  }
793 
794  // Brem Entrance PS Layer 2
795 
796  theBremParticle.propagateToPreshowerLayer2(false);
797  if(theBremParticle.getSuccess()!=0){
798  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS2,
799  math::XYZPoint(theBremParticle.vertex()),
800  math::XYZTLorentzVector(theBremParticle.momentum())));
801  isBelowPS=true;
802  } else {
803  PFTrajectoryPoint dummyPS2;
804  brem.addPoint(dummyPS2);
805  }
806 
807  theBremParticle.propagateToEcalEntrance(false);
808 
809  if(theBremParticle.getSuccess()!=0){
810  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALEntrance,
811  math::XYZPoint(theBremParticle.vertex()),
812  math::XYZTLorentzVector(theBremParticle.momentum())));
813 
814  double ecalShowerDepth
815  = PFCluster::getDepthCorrection(theBremParticle.momentum().E(),
816  isBelowPS,
817  false);
818 
819  math::XYZPoint meanShower=math::XYZPoint(theBremParticle.vertex())+
820  math::XYZTLorentzVector(theBremParticle.momentum()).Vect().Unit()*ecalShowerDepth;
821 
822  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALShowerMax,
823  meanShower,
824  math::XYZTLorentzVector(theBremParticle.momentum())));}
825  else {
826  if ((dp_tang>5.) && ((dp_tang/sdp_tang)>3))
827  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE ECAL HAS FAILED";
828  PFTrajectoryPoint dummyECAL;
829  brem.addPoint(dummyECAL);
830  PFTrajectoryPoint dummyMaxSh;
831  brem.addPoint(dummyMaxSh);
832  }
833 
834 
835 
836  //HCAL entrance
837  theBremParticle.propagateToHcalEntrance(false);
838  if(theBremParticle.getSuccess()!=0)
839  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALEntrance,
840  math::XYZPoint(theBremParticle.vertex()),
841  math::XYZTLorentzVector(theBremParticle.momentum())));
842  else{
843  if ((dp_tang>5.) && ((dp_tang/sdp_tang)>3))
844  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE HCAL ENTRANCE HAS FAILED";
845  PFTrajectoryPoint dummyHCALentrance;
846  brem.addPoint(dummyHCALentrance);
847  }
848 
849  //HCAL exit
850  theBremParticle.propagateToHcalExit(false);
851  if(theBremParticle.getSuccess()!=0)
852  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALExit,
853  math::XYZPoint(theBremParticle.vertex()),
854  math::XYZTLorentzVector(theBremParticle.momentum())));
855  else{
856  if ((dp_tang>5.) && ((dp_tang/sdp_tang)>3))
857  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE HCAL EXIT HAS FAILED";
858  PFTrajectoryPoint dummyHCALexit;
859  brem.addPoint(dummyHCALexit);
860  }
861 
862  brem.calculatePositionREP();
863  pftrack.addBrem(brem);
864  iTrajPos++;
865  }
866 
867 
868 
869 
870  // ***************************** OUTER State *************************************
871 
872  if(outTSOS.isValid()) {
873  GlobalVector OutMom;
874  GlobalPoint OutPos;
875 
876  // DANIELE ????? if the out is not valid maybe take the last tangent?
877  // From Wolfgang. It should be always valid
878 
879  mtsMode_->momentumFromModeCartesian(outTSOS,OutMom);
880  mtsMode_->positionFromModeCartesian(outTSOS,OutPos);
881 
882 
883 
884  // last tjpoint
885  pftrack.addPoint(PFTrajectoryPoint(iid,-1,
886  math::XYZPoint(OutPos.x(),OutPos.y(), OutPos.z()),
887  math::XYZTLorentzVector(OutMom.x(),OutMom.y(),OutMom.z(),OutMom.mag())));
888 
889 
890  float ptot_out = sqrt((OutMom.x()*OutMom.x())+(OutMom.y()*OutMom.y())+(OutMom.z()*OutMom.z()));
891  float pfenergy_out =sqrt((pfmass*pfmass)+(ptot_out *ptot_out));
892  BaseParticlePropagator theOutParticle =
894  OutMom.y(),
895  OutMom.z(),
896  pfenergy_out),
897  XYZTLorentzVector(OutPos.x(),
898  OutPos.y(),
899  OutPos.z(),
900  0.)),
901  0.,0.,B_.z());
902  theOutParticle.setCharge(track.charge());
903  isBelowPS=false;
904  theOutParticle.propagateToPreshowerLayer1(false);
905  if(theOutParticle.getSuccess()!=0)
906  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS1,
907  math::XYZPoint(theOutParticle.vertex()),
908  math::XYZTLorentzVector(theOutParticle.momentum())));
909  else {
910  PFTrajectoryPoint dummyPS1;
911  pftrack.addPoint(dummyPS1);
912  }
913 
914 
915  theOutParticle.propagateToPreshowerLayer2(false);
916  if(theOutParticle.getSuccess()!=0){
917  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS2,
918  math::XYZPoint(theOutParticle.vertex()),
919  math::XYZTLorentzVector(theOutParticle.momentum())));
920  isBelowPS=true;
921  } else {
922  PFTrajectoryPoint dummyPS2;
923  pftrack.addPoint(dummyPS2);
924  }
925 
926  theOutParticle.propagateToEcalEntrance(false);
927 
928  if(theOutParticle.getSuccess()!=0){
929  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALEntrance,
930  math::XYZPoint(theOutParticle.vertex()),
931  math::XYZTLorentzVector(theOutParticle.momentum())));
932  double EDepthCorr = 0.01;
933  double ecalShowerDepth
934  = PFCluster::getDepthCorrection(EDepthCorr,
935  isBelowPS,
936  false);
937 
938  math::XYZPoint meanShower=math::XYZPoint(theOutParticle.vertex())+
939  math::XYZTLorentzVector(theOutParticle.momentum()).Vect().Unit()*ecalShowerDepth;
940 
941  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALShowerMax,
942  meanShower,
943  math::XYZTLorentzVector(theOutParticle.momentum())));}
944  else {
945  if (PT>5.)
946  LogWarning("PFTrackTransformer")<<"GSF TRACK "<<pftrack<< " PROPAGATION TO THE ECAL HAS FAILED";
947  PFTrajectoryPoint dummyECAL;
948  pftrack.addPoint(dummyECAL);
949  PFTrajectoryPoint dummyMaxSh;
950  pftrack.addPoint(dummyMaxSh);
951  }
952 
953 
954 
955  //HCAL entrance
956  theOutParticle.propagateToHcalEntrance(false);
957  if(theOutParticle.getSuccess()!=0)
958  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALEntrance,
959  math::XYZPoint(theOutParticle.vertex()),
960  math::XYZTLorentzVector(theOutParticle.momentum())));
961  else{
962  if (PT>5.)
963  LogWarning("PFTrackTransformer")<<"GSF TRACK "<<pftrack<< " PROPAGATION TO THE HCAL ENTRANCE HAS FAILED";
964  PFTrajectoryPoint dummyHCALentrance;
965  pftrack.addPoint(dummyHCALentrance);
966  }
967 
968  //HCAL exit
969  theOutParticle.propagateToHcalExit(false);
970  if(theOutParticle.getSuccess()!=0)
971  pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALExit,
972  math::XYZPoint(theOutParticle.vertex()),
973  math::XYZTLorentzVector(theOutParticle.momentum())));
974  else{
975  if (PT>5.)
976  LogWarning("PFTrackTransformer")<<"GSF TRACK "<<pftrack<< " PROPAGATION TO THE HCAL EXIT HAS FAILED";
977  PFTrajectoryPoint dummyHCALexit;
978  pftrack.addPoint(dummyHCALexit);
979  }
980 
981 
982 
983 
984  //######## Photon at the OUTER State ##########
985 
986  dp_tang = OutMom.mag();
987  // for the moment same inner error just for semplicity
988  sdp_tang = track.ptModeError()*(track.pMode()/track.ptMode());
989  iTrajPoint = iTrajPos + 2;
990  PFBrem brem(dp_tang,sdp_tang,iTrajPoint);
991 
992  theBremParticle =
994  OutMom.y(),
995  OutMom.z(),
996  dp_tang),
997  XYZTLorentzVector(OutPos.x(),
998  OutPos.y(),
999  OutPos.z(),
1000  0.)),
1001  0.,0.,B_.z());
1002  theBremParticle.setCharge(gamma_charge);
1003 
1004 
1005  PFTrajectoryPoint dummyClosest; // Added just to have the right number order in PFTrack.cc
1006  brem.addPoint(dummyClosest);
1007 
1008 
1009  PFTrajectoryPoint dummyBeamPipe; // Added just to have the right number order in PFTrack.cc
1010  brem.addPoint(dummyBeamPipe);
1011 
1012 
1013 
1014  isBelowPS=false;
1015  theBremParticle.propagateToPreshowerLayer1(false);
1016  if(theBremParticle.getSuccess()!=0)
1017  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS1,
1018  math::XYZPoint(theBremParticle.vertex()),
1019  math::XYZTLorentzVector(theBremParticle.momentum())));
1020  else {
1021  PFTrajectoryPoint dummyPS1;
1022  brem.addPoint(dummyPS1);
1023  }
1024 
1025  // Brem Entrance PS Layer 2
1026 
1027  theBremParticle.propagateToPreshowerLayer2(false);
1028  if(theBremParticle.getSuccess()!=0){
1029  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS2,
1030  math::XYZPoint(theBremParticle.vertex()),
1031  math::XYZTLorentzVector(theBremParticle.momentum())));
1032  isBelowPS=true;
1033  } else {
1034  PFTrajectoryPoint dummyPS2;
1035  brem.addPoint(dummyPS2);
1036  }
1037 
1038  theBremParticle.propagateToEcalEntrance(false);
1039 
1040  if(theBremParticle.getSuccess()!=0){
1041  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALEntrance,
1042  math::XYZPoint(theBremParticle.vertex()),
1043  math::XYZTLorentzVector(theBremParticle.momentum())));
1044  double ecalShowerDepth
1045  = PFCluster::getDepthCorrection(theBremParticle.momentum().E(),
1046  isBelowPS,
1047  false);
1048 
1049  math::XYZPoint meanShower=math::XYZPoint(theBremParticle.vertex())+
1050  math::XYZTLorentzVector(theBremParticle.momentum()).Vect().Unit()*ecalShowerDepth;
1051 
1052  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALShowerMax,
1053  meanShower,
1054  math::XYZTLorentzVector(theBremParticle.momentum())));}
1055  else {
1056  if ((dp_tang>5.) && ((dp_tang/sdp_tang)>3))
1057  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE ECAL HAS FAILED";
1058  PFTrajectoryPoint dummyECAL;
1059  brem.addPoint(dummyECAL);
1060  PFTrajectoryPoint dummyMaxSh;
1061  brem.addPoint(dummyMaxSh);
1062  }
1063 
1064 
1065 
1066  //HCAL entrance
1067  theBremParticle.propagateToHcalEntrance(false);
1068  if(theBremParticle.getSuccess()!=0)
1069  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALEntrance,
1070  math::XYZPoint(theBremParticle.vertex()),
1071  math::XYZTLorentzVector(theBremParticle.momentum())));
1072  else{
1073  if ((dp_tang>5.) && ((dp_tang/sdp_tang)>3))
1074  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE HCAL ENTRANCE HAS FAILED";
1075  PFTrajectoryPoint dummyHCALentrance;
1076  brem.addPoint(dummyHCALentrance);
1077  }
1078 
1079  //HCAL exit
1080  theBremParticle.propagateToHcalExit(false);
1081  if(theBremParticle.getSuccess()!=0)
1082  brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALExit,
1083  math::XYZPoint(theBremParticle.vertex()),
1084  math::XYZTLorentzVector(theBremParticle.momentum())));
1085  else{
1086  if ((dp_tang>5.) && ((dp_tang/sdp_tang)>3))
1087  LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE HCAL EXIT HAS FAILED";
1088  PFTrajectoryPoint dummyHCALexit;
1089  brem.addPoint(dummyHCALexit);
1090  }
1091 
1092  brem.calculatePositionREP();
1093  pftrack.addBrem(brem);
1094  iTrajPos++;
1095  }
1096 
1097  return true;
1098 }
#define LogDebug(id)
double p() const
momentum vector magnitude
Definition: TrackBase.h:128
void setCharge(float q)
set the MEASURED charge
Definition: RawParticle.cc:138
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:22
double pyMode() const
y coordinate of momentum vector from mode
Definition: GsfTrack.h:52
bool positionFromModeCartesian(const TrajectoryStateOnSurface tsos, GlobalPoint &position) const
bool propagateToPreshowerLayer1(bool first=true)
T perp() const
Definition: PV3DBase.h:66
double pMode() const
momentum vector magnitude from mode
Definition: GsfTrack.h:46
#define PT
math::XYZVector B_
B field.
bool addPointsAndBrems(reco::GsfPFRecTrack &pftrack, const reco::Track &track, const Trajectory &traj, const bool &GetMode) const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:57
int getSuccess() const
Has propagation been performed and was barrel or endcap reached ?
General CMS geometry parameters used during Particle Flow reconstruction or drawing. All methods and members are static.
Definition: PFGeometry.h:23
#define abs(x)
Definition: mlp_lapack.h:159
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:132
bool momentumFromModeCartesian(const TrajectoryStateOnSurface tsos, GlobalVector &momentum) const
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:47
double sigmafittedDP() const
Definition: PFGsfHelper.cc:140
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
TrajectoryStateOnSurface outerStateOnSurface(const reco::GsfTrack &tk) const
const MultiTrajectoryStateMode * mtsMode_
double ptModeError() const
error on Pt (set to 1000 TeV if charge==0 for safety) from mode
Definition: GsfTrack.h:80
PropagationDirection const & direction() const
Definition: Trajectory.cc:195
void addBrem(const reco::PFBrem &brem)
add a Bremsstrahlung photon
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
static const float outerZ(PFGeometry::Layers_t layer)
return outer position along z axis of a given layer
Definition: PFGeometry.h:65
DataContainer const & measurements() const
Definition: Trajectory.h:169
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:285
T mag() const
Definition: PV3DBase.h:61
Definition: DDAxes.h:10
bool addPoints(reco::PFRecTrack &pftrack, const reco::Track &track, const Trajectory &traj) const
Add points to a PFTrack. return false if a TSOS is invalid.
T sqrt(T t)
Definition: SSEVec.h:28
double pt() const
track transverse momentum
Definition: TrackBase.h:130
bool propagateToHcalExit(bool first=true)
T z() const
Definition: PV3DBase.h:58
double pzMode() const
z coordinate of momentum vector from mode
Definition: GsfTrack.h:54
void addPoint(const reco::PFTrajectoryPoint &trajPt)
Definition: PFTrack.cc:42
bool propagateToEcalEntrance(bool first=true)
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:136
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:155
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:284
const GsfTrackExtraRef & gsfExtra() const
reference to &quot;extra&quot; object
Definition: GsfTrack.h:31
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:49
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
static const float outerRadius(PFGeometry::Layers_t layer)
return outer radius of a given layer
Definition: PFGeometry.h:57
void calculatePositionREP()
Definition: PFTrack.cc:68
double pxMode() const
x coordinate of momentum vector from mode
Definition: GsfTrack.h:50
double fittedDP() const
Definition: PFGsfHelper.cc:136
TrajectoryStateOnSurface innerStateOnSurface(const reco::GsfTrack &tk) const
bool propagateToHcalEntrance(bool first=true)
unsigned int algoType() const
Definition: PFRecTrack.h:47
int charge() const
track electric charge
Definition: TrackBase.h:112
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
GlobalVector computeP(bool ComputeMode) const
Definition: PFGsfHelper.cc:130
T x() const
Definition: PV3DBase.h:56
double ptMode() const
track transverse momentum from mode
Definition: GsfTrack.h:48
bool propagateToPreshowerLayer2(bool first=true)
mathSSE::Vec4< T > v
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:134
Global3DVector GlobalVector
Definition: GlobalVector.h:10
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
PFTrackTransformer(math::XYZVector)