CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/TrackPropagation/SteppingHelixPropagator/src/SteppingHelixPropagator.cc

Go to the documentation of this file.
00001 
00013 //
00014 // Original Author:  Vyacheslav Krutelyov
00015 //         Created:  Fri Mar  3 16:01:24 CST 2006
00016 // $Id: SteppingHelixPropagator.cc,v 1.75 2011/02/12 17:55:25 slava77 Exp $
00017 //
00018 //
00019 
00020 
00021 #include "MagneticField/Engine/interface/MagneticField.h"
00022 #include "MagneticField/VolumeBasedEngine/interface/VolumeBasedMagneticField.h"
00023 #include "MagneticField/VolumeGeometry/interface/MagVolume.h"
00024 #include "MagneticField/Interpolation/interface/MFGrid.h"
00025 
00026 #include "DataFormats/GeometrySurface/interface/Cylinder.h"
00027 #include "DataFormats/GeometrySurface/interface/Plane.h"
00028 #include "DataFormats/GeometrySurface/interface/Cone.h"
00029 
00030 #include "TrackingTools/GeomPropagators/interface/PropagationExceptions.h"
00031 
00032 #include "TrackingTools/AnalyticalJacobians/interface/JacobianCartesianToCurvilinear.h"
00033 #include "TrackingTools/AnalyticalJacobians/interface/JacobianCurvilinearToCartesian.h"
00034 
00035 #include "TrackPropagation/SteppingHelixPropagator/interface/SteppingHelixPropagator.h"
00036 
00037 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00038 
00039 #include <sstream>
00040 #include <typeinfo>
00041 
00042 SteppingHelixPropagator::SteppingHelixPropagator() :
00043   Propagator(anyDirection)
00044 {
00045   field_ = 0;
00046 }
00047 
00048 SteppingHelixPropagator::SteppingHelixPropagator(const MagneticField* field, 
00049                                                  PropagationDirection dir):
00050   Propagator(dir),
00051   unit55_(AlgebraicMatrixID())
00052 {
00053   field_ = field;
00054   vbField_ = dynamic_cast<const VolumeBasedMagneticField*>(field_);
00055   covCurvRot_ = AlgebraicMatrix55();
00056   dCCurvTransform_ = unit55_;
00057   debug_ = false;
00058   noMaterialMode_ = false;
00059   noErrorPropagation_ = false;
00060   applyRadX0Correction_ = true;
00061   useMagVolumes_ = true;
00062   useIsYokeFlag_ = true;
00063   useMatVolumes_ = true;
00064   useInTeslaFromMagField_ = false; //overrides behavior only if true
00065   returnTangentPlane_ = true;
00066   sendLogWarning_ = false;
00067   useTuningForL2Speed_ = false;
00068   for (int i = 0; i <= MAX_POINTS; i++){
00069     svBuf_[i].p3 = Vector(0,0,0);
00070     svBuf_[i].r3 = Point(0,0,0);
00071     svBuf_[i].bf = Vector(0,0,0);
00072     svBuf_[i].bfGradLoc =  Vector(0,0,0);
00073     svBuf_[i].covCurv = AlgebraicSymMatrix55();
00074     svBuf_[i].matDCovCurv = AlgebraicSymMatrix55();
00075     svBuf_[i].isComplete = true;
00076     svBuf_[i].isValid_ = true;
00077     svBuf_[i].hasErrorPropagated_ = !noErrorPropagation_;
00078   }
00079   defaultStep_ = 5.;
00080 
00081   ecShiftPos_ = 0;
00082   ecShiftNeg_ = 0;
00083 
00084 }
00085 
00086 TrajectoryStateOnSurface 
00087 SteppingHelixPropagator::propagate(const FreeTrajectoryState& ftsStart, const Plane& pDest) const {
00088   return propagateWithPath(ftsStart, pDest).first;
00089 }
00090 
00091 TrajectoryStateOnSurface 
00092 SteppingHelixPropagator::propagate(const FreeTrajectoryState& ftsStart, const Cylinder& cDest) const
00093 {
00094   return propagateWithPath(ftsStart, cDest).first;
00095 }
00096 
00097 FreeTrajectoryState
00098 SteppingHelixPropagator::propagate(const FreeTrajectoryState& ftsStart, const GlobalPoint& pDest) const
00099 {
00100   return propagateWithPath(ftsStart, pDest).first;
00101 }
00102 
00103 FreeTrajectoryState
00104 SteppingHelixPropagator::propagate(const FreeTrajectoryState& ftsStart, 
00105                                    const GlobalPoint& pDest1, const GlobalPoint& pDest2) const
00106 {
00107   return propagateWithPath(ftsStart, pDest1, pDest2).first;
00108 }
00109 
00110 FreeTrajectoryState
00111 SteppingHelixPropagator::propagate(const FreeTrajectoryState& ftsStart, 
00112                                    const reco::BeamSpot& beamSpot) const
00113 {
00114   return propagateWithPath(ftsStart, beamSpot).first;
00115 }
00116 
00117 
00118 std::pair<TrajectoryStateOnSurface, double> 
00119 SteppingHelixPropagator::propagateWithPath(const FreeTrajectoryState& ftsStart, 
00120                                            const Plane& pDest) const {
00121 
00122   setIState(SteppingHelixStateInfo(ftsStart));
00123 
00124   const StateInfo& svCurrent = propagate(svBuf_[0], pDest);
00125 
00126   return TsosPP(svCurrent.getStateOnSurface(pDest), svCurrent.path());
00127 }
00128 
00129 std::pair<TrajectoryStateOnSurface, double> 
00130 SteppingHelixPropagator::propagateWithPath(const FreeTrajectoryState& ftsStart, 
00131                                            const Cylinder& cDest) const {
00132 
00133   setIState(SteppingHelixStateInfo(ftsStart));
00134 
00135   const StateInfo& svCurrent = propagate(svBuf_[0], cDest);
00136 
00137   return TsosPP(svCurrent.getStateOnSurface(cDest, returnTangentPlane_), svCurrent.path());
00138 }
00139 
00140 
00141 std::pair<FreeTrajectoryState, double> 
00142 SteppingHelixPropagator::propagateWithPath(const FreeTrajectoryState& ftsStart, 
00143                                            const GlobalPoint& pDest) const {
00144   setIState(SteppingHelixStateInfo(ftsStart));
00145 
00146   const StateInfo& svCurrent = propagate(svBuf_[0], pDest);
00147 
00148   FreeTrajectoryState ftsDest;
00149   svCurrent.getFreeState(ftsDest);
00150 
00151   return FtsPP(ftsDest, svCurrent.path());
00152 }
00153 
00154 std::pair<FreeTrajectoryState, double> 
00155 SteppingHelixPropagator::propagateWithPath(const FreeTrajectoryState& ftsStart, 
00156                                            const GlobalPoint& pDest1, const GlobalPoint& pDest2) const {
00157 
00158   if ((pDest1-pDest2).mag() < 1e-10){
00159     if (sendLogWarning_){
00160       edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: the points should be at a bigger distance"
00161                                                 <<std::endl;
00162     }
00163     return FtsPP();
00164   }
00165   setIState(SteppingHelixStateInfo(ftsStart));
00166   
00167   const StateInfo& svCurrent = propagate(svBuf_[0], pDest1, pDest2);
00168 
00169   FreeTrajectoryState ftsDest;
00170   svCurrent.getFreeState(ftsDest);
00171 
00172   return FtsPP(ftsDest, svCurrent.path());
00173 }
00174 
00175 
00176 std::pair<FreeTrajectoryState, double>  
00177 SteppingHelixPropagator::propagateWithPath(const FreeTrajectoryState& ftsStart,  
00178                                            const reco::BeamSpot& beamSpot) const {
00179   GlobalPoint pDest1(beamSpot.x0(), beamSpot.y0(), beamSpot.z0());
00180   GlobalPoint pDest2(pDest1.x() + beamSpot.dxdz()*1e3, 
00181                      pDest1.y() + beamSpot.dydz()*1e3,
00182                      pDest1.z() + 1e3);
00183   return propagateWithPath(ftsStart, pDest1, pDest2);
00184 }
00185 
00186 const SteppingHelixStateInfo&
00187 SteppingHelixPropagator::propagate(const SteppingHelixStateInfo& sStart, 
00188                                    const Surface& sDest) const {
00189   
00190   if (! sStart.isValid()){
00191     if (sendLogWarning_){
00192       edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: invalid input state"
00193                                                 <<std::endl;
00194     }
00195     return invalidState_;
00196   }
00197 
00198   const Plane* pDest = dynamic_cast<const Plane*>(&sDest);
00199   if (pDest != 0) return propagate(sStart, *pDest);
00200 
00201   const Cylinder* cDest = dynamic_cast<const Cylinder*>(&sDest);
00202   if (cDest != 0) return propagate(sStart, *cDest);
00203 
00204   throw PropagationException("The surface is neither Cylinder nor Plane");
00205 
00206 }
00207 
00208 const SteppingHelixStateInfo&
00209 SteppingHelixPropagator::propagate(const SteppingHelixStateInfo& sStart, 
00210                                    const Plane& pDest) const {
00211   
00212   if (! sStart.isValid()){
00213     if (sendLogWarning_){
00214       edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: invalid input state"
00215                                                 <<std::endl;
00216     }    
00217     return invalidState_;
00218   }
00219   setIState(sStart);
00220   
00221   Point rPlane(pDest.position().x(), pDest.position().y(), pDest.position().z());
00222   Vector nPlane(pDest.rotation().zx(), pDest.rotation().zy(), pDest.rotation().zz()); nPlane /= nPlane.mag();
00223 
00224   double pars[6] = { rPlane.x(), rPlane.y(), rPlane.z(),
00225                      nPlane.x(), nPlane.y(), nPlane.z() };
00226   
00227   propagate(PLANE_DT, pars);
00228   
00229   //(re)set it before leaving: dir =1 (-1) if path increased (decreased) and 0 if it didn't change
00230   //need to implement this somewhere else as a separate function
00231   double lDir = 0;
00232   if (sStart.path() < svBuf_[cIndex_(nPoints_-1)].path()) lDir = 1.;
00233   if (sStart.path() > svBuf_[cIndex_(nPoints_-1)].path()) lDir = -1.;
00234   svBuf_[cIndex_(nPoints_-1)].dir = lDir;
00235   return svBuf_[cIndex_(nPoints_-1)];
00236 }
00237 
00238 const SteppingHelixStateInfo&
00239 SteppingHelixPropagator::propagate(const SteppingHelixStateInfo& sStart, 
00240                                    const Cylinder& cDest) const {
00241   
00242   if (! sStart.isValid()){
00243     if (sendLogWarning_){
00244       edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: invalid input state"
00245                                                 <<std::endl;
00246     }    
00247     return invalidState_;
00248   }
00249   setIState(sStart);
00250   
00251   double pars[6] = {0,0,0,0,0,0};
00252   pars[RADIUS_P] = cDest.radius();
00253 
00254   
00255   propagate(RADIUS_DT, pars);
00256   
00257   //(re)set it before leaving: dir =1 (-1) if path increased (decreased) and 0 if it didn't change
00258   //need to implement this somewhere else as a separate function
00259   double lDir = 0;
00260   if (sStart.path() < svBuf_[cIndex_(nPoints_-1)].path()) lDir = 1.;
00261   if (sStart.path() > svBuf_[cIndex_(nPoints_-1)].path()) lDir = -1.;
00262   svBuf_[cIndex_(nPoints_-1)].dir = lDir;
00263   return svBuf_[cIndex_(nPoints_-1)];
00264 }
00265 
00266 const SteppingHelixStateInfo&
00267 SteppingHelixPropagator::propagate(const SteppingHelixStateInfo& sStart, 
00268                                    const GlobalPoint& pDest) const {
00269   
00270   if (! sStart.isValid()){
00271     if (sendLogWarning_){
00272       edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: invalid input state"
00273                                                 <<std::endl;
00274     }    
00275     return invalidState_;
00276   }
00277   setIState(sStart);
00278   
00279   double pars[6] = {pDest.x(), pDest.y(), pDest.z(), 0, 0, 0};
00280 
00281   
00282   propagate(POINT_PCA_DT, pars);
00283   
00284   return svBuf_[cIndex_(nPoints_-1)];
00285 }
00286 
00287 const SteppingHelixStateInfo&
00288 SteppingHelixPropagator::propagate(const SteppingHelixStateInfo& sStart, 
00289                                    const GlobalPoint& pDest1, const GlobalPoint& pDest2) const {
00290   
00291   if ((pDest1-pDest2).mag() < 1e-10 || !sStart.isValid()){
00292     if (sendLogWarning_){
00293       if ((pDest1-pDest2).mag() < 1e-10)
00294         edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: points are too close"
00295                                                   <<std::endl;
00296       if (!sStart.isValid())
00297         edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: invalid input state"
00298                                                   <<std::endl;
00299     }
00300     return invalidState_;
00301   }
00302   setIState(sStart);
00303   
00304   double pars[6] = {pDest1.x(), pDest1.y(), pDest1.z(),
00305                     pDest2.x(), pDest2.y(), pDest2.z()};
00306   
00307   propagate(LINE_PCA_DT, pars);
00308   
00309   return svBuf_[cIndex_(nPoints_-1)];
00310 }
00311 
00312 void SteppingHelixPropagator::setIState(const SteppingHelixStateInfo& sStart) const {
00313   nPoints_ = 0;
00314   svBuf_[cIndex_(nPoints_)] = sStart; //do this anyways to have a fresh start
00315   if (sStart.isComplete ) {
00316     svBuf_[cIndex_(nPoints_)] = sStart;
00317     nPoints_++;
00318   } else {
00319     loadState(svBuf_[cIndex_(nPoints_)], sStart.p3, sStart.r3, sStart.q,
00320               propagationDirection(), sStart.covCurv);
00321     nPoints_++;
00322   }
00323   svBuf_[cIndex_(0)].hasErrorPropagated_ = sStart.hasErrorPropagated_ & !noErrorPropagation_;
00324 }
00325 
00326 SteppingHelixPropagator::Result 
00327 SteppingHelixPropagator::propagate(SteppingHelixPropagator::DestType type, 
00328                                    const double pars[6], double epsilon)  const{
00329 
00330   static const std::string metname = "SteppingHelixPropagator";
00331   StateInfo* svCurrent = &svBuf_[cIndex_(nPoints_-1)];
00332 
00333   //check if it's going to work at all
00334   double tanDist = 0;
00335   double dist = 0;
00336   PropagationDirection refDirection = anyDirection;
00337   Result result = refToDest(type, (*svCurrent), pars, dist, tanDist, refDirection);
00338 
00339   if (result != SteppingHelixStateInfo::OK || fabs(dist) > 1e6){
00340     svCurrent->status_ = result;
00341     if (fabs(dist) > 1e6) svCurrent->status_ = SteppingHelixStateInfo::INACC;
00342     svCurrent->isValid_ = false;
00343     svCurrent->field = field_;
00344     if (sendLogWarning_){
00345       edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Failed after first refToDest check with status "
00346                               <<SteppingHelixStateInfo::ResultName[result]
00347                               <<std::endl;
00348     }
00349     return result;
00350   }
00351 
00352   result = SteppingHelixStateInfo::UNDEFINED;
00353   bool makeNextStep = true;
00354   double dStep = defaultStep_;
00355   PropagationDirection dir,oldDir;
00356   dir = propagationDirection(); 
00357   oldDir = dir;
00358   int nOsc = 0;
00359 
00360   double distMag = 1e12;
00361   double tanDistMag = 1e12;
00362 
00363   double distMat = 1e12;
00364   double tanDistMat = 1e12;
00365 
00366   double tanDistNextCheck = -0.1;//just need a negative start val
00367   double tanDistMagNextCheck = -0.1;
00368   double tanDistMatNextCheck = -0.1;
00369   double oldDStep = 0;
00370   PropagationDirection oldRefDirection = propagationDirection();
00371 
00372   Result resultToMat = SteppingHelixStateInfo::UNDEFINED;
00373   Result resultToMag = SteppingHelixStateInfo::UNDEFINED;
00374 
00375   bool isFirstStep = true;
00376   bool expectNewMagVolume = false;
00377 
00378   int loopCount = 0;
00379   while (makeNextStep){
00380     dStep = defaultStep_;
00381     svCurrent = &svBuf_[cIndex_(nPoints_-1)];
00382     double curZ = svCurrent->r3.z();
00383     double curR = svCurrent->r3.perp();
00384     if ( fabs(curZ) < 440 && curR < 260) dStep = defaultStep_*2;
00385 
00386     //more such ifs might be scattered around
00387     //even though dStep is large, it will still make stops at each volume boundary
00388     if (useTuningForL2Speed_){
00389       dStep = 100.;
00390       if (! useIsYokeFlag_ && fabs(curZ) < 667 && curR > 380 && curR < 850){
00391         dStep = 5*(1+0.2*svCurrent->p3.mag());
00392       }
00393     }
00394 
00395     //    refDirection = propagationDirection();
00396 
00397     tanDistNextCheck -= oldDStep;
00398     tanDistMagNextCheck -= oldDStep;
00399     tanDistMatNextCheck -= oldDStep;
00400     
00401     if (tanDistNextCheck < 0){
00402       //use pre-computed values if it's the first step
00403       if (! isFirstStep) refToDest(type, (*svCurrent), pars, dist, tanDist, refDirection);
00404       // constrain allowed path for a tangential approach
00405       if (fabs(tanDist) > 4.*(fabs(dist)) ) tanDist *= tanDist == 0 ? 0 :fabs(dist/tanDist*4.);
00406 
00407       tanDistNextCheck = fabs(tanDist)*0.5 - 0.5; //need a better guess (to-do)
00408       //reasonable limit
00409       if (tanDistNextCheck >  defaultStep_*20. ) tanDistNextCheck = defaultStep_*20.;
00410       oldRefDirection = refDirection;
00411     } else {
00412       tanDist  = tanDist > 0. ? tanDist - oldDStep : tanDist + oldDStep; 
00413       refDirection = oldRefDirection;
00414       if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Skipped refToDest: guess tanDist = "<<tanDist
00415                                    <<" next check at "<<tanDistNextCheck<<std::endl;
00416     }
00418     double fastSkipDist = fabs(dist) > fabs(tanDist) ? tanDist : dist;
00419 
00420     if (propagationDirection() == anyDirection){
00421       dir = refDirection;
00422     } else {
00423       dir = propagationDirection();
00424       if (fabs(tanDist)<0.1 && refDirection != dir ){
00425         //how did it get here?  nOsc++;
00426         dir = refDirection;
00427         if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"NOTE: overstepped last time: switch direction (can do it if within 1 mm)"<<std::endl;
00428       }
00429     }    
00430 
00431     if (useMagVolumes_ && ! (fabs(dist) < fabs(epsilon))){//need to know the general direction
00432       if (tanDistMagNextCheck < 0){
00433         resultToMag = refToMagVolume((*svCurrent), dir, distMag, tanDistMag, fabs(fastSkipDist), expectNewMagVolume, fabs(tanDist));
00434         // constrain allowed path for a tangential approach
00435         if (fabs(tanDistMag) > 4.*(fabs(distMag)) ) tanDistMag *= tanDistMag == 0 ? 0 : fabs(distMag/tanDistMag*4.);
00436 
00437         tanDistMagNextCheck = fabs(tanDistMag)*0.5-0.5; //need a better guess (to-do)
00438         //reasonable limit; "turn off" checking if bounds are further than the destination
00439         if (tanDistMagNextCheck >  defaultStep_*20. 
00440             || fabs(dist) < fabs(distMag)
00441             || resultToMag ==SteppingHelixStateInfo::INACC) 
00442           tanDistMagNextCheck  = defaultStep_*20 > fabs(fastSkipDist) ? fabs(fastSkipDist) : defaultStep_*20;   
00443         if (resultToMag != SteppingHelixStateInfo::INACC 
00444             && resultToMag != SteppingHelixStateInfo::OK) tanDistMagNextCheck = -1;
00445       } else {
00446         //      resultToMag = SteppingHelixStateInfo::OK;
00447         tanDistMag  = tanDistMag > 0. ? tanDistMag - oldDStep : tanDistMag + oldDStep; 
00448         if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Skipped refToMag: guess tanDistMag = "<<tanDistMag
00449                                      <<" next check at "<<tanDistMagNextCheck;
00450       }
00451     }
00452 
00453     if (useMatVolumes_ && ! (fabs(dist) < fabs(epsilon))){//need to know the general direction
00454       if (tanDistMatNextCheck < 0){
00455         resultToMat = refToMatVolume((*svCurrent), dir, distMat, tanDistMat, fabs(fastSkipDist));
00456         // constrain allowed path for a tangential approach
00457         if (fabs(tanDistMat) > 4.*(fabs(distMat)) ) tanDistMat *= tanDistMat == 0 ? 0 : fabs(distMat/tanDistMat*4.);
00458 
00459         tanDistMatNextCheck = fabs(tanDistMat)*0.5-0.5; //need a better guess (to-do)
00460         //reasonable limit; "turn off" checking if bounds are further than the destination
00461         if (tanDistMatNextCheck >  defaultStep_*20. 
00462             || fabs(dist) < fabs(distMat)
00463             || resultToMat ==SteppingHelixStateInfo::INACC ) 
00464           tanDistMatNextCheck = defaultStep_*20 > fabs(fastSkipDist) ? fabs(fastSkipDist) : defaultStep_*20;
00465         if (resultToMat != SteppingHelixStateInfo::INACC 
00466             && resultToMat != SteppingHelixStateInfo::OK) tanDistMatNextCheck = -1;
00467       } else {
00468         //      resultToMat = SteppingHelixStateInfo::OK;
00469         tanDistMat  = tanDistMat > 0. ? tanDistMat - oldDStep : tanDistMat + oldDStep; 
00470         if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Skipped refToMat: guess tanDistMat = "<<tanDistMat
00471                                      <<" next check at "<<tanDistMatNextCheck;
00472       }
00473     }
00474 
00475     double rDotP = svCurrent->r3.dot(svCurrent->p3);
00476     if ((fabs(curZ) > 1.5e3 || curR >800.) 
00477         && ((dir == alongMomentum && rDotP > 0) 
00478             || (dir == oppositeToMomentum && rDotP < 0) )
00479         ){
00480       dStep = fabs(tanDist) -1e-12;
00481     }
00482     double tanDistMin = fabs(tanDist);
00483     if (tanDistMin > fabs(tanDistMag)+0.05 && 
00484         (resultToMag == SteppingHelixStateInfo::OK || resultToMag == SteppingHelixStateInfo::WRONG_VOLUME)){
00485       tanDistMin = fabs(tanDistMag)+0.05;     //try to step into the next volume
00486       expectNewMagVolume = true;
00487     } else expectNewMagVolume = false;
00488 
00489     if (tanDistMin > fabs(tanDistMat)+0.05 && resultToMat == SteppingHelixStateInfo::OK){
00490       tanDistMin = fabs(tanDistMat)+0.05;     //try to step into the next volume
00491       if (expectNewMagVolume) expectNewMagVolume = false;
00492     }
00493 
00494     if (tanDistMin*fabs(svCurrent->dEdx) > 0.5*svCurrent->p3.mag()){
00495       tanDistMin = 0.5*svCurrent->p3.mag()/fabs(svCurrent->dEdx);
00496       if (expectNewMagVolume) expectNewMagVolume = false;
00497     }
00498 
00499 
00500 
00501     double tanDistMinLazy = fabs(tanDistMin);
00502     if ((type == POINT_PCA_DT || type == LINE_PCA_DT)
00503         && fabs(tanDist) < 2.*fabs(tanDistMin) ){
00504       //being lazy here; the best is to take into account the curvature
00505       tanDistMinLazy = fabs(tanDistMin)*0.5;
00506     }
00507  
00508     if (fabs(tanDistMinLazy) < dStep){
00509       dStep = fabs(tanDistMinLazy); 
00510     }
00511 
00512     //keep this path length for the next step
00513     oldDStep = dStep;
00514 
00515     if (dStep > 1e-10 && ! (fabs(dist) < fabs(epsilon))){
00516       StateInfo* svNext = &svBuf_[cIndex_(nPoints_)];
00517       makeAtomStep((*svCurrent), (*svNext), dStep, dir, HEL_AS_F);
00518 //       if (useMatVolumes_ && expectNewMagVolume 
00519 //        && svCurrent->magVol == svNext->magVol){
00520 //      double tmpDist=0;
00521 //      double tmpDistMag = 0;
00522 //      if (refToMagVolume((*svNext), dir, tmpDist, tmpDistMag, fabs(dist)) != SteppingHelixStateInfo::OK){
00523 //      //the point appears to be outside, but findVolume claims the opposite
00524 //        dStep += 0.05*fabs(tanDistMag/distMag); oldDStep = dStep; //do it again with a bigger step
00525 //        if (debug_) LogTrace(metname)
00526 //          <<"Failed to get into new mag volume: will try with new bigger step "<<dStep<<std::endl;
00527 //        makeAtomStep((*svCurrent), (*svNext), dStep, dir, HEL_AS_F);    
00528 //      }
00529 //       }
00530       nPoints_++;    svCurrent = &svBuf_[cIndex_(nPoints_-1)];
00531       if (oldDir != dir){
00532         nOsc++;
00533         tanDistNextCheck = -1;//check dist after osc
00534         tanDistMagNextCheck = -1;
00535         tanDistMatNextCheck = -1;
00536       }
00537       oldDir = dir;
00538     }
00539 
00540     if (nOsc>1 && fabs(dStep)>epsilon){
00541       if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Ooops"<<std::endl;
00542     }
00543 
00544     if (fabs(dist) < fabs(epsilon)  ) result = SteppingHelixStateInfo::OK;
00545 
00546     if ((type == POINT_PCA_DT || type == LINE_PCA_DT )
00547         && fabs(dStep) < fabs(epsilon)  ){
00548       //now check if it's not a branch point (peek ahead at 1 cm)
00549       double nextDist = 0;
00550       double nextTanDist = 0;
00551       PropagationDirection nextRefDirection = anyDirection;
00552       StateInfo* svNext = &svBuf_[cIndex_(nPoints_)];
00553       makeAtomStep((*svCurrent), (*svNext), 1., dir, HEL_AS_F);
00554       nPoints_++;     svCurrent = &svBuf_[cIndex_(nPoints_-1)];
00555       refToDest(type, (*svCurrent), pars, nextDist, nextTanDist, nextRefDirection);
00556       if ( fabs(nextDist) > fabs(dist)){
00557         nPoints_--;      svCurrent = &svBuf_[cIndex_(nPoints_-1)];
00558         result = SteppingHelixStateInfo::OK;
00559         if (debug_){
00560           LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Found real local minimum in PCA"<<std::endl;
00561         }
00562       } else {
00563         //keep this trial point and continue
00564         dStep = defaultStep_;
00565         if (debug_){
00566           LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Found branch point in PCA"<<std::endl;
00567         }
00568       }
00569     }
00570 
00571     if (nPoints_ > MAX_STEPS*1./defaultStep_  || loopCount > MAX_STEPS*100
00572         || nOsc > 6 ) result = SteppingHelixStateInfo::FAULT;
00573 
00574     if (svCurrent->p3.mag() < 0.1 ) result = SteppingHelixStateInfo::RANGEOUT;
00575 
00576     curZ = svCurrent->r3.z();
00577     curR = svCurrent->r3.perp();
00578     if ( curR > 20000 || fabs(curZ) > 20000 ) result = SteppingHelixStateInfo::INACC;
00579 
00580     makeNextStep = result == SteppingHelixStateInfo::UNDEFINED;
00581     svCurrent->status_ = result;
00582     svCurrent->isValid_ = (result == SteppingHelixStateInfo::OK || makeNextStep );
00583     svCurrent->field = field_;
00584 
00585     isFirstStep = false;
00586     loopCount++;
00587   }
00588 
00589   if (sendLogWarning_ && result != SteppingHelixStateInfo::OK){
00590     edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Propagation failed with status "
00591                             <<SteppingHelixStateInfo::ResultName[result]
00592                             <<std::endl;
00593     if (result == SteppingHelixStateInfo::RANGEOUT)
00594       edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Momentum at last point is too low (<0.1) p_last = "
00595                               <<svCurrent->p3.mag()
00596                               <<std::endl;
00597     if (result == SteppingHelixStateInfo::INACC)
00598       edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Went too far: the last point is at "<<svCurrent->r3
00599                               <<std::endl;
00600     if (result == SteppingHelixStateInfo::FAULT && nOsc > 6)
00601       edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Infinite loop condidtion detected: going in cycles. Break after 6 cycles"
00602                               <<std::endl;
00603     if (result == SteppingHelixStateInfo::FAULT && nPoints_ > MAX_STEPS*1./defaultStep_)
00604       edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Tired to go farther. Made too many steps: more than "
00605                               <<MAX_STEPS*1./defaultStep_
00606                               <<std::endl;
00607     
00608   }
00609 
00610   if (debug_){
00611     switch (type) {
00612     case RADIUS_DT:
00613       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to radius "<<pars[RADIUS_P]<<std::endl;
00614       break;
00615     case Z_DT:
00616       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to z "<<pars[Z_P]<<std::endl;
00617       break;
00618     case PATHL_DT:
00619       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to pathL "<<pars[PATHL_P]<<std::endl;
00620       break;
00621     case PLANE_DT:
00622       {
00623         Point rPlane(pars[0], pars[1], pars[2]);
00624         Vector nPlane(pars[3], pars[4], pars[5]);
00625         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to plane r0:"<<rPlane<<" n:"<<nPlane<<std::endl;
00626       }
00627       break;
00628     case POINT_PCA_DT:
00629       {
00630         Point rDest(pars[0], pars[1], pars[2]);
00631         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to PCA to point "<<rDest<<std::endl;
00632       }
00633       break;
00634     case LINE_PCA_DT:
00635       {
00636         Point rDest1(pars[0], pars[1], pars[2]);
00637         Point rDest2(pars[3], pars[4], pars[5]);
00638         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to PCA to line "<<rDest1<<" - "<<rDest2<<std::endl;
00639       }
00640       break;
00641     default:
00642       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to NOT IMPLEMENTED"<<std::endl;
00643       break;
00644     }
00645     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Made "<<nPoints_-1<<" steps and stopped at(cur step) "<<svCurrent->r3<<" nOsc "<<nOsc<<std::endl;
00646   }
00647   
00648   return result;
00649 }
00650   
00651 void SteppingHelixPropagator::loadState(SteppingHelixPropagator::StateInfo& svCurrent, 
00652                                         const SteppingHelixPropagator::Vector& p3, 
00653                                         const SteppingHelixPropagator::Point& r3, int charge,
00654                                         PropagationDirection dir,
00655                                         const AlgebraicSymMatrix55& covCurv) const{
00656   static const std::string metname = "SteppingHelixPropagator";
00657 
00658   svCurrent.q = charge;
00659   svCurrent.p3 = p3;
00660   svCurrent.r3 = r3;
00661   svCurrent.dir = dir == alongMomentum ? 1.: -1.;
00662 
00663   svCurrent.path_ = 0; // this could've held the initial path
00664   svCurrent.radPath_ = 0;
00665 
00666   GlobalPoint gPointNegZ(svCurrent.r3.x(), svCurrent.r3.y(), -fabs(svCurrent.r3.z()));
00667   GlobalPoint gPointNorZ(svCurrent.r3.x(), svCurrent.r3.y(), svCurrent.r3.z());
00668 
00669   float gpmag = gPointNegZ.mag2();
00670   if (gpmag > 1e20f ) {
00671     LogTrace(metname)<<"Initial point is too far";
00672     svCurrent.isValid_ = false;
00673     return;
00674   }
00675   if (! (gpmag == gpmag) ) {
00676     LogTrace(metname)<<"Initial point is a nan";
00677     svCurrent.isValid_ = false;
00678     return;
00679   }
00680 
00681   GlobalVector bf(0,0,0);
00682   // = field_->inTesla(gPoint);
00683   if (useMagVolumes_){
00684     if (vbField_ ){
00685       if (vbField_->isZSymmetric()){
00686         svCurrent.magVol = vbField_->findVolume(gPointNegZ);
00687       } else {
00688         svCurrent.magVol = vbField_->findVolume(gPointNorZ);
00689       }
00690       if (useIsYokeFlag_){
00691         double curRad = svCurrent.r3.perp();
00692         if (curRad > 380 && curRad < 850 && fabs(svCurrent.r3.z()) < 667){
00693           svCurrent.isYokeVol = isYokeVolume(svCurrent.magVol);
00694         } else {
00695           svCurrent.isYokeVol = false;
00696         }
00697       }
00698     } else {
00699       edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Failed to cast into VolumeBasedMagneticField: fall back to the default behavior"<<std::endl;
00700       svCurrent.magVol = 0;
00701     }
00702     if (debug_){
00703       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got volume at "<<svCurrent.magVol<<std::endl;
00704     }
00705   }
00706   
00707   if (useMagVolumes_ && svCurrent.magVol != 0 && ! useInTeslaFromMagField_){
00708     if (vbField_->isZSymmetric()){
00709       bf = svCurrent.magVol->inTesla(gPointNegZ);
00710       if (debug_){
00711         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Loaded bfield  float: "<<bf
00712                          <<" at global float "<< gPointNegZ<<" double "<< svCurrent.r3<<std::endl;
00713         LocalPoint lPoint(svCurrent.magVol->toLocal(gPointNegZ));
00714         LocalVector lbf = svCurrent.magVol->fieldInTesla(lPoint);
00715         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"\t cf in local locF: "<<lbf<<" at "<<lPoint<<std::endl;
00716       }
00717     } else {
00718       bf = svCurrent.magVol->inTesla(gPointNorZ);
00719       if (debug_){
00720         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Loaded bfield float: "<<bf
00721                          <<" at global float "<< gPointNorZ<<" double "<< svCurrent.r3<<std::endl;
00722         LocalPoint lPoint(svCurrent.magVol->toLocal(gPointNorZ));
00723         LocalVector lbf = svCurrent.magVol->fieldInTesla(lPoint);
00724         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"\t cf in local locF: "<<lbf<<" at "<<lPoint<<std::endl;
00725       }
00726     }
00727     if (r3.z() > 0 && vbField_->isZSymmetric() ){
00728       svCurrent.bf.set(-bf.x(), -bf.y(), bf.z());
00729     } else {
00730       svCurrent.bf.set(bf.x(), bf.y(), bf.z());
00731     }
00732   } else {
00733     GlobalPoint gPoint(r3.x(), r3.y(), r3.z());
00734     bf = field_->inTesla(gPoint);
00735     svCurrent.bf.set(bf.x(), bf.y(), bf.z());
00736   }
00737   if (svCurrent.bf.mag2() < 1e-32) svCurrent.bf.set(0., 0., 1e-16);
00738   if (debug_){
00739     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific
00740                      <<"Loaded bfield double: "<<svCurrent.bf<<"  from float: "<<bf
00741                      <<" at float "<< gPointNorZ<<" "<<gPointNegZ<<" double "<< svCurrent.r3<<std::endl;
00742   }
00743 
00744 
00745 
00746   double dEdXPrime = 0;
00747   double dEdx = 0;
00748   double radX0 = 0;
00749   MatBounds rzLims;
00750   dEdx = getDeDx(svCurrent, dEdXPrime, radX0, rzLims);
00751   svCurrent.dEdx = dEdx;    svCurrent.dEdXPrime = dEdXPrime;
00752   svCurrent.radX0 = radX0;
00753   svCurrent.rzLims = rzLims;
00754 
00755   svCurrent.covCurv =covCurv;
00756 
00757   svCurrent.isComplete = true;
00758   svCurrent.isValid_ = true;
00759 
00760   if (debug_){
00761     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Loaded at  path: "<<svCurrent.path()<<" radPath: "<<svCurrent.radPath()
00762                      <<" p3 "<<" pt: "<<svCurrent.p3.perp()<<" phi: "<<svCurrent.p3.phi()
00763                      <<" eta: "<<svCurrent.p3.eta()
00764                      <<" "<<svCurrent.p3
00765                      <<" r3: "<<svCurrent.r3
00766                      <<" bField: "<<svCurrent.bf.mag()
00767                      <<std::endl;
00768     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Input Covariance in Curvilinear RF "<<covCurv<<std::endl;
00769   }
00770 }
00771 
00772 void SteppingHelixPropagator::getNextState(const SteppingHelixPropagator::StateInfo& svPrevious, 
00773                                            SteppingHelixPropagator::StateInfo& svNext,
00774                                            double dP, const SteppingHelixPropagator::Vector& tau,
00775                                            const SteppingHelixPropagator::Vector& drVec, double dS, double dX0,
00776                                            const AlgebraicMatrix55& dCovCurvTransform) const{
00777   static const std::string metname = "SteppingHelixPropagator";
00778   svNext.q = svPrevious.q;
00779   svNext.dir = dS > 0.0 ? 1.: -1.; 
00780   svNext.p3 = tau;  svNext.p3*=(svPrevious.p3.mag() - svNext.dir*fabs(dP));
00781 
00782   svNext.r3 = svPrevious.r3; svNext.r3 += drVec;
00783 
00784   svNext.path_ = svPrevious.path() + dS;
00785   svNext.radPath_ = svPrevious.radPath() + dX0;
00786 
00787   GlobalPoint gPointNegZ(svNext.r3.x(), svNext.r3.y(), -fabs(svNext.r3.z()));
00788   GlobalPoint gPointNorZ(svNext.r3.x(), svNext.r3.y(), svNext.r3.z());
00789 
00790   GlobalVector bf(0,0,0); 
00791 
00792   if (useMagVolumes_){
00793     if (vbField_ != 0){
00794        if (vbField_->isZSymmetric()){
00795          svNext.magVol = vbField_->findVolume(gPointNegZ);
00796        } else {
00797          svNext.magVol = vbField_->findVolume(gPointNorZ);
00798        }
00799       if (useIsYokeFlag_){
00800         double curRad = svNext.r3.perp();
00801         if (curRad > 380 && curRad < 850 && fabs(svNext.r3.z()) < 667){
00802           svNext.isYokeVol = isYokeVolume(svNext.magVol);
00803         } else {
00804           svNext.isYokeVol = false;
00805         }
00806       }
00807     } else {
00808       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Failed to cast into VolumeBasedMagneticField"<<std::endl;
00809       svNext.magVol = 0;
00810     }
00811     if (debug_){
00812       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got volume at "<<svNext.magVol<<std::endl;
00813     }
00814   }
00815 
00816   if (useMagVolumes_ && svNext.magVol != 0 && ! useInTeslaFromMagField_){
00817     if (vbField_->isZSymmetric()){
00818       bf = svNext.magVol->inTesla(gPointNegZ);
00819     } else {
00820       bf = svNext.magVol->inTesla(gPointNorZ);
00821     }
00822     if (svNext.r3.z() > 0  && vbField_->isZSymmetric() ){
00823       svNext.bf.set(-bf.x(), -bf.y(), bf.z());
00824     } else {
00825       svNext.bf.set(bf.x(), bf.y(), bf.z());
00826     }
00827   } else {
00828     GlobalPoint gPoint(svNext.r3.x(), svNext.r3.y(), svNext.r3.z());
00829     bf = field_->inTesla(gPoint);
00830     svNext.bf.set(bf.x(), bf.y(), bf.z());
00831   }
00832   if (svNext.bf.mag2() < 1e-32) svNext.bf.set(0., 0., 1e-16);
00833   
00834   
00835   double dEdXPrime = 0;
00836   double dEdx = 0;
00837   double radX0 = 0;
00838   MatBounds rzLims;
00839   dEdx = getDeDx(svNext, dEdXPrime, radX0, rzLims);
00840   svNext.dEdx = dEdx;    svNext.dEdXPrime = dEdXPrime;
00841   svNext.radX0 = radX0;
00842   svNext.rzLims = rzLims;
00843 
00844   //update Emat only if it's valid
00845   svNext.hasErrorPropagated_ = svPrevious.hasErrorPropagated_;
00846   if (svPrevious.hasErrorPropagated_){
00847     {
00848       AlgebraicMatrix55 tmp = dCovCurvTransform*svPrevious.covCurv;
00849       ROOT::Math::AssignSym::Evaluate(svNext.covCurv, tmp*ROOT::Math::Transpose(dCovCurvTransform));
00850       
00851       svNext.covCurv += svPrevious.matDCovCurv;
00852     }
00853   } else {
00854     //could skip dragging along the unprop. cov later .. now
00855     // svNext.cov = svPrevious.cov;
00856   }
00857 
00858   if (debug_){
00859     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Now at  path: "<<svNext.path()<<" radPath: "<<svNext.radPath()
00860                      <<" p3 "<<" pt: "<<svNext.p3.perp()<<" phi: "<<svNext.p3.phi()
00861                      <<" eta: "<<svNext.p3.eta()
00862                      <<" "<<svNext.p3
00863                      <<" r3: "<<svNext.r3
00864                      <<" dPhi: "<<acos(svNext.p3.unit().dot(svPrevious.p3.unit()))
00865                      <<" bField: "<<svNext.bf.mag()
00866                      <<" dedx: "<<svNext.dEdx
00867                      <<std::endl;
00868     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"New Covariance "<<svNext.covCurv<<std::endl;
00869     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Transf by dCovTransform "<<dCovCurvTransform<<std::endl;
00870   }
00871 }
00872 
00873 void SteppingHelixPropagator::setRep(SteppingHelixPropagator::Basis& rep, 
00874                                      const SteppingHelixPropagator::Vector& tau) const{
00875   Vector zRep(0., 0., 1.);
00876   rep.lX = tau;
00877   rep.lY = zRep.cross(tau); rep.lY *= 1./tau.perp();
00878   rep.lZ = rep.lX.cross(rep.lY);
00879 }
00880 
00881 bool SteppingHelixPropagator::makeAtomStep(SteppingHelixPropagator::StateInfo& svCurrent,
00882                                            SteppingHelixPropagator::StateInfo& svNext,
00883                                            double dS, 
00884                                            PropagationDirection dir, 
00885                                            SteppingHelixPropagator::Fancy fancy) const{
00886   static const std::string metname = "SteppingHelixPropagator";
00887   if (debug_){
00888     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Make atom step "<<svCurrent.path()<<" with step "<<dS<<" in direction "<<dir<<std::endl;
00889   }
00890 
00891   double dP = 0;
00892   double curP = svCurrent.p3.mag();
00893   Vector tau = svCurrent.p3; tau *= 1./curP;
00894   Vector tauNext(tau);
00895   Vector drVec(0,0,0);
00896 
00897   dS = dir == alongMomentum ? fabs(dS) : -fabs(dS);
00898 
00899 
00900   double radX0 = 1e24;
00901 
00902   switch (fancy){
00903   case HEL_AS_F:
00904   case HEL_ALL_F:{
00905     double p0 = curP;// see above = svCurrent.p3.mag();
00906     double p0Inv = 1./p0;
00907     double b0 = svCurrent.bf.mag();
00908 
00909     //get to the mid-point first
00910     double phi = (2.99792458e-3*svCurrent.q*b0*p0Inv)*dS/2.;
00911     bool phiSmall = fabs(phi) < 1e-4;
00912 
00913     double cosPhi = 0;
00914     double sinPhi = 0;
00915 
00916     double oneLessCosPhi=0;
00917     double oneLessCosPhiOPhi=0;
00918     double sinPhiOPhi=0;
00919     double phiLessSinPhiOPhi=0;
00920 
00921     if (phiSmall){
00922       double phi2 = phi*phi;
00923       double phi3 = phi2*phi;
00924       double phi4 = phi3*phi;
00925       sinPhi = phi - phi3/6. + phi4*phi/120.;
00926       cosPhi = 1. -phi2/2. + phi4/24.;
00927       oneLessCosPhi = phi2/2. - phi4/24. + phi2*phi4/720.; // 0.5*phi*phi;//*(1.- phi*phi/12.);
00928       oneLessCosPhiOPhi = 0.5*phi - phi3/24. + phi2*phi3/720.;//*(1.- phi*phi/12.);
00929       sinPhiOPhi = 1. - phi*phi/6. + phi4/120.;
00930       phiLessSinPhiOPhi = phi*phi/6. - phi4/120. + phi4*phi2/5040.;//*(1. - phi*phi/20.);
00931     } else {
00932       cosPhi = cos(phi);
00933       sinPhi = sin(phi);
00934       oneLessCosPhi = 1.-cosPhi;
00935       oneLessCosPhiOPhi = oneLessCosPhi/phi;
00936       sinPhiOPhi = sinPhi/phi;
00937       phiLessSinPhiOPhi = 1 - sinPhiOPhi;
00938     }
00939 
00940     Vector bHat = svCurrent.bf; bHat *= 1./b0; //bHat.mag();
00941     //    bool bAlongZ = fabs(bHat.z()) > 0.9999;
00942 
00943     Vector btVec(bHat.cross(tau)); // for balong z btVec.z()==0
00944     double tauB =  tau.dot(bHat);
00945     Vector bbtVec(bHat*tauB - tau); // (-tau.x(), -tau.y(), 0)
00946 
00947     //don't need it here    tauNext = tau + bbtVec*oneLessCosPhi - btVec*sinPhi;
00948     drVec = bbtVec*phiLessSinPhiOPhi; drVec -= btVec*oneLessCosPhiOPhi; drVec += tau; 
00949     drVec *= dS/2.;
00950 
00951     double dEdx = svCurrent.dEdx;
00952     double dEdXPrime = svCurrent.dEdXPrime;
00953     radX0 = svCurrent.radX0;
00954     dP = dEdx*dS;
00955 
00956     //improve with above values:
00957     drVec += svCurrent.r3;
00958     GlobalVector bfGV(0,0,0);
00959     Vector bf(0,0,0); //(bfGV.x(), bfGV.y(), bfGV.z());
00960     // = svCurrent.magVol->inTesla(GlobalPoint(drVec.x(), drVec.y(), -fabs(drVec.z())));
00961     if (useMagVolumes_ && svCurrent.magVol != 0 && ! useInTeslaFromMagField_){
00962       // this negative-z business will break at some point
00963       if (vbField_->isZSymmetric()){
00964         bfGV = svCurrent.magVol->inTesla(GlobalPoint(drVec.x(), drVec.y(), -fabs(drVec.z())));
00965       } else {
00966         bfGV = svCurrent.magVol->inTesla(GlobalPoint(drVec.x(), drVec.y(), drVec.z()));
00967       }
00968       if (drVec.z() > 0 && vbField_->isZSymmetric()){
00969         bf.set(-bfGV.x(), -bfGV.y(), bfGV.z());
00970       } else {
00971         bf.set(bfGV.x(), bfGV.y(), bfGV.z());
00972       }
00973     } else {
00974       bfGV = field_->inTesla(GlobalPoint(drVec.x(), drVec.y(), drVec.z()));
00975       bf.set(bfGV.x(), bfGV.y(), bfGV.z());
00976     }
00977     double b0Init = b0;
00978     b0 = bf.mag();
00979     if (b0 < 1e-16) {
00980       b0 = 1e-16;
00981       bf.set(0., 0., 1e-16);
00982     }
00983     if (debug_){
00984       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Improved b "<<b0
00985                        <<" at r3 "<<drVec<<std::endl;
00986     }
00987 
00988     if (fabs((b0-b0Init)*dS) > 1){
00989       //missed the mag volume boundary?
00990       if (debug_){
00991         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Large bf*dS change "<<fabs((b0-svCurrent.bf.mag())*dS)
00992                          <<" --> recalc dedx"<<std::endl;
00993       }
00994       svNext.r3 = drVec;
00995       svNext.bf = bf;
00996       svNext.p3 = svCurrent.p3;
00997       svNext.isYokeVol = svCurrent.isYokeVol;
00998       MatBounds rzTmp;
00999       dEdx = getDeDx(svNext, dEdXPrime, radX0, rzTmp);
01000       dP = dEdx*dS;      
01001       if (debug_){
01002         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"New dEdX= "<<dEdx
01003                          <<" dP= "<<dP
01004                          <<" for p0 "<<p0<<std::endl;
01005       }
01006     }
01007     //p0 is mid-way and b0 from mid-point
01008     p0 += dP/2.; if (p0 < 1e-2) p0 = 1e-2;
01009     p0Inv = 1./p0;
01010 
01011     phi = (2.99792458e-3*svCurrent.q*b0*p0Inv)*dS;
01012     phiSmall = fabs(phi) < 1e-4;
01013 
01014     if (phiSmall){
01015       double phi2 = phi*phi;
01016       double phi3 = phi2*phi;
01017       double phi4 = phi3*phi;
01018       sinPhi = phi - phi3/6. + phi4*phi/120.;
01019       cosPhi = 1. -phi2/2. + phi4/24.;
01020       oneLessCosPhi = phi2/2. - phi4/24. + phi2*phi4/720.; // 0.5*phi*phi;//*(1.- phi*phi/12.);
01021       oneLessCosPhiOPhi = 0.5*phi - phi3/24. + phi2*phi3/720.;//*(1.- phi*phi/12.);
01022       sinPhiOPhi = 1. - phi*phi/6. + phi4/120.;
01023       phiLessSinPhiOPhi = phi*phi/6. - phi4/120. + phi4*phi2/5040.;//*(1. - phi*phi/20.);
01024     }else {
01025       cosPhi = cos(phi); 
01026       sinPhi = sin(phi);
01027       oneLessCosPhi = 1.-cosPhi;
01028       oneLessCosPhiOPhi = oneLessCosPhi/phi;
01029       sinPhiOPhi = sinPhi/phi;
01030       phiLessSinPhiOPhi = 1. - sinPhiOPhi;
01031     }
01032 
01033     bHat = bf; bHat *= 1./b0;// as above =1./bHat.mag();
01034     //    bAlongZ = fabs(bHat.z()) > 0.9999;
01035     btVec = bHat.cross(tau); // for b||z (-tau.y(), tau.x() ,0)
01036     tauB = tau.dot(bHat);
01037     bbtVec = bHat*tauB - tau; //bHat.cross(btVec); for b||z: (-tau.x(), -tau.y(), 0)
01038 
01039     tauNext = bbtVec*oneLessCosPhi; tauNext -= btVec*sinPhi;     tauNext += tau; //for b||z tauNext.z() == tau.z()
01040     double tauNextPerpInv = 1./tauNext.perp();
01041     drVec   = bbtVec*phiLessSinPhiOPhi; drVec -= btVec*oneLessCosPhiOPhi;  drVec += tau;
01042     drVec *= dS;
01043     
01044     
01045     if (svCurrent.hasErrorPropagated_){
01046       double theta02 = 0;
01047       double dX0 = fabs(dS)/radX0;
01048 
01049       if (applyRadX0Correction_){
01050         // this provides the integrand for theta^2
01051         // if summed up along the path, should result in 
01052         // theta_total^2 = Int_0^x0{ f(x)dX} = (13.6/p0)^2*x0*(1+0.036*ln(x0+1))
01053         // x0+1 above is to make the result infrared safe.
01054         double x0 = fabs(svCurrent.radPath());
01055         double alphaX0 = 13.6e-3*p0Inv; alphaX0 *= alphaX0;
01056         double betaX0 = 0.038;
01057         double logx0p1 = log(x0+1);
01058         theta02 = dX0*alphaX0*(1+betaX0*logx0p1)*(1 + betaX0*logx0p1 + 2.*betaX0*x0/(x0+1) );
01059       } else {
01060         theta02 = 196e-6* p0Inv * p0Inv * dX0; //14.e-3/p0*sqrt(fabs(dS)/radX0); // .. drop log term (this is non-additive)
01061       }
01062       
01063       double epsilonP0 = 1.+ dP/(p0-0.5*dP);
01064       //      double omegaP0 = -dP/(p0-0.5*dP) + dS*dEdXPrime;      
01065       //      double dsp = dS/(p0-0.5*dP); //use the initial p0 (not the mid-point) to keep the transport properly additive
01066       
01067       Vector tbtVec(bHat - tauB*tau); // for b||z tau.z()*(-tau.x(), -tau.y(), 1.-tau.z())
01068       
01069       dCCurvTransform_ = unit55_;
01070       {
01071         //Slightly modified copy of the curvilinear jacobian (don't use the original just because it's in float precision
01072         // and seems to have some assumptions about the field values
01073         // notation changes: p1--> tau, p2-->tauNext
01074         // theta --> phi
01075         //      Vector p1 = tau;
01076         //      Vector p2 = tauNext;
01077         Point xStart = svCurrent.r3;
01078         Vector dx = drVec;
01079         //GlobalVector h  = MagneticField::inInverseGeV(xStart);
01080         // Martijn: field is now given as parameter.. GlobalVector h  = globalParameters.magneticFieldInInverseGeV(xStart);
01081 
01082         //double qbp = fts.signedInverseMomentum();
01083         double qbp = svCurrent.q*p0Inv;
01084         //      double absS = dS;
01085   
01086         // calculate transport matrix
01087         // Origin: TRPRFN
01088         double t11 = tau.x(); double t12 = tau.y(); double t13 = tau.z();
01089         double t21 = tauNext.x(); double t22 = tauNext.y(); double t23 = tauNext.z();
01090         double cosl0 = tau.perp(); 
01091         //      double cosl1 = 1./tauNext.perp(); //not quite a cos .. it's a cosec--> change to cosecl1 below
01092         double cosecl1 = tauNextPerpInv;
01093         //AlgebraicMatrix a(5,5,1);
01094         // define average magnetic field and gradient 
01095         // at initial point - inlike TRPRFN
01096         Vector hn = bHat;
01097         //      double qp = -2.99792458e-3*b0;
01098         //   double q = -h.mag()*qbp;
01099 
01100         double q = -phi/dS; //qp*qbp; // -phi/dS
01101         //      double theta = -phi; 
01102         double sint = -sinPhi; double cost = cosPhi;
01103         double hn1 = hn.x(); double hn2 = hn.y(); double hn3 = hn.z();
01104         double dx1 = dx.x(); double dx2 = dx.y(); double dx3 = dx.z();
01105         //      double hnDt1 = hn1*t11 + hn2*t12 + hn3*t13;
01106 
01107         double gamma =  hn1*t21 + hn2*t22 + hn3*t23;
01108         double an1 =  hn2*t23 - hn3*t22;
01109         double an2 =  hn3*t21 - hn1*t23;
01110         double an3 =  hn1*t22 - hn2*t21;
01111         //        double auInv = sqrt(1.- t13*t13); double au = auInv>0 ? 1./auInv : 1e24;
01112         double auInv = cosl0; double au = auInv>0 ? 1./auInv : 1e24;
01113         //        double auInv = sqrt(t11*t11 + t12*t12); double au = auInv>0 ? 1./auInv : 1e24;
01114         double u11 = -au*t12; double u12 = au*t11;
01115         double v11 = -t13*u12; double v12 = t13*u11; double v13 = auInv;//t11*u12 - t12*u11;
01116         auInv = sqrt(1. - t23*t23); au = auInv>0 ? 1./auInv : 1e24;
01117         //        auInv = sqrt(t21*t21 + t22*t22); au = auInv>0 ? 1./auInv : 1e24;
01118         double u21 = -au*t22; double u22 = au*t21;
01119         double v21 = -t23*u22; double v22 = t23*u21; double v23 = auInv;//t21*u22 - t22*u21;
01120         // now prepare the transport matrix
01121         // pp only needed in high-p case (WA)
01122         //   double pp = 1./qbp;
01124         // moved up (where -h.mag() is needed()
01125         //   double qp = q*pp;
01126         double anv =  -(hn1*u21 + hn2*u22          );
01127         double anu =   (hn1*v21 + hn2*v22 + hn3*v23); 
01128         double omcost = oneLessCosPhi; double tmsint = -phi*phiLessSinPhiOPhi;
01129         
01130         double hu1 =         - hn3*u12;
01131         double hu2 =  hn3*u11;
01132         double hu3 =  hn1*u12 - hn2*u11;
01133         
01134         double hv1 =  hn2*v13 - hn3*v12;
01135         double hv2 =  hn3*v11 - hn1*v13;
01136         double hv3 =  hn1*v12 - hn2*v11;
01137         
01138         //   1/p - doesn't change since |tau| = |tauNext| ... not. It changes now
01139         dCCurvTransform_(0,0) = 1./(epsilonP0*epsilonP0)*(1. + dS*dEdXPrime);
01140         
01141         //   lambda
01142         
01143         dCCurvTransform_(1,0) = phi*p0/svCurrent.q*cosecl1*
01144           (sinPhi*bbtVec.z() - cosPhi*btVec.z());
01145         //was dCCurvTransform_(1,0) = -qp*anv*(t21*dx1 + t22*dx2 + t23*dx3); //NOTE (SK) this was found to have an opposite sign
01146         //from independent re-calculation ... in fact the tauNext.dot.dR piece isnt reproduced 
01147         
01148         dCCurvTransform_(1,1) = cost*(v11*v21 + v12*v22 + v13*v23) +
01149           sint*(hv1*v21 + hv2*v22 + hv3*v23) +
01150           omcost*(hn1*v11 + hn2*v12 + hn3*v13) * (hn1*v21 + hn2*v22 + hn3*v23) +
01151           anv*(-sint*(v11*t21 + v12*t22 + v13*t23) +
01152                omcost*(v11*an1 + v12*an2 + v13*an3) -
01153                tmsint*gamma*(hn1*v11 + hn2*v12 + hn3*v13) );
01154         
01155         dCCurvTransform_(1,2) = cost*(u11*v21 + u12*v22          ) +
01156           sint*(hu1*v21 + hu2*v22 + hu3*v23) +
01157           omcost*(hn1*u11 + hn2*u12          ) * (hn1*v21 + hn2*v22 + hn3*v23) +
01158           anv*(-sint*(u11*t21 + u12*t22          ) +
01159                omcost*(u11*an1 + u12*an2          ) -
01160                tmsint*gamma*(hn1*u11 + hn2*u12          ) );
01161         dCCurvTransform_(1,2) *= cosl0;
01162         
01163         // Commented out in part for reproducibility purposes: these terms are zero in cart->curv 
01164         //      dCCurvTransform_(1,3) = -q*anv*(u11*t21 + u12*t22          ); //don't show up in cartesian setup-->curv
01165         //why would lambdaNext depend explicitely on initial position ? any arbitrary init point can be chosen not 
01166         // affecting the final state's momentum direction ... is this the field gradient in curvilinear coord?
01167         //      dCCurvTransform_(1,4) = -q*anv*(v11*t21 + v12*t22 + v13*t23); //don't show up in cartesian setup-->curv
01168         
01169         //   phi
01170         
01171         dCCurvTransform_(2,0) = - phi*p0/svCurrent.q*cosecl1*cosecl1*
01172           (oneLessCosPhi*bHat.z()*btVec.mag2() + sinPhi*btVec.z() + cosPhi*tbtVec.z()) ;
01173         //was   dCCurvTransform_(2,0) = -qp*anu*(t21*dx1 + t22*dx2 + t23*dx3)*cosecl1;
01174         
01175         dCCurvTransform_(2,1) = cost*(v11*u21 + v12*u22          ) +
01176           sint*(hv1*u21 + hv2*u22          ) +
01177           omcost*(hn1*v11 + hn2*v12 + hn3*v13) *
01178           (hn1*u21 + hn2*u22          ) +
01179           anu*(-sint*(v11*t21 + v12*t22 + v13*t23) +
01180                omcost*(v11*an1 + v12*an2 + v13*an3) -
01181                tmsint*gamma*(hn1*v11 + hn2*v12 + hn3*v13) );
01182         dCCurvTransform_(2,1) *= cosecl1;
01183         
01184         dCCurvTransform_(2,2) = cost*(u11*u21 + u12*u22          ) +
01185           sint*(hu1*u21 + hu2*u22          ) +
01186           omcost*(hn1*u11 + hn2*u12          ) *
01187           (hn1*u21 + hn2*u22          ) +
01188           anu*(-sint*(u11*t21 + u12*t22          ) +
01189                omcost*(u11*an1 + u12*an2          ) -
01190                tmsint*gamma*(hn1*u11 + hn2*u12          ) );
01191         dCCurvTransform_(2,2) *= cosecl1*cosl0;
01192         
01193         // Commented out in part for reproducibility purposes: these terms are zero in cart->curv 
01194         // dCCurvTransform_(2,3) = -q*anu*(u11*t21 + u12*t22          )*cosecl1;
01195         //why would lambdaNext depend explicitely on initial position ? any arbitrary init point can be chosen not 
01196         // affecting the final state's momentum direction ... is this the field gradient in curvilinear coord?
01197         // dCCurvTransform_(2,4) = -q*anu*(v11*t21 + v12*t22 + v13*t23)*cosecl1;
01198         
01199         //   yt
01200         
01201         double pp = 1./qbp;
01202         // (SK) these terms seem to consistently have a sign opp from private derivation
01203         dCCurvTransform_(3,0) = - pp*(u21*dx1 + u22*dx2            ); //NB: modified from the original: changed the sign
01204         dCCurvTransform_(4,0) = - pp*(v21*dx1 + v22*dx2 + v23*dx3);  
01205         
01206         
01207         dCCurvTransform_(3,1) = (sint*(v11*u21 + v12*u22          ) +
01208                                  omcost*(hv1*u21 + hv2*u22          ) +
01209                                  tmsint*(hn1*u21 + hn2*u22          ) *
01210                                  (hn1*v11 + hn2*v12 + hn3*v13))/q;
01211         
01212         dCCurvTransform_(3,2) = (sint*(u11*u21 + u12*u22          ) +
01213                                  omcost*(hu1*u21 + hu2*u22          ) +
01214                                  tmsint*(hn1*u21 + hn2*u22          ) *
01215                                  (hn1*u11 + hn2*u12          ))*cosl0/q;
01216         
01217         dCCurvTransform_(3,3) = (u11*u21 + u12*u22          );
01218         
01219         dCCurvTransform_(3,4) = (v11*u21 + v12*u22          );
01220         
01221         //   zt
01222         
01223         dCCurvTransform_(4,1) = (sint*(v11*v21 + v12*v22 + v13*v23) +
01224                                  omcost*(hv1*v21 + hv2*v22 + hv3*v23) +
01225                                  tmsint*(hn1*v21 + hn2*v22 + hn3*v23) *
01226                                  (hn1*v11 + hn2*v12 + hn3*v13))/q;
01227         
01228         dCCurvTransform_(4,2) = (sint*(u11*v21 + u12*v22          ) +
01229                                  omcost*(hu1*v21 + hu2*v22 + hu3*v23) +
01230                                  tmsint*(hn1*v21 + hn2*v22 + hn3*v23) *
01231                                  (hn1*u11 + hn2*u12          ))*cosl0/q;
01232         
01233         dCCurvTransform_(4,3) = (u11*v21 + u12*v22          );
01234         
01235         dCCurvTransform_(4,4) = (v11*v21 + v12*v22 + v13*v23);
01236         // end of TRPRFN
01237       }
01238     
01239       if (debug_){
01240         Basis rep; setRep(rep, tauNext);
01241         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"rep X: "<<rep.lX<<" "<<rep.lX.mag()
01242                          <<"\t Y: "<<rep.lY<<" "<<rep.lY.mag()
01243                          <<"\t Z: "<<rep.lZ<<" "<<rep.lZ.mag();
01244       }
01245       //mind the sign of dS and dP (dS*dP < 0 allways)
01246       //covariance should grow no matter which direction you propagate
01247       //==> take abs values.
01248       //reset not needed: fill all below  svCurrent.matDCov *= 0.;
01249       double mulRR = theta02*dS*dS/3.;
01250       double mulRP = theta02*fabs(dS)*p0/2.;
01251       double mulPP = theta02*p0*p0;
01252       double losPP = dP*dP*1.6/fabs(dS)*(1.0 + p0*1e-3);
01253       //another guess .. makes sense for 1 cm steps 2./dS == 2 [cm] / dS [cm] at low pt
01254       //double it by 1TeV
01255       //not gaussian anyways
01256       // derived from the fact that sigma_p/eLoss ~ 0.08 after ~ 200 steps
01257 
01258       //curvilinear
01259       double sinThetaInv = tauNextPerpInv;
01260       double p0Mat = p0+ 0.5*dP; // FIXME change this to p0 after it's clear that there's agreement in everything else
01261       double p0Mat2 = p0Mat*p0Mat;
01262       // with 6x6 formulation
01263       svCurrent.matDCovCurv*=0;
01264       
01265       svCurrent.matDCovCurv(0,0) = losPP/(p0Mat2*p0Mat2);
01266       //      svCurrent.matDCovCurv(0,1) = 0;
01267       //      svCurrent.matDCovCurv(0,2) = 0;
01268       //      svCurrent.matDCovCurv(0,3) = 0;
01269       //      svCurrent.matDCovCurv(0,4) = 0;
01270 
01271       //      svCurrent.matDCovCurv(1,0) = 0;
01272       svCurrent.matDCovCurv(1,1) = mulPP/p0Mat2;
01273       //      svCurrent.matDCovCurv(1,2) = 0;
01274       //      svCurrent.matDCovCurv(1,3) = 0;
01275       svCurrent.matDCovCurv(1,4) = mulRP/p0Mat;
01276 
01277       //      svCurrent.matDCovCurv(2,0) = 0;
01278       //      svCurrent.matDCovCurv(2,1) = 0;
01279       svCurrent.matDCovCurv(2,2) = mulPP/p0Mat2*(sinThetaInv*sinThetaInv);
01280       svCurrent.matDCovCurv(2,3) = mulRP/p0Mat*sinThetaInv;
01281       //      svCurrent.matDCovCurv(2,4) = 0;
01282 
01283       //      svCurrent.matDCovCurv(3,0) = 0;
01284       //      svCurrent.matDCovCurv(3,1) = 0;
01285       svCurrent.matDCovCurv(3,2) = mulRP/p0Mat*sinThetaInv;
01286       //      svCurrent.matDCovCurv(3,0) = 0;
01287       svCurrent.matDCovCurv(3,3) = mulRR;
01288       //      svCurrent.matDCovCurv(3,4) = 0;
01289 
01290       //      svCurrent.matDCovCurv(4,0) = 0;
01291       svCurrent.matDCovCurv(4,1) = mulRP/p0Mat;
01292       //      svCurrent.matDCovCurv(4,2) = 0;
01293       //      svCurrent.matDCovCurv(4,3) = 0;
01294       svCurrent.matDCovCurv(4,4) = mulRR;
01295     }
01296     break;
01297   }
01298   default:
01299     break;
01300   }
01301   
01302   double pMag = curP;//svCurrent.p3.mag();
01303   
01304   if (dir == oppositeToMomentum) dP = fabs(dP);
01305   else if( dP < 0) { //the case of negative dP
01306     dP = -dP > pMag ? -pMag+1e-5 : dP;
01307   }
01308   
01309   getNextState(svCurrent, svNext, dP, tauNext, drVec, dS, dS/radX0,
01310                dCCurvTransform_);
01311   return true;
01312 }
01313 
01314 double SteppingHelixPropagator::getDeDx(const SteppingHelixPropagator::StateInfo& sv, 
01315                                         double& dEdXPrime, double& radX0,
01316                                         MatBounds& rzLims) const{
01317   radX0 = 1.e24;
01318   dEdXPrime = 0.;
01319   rzLims = MatBounds();
01320   if (noMaterialMode_) return 0;
01321 
01322   double dEdx = 0.;
01323 
01324   double lR = sv.r3.perp();
01325   double lZ = fabs(sv.r3.z());
01326 
01327   //assume "Iron" .. seems to be quite the same for brass/iron/PbW04
01328   //good for Fe within 3% for 0.2 GeV to 10PeV
01329 
01330   double dEdX_HCal = 0.95; //extracted from sim
01331   double dEdX_ECal = 0.45;
01332   double dEdX_coil = 0.35; //extracted from sim .. closer to 40% in fact
01333   double dEdX_Fe =   1;
01334   double dEdX_MCh =  0.053; //chambers on average
01335   double dEdX_Trk =  0.0114;
01336   double dEdX_Air =  2E-4;
01337   double dEdX_Vac =  0.0;
01338 
01339   double radX0_HCal = 1.44/0.8; //guessing
01340   double radX0_ECal = 0.89/0.7;
01341   double radX0_coil = 4.; //
01342   double radX0_Fe =   1.76;
01343   double radX0_MCh =  1e3; //
01344   double radX0_Trk =  320.;
01345   double radX0_Air =  3.e4;
01346   double radX0_Vac =  3.e9; //"big" number for vacuum
01347 
01348 
01349   //not all the boundaries are set below: this will be a default
01350   if (! (lR < 380 && lZ < 785)){
01351     if (lZ > 785 ) rzLims = MatBounds(0, 1e4, 785, 1e4);
01352     if (lZ < 785 ) rzLims = MatBounds(380, 1e4, 0, 785);
01353   }
01354 
01355   //this def makes sense assuming we don't jump into endcap volume from the other z-side in one step
01356   //also, it is a positive shift only (at least for now): can't move ec further into the detector
01357   double ecShift = sv.r3.z() > 0 ? fabs(ecShiftPos_) : fabs(ecShiftNeg_); 
01358 
01359   //this should roughly figure out where things are 
01360   //(numbers taken from Fig1.1.2 TDR and from geom xmls)
01361   if (lR < 2.9){ //inside beampipe
01362     dEdx = dEdX_Vac; radX0 = radX0_Vac;
01363     rzLims = MatBounds(0, 2.9, 0, 1E4);
01364   }
01365   else if (lR < 129){
01366     if (lZ < 294){ 
01367       rzLims = MatBounds(2.9, 129, 0, 294); //tracker boundaries
01368       dEdx = dEdX_Trk; radX0 = radX0_Trk; 
01369       //somewhat empirical formula that ~ matches the average if going from 0,0,0
01370       //assuming "uniform" tracker material
01371       //doesn't really track material layer to layer
01372       double lEtaDet = fabs(sv.r3.eta());
01373       double scaleRadX = lEtaDet > 1.5 ? 0.7724 : 1./cosh(0.5*lEtaDet);//sin(2.*atan(exp(-0.5*lEtaDet)));
01374       scaleRadX *= scaleRadX;
01375       if (lEtaDet > 2 && lZ > 20) scaleRadX *= (lEtaDet-1.);
01376       if (lEtaDet > 2.5 && lZ > 20) scaleRadX *= (lEtaDet-1.);
01377       radX0 *= scaleRadX;
01378     }
01379     //endcap part begins here
01380     else if ( lZ < 294 + ecShift ){
01381       //gap in front of EE here, piece inside 2.9<R<129
01382       rzLims = MatBounds(2.9, 129, 294, 294 + ecShift); 
01383       dEdx = dEdX_Air; radX0 = radX0_Air;
01384     }
01385     else if (lZ < 372 + ecShift){ 
01386       rzLims = MatBounds(2.9, 129, 294 + ecShift, 372 + ecShift); //EE here, piece inside 2.9<R<129
01387       dEdx = dEdX_ECal; radX0 = radX0_ECal; 
01388     }//EE averaged out over a larger space
01389     else if (lZ < 398 + ecShift){
01390       rzLims = MatBounds(2.9, 129, 372 + ecShift, 398 + ecShift); //whatever goes behind EE 2.9<R<129 is air up to Z=398
01391       dEdx = dEdX_HCal*0.05; radX0 = radX0_Air; 
01392     }//betw EE and HE
01393     else if (lZ < 555 + ecShift){ 
01394       rzLims = MatBounds(2.9, 129, 398 + ecShift, 555 + ecShift); //HE piece 2.9<R<129; 
01395       dEdx = dEdX_HCal*0.96; radX0 = radX0_HCal/0.96; 
01396     } //HE calor abit less dense
01397     else {
01398       //      rzLims = MatBounds(2.9, 129, 555, 785);
01399       // set the boundaries first: they serve as stop-points here
01400       // the material is set below
01401       if (lZ < 568  + ecShift) rzLims = MatBounds(2.9, 129, 555 + ecShift, 568 + ecShift); //a piece of HE support R<129, 555<Z<568
01402       else if (lZ < 625 + ecShift){
01403         if (lR < 85 + ecShift) rzLims = MatBounds(2.9, 85, 568 + ecShift, 625 + ecShift); 
01404         else rzLims = MatBounds(85, 129, 568 + ecShift, 625 + ecShift);
01405       } else if (lZ < 785 + ecShift) rzLims = MatBounds(2.9, 129, 625 + ecShift, 785 + ecShift);
01406       else if (lZ < 1500 + ecShift)  rzLims = MatBounds(2.9, 129, 785 + ecShift, 1500 + ecShift);
01407       else rzLims = MatBounds(2.9, 129, 1500 + ecShift, 1E4);
01408 
01409       //iron .. don't care about no material in front of HF (too forward)
01410       if (! (lZ > 568 + ecShift && lZ < 625 + ecShift && lR > 85 ) // HE support 
01411           && ! (lZ > 785 + ecShift && lZ < 850 + ecShift && lR > 118)) {dEdx = dEdX_Fe; radX0 = radX0_Fe; }
01412       else  { dEdx = dEdX_MCh; radX0 = radX0_MCh; } //ME at eta > 2.2
01413     }
01414   }
01415   else if (lR < 287){
01416     if (lZ < 372 + ecShift && lR < 177){ // 129<<R<177
01417       if (lZ < 304) rzLims = MatBounds(129, 177, 0, 304); //EB  129<<R<177 0<Z<304
01418       else if (lZ < 343){ // 129<<R<177 304<Z<343
01419         if (lR < 135 ) rzLims = MatBounds(129, 135, 304, 343);// tk piece 129<<R<135 304<Z<343
01420         else if (lR < 172 ){ //
01421           if (lZ < 311 ) rzLims = MatBounds(135, 172, 304, 311);
01422           else rzLims = MatBounds(135, 172, 311, 343);
01423         } else {
01424           if (lZ < 328) rzLims = MatBounds(172, 177, 304, 328);
01425           else rzLims = MatBounds(172, 177, 328, 343);
01426         }
01427       }
01428       else if ( lZ < 343 + ecShift){
01429         rzLims = MatBounds(129, 177, 343, 343 + ecShift); //gap
01430       }
01431       else {
01432         if (lR < 156 ) rzLims = MatBounds(129, 156, 343 + ecShift, 372 + ecShift);
01433         else if ( (lZ - 343 - ecShift) > (lR - 156)*1.38 ) 
01434           rzLims = MatBounds(156, 177, 127.72 + ecShift, 372 + ecShift, atan(1.38), 0);
01435         else rzLims = MatBounds(156, 177, 343 + ecShift, 127.72 + ecShift, 0, atan(1.38));
01436       }
01437 
01438       if (!(lR > 135 && lZ <343 + ecShift && lZ > 304 )
01439           && ! (lR > 156 && lZ < 372 + ecShift  && lZ > 343 + ecShift  && ((lZ-343. - ecShift )< (lR-156.)*1.38)))
01440         {
01441           //the crystals are the same length, but they are not 100% of material
01442           double cosThetaEquiv = 0.8/sqrt(1.+lZ*lZ/lR/lR) + 0.2;
01443           if (lZ > 343) cosThetaEquiv = 1.;
01444           dEdx = dEdX_ECal*cosThetaEquiv; radX0 = radX0_ECal/cosThetaEquiv; 
01445         } //EB
01446       else { 
01447         if ( (lZ > 304 && lZ < 328 && lR < 177 && lR > 135) 
01448              && ! (lZ > 311 && lR < 172) ) {dEdx = dEdX_Fe; radX0 = radX0_Fe; } //Tk_Support
01449         else if ( lZ > 343 && lZ < 343 + ecShift) { dEdx = dEdX_Air; radX0 = radX0_Air; }
01450         else {dEdx = dEdX_ECal*0.2; radX0 = radX0_Air;} //cables go here <-- will be abit too dense for ecShift > 0
01451       }
01452     }
01453     else if (lZ < 554 + ecShift){ // 129<R<177 372<Z<554  AND 177<R<287 0<Z<554
01454       if (lR < 177){ //  129<R<177 372<Z<554
01455         if ( lZ > 372 + ecShift && lZ < 398 + ecShift )rzLims = MatBounds(129, 177, 372 + ecShift, 398 + ecShift); // EE 129<R<177 372<Z<398
01456         else if (lZ < 548 + ecShift) rzLims = MatBounds(129, 177, 398 + ecShift, 548 + ecShift); // HE 129<R<177 398<Z<548
01457         else rzLims = MatBounds(129, 177, 548 + ecShift, 554 + ecShift); // HE gap 129<R<177 548<Z<554
01458       }
01459       else if (lR < 193){ // 177<R<193 0<Z<554
01460         if ((lZ - 307) < (lR - 177.)*1.739) rzLims = MatBounds(177, 193, 0, -0.803, 0, atan(1.739));
01461         else if ( lZ < 389)  rzLims = MatBounds(177, 193, -0.803, 389, atan(1.739), 0.);
01462         else if ( lZ < 389 + ecShift)  rzLims = MatBounds(177, 193, 389, 389 + ecShift); // air insert
01463         else if ( lZ < 548 + ecShift ) rzLims = MatBounds(177, 193, 389 + ecShift, 548 + ecShift);// HE 177<R<193 389<Z<548
01464         else rzLims = MatBounds(177, 193, 548 + ecShift, 554 + ecShift); // HE gap 177<R<193 548<Z<554
01465       }
01466       else if (lR < 264){ // 193<R<264 0<Z<554
01467         double anApex = 375.7278 - 193./1.327; // 230.28695599096
01468         if ( (lZ - 375.7278) < (lR - 193.)/1.327) rzLims = MatBounds(193, 264, 0, anApex, 0, atan(1./1.327)); //HB
01469         else if ( (lZ - 392.7278 ) < (lR - 193.)/1.327) 
01470           rzLims = MatBounds(193, 264, anApex, anApex+17., atan(1./1.327), atan(1./1.327)); // HB-HE gap
01471         else if ( (lZ - 392.7278 - ecShift ) < (lR - 193.)/1.327) 
01472           rzLims = MatBounds(193, 264, anApex+17., anApex+17. + ecShift, atan(1./1.327), atan(1./1.327)); // HB-HE gap air insert
01473         // HE (372,193)-(517,193)-(517,264)-(417.8,264)
01474         else if ( lZ < 517 + ecShift ) rzLims = MatBounds(193, 264, anApex+17. + ecShift, 517 + ecShift, atan(1./1.327), 0);
01475         else if (lZ < 548 + ecShift){ // HE+gap 193<R<264 517<Z<548
01476           if (lR < 246 ) rzLims = MatBounds(193, 246, 517 + ecShift, 548 + ecShift); // HE 193<R<246 517<Z<548
01477           else rzLims = MatBounds(246, 264, 517 + ecShift, 548 + ecShift); // HE gap 246<R<264 517<Z<548
01478         } 
01479         else rzLims = MatBounds(193, 264, 548 + ecShift, 554 + ecShift); // HE gap  193<R<246 548<Z<554
01480       }
01481       else if ( lR < 275){ // 264<R<275 0<Z<554
01482         if (lZ < 433) rzLims = MatBounds(264, 275, 0, 433); //HB
01483         else if (lZ < 554 ) rzLims = MatBounds(264, 275, 433, 554); // HE gap 264<R<275 433<Z<554
01484         else rzLims = MatBounds(264, 275, 554, 554 + ecShift); // HE gap 264<R<275 554<Z<554 air insert
01485       }
01486       else { // 275<R<287 0<Z<554
01487         if (lZ < 402) rzLims = MatBounds(275, 287, 0, 402);// HB 275<R<287 0<Z<402
01488         else if (lZ < 554) rzLims = MatBounds(275, 287, 402, 554);//  //HE gap 275<R<287 402<Z<554
01489         else rzLims = MatBounds(275, 287, 554, 554 + ecShift); //HE gap 275<R<287 554<Z<554 air insert
01490       }
01491 
01492       if ((lZ < 433 || lR < 264) && (lZ < 402 || lR < 275) && (lZ < 517 + ecShift || lR < 246) //notches
01493           //I should've made HE and HF different .. now need to shorten HE to match
01494           && lZ < 548 + ecShift
01495           && ! (lZ < 389 + ecShift && lZ > 335 && lR < 193 ) //not a gap EB-EE 129<R<193
01496           && ! (lZ > 307 && lZ < 335 && lR < 193 && ((lZ - 307) > (lR - 177.)*1.739)) //not a gap 
01497           && ! (lR < 177 && lZ < 398 + ecShift) //under the HE nose
01498           && ! (lR < 264 && lR > 193 && fabs(441.5+0.5*ecShift - lZ + (lR - 269.)/1.327) < (8.5 + ecShift*0.5)) ) //not a gap
01499         { dEdx = dEdX_HCal; radX0 = radX0_HCal; //hcal
01500         }
01501       else if( (lR < 193 && lZ > 389 && lZ < 389 + ecShift ) || (lR > 264 && lR < 287 && lZ > 554 && lZ < 554 + ecShift)
01502                ||(lR < 264 && lR > 193 && fabs(441.5+8.5+0.5*ecShift - lZ + (lR - 269.)/1.327) < ecShift*0.5) )  {
01503         dEdx = dEdX_Air; radX0 = radX0_Air; 
01504       }
01505       else {dEdx = dEdX_HCal*0.05; radX0 = radX0_Air; }//endcap gap
01506     }
01507     //the rest is a tube of endcap volume  129<R<251 554<Z<largeValue
01508     else if (lZ < 564 + ecShift){ // 129<R<287  554<Z<564
01509       if (lR < 251) {
01510         rzLims = MatBounds(129, 251, 554 + ecShift, 564 + ecShift);  // HE support 129<R<251  554<Z<564    
01511         dEdx = dEdX_Fe; radX0 = radX0_Fe; 
01512       }//HE support
01513       else { 
01514         rzLims = MatBounds(251, 287, 554 + ecShift, 564 + ecShift); //HE/ME gap 251<R<287  554<Z<564    
01515         dEdx = dEdX_MCh; radX0 = radX0_MCh; 
01516       }
01517     }
01518     else if (lZ < 625 + ecShift){ // ME/1/1 129<R<287  564<Z<625
01519       rzLims = MatBounds(129, 287, 564 + ecShift, 625 + ecShift);      
01520       dEdx = dEdX_MCh; radX0 = radX0_MCh; 
01521     }
01522     else if (lZ < 785 + ecShift){ //129<R<287  625<Z<785
01523       if (lR < 275) rzLims = MatBounds(129, 275, 625 + ecShift, 785 + ecShift); //YE/1 129<R<275 625<Z<785
01524       else { // 275<R<287  625<Z<785
01525         if (lZ < 720 + ecShift) rzLims = MatBounds(275, 287, 625 + ecShift, 720 + ecShift); // ME/1/2 275<R<287  625<Z<720
01526         else rzLims = MatBounds(275, 287, 720 + ecShift, 785 + ecShift); // YE/1 275<R<287  720<Z<785
01527       }
01528       if (! (lR > 275 && lZ < 720 + ecShift)) { dEdx = dEdX_Fe; radX0 = radX0_Fe; }//iron
01529       else { dEdx = dEdX_MCh; radX0 = radX0_MCh; }
01530     }
01531     else if (lZ < 850 + ecShift){
01532       rzLims = MatBounds(129, 287, 785 + ecShift, 850 + ecShift);
01533       dEdx = dEdX_MCh; radX0 = radX0_MCh; 
01534     }
01535     else if (lZ < 910 + ecShift){
01536       rzLims = MatBounds(129, 287, 850 + ecShift, 910 + ecShift);
01537       dEdx = dEdX_Fe; radX0 = radX0_Fe; 
01538     }//iron
01539     else if (lZ < 975 + ecShift){
01540       rzLims = MatBounds(129, 287, 910 + ecShift, 975 + ecShift);
01541       dEdx = dEdX_MCh; radX0 = radX0_MCh; 
01542     }
01543     else if (lZ < 1000 + ecShift){
01544       rzLims = MatBounds(129, 287, 975 + ecShift, 1000 + ecShift);
01545       dEdx = dEdX_Fe; radX0 = radX0_Fe; 
01546     }//iron
01547     else if (lZ < 1063 + ecShift){
01548       rzLims = MatBounds(129, 287, 1000 + ecShift, 1063 + ecShift);
01549       dEdx = dEdX_MCh; radX0 = radX0_MCh; 
01550     }
01551     else if ( lZ < 1073 + ecShift){
01552       rzLims = MatBounds(129, 287, 1063 + ecShift, 1073 + ecShift);
01553       dEdx = dEdX_Fe; radX0 = radX0_Fe; 
01554     }
01555     else if (lZ < 1E4 )  { 
01556       rzLims = MatBounds(129, 287, 1073 + ecShift, 1E4);
01557       dEdx = dEdX_Air; radX0 = radX0_Air;
01558     }
01559     else { 
01560       
01561       dEdx = dEdX_Air; radX0 = radX0_Air;
01562     }
01563   }
01564   else if (lR <380 && lZ < 667){
01565     if (lZ < 630){
01566       if (lR < 315) rzLims = MatBounds(287, 315, 0, 630); 
01567       else if (lR < 341 ) rzLims = MatBounds(315, 341, 0, 630); //b-field ~linear rapid fall here
01568       else rzLims = MatBounds(341, 380, 0, 630);      
01569     } else rzLims = MatBounds(287, 380, 630, 667);  
01570 
01571     if (lZ < 630) { dEdx = dEdX_coil; radX0 = radX0_coil; }//a guess for the solenoid average
01572     else {dEdx = dEdX_Air; radX0 = radX0_Air; }//endcap gap
01573   }
01574   else {
01575     if (lZ < 667) {
01576       if (lR < 850){
01577         bool isIron = false;
01578         //sanity check in addition to flags
01579         if (useIsYokeFlag_ && useMagVolumes_ && sv.magVol != 0){
01580           isIron = sv.isYokeVol;
01581         } else {
01582           double bMag = sv.bf.mag();
01583           isIron = (bMag > 0.75 && ! (lZ > 500 && lR <500 && bMag < 1.15)
01584                     && ! (lZ < 450 && lR > 420 && bMag < 1.15 ) );
01585         }
01586         //tell the material stepper where mat bounds are
01587         rzLims = MatBounds(380, 850, 0, 667);
01588         if (isIron) { dEdx = dEdX_Fe; radX0 = radX0_Fe; }//iron
01589         else { dEdx = dEdX_MCh; radX0 = radX0_MCh; }
01590       } else {
01591         rzLims = MatBounds(850, 1E4, 0, 667);
01592         dEdx = dEdX_Air; radX0 = radX0_Air; 
01593       }
01594     } 
01595     else if (lR > 750 ){
01596       rzLims = MatBounds(750, 1E4, 667, 1E4);
01597       dEdx = dEdX_Air; radX0 = radX0_Air; 
01598     }
01599     else if (lZ < 667 + ecShift){
01600       rzLims = MatBounds(287, 750, 667, 667 + ecShift);
01601       dEdx = dEdX_Air; radX0 = radX0_Air;       
01602     }
01603     //the rest is endcap piece with 287<R<750 Z>667
01604     else if (lZ < 724 + ecShift){
01605       if (lR < 380 ) rzLims = MatBounds(287, 380, 667 + ecShift, 724 + ecShift); 
01606       else rzLims = MatBounds(380, 750, 667 + ecShift, 724 + ecShift); 
01607       dEdx = dEdX_MCh; radX0 = radX0_MCh; 
01608     }
01609     else if (lZ < 785 + ecShift){
01610       if (lR < 380 ) rzLims = MatBounds(287, 380, 724 + ecShift, 785 + ecShift); 
01611       else rzLims = MatBounds(380, 750, 724 + ecShift, 785 + ecShift); 
01612       dEdx = dEdX_Fe; radX0 = radX0_Fe; 
01613     }//iron
01614     else if (lZ < 850 + ecShift){
01615       rzLims = MatBounds(287, 750, 785 + ecShift, 850 + ecShift); 
01616       dEdx = dEdX_MCh; radX0 = radX0_MCh; 
01617     }
01618     else if (lZ < 910 + ecShift){
01619       rzLims = MatBounds(287, 750, 850 + ecShift, 910 + ecShift); 
01620       dEdx = dEdX_Fe; radX0 = radX0_Fe; 
01621     }//iron
01622     else if (lZ < 975 + ecShift){
01623       rzLims = MatBounds(287, 750, 910 + ecShift, 975 + ecShift); 
01624       dEdx = dEdX_MCh; radX0 = radX0_MCh; 
01625     }
01626     else if (lZ < 1000 + ecShift){
01627       rzLims = MatBounds(287, 750, 975 + ecShift, 1000 + ecShift); 
01628       dEdx = dEdX_Fe; radX0 = radX0_Fe; 
01629     }//iron
01630     else if (lZ < 1063 + ecShift){
01631       if (lR < 360){
01632         rzLims = MatBounds(287, 360, 1000 + ecShift, 1063 + ecShift);
01633         dEdx = dEdX_MCh; radX0 = radX0_MCh; 
01634       } 
01635       //put empty air where me4/2 should be (future)
01636       else {
01637         rzLims = MatBounds(360, 750, 1000 + ecShift, 1063 + ecShift);
01638         dEdx = dEdX_Air; radX0 = radX0_Air;
01639       }
01640     }
01641     else if ( lZ < 1073 + ecShift){
01642       rzLims = MatBounds(287, 750, 1063 + ecShift, 1073 + ecShift);
01643       //this plate does not exist: air
01644       dEdx = dEdX_Air; radX0 = radX0_Air; 
01645     }
01646     else if (lZ < 1E4 )  { 
01647       rzLims = MatBounds(287, 750, 1073 + ecShift, 1E4);
01648       dEdx = dEdX_Air; radX0 = radX0_Air;
01649     }
01650     else {dEdx = dEdX_Air; radX0 = radX0_Air; }//air
01651   }
01652   
01653   //dEdx so far is a relative number (relative to iron)
01654   //scale by what's expected for iron (the function was fit from pdg table)
01655   //0.065 (PDG) --> 0.044 to better match with MPV
01656   double p0 = sv.p3.mag();
01657   double logp0 = log(p0);
01658   double p0powN33 = 0; 
01659   if (p0>3.) {
01660     // p0powN33 = exp(-0.33*logp0); //calculate for p>3GeV
01661     double xx=1./p0; xx=sqrt(sqrt(sqrt(sqrt(xx)))); xx=xx*xx*xx*xx*xx; // this is (p0)**(-5/16), close enough to -0.33
01662     p0powN33 = xx;
01663   }
01664   double dEdX_mat = -(11.4 + 0.96*fabs(logp0+log(2.8)) + 0.033*p0*(1.0 - p0powN33) )*1e-3; 
01665   //in GeV/cm .. 0.8 to get closer to the median or MPV
01666 
01667   dEdXPrime = dEdx == 0 ? 0 : -dEdx*(0.96/p0 + 0.033 - 0.022*p0powN33)*1e-3; //== d(dEdX)/dp
01668   dEdx *= dEdX_mat;
01669 
01670   return dEdx;
01671 }
01672 
01673 
01674 int SteppingHelixPropagator::cIndex_(int ind) const{
01675   int result = ind%MAX_POINTS;  
01676   if (ind != 0 && result == 0){
01677     result = MAX_POINTS;
01678   }
01679   return result;
01680 }
01681 
01682 SteppingHelixPropagator::Result
01683 SteppingHelixPropagator::refToDest(SteppingHelixPropagator::DestType dest, 
01684                                    const SteppingHelixPropagator::StateInfo& sv,
01685                                    const double pars[6], 
01686                                    double& dist, double& tanDist, 
01687                                    PropagationDirection& refDirection,
01688                                    double fastSkipDist) const{
01689   static const std::string metname = "SteppingHelixPropagator";
01690   Result result = SteppingHelixStateInfo::NOT_IMPLEMENTED;
01691 
01692   switch (dest){
01693   case RADIUS_DT:
01694     {
01695       double curR = sv.r3.perp();
01696       dist = pars[RADIUS_P] - curR;
01697       if (fabs(dist) > fastSkipDist){
01698         result = SteppingHelixStateInfo::INACC;
01699         break;
01700       }
01701       double curP2 = sv.p3.mag2();
01702       double curPtPos2 = sv.p3.perp2(); if(curPtPos2< 1e-16) curPtPos2=1e-16;
01703 
01704       double cosDPhiPR = (sv.r3.x()*sv.p3.x()+sv.r3.y()*sv.p3.y());//only the sign is needed cos((sv.r3.deltaPhi(sv.p3)));
01705       refDirection = dist*cosDPhiPR > 0 ?
01706         alongMomentum : oppositeToMomentum;
01707       tanDist = dist*sqrt(curP2/curPtPos2);      
01708       result = SteppingHelixStateInfo::OK;
01709     }
01710     break;
01711   case Z_DT:
01712     {
01713       double curZ = sv.r3.z();
01714       dist = pars[Z_P] - curZ;
01715       if (fabs(dist) > fastSkipDist){
01716         result = SteppingHelixStateInfo::INACC;
01717         break;
01718       }
01719       double curP = sv.p3.mag();
01720       refDirection = sv.p3.z()*dist > 0. ?
01721         alongMomentum : oppositeToMomentum;
01722       tanDist = dist/sv.p3.z()*curP;
01723       result = SteppingHelixStateInfo::OK;
01724     }
01725     break;
01726   case PLANE_DT:
01727     {
01728       Point rPlane(pars[0], pars[1], pars[2]);
01729       Vector nPlane(pars[3], pars[4], pars[5]);
01730       
01731 
01732       // unfortunately this doesn't work: the numbers are too large
01733       //      bool pVertical = fabs(pars[5])>0.9999;
01734       //      double dRDotN = pVertical? (sv.r3.z() - rPlane.z())*nPlane.z() :(sv.r3 - rPlane).dot(nPlane);
01735       double dRDotN = (sv.r3.x()-rPlane.x())*nPlane.x() + (sv.r3.y()-rPlane.y())*nPlane.y() + (sv.r3.z()-rPlane.z())*nPlane.z();//(sv.r3 - rPlane).dot(nPlane);
01736 
01737       dist = fabs(dRDotN);
01738       if (dist > fastSkipDist){
01739         result = SteppingHelixStateInfo::INACC;
01740         break;
01741       }
01742       double curP = sv.p3.mag();
01743       double p0 = curP;
01744       double p0Inv = 1./p0;
01745       Vector tau(sv.p3); tau *=p0Inv; 
01746       double tN =  tau.dot(nPlane);
01747       refDirection = tN*dRDotN < 0. ?
01748         alongMomentum : oppositeToMomentum;
01749       double b0 = sv.bf.mag();
01750       if (fabs(tN)>1e-24){
01751         tanDist = -dRDotN/tN;
01752       } else {
01753         tN = 1e-24;
01754         if (fabs(dRDotN)>1e-24) tanDist = 1e6;
01755         else tanDist = 1;
01756       }
01757       if (fabs(tanDist) > 1e4) tanDist = 1e4;
01758       if (b0>1.5e-6){
01759         double b0Inv = 1./b0;
01760         double tNInv = 1./tN;
01761         Vector bHat(sv.bf); bHat *=b0Inv;
01762         double bHatN = bHat.dot(nPlane);
01763         double cosPB = bHat.dot(tau);
01764         double kVal = 2.99792458e-3*sv.q*p0Inv*b0;
01765         double aVal = tanDist*kVal;
01766         Vector lVec = bHat.cross(tau);
01767         double bVal = lVec.dot(nPlane)*tNInv;
01768         if (fabs(aVal*bVal)< 0.3){
01769           double cVal = 1. - bHatN*cosPB*tNInv; // - sv.bf.cross(lVec).dot(nPlane)*b0Inv*tNInv; //1- bHat_n*bHat_tau/tau_n;
01770           double aacVal = cVal*aVal*aVal;
01771           if (fabs(aacVal)<1){
01772             double tanDCorr = bVal/2. + (bVal*bVal/2. + cVal/6)*aVal;
01773             tanDCorr *= aVal;
01774             //+ (-bVal/24. + 0.625*bVal*bVal*bVal + 5./12.*bVal*cVal)*aVal*aVal*aVal
01775             if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<tanDist<<" vs "
01776                                          <<tanDist*(1.+tanDCorr)<<" corr "<<tanDist*tanDCorr<<std::endl;
01777             tanDist *= (1.+tanDCorr);
01778           } else {
01779             if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"AACVal "<< fabs(aacVal)
01780                                          <<" = "<<aVal<<"**2 * "<<cVal<<" too large:: will not converge"<<std::endl;
01781           }
01782         } else {
01783           if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"ABVal "<< fabs(aVal*bVal)
01784                                        <<" = "<<aVal<<" * "<<bVal<<" too large:: will not converge"<<std::endl;
01785         }
01786       }
01787       result = SteppingHelixStateInfo::OK;
01788     }
01789     break;
01790   case CONE_DT:
01791     {
01792       Point cVertex(pars[0], pars[1], pars[2]);
01793       if (cVertex.perp2() < 1e-10){
01794         //assumes the cone axis/vertex is along z
01795         Vector relV3 = sv.r3 - cVertex;
01796         double relV3mag = relV3.mag();
01797         //      double relV3Theta = relV3.theta();
01798         double theta(pars[3]);
01799         //      double dTheta = theta-relV3Theta;
01800         double sinTheta = sin(theta);
01801         double cosTheta = cos(theta);
01802         double cosV3Theta = relV3.z()/relV3mag;
01803         if (cosV3Theta>1) cosV3Theta=1;
01804         if (cosV3Theta<-1) cosV3Theta=-1;
01805         double sinV3Theta = sqrt(1.-cosV3Theta*cosV3Theta);
01806 
01807         double sinDTheta = sinTheta*cosV3Theta - cosTheta*sinV3Theta;//sin(dTheta);
01808         double cosDTheta = cosTheta*cosV3Theta + sinTheta*sinV3Theta;//cos(dTheta);
01809         bool isInside = sinTheta > sinV3Theta  && cosTheta*cosV3Theta > 0;
01810         dist = isInside || cosDTheta > 0 ? relV3mag*sinDTheta : relV3mag;
01811         if (fabs(dist) > fastSkipDist){
01812           result = SteppingHelixStateInfo::INACC;
01813           break;
01814         }
01815 
01816         double relV3phi=relV3.phi();
01817         double normPhi = isInside ? 
01818           Geom::pi() + relV3phi : relV3phi;
01819         double normTheta = theta > Geom::pi()/2. ?
01820           ( isInside ? 1.5*Geom::pi()  - theta : theta - Geom::pi()/2. )
01821           : ( isInside ? Geom::pi()/2. - theta : theta + Geom::pi()/2 );
01822         //this is a normVector from the cone to the point
01823         Vector norm; norm.setRThetaPhi(fabs(dist), normTheta, normPhi);
01824         double curP = sv.p3.mag(); double cosp3theta = sv.p3.z()/curP; 
01825         if (cosp3theta>1) cosp3theta=1;
01826         if (cosp3theta<-1) cosp3theta=-1;
01827         double sineConeP = sinTheta*cosp3theta - cosTheta*sqrt(1.-cosp3theta*cosp3theta);
01828 
01829         double sinSolid = norm.dot(sv.p3)/(fabs(dist)*curP);
01830         tanDist = fabs(sinSolid) > fabs(sineConeP) ? dist/fabs(sinSolid) : dist/fabs(sineConeP);
01831 
01832 
01833         refDirection = sinSolid > 0 ?
01834           oppositeToMomentum : alongMomentum;
01835         if (debug_){
01836           LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Cone pars: theta "<<theta
01837                            <<" normTheta "<<norm.theta()
01838                            <<" p3Theta "<<sv.p3.theta()
01839                            <<" sinD: "<< sineConeP
01840                            <<" sinSolid "<<sinSolid;
01841         }
01842         if (debug_){
01843           LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"refToDest:toCone the point is "
01844                            <<(isInside? "in" : "out")<<"side the cone"
01845                            <<std::endl;
01846         }
01847       }
01848       result = SteppingHelixStateInfo::OK;
01849     }
01850     break;
01851     //   case CYLINDER_DT:
01852     //     break;
01853   case PATHL_DT:
01854     {
01855       double curS = fabs(sv.path());
01856       dist = pars[PATHL_P] - curS;
01857       if (fabs(dist) > fastSkipDist){
01858         result = SteppingHelixStateInfo::INACC;
01859         break;
01860       }
01861       refDirection = pars[PATHL_P] > 0 ? 
01862         alongMomentum : oppositeToMomentum;
01863       tanDist = dist;
01864       result = SteppingHelixStateInfo::OK;
01865     }
01866     break;
01867   case POINT_PCA_DT:
01868     {
01869       Point pDest(pars[0], pars[1], pars[2]);
01870       double curP = sv.p3.mag();
01871       dist = (sv.r3 - pDest).mag()+ 1e-24;//add a small number to avoid 1/0
01872       tanDist = (sv.r3.dot(sv.p3) - pDest.dot(sv.p3))/curP;
01873       //account for bending in magnetic field (quite approximate)
01874       double b0 = sv.bf.mag();
01875       if (b0>1.5e-6){
01876         double p0 = curP;
01877         double kVal = 2.99792458e-3*sv.q/p0*b0;
01878         double aVal = fabs(dist*kVal);
01879         tanDist *= 1./(1.+ aVal);
01880         if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"corrected by aVal "<<aVal<<" to "<<tanDist;
01881       }
01882       refDirection = tanDist < 0 ?
01883         alongMomentum : oppositeToMomentum;
01884       result = SteppingHelixStateInfo::OK;
01885     }
01886     break;
01887   case LINE_PCA_DT:
01888     {
01889       Point rLine(pars[0], pars[1], pars[2]);
01890       Vector dLine(pars[3], pars[4], pars[5]);
01891       dLine = (dLine - rLine);
01892       dLine *= 1./dLine.mag(); if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"dLine "<<dLine;
01893 
01894       Vector dR = sv.r3 - rLine;
01895       double curP = sv.p3.mag();
01896       Vector dRPerp = dR - dLine*(dR.dot(dLine)); 
01897       if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"dRperp "<<dRPerp;
01898 
01899       dist = dRPerp.mag() + 1e-24;//add a small number to avoid 1/0
01900       tanDist = dRPerp.dot(sv.p3)/curP;
01901       //angle wrt line
01902       double cosAlpha2 = dLine.dot(sv.p3)/curP; cosAlpha2 *= cosAlpha2;
01903       tanDist *= 1./sqrt(fabs(1.-cosAlpha2)+1e-96);
01904       //correct for dPhi in magnetic field: this isn't made quite right here 
01905       //(the angle between the line and the trajectory plane is neglected .. conservative)
01906       double b0 = sv.bf.mag();
01907       if (b0>1.5e-6){
01908         double p0 = curP;
01909         double kVal = 2.99792458e-3*sv.q/p0*b0;
01910         double aVal = fabs(dist*kVal);
01911         tanDist *= 1./(1.+ aVal);
01912         if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"corrected by aVal "<<aVal<<" to "<<tanDist;
01913       }
01914       refDirection = tanDist < 0 ?
01915         alongMomentum : oppositeToMomentum;
01916       result = SteppingHelixStateInfo::OK;
01917     }
01918     break;
01919   default:
01920     {
01921       //some large number
01922       dist = 1e12;
01923       tanDist = 1e12;
01924       refDirection = anyDirection;
01925       result = SteppingHelixStateInfo::NOT_IMPLEMENTED;
01926     }
01927     break;
01928   }
01929 
01930   double tanDistConstrained = tanDist;
01931   if (fabs(tanDist) > 4.*fabs(dist) ) tanDistConstrained *= tanDist == 0 ? 0 : fabs(dist/tanDist*4.);
01932 
01933   if (debug_){
01934     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"refToDest input: dest"<<dest<<" pars[]: ";
01935     for (int i = 0; i < 6; i++){
01936       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<", "<<i<<" "<<pars[i];
01937     }
01938     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<std::endl;
01939     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"refToDest output: "
01940                      <<"\t dist" << dist
01941                      <<"\t tanDist "<< tanDist      
01942                      <<"\t tanDistConstr "<< tanDistConstrained      
01943                      <<"\t refDirection "<< refDirection
01944                      <<std::endl;
01945   }
01946   tanDist = tanDistConstrained;
01947 
01948   return result;
01949 }
01950 
01951 SteppingHelixPropagator::Result
01952 SteppingHelixPropagator::refToMagVolume(const SteppingHelixPropagator::StateInfo& sv,
01953                                         PropagationDirection dir,
01954                                         double& dist, double& tanDist,
01955                                         double fastSkipDist, bool expectNewMagVolume, double maxStep) const{
01956 
01957   static const std::string metname = "SteppingHelixPropagator";
01958   Result result = SteppingHelixStateInfo::NOT_IMPLEMENTED;
01959   const MagVolume* cVol = sv.magVol;
01960 
01961   if (cVol == 0) return result;
01962   const std::vector<VolumeSide>& cVolFaces(cVol->faces());
01963 
01964   double distToFace[6] = {0,0,0,0,0,0};
01965   double tanDistToFace[6] = {0,0,0,0,0,0};
01966   PropagationDirection refDirectionToFace[6] = {anyDirection, anyDirection, anyDirection, anyDirection, anyDirection, anyDirection};
01967   Result resultToFace[6] = {result, result, result, result, result, result};
01968   int iFDest = -1;
01969   int iDistMin = -1;
01970   
01971   unsigned int iFDestSorted[6] = {0,0,0,0,0,0};
01972   unsigned int nDestSorted =0;
01973   unsigned int nearParallels = 0;
01974 
01975   double curP = sv.p3.mag();
01976 
01977   if (debug_){
01978     LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Trying volume "<<DDSolidShapesName::name(cVol->shapeType())
01979                      <<" with "<<cVolFaces.size()<<" faces"<<std::endl;
01980   }
01981 
01982   unsigned int nFaces = cVolFaces.size();
01983   for (unsigned int iFace = 0; iFace < nFaces; ++iFace){
01984     if (iFace > 5){
01985       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Too many faces"<<std::endl;
01986     }
01987     if (debug_){
01988       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Start with face "<<iFace<<std::endl;
01989     }
01990 //     const Plane* cPlane = dynamic_cast<const Plane*>(&cVolFaces[iFace].surface());
01991 //     const Cylinder* cCyl = dynamic_cast<const Cylinder*>(&cVolFaces[iFace].surface());
01992 //     const Cone* cCone = dynamic_cast<const Cone*>(&cVolFaces[iFace].surface());
01993     const Surface* cPlane = 0; //only need to know the loc->glob transform
01994     const Cylinder* cCyl = 0;
01995     const Cone* cCone = 0;
01996     if (typeid(cVolFaces[iFace].surface()) == typeid(const Plane&)){
01997       cPlane = &cVolFaces[iFace].surface();
01998     } else if (typeid(cVolFaces[iFace].surface()) == typeid(const Cylinder&)){
01999       cCyl = dynamic_cast<const Cylinder*>(&cVolFaces[iFace].surface());
02000     } else if (typeid(cVolFaces[iFace].surface()) == typeid(const Cone&)){
02001       cCone = dynamic_cast<const Cone*>(&cVolFaces[iFace].surface());
02002     } else {
02003       edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Could not cast a volume side surface to a known type"<<std::endl;
02004     }
02005     
02006     if (debug_){
02007       if (cPlane!=0) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The face is a plane at "<<cPlane<<std::endl;
02008       if (cCyl!=0) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The face is a cylinder at "<<cCyl<<std::endl;
02009     }
02010 
02011     double pars[6] = {0,0,0,0,0,0};
02012     DestType dType = UNDEFINED_DT;
02013     if (cPlane != 0){
02014       Point rPlane(cPlane->position().x(),cPlane->position().y(),cPlane->position().z());
02015       // = cPlane->toGlobal(LocalVector(0,0,1.)); nPlane = nPlane.unit();
02016       Vector nPlane(cPlane->rotation().zx(), cPlane->rotation().zy(), cPlane->rotation().zz()); nPlane /= nPlane.mag();
02017       
02018       if (sv.r3.z() < 0 || !vbField_->isZSymmetric() ){
02019         pars[0] = rPlane.x(); pars[1] = rPlane.y(); pars[2] = rPlane.z();
02020         pars[3] = nPlane.x(); pars[4] = nPlane.y(); pars[5] = nPlane.z();
02021       } else {
02022         pars[0] = rPlane.x(); pars[1] = rPlane.y(); pars[2] = -rPlane.z();
02023         pars[3] = nPlane.x(); pars[4] = nPlane.y(); pars[5] = -nPlane.z();
02024       }
02025       dType = PLANE_DT;
02026     } else if (cCyl != 0){
02027       if (debug_){
02028         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Cylinder at "<<cCyl->position()
02029                          <<" rorated by "<<cCyl->rotation()
02030                          <<std::endl;
02031       }
02032       pars[RADIUS_P] = cCyl->radius();
02033       dType = RADIUS_DT;
02034     } else if (cCone != 0){
02035       if (debug_){
02036         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Cone at "<<cCone->position()
02037                          <<" rorated by "<<cCone->rotation()
02038                          <<" vertex at "<<cCone->vertex()
02039                          <<" angle of "<<cCone->openingAngle()
02040                          <<std::endl;
02041       }
02042       if (sv.r3.z() < 0 || !vbField_->isZSymmetric()){
02043         pars[0] = cCone->vertex().x(); pars[1] = cCone->vertex().y(); 
02044         pars[2] = cCone->vertex().z();
02045         pars[3] = cCone->openingAngle();
02046       } else {
02047         pars[0] = cCone->vertex().x(); pars[1] = cCone->vertex().y(); 
02048         pars[2] = -cCone->vertex().z();
02049         pars[3] = Geom::pi() - cCone->openingAngle();
02050       }
02051       dType = CONE_DT;
02052     } else {
02053       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Unknown surface"<<std::endl;
02054       resultToFace[iFace] = SteppingHelixStateInfo::UNDEFINED;
02055       continue;
02056     }
02057     resultToFace[iFace] = 
02058       refToDest(dType, sv, pars, 
02059                 distToFace[iFace], tanDistToFace[iFace], refDirectionToFace[iFace], fastSkipDist);    
02060     
02061     
02062     if (resultToFace[iFace] != SteppingHelixStateInfo::OK){
02063       if (resultToFace[iFace] == SteppingHelixStateInfo::INACC) result = SteppingHelixStateInfo::INACC;
02064     }
02065 
02066 
02067       
02068     //keep those in right direction for later use
02069     if (resultToFace[iFace] == SteppingHelixStateInfo::OK){
02070       double invDTFPosiF = 1./(1e-32+fabs(tanDistToFace[iFace]));
02071       double dSlope = fabs(distToFace[iFace])*invDTFPosiF;
02072       double maxStepL = maxStep> 100 ? 100 : maxStep; if (maxStepL < 10) maxStepL = 10.;
02073       bool isNearParallel = fabs(tanDistToFace[iFace]) + 100.*curP*dSlope < maxStepL //
02074         //a better choice is to use distance to next check of mag volume instead of 100cm; the last is for ~1.5arcLength(4T)+tandistance< maxStep
02075         && dSlope < 0.15 ; //
02076       if (refDirectionToFace[iFace] == dir || isNearParallel){
02077         if (isNearParallel) nearParallels++;
02078         iFDestSorted[nDestSorted] = iFace;
02079         nDestSorted++;
02080       }
02081     }
02082     
02083     //pick a shortest distance here (right dir only for now)
02084     if (refDirectionToFace[iFace] == dir){
02085       if (iDistMin == -1) iDistMin = iFace;
02086       else if (fabs(distToFace[iFace]) < fabs(distToFace[iDistMin])) iDistMin = iFace;
02087     }
02088     if (debug_) 
02089       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<cVol<<" "<<iFace<<" "
02090                        <<tanDistToFace[iFace]<<" "<<distToFace[iFace]<<" "<<refDirectionToFace[iFace]<<" "<<dir<<std::endl;
02091     
02092   }
02093   
02094   for (unsigned int i = 0;i<nDestSorted; ++i){
02095     unsigned int iMax = nDestSorted-i-1;
02096     for (unsigned int j=0;j<nDestSorted-i; ++j){
02097       if (fabs(tanDistToFace[iFDestSorted[j]]) > fabs(tanDistToFace[iFDestSorted[iMax]]) ){
02098         iMax = j;
02099       }
02100     }
02101     unsigned int iTmp = iFDestSorted[nDestSorted-i-1];
02102     iFDestSorted[nDestSorted-i-1] = iFDestSorted[iMax];
02103     iFDestSorted[iMax] = iTmp;
02104   }
02105 
02106   if (debug_){
02107     for (unsigned int i=0;i<nDestSorted;++i){
02108       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<cVol<<" "<<i<<" "<<iFDestSorted[i]<<" "<<tanDistToFace[iFDestSorted[i]]<<std::endl;
02109     }
02110   }
02111 
02112   //now go from the shortest to the largest distance hoping to get a point in the volume.
02113   //other than in case of a near-parallel travel this should stop after the first try
02114   
02115   for (unsigned int i=0; i<nDestSorted;++i){
02116     iFDest = iFDestSorted[i];
02117 
02118     double sign = dir == alongMomentum ? 1. : -1.;
02119     Point gPointEst(sv.r3);
02120     Vector lDelta(sv.p3); lDelta *= sign/curP*fabs(distToFace[iFDest]);
02121     gPointEst += lDelta;
02122     if (debug_){
02123       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Linear est point "<<gPointEst
02124                        <<" for iFace "<<iFDest<<std::endl;
02125     }
02126     GlobalPoint gPointEstNegZ(gPointEst.x(), gPointEst.y(), -fabs(gPointEst.z()));
02127     GlobalPoint gPointEstNorZ(gPointEst.x(), gPointEst.y(), gPointEst.z() );
02128     if (  (vbField_->isZSymmetric() && cVol->inside(gPointEstNegZ))
02129           || ( !vbField_->isZSymmetric() && cVol->inside(gPointEstNorZ) )  ){
02130       if (debug_){
02131         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is inside the volume"<<std::endl;
02132       }
02133       
02134       result = SteppingHelixStateInfo::OK;
02135       dist = distToFace[iFDest];
02136       tanDist = tanDistToFace[iFDest];
02137       if (debug_){
02138         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got a point near closest boundary -- face "<<iFDest<<std::endl;
02139       }
02140       break;
02141     }
02142   }
02143   
02144   if (result != SteppingHelixStateInfo::OK && expectNewMagVolume){
02145     double sign = dir == alongMomentum ? 1. : -1.;
02146 
02147     //check if it's a wrong volume situation
02148     if (nDestSorted-nearParallels > 0 ) result = SteppingHelixStateInfo::WRONG_VOLUME;
02149     else {
02150       //get here if all faces in the corr direction were skipped
02151       Point gPointEst(sv.r3);
02152       double lDist = fabs(distToFace[iDistMin]);
02153       if (lDist > fastSkipDist) lDist = fastSkipDist;
02154       Vector lDelta(sv.p3); lDelta *= sign/curP*lDist;
02155       gPointEst += lDelta;
02156       if (debug_){
02157         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Linear est point to shortest dist "<<gPointEst
02158                          <<" for iFace "<<iDistMin<<" at distance "<<lDist*sign<<std::endl;
02159       }
02160       GlobalPoint gPointEstNegZ(gPointEst.x(), gPointEst.y(), -fabs(gPointEst.z()));
02161       GlobalPoint gPointEstNorZ(gPointEst.x(), gPointEst.y(), gPointEst.z() );
02162       if (  (vbField_->isZSymmetric() && cVol->inside(gPointEstNegZ))
02163           || ( !vbField_->isZSymmetric() && cVol->inside(gPointEstNorZ) ) ){
02164         if (debug_){
02165           LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is inside the volume"<<std::endl;
02166         }
02167         
02168       }else {
02169         result = SteppingHelixStateInfo::WRONG_VOLUME;
02170       }
02171     }
02172     
02173     if (result == SteppingHelixStateInfo::WRONG_VOLUME){
02174       dist = sign*0.05;
02175       tanDist = dist*1.01;
02176       if( debug_){
02177         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Wrong volume located: return small dist, tandist"<<std::endl;
02178       }
02179     }
02180   }
02181 
02182   if (result == SteppingHelixStateInfo::INACC){
02183     if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"All faces are too far"<<std::endl;
02184   } else if (result == SteppingHelixStateInfo::WRONG_VOLUME){
02185     if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Appear to be in a wrong volume"<<std::endl;
02186   } else if (result != SteppingHelixStateInfo::OK){
02187     if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Something else went wrong"<<std::endl;
02188   }
02189 
02190   
02191   return result;
02192 }
02193 
02194 
02195 SteppingHelixPropagator::Result
02196 SteppingHelixPropagator::refToMatVolume(const SteppingHelixPropagator::StateInfo& sv,
02197                                         PropagationDirection dir,
02198                                         double& dist, double& tanDist,
02199                                         double fastSkipDist) const{
02200 
02201   static const std::string metname = "SteppingHelixPropagator";
02202   Result result = SteppingHelixStateInfo::NOT_IMPLEMENTED;
02203 
02204   double parLim[6] = {sv.rzLims.rMin, sv.rzLims.rMax, 
02205                       sv.rzLims.zMin, sv.rzLims.zMax, 
02206                       sv.rzLims.th1, sv.rzLims.th2 };
02207 
02208   double distToFace[4] = {0,0,0,0};
02209   double tanDistToFace[4] = {0,0,0,0};
02210   PropagationDirection refDirectionToFace[4] = {anyDirection, anyDirection, anyDirection, anyDirection};
02211   Result resultToFace[4] = {result, result, result, result};
02212   int iFDest = -1;
02213   
02214   double curP = sv.p3.mag();
02215 
02216   for (unsigned int iFace = 0; iFace < 4; iFace++){
02217     if (debug_){
02218       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Start with mat face "<<iFace<<std::endl;
02219     }
02220 
02221     double pars[6] = {0,0,0,0,0,0};
02222     DestType dType = UNDEFINED_DT;
02223     if (iFace > 1){
02224       if (fabs(parLim[iFace+2])< 1e-6){//plane
02225         if (sv.r3.z() < 0){
02226           pars[0] = 0; pars[1] = 0; pars[2] = -parLim[iFace];
02227           pars[3] = 0; pars[4] = 0; pars[5] = 1;
02228         } else {
02229           pars[0] = 0; pars[1] = 0; pars[2] = parLim[iFace];
02230           pars[3] = 0; pars[4] = 0; pars[5] = 1;
02231         }
02232         dType = PLANE_DT;
02233       } else {
02234         if (sv.r3.z() > 0){
02235           pars[0] = 0; pars[1] = 0; 
02236           pars[2] = parLim[iFace];
02237           pars[3] = Geom::pi()/2. - parLim[iFace+2];
02238         } else {
02239           pars[0] = 0; pars[1] = 0; 
02240           pars[2] = - parLim[iFace];
02241           pars[3] = Geom::pi()/2. + parLim[iFace+2];
02242         }
02243         dType = CONE_DT;        
02244       }
02245     } else {
02246       pars[RADIUS_P] = parLim[iFace];
02247       dType = RADIUS_DT;
02248     }
02249 
02250     resultToFace[iFace] = 
02251       refToDest(dType, sv, pars, 
02252                 distToFace[iFace], tanDistToFace[iFace], refDirectionToFace[iFace], fastSkipDist);
02253     
02254     if (resultToFace[iFace] != SteppingHelixStateInfo::OK){
02255       if (resultToFace[iFace] == SteppingHelixStateInfo::INACC) result = SteppingHelixStateInfo::INACC;
02256       continue;
02257     }
02258     if (refDirectionToFace[iFace] == dir || fabs(distToFace[iFace]) < 2e-2*fabs(tanDistToFace[iFace]) ){
02259       double sign = dir == alongMomentum ? 1. : -1.;
02260       Point gPointEst(sv.r3);
02261       Vector lDelta(sv.p3); lDelta *= sign*fabs(distToFace[iFace])/curP;
02262       gPointEst += lDelta;
02263       if (debug_){
02264         LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Linear est point "<<gPointEst
02265                          <<std::endl;
02266       }
02267       double lZ = fabs(gPointEst.z());
02268       double lR = gPointEst.perp();
02269       double tan4 = parLim[4] == 0 ? 0 : tan(parLim[4]);
02270       double tan5 = parLim[5] == 0 ? 0 : tan(parLim[5]);
02271       if ( (lZ - parLim[2]) > lR*tan4 
02272            && (lZ - parLim[3]) < lR*tan5  
02273            && lR > parLim[0] && lR < parLim[1]){
02274         if (debug_){
02275           LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is inside the volume"<<std::endl;
02276         }
02277         //OK, guessed a point still inside the volume
02278         if (iFDest == -1){
02279           iFDest = iFace;
02280         } else {
02281           if (fabs(tanDistToFace[iFDest]) > fabs(tanDistToFace[iFace])){
02282             iFDest = iFace;
02283           }
02284         }
02285       } else {
02286         if (debug_){
02287           LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is NOT inside the volume"<<std::endl;
02288         }
02289       }
02290     }
02291 
02292   }
02293   if (iFDest != -1){
02294     result = SteppingHelixStateInfo::OK;
02295     dist = distToFace[iFDest];
02296     tanDist = tanDistToFace[iFDest];
02297     if (debug_){
02298       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got a point near closest boundary -- face "<<iFDest<<std::endl;
02299     }
02300   } else {
02301     if (debug_){
02302       LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Failed to find a dest point inside the volume"<<std::endl;
02303     }
02304   }
02305 
02306   return result;
02307 }
02308 
02309 
02310 bool SteppingHelixPropagator::isYokeVolume(const MagVolume* vol) const {
02311   static const std::string metname = "SteppingHelixPropagator";
02312   if (vol == 0) return false;
02313   /*
02314   const MFGrid* mGrid = reinterpret_cast<const MFGrid*>(vol->provider());
02315   std::vector<int> dims(mGrid->dimensions());
02316   
02317   LocalVector lVCen(mGrid->nodeValue(dims[0]/2, dims[1]/2, dims[2]/2));
02318   LocalVector lVZLeft(mGrid->nodeValue(dims[0]/2, dims[1]/2, dims[2]/5));
02319   LocalVector lVZRight(mGrid->nodeValue(dims[0]/2, dims[1]/2, (dims[2]*4)/5));
02320 
02321   double mag2VCen = lVCen.mag2();
02322   double mag2VZLeft = lVZLeft.mag2();
02323   double mag2VZRight = lVZRight.mag2();
02324 
02325   bool result = false;
02326   if (mag2VCen > 0.6 && mag2VZLeft > 0.6 && mag2VZRight > 0.6){
02327     if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Volume is magnetic, located at "<<vol->position()<<std::endl;    
02328     result = true;
02329   } else {
02330     if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Volume is not magnetic, located at "<<vol->position()<<std::endl;
02331   }
02332 
02333   */
02334   bool result = vol->isIron();
02335   if (result){
02336     if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Volume is magnetic, located at "<<vol->position()<<std::endl;
02337   } else {
02338     if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Volume is not magnetic, located at "<<vol->position()<<std::endl;
02339   }
02340 
02341   return result;
02342 }