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 if (theVolume != 0) {
00104 LogDebug("RKPropagatorInS") << "RKPropagatorInS: starting prop to plane in volume with pos " << theVolume->position()
00105 << " Z axis " << theVolume->toGlobal( LocalVector(0,0,1)) ;
00106
00107 LogDebug("RKPropagatorInS") << "The starting position is " << ts.position() << " (global) "
00108 << theVolume->toLocal(ts.position()) << " (local) " ;
00109 FrameChanger changer;
00110 FrameChanger::PlanePtr localPlane = changer.transformPlane( plane, *theVolume);
00111 LogDebug("RKPropagatorInS") << "The plane position is " << plane.position() << " (global) "
00112 << localPlane->position() << " (local) " ;
00113
00114 LogDebug("RKPropagatorInS") << "The initial distance to plane is " << plane.localZ( ts.position()) ;
00115
00116 StraightLinePlaneCrossing cross( ts.position().basicVector(), ts.momentum().basicVector());
00117 std::pair<bool,double> res3 = cross.pathLength(plane);
00118 LogDebug("RKPropagatorInS") << "straight line distance " << res3.first << " " << res3.second ;
00119 }
00120
00121 typedef RKAdaptiveSolver<double,RKOneCashKarpStep, 6> Solver;
00122 typedef Solver::Vector RKVector;
00123
00124 RKLocalFieldProvider field( fieldProvider());
00125 PathToPlane2Order pathLength( field, &field.frame());
00126 CartesianLorentzForce deriv(field, ts.charge());
00127
00128 RKCartesianDistance dist;
00129 double eps = theTolerance;
00130 Solver solver;
00131 double stot = 0;
00132 PropagationDirection currentDirection = propagationDirection();
00133
00134
00135 RKVector start( CartesianStateAdaptor::rkstate( rkPosition(gpos), rkMomentum(gmom)));
00136 while (true) {
00137 CartesianStateAdaptor startState(start);
00138
00139 std::pair<bool,double> path = pathLength( plane, startState.position(),
00140 startState.momentum(),
00141 (double) ts.charge(), currentDirection);
00142 if (!path.first) {
00143 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Path length calculation to plane failed!"
00144 << "...distance to plane " << plane.localZ( globalPosition(startState.position()))
00145 << "...Local starting position in volume " << startState.position()
00146 << "...Magnetic field " << field.inTesla( startState.position()) ;
00147
00148
00149 return GlobalParametersWithPath();
00150 }
00151 else {
00152 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Path lenght to plane is " << path.second ;
00153 }
00154
00155 double sstep = path.second;
00156 if ( std::abs(sstep) < eps) {
00157 LogDebug("RKPropagatorInS") << "On-surface accuracy not reached, but pathLength calculation says we are there! "
00158 << "path " << path.second << " distance to plane is " << startZ ;
00159 GlobalTrajectoryParameters res( gtpFromVolumeLocal( startState, ts.charge()));
00160 return GlobalParametersWithPath( res, stot);
00161 }
00162
00163 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Solving for " << sstep
00164 << " current distance to plane is " << startZ ;
00165
00166 RKVector rkresult = solver( 0, start, sstep, deriv, dist, eps);
00167 stot += sstep;
00168 CartesianStateAdaptor cur( rkresult);
00169 double remainingZ = plane.localZ( globalPosition(cur.position()));
00170 if ( fabs(remainingZ) < eps) {
00171 LogDebug("RKPropagatorInS") << "On-surface accuracy reached! " << remainingZ ;
00172 GlobalTrajectoryParameters res( gtpFromVolumeLocal( cur, ts.charge()));
00173 return GlobalParametersWithPath( res, stot);
00174 }
00175
00176 start = rkresult;
00177 if (remainingZ * startZ > 0) {
00178 LogDebug("RKPropagatorInS") << "Accuracy not reached yet, trying in same direction again "
00179 << remainingZ ;
00180 }
00181 else {
00182 LogDebug("RKPropagatorInS") << "Accuracy not reached yet, trying in opposite direction "
00183 << remainingZ ;
00184 currentDirection = invertDirection( currentDirection);
00185 }
00186 startZ = remainingZ;
00187 }
00188 }
00189
00190 GlobalParametersWithPath
00191 RKPropagatorInS::propagateParametersOnCylinder( const FreeTrajectoryState& ts,
00192 const Cylinder& cyl) const
00193 {
00194 typedef RKAdaptiveSolver<double,RKOneCashKarpStep, 6> Solver;
00195 typedef Solver::Vector RKVector;
00196
00197
00198 GlobalPoint sp = cyl.toGlobal(LocalPoint(0.,0.));
00199 if (sp.x()!=0. || sp.y()!=0.) {
00200 throw PropagationException("Cannot propagate to an arbitrary cylinder");
00201 }
00202
00203 GlobalPoint gpos( ts.position());
00204 GlobalVector gmom( ts.momentum());
00205 LocalPoint pos(cyl.toLocal(gpos));
00206 LocalVector mom(cyl.toLocal(gmom));
00207 double startR = cyl.radius() - pos.perp();
00208
00209
00210
00211
00212
00213 double rho = ts.transverseCurvature();
00214
00215
00216
00217
00218 if( fabs(rho)<1.e-10 ) {
00219
00220
00221
00222
00223
00224
00225
00226 StraightLineBarrelCylinderCrossing cylCrossing(gpos,gmom,propagationDirection());
00227
00228
00229
00230
00231 std::pair<bool,double> propResult = cylCrossing.pathLength(cyl);
00232 if ( propResult.first && theVolume !=0) {
00233 double s = propResult.second;
00234
00235 GlobalPoint x (cylCrossing.position(s));
00236 GlobalTrajectoryParameters res( x, gmom, ts.charge(), theVolume);
00237 LogDebug("RKPropagatorInS") << "Straight line propagation to cylinder succeeded !!";
00238 return GlobalParametersWithPath( res, s);
00239 } else {
00240
00241 edm::LogError("RKPropagatorInS") <<"Straight line propagation to cylinder failed !!";
00242 return GlobalParametersWithPath( );
00243 }
00244 }
00245
00246
00247 RKLocalFieldProvider field( fieldProvider(cyl));
00248
00249 CartesianLorentzForce deriv(field, ts.charge());
00250
00251 RKCartesianDistance dist;
00252 double eps = theTolerance;
00253 Solver solver;
00254 double stot = 0;
00255 PropagationDirection currentDirection = propagationDirection();
00256
00257 RKVector start( CartesianStateAdaptor::rkstate( pos.basicVector(), mom.basicVector()));
00258 while (true) {
00259 CartesianStateAdaptor startState(start);
00260 StraightLineCylinderCrossing pathLength( LocalPoint(startState.position()),
00261 LocalVector(startState.momentum()),
00262 currentDirection, eps);
00263
00264 std::pair<bool,double> path = pathLength.pathLength( cyl);
00265 if (!path.first) {
00266 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Path length calculation to cylinder failed!"
00267 << "Radius " << cyl.radius() << " pos.perp() " << LocalPoint(startState.position()).perp() ;
00268 return GlobalParametersWithPath();
00269 }
00270 else {
00271 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Path lenght to cylinder is " << path.second
00272 << " from point (R,z) " << startState.position().perp()
00273 << ", " << startState.position().z()
00274 << " to R " << cyl.radius()
00275 ;
00276 }
00277
00278 double sstep = path.second;
00279 if ( std::abs(sstep) < eps) {
00280 LogDebug("RKPropagatorInS") << "accuracy not reached, but pathLength calculation says we are there! "
00281 << path.second ;
00282
00283 GlobalTrajectoryParameters res( gtpFromLocal( startState.position(),
00284 startState.momentum(),
00285 ts.charge(), cyl));
00286 return GlobalParametersWithPath( res, stot);
00287 }
00288
00289 LogDebug("RKPropagatorInS") << "RKPropagatorInS: Solving for " << sstep
00290 << " current distance to cylinder is " << startR ;
00291
00292 RKVector rkresult = solver( 0, start, sstep, deriv, dist, eps);
00293 stot += sstep;
00294 CartesianStateAdaptor cur( rkresult);
00295 double remainingR = cyl.radius() - cur.position().perp();
00296 if ( fabs(remainingR) < eps) {
00297 LogDebug("RKPropagatorInS") << "Accuracy reached! " << remainingR ;
00298 GlobalTrajectoryParameters res( gtpFromLocal( cur.position(),
00299 cur.momentum(),
00300 ts.charge(), cyl));
00301 return GlobalParametersWithPath( res, stot);
00302 }
00303
00304 start = rkresult;
00305 if (remainingR * startR > 0) {
00306 LogDebug("RKPropagatorInS") << "Accuracy not reached yet, trying in same direction again "
00307 << remainingR ;
00308 }
00309 else {
00310 LogDebug("RKPropagatorInS") << "Accuracy not reached yet, trying in opposite direction "
00311 << remainingR ;
00312 currentDirection = invertDirection( currentDirection);
00313 }
00314 startR = remainingR;
00315 }
00316 }
00317
00318 TrajectoryStateOnSurface
00319 RKPropagatorInS::propagate(const FreeTrajectoryState& fts, const Plane& plane) const
00320 {
00321 return propagateWithPath( fts, plane).first;
00322 }
00323
00324 TrajectoryStateOnSurface
00325 RKPropagatorInS::propagate( const FreeTrajectoryState& fts, const Cylinder& cyl) const
00326 {
00327 return propagateWithPath( fts, cyl).first;
00328 }
00329
00330 Propagator * RKPropagatorInS::clone() const
00331 {
00332 return new RKPropagatorInS(*this);
00333 }
00334
00335 GlobalTrajectoryParameters RKPropagatorInS::gtpFromLocal( const Basic3DVector<double>& lpos,
00336 const Basic3DVector<double>& lmom,
00337 TrackCharge ch, const Surface& surf) const
00338 {
00339 return GlobalTrajectoryParameters( surf.toGlobal( LocalPoint( lpos)),
00340 surf.toGlobal( LocalVector( lmom)), ch, theVolume);
00341 }
00342
00343 RKLocalFieldProvider RKPropagatorInS::fieldProvider() const
00344 {
00345 return RKLocalFieldProvider( *theVolume);
00346 }
00347
00348 RKLocalFieldProvider RKPropagatorInS::fieldProvider( const Cylinder& cyl) const
00349 {
00350 return RKLocalFieldProvider( *theVolume, cyl);
00351 }
00352
00353 PropagationDirection RKPropagatorInS::invertDirection( PropagationDirection dir) const
00354 {
00355 if (dir == anyDirection) return dir;
00356 return ( dir == alongMomentum ? oppositeToMomentum : alongMomentum);
00357 }
00358
00359 Basic3DVector<double> RKPropagatorInS::rkPosition( const GlobalPoint& pos) const
00360 {
00361 if (theVolume != 0) return theVolume->toLocal( pos).basicVector();
00362 else return pos.basicVector();
00363 }
00364
00365 Basic3DVector<double> RKPropagatorInS::rkMomentum( const GlobalVector& mom) const
00366 {
00367 if (theVolume != 0) return theVolume->toLocal( mom).basicVector();
00368 else return mom.basicVector();
00369 }
00370
00371 GlobalPoint RKPropagatorInS::globalPosition( const Basic3DVector<double>& pos) const
00372 {
00373 if (theVolume != 0) return theVolume->toGlobal( LocalPoint(pos));
00374 else return GlobalPoint(pos);
00375 }
00376
00377 GlobalVector RKPropagatorInS::globalMomentum( const Basic3DVector<double>& mom) const
00378
00379 {
00380 if (theVolume != 0) return theVolume->toGlobal( LocalVector(mom));
00381 else return GlobalVector(mom);
00382 }
00383
00384 GlobalTrajectoryParameters
00385 RKPropagatorInS::gtpFromVolumeLocal( const CartesianStateAdaptor& state,
00386 TrackCharge charge) const
00387 {
00388 return GlobalTrajectoryParameters( globalPosition(state.position()),
00389 globalMomentum(state.momentum()),
00390 charge, theVolume);
00391 }