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