00001 #include "TrackPropagation/RungeKutta/interface/RKPropagatorInS.h"
00002 #include "RKCartesianDistance.h"
00003 #include "CartesianLorentzForce.h"
00004 #include "RKLocalFieldProvider.h"
00005 #include "DataFormats/GeometrySurface/interface/Plane.h"
00006 #include "DataFormats/GeometrySurface/interface/Cylinder.h"
00007 #include "RKAdaptiveSolver.h"
00008 #include "RKOne4OrderStep.h"
00009 #include "RKOneCashKarpStep.h"
00010 #include "PathToPlane2Order.h"
00011 #include "CartesianStateAdaptor.h"
00012 #include "TrackingTools/GeomPropagators/interface/StraightLineCylinderCrossing.h"
00013 #include "TrackingTools/GeomPropagators/interface/StraightLineBarrelCylinderCrossing.h"
00014 #include "MagneticField/VolumeGeometry/interface/MagVolume.h"
00015 #include "TrackingTools/GeomPropagators/interface/PropagationDirectionFromPath.h"
00016 #include "FrameChanger.h"
00017 #include "TrackingTools/GeomPropagators/interface/StraightLinePlaneCrossing.h"
00018 #include "AnalyticalErrorPropagation.h"
00019 #include "GlobalParametersWithPath.h"
00020 #include "TrackingTools/GeomPropagators/interface/PropagationExceptions.h"
00021
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023 #include "FWCore/Utilities/interface/Likely.h"
00024
00025
00026 std::pair<TrajectoryStateOnSurface,double>
00027 RKPropagatorInS::propagateWithPath(const FreeTrajectoryState& fts,
00028 const Plane& plane) const
00029 {
00030 GlobalParametersWithPath gp = propagateParametersOnPlane( fts, plane);
00031 if unlikely(!gp) return TsosWP(TrajectoryStateOnSurface(),0.);
00032
00033 SurfaceSideDefinition::SurfaceSide side = PropagationDirectionFromPath()(gp.s(),propagationDirection())==alongMomentum
00034 ? SurfaceSideDefinition::beforeSurface : SurfaceSideDefinition::afterSurface;
00035 AnalyticalErrorPropagation errorprop;
00036 return errorprop( fts, plane, side, gp.parameters(),gp.s());
00037 }
00038
00039 std::pair< TrajectoryStateOnSurface, double>
00040 RKPropagatorInS::propagateWithPath (const FreeTrajectoryState& fts, const Cylinder& cyl) const
00041 {
00042 GlobalParametersWithPath gp = propagateParametersOnCylinder( fts, cyl);
00043 if unlikely(!gp) return TsosWP(TrajectoryStateOnSurface(),0.);
00044
00045 SurfaceSideDefinition::SurfaceSide side = PropagationDirectionFromPath()(gp.s(),propagationDirection())==alongMomentum
00046 ? SurfaceSideDefinition::beforeSurface : SurfaceSideDefinition::afterSurface;
00047 AnalyticalErrorPropagation errorprop;
00048 return errorprop( fts, cyl, side, gp.parameters(),gp.s());
00049
00050 }
00051
00052 GlobalParametersWithPath
00053 RKPropagatorInS::propagateParametersOnPlane( const FreeTrajectoryState& ts,
00054 const Plane& plane) const
00055 {
00056
00057 GlobalPoint gpos( ts.position());
00058 GlobalVector gmom( ts.momentum());
00059 double startZ = plane.localZ(gpos);
00060
00061 double rho = ts.transverseCurvature();
00062
00063
00064
00065
00066 if unlikely( fabs(rho)<1.e-10 ) {
00067
00068
00069
00070
00071
00072
00073 LogDebug("RKPropagatorInS")<<" startZ = "<<startZ;
00074
00075 if unlikely(fabs(startZ) < 1e-5){
00076 LogDebug("RKPropagatorInS")<< "Propagation is not performed: state is already on final surface.";
00077 GlobalTrajectoryParameters res( gpos, gmom, ts.charge(),
00078 theVolume);
00079 return GlobalParametersWithPath( res, 0.0);
00080 }
00081
00082 StraightLinePlaneCrossing::PositionType pos(gpos);
00083 StraightLinePlaneCrossing::DirectionType dir(gmom);
00084 StraightLinePlaneCrossing planeCrossing(pos,dir, propagationDirection());
00085
00086
00087
00088 std::pair<bool,double> propResult = planeCrossing.pathLength(plane);
00089 if likely( propResult.first && theVolume !=0) {
00090 double s = propResult.second;
00091
00092 GlobalPoint x (planeCrossing.position(s));
00093 GlobalTrajectoryParameters res( x, gmom, ts.charge(), theVolume);
00094 return GlobalParametersWithPath( res, s);
00095 }
00096
00097 LogDebug("RKPropagatorInS")<< "Straight line propgation to plane failed !!";
00098 return GlobalParametersWithPath( );
00099
00100 }
00101
00102
00103 #ifdef EDM_LM_DEBUG
00104 if (theVolume != 0) {
00105 LogDebug("RKPropagatorInS") << "RKPropagatorInS: starting prop to plane in volume with pos " << theVolume->position()
00106 << " Z axis " << theVolume->toGlobal( LocalVector(0,0,1)) ;
00107
00108 LogDebug("RKPropagatorInS") << "The starting position is " << ts.position() << " (global) "
00109 << theVolume->toLocal(ts.position()) << " (local) " ;
00110
00111 FrameChanger changer;
00112 FrameChanger::PlanePtr localPlane = changer.transformPlane( plane, *theVolume);
00113 LogDebug("RKPropagatorInS") << "The plane position is " << plane.position() << " (global) "
00114 << localPlane->position() << " (local) " ;
00115
00116 LogDebug("RKPropagatorInS") << "The initial distance to plane is " << plane.localZ( ts.position()) ;
00117
00118 StraightLinePlaneCrossing cross( ts.position().basicVector(), ts.momentum().basicVector());
00119 std::pair<bool,double> res3 = cross.pathLength(plane);
00120 LogDebug("RKPropagatorInS") << "straight line distance " << res3.first << " " << res3.second ;
00121 }
00122 #endif
00123
00124 typedef RKAdaptiveSolver<double,RKOneCashKarpStep, 6> Solver;
00125 typedef Solver::Vector RKVector;
00126
00127 RKLocalFieldProvider field( fieldProvider());
00128 PathToPlane2Order pathLength( field, &field.frame());
00129 CartesianLorentzForce deriv(field, ts.charge());
00130
00131 RKCartesianDistance dist;
00132 double eps = theTolerance;
00133 Solver solver;
00134 double stot = 0;
00135 PropagationDirection currentDirection = propagationDirection();
00136
00137
00138 RKVector start( CartesianStateAdaptor::rkstate( rkPosition(gpos), rkMomentum(gmom)));
00139 int safeGuard = 0;
00140 while (safeGuard++<100) {
00141 CartesianStateAdaptor startState(start);
00142
00143 std::pair<bool,double> path = pathLength( plane, startState.position(),
00144 startState.momentum(),
00145 (double) ts.charge(), currentDirection);
00146 if unlikely(!path.first) {
00147 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Path length calculation to plane failed!"
00148 << "...distance to plane " << plane.localZ( globalPosition(startState.position()))
00149 << "...Local starting position in volume " << startState.position()
00150 << "...Magnetic field " << field.inTesla( startState.position()) ;
00151
00152
00153 return GlobalParametersWithPath();
00154 }
00155
00156 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Path lenght to plane is " << path.second ;
00157
00158
00159 double sstep = path.second;
00160 if unlikely( std::abs(sstep) < eps) {
00161 LogDebug("RKPropagatorInS") << "On-surface accuracy not reached, but pathLength calculation says we are there! "
00162 << "path " << path.second << " distance to plane is " << startZ ;
00163 GlobalTrajectoryParameters res( gtpFromVolumeLocal( startState, ts.charge()));
00164 return GlobalParametersWithPath( res, stot);
00165 }
00166
00167 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Solving for " << sstep
00168 << " current distance to plane is " << startZ ;
00169
00170 RKVector rkresult = solver( 0, start, sstep, deriv, dist, eps);
00171 stot += sstep;
00172 CartesianStateAdaptor cur( rkresult);
00173 double remainingZ = plane.localZ( globalPosition(cur.position()));
00174
00175 if ( fabs(remainingZ) < eps) {
00176 LogDebug("RKPropagatorInS") << "On-surface accuracy reached! " << remainingZ ;
00177 GlobalTrajectoryParameters res( gtpFromVolumeLocal( cur, ts.charge()));
00178 return GlobalParametersWithPath( res, stot);
00179 }
00180
00181 start = rkresult;
00182
00183 if (remainingZ * startZ > 0) {
00184 LogDebug("RKPropagatorInS") << "Accuracy not reached yet, trying in same direction again "
00185 << remainingZ ;
00186 }
00187 else {
00188 LogDebug("RKPropagatorInS") << "Accuracy not reached yet, trying in opposite direction "
00189 << remainingZ ;
00190 currentDirection = invertDirection( currentDirection);
00191 }
00192 startZ = remainingZ;
00193 }
00194
00195 edm::LogError("FailedPropagation") << " too many iterations trying to reach plane ";
00196 return GlobalParametersWithPath();
00197 }
00198
00199 GlobalParametersWithPath
00200 RKPropagatorInS::propagateParametersOnCylinder( const FreeTrajectoryState& ts,
00201 const Cylinder& cyl) const {
00202
00203 typedef RKAdaptiveSolver<double,RKOneCashKarpStep, 6> Solver;
00204 typedef Solver::Vector RKVector;
00205
00206
00207 GlobalPoint sp = cyl.position();
00208 if unlikely(sp.x()!=0. || sp.y()!=0.) {
00209 throw PropagationException("Cannot propagate to an arbitrary cylinder");
00210 }
00211
00212 GlobalPoint gpos( ts.position());
00213 GlobalVector gmom( ts.momentum());
00214 LocalPoint pos(cyl.toLocal(gpos));
00215 LocalVector mom(cyl.toLocal(gmom));
00216 double startR = cyl.radius() - pos.perp();
00217
00218
00219
00220
00221
00222 double rho = ts.transverseCurvature();
00223
00224
00225
00226
00227 if unlikely( fabs(rho)<1.e-10 ) {
00228
00229
00230
00231
00232
00233
00234
00235 StraightLineBarrelCylinderCrossing cylCrossing(gpos,gmom,propagationDirection());
00236
00237
00238
00239
00240 std::pair<bool,double> propResult = cylCrossing.pathLength(cyl);
00241 if likely( propResult.first && theVolume !=0) {
00242 double s = propResult.second;
00243
00244 GlobalPoint x (cylCrossing.position(s));
00245 GlobalTrajectoryParameters res( x, gmom, ts.charge(), theVolume);
00246 LogDebug("RKPropagatorInS") << "Straight line propagation to cylinder succeeded !!";
00247 return GlobalParametersWithPath( res, s);
00248 }
00249
00250
00251 edm::LogError("RKPropagatorInS") <<"Straight line propagation to cylinder failed !!";
00252 return GlobalParametersWithPath( );
00253
00254 }
00255
00256
00257 RKLocalFieldProvider field( fieldProvider(cyl));
00258
00259 CartesianLorentzForce deriv(field, ts.charge());
00260
00261 RKCartesianDistance dist;
00262 double eps = theTolerance;
00263 Solver solver;
00264 double stot = 0;
00265 PropagationDirection currentDirection = propagationDirection();
00266
00267 RKVector start( CartesianStateAdaptor::rkstate( pos.basicVector(), mom.basicVector()));
00268 int safeGuard = 0;
00269 while (safeGuard++<100) {
00270 CartesianStateAdaptor startState(start);
00271 StraightLineCylinderCrossing pathLength( LocalPoint(startState.position()),
00272 LocalVector(startState.momentum()),
00273 currentDirection, eps);
00274
00275 std::pair<bool,double> path = pathLength.pathLength( cyl);
00276 if unlikely(!path.first) {
00277 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Path length calculation to cylinder failed!"
00278 << "Radius " << cyl.radius() << " pos.perp() " << LocalPoint(startState.position()).perp() ;
00279 return GlobalParametersWithPath();
00280 }
00281
00282 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Path lenght to cylinder is " << path.second
00283 << " from point (R,z) " << startState.position().perp()
00284 << ", " << startState.position().z()
00285 << " to R " << cyl.radius()
00286 ;
00287
00288
00289 double sstep = path.second;
00290 if unlikely( std::abs(sstep) < eps) {
00291 LogDebug("RKPropagatorInS") << "accuracy not reached, but pathLength calculation says we are there! "
00292 << path.second ;
00293
00294 GlobalTrajectoryParameters res( gtpFromLocal( startState.position(),
00295 startState.momentum(),
00296 ts.charge(), cyl));
00297 return GlobalParametersWithPath( res, stot);
00298 }
00299
00300 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Solving for " << sstep
00301 << " current distance to cylinder is " << startR ;
00302
00303 RKVector rkresult = solver( 0, start, sstep, deriv, dist, eps);
00304 stot += sstep;
00305 CartesianStateAdaptor cur( rkresult);
00306 double remainingR = cyl.radius() - cur.position().perp();
00307
00308 if ( fabs(remainingR) < eps) {
00309 LogDebug("RKPropagatorInS") << "Accuracy reached! " << remainingR ;
00310 GlobalTrajectoryParameters res( gtpFromLocal( cur.position(),
00311 cur.momentum(),
00312 ts.charge(), cyl));
00313 return GlobalParametersWithPath( res, stot);
00314 }
00315
00316 start = rkresult;
00317 if (remainingR * startR > 0) {
00318 LogDebug("RKPropagatorInS") << "Accuracy not reached yet, trying in same direction again "
00319 << remainingR ;
00320 }
00321 else {
00322 LogDebug("RKPropagatorInS") << "Accuracy not reached yet, trying in opposite direction "
00323 << remainingR ;
00324 currentDirection = invertDirection( currentDirection);
00325 }
00326 startR = remainingR;
00327 }
00328
00329 edm::LogError("FailedPropagation") << " too many iterations trying to reach cylinder ";
00330 return GlobalParametersWithPath();
00331 }
00332
00333 TrajectoryStateOnSurface
00334 RKPropagatorInS::propagate(const FreeTrajectoryState& fts, const Plane& plane) const
00335 {
00336 return propagateWithPath( fts, plane).first;
00337 }
00338
00339 TrajectoryStateOnSurface
00340 RKPropagatorInS::propagate( const FreeTrajectoryState& fts, const Cylinder& cyl) const
00341 {
00342 return propagateWithPath( fts, cyl).first;
00343 }
00344
00345 Propagator * RKPropagatorInS::clone() const
00346 {
00347 return new RKPropagatorInS(*this);
00348 }
00349
00350 GlobalTrajectoryParameters RKPropagatorInS::gtpFromLocal( const Basic3DVector<double>& lpos,
00351 const Basic3DVector<double>& lmom,
00352 TrackCharge ch, const Surface& surf) const
00353 {
00354 return GlobalTrajectoryParameters( surf.toGlobal( LocalPoint( lpos)),
00355 surf.toGlobal( LocalVector( lmom)), ch, theVolume);
00356 }
00357
00358 RKLocalFieldProvider RKPropagatorInS::fieldProvider() const
00359 {
00360 return RKLocalFieldProvider( *theVolume);
00361 }
00362
00363 RKLocalFieldProvider RKPropagatorInS::fieldProvider( const Cylinder& cyl) const
00364 {
00365 return RKLocalFieldProvider( *theVolume, cyl);
00366 }
00367
00368 PropagationDirection RKPropagatorInS::invertDirection( PropagationDirection dir) const
00369 {
00370 if (dir == anyDirection) return dir;
00371 return ( dir == alongMomentum ? oppositeToMomentum : alongMomentum);
00372 }
00373
00374 Basic3DVector<double> RKPropagatorInS::rkPosition( const GlobalPoint& pos) const
00375 {
00376 if (theVolume != 0) return theVolume->toLocal( pos).basicVector();
00377 else return pos.basicVector();
00378 }
00379
00380 Basic3DVector<double> RKPropagatorInS::rkMomentum( const GlobalVector& mom) const
00381 {
00382 if (theVolume != 0) return theVolume->toLocal( mom).basicVector();
00383 else return mom.basicVector();
00384 }
00385
00386 GlobalPoint RKPropagatorInS::globalPosition( const Basic3DVector<double>& pos) const
00387 {
00388 if (theVolume != 0) return theVolume->toGlobal( LocalPoint(pos));
00389 else return GlobalPoint(pos);
00390 }
00391
00392 GlobalVector RKPropagatorInS::globalMomentum( const Basic3DVector<double>& mom) const
00393
00394 {
00395 if (theVolume != 0) return theVolume->toGlobal( LocalVector(mom));
00396 else return GlobalVector(mom);
00397 }
00398
00399 GlobalTrajectoryParameters
00400 RKPropagatorInS::gtpFromVolumeLocal( const CartesianStateAdaptor& state,
00401 TrackCharge charge) const
00402 {
00403 return GlobalTrajectoryParameters( globalPosition(state.position()),
00404 globalMomentum(state.momentum()),
00405 charge, theVolume);
00406 }