8 : particle_(myPart), rCyl(
R), rCyl2(
R *
R), zCyl(
Z),
bField(
B), properDecayTime(
t) {
13 : particle_(myPart), rCyl(
R), rCyl2(
R *
R), zCyl(
Z),
bField(
B), properDecayTime(1E99) {
77 double solution = tminus < 0 ? tplus /
pT2 : tminus /
pT2;
84 if (fabs(zProp) >
zCyl) {
87 solution = tminus < 0 ? tplus : tminus;
119 if (
success == 2 && xProp * xProp + yProp * yProp >
rCyl2) {
159 double phiProp, zProp;
166 double sinPhiProp = (rProp * rProp -
radius *
radius - dist * dist) / (2. * dist *
radius);
172 if (1. - fabs(sinPhiProp) < 1E-9) {
178 double delta = r0 * r0 + rcyl2 * (1. -
q0);
189 phiProp = std::asin(sinPhiProp);
200 if (fabs(phiProp - phi0) > 2. *
M_PI)
201 phiProp += (phiProp > phi0 ? -2. *
M_PI : +2. *
M_PI);
207 if ((phiProp - phi0) *
radius < 0.0)
219 if (fabs(zProp) >
zCyl) {
220 zProp =
zCyl * fabs(pZ) / pZ;
234 zProp =
zCyl * fabs(pZ) / pZ;
262 phiProp = phi0 + (phiProp - phi0) *
factor;
266 double xProp = xC +
radius * sProp;
267 double yProp = yC -
radius * cProp;
269 double pxProp =
pT * cProp;
270 double pyProp =
pT * sProp;
276 if (
success == 2 && xProp * xProp + yProp * yProp >
rCyl2) {
343 rz = fabs(fabs(
radius) - distz) +
std::sqrt(x0 * x0 + y0 * y0) + 0.0000001;
344 r0 = fabs(fabs(
radius) - dist0) + 0.0000001;
362 if (fabs(rz - r0) < 1E-10)
522 if (
done && (c2teta < 0.99014 || c2teta > 0.9998755))
604 double Dxy2 = Dxy * Dxy;
605 double Sxy2 = Sxy * Sxy;
608 double ATxy = std::atan2(
dy,
dx);
610 double a =
r2 - Dxy2 / SDxy;
611 double b =
r * (
r2 - Sxy);
612 double c = r4 - 2. * Sxy *
r2 + Sxy2 + Dxy2;
616 std::array<double, 4> helixRs = {{0.}};
619 helixRs[2] = -helixRs[0];
620 helixRs[3] = -helixRs[1];
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];
628 double solution1 = 0.;
629 double solution2 = 0.;
633 for (
unsigned int i = 0;
i < 4; ++
i) {
634 helixPhis[
i] = ATxy + helixPhis[
i];
647 if (fabs(fabs(fabs(helixRs[
i]) - dist) - fabs(
radius)) < 1
e-6) {
649 if (solution1 == 0.) {
650 solution1 = helixRs[
i];
653 }
else if (solution2 == 0.) {
654 if (helixRs[
i] < solution1) {
655 solution2 = solution1;
656 solution1 = helixRs[
i];
660 solution2 = helixRs[
i];
665 std::cout <<
"warning !!! More than two solutions for this track !!! " << std::endl;
672 double helixPhi = 0.;
675 if (solution1 * solution2 == 0.) {
676 std::cout <<
"warning !!! Less than two solution for this track! " << solution1 <<
" " << solution2 << std::endl;
679 }
else if (solution1 * solution2 < 0.) {
681 double norm =
pT / SSDxy;
687 if (solution1 < 0.) {
697 double norm =
pT / SSDxy;
718 ip = distz - fabs(
radius);
double R2() const
vertex radius**2
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
void setMomentum(const XYZTLorentzVector &vtx)
set the momentum
void setCharge(float q)
set the MEASURED charge
double helixCentrePhi() const
The azimuth if the vector joining the cylinder and the helix axes.
double c_light() const
The speed of light in mm/ns (!) without clhep (yeaaahhh!)
bool propagateToPreshowerLayer1(bool first=true)
double Pz() const
z of the momentum
double rCyl
Simulated particle that is to be resp has been propagated.
bool propagateToNominalVertex(const XYZTLorentzVector &hit2=XYZTLorentzVector(0., 0., 0., 0.))
double helixStartPhi() const
The azimuth of the momentum at the vertex.
void init()
Initialize internal switches and quantities.
bool firstLoop
Do only the first half-loop.
double Z() const
z of vertex
int success
0:propagation still be done, 1:reached 'barrel', 2:reached 'endcaps'
double properTime
The proper time of the particle.
Sin< T >::type sin(const T &t)
BaseParticlePropagator()
Default c'tor.
double xyImpactParameter(double x0=0., double y0=0.) const
Transverse impact parameter.
BaseParticlePropagator propagated() const
double helixRadius() const
The helix Radius.
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 E() const
energy of the momentum
const XYZTLorentzVector & momentum() const
the momentum fourvector
double charge() const
get the MEASURED charge
double pt() const
transverse momentum
bool propagateToEcal(bool first=true)
double helixCentreDistToAxis() const
The distance between the cylinder and the helix axes.
double Pt() const
transverse momentum
bool propagateToClosestApproach(double x0=0., double y0=0, bool first=true)
bool propagateToVFcalEntrance(bool first=true)
double cos2ThetaV() const
double properDecayTime
The proper decay time of the particle.
double Py() const
y of the momentum
bool propagateToHcalExit(bool first=true)
Cos< T >::type cos(const T &t)
Abs< T >::type abs(const T &t)
bool onEndcap() const
Is the vertex already on the cylinder endcap ?
double helixCentreX() const
The x coordinate of the helix axis.
int propDir
The propagation direction.
bool propagateToEcalEntrance(bool first=true)
double Perp2() const
perpendicular momentum squared
bool decayed
The particle decayed while propagated !
double T() const
vertex time
bool onSurface() const
Is the vertex already on the cylinder surface ?
double mass() const
get the MEASURED mass
bool onBarrel() const
Is the vertex already on the cylinder barrel ?
bool debug
The debug level.
double Y() const
y of vertex
bool propagateToHcalEntrance(bool first=true)
double X() const
x of vertex
double bField
Magnetic field in the cylinder, oriented along the Z axis.
bool propagateToHOLayer(bool first=true)
bool propagateToPreshowerLayer2(bool first=true)
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
bool inside() const
Is the vertex inside the cylinder ? (stricly inside : true)
bool fiducial
The particle traverses some real material.
double helixCentreY() const
The y coordinate of the helix axis.
math::XYZTLorentzVector XYZTLorentzVector
const XYZTLorentzVector & vertex() const
the vertex fourvector
double Px() const
x of the momentum