CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/FastSimulation/BaseParticlePropagator/src/BaseParticlePropagator.cc

Go to the documentation of this file.
00001 //FAMOS headers
00002 #include "FastSimulation/BaseParticlePropagator/interface/BaseParticlePropagator.h"
00003 
00004 BaseParticlePropagator::BaseParticlePropagator() :
00005   RawParticle(), rCyl(0.), rCyl2(0.), zCyl(0.), bField(0), properDecayTime(1E99)
00006 {;}
00007 
00008 BaseParticlePropagator::BaseParticlePropagator( 
00009          const RawParticle& myPart,double R, double Z, double B, double t ) :
00010   RawParticle(myPart), rCyl(R), rCyl2(R*R), zCyl(Z), bField(B), properDecayTime(t) 
00011 {
00012   init();
00013 }
00014 
00015 BaseParticlePropagator::BaseParticlePropagator( 
00016          const RawParticle& myPart,double R, double Z, double B) :
00017   RawParticle(myPart), rCyl(R), rCyl2(R*R), zCyl(Z), bField(B), properDecayTime(1E99) 
00018 {
00019   init();
00020 }
00021 
00022 void 
00023 BaseParticlePropagator::init() {
00024   
00025   // The status of the propagation
00026   success = 0;
00027   // Propagate only for the first half-loop
00028   firstLoop = true;
00029   //
00030   fiducial = true;
00031   // The particle has not yet decayed
00032   decayed = false;
00033   // The proper time is set to zero
00034   properTime = 0.;
00035   // The propagation direction is along the momentum
00036   propDir = 1;
00037   //
00038   debug = false;
00039 
00040 }
00041 
00042 bool 
00043 BaseParticlePropagator::propagate() { 
00044 
00045 
00046   //
00047   // Check that the particle is not already on the cylinder surface
00048   //
00049   double rPos2 = R2();
00050 
00051   if ( onBarrel(rPos2) ) { 
00052     success = 1;
00053     return true;
00054   }
00055   //
00056   if ( onEndcap(rPos2) ) { 
00057     success = 2;
00058     return true;
00059   }
00060   //
00061   // Treat neutral particles (or no magnetic field) first 
00062   //
00063   if ( fabs(charge()) < 1E-12 || bField < 1E-12 ) {
00064 
00065     //
00066     // First check that the particle crosses the cylinder
00067     //
00068     double pT2 = Perp2();
00069     //    double pT2 = pT*pT;
00070     //    double b2 = rCyl * pT;
00071     double b2 = rCyl2 * pT2;
00072     double ac = fabs ( X()*Py() - Y()*Px() );
00073     double ac2 = ac*ac;
00074     //
00075     // The particle never crosses (or never crossed) the cylinder
00076     //
00077     //    if ( ac > b2 ) return false;
00078     if ( ac2 > b2 ) return false;
00079 
00080     //    double delta  = std::sqrt(b2*b2 - ac*ac);
00081     double delta  = std::sqrt(b2 - ac2);
00082     double tplus  = -( X()*Px()+Y()*Py() ) + delta;
00083     double tminus = -( X()*Px()+Y()*Py() ) - delta;
00084     //
00085     // Find the first (time-wise) intersection point with the cylinder
00086     //
00087 
00088     double solution = tminus < 0 ? tplus/pT2 : tminus/pT2;
00089     if ( solution < 0. ) return false;
00090     //
00091     // Check that the particle does (did) not exit from the endcaps first.
00092     //
00093     double zProp = Z() + Pz() * solution;
00094     if ( fabs(zProp) > zCyl ) {
00095       tplus  = ( zCyl - Z() ) / Pz();
00096       tminus = (-zCyl - Z() ) / Pz();
00097       solution = tminus < 0 ? tplus : tminus;
00098       if ( solution < 0. ) return false;
00099       zProp = Z() + Pz() * solution;
00100       success = 2;
00101     }
00102     else {
00103       success = 1;
00104     }
00105 
00106     //
00107     // Check the decay time
00108     //
00109     double delTime = propDir * mass() * solution;
00110     double factor = 1.;
00111     properTime += delTime;
00112     if ( properTime > properDecayTime ) {
00113       factor = 1.-(properTime-properDecayTime)/delTime;
00114       properTime = properDecayTime;
00115       decayed = true;
00116     }
00117 
00118     //
00119     // Compute the coordinates of the RawParticle after propagation
00120     //
00121     double xProp = X() + Px() * solution * factor;
00122     double yProp = Y() + Py() * solution * factor;
00123            zProp = Z() + Pz() * solution * factor;
00124     double tProp = T() + E()  * solution * factor;
00125 
00126     //
00127     // Last check : Although propagated to the endcaps, R could still be
00128     // larger than rCyl because the cylinder was too short...
00129     //
00130     if ( success == 2 && xProp*xProp+yProp*yProp > rCyl2 ) { 
00131       success = -1;
00132       return true;
00133     }
00134 
00135     //
00136     // ... and update the particle with its propagated coordinates
00137     //
00138     setVertex( XYZTLorentzVector(xProp,yProp,zProp,tProp) );
00139     
00140     return true;
00141   }
00142   //
00143   // Treat charged particles with magnetic field next.
00144   //
00145   else {
00146     //
00147     // First check that the particle can cross the cylinder
00148     //
00149     double pT = Pt();
00150     double radius = helixRadius(pT);
00151     double phi0 = helixStartPhi();
00152     double xC = helixCentreX(radius,phi0);
00153     double yC = helixCentreY(radius,phi0);
00154     double dist = helixCentreDistToAxis(xC,yC);
00155     //
00156     // The helix either is outside or includes the cylinder -> no intersection
00157     //
00158     if ( fabs ( fabs(radius) - dist ) > rCyl ) return false;
00159     //
00160     // The particle is already away from the endcaps, and will never come back
00161     // Could be back-propagated, so warn the user that it is so.
00162     // 
00163     if ( Z() * Pz() > zCyl * fabs(Pz()) ) { 
00164       success = -2;
00165       return true;
00166     }
00167 
00168     double pZ = Pz();
00169     double phiProp, zProp;
00170     //
00171     // Does the helix cross the cylinder barrel ? If not, do the first half-loop
00172     //
00173     double rProp = std::min(fabs(radius)+dist-0.000001, rCyl);
00174 
00175     // Check for rounding errors in the ArcSin.
00176     double sinPhiProp = 
00177       (rProp*rProp-radius*radius-dist*dist)/( 2.*dist*radius);
00178     
00179     //
00180     double deltaPhi = 1E99;
00181 
00182     // Taylor development up to third order for large momenta 
00183     if ( 1.-fabs(sinPhiProp) < 1E-9 ) { 
00184 
00185       double cphi0 = std::cos(phi0);
00186       double sphi0 = std::sin(phi0);
00187       double r0 = (X()*cphi0 + Y()*sphi0)/radius;
00188       double q0 = (X()*sphi0 - Y()*cphi0)/radius;
00189       double rcyl2 = (rCyl2 - X()*X() - Y()*Y())/(radius*radius);
00190       double delta = r0*r0 + rcyl2*(1.-q0);
00191 
00192       // This is a solution of a second order equation, assuming phi = phi0 + epsilon
00193       // and epsilon is small. 
00194       deltaPhi = radius > 0 ? 
00195         ( -r0 + std::sqrt(delta ) ) / (1.-q0) :   
00196         ( -r0 - std::sqrt(delta ) ) / (1.-q0); 
00197       
00198     }
00199 
00200     // Use complete calculation otherwise, or if the delta phi is large anyway.
00201     if ( fabs(deltaPhi) > 1E-3 ) { 
00202 
00203     //    phiProp =  fabs(sinPhiProp) < 0.99999999  ? 
00204     //      std::asin(sinPhiProp) : M_PI/2. * sinPhiProp;
00205       phiProp =  std::asin(sinPhiProp);
00206 
00207       //
00208       // Compute phi after propagation (two solutions, but the asin definition
00209       // (between -pi/2 and +pi/2) takes care of the wrong one.
00210       //
00211       phiProp = helixCentrePhi(xC,yC) + 
00212         ( inside(rPos2) || onSurface(rPos2) ? phiProp : M_PI-phiProp ); 
00213       
00214       //
00215       // Solve the obvious two-pi ambiguities: more than one turn!
00216       //
00217       if ( fabs(phiProp - phi0) > 2.*M_PI )
00218         phiProp += ( phiProp > phi0 ? -2.*M_PI : +2.*M_PI );
00219       
00220       //
00221       // Check that the phi difference has the right sign, according to Q*B
00222       // (Another two pi ambuiguity, slightly more subtle)
00223       //
00224       if ( (phiProp-phi0)*radius < 0.0 )
00225         radius > 0.0 ? phiProp += 2 * M_PI : phiProp -= 2 * M_PI;
00226       
00227     } else { 
00228 
00229       // Use Taylor
00230       phiProp = phi0 + deltaPhi;
00231 
00232     }
00233 
00234     //
00235     // Check that the particle does not exit from the endcaps first.
00236     //
00237     zProp = Z() + ( phiProp - phi0 ) * pZ * radius / pT ;
00238     if ( fabs(zProp) > zCyl ) {
00239 
00240       zProp = zCyl * fabs(pZ)/pZ;
00241       phiProp = phi0 + (zProp - Z()) / radius * pT / pZ;
00242       success = 2;
00243 
00244     } else {
00245 
00246       //
00247       // If the particle does not cross the barrel, either process the
00248       // first-half loop only or propagate all the way down to the endcaps
00249       //
00250       if ( rProp < rCyl ) {
00251 
00252         if ( firstLoop || fabs(pZ)/pT < 1E-10 ) {
00253 
00254           success = 0;
00255 
00256         } else {
00257 
00258           zProp = zCyl * fabs(pZ)/pZ;
00259           phiProp = phi0 + (zProp - Z()) / radius * pT / pZ;
00260           success = 2;
00261 
00262         }
00263       //
00264       // The particle crossed the barrel 
00265       // 
00266       } else {
00267 
00268         success = 1;
00269 
00270       }
00271     }
00272     //
00273     // Compute the coordinates of the RawParticle after propagation
00274     //
00275     //
00276     // Check the decay time
00277     //
00278     double delTime = propDir * (phiProp-phi0)*radius*mass() / pT;
00279     double factor = 1.;
00280     properTime += delTime;
00281     if ( properTime > properDecayTime ) {
00282       factor = 1.-(properTime-properDecayTime)/delTime;
00283       properTime = properDecayTime;
00284       decayed = true;
00285     }
00286 
00287     zProp = Z() + (zProp-Z())*factor;
00288 
00289     phiProp = phi0 + (phiProp-phi0)*factor;
00290 
00291     double sProp = std::sin(phiProp);
00292     double cProp = std::cos(phiProp);
00293     double xProp = xC + radius * sProp;
00294     double yProp = yC - radius * cProp;
00295     double tProp = T() + (phiProp-phi0)*radius*E()/pT;
00296     double pxProp = pT * cProp;
00297     double pyProp = pT * sProp;
00298 
00299     //
00300     // Last check : Although propagated to the endcaps, R could still be
00301     // larger than rCyl because the cylinder was too short...
00302     //
00303     if ( success == 2 && xProp*xProp+yProp*yProp > rCyl2 ) { 
00304       success = -1;
00305       return true;
00306     }
00307     //
00308     // ... and update the particle vertex and momentum
00309     //
00310     setVertex( XYZTLorentzVector(xProp,yProp,zProp,tProp) );
00311     SetXYZT(pxProp,pyProp,pZ,E());
00312     //    SetPx(pxProp);
00313     //    SetPy(pyProp);
00314     return true;
00315  
00316   }
00317 }
00318 
00319 bool
00320 BaseParticlePropagator::backPropagate() {
00321 
00322   // Backpropagate
00323   SetXYZT(-Px(),-Py(),-Pz(),E());
00324   setCharge(-charge());
00325   propDir = -1;
00326   bool done = propagate();
00327   SetXYZT(-Px(),-Py(),-Pz(),E());
00328   setCharge(-charge());
00329   propDir = +1;
00330 
00331   return done;
00332 }
00333 
00334 BaseParticlePropagator
00335 BaseParticlePropagator::propagated() const {
00336   // Copy the input particle in a new instance
00337   BaseParticlePropagator myPropPart(*this);
00338   // Allow for many loops
00339   myPropPart.firstLoop=false;
00340   // Propagate the particle ! 
00341   myPropPart.propagate();
00342   // Back propagate if the forward propagation was not a success
00343   if (myPropPart.success < 0 )
00344     myPropPart.backPropagate();
00345   // Return the new instance
00346   return myPropPart;
00347 }
00348 
00349 void 
00350 BaseParticlePropagator::setPropagationConditions(double R, double Z, bool f) { 
00351   rCyl = R;
00352   rCyl2 = R*R;
00353   zCyl = Z;
00354   firstLoop = f;
00355 }
00356 
00357 bool
00358 BaseParticlePropagator::propagateToClosestApproach(double x0, double y0, bool first) {
00359 
00360   // Pre-computed quantities
00361   double pT = Pt();
00362   double radius = helixRadius(pT);
00363   double phi0 = helixStartPhi();
00364 
00365   // The Helix centre coordinates are computed wrt to the z axis
00366   // Actually the z axis might be centred on 0,0: it's the beam spot position (x0,y0)!
00367   double xC = helixCentreX(radius,phi0);
00368   double yC = helixCentreY(radius,phi0);
00369   double distz = helixCentreDistToAxis(xC-x0,yC-y0);
00370   double dist0 = helixCentreDistToAxis(xC,yC);
00371 
00372   //
00373   // Propagation to point of clostest approach to z axis
00374   //
00375   double rz,r0,z;
00376   if ( charge() != 0.0 && bField != 0.0 ) {
00377     rz = fabs ( fabs(radius) - distz ) + std::sqrt(x0*x0+y0*y0) + 0.0000001;
00378     r0 = fabs ( fabs(radius) - dist0 ) + 0.0000001;
00379   } else {
00380     rz = fabs( Px() * (Y()-y0) - Py() * (X()-x0) ) / Pt(); 
00381     r0 = fabs( Px() *  Y()     - Py() *  X() ) / Pt(); 
00382   }
00383 
00384   z = 999999.;
00385 
00386   // Propagate to the first interesection point 
00387   // with cylinder of radius sqrt(x0**2+y0**2)
00388   setPropagationConditions(rz , z, first);
00389   bool done = backPropagate();
00390 
00391   // Unsuccessful propagation - should not happen!
00392   if ( !done ) return done; 
00393 
00394   // The z axis is (0,0) - no need to go further
00395   if ( fabs(rz-r0) < 1E-10 ) return done;
00396   double dist1 = (X()-x0)*(X()-x0) + (Y()-y0)*(Y()-y0);
00397 
00398   // We are already at closest approach - no need to go further
00399   if ( dist1 < 1E-10 ) return done;
00400 
00401   // Keep for later if it happens to be the right solution
00402   XYZTLorentzVector vertex1 = vertex();
00403   XYZTLorentzVector momentum1 = momentum();
00404   
00405   // Propagate to the distance of closest approach to (0,0)
00406   setPropagationConditions(r0 , z, first);
00407   done = backPropagate();
00408   if ( !done ) return done; 
00409 
00410   // Propagate to the first interesection point 
00411   // with cylinder of radius sqrt(x0**2+y0**2)
00412   setPropagationConditions(rz , z, first);
00413   done = backPropagate();
00414   if ( !done ) return done; 
00415   double dist2 = (X()-x0)*(X()-x0) + (Y()-y0)*(Y()-y0);
00416 
00417   // Keep the good solution.
00418   if ( dist2 > dist1 ) { 
00419     setVertex(vertex1);
00420     SetXYZT(momentum1.X(),momentum1.Y(),momentum1.Z(),momentum1.E());
00421     dist2 = dist1;
00422   }
00423 
00424   // Done
00425   return done;
00426 
00427 }
00428 
00429 bool
00430 BaseParticlePropagator::propagateToEcal(bool first) {
00431   //
00432   // Propagation to Ecal (including preshower) after the
00433   // last Tracker layer
00434   // TODO: include proper geometry
00435   // Geometry taken from CMS ECAL TDR
00436   //
00437   // First propagate to global barrel / endcap cylinder 
00438   //  setPropagationConditions(1290. , 3045 , first);
00439   //  New position of the preshower as of CMSSW_1_7_X
00440   setPropagationConditions(129.0 , 303.353 , first);
00441   return propagate();
00442 
00443 }
00444 
00445 bool
00446 BaseParticlePropagator::propagateToPreshowerLayer1(bool first) {
00447   //
00448   // Propagation to Preshower Layer 1
00449   // TODO: include proper geometry
00450   // Geometry taken from CMS ECAL TDR
00451   //
00452   // First propagate to global barrel / endcap cylinder 
00453   //  setPropagationConditions(1290., 3045 , first);
00454   //  New position of the preshower as of CMSSW_1_7_X
00455   setPropagationConditions(129.0, 303.353 , first);
00456   bool done = propagate();
00457 
00458   // Check that were are on the Layer 1 
00459   if ( done && (R2() > 125.0*125.0 || R2() < 45.0*45.0) ) 
00460     success = 0;
00461   
00462   return done;
00463 }
00464 
00465 bool
00466 BaseParticlePropagator::propagateToPreshowerLayer2(bool first) {
00467   //
00468   // Propagation to Preshower Layer 2
00469   // TODO: include proper geometry
00470   // Geometry taken from CMS ECAL TDR
00471   //
00472   // First propagate to global barrel / endcap cylinder 
00473   //  setPropagationConditions(1290. , 3090 , first);
00474   //  New position of the preshower as of CMSSW_1_7_X
00475   setPropagationConditions(129.0 , 307.838 , first);
00476   bool done = propagate();
00477 
00478   // Check that we are on Layer 2 
00479   if ( done && (R2() > 125.0*125.0 || R2() < 45.0*45.0 ) )
00480     success = 0;
00481 
00482   return done;
00483 
00484 }
00485 
00486 bool
00487 BaseParticlePropagator::propagateToEcalEntrance(bool first) {
00488   //
00489   // Propagation to ECAL entrance
00490   // TODO: include proper geometry
00491   // Geometry taken from CMS ECAL TDR
00492   //
00493   // First propagate to global barrel / endcap cylinder 
00494   setPropagationConditions(129.0 , 320.9,first);
00495   bool done = propagate();
00496 
00497   // Go to endcap cylinder in the "barrel cut corner" 
00498   // eta = 1.479 -> cos^2(theta) = 0.81230
00499   //  if ( done && eta > 1.479 && success == 1 ) {
00500   if ( done && cos2ThetaV() > 0.81230 && success == 1 ) {
00501     setPropagationConditions(152.6 , 320.9, first);
00502     done = propagate();
00503   }
00504 
00505   // We are not in the ECAL acceptance
00506   // eta = 3.0 -> cos^2(theta) = 0.99013
00507   if ( cos2ThetaV() > 0.99014 ) success = 0;
00508 
00509   return done;
00510 }
00511 
00512 bool
00513 BaseParticlePropagator::propagateToHcalEntrance(bool first) {
00514   //
00515   // Propagation to HCAL entrance
00516   // TODO: include proper geometry
00517   // Geometry taken from CMSSW_3_1_X xml (Sunanda)
00518   //
00519 
00520   // First propagate to global barrel / endcap cylinder 
00521   setPropagationConditions(177.5 , 335.0, first);
00522   propDir = 0;
00523   bool done = propagate();
00524   propDir = 1;
00525 
00526   // If went through the bottom of HB cylinder -> re-propagate to HE surface
00527   if (done && success == 2) {
00528     setPropagationConditions(300.0, 400.458, first);
00529     propDir = 0;
00530     done = propagate();
00531     propDir = 1;
00532   }
00533 
00534 
00535   // out of the HB/HE acceptance
00536   // eta = 3.0 -> cos^2(theta) = 0.99014
00537   if ( done && cos2ThetaV() > 0.99014 ) success = 0;
00538 
00539   return done;
00540 }
00541 
00542 bool
00543 BaseParticlePropagator::propagateToVFcalEntrance(bool first) {
00544   //
00545   // Propagation to VFCAL entrance
00546   // TODO: include proper geometry
00547   // Geometry taken from DAQ TDR Chapter 13
00548 
00549   setPropagationConditions(400.0 , 1110.0, first);
00550   propDir = 0;
00551   bool done = propagate();
00552   propDir = 1;
00553 
00554   if (!done) success = 0;
00555 
00556   // We are not in the VFCAL acceptance
00557   // eta = 3.0  -> cos^2(theta) = 0.99014
00558   // eta = 5.2  -> cos^2(theta) = 0.9998755
00559   double c2teta = cos2ThetaV();
00560   if ( done && ( c2teta < 0.99014 || c2teta > 0.9998755 ) ) success = 0;
00561 
00562   return done;
00563 }
00564 
00565 bool
00566 BaseParticlePropagator::propagateToHcalExit(bool first) {
00567   //
00568   // Propagation to HCAL exit
00569   // TODO: include proper geometry
00570   // Geometry taken from CMSSW_3_1_X xml (Sunanda)
00571   //
00572 
00573   // Approximate it to a single cylinder as it is not that crucial.
00574   setPropagationConditions(263.9 , 554.1, first);
00575   //  this->rawPart().setCharge(0.0); ?? Shower Propagation ??
00576   propDir = 0;
00577   bool done = propagate();
00578   propDir = 1;
00579 
00580   return done;
00581 }
00582 
00583 bool
00584 BaseParticlePropagator::propagateToHOLayer(bool first) {
00585   //
00586   // Propagation to Layer0 of HO (Ring-0)
00587   // TODO: include proper geometry
00588   // Geometry taken from CMSSW_3_1_X xml (Sunanda)
00589   //
00590 
00591   // Approximate it to a single cylinder as it is not that crucial.
00592   setPropagationConditions(387.6, 1000.0, first); //Dia is the middle of HO \sqrt{384.8**2+46.2**2}
00593   //  setPropagationConditions(410.2, 1000.0, first); //sqrt{406.6**2+54.2**2} //for outer layer
00594   //  this->rawPart().setCharge(0.0); ?? Shower Propagation ??
00595   propDir = 0;
00596   bool done = propagate();
00597   propDir = 1;
00598   
00599   if ( done && abs(Z()) > 700.25) success = 0;
00600   
00601   return done;
00602 }
00603 
00604 
00605 bool
00606 BaseParticlePropagator::propagateToNominalVertex(const XYZTLorentzVector& v) 
00607 {
00608 
00609 
00610   // Not implemented for neutrals (used for electrons only)
00611   if ( charge() == 0. || bField == 0.) return false;
00612 
00613   // Define the proper pT direction to meet the point (vx,vy) in (x,y)
00614   double dx = X()-v.X();
00615   double dy = Y()-v.Y();
00616   double phi = std::atan2(dy,dx) + std::asin ( std::sqrt(dx*dx+dy*dy)/(2.*helixRadius()) );
00617 
00618   // The absolute pT (kept to the original value)
00619   double pT = pt();
00620 
00621   // Set the new pT
00622   SetPx(pT*std::cos(phi));
00623   SetPy(pT*std::sin(phi));
00624 
00625   return propagateToClosestApproach(v.X(),v.Y());
00626   
00627 }
00628 
00629 bool
00630 BaseParticlePropagator::propagateToBeamCylinder(const XYZTLorentzVector& v, double radius) 
00631 {
00632 
00633   // For neutral (or BField = 0, simply check that the track passes through the cylinder
00634   if ( charge() == 0. || bField == 0.) 
00635     return fabs( Px() * Y() - Py() * X() ) / Pt() < radius; 
00636 
00637   // Now go to the charged particles
00638 
00639   // A few needed intermediate variables (to make the code understandable)
00640   // The inner cylinder
00641   double r = radius;
00642   double r2 = r*r;
00643   double r4 = r2*r2;
00644   // The two hits
00645   double dx = X()-v.X();
00646   double dy = Y()-v.Y();
00647   double dz = Z()-v.Z();
00648   double Sxy = X()*v.X() + Y()*v.Y();
00649   double Dxy = Y()*v.X() - X()*v.Y();
00650   double Dxy2 = Dxy*Dxy;
00651   double Sxy2 = Sxy*Sxy;
00652   double SDxy = dx*dx + dy*dy;
00653   double SSDxy = std::sqrt(SDxy);
00654   double ATxy = std::atan2(dy,dx);
00655   // The coefficients of the second order equation for the trajectory radius
00656   double a = r2 - Dxy2/SDxy;
00657   double b = r * (r2 - Sxy);
00658   double c = r4 - 2.*Sxy*r2 + Sxy2 + Dxy2;
00659 
00660   // Here are the four possible solutions for
00661   // 1) The trajectory radius
00662   std::vector<double> helixRs(4,static_cast<double>(0.));
00663   helixRs[0] = (b - std::sqrt(b*b - a*c))/(2.*a); 
00664   helixRs[1] = (b + std::sqrt(b*b - a*c))/(2.*a); 
00665   helixRs[2] = -helixRs[0]; 
00666   helixRs[3] = -helixRs[1]; 
00667   // 2) The azimuthal direction at the second point
00668   std::vector<double> helixPhis(4,static_cast<double>(0.));
00669   helixPhis[0] = std::asin ( SSDxy/(2.*helixRs[0]) );
00670   helixPhis[1] = std::asin ( SSDxy/(2.*helixRs[1]) );
00671   helixPhis[2] = -helixPhis[0];
00672   helixPhis[3] = -helixPhis[1];
00673   // Only two solutions are valid though
00674   double solution1=0.;
00675   double solution2=0.;
00676   double phi1=0.;
00677   double phi2=0.;
00678   // Loop over the four possibilities
00679   for ( unsigned int i=0; i<4; ++i ) {
00680     helixPhis[i] = ATxy + helixPhis[i];
00681     // The centre of the helix
00682     double xC = helixCentreX(helixRs[i],helixPhis[i]);
00683     double yC = helixCentreY(helixRs[i],helixPhis[i]);
00684     // The distance of that centre to the origin
00685     double dist = helixCentreDistToAxis(xC, yC);
00686     /*
00687     std::cout << "Solution "<< i << " = " << helixRs[i] << " accepted ? " 
00688               << fabs(fabs(helixRs[i]) - dist ) << " - " << radius << " = " 
00689               << fabs(fabs(helixRs[i]) - dist ) - fabs(radius) << " " 
00690               << ( fabs( fabs(fabs(helixRs[i]) - dist) -fabs(radius) ) < 1e-6)  << std::endl;
00691     */
00692     // Check that the propagation will lead to a trajectroy tangent to the inner cylinder
00693     if ( fabs( fabs(fabs(helixRs[i]) - dist) -fabs(radius) ) < 1e-6 ) { 
00694       // Fill first solution
00695       if ( solution1 == 0. ) { 
00696         solution1 = helixRs[i];
00697         phi1 = helixPhis[i];
00698       // Fill second solution (ordered)
00699       } else if ( solution2 == 0. ) {
00700         if ( helixRs[i] < solution1 ) {
00701           solution2 = solution1;
00702           solution1 = helixRs[i];
00703           phi2 = phi1;
00704           phi1 = helixPhis[i];
00705         } else {
00706           solution2 = helixRs[i];
00707           phi2 = helixPhis[i];
00708         }
00709       // Must not happen! 
00710       } else { 
00711         std::cout << "warning !!! More than two solutions for this track !!! " << std::endl;
00712       }
00713     }
00714   }
00715 
00716   // Find the largest possible pT compatible with coming from within the inne cylinder
00717   double pT = 0.;
00718   double helixPhi = 0.;
00719   double helixR = 0.;
00720   // This should not happen
00721   if ( solution1*solution2 == 0. ) { 
00722     std::cout << "warning !!! Less than two solution for this track! " 
00723               << solution1 << " " << solution2 << std::endl;
00724     return false;
00725   // If the two solutions have opposite sign, it means that pT = infinity is a possibility.
00726   } else if ( solution1*solution2 < 0. ) { 
00727     pT = 1000.;
00728     double norm = pT/SSDxy;
00729     setCharge(+1.);
00730     SetXYZT(dx*norm,dy*norm,dz*norm,0.);
00731     SetE(std::sqrt(Vect().Mag2()));
00732   // Otherwise take the solution that gives the largest transverse momentum 
00733   } else { 
00734     if (solution1<0.) { 
00735       helixR   = solution1;
00736       helixPhi = phi1;
00737       setCharge(+1.);
00738     } else {
00739       helixR = solution2;
00740       helixPhi = phi2;
00741       setCharge(-1.);
00742     }
00743     pT = fabs(helixR) * 1e-5 * c_light() *bField;
00744     double norm = pT/SSDxy;
00745     SetXYZT(pT*std::cos(helixPhi),pT*std::sin(helixPhi),dz*norm,0.);
00746     SetE(std::sqrt(Vect().Mag2()));      
00747   }
00748     
00749   // Propagate to closest approach to get the Z value (a bit of an overkill)
00750   return propagateToClosestApproach();
00751   
00752 }
00753 
00754 double 
00755 BaseParticlePropagator::xyImpactParameter(double x0, double y0) const {
00756 
00757   double ip=0.;
00758   double pT = Pt();
00759 
00760   if ( charge() != 0.0 && bField != 0.0 ) {
00761     double radius = helixRadius(pT);
00762     double phi0 = helixStartPhi();
00763     
00764     // The Helix centre coordinates are computed wrt to the z axis
00765     double xC = helixCentreX(radius,phi0);
00766     double yC = helixCentreY(radius,phi0);
00767     double distz = helixCentreDistToAxis(xC-x0,yC-y0);
00768     ip = distz - fabs(radius);
00769   } else {
00770     ip = fabs( Px() * (Y()-y0) - Py() * (X()-x0) ) / pT; 
00771   }
00772 
00773   return ip;
00774 
00775 }