CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BaseParticlePropagator.cc
Go to the documentation of this file.
1 //FAMOS headers
3 
5  RawParticle(), rCyl(0.), rCyl2(0.), zCyl(0.), bField(0), properDecayTime(1E99)
6 {;}
7 
9  const RawParticle& myPart,double R, double Z, double B, double t ) :
10  RawParticle(myPart), rCyl(R), rCyl2(R*R), zCyl(Z), bField(B), properDecayTime(t)
11 {
12  init();
13 }
14 
16  const RawParticle& myPart,double R, double Z, double B) :
17  RawParticle(myPart), rCyl(R), rCyl2(R*R), zCyl(Z), bField(B), properDecayTime(1E99)
18 {
19  init();
20 }
21 
22 void
24 
25  // The status of the propagation
26  success = 0;
27  // Propagate only for the first half-loop
28  firstLoop = true;
29  //
30  fiducial = true;
31  // The particle has not yet decayed
32  decayed = false;
33  // The proper time is set to zero
34  properTime = 0.;
35  // The propagation direction is along the momentum
36  propDir = 1;
37  //
38  debug = false;
39 
40 }
41 
42 bool
44 
45 
46  //
47  // Check that the particle is not already on the cylinder surface
48  //
49  double rPos2 = R2();
50 
51  if ( onBarrel(rPos2) ) {
52  success = 1;
53  return true;
54  }
55  //
56  if ( onEndcap(rPos2) ) {
57  success = 2;
58  return true;
59  }
60  //
61  // Treat neutral particles (or no magnetic field) first
62  //
63  if ( fabs(charge()) < 1E-12 || bField < 1E-12 ) {
64 
65  //
66  // First check that the particle crosses the cylinder
67  //
68  double pT2 = Perp2();
69  // double pT2 = pT*pT;
70  // double b2 = rCyl * pT;
71  double b2 = rCyl2 * pT2;
72  double ac = fabs ( X()*Py() - Y()*Px() );
73  double ac2 = ac*ac;
74  //
75  // The particle never crosses (or never crossed) the cylinder
76  //
77  // if ( ac > b2 ) return false;
78  if ( ac2 > b2 ) return false;
79 
80  // double delta = std::sqrt(b2*b2 - ac*ac);
81  double delta = std::sqrt(b2 - ac2);
82  double tplus = -( X()*Px()+Y()*Py() ) + delta;
83  double tminus = -( X()*Px()+Y()*Py() ) - delta;
84  //
85  // Find the first (time-wise) intersection point with the cylinder
86  //
87 
88  double solution = tminus < 0 ? tplus/pT2 : tminus/pT2;
89  if ( solution < 0. ) return false;
90  //
91  // Check that the particle does (did) not exit from the endcaps first.
92  //
93  double zProp = Z() + Pz() * solution;
94  if ( fabs(zProp) > zCyl ) {
95  tplus = ( zCyl - Z() ) / Pz();
96  tminus = (-zCyl - Z() ) / Pz();
97  solution = tminus < 0 ? tplus : tminus;
98  if ( solution < 0. ) return false;
99  zProp = Z() + Pz() * solution;
100  success = 2;
101  }
102  else {
103  success = 1;
104  }
105 
106  //
107  // Check the decay time
108  //
109  double delTime = propDir * mass() * solution;
110  double factor = 1.;
111  properTime += delTime;
112  if ( properTime > properDecayTime ) {
113  factor = 1.-(properTime-properDecayTime)/delTime;
115  decayed = true;
116  }
117 
118  //
119  // Compute the coordinates of the RawParticle after propagation
120  //
121  double xProp = X() + Px() * solution * factor;
122  double yProp = Y() + Py() * solution * factor;
123  zProp = Z() + Pz() * solution * factor;
124  double tProp = T() + E() * solution * factor;
125 
126  //
127  // Last check : Although propagated to the endcaps, R could still be
128  // larger than rCyl because the cylinder was too short...
129  //
130  if ( success == 2 && xProp*xProp+yProp*yProp > rCyl2 ) {
131  success = -1;
132  return true;
133  }
134 
135  //
136  // ... and update the particle with its propagated coordinates
137  //
138  setVertex( XYZTLorentzVector(xProp,yProp,zProp,tProp) );
139 
140  return true;
141  }
142  //
143  // Treat charged particles with magnetic field next.
144  //
145  else {
146  //
147  // First check that the particle can cross the cylinder
148  //
149  double pT = Pt();
150  double radius = helixRadius(pT);
151  double phi0 = helixStartPhi();
152  double xC = helixCentreX(radius,phi0);
153  double yC = helixCentreY(radius,phi0);
154  double dist = helixCentreDistToAxis(xC,yC);
155  //
156  // The helix either is outside or includes the cylinder -> no intersection
157  //
158  if ( fabs ( fabs(radius) - dist ) > rCyl ) return false;
159  //
160  // The particle is already away from the endcaps, and will never come back
161  // Could be back-propagated, so warn the user that it is so.
162  //
163  if ( Z() * Pz() > zCyl * fabs(Pz()) ) {
164  success = -2;
165  return true;
166  }
167 
168  double pZ = Pz();
169  double phiProp, zProp;
170  //
171  // Does the helix cross the cylinder barrel ? If not, do the first half-loop
172  //
173  double rProp = std::min(fabs(radius)+dist-0.000001, rCyl);
174 
175  // Check for rounding errors in the ArcSin.
176  double sinPhiProp =
177  (rProp*rProp-radius*radius-dist*dist)/( 2.*dist*radius);
178 
179  //
180  double deltaPhi = 1E99;
181 
182  // Taylor development up to third order for large momenta
183  if ( 1.-fabs(sinPhiProp) < 1E-9 ) {
184 
185  double cphi0 = std::cos(phi0);
186  double sphi0 = std::sin(phi0);
187  double r0 = (X()*cphi0 + Y()*sphi0)/radius;
188  double q0 = (X()*sphi0 - Y()*cphi0)/radius;
189  double rcyl2 = (rCyl2 - X()*X() - Y()*Y())/(radius*radius);
190  double delta = r0*r0 + rcyl2*(1.-q0);
191 
192  // This is a solution of a second order equation, assuming phi = phi0 + epsilon
193  // and epsilon is small.
194  deltaPhi = radius > 0 ?
195  ( -r0 + std::sqrt(delta ) ) / (1.-q0) :
196  ( -r0 - std::sqrt(delta ) ) / (1.-q0);
197 
198  }
199 
200  // Use complete calculation otherwise, or if the delta phi is large anyway.
201  if ( fabs(deltaPhi) > 1E-3 ) {
202 
203  // phiProp = fabs(sinPhiProp) < 0.99999999 ?
204  // std::asin(sinPhiProp) : M_PI/2. * sinPhiProp;
205  phiProp = std::asin(sinPhiProp);
206 
207  //
208  // Compute phi after propagation (two solutions, but the asin definition
209  // (between -pi/2 and +pi/2) takes care of the wrong one.
210  //
211  phiProp = helixCentrePhi(xC,yC) +
212  ( inside(rPos2) || onSurface(rPos2) ? phiProp : M_PI-phiProp );
213 
214  //
215  // Solve the obvious two-pi ambiguities: more than one turn!
216  //
217  if ( fabs(phiProp - phi0) > 2.*M_PI )
218  phiProp += ( phiProp > phi0 ? -2.*M_PI : +2.*M_PI );
219 
220  //
221  // Check that the phi difference has the right sign, according to Q*B
222  // (Another two pi ambuiguity, slightly more subtle)
223  //
224  if ( (phiProp-phi0)*radius < 0.0 )
225  radius > 0.0 ? phiProp += 2 * M_PI : phiProp -= 2 * M_PI;
226 
227  } else {
228 
229  // Use Taylor
230  phiProp = phi0 + deltaPhi;
231 
232  }
233 
234  //
235  // Check that the particle does not exit from the endcaps first.
236  //
237  zProp = Z() + ( phiProp - phi0 ) * pZ * radius / pT ;
238  if ( fabs(zProp) > zCyl ) {
239 
240  zProp = zCyl * fabs(pZ)/pZ;
241  phiProp = phi0 + (zProp - Z()) / radius * pT / pZ;
242  success = 2;
243 
244  } else {
245 
246  //
247  // If the particle does not cross the barrel, either process the
248  // first-half loop only or propagate all the way down to the endcaps
249  //
250  if ( rProp < rCyl ) {
251 
252  if ( firstLoop || fabs(pZ)/pT < 1E-10 ) {
253 
254  success = 0;
255 
256  } else {
257 
258  zProp = zCyl * fabs(pZ)/pZ;
259  phiProp = phi0 + (zProp - Z()) / radius * pT / pZ;
260  success = 2;
261 
262  }
263  //
264  // The particle crossed the barrel
265  //
266  } else {
267 
268  success = 1;
269 
270  }
271  }
272  //
273  // Compute the coordinates of the RawParticle after propagation
274  //
275  //
276  // Check the decay time
277  //
278  double delTime = propDir * (phiProp-phi0)*radius*mass() / pT;
279  double factor = 1.;
280  properTime += delTime;
281  if ( properTime > properDecayTime ) {
282  factor = 1.-(properTime-properDecayTime)/delTime;
284  decayed = true;
285  }
286 
287  zProp = Z() + (zProp-Z())*factor;
288 
289  phiProp = phi0 + (phiProp-phi0)*factor;
290 
291  double sProp = std::sin(phiProp);
292  double cProp = std::cos(phiProp);
293  double xProp = xC + radius * sProp;
294  double yProp = yC - radius * cProp;
295  double tProp = T() + (phiProp-phi0)*radius*E()/pT;
296  double pxProp = pT * cProp;
297  double pyProp = pT * sProp;
298 
299  //
300  // Last check : Although propagated to the endcaps, R could still be
301  // larger than rCyl because the cylinder was too short...
302  //
303  if ( success == 2 && xProp*xProp+yProp*yProp > rCyl2 ) {
304  success = -1;
305  return true;
306  }
307  //
308  // ... and update the particle vertex and momentum
309  //
310  setVertex( XYZTLorentzVector(xProp,yProp,zProp,tProp) );
311  SetXYZT(pxProp,pyProp,pZ,E());
312  // SetPx(pxProp);
313  // SetPy(pyProp);
314  return true;
315 
316  }
317 }
318 
319 bool
321 
322  // Backpropagate
323  SetXYZT(-Px(),-Py(),-Pz(),E());
324  setCharge(-charge());
325  propDir = -1;
326  bool done = propagate();
327  SetXYZT(-Px(),-Py(),-Pz(),E());
328  setCharge(-charge());
329  propDir = +1;
330 
331  return done;
332 }
333 
336  // Copy the input particle in a new instance
337  BaseParticlePropagator myPropPart(*this);
338  // Allow for many loops
339  myPropPart.firstLoop=false;
340  // Propagate the particle !
341  myPropPart.propagate();
342  // Back propagate if the forward propagation was not a success
343  if (myPropPart.success < 0 )
344  myPropPart.backPropagate();
345  // Return the new instance
346  return myPropPart;
347 }
348 
349 void
351  rCyl = R;
352  rCyl2 = R*R;
353  zCyl = Z;
354  firstLoop = f;
355 }
356 
357 bool
359 
360  // Pre-computed quantities
361  double pT = Pt();
362  double radius = helixRadius(pT);
363  double phi0 = helixStartPhi();
364 
365  // The Helix centre coordinates are computed wrt to the z axis
366  // Actually the z axis might be centred on 0,0: it's the beam spot position (x0,y0)!
367  double xC = helixCentreX(radius,phi0);
368  double yC = helixCentreY(radius,phi0);
369  double distz = helixCentreDistToAxis(xC-x0,yC-y0);
370  double dist0 = helixCentreDistToAxis(xC,yC);
371 
372  //
373  // Propagation to point of clostest approach to z axis
374  //
375  double rz,r0,z;
376  if ( charge() != 0.0 && bField != 0.0 ) {
377  rz = fabs ( fabs(radius) - distz ) + std::sqrt(x0*x0+y0*y0) + 0.0000001;
378  r0 = fabs ( fabs(radius) - dist0 ) + 0.0000001;
379  } else {
380  rz = fabs( Px() * (Y()-y0) - Py() * (X()-x0) ) / Pt();
381  r0 = fabs( Px() * Y() - Py() * X() ) / Pt();
382  }
383 
384  z = 999999.;
385 
386  // Propagate to the first interesection point
387  // with cylinder of radius sqrt(x0**2+y0**2)
388  setPropagationConditions(rz , z, first);
389  bool done = backPropagate();
390 
391  // Unsuccessful propagation - should not happen!
392  if ( !done ) return done;
393 
394  // The z axis is (0,0) - no need to go further
395  if ( fabs(rz-r0) < 1E-10 ) return done;
396  double dist1 = (X()-x0)*(X()-x0) + (Y()-y0)*(Y()-y0);
397 
398  // We are already at closest approach - no need to go further
399  if ( dist1 < 1E-10 ) return done;
400 
401  // Keep for later if it happens to be the right solution
402  XYZTLorentzVector vertex1 = vertex();
403  XYZTLorentzVector momentum1 = momentum();
404 
405  // Propagate to the distance of closest approach to (0,0)
406  setPropagationConditions(r0 , z, first);
407  done = backPropagate();
408  if ( !done ) return done;
409 
410  // Propagate to the first interesection point
411  // with cylinder of radius sqrt(x0**2+y0**2)
412  setPropagationConditions(rz , z, first);
413  done = backPropagate();
414  if ( !done ) return done;
415  double dist2 = (X()-x0)*(X()-x0) + (Y()-y0)*(Y()-y0);
416 
417  // Keep the good solution.
418  if ( dist2 > dist1 ) {
419  setVertex(vertex1);
420  SetXYZT(momentum1.X(),momentum1.Y(),momentum1.Z(),momentum1.E());
421  dist2 = dist1;
422  }
423 
424  // Done
425  return done;
426 
427 }
428 
429 bool
431  //
432  // Propagation to Ecal (including preshower) after the
433  // last Tracker layer
434  // TODO: include proper geometry
435  // Geometry taken from CMS ECAL TDR
436  //
437  // First propagate to global barrel / endcap cylinder
438  // setPropagationConditions(1290. , 3045 , first);
439  // New position of the preshower as of CMSSW_1_7_X
440  setPropagationConditions(129.0 , 303.353 , first);
441  return propagate();
442 
443 }
444 
445 bool
447  //
448  // Propagation to Preshower Layer 1
449  // TODO: include proper geometry
450  // Geometry taken from CMS ECAL TDR
451  //
452  // First propagate to global barrel / endcap cylinder
453  // setPropagationConditions(1290., 3045 , first);
454  // New position of the preshower as of CMSSW_1_7_X
455  setPropagationConditions(129.0, 303.353 , first);
456  bool done = propagate();
457 
458  // Check that were are on the Layer 1
459  if ( done && (R2() > 125.0*125.0 || R2() < 45.0*45.0) )
460  success = 0;
461 
462  return done;
463 }
464 
465 bool
467  //
468  // Propagation to Preshower Layer 2
469  // TODO: include proper geometry
470  // Geometry taken from CMS ECAL TDR
471  //
472  // First propagate to global barrel / endcap cylinder
473  // setPropagationConditions(1290. , 3090 , first);
474  // New position of the preshower as of CMSSW_1_7_X
475  setPropagationConditions(129.0 , 307.838 , first);
476  bool done = propagate();
477 
478  // Check that we are on Layer 2
479  if ( done && (R2() > 125.0*125.0 || R2() < 45.0*45.0 ) )
480  success = 0;
481 
482  return done;
483 
484 }
485 
486 bool
488  //
489  // Propagation to ECAL entrance
490  // TODO: include proper geometry
491  // Geometry taken from CMS ECAL TDR
492  //
493  // First propagate to global barrel / endcap cylinder
494  setPropagationConditions(129.0 , 320.9,first);
495  bool done = propagate();
496 
497  // Go to endcap cylinder in the "barrel cut corner"
498  // eta = 1.479 -> cos^2(theta) = 0.81230
499  // if ( done && eta > 1.479 && success == 1 ) {
500  if ( done && cos2ThetaV() > 0.81230 && success == 1 ) {
501  setPropagationConditions(152.6 , 320.9, first);
502  done = propagate();
503  }
504 
505  // We are not in the ECAL acceptance
506  // eta = 3.0 -> cos^2(theta) = 0.99013
507  if ( cos2ThetaV() > 0.99014 ) success = 0;
508 
509  return done;
510 }
511 
512 bool
514  //
515  // Propagation to HCAL entrance
516  // TODO: include proper geometry
517  // Geometry taken from CMSSW_3_1_X xml (Sunanda)
518  //
519 
520  // First propagate to global barrel / endcap cylinder
521  setPropagationConditions(177.5 , 335.0, first);
522  propDir = 0;
523  bool done = propagate();
524  propDir = 1;
525 
526  // If went through the bottom of HB cylinder -> re-propagate to HE surface
527  if (done && success == 2) {
528  setPropagationConditions(300.0, 400.458, first);
529  propDir = 0;
530  done = propagate();
531  propDir = 1;
532  }
533 
534 
535  // out of the HB/HE acceptance
536  // eta = 3.0 -> cos^2(theta) = 0.99014
537  if ( done && cos2ThetaV() > 0.99014 ) success = 0;
538 
539  return done;
540 }
541 
542 bool
544  //
545  // Propagation to VFCAL entrance
546  // TODO: include proper geometry
547  // Geometry taken from DAQ TDR Chapter 13
548 
549  setPropagationConditions(400.0 , 1110.0, first);
550  propDir = 0;
551  bool done = propagate();
552  propDir = 1;
553 
554  if (!done) success = 0;
555 
556  // We are not in the VFCAL acceptance
557  // eta = 3.0 -> cos^2(theta) = 0.99014
558  // eta = 5.2 -> cos^2(theta) = 0.9998755
559  double c2teta = cos2ThetaV();
560  if ( done && ( c2teta < 0.99014 || c2teta > 0.9998755 ) ) success = 0;
561 
562  return done;
563 }
564 
565 bool
567  //
568  // Propagation to HCAL exit
569  // TODO: include proper geometry
570  // Geometry taken from CMSSW_3_1_X xml (Sunanda)
571  //
572 
573  // Approximate it to a single cylinder as it is not that crucial.
574  setPropagationConditions(263.9 , 554.1, first);
575  // this->rawPart().setCharge(0.0); ?? Shower Propagation ??
576  propDir = 0;
577  bool done = propagate();
578  propDir = 1;
579 
580  return done;
581 }
582 
583 bool
585  //
586  // Propagation to Layer0 of HO (Ring-0)
587  // TODO: include proper geometry
588  // Geometry taken from CMSSW_3_1_X xml (Sunanda)
589  //
590 
591  // Approximate it to a single cylinder as it is not that crucial.
592  setPropagationConditions(387.6, 1000.0, first); //Dia is the middle of HO \sqrt{384.8**2+46.2**2}
593  // setPropagationConditions(410.2, 1000.0, first); //sqrt{406.6**2+54.2**2} //for outer layer
594  // this->rawPart().setCharge(0.0); ?? Shower Propagation ??
595  propDir = 0;
596  bool done = propagate();
597  propDir = 1;
598 
599  if ( done && abs(Z()) > 700.25) success = 0;
600 
601  return done;
602 }
603 
604 
605 bool
607 {
608 
609 
610  // Not implemented for neutrals (used for electrons only)
611  if ( charge() == 0. || bField == 0.) return false;
612 
613  // Define the proper pT direction to meet the point (vx,vy) in (x,y)
614  double dx = X()-v.X();
615  double dy = Y()-v.Y();
616  double phi = std::atan2(dy,dx) + std::asin ( std::sqrt(dx*dx+dy*dy)/(2.*helixRadius()) );
617 
618  // The absolute pT (kept to the original value)
619  double pT = pt();
620 
621  // Set the new pT
622  SetPx(pT*std::cos(phi));
623  SetPy(pT*std::sin(phi));
624 
625  return propagateToClosestApproach(v.X(),v.Y());
626 
627 }
628 
629 bool
631 {
632 
633  // For neutral (or BField = 0, simply check that the track passes through the cylinder
634  if ( charge() == 0. || bField == 0.)
635  return fabs( Px() * Y() - Py() * X() ) / Pt() < radius;
636 
637  // Now go to the charged particles
638 
639  // A few needed intermediate variables (to make the code understandable)
640  // The inner cylinder
641  double r = radius;
642  double r2 = r*r;
643  double r4 = r2*r2;
644  // The two hits
645  double dx = X()-v.X();
646  double dy = Y()-v.Y();
647  double dz = Z()-v.Z();
648  double Sxy = X()*v.X() + Y()*v.Y();
649  double Dxy = Y()*v.X() - X()*v.Y();
650  double Dxy2 = Dxy*Dxy;
651  double Sxy2 = Sxy*Sxy;
652  double SDxy = dx*dx + dy*dy;
653  double SSDxy = std::sqrt(SDxy);
654  double ATxy = std::atan2(dy,dx);
655  // The coefficients of the second order equation for the trajectory radius
656  double a = r2 - Dxy2/SDxy;
657  double b = r * (r2 - Sxy);
658  double c = r4 - 2.*Sxy*r2 + Sxy2 + Dxy2;
659 
660  // Here are the four possible solutions for
661  // 1) The trajectory radius
662  std::vector<double> helixRs(4,static_cast<double>(0.));
663  helixRs[0] = (b - std::sqrt(b*b - a*c))/(2.*a);
664  helixRs[1] = (b + std::sqrt(b*b - a*c))/(2.*a);
665  helixRs[2] = -helixRs[0];
666  helixRs[3] = -helixRs[1];
667  // 2) The azimuthal direction at the second point
668  std::vector<double> helixPhis(4,static_cast<double>(0.));
669  helixPhis[0] = std::asin ( SSDxy/(2.*helixRs[0]) );
670  helixPhis[1] = std::asin ( SSDxy/(2.*helixRs[1]) );
671  helixPhis[2] = -helixPhis[0];
672  helixPhis[3] = -helixPhis[1];
673  // Only two solutions are valid though
674  double solution1=0.;
675  double solution2=0.;
676  double phi1=0.;
677  double phi2=0.;
678  // Loop over the four possibilities
679  for ( unsigned int i=0; i<4; ++i ) {
680  helixPhis[i] = ATxy + helixPhis[i];
681  // The centre of the helix
682  double xC = helixCentreX(helixRs[i],helixPhis[i]);
683  double yC = helixCentreY(helixRs[i],helixPhis[i]);
684  // The distance of that centre to the origin
685  double dist = helixCentreDistToAxis(xC, yC);
686  /*
687  std::cout << "Solution "<< i << " = " << helixRs[i] << " accepted ? "
688  << fabs(fabs(helixRs[i]) - dist ) << " - " << radius << " = "
689  << fabs(fabs(helixRs[i]) - dist ) - fabs(radius) << " "
690  << ( fabs( fabs(fabs(helixRs[i]) - dist) -fabs(radius) ) < 1e-6) << std::endl;
691  */
692  // Check that the propagation will lead to a trajectroy tangent to the inner cylinder
693  if ( fabs( fabs(fabs(helixRs[i]) - dist) -fabs(radius) ) < 1e-6 ) {
694  // Fill first solution
695  if ( solution1 == 0. ) {
696  solution1 = helixRs[i];
697  phi1 = helixPhis[i];
698  // Fill second solution (ordered)
699  } else if ( solution2 == 0. ) {
700  if ( helixRs[i] < solution1 ) {
701  solution2 = solution1;
702  solution1 = helixRs[i];
703  phi2 = phi1;
704  phi1 = helixPhis[i];
705  } else {
706  solution2 = helixRs[i];
707  phi2 = helixPhis[i];
708  }
709  // Must not happen!
710  } else {
711  std::cout << "warning !!! More than two solutions for this track !!! " << std::endl;
712  }
713  }
714  }
715 
716  // Find the largest possible pT compatible with coming from within the inne cylinder
717  double pT = 0.;
718  double helixPhi = 0.;
719  double helixR = 0.;
720  // This should not happen
721  if ( solution1*solution2 == 0. ) {
722  std::cout << "warning !!! Less than two solution for this track! "
723  << solution1 << " " << solution2 << std::endl;
724  return false;
725  // If the two solutions have opposite sign, it means that pT = infinity is a possibility.
726  } else if ( solution1*solution2 < 0. ) {
727  pT = 1000.;
728  double norm = pT/SSDxy;
729  setCharge(+1.);
730  SetXYZT(dx*norm,dy*norm,dz*norm,0.);
731  SetE(std::sqrt(Vect().Mag2()));
732  // Otherwise take the solution that gives the largest transverse momentum
733  } else {
734  if (solution1<0.) {
735  helixR = solution1;
736  helixPhi = phi1;
737  setCharge(+1.);
738  } else {
739  helixR = solution2;
740  helixPhi = phi2;
741  setCharge(-1.);
742  }
743  pT = fabs(helixR) * 1e-5 * c_light() *bField;
744  double norm = pT/SSDxy;
745  SetXYZT(pT*std::cos(helixPhi),pT*std::sin(helixPhi),dz*norm,0.);
746  SetE(std::sqrt(Vect().Mag2()));
747  }
748 
749  // Propagate to closest approach to get the Z value (a bit of an overkill)
751 
752 }
753 
754 double
755 BaseParticlePropagator::xyImpactParameter(double x0, double y0) const {
756 
757  double ip=0.;
758  double pT = Pt();
759 
760  if ( charge() != 0.0 && bField != 0.0 ) {
761  double radius = helixRadius(pT);
762  double phi0 = helixStartPhi();
763 
764  // The Helix centre coordinates are computed wrt to the z axis
765  double xC = helixCentreX(radius,phi0);
766  double yC = helixCentreY(radius,phi0);
767  double distz = helixCentreDistToAxis(xC-x0,yC-y0);
768  ip = distz - fabs(radius);
769  } else {
770  ip = fabs( Px() * (Y()-y0) - Py() * (X()-x0) ) / pT;
771  }
772 
773  return ip;
774 
775 }
const double Z[kNumberCalorimeter]
dbl * delta
Definition: mlp_gen.cc:36
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
void setCharge(float q)
set the MEASURED charge
Definition: RawParticle.cc:139
int i
Definition: DBlmapReader.cc:9
double helixRadius() const
The helix Radius.
double helixCentrePhi() const
The azimuth if the vector joining the cylinder and the helix axes.
bool propagateToPreshowerLayer1(bool first=true)
double rCyl
Simulated particle that is to be resp has been propagated.
bool propagateToNominalVertex(const XYZTLorentzVector &hit2=XYZTLorentzVector(0., 0., 0., 0.))
void init()
Initialize internal switches and quantities.
double z() const
z of vertex
Definition: RawParticle.h:272
bool firstLoop
Do only the first half-loop.
int success
0:propagation still be done, 1:reached &#39;barrel&#39;, 2:reached &#39;endcaps&#39;
double properTime
The proper time of the particle.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
BaseParticlePropagator()
Default c&#39;tor.
bool propagateToBeamCylinder(const XYZTLorentzVector &v, double radius=0.)
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)
double mass() const
get the MEASURED mass
Definition: RawParticle.h:283
double c_light() const
The speed of light in mm/ns (!) without clhep (yeaaahhh!)
double R() const
vertex radius
Definition: RawParticle.h:278
double R2() const
vertex radius**2
Definition: RawParticle.h:279
bool propagateToEcal(bool first=true)
bool propagateToClosestApproach(double x0=0., double y0=0, bool first=true)
double helixCentreX() const
The x coordinate of the helix axis.
bool propagateToVFcalEntrance(bool first=true)
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:286
double helixCentreDistToAxis() const
The distance between the cylinder and the helix axes.
double r2() const
vertex radius**2
Definition: RawParticle.h:281
double properDecayTime
The proper decay time of the particle.
T sqrt(T t)
Definition: SSEVec.h:48
bool propagateToHcalExit(bool first=true)
double cos2ThetaV() const
Definition: RawParticle.h:269
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double Y() const
y of vertex
Definition: RawParticle.h:275
double Z() const
z of vertex
Definition: RawParticle.h:276
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
double charge() const
get the MEASURED charge
Definition: RawParticle.h:282
#define M_PI
int propDir
The propagation direction.
bool propagateToEcalEntrance(bool first=true)
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:285
bool onEndcap() const
Is the vertex already on the cylinder endcap ?
bool decayed
The particle decayed while propagated !
BaseParticlePropagator propagated() const
double helixCentreY() const
The y coordinate of the helix axis.
bool inside() const
Is the vertex inside the cylinder ? (stricly inside : true)
double b
Definition: hdecay.h:120
bool onSurface() const
Is the vertex already on the cylinder surface ?
double r() const
vertex radius
Definition: RawParticle.h:280
double X() const
x of vertex
Definition: RawParticle.h:274
bool debug
The debug level.
double helixStartPhi() const
The azimuth of the momentum at the vertex.
bool propagateToHcalEntrance(bool first=true)
double a
Definition: hdecay.h:121
tuple cout
Definition: gather_cfg.py:121
double bField
Magnetic field in the cylinder, oriented along the Z axis.
double T() const
vertex time
Definition: RawParticle.h:277
bool propagateToHOLayer(bool first=true)
bool onBarrel() const
Is the vertex already on the cylinder barrel ?
bool propagateToPreshowerLayer2(bool first=true)
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
Definition: RawParticle.h:288
double xyImpactParameter(double x0=0., double y0=0.) const
Transverse impact parameter.
bool fiducial
The particle traverses some real material.
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15