00001
00013
00014
00015
00016
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;
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
00230
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
00258
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;
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
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;
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
00387
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
00396
00397 tanDistNextCheck -= oldDStep;
00398 tanDistMagNextCheck -= oldDStep;
00399 tanDistMatNextCheck -= oldDStep;
00400
00401 if (tanDistNextCheck < 0){
00402
00403 if (! isFirstStep) refToDest(type, (*svCurrent), pars, dist, tanDist, refDirection);
00404
00405 if (fabs(tanDist) > 4.*(fabs(dist)) ) tanDist *= tanDist == 0 ? 0 :fabs(dist/tanDist*4.);
00406
00407 tanDistNextCheck = fabs(tanDist)*0.5 - 0.5;
00408
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
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))){
00432 if (tanDistMagNextCheck < 0){
00433 resultToMag = refToMagVolume((*svCurrent), dir, distMag, tanDistMag, fabs(fastSkipDist), expectNewMagVolume, fabs(tanDist));
00434
00435 if (fabs(tanDistMag) > 4.*(fabs(distMag)) ) tanDistMag *= tanDistMag == 0 ? 0 : fabs(distMag/tanDistMag*4.);
00436
00437 tanDistMagNextCheck = fabs(tanDistMag)*0.5-0.5;
00438
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
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))){
00454 if (tanDistMatNextCheck < 0){
00455 resultToMat = refToMatVolume((*svCurrent), dir, distMat, tanDistMat, fabs(fastSkipDist));
00456
00457 if (fabs(tanDistMat) > 4.*(fabs(distMat)) ) tanDistMat *= tanDistMat == 0 ? 0 : fabs(distMat/tanDistMat*4.);
00458
00459 tanDistMatNextCheck = fabs(tanDistMat)*0.5-0.5;
00460
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
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;
00486 expectNewMagVolume = true;
00487 } else expectNewMagVolume = false;
00488
00489 if (tanDistMin > fabs(tanDistMat)+0.05 && resultToMat == SteppingHelixStateInfo::OK){
00490 tanDistMin = fabs(tanDistMat)+0.05;
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
00505 tanDistMinLazy = fabs(tanDistMin)*0.5;
00506 }
00507
00508 if (fabs(tanDistMinLazy) < dStep){
00509 dStep = fabs(tanDistMinLazy);
00510 }
00511
00512
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
00519
00520
00521
00522
00523
00524
00525
00526
00527
00528
00529
00530 nPoints_++; svCurrent = &svBuf_[cIndex_(nPoints_-1)];
00531 if (oldDir != dir){
00532 nOsc++;
00533 tanDistNextCheck = -1;
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
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
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;
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 float pmag2 = p3.mag2();
00671 if (gpmag > 1e20f ) {
00672 LogTrace(metname)<<"Initial point is too far";
00673 svCurrent.isValid_ = false;
00674 return;
00675 }
00676 if (! (gpmag == gpmag) ) {
00677 LogTrace(metname)<<"Initial point is a nan";
00678 edm::LogWarning("SteppingHelixPropagatorNAN")<<"Initial point is a nan";
00679 svCurrent.isValid_ = false;
00680 return;
00681 }
00682 if (! (pmag2 == pmag2) ) {
00683 LogTrace(metname)<<"Initial momentum is a nan";
00684 edm::LogWarning("SteppingHelixPropagatorNAN")<<"Initial momentum is a nan";
00685 svCurrent.isValid_ = false;
00686 return;
00687 }
00688
00689 GlobalVector bf(0,0,0);
00690
00691 if (useMagVolumes_){
00692 if (vbField_ ){
00693 if (vbField_->isZSymmetric()){
00694 svCurrent.magVol = vbField_->findVolume(gPointNegZ);
00695 } else {
00696 svCurrent.magVol = vbField_->findVolume(gPointNorZ);
00697 }
00698 if (useIsYokeFlag_){
00699 double curRad = svCurrent.r3.perp();
00700 if (curRad > 380 && curRad < 850 && fabs(svCurrent.r3.z()) < 667){
00701 svCurrent.isYokeVol = isYokeVolume(svCurrent.magVol);
00702 } else {
00703 svCurrent.isYokeVol = false;
00704 }
00705 }
00706 } else {
00707 edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Failed to cast into VolumeBasedMagneticField: fall back to the default behavior"<<std::endl;
00708 svCurrent.magVol = 0;
00709 }
00710 if (debug_){
00711 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got volume at "<<svCurrent.magVol<<std::endl;
00712 }
00713 }
00714
00715 if (useMagVolumes_ && svCurrent.magVol != 0 && ! useInTeslaFromMagField_){
00716 if (vbField_->isZSymmetric()){
00717 bf = svCurrent.magVol->inTesla(gPointNegZ);
00718 if (debug_){
00719 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Loaded bfield float: "<<bf
00720 <<" at global float "<< gPointNegZ<<" double "<< svCurrent.r3<<std::endl;
00721 LocalPoint lPoint(svCurrent.magVol->toLocal(gPointNegZ));
00722 LocalVector lbf = svCurrent.magVol->fieldInTesla(lPoint);
00723 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"\t cf in local locF: "<<lbf<<" at "<<lPoint<<std::endl;
00724 }
00725 } else {
00726 bf = svCurrent.magVol->inTesla(gPointNorZ);
00727 if (debug_){
00728 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Loaded bfield float: "<<bf
00729 <<" at global float "<< gPointNorZ<<" double "<< svCurrent.r3<<std::endl;
00730 LocalPoint lPoint(svCurrent.magVol->toLocal(gPointNorZ));
00731 LocalVector lbf = svCurrent.magVol->fieldInTesla(lPoint);
00732 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"\t cf in local locF: "<<lbf<<" at "<<lPoint<<std::endl;
00733 }
00734 }
00735 if (r3.z() > 0 && vbField_->isZSymmetric() ){
00736 svCurrent.bf.set(-bf.x(), -bf.y(), bf.z());
00737 } else {
00738 svCurrent.bf.set(bf.x(), bf.y(), bf.z());
00739 }
00740 } else {
00741 GlobalPoint gPoint(r3.x(), r3.y(), r3.z());
00742 bf = field_->inTesla(gPoint);
00743 svCurrent.bf.set(bf.x(), bf.y(), bf.z());
00744 }
00745 if (svCurrent.bf.mag2() < 1e-32) svCurrent.bf.set(0., 0., 1e-16);
00746 if (debug_){
00747 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific
00748 <<"Loaded bfield double: "<<svCurrent.bf<<" from float: "<<bf
00749 <<" at float "<< gPointNorZ<<" "<<gPointNegZ<<" double "<< svCurrent.r3<<std::endl;
00750 }
00751
00752
00753
00754 double dEdXPrime = 0;
00755 double dEdx = 0;
00756 double radX0 = 0;
00757 MatBounds rzLims;
00758 dEdx = getDeDx(svCurrent, dEdXPrime, radX0, rzLims);
00759 svCurrent.dEdx = dEdx; svCurrent.dEdXPrime = dEdXPrime;
00760 svCurrent.radX0 = radX0;
00761 svCurrent.rzLims = rzLims;
00762
00763 svCurrent.covCurv =covCurv;
00764
00765 svCurrent.isComplete = true;
00766 svCurrent.isValid_ = true;
00767
00768 if (debug_){
00769 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Loaded at path: "<<svCurrent.path()<<" radPath: "<<svCurrent.radPath()
00770 <<" p3 "<<" pt: "<<svCurrent.p3.perp()<<" phi: "<<svCurrent.p3.phi()
00771 <<" eta: "<<svCurrent.p3.eta()
00772 <<" "<<svCurrent.p3
00773 <<" r3: "<<svCurrent.r3
00774 <<" bField: "<<svCurrent.bf.mag()
00775 <<std::endl;
00776 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Input Covariance in Curvilinear RF "<<covCurv<<std::endl;
00777 }
00778 }
00779
00780 void SteppingHelixPropagator::getNextState(const SteppingHelixPropagator::StateInfo& svPrevious,
00781 SteppingHelixPropagator::StateInfo& svNext,
00782 double dP, const SteppingHelixPropagator::Vector& tau,
00783 const SteppingHelixPropagator::Vector& drVec, double dS, double dX0,
00784 const AlgebraicMatrix55& dCovCurvTransform) const{
00785 static const std::string metname = "SteppingHelixPropagator";
00786 svNext.q = svPrevious.q;
00787 svNext.dir = dS > 0.0 ? 1.: -1.;
00788 svNext.p3 = tau; svNext.p3*=(svPrevious.p3.mag() - svNext.dir*fabs(dP));
00789
00790 svNext.r3 = svPrevious.r3; svNext.r3 += drVec;
00791
00792 svNext.path_ = svPrevious.path() + dS;
00793 svNext.radPath_ = svPrevious.radPath() + dX0;
00794
00795 GlobalPoint gPointNegZ(svNext.r3.x(), svNext.r3.y(), -fabs(svNext.r3.z()));
00796 GlobalPoint gPointNorZ(svNext.r3.x(), svNext.r3.y(), svNext.r3.z());
00797
00798 GlobalVector bf(0,0,0);
00799
00800 if (useMagVolumes_){
00801 if (vbField_ != 0){
00802 if (vbField_->isZSymmetric()){
00803 svNext.magVol = vbField_->findVolume(gPointNegZ);
00804 } else {
00805 svNext.magVol = vbField_->findVolume(gPointNorZ);
00806 }
00807 if (useIsYokeFlag_){
00808 double curRad = svNext.r3.perp();
00809 if (curRad > 380 && curRad < 850 && fabs(svNext.r3.z()) < 667){
00810 svNext.isYokeVol = isYokeVolume(svNext.magVol);
00811 } else {
00812 svNext.isYokeVol = false;
00813 }
00814 }
00815 } else {
00816 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Failed to cast into VolumeBasedMagneticField"<<std::endl;
00817 svNext.magVol = 0;
00818 }
00819 if (debug_){
00820 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got volume at "<<svNext.magVol<<std::endl;
00821 }
00822 }
00823
00824 if (useMagVolumes_ && svNext.magVol != 0 && ! useInTeslaFromMagField_){
00825 if (vbField_->isZSymmetric()){
00826 bf = svNext.magVol->inTesla(gPointNegZ);
00827 } else {
00828 bf = svNext.magVol->inTesla(gPointNorZ);
00829 }
00830 if (svNext.r3.z() > 0 && vbField_->isZSymmetric() ){
00831 svNext.bf.set(-bf.x(), -bf.y(), bf.z());
00832 } else {
00833 svNext.bf.set(bf.x(), bf.y(), bf.z());
00834 }
00835 } else {
00836 GlobalPoint gPoint(svNext.r3.x(), svNext.r3.y(), svNext.r3.z());
00837 bf = field_->inTesla(gPoint);
00838 svNext.bf.set(bf.x(), bf.y(), bf.z());
00839 }
00840 if (svNext.bf.mag2() < 1e-32) svNext.bf.set(0., 0., 1e-16);
00841
00842
00843 double dEdXPrime = 0;
00844 double dEdx = 0;
00845 double radX0 = 0;
00846 MatBounds rzLims;
00847 dEdx = getDeDx(svNext, dEdXPrime, radX0, rzLims);
00848 svNext.dEdx = dEdx; svNext.dEdXPrime = dEdXPrime;
00849 svNext.radX0 = radX0;
00850 svNext.rzLims = rzLims;
00851
00852
00853 svNext.hasErrorPropagated_ = svPrevious.hasErrorPropagated_;
00854 if (svPrevious.hasErrorPropagated_){
00855 {
00856 AlgebraicMatrix55 tmp = dCovCurvTransform*svPrevious.covCurv;
00857 ROOT::Math::AssignSym::Evaluate(svNext.covCurv, tmp*ROOT::Math::Transpose(dCovCurvTransform));
00858
00859 svNext.covCurv += svPrevious.matDCovCurv;
00860 }
00861 } else {
00862
00863
00864 }
00865
00866 if (debug_){
00867 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Now at path: "<<svNext.path()<<" radPath: "<<svNext.radPath()
00868 <<" p3 "<<" pt: "<<svNext.p3.perp()<<" phi: "<<svNext.p3.phi()
00869 <<" eta: "<<svNext.p3.eta()
00870 <<" "<<svNext.p3
00871 <<" r3: "<<svNext.r3
00872 <<" dPhi: "<<acos(svNext.p3.unit().dot(svPrevious.p3.unit()))
00873 <<" bField: "<<svNext.bf.mag()
00874 <<" dedx: "<<svNext.dEdx
00875 <<std::endl;
00876 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"New Covariance "<<svNext.covCurv<<std::endl;
00877 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Transf by dCovTransform "<<dCovCurvTransform<<std::endl;
00878 }
00879 }
00880
00881 void SteppingHelixPropagator::setRep(SteppingHelixPropagator::Basis& rep,
00882 const SteppingHelixPropagator::Vector& tau) const{
00883 Vector zRep(0., 0., 1.);
00884 rep.lX = tau;
00885 rep.lY = zRep.cross(tau); rep.lY *= 1./tau.perp();
00886 rep.lZ = rep.lX.cross(rep.lY);
00887 }
00888
00889 bool SteppingHelixPropagator::makeAtomStep(SteppingHelixPropagator::StateInfo& svCurrent,
00890 SteppingHelixPropagator::StateInfo& svNext,
00891 double dS,
00892 PropagationDirection dir,
00893 SteppingHelixPropagator::Fancy fancy) const{
00894 static const std::string metname = "SteppingHelixPropagator";
00895 if (debug_){
00896 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Make atom step "<<svCurrent.path()<<" with step "<<dS<<" in direction "<<dir<<std::endl;
00897 }
00898
00899 double dP = 0;
00900 double curP = svCurrent.p3.mag();
00901 Vector tau = svCurrent.p3; tau *= 1./curP;
00902 Vector tauNext(tau);
00903 Vector drVec(0,0,0);
00904
00905 dS = dir == alongMomentum ? fabs(dS) : -fabs(dS);
00906
00907
00908 double radX0 = 1e24;
00909
00910 switch (fancy){
00911 case HEL_AS_F:
00912 case HEL_ALL_F:{
00913 double p0 = curP;
00914 double p0Inv = 1./p0;
00915 double b0 = svCurrent.bf.mag();
00916
00917
00918 double phi = (2.99792458e-3*svCurrent.q*b0*p0Inv)*dS/2.;
00919 bool phiSmall = fabs(phi) < 1e-4;
00920
00921 double cosPhi = 0;
00922 double sinPhi = 0;
00923
00924 double oneLessCosPhi=0;
00925 double oneLessCosPhiOPhi=0;
00926 double sinPhiOPhi=0;
00927 double phiLessSinPhiOPhi=0;
00928
00929 if (phiSmall){
00930 double phi2 = phi*phi;
00931 double phi3 = phi2*phi;
00932 double phi4 = phi3*phi;
00933 sinPhi = phi - phi3/6. + phi4*phi/120.;
00934 cosPhi = 1. -phi2/2. + phi4/24.;
00935 oneLessCosPhi = phi2/2. - phi4/24. + phi2*phi4/720.;
00936 oneLessCosPhiOPhi = 0.5*phi - phi3/24. + phi2*phi3/720.;
00937 sinPhiOPhi = 1. - phi*phi/6. + phi4/120.;
00938 phiLessSinPhiOPhi = phi*phi/6. - phi4/120. + phi4*phi2/5040.;
00939 } else {
00940 cosPhi = cos(phi);
00941 sinPhi = sin(phi);
00942 oneLessCosPhi = 1.-cosPhi;
00943 oneLessCosPhiOPhi = oneLessCosPhi/phi;
00944 sinPhiOPhi = sinPhi/phi;
00945 phiLessSinPhiOPhi = 1 - sinPhiOPhi;
00946 }
00947
00948 Vector bHat = svCurrent.bf; bHat *= 1./b0;
00949
00950
00951 Vector btVec(bHat.cross(tau));
00952 double tauB = tau.dot(bHat);
00953 Vector bbtVec(bHat*tauB - tau);
00954
00955
00956 drVec = bbtVec*phiLessSinPhiOPhi; drVec -= btVec*oneLessCosPhiOPhi; drVec += tau;
00957 drVec *= dS/2.;
00958
00959 double dEdx = svCurrent.dEdx;
00960 double dEdXPrime = svCurrent.dEdXPrime;
00961 radX0 = svCurrent.radX0;
00962 dP = dEdx*dS;
00963
00964
00965 drVec += svCurrent.r3;
00966 GlobalVector bfGV(0,0,0);
00967 Vector bf(0,0,0);
00968
00969 if (useMagVolumes_ && svCurrent.magVol != 0 && ! useInTeslaFromMagField_){
00970
00971 if (vbField_->isZSymmetric()){
00972 bfGV = svCurrent.magVol->inTesla(GlobalPoint(drVec.x(), drVec.y(), -fabs(drVec.z())));
00973 } else {
00974 bfGV = svCurrent.magVol->inTesla(GlobalPoint(drVec.x(), drVec.y(), drVec.z()));
00975 }
00976 if (drVec.z() > 0 && vbField_->isZSymmetric()){
00977 bf.set(-bfGV.x(), -bfGV.y(), bfGV.z());
00978 } else {
00979 bf.set(bfGV.x(), bfGV.y(), bfGV.z());
00980 }
00981 } else {
00982 bfGV = field_->inTesla(GlobalPoint(drVec.x(), drVec.y(), drVec.z()));
00983 bf.set(bfGV.x(), bfGV.y(), bfGV.z());
00984 }
00985 double b0Init = b0;
00986 b0 = bf.mag();
00987 if (b0 < 1e-16) {
00988 b0 = 1e-16;
00989 bf.set(0., 0., 1e-16);
00990 }
00991 if (debug_){
00992 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Improved b "<<b0
00993 <<" at r3 "<<drVec<<std::endl;
00994 }
00995
00996 if (fabs((b0-b0Init)*dS) > 1){
00997
00998 if (debug_){
00999 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Large bf*dS change "<<fabs((b0-svCurrent.bf.mag())*dS)
01000 <<" --> recalc dedx"<<std::endl;
01001 }
01002 svNext.r3 = drVec;
01003 svNext.bf = bf;
01004 svNext.p3 = svCurrent.p3;
01005 svNext.isYokeVol = svCurrent.isYokeVol;
01006 MatBounds rzTmp;
01007 dEdx = getDeDx(svNext, dEdXPrime, radX0, rzTmp);
01008 dP = dEdx*dS;
01009 if (debug_){
01010 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"New dEdX= "<<dEdx
01011 <<" dP= "<<dP
01012 <<" for p0 "<<p0<<std::endl;
01013 }
01014 }
01015
01016 p0 += dP/2.; if (p0 < 1e-2) p0 = 1e-2;
01017 p0Inv = 1./p0;
01018
01019 phi = (2.99792458e-3*svCurrent.q*b0*p0Inv)*dS;
01020 phiSmall = fabs(phi) < 1e-4;
01021
01022 if (phiSmall){
01023 double phi2 = phi*phi;
01024 double phi3 = phi2*phi;
01025 double phi4 = phi3*phi;
01026 sinPhi = phi - phi3/6. + phi4*phi/120.;
01027 cosPhi = 1. -phi2/2. + phi4/24.;
01028 oneLessCosPhi = phi2/2. - phi4/24. + phi2*phi4/720.;
01029 oneLessCosPhiOPhi = 0.5*phi - phi3/24. + phi2*phi3/720.;
01030 sinPhiOPhi = 1. - phi*phi/6. + phi4/120.;
01031 phiLessSinPhiOPhi = phi*phi/6. - phi4/120. + phi4*phi2/5040.;
01032 }else {
01033 cosPhi = cos(phi);
01034 sinPhi = sin(phi);
01035 oneLessCosPhi = 1.-cosPhi;
01036 oneLessCosPhiOPhi = oneLessCosPhi/phi;
01037 sinPhiOPhi = sinPhi/phi;
01038 phiLessSinPhiOPhi = 1. - sinPhiOPhi;
01039 }
01040
01041 bHat = bf; bHat *= 1./b0;
01042
01043 btVec = bHat.cross(tau);
01044 tauB = tau.dot(bHat);
01045 bbtVec = bHat*tauB - tau;
01046
01047 tauNext = bbtVec*oneLessCosPhi; tauNext -= btVec*sinPhi; tauNext += tau;
01048 double tauNextPerpInv = 1./tauNext.perp();
01049 drVec = bbtVec*phiLessSinPhiOPhi; drVec -= btVec*oneLessCosPhiOPhi; drVec += tau;
01050 drVec *= dS;
01051
01052
01053 if (svCurrent.hasErrorPropagated_){
01054 double theta02 = 0;
01055 double dX0 = fabs(dS)/radX0;
01056
01057 if (applyRadX0Correction_){
01058
01059
01060
01061
01062 double x0 = fabs(svCurrent.radPath());
01063 double alphaX0 = 13.6e-3*p0Inv; alphaX0 *= alphaX0;
01064 double betaX0 = 0.038;
01065 double logx0p1 = log(x0+1);
01066 theta02 = dX0*alphaX0*(1+betaX0*logx0p1)*(1 + betaX0*logx0p1 + 2.*betaX0*x0/(x0+1) );
01067 } else {
01068 theta02 = 196e-6* p0Inv * p0Inv * dX0;
01069 }
01070
01071 double epsilonP0 = 1.+ dP/(p0-0.5*dP);
01072
01073
01074
01075 Vector tbtVec(bHat - tauB*tau);
01076
01077 dCCurvTransform_ = unit55_;
01078 {
01079
01080
01081
01082
01083
01084
01085 Point xStart = svCurrent.r3;
01086 Vector dx = drVec;
01087
01088
01089
01090
01091 double qbp = svCurrent.q*p0Inv;
01092
01093
01094
01095
01096 double t11 = tau.x(); double t12 = tau.y(); double t13 = tau.z();
01097 double t21 = tauNext.x(); double t22 = tauNext.y(); double t23 = tauNext.z();
01098 double cosl0 = tau.perp();
01099
01100 double cosecl1 = tauNextPerpInv;
01101
01102
01103
01104 Vector hn = bHat;
01105
01106
01107
01108 double q = -phi/dS;
01109
01110 double sint = -sinPhi; double cost = cosPhi;
01111 double hn1 = hn.x(); double hn2 = hn.y(); double hn3 = hn.z();
01112 double dx1 = dx.x(); double dx2 = dx.y(); double dx3 = dx.z();
01113
01114
01115 double gamma = hn1*t21 + hn2*t22 + hn3*t23;
01116 double an1 = hn2*t23 - hn3*t22;
01117 double an2 = hn3*t21 - hn1*t23;
01118 double an3 = hn1*t22 - hn2*t21;
01119
01120 double auInv = cosl0; double au = auInv>0 ? 1./auInv : 1e24;
01121
01122 double u11 = -au*t12; double u12 = au*t11;
01123 double v11 = -t13*u12; double v12 = t13*u11; double v13 = auInv;
01124 auInv = sqrt(1. - t23*t23); au = auInv>0 ? 1./auInv : 1e24;
01125
01126 double u21 = -au*t22; double u22 = au*t21;
01127 double v21 = -t23*u22; double v22 = t23*u21; double v23 = auInv;
01128
01129
01130
01132
01133
01134 double anv = -(hn1*u21 + hn2*u22 );
01135 double anu = (hn1*v21 + hn2*v22 + hn3*v23);
01136 double omcost = oneLessCosPhi; double tmsint = -phi*phiLessSinPhiOPhi;
01137
01138 double hu1 = - hn3*u12;
01139 double hu2 = hn3*u11;
01140 double hu3 = hn1*u12 - hn2*u11;
01141
01142 double hv1 = hn2*v13 - hn3*v12;
01143 double hv2 = hn3*v11 - hn1*v13;
01144 double hv3 = hn1*v12 - hn2*v11;
01145
01146
01147 dCCurvTransform_(0,0) = 1./(epsilonP0*epsilonP0)*(1. + dS*dEdXPrime);
01148
01149
01150
01151 dCCurvTransform_(1,0) = phi*p0/svCurrent.q*cosecl1*
01152 (sinPhi*bbtVec.z() - cosPhi*btVec.z());
01153
01154
01155
01156 dCCurvTransform_(1,1) = cost*(v11*v21 + v12*v22 + v13*v23) +
01157 sint*(hv1*v21 + hv2*v22 + hv3*v23) +
01158 omcost*(hn1*v11 + hn2*v12 + hn3*v13) * (hn1*v21 + hn2*v22 + hn3*v23) +
01159 anv*(-sint*(v11*t21 + v12*t22 + v13*t23) +
01160 omcost*(v11*an1 + v12*an2 + v13*an3) -
01161 tmsint*gamma*(hn1*v11 + hn2*v12 + hn3*v13) );
01162
01163 dCCurvTransform_(1,2) = cost*(u11*v21 + u12*v22 ) +
01164 sint*(hu1*v21 + hu2*v22 + hu3*v23) +
01165 omcost*(hn1*u11 + hn2*u12 ) * (hn1*v21 + hn2*v22 + hn3*v23) +
01166 anv*(-sint*(u11*t21 + u12*t22 ) +
01167 omcost*(u11*an1 + u12*an2 ) -
01168 tmsint*gamma*(hn1*u11 + hn2*u12 ) );
01169 dCCurvTransform_(1,2) *= cosl0;
01170
01171
01172
01173
01174
01175
01176
01177
01178
01179 dCCurvTransform_(2,0) = - phi*p0/svCurrent.q*cosecl1*cosecl1*
01180 (oneLessCosPhi*bHat.z()*btVec.mag2() + sinPhi*btVec.z() + cosPhi*tbtVec.z()) ;
01181
01182
01183 dCCurvTransform_(2,1) = cost*(v11*u21 + v12*u22 ) +
01184 sint*(hv1*u21 + hv2*u22 ) +
01185 omcost*(hn1*v11 + hn2*v12 + hn3*v13) *
01186 (hn1*u21 + hn2*u22 ) +
01187 anu*(-sint*(v11*t21 + v12*t22 + v13*t23) +
01188 omcost*(v11*an1 + v12*an2 + v13*an3) -
01189 tmsint*gamma*(hn1*v11 + hn2*v12 + hn3*v13) );
01190 dCCurvTransform_(2,1) *= cosecl1;
01191
01192 dCCurvTransform_(2,2) = cost*(u11*u21 + u12*u22 ) +
01193 sint*(hu1*u21 + hu2*u22 ) +
01194 omcost*(hn1*u11 + hn2*u12 ) *
01195 (hn1*u21 + hn2*u22 ) +
01196 anu*(-sint*(u11*t21 + u12*t22 ) +
01197 omcost*(u11*an1 + u12*an2 ) -
01198 tmsint*gamma*(hn1*u11 + hn2*u12 ) );
01199 dCCurvTransform_(2,2) *= cosecl1*cosl0;
01200
01201
01202
01203
01204
01205
01206
01207
01208
01209 double pp = 1./qbp;
01210
01211 dCCurvTransform_(3,0) = - pp*(u21*dx1 + u22*dx2 );
01212 dCCurvTransform_(4,0) = - pp*(v21*dx1 + v22*dx2 + v23*dx3);
01213
01214
01215 dCCurvTransform_(3,1) = (sint*(v11*u21 + v12*u22 ) +
01216 omcost*(hv1*u21 + hv2*u22 ) +
01217 tmsint*(hn1*u21 + hn2*u22 ) *
01218 (hn1*v11 + hn2*v12 + hn3*v13))/q;
01219
01220 dCCurvTransform_(3,2) = (sint*(u11*u21 + u12*u22 ) +
01221 omcost*(hu1*u21 + hu2*u22 ) +
01222 tmsint*(hn1*u21 + hn2*u22 ) *
01223 (hn1*u11 + hn2*u12 ))*cosl0/q;
01224
01225 dCCurvTransform_(3,3) = (u11*u21 + u12*u22 );
01226
01227 dCCurvTransform_(3,4) = (v11*u21 + v12*u22 );
01228
01229
01230
01231 dCCurvTransform_(4,1) = (sint*(v11*v21 + v12*v22 + v13*v23) +
01232 omcost*(hv1*v21 + hv2*v22 + hv3*v23) +
01233 tmsint*(hn1*v21 + hn2*v22 + hn3*v23) *
01234 (hn1*v11 + hn2*v12 + hn3*v13))/q;
01235
01236 dCCurvTransform_(4,2) = (sint*(u11*v21 + u12*v22 ) +
01237 omcost*(hu1*v21 + hu2*v22 + hu3*v23) +
01238 tmsint*(hn1*v21 + hn2*v22 + hn3*v23) *
01239 (hn1*u11 + hn2*u12 ))*cosl0/q;
01240
01241 dCCurvTransform_(4,3) = (u11*v21 + u12*v22 );
01242
01243 dCCurvTransform_(4,4) = (v11*v21 + v12*v22 + v13*v23);
01244
01245 }
01246
01247 if (debug_){
01248 Basis rep; setRep(rep, tauNext);
01249 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"rep X: "<<rep.lX<<" "<<rep.lX.mag()
01250 <<"\t Y: "<<rep.lY<<" "<<rep.lY.mag()
01251 <<"\t Z: "<<rep.lZ<<" "<<rep.lZ.mag();
01252 }
01253
01254
01255
01256
01257 double mulRR = theta02*dS*dS/3.;
01258 double mulRP = theta02*fabs(dS)*p0/2.;
01259 double mulPP = theta02*p0*p0;
01260 double losPP = dP*dP*1.6/fabs(dS)*(1.0 + p0*1e-3);
01261
01262
01263
01264
01265
01266
01267 double sinThetaInv = tauNextPerpInv;
01268 double p0Mat = p0+ 0.5*dP;
01269 double p0Mat2 = p0Mat*p0Mat;
01270
01271 svCurrent.matDCovCurv*=0;
01272
01273 svCurrent.matDCovCurv(0,0) = losPP/(p0Mat2*p0Mat2);
01274
01275
01276
01277
01278
01279
01280 svCurrent.matDCovCurv(1,1) = mulPP/p0Mat2;
01281
01282
01283 svCurrent.matDCovCurv(1,4) = mulRP/p0Mat;
01284
01285
01286
01287 svCurrent.matDCovCurv(2,2) = mulPP/p0Mat2*(sinThetaInv*sinThetaInv);
01288 svCurrent.matDCovCurv(2,3) = mulRP/p0Mat*sinThetaInv;
01289
01290
01291
01292
01293 svCurrent.matDCovCurv(3,2) = mulRP/p0Mat*sinThetaInv;
01294
01295 svCurrent.matDCovCurv(3,3) = mulRR;
01296
01297
01298
01299 svCurrent.matDCovCurv(4,1) = mulRP/p0Mat;
01300
01301
01302 svCurrent.matDCovCurv(4,4) = mulRR;
01303 }
01304 break;
01305 }
01306 default:
01307 break;
01308 }
01309
01310 double pMag = curP;
01311
01312 if (dir == oppositeToMomentum) dP = fabs(dP);
01313 else if( dP < 0) {
01314 dP = -dP > pMag ? -pMag+1e-5 : dP;
01315 }
01316
01317 getNextState(svCurrent, svNext, dP, tauNext, drVec, dS, dS/radX0,
01318 dCCurvTransform_);
01319 return true;
01320 }
01321
01322 double SteppingHelixPropagator::getDeDx(const SteppingHelixPropagator::StateInfo& sv,
01323 double& dEdXPrime, double& radX0,
01324 MatBounds& rzLims) const{
01325 radX0 = 1.e24;
01326 dEdXPrime = 0.;
01327 rzLims = MatBounds();
01328 if (noMaterialMode_) return 0;
01329
01330 double dEdx = 0.;
01331
01332 double lR = sv.r3.perp();
01333 double lZ = fabs(sv.r3.z());
01334
01335
01336
01337
01338 double dEdX_HCal = 0.95;
01339 double dEdX_ECal = 0.45;
01340 double dEdX_coil = 0.35;
01341 double dEdX_Fe = 1;
01342 double dEdX_MCh = 0.053;
01343 double dEdX_Trk = 0.0114;
01344 double dEdX_Air = 2E-4;
01345 double dEdX_Vac = 0.0;
01346
01347 double radX0_HCal = 1.44/0.8;
01348 double radX0_ECal = 0.89/0.7;
01349 double radX0_coil = 4.;
01350 double radX0_Fe = 1.76;
01351 double radX0_MCh = 1e3;
01352 double radX0_Trk = 320.;
01353 double radX0_Air = 3.e4;
01354 double radX0_Vac = 3.e9;
01355
01356
01357
01358 if (! (lR < 380 && lZ < 785)){
01359 if (lZ > 785 ) rzLims = MatBounds(0, 1e4, 785, 1e4);
01360 if (lZ < 785 ) rzLims = MatBounds(380, 1e4, 0, 785);
01361 }
01362
01363
01364
01365 double ecShift = sv.r3.z() > 0 ? fabs(ecShiftPos_) : fabs(ecShiftNeg_);
01366
01367
01368
01369 if (lR < 2.9){
01370 dEdx = dEdX_Vac; radX0 = radX0_Vac;
01371 rzLims = MatBounds(0, 2.9, 0, 1E4);
01372 }
01373 else if (lR < 129){
01374 if (lZ < 294){
01375 rzLims = MatBounds(2.9, 129, 0, 294);
01376 dEdx = dEdX_Trk; radX0 = radX0_Trk;
01377
01378
01379
01380 double lEtaDet = fabs(sv.r3.eta());
01381 double scaleRadX = lEtaDet > 1.5 ? 0.7724 : 1./cosh(0.5*lEtaDet);
01382 scaleRadX *= scaleRadX;
01383 if (lEtaDet > 2 && lZ > 20) scaleRadX *= (lEtaDet-1.);
01384 if (lEtaDet > 2.5 && lZ > 20) scaleRadX *= (lEtaDet-1.);
01385 radX0 *= scaleRadX;
01386 }
01387
01388 else if ( lZ < 294 + ecShift ){
01389
01390 rzLims = MatBounds(2.9, 129, 294, 294 + ecShift);
01391 dEdx = dEdX_Air; radX0 = radX0_Air;
01392 }
01393 else if (lZ < 372 + ecShift){
01394 rzLims = MatBounds(2.9, 129, 294 + ecShift, 372 + ecShift);
01395 dEdx = dEdX_ECal; radX0 = radX0_ECal;
01396 }
01397 else if (lZ < 398 + ecShift){
01398 rzLims = MatBounds(2.9, 129, 372 + ecShift, 398 + ecShift);
01399 dEdx = dEdX_HCal*0.05; radX0 = radX0_Air;
01400 }
01401 else if (lZ < 555 + ecShift){
01402 rzLims = MatBounds(2.9, 129, 398 + ecShift, 555 + ecShift);
01403 dEdx = dEdX_HCal*0.96; radX0 = radX0_HCal/0.96;
01404 }
01405 else {
01406
01407
01408
01409 if (lZ < 568 + ecShift) rzLims = MatBounds(2.9, 129, 555 + ecShift, 568 + ecShift);
01410 else if (lZ < 625 + ecShift){
01411 if (lR < 85 + ecShift) rzLims = MatBounds(2.9, 85, 568 + ecShift, 625 + ecShift);
01412 else rzLims = MatBounds(85, 129, 568 + ecShift, 625 + ecShift);
01413 } else if (lZ < 785 + ecShift) rzLims = MatBounds(2.9, 129, 625 + ecShift, 785 + ecShift);
01414 else if (lZ < 1500 + ecShift) rzLims = MatBounds(2.9, 129, 785 + ecShift, 1500 + ecShift);
01415 else rzLims = MatBounds(2.9, 129, 1500 + ecShift, 1E4);
01416
01417
01418 if (! (lZ > 568 + ecShift && lZ < 625 + ecShift && lR > 85 )
01419 && ! (lZ > 785 + ecShift && lZ < 850 + ecShift && lR > 118)) {dEdx = dEdX_Fe; radX0 = radX0_Fe; }
01420 else { dEdx = dEdX_MCh; radX0 = radX0_MCh; }
01421 }
01422 }
01423 else if (lR < 287){
01424 if (lZ < 372 + ecShift && lR < 177){
01425 if (lZ < 304) rzLims = MatBounds(129, 177, 0, 304);
01426 else if (lZ < 343){
01427 if (lR < 135 ) rzLims = MatBounds(129, 135, 304, 343);
01428 else if (lR < 172 ){
01429 if (lZ < 311 ) rzLims = MatBounds(135, 172, 304, 311);
01430 else rzLims = MatBounds(135, 172, 311, 343);
01431 } else {
01432 if (lZ < 328) rzLims = MatBounds(172, 177, 304, 328);
01433 else rzLims = MatBounds(172, 177, 328, 343);
01434 }
01435 }
01436 else if ( lZ < 343 + ecShift){
01437 rzLims = MatBounds(129, 177, 343, 343 + ecShift);
01438 }
01439 else {
01440 if (lR < 156 ) rzLims = MatBounds(129, 156, 343 + ecShift, 372 + ecShift);
01441 else if ( (lZ - 343 - ecShift) > (lR - 156)*1.38 )
01442 rzLims = MatBounds(156, 177, 127.72 + ecShift, 372 + ecShift, atan(1.38), 0);
01443 else rzLims = MatBounds(156, 177, 343 + ecShift, 127.72 + ecShift, 0, atan(1.38));
01444 }
01445
01446 if (!(lR > 135 && lZ <343 + ecShift && lZ > 304 )
01447 && ! (lR > 156 && lZ < 372 + ecShift && lZ > 343 + ecShift && ((lZ-343. - ecShift )< (lR-156.)*1.38)))
01448 {
01449
01450 double cosThetaEquiv = 0.8/sqrt(1.+lZ*lZ/lR/lR) + 0.2;
01451 if (lZ > 343) cosThetaEquiv = 1.;
01452 dEdx = dEdX_ECal*cosThetaEquiv; radX0 = radX0_ECal/cosThetaEquiv;
01453 }
01454 else {
01455 if ( (lZ > 304 && lZ < 328 && lR < 177 && lR > 135)
01456 && ! (lZ > 311 && lR < 172) ) {dEdx = dEdX_Fe; radX0 = radX0_Fe; }
01457 else if ( lZ > 343 && lZ < 343 + ecShift) { dEdx = dEdX_Air; radX0 = radX0_Air; }
01458 else {dEdx = dEdX_ECal*0.2; radX0 = radX0_Air;}
01459 }
01460 }
01461 else if (lZ < 554 + ecShift){
01462 if (lR < 177){
01463 if ( lZ > 372 + ecShift && lZ < 398 + ecShift )rzLims = MatBounds(129, 177, 372 + ecShift, 398 + ecShift);
01464 else if (lZ < 548 + ecShift) rzLims = MatBounds(129, 177, 398 + ecShift, 548 + ecShift);
01465 else rzLims = MatBounds(129, 177, 548 + ecShift, 554 + ecShift);
01466 }
01467 else if (lR < 193){
01468 if ((lZ - 307) < (lR - 177.)*1.739) rzLims = MatBounds(177, 193, 0, -0.803, 0, atan(1.739));
01469 else if ( lZ < 389) rzLims = MatBounds(177, 193, -0.803, 389, atan(1.739), 0.);
01470 else if ( lZ < 389 + ecShift) rzLims = MatBounds(177, 193, 389, 389 + ecShift);
01471 else if ( lZ < 548 + ecShift ) rzLims = MatBounds(177, 193, 389 + ecShift, 548 + ecShift);
01472 else rzLims = MatBounds(177, 193, 548 + ecShift, 554 + ecShift);
01473 }
01474 else if (lR < 264){
01475 double anApex = 375.7278 - 193./1.327;
01476 if ( (lZ - 375.7278) < (lR - 193.)/1.327) rzLims = MatBounds(193, 264, 0, anApex, 0, atan(1./1.327));
01477 else if ( (lZ - 392.7278 ) < (lR - 193.)/1.327)
01478 rzLims = MatBounds(193, 264, anApex, anApex+17., atan(1./1.327), atan(1./1.327));
01479 else if ( (lZ - 392.7278 - ecShift ) < (lR - 193.)/1.327)
01480 rzLims = MatBounds(193, 264, anApex+17., anApex+17. + ecShift, atan(1./1.327), atan(1./1.327));
01481
01482 else if ( lZ < 517 + ecShift ) rzLims = MatBounds(193, 264, anApex+17. + ecShift, 517 + ecShift, atan(1./1.327), 0);
01483 else if (lZ < 548 + ecShift){
01484 if (lR < 246 ) rzLims = MatBounds(193, 246, 517 + ecShift, 548 + ecShift);
01485 else rzLims = MatBounds(246, 264, 517 + ecShift, 548 + ecShift);
01486 }
01487 else rzLims = MatBounds(193, 264, 548 + ecShift, 554 + ecShift);
01488 }
01489 else if ( lR < 275){
01490 if (lZ < 433) rzLims = MatBounds(264, 275, 0, 433);
01491 else if (lZ < 554 ) rzLims = MatBounds(264, 275, 433, 554);
01492 else rzLims = MatBounds(264, 275, 554, 554 + ecShift);
01493 }
01494 else {
01495 if (lZ < 402) rzLims = MatBounds(275, 287, 0, 402);
01496 else if (lZ < 554) rzLims = MatBounds(275, 287, 402, 554);
01497 else rzLims = MatBounds(275, 287, 554, 554 + ecShift);
01498 }
01499
01500 if ((lZ < 433 || lR < 264) && (lZ < 402 || lR < 275) && (lZ < 517 + ecShift || lR < 246)
01501
01502 && lZ < 548 + ecShift
01503 && ! (lZ < 389 + ecShift && lZ > 335 && lR < 193 )
01504 && ! (lZ > 307 && lZ < 335 && lR < 193 && ((lZ - 307) > (lR - 177.)*1.739))
01505 && ! (lR < 177 && lZ < 398 + ecShift)
01506 && ! (lR < 264 && lR > 193 && fabs(441.5+0.5*ecShift - lZ + (lR - 269.)/1.327) < (8.5 + ecShift*0.5)) )
01507 { dEdx = dEdX_HCal; radX0 = radX0_HCal;
01508 }
01509 else if( (lR < 193 && lZ > 389 && lZ < 389 + ecShift ) || (lR > 264 && lR < 287 && lZ > 554 && lZ < 554 + ecShift)
01510 ||(lR < 264 && lR > 193 && fabs(441.5+8.5+0.5*ecShift - lZ + (lR - 269.)/1.327) < ecShift*0.5) ) {
01511 dEdx = dEdX_Air; radX0 = radX0_Air;
01512 }
01513 else {dEdx = dEdX_HCal*0.05; radX0 = radX0_Air; }
01514 }
01515
01516 else if (lZ < 564 + ecShift){
01517 if (lR < 251) {
01518 rzLims = MatBounds(129, 251, 554 + ecShift, 564 + ecShift);
01519 dEdx = dEdX_Fe; radX0 = radX0_Fe;
01520 }
01521 else {
01522 rzLims = MatBounds(251, 287, 554 + ecShift, 564 + ecShift);
01523 dEdx = dEdX_MCh; radX0 = radX0_MCh;
01524 }
01525 }
01526 else if (lZ < 625 + ecShift){
01527 rzLims = MatBounds(129, 287, 564 + ecShift, 625 + ecShift);
01528 dEdx = dEdX_MCh; radX0 = radX0_MCh;
01529 }
01530 else if (lZ < 785 + ecShift){
01531 if (lR < 275) rzLims = MatBounds(129, 275, 625 + ecShift, 785 + ecShift);
01532 else {
01533 if (lZ < 720 + ecShift) rzLims = MatBounds(275, 287, 625 + ecShift, 720 + ecShift);
01534 else rzLims = MatBounds(275, 287, 720 + ecShift, 785 + ecShift);
01535 }
01536 if (! (lR > 275 && lZ < 720 + ecShift)) { dEdx = dEdX_Fe; radX0 = radX0_Fe; }
01537 else { dEdx = dEdX_MCh; radX0 = radX0_MCh; }
01538 }
01539 else if (lZ < 850 + ecShift){
01540 rzLims = MatBounds(129, 287, 785 + ecShift, 850 + ecShift);
01541 dEdx = dEdX_MCh; radX0 = radX0_MCh;
01542 }
01543 else if (lZ < 910 + ecShift){
01544 rzLims = MatBounds(129, 287, 850 + ecShift, 910 + ecShift);
01545 dEdx = dEdX_Fe; radX0 = radX0_Fe;
01546 }
01547 else if (lZ < 975 + ecShift){
01548 rzLims = MatBounds(129, 287, 910 + ecShift, 975 + ecShift);
01549 dEdx = dEdX_MCh; radX0 = radX0_MCh;
01550 }
01551 else if (lZ < 1000 + ecShift){
01552 rzLims = MatBounds(129, 287, 975 + ecShift, 1000 + ecShift);
01553 dEdx = dEdX_Fe; radX0 = radX0_Fe;
01554 }
01555 else if (lZ < 1063 + ecShift){
01556 rzLims = MatBounds(129, 287, 1000 + ecShift, 1063 + ecShift);
01557 dEdx = dEdX_MCh; radX0 = radX0_MCh;
01558 }
01559 else if ( lZ < 1073 + ecShift){
01560 rzLims = MatBounds(129, 287, 1063 + ecShift, 1073 + ecShift);
01561 dEdx = dEdX_Fe; radX0 = radX0_Fe;
01562 }
01563 else if (lZ < 1E4 ) {
01564 rzLims = MatBounds(129, 287, 1073 + ecShift, 1E4);
01565 dEdx = dEdX_Air; radX0 = radX0_Air;
01566 }
01567 else {
01568
01569 dEdx = dEdX_Air; radX0 = radX0_Air;
01570 }
01571 }
01572 else if (lR <380 && lZ < 667){
01573 if (lZ < 630){
01574 if (lR < 315) rzLims = MatBounds(287, 315, 0, 630);
01575 else if (lR < 341 ) rzLims = MatBounds(315, 341, 0, 630);
01576 else rzLims = MatBounds(341, 380, 0, 630);
01577 } else rzLims = MatBounds(287, 380, 630, 667);
01578
01579 if (lZ < 630) { dEdx = dEdX_coil; radX0 = radX0_coil; }
01580 else {dEdx = dEdX_Air; radX0 = radX0_Air; }
01581 }
01582 else {
01583 if (lZ < 667) {
01584 if (lR < 850){
01585 bool isIron = false;
01586
01587 if (useIsYokeFlag_ && useMagVolumes_ && sv.magVol != 0){
01588 isIron = sv.isYokeVol;
01589 } else {
01590 double bMag = sv.bf.mag();
01591 isIron = (bMag > 0.75 && ! (lZ > 500 && lR <500 && bMag < 1.15)
01592 && ! (lZ < 450 && lR > 420 && bMag < 1.15 ) );
01593 }
01594
01595 rzLims = MatBounds(380, 850, 0, 667);
01596 if (isIron) { dEdx = dEdX_Fe; radX0 = radX0_Fe; }
01597 else { dEdx = dEdX_MCh; radX0 = radX0_MCh; }
01598 } else {
01599 rzLims = MatBounds(850, 1E4, 0, 667);
01600 dEdx = dEdX_Air; radX0 = radX0_Air;
01601 }
01602 }
01603 else if (lR > 750 ){
01604 rzLims = MatBounds(750, 1E4, 667, 1E4);
01605 dEdx = dEdX_Air; radX0 = radX0_Air;
01606 }
01607 else if (lZ < 667 + ecShift){
01608 rzLims = MatBounds(287, 750, 667, 667 + ecShift);
01609 dEdx = dEdX_Air; radX0 = radX0_Air;
01610 }
01611
01612 else if (lZ < 724 + ecShift){
01613 if (lR < 380 ) rzLims = MatBounds(287, 380, 667 + ecShift, 724 + ecShift);
01614 else rzLims = MatBounds(380, 750, 667 + ecShift, 724 + ecShift);
01615 dEdx = dEdX_MCh; radX0 = radX0_MCh;
01616 }
01617 else if (lZ < 785 + ecShift){
01618 if (lR < 380 ) rzLims = MatBounds(287, 380, 724 + ecShift, 785 + ecShift);
01619 else rzLims = MatBounds(380, 750, 724 + ecShift, 785 + ecShift);
01620 dEdx = dEdX_Fe; radX0 = radX0_Fe;
01621 }
01622 else if (lZ < 850 + ecShift){
01623 rzLims = MatBounds(287, 750, 785 + ecShift, 850 + ecShift);
01624 dEdx = dEdX_MCh; radX0 = radX0_MCh;
01625 }
01626 else if (lZ < 910 + ecShift){
01627 rzLims = MatBounds(287, 750, 850 + ecShift, 910 + ecShift);
01628 dEdx = dEdX_Fe; radX0 = radX0_Fe;
01629 }
01630 else if (lZ < 975 + ecShift){
01631 rzLims = MatBounds(287, 750, 910 + ecShift, 975 + ecShift);
01632 dEdx = dEdX_MCh; radX0 = radX0_MCh;
01633 }
01634 else if (lZ < 1000 + ecShift){
01635 rzLims = MatBounds(287, 750, 975 + ecShift, 1000 + ecShift);
01636 dEdx = dEdX_Fe; radX0 = radX0_Fe;
01637 }
01638 else if (lZ < 1063 + ecShift){
01639 if (lR < 360){
01640 rzLims = MatBounds(287, 360, 1000 + ecShift, 1063 + ecShift);
01641 dEdx = dEdX_MCh; radX0 = radX0_MCh;
01642 }
01643
01644 else {
01645 rzLims = MatBounds(360, 750, 1000 + ecShift, 1063 + ecShift);
01646 dEdx = dEdX_Air; radX0 = radX0_Air;
01647 }
01648 }
01649 else if ( lZ < 1073 + ecShift){
01650 rzLims = MatBounds(287, 750, 1063 + ecShift, 1073 + ecShift);
01651
01652 dEdx = dEdX_Air; radX0 = radX0_Air;
01653 }
01654 else if (lZ < 1E4 ) {
01655 rzLims = MatBounds(287, 750, 1073 + ecShift, 1E4);
01656 dEdx = dEdX_Air; radX0 = radX0_Air;
01657 }
01658 else {dEdx = dEdX_Air; radX0 = radX0_Air; }
01659 }
01660
01661
01662
01663
01664 double p0 = sv.p3.mag();
01665 double logp0 = log(p0);
01666 double p0powN33 = 0;
01667 if (p0>3.) {
01668
01669 double xx=1./p0; xx=sqrt(sqrt(sqrt(sqrt(xx)))); xx=xx*xx*xx*xx*xx;
01670 p0powN33 = xx;
01671 }
01672 double dEdX_mat = -(11.4 + 0.96*fabs(logp0+log(2.8)) + 0.033*p0*(1.0 - p0powN33) )*1e-3;
01673
01674
01675 dEdXPrime = dEdx == 0 ? 0 : -dEdx*(0.96/p0 + 0.033 - 0.022*p0powN33)*1e-3;
01676 dEdx *= dEdX_mat;
01677
01678 return dEdx;
01679 }
01680
01681
01682 int SteppingHelixPropagator::cIndex_(int ind) const{
01683 int result = ind%MAX_POINTS;
01684 if (ind != 0 && result == 0){
01685 result = MAX_POINTS;
01686 }
01687 return result;
01688 }
01689
01690 SteppingHelixPropagator::Result
01691 SteppingHelixPropagator::refToDest(SteppingHelixPropagator::DestType dest,
01692 const SteppingHelixPropagator::StateInfo& sv,
01693 const double pars[6],
01694 double& dist, double& tanDist,
01695 PropagationDirection& refDirection,
01696 double fastSkipDist) const{
01697 static const std::string metname = "SteppingHelixPropagator";
01698 Result result = SteppingHelixStateInfo::NOT_IMPLEMENTED;
01699
01700 switch (dest){
01701 case RADIUS_DT:
01702 {
01703 double curR = sv.r3.perp();
01704 dist = pars[RADIUS_P] - curR;
01705 if (fabs(dist) > fastSkipDist){
01706 result = SteppingHelixStateInfo::INACC;
01707 break;
01708 }
01709 double curP2 = sv.p3.mag2();
01710 double curPtPos2 = sv.p3.perp2(); if(curPtPos2< 1e-16) curPtPos2=1e-16;
01711
01712 double cosDPhiPR = (sv.r3.x()*sv.p3.x()+sv.r3.y()*sv.p3.y());
01713 refDirection = dist*cosDPhiPR > 0 ?
01714 alongMomentum : oppositeToMomentum;
01715 tanDist = dist*sqrt(curP2/curPtPos2);
01716 result = SteppingHelixStateInfo::OK;
01717 }
01718 break;
01719 case Z_DT:
01720 {
01721 double curZ = sv.r3.z();
01722 dist = pars[Z_P] - curZ;
01723 if (fabs(dist) > fastSkipDist){
01724 result = SteppingHelixStateInfo::INACC;
01725 break;
01726 }
01727 double curP = sv.p3.mag();
01728 refDirection = sv.p3.z()*dist > 0. ?
01729 alongMomentum : oppositeToMomentum;
01730 tanDist = dist/sv.p3.z()*curP;
01731 result = SteppingHelixStateInfo::OK;
01732 }
01733 break;
01734 case PLANE_DT:
01735 {
01736 Point rPlane(pars[0], pars[1], pars[2]);
01737 Vector nPlane(pars[3], pars[4], pars[5]);
01738
01739
01740
01741
01742
01743 double dRDotN = (sv.r3.x()-rPlane.x())*nPlane.x() + (sv.r3.y()-rPlane.y())*nPlane.y() + (sv.r3.z()-rPlane.z())*nPlane.z();
01744
01745 dist = fabs(dRDotN);
01746 if (dist > fastSkipDist){
01747 result = SteppingHelixStateInfo::INACC;
01748 break;
01749 }
01750 double curP = sv.p3.mag();
01751 double p0 = curP;
01752 double p0Inv = 1./p0;
01753 Vector tau(sv.p3); tau *=p0Inv;
01754 double tN = tau.dot(nPlane);
01755 refDirection = tN*dRDotN < 0. ?
01756 alongMomentum : oppositeToMomentum;
01757 double b0 = sv.bf.mag();
01758 if (fabs(tN)>1e-24){
01759 tanDist = -dRDotN/tN;
01760 } else {
01761 tN = 1e-24;
01762 if (fabs(dRDotN)>1e-24) tanDist = 1e6;
01763 else tanDist = 1;
01764 }
01765 if (fabs(tanDist) > 1e4) tanDist = 1e4;
01766 if (b0>1.5e-6){
01767 double b0Inv = 1./b0;
01768 double tNInv = 1./tN;
01769 Vector bHat(sv.bf); bHat *=b0Inv;
01770 double bHatN = bHat.dot(nPlane);
01771 double cosPB = bHat.dot(tau);
01772 double kVal = 2.99792458e-3*sv.q*p0Inv*b0;
01773 double aVal = tanDist*kVal;
01774 Vector lVec = bHat.cross(tau);
01775 double bVal = lVec.dot(nPlane)*tNInv;
01776 if (fabs(aVal*bVal)< 0.3){
01777 double cVal = 1. - bHatN*cosPB*tNInv;
01778 double aacVal = cVal*aVal*aVal;
01779 if (fabs(aacVal)<1){
01780 double tanDCorr = bVal/2. + (bVal*bVal/2. + cVal/6)*aVal;
01781 tanDCorr *= aVal;
01782
01783 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<tanDist<<" vs "
01784 <<tanDist*(1.+tanDCorr)<<" corr "<<tanDist*tanDCorr<<std::endl;
01785 tanDist *= (1.+tanDCorr);
01786 } else {
01787 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"AACVal "<< fabs(aacVal)
01788 <<" = "<<aVal<<"**2 * "<<cVal<<" too large:: will not converge"<<std::endl;
01789 }
01790 } else {
01791 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"ABVal "<< fabs(aVal*bVal)
01792 <<" = "<<aVal<<" * "<<bVal<<" too large:: will not converge"<<std::endl;
01793 }
01794 }
01795 result = SteppingHelixStateInfo::OK;
01796 }
01797 break;
01798 case CONE_DT:
01799 {
01800 Point cVertex(pars[0], pars[1], pars[2]);
01801 if (cVertex.perp2() < 1e-10){
01802
01803 Vector relV3 = sv.r3 - cVertex;
01804 double relV3mag = relV3.mag();
01805
01806 double theta(pars[3]);
01807
01808 double sinTheta = sin(theta);
01809 double cosTheta = cos(theta);
01810 double cosV3Theta = relV3.z()/relV3mag;
01811 if (cosV3Theta>1) cosV3Theta=1;
01812 if (cosV3Theta<-1) cosV3Theta=-1;
01813 double sinV3Theta = sqrt(1.-cosV3Theta*cosV3Theta);
01814
01815 double sinDTheta = sinTheta*cosV3Theta - cosTheta*sinV3Theta;
01816 double cosDTheta = cosTheta*cosV3Theta + sinTheta*sinV3Theta;
01817 bool isInside = sinTheta > sinV3Theta && cosTheta*cosV3Theta > 0;
01818 dist = isInside || cosDTheta > 0 ? relV3mag*sinDTheta : relV3mag;
01819 if (fabs(dist) > fastSkipDist){
01820 result = SteppingHelixStateInfo::INACC;
01821 break;
01822 }
01823
01824 double relV3phi=relV3.phi();
01825 double normPhi = isInside ?
01826 Geom::pi() + relV3phi : relV3phi;
01827 double normTheta = theta > Geom::pi()/2. ?
01828 ( isInside ? 1.5*Geom::pi() - theta : theta - Geom::pi()/2. )
01829 : ( isInside ? Geom::pi()/2. - theta : theta + Geom::pi()/2 );
01830
01831 Vector norm; norm.setRThetaPhi(fabs(dist), normTheta, normPhi);
01832 double curP = sv.p3.mag(); double cosp3theta = sv.p3.z()/curP;
01833 if (cosp3theta>1) cosp3theta=1;
01834 if (cosp3theta<-1) cosp3theta=-1;
01835 double sineConeP = sinTheta*cosp3theta - cosTheta*sqrt(1.-cosp3theta*cosp3theta);
01836
01837 double sinSolid = norm.dot(sv.p3)/(fabs(dist)*curP);
01838 tanDist = fabs(sinSolid) > fabs(sineConeP) ? dist/fabs(sinSolid) : dist/fabs(sineConeP);
01839
01840
01841 refDirection = sinSolid > 0 ?
01842 oppositeToMomentum : alongMomentum;
01843 if (debug_){
01844 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Cone pars: theta "<<theta
01845 <<" normTheta "<<norm.theta()
01846 <<" p3Theta "<<sv.p3.theta()
01847 <<" sinD: "<< sineConeP
01848 <<" sinSolid "<<sinSolid;
01849 }
01850 if (debug_){
01851 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"refToDest:toCone the point is "
01852 <<(isInside? "in" : "out")<<"side the cone"
01853 <<std::endl;
01854 }
01855 }
01856 result = SteppingHelixStateInfo::OK;
01857 }
01858 break;
01859
01860
01861 case PATHL_DT:
01862 {
01863 double curS = fabs(sv.path());
01864 dist = pars[PATHL_P] - curS;
01865 if (fabs(dist) > fastSkipDist){
01866 result = SteppingHelixStateInfo::INACC;
01867 break;
01868 }
01869 refDirection = pars[PATHL_P] > 0 ?
01870 alongMomentum : oppositeToMomentum;
01871 tanDist = dist;
01872 result = SteppingHelixStateInfo::OK;
01873 }
01874 break;
01875 case POINT_PCA_DT:
01876 {
01877 Point pDest(pars[0], pars[1], pars[2]);
01878 double curP = sv.p3.mag();
01879 dist = (sv.r3 - pDest).mag()+ 1e-24;
01880 tanDist = (sv.r3.dot(sv.p3) - pDest.dot(sv.p3))/curP;
01881
01882 double b0 = sv.bf.mag();
01883 if (b0>1.5e-6){
01884 double p0 = curP;
01885 double kVal = 2.99792458e-3*sv.q/p0*b0;
01886 double aVal = fabs(dist*kVal);
01887 tanDist *= 1./(1.+ aVal);
01888 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"corrected by aVal "<<aVal<<" to "<<tanDist;
01889 }
01890 refDirection = tanDist < 0 ?
01891 alongMomentum : oppositeToMomentum;
01892 result = SteppingHelixStateInfo::OK;
01893 }
01894 break;
01895 case LINE_PCA_DT:
01896 {
01897 Point rLine(pars[0], pars[1], pars[2]);
01898 Vector dLine(pars[3], pars[4], pars[5]);
01899 dLine = (dLine - rLine);
01900 dLine *= 1./dLine.mag(); if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"dLine "<<dLine;
01901
01902 Vector dR = sv.r3 - rLine;
01903 double curP = sv.p3.mag();
01904 Vector dRPerp = dR - dLine*(dR.dot(dLine));
01905 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"dRperp "<<dRPerp;
01906
01907 dist = dRPerp.mag() + 1e-24;
01908 tanDist = dRPerp.dot(sv.p3)/curP;
01909
01910 double cosAlpha2 = dLine.dot(sv.p3)/curP; cosAlpha2 *= cosAlpha2;
01911 tanDist *= 1./sqrt(fabs(1.-cosAlpha2)+1e-96);
01912
01913
01914 double b0 = sv.bf.mag();
01915 if (b0>1.5e-6){
01916 double p0 = curP;
01917 double kVal = 2.99792458e-3*sv.q/p0*b0;
01918 double aVal = fabs(dist*kVal);
01919 tanDist *= 1./(1.+ aVal);
01920 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"corrected by aVal "<<aVal<<" to "<<tanDist;
01921 }
01922 refDirection = tanDist < 0 ?
01923 alongMomentum : oppositeToMomentum;
01924 result = SteppingHelixStateInfo::OK;
01925 }
01926 break;
01927 default:
01928 {
01929
01930 dist = 1e12;
01931 tanDist = 1e12;
01932 refDirection = anyDirection;
01933 result = SteppingHelixStateInfo::NOT_IMPLEMENTED;
01934 }
01935 break;
01936 }
01937
01938 double tanDistConstrained = tanDist;
01939 if (fabs(tanDist) > 4.*fabs(dist) ) tanDistConstrained *= tanDist == 0 ? 0 : fabs(dist/tanDist*4.);
01940
01941 if (debug_){
01942 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"refToDest input: dest"<<dest<<" pars[]: ";
01943 for (int i = 0; i < 6; i++){
01944 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<", "<<i<<" "<<pars[i];
01945 }
01946 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<std::endl;
01947 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"refToDest output: "
01948 <<"\t dist" << dist
01949 <<"\t tanDist "<< tanDist
01950 <<"\t tanDistConstr "<< tanDistConstrained
01951 <<"\t refDirection "<< refDirection
01952 <<std::endl;
01953 }
01954 tanDist = tanDistConstrained;
01955
01956 return result;
01957 }
01958
01959 SteppingHelixPropagator::Result
01960 SteppingHelixPropagator::refToMagVolume(const SteppingHelixPropagator::StateInfo& sv,
01961 PropagationDirection dir,
01962 double& dist, double& tanDist,
01963 double fastSkipDist, bool expectNewMagVolume, double maxStep) const{
01964
01965 static const std::string metname = "SteppingHelixPropagator";
01966 Result result = SteppingHelixStateInfo::NOT_IMPLEMENTED;
01967 const MagVolume* cVol = sv.magVol;
01968
01969 if (cVol == 0) return result;
01970 const std::vector<VolumeSide>& cVolFaces(cVol->faces());
01971
01972 double distToFace[6] = {0,0,0,0,0,0};
01973 double tanDistToFace[6] = {0,0,0,0,0,0};
01974 PropagationDirection refDirectionToFace[6] = {anyDirection, anyDirection, anyDirection, anyDirection, anyDirection, anyDirection};
01975 Result resultToFace[6] = {result, result, result, result, result, result};
01976 int iFDest = -1;
01977 int iDistMin = -1;
01978
01979 unsigned int iFDestSorted[6] = {0,0,0,0,0,0};
01980 unsigned int nDestSorted =0;
01981 unsigned int nearParallels = 0;
01982
01983 double curP = sv.p3.mag();
01984
01985 if (debug_){
01986 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Trying volume "<<DDSolidShapesName::name(cVol->shapeType())
01987 <<" with "<<cVolFaces.size()<<" faces"<<std::endl;
01988 }
01989
01990 unsigned int nFaces = cVolFaces.size();
01991 for (unsigned int iFace = 0; iFace < nFaces; ++iFace){
01992 if (iFace > 5){
01993 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Too many faces"<<std::endl;
01994 }
01995 if (debug_){
01996 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Start with face "<<iFace<<std::endl;
01997 }
01998
01999
02000
02001 const Surface* cPlane = 0;
02002 const Cylinder* cCyl = 0;
02003 const Cone* cCone = 0;
02004 if (typeid(cVolFaces[iFace].surface()) == typeid(const Plane&)){
02005 cPlane = &cVolFaces[iFace].surface();
02006 } else if (typeid(cVolFaces[iFace].surface()) == typeid(const Cylinder&)){
02007 cCyl = dynamic_cast<const Cylinder*>(&cVolFaces[iFace].surface());
02008 } else if (typeid(cVolFaces[iFace].surface()) == typeid(const Cone&)){
02009 cCone = dynamic_cast<const Cone*>(&cVolFaces[iFace].surface());
02010 } else {
02011 edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Could not cast a volume side surface to a known type"<<std::endl;
02012 }
02013
02014 if (debug_){
02015 if (cPlane!=0) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The face is a plane at "<<cPlane<<std::endl;
02016 if (cCyl!=0) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The face is a cylinder at "<<cCyl<<std::endl;
02017 }
02018
02019 double pars[6] = {0,0,0,0,0,0};
02020 DestType dType = UNDEFINED_DT;
02021 if (cPlane != 0){
02022 Point rPlane(cPlane->position().x(),cPlane->position().y(),cPlane->position().z());
02023
02024 Vector nPlane(cPlane->rotation().zx(), cPlane->rotation().zy(), cPlane->rotation().zz()); nPlane /= nPlane.mag();
02025
02026 if (sv.r3.z() < 0 || !vbField_->isZSymmetric() ){
02027 pars[0] = rPlane.x(); pars[1] = rPlane.y(); pars[2] = rPlane.z();
02028 pars[3] = nPlane.x(); pars[4] = nPlane.y(); pars[5] = nPlane.z();
02029 } else {
02030 pars[0] = rPlane.x(); pars[1] = rPlane.y(); pars[2] = -rPlane.z();
02031 pars[3] = nPlane.x(); pars[4] = nPlane.y(); pars[5] = -nPlane.z();
02032 }
02033 dType = PLANE_DT;
02034 } else if (cCyl != 0){
02035 if (debug_){
02036 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Cylinder at "<<cCyl->position()
02037 <<" rorated by "<<cCyl->rotation()
02038 <<std::endl;
02039 }
02040 pars[RADIUS_P] = cCyl->radius();
02041 dType = RADIUS_DT;
02042 } else if (cCone != 0){
02043 if (debug_){
02044 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Cone at "<<cCone->position()
02045 <<" rorated by "<<cCone->rotation()
02046 <<" vertex at "<<cCone->vertex()
02047 <<" angle of "<<cCone->openingAngle()
02048 <<std::endl;
02049 }
02050 if (sv.r3.z() < 0 || !vbField_->isZSymmetric()){
02051 pars[0] = cCone->vertex().x(); pars[1] = cCone->vertex().y();
02052 pars[2] = cCone->vertex().z();
02053 pars[3] = cCone->openingAngle();
02054 } else {
02055 pars[0] = cCone->vertex().x(); pars[1] = cCone->vertex().y();
02056 pars[2] = -cCone->vertex().z();
02057 pars[3] = Geom::pi() - cCone->openingAngle();
02058 }
02059 dType = CONE_DT;
02060 } else {
02061 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Unknown surface"<<std::endl;
02062 resultToFace[iFace] = SteppingHelixStateInfo::UNDEFINED;
02063 continue;
02064 }
02065 resultToFace[iFace] =
02066 refToDest(dType, sv, pars,
02067 distToFace[iFace], tanDistToFace[iFace], refDirectionToFace[iFace], fastSkipDist);
02068
02069
02070 if (resultToFace[iFace] != SteppingHelixStateInfo::OK){
02071 if (resultToFace[iFace] == SteppingHelixStateInfo::INACC) result = SteppingHelixStateInfo::INACC;
02072 }
02073
02074
02075
02076
02077 if (resultToFace[iFace] == SteppingHelixStateInfo::OK){
02078 double invDTFPosiF = 1./(1e-32+fabs(tanDistToFace[iFace]));
02079 double dSlope = fabs(distToFace[iFace])*invDTFPosiF;
02080 double maxStepL = maxStep> 100 ? 100 : maxStep; if (maxStepL < 10) maxStepL = 10.;
02081 bool isNearParallel = fabs(tanDistToFace[iFace]) + 100.*curP*dSlope < maxStepL
02082
02083 && dSlope < 0.15 ;
02084 if (refDirectionToFace[iFace] == dir || isNearParallel){
02085 if (isNearParallel) nearParallels++;
02086 iFDestSorted[nDestSorted] = iFace;
02087 nDestSorted++;
02088 }
02089 }
02090
02091
02092 if (refDirectionToFace[iFace] == dir){
02093 if (iDistMin == -1) iDistMin = iFace;
02094 else if (fabs(distToFace[iFace]) < fabs(distToFace[iDistMin])) iDistMin = iFace;
02095 }
02096 if (debug_)
02097 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<cVol<<" "<<iFace<<" "
02098 <<tanDistToFace[iFace]<<" "<<distToFace[iFace]<<" "<<refDirectionToFace[iFace]<<" "<<dir<<std::endl;
02099
02100 }
02101
02102 for (unsigned int i = 0;i<nDestSorted; ++i){
02103 unsigned int iMax = nDestSorted-i-1;
02104 for (unsigned int j=0;j<nDestSorted-i; ++j){
02105 if (fabs(tanDistToFace[iFDestSorted[j]]) > fabs(tanDistToFace[iFDestSorted[iMax]]) ){
02106 iMax = j;
02107 }
02108 }
02109 unsigned int iTmp = iFDestSorted[nDestSorted-i-1];
02110 iFDestSorted[nDestSorted-i-1] = iFDestSorted[iMax];
02111 iFDestSorted[iMax] = iTmp;
02112 }
02113
02114 if (debug_){
02115 for (unsigned int i=0;i<nDestSorted;++i){
02116 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<cVol<<" "<<i<<" "<<iFDestSorted[i]<<" "<<tanDistToFace[iFDestSorted[i]]<<std::endl;
02117 }
02118 }
02119
02120
02121
02122
02123 for (unsigned int i=0; i<nDestSorted;++i){
02124 iFDest = iFDestSorted[i];
02125
02126 double sign = dir == alongMomentum ? 1. : -1.;
02127 Point gPointEst(sv.r3);
02128 Vector lDelta(sv.p3); lDelta *= sign/curP*fabs(distToFace[iFDest]);
02129 gPointEst += lDelta;
02130 if (debug_){
02131 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Linear est point "<<gPointEst
02132 <<" for iFace "<<iFDest<<std::endl;
02133 }
02134 GlobalPoint gPointEstNegZ(gPointEst.x(), gPointEst.y(), -fabs(gPointEst.z()));
02135 GlobalPoint gPointEstNorZ(gPointEst.x(), gPointEst.y(), gPointEst.z() );
02136 if ( (vbField_->isZSymmetric() && cVol->inside(gPointEstNegZ))
02137 || ( !vbField_->isZSymmetric() && cVol->inside(gPointEstNorZ) ) ){
02138 if (debug_){
02139 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is inside the volume"<<std::endl;
02140 }
02141
02142 result = SteppingHelixStateInfo::OK;
02143 dist = distToFace[iFDest];
02144 tanDist = tanDistToFace[iFDest];
02145 if (debug_){
02146 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got a point near closest boundary -- face "<<iFDest<<std::endl;
02147 }
02148 break;
02149 }
02150 }
02151
02152 if (result != SteppingHelixStateInfo::OK && expectNewMagVolume){
02153 double sign = dir == alongMomentum ? 1. : -1.;
02154
02155
02156 if (nDestSorted-nearParallels > 0 ) result = SteppingHelixStateInfo::WRONG_VOLUME;
02157 else {
02158
02159 Point gPointEst(sv.r3);
02160 double lDist = iDistMin == -1 ? fastSkipDist : fabs(distToFace[iDistMin]);
02161 if (lDist > fastSkipDist) lDist = fastSkipDist;
02162 Vector lDelta(sv.p3); lDelta *= sign/curP*lDist;
02163 gPointEst += lDelta;
02164 if (debug_){
02165 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Linear est point to shortest dist "<<gPointEst
02166 <<" for iFace "<<iDistMin<<" at distance "<<lDist*sign<<std::endl;
02167 }
02168 GlobalPoint gPointEstNegZ(gPointEst.x(), gPointEst.y(), -fabs(gPointEst.z()));
02169 GlobalPoint gPointEstNorZ(gPointEst.x(), gPointEst.y(), gPointEst.z() );
02170 if ( (vbField_->isZSymmetric() && cVol->inside(gPointEstNegZ))
02171 || ( !vbField_->isZSymmetric() && cVol->inside(gPointEstNorZ) ) ){
02172 if (debug_){
02173 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is inside the volume"<<std::endl;
02174 }
02175
02176 }else {
02177 result = SteppingHelixStateInfo::WRONG_VOLUME;
02178 }
02179 }
02180
02181 if (result == SteppingHelixStateInfo::WRONG_VOLUME){
02182 dist = sign*0.05;
02183 tanDist = dist*1.01;
02184 if( debug_){
02185 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Wrong volume located: return small dist, tandist"<<std::endl;
02186 }
02187 }
02188 }
02189
02190 if (result == SteppingHelixStateInfo::INACC){
02191 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"All faces are too far"<<std::endl;
02192 } else if (result == SteppingHelixStateInfo::WRONG_VOLUME){
02193 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Appear to be in a wrong volume"<<std::endl;
02194 } else if (result != SteppingHelixStateInfo::OK){
02195 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Something else went wrong"<<std::endl;
02196 }
02197
02198
02199 return result;
02200 }
02201
02202
02203 SteppingHelixPropagator::Result
02204 SteppingHelixPropagator::refToMatVolume(const SteppingHelixPropagator::StateInfo& sv,
02205 PropagationDirection dir,
02206 double& dist, double& tanDist,
02207 double fastSkipDist) const{
02208
02209 static const std::string metname = "SteppingHelixPropagator";
02210 Result result = SteppingHelixStateInfo::NOT_IMPLEMENTED;
02211
02212 double parLim[6] = {sv.rzLims.rMin, sv.rzLims.rMax,
02213 sv.rzLims.zMin, sv.rzLims.zMax,
02214 sv.rzLims.th1, sv.rzLims.th2 };
02215
02216 double distToFace[4] = {0,0,0,0};
02217 double tanDistToFace[4] = {0,0,0,0};
02218 PropagationDirection refDirectionToFace[4] = {anyDirection, anyDirection, anyDirection, anyDirection};
02219 Result resultToFace[4] = {result, result, result, result};
02220 int iFDest = -1;
02221
02222 double curP = sv.p3.mag();
02223
02224 for (unsigned int iFace = 0; iFace < 4; iFace++){
02225 if (debug_){
02226 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Start with mat face "<<iFace<<std::endl;
02227 }
02228
02229 double pars[6] = {0,0,0,0,0,0};
02230 DestType dType = UNDEFINED_DT;
02231 if (iFace > 1){
02232 if (fabs(parLim[iFace+2])< 1e-6){
02233 if (sv.r3.z() < 0){
02234 pars[0] = 0; pars[1] = 0; pars[2] = -parLim[iFace];
02235 pars[3] = 0; pars[4] = 0; pars[5] = 1;
02236 } else {
02237 pars[0] = 0; pars[1] = 0; pars[2] = parLim[iFace];
02238 pars[3] = 0; pars[4] = 0; pars[5] = 1;
02239 }
02240 dType = PLANE_DT;
02241 } else {
02242 if (sv.r3.z() > 0){
02243 pars[0] = 0; pars[1] = 0;
02244 pars[2] = parLim[iFace];
02245 pars[3] = Geom::pi()/2. - parLim[iFace+2];
02246 } else {
02247 pars[0] = 0; pars[1] = 0;
02248 pars[2] = - parLim[iFace];
02249 pars[3] = Geom::pi()/2. + parLim[iFace+2];
02250 }
02251 dType = CONE_DT;
02252 }
02253 } else {
02254 pars[RADIUS_P] = parLim[iFace];
02255 dType = RADIUS_DT;
02256 }
02257
02258 resultToFace[iFace] =
02259 refToDest(dType, sv, pars,
02260 distToFace[iFace], tanDistToFace[iFace], refDirectionToFace[iFace], fastSkipDist);
02261
02262 if (resultToFace[iFace] != SteppingHelixStateInfo::OK){
02263 if (resultToFace[iFace] == SteppingHelixStateInfo::INACC) result = SteppingHelixStateInfo::INACC;
02264 continue;
02265 }
02266 if (refDirectionToFace[iFace] == dir || fabs(distToFace[iFace]) < 2e-2*fabs(tanDistToFace[iFace]) ){
02267 double sign = dir == alongMomentum ? 1. : -1.;
02268 Point gPointEst(sv.r3);
02269 Vector lDelta(sv.p3); lDelta *= sign*fabs(distToFace[iFace])/curP;
02270 gPointEst += lDelta;
02271 if (debug_){
02272 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Linear est point "<<gPointEst
02273 <<std::endl;
02274 }
02275 double lZ = fabs(gPointEst.z());
02276 double lR = gPointEst.perp();
02277 double tan4 = parLim[4] == 0 ? 0 : tan(parLim[4]);
02278 double tan5 = parLim[5] == 0 ? 0 : tan(parLim[5]);
02279 if ( (lZ - parLim[2]) > lR*tan4
02280 && (lZ - parLim[3]) < lR*tan5
02281 && lR > parLim[0] && lR < parLim[1]){
02282 if (debug_){
02283 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is inside the volume"<<std::endl;
02284 }
02285
02286 if (iFDest == -1){
02287 iFDest = iFace;
02288 } else {
02289 if (fabs(tanDistToFace[iFDest]) > fabs(tanDistToFace[iFace])){
02290 iFDest = iFace;
02291 }
02292 }
02293 } else {
02294 if (debug_){
02295 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is NOT inside the volume"<<std::endl;
02296 }
02297 }
02298 }
02299
02300 }
02301 if (iFDest != -1){
02302 result = SteppingHelixStateInfo::OK;
02303 dist = distToFace[iFDest];
02304 tanDist = tanDistToFace[iFDest];
02305 if (debug_){
02306 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got a point near closest boundary -- face "<<iFDest<<std::endl;
02307 }
02308 } else {
02309 if (debug_){
02310 LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Failed to find a dest point inside the volume"<<std::endl;
02311 }
02312 }
02313
02314 return result;
02315 }
02316
02317
02318 bool SteppingHelixPropagator::isYokeVolume(const MagVolume* vol) const {
02319 static const std::string metname = "SteppingHelixPropagator";
02320 if (vol == 0) return false;
02321
02322
02323
02324
02325
02326
02327
02328
02329
02330
02331
02332
02333
02334
02335
02336
02337
02338
02339
02340
02341
02342 bool result = vol->isIron();
02343 if (result){
02344 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Volume is magnetic, located at "<<vol->position()<<std::endl;
02345 } else {
02346 if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Volume is not magnetic, located at "<<vol->position()<<std::endl;
02347 }
02348
02349 return result;
02350 }