CMS 3D CMS Logo

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