CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SteppingHelixPropagator.cc
Go to the documentation of this file.
1 
13 //
14 // Original Author: Vyacheslav Krutelyov
15 // Created: Fri Mar 3 16:01:24 CST 2006
16 // $Id: SteppingHelixPropagator.cc,v 1.77 2012/01/18 21:38:09 slava77 Exp $
17 //
18 //
19 
20 
25 
29 
31 
34 
36 
38 
39 #include <sstream>
40 #include <typeinfo>
41 
44 {
45  field_ = 0;
46 }
47 
50  Propagator(dir),
51  unit55_(AlgebraicMatrixID())
52 {
53  field_ = field;
54  vbField_ = dynamic_cast<const VolumeBasedMagneticField*>(field_);
57  debug_ = false;
58  noMaterialMode_ = false;
59  noErrorPropagation_ = false;
60  applyRadX0Correction_ = true;
61  useMagVolumes_ = true;
62  useIsYokeFlag_ = true;
63  useMatVolumes_ = true;
64  useInTeslaFromMagField_ = false; //overrides behavior only if true
65  returnTangentPlane_ = true;
66  sendLogWarning_ = false;
67  useTuningForL2Speed_ = false;
68  for (int i = 0; i <= MAX_POINTS; i++){
69  svBuf_[i].p3 = Vector(0,0,0);
70  svBuf_[i].r3 = Point(0,0,0);
71  svBuf_[i].bf = Vector(0,0,0);
72  svBuf_[i].bfGradLoc = Vector(0,0,0);
75  svBuf_[i].isComplete = true;
76  svBuf_[i].isValid_ = true;
78  }
79  defaultStep_ = 5.;
80 
81  ecShiftPos_ = 0;
82  ecShiftNeg_ = 0;
83 
84 }
85 
87 SteppingHelixPropagator::propagate(const FreeTrajectoryState& ftsStart, const Plane& pDest) const {
88  return propagateWithPath(ftsStart, pDest).first;
89 }
90 
93 {
94  return propagateWithPath(ftsStart, cDest).first;
95 }
96 
99 {
100  return propagateWithPath(ftsStart, pDest).first;
101 }
102 
105  const GlobalPoint& pDest1, const GlobalPoint& pDest2) const
106 {
107  return propagateWithPath(ftsStart, pDest1, pDest2).first;
108 }
109 
112  const reco::BeamSpot& beamSpot) const
113 {
114  return propagateWithPath(ftsStart, beamSpot).first;
115 }
116 
117 
118 std::pair<TrajectoryStateOnSurface, double>
120  const Plane& pDest) const {
121 
123 
124  const StateInfo& svCurrent = propagate(svBuf_[0], pDest);
125 
126  return TsosPP(svCurrent.getStateOnSurface(pDest), svCurrent.path());
127 }
128 
129 std::pair<TrajectoryStateOnSurface, double>
131  const Cylinder& cDest) const {
132 
134 
135  const StateInfo& svCurrent = propagate(svBuf_[0], cDest);
136 
137  return TsosPP(svCurrent.getStateOnSurface(cDest, returnTangentPlane_), svCurrent.path());
138 }
139 
140 
141 std::pair<FreeTrajectoryState, double>
143  const GlobalPoint& pDest) const {
145 
146  const StateInfo& svCurrent = propagate(svBuf_[0], pDest);
147 
148  FreeTrajectoryState ftsDest;
149  svCurrent.getFreeState(ftsDest);
150 
151  return FtsPP(ftsDest, svCurrent.path());
152 }
153 
154 std::pair<FreeTrajectoryState, double>
156  const GlobalPoint& pDest1, const GlobalPoint& pDest2) const {
157 
158  if ((pDest1-pDest2).mag() < 1e-10){
159  if (sendLogWarning_){
160  edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: the points should be at a bigger distance"
161  <<std::endl;
162  }
163  return FtsPP();
164  }
166 
167  const StateInfo& svCurrent = propagate(svBuf_[0], pDest1, pDest2);
168 
169  FreeTrajectoryState ftsDest;
170  svCurrent.getFreeState(ftsDest);
171 
172  return FtsPP(ftsDest, svCurrent.path());
173 }
174 
175 
176 std::pair<FreeTrajectoryState, double>
178  const reco::BeamSpot& beamSpot) const {
179  GlobalPoint pDest1(beamSpot.x0(), beamSpot.y0(), beamSpot.z0());
180  GlobalPoint pDest2(pDest1.x() + beamSpot.dxdz()*1e3,
181  pDest1.y() + beamSpot.dydz()*1e3,
182  pDest1.z() + 1e3);
183  return propagateWithPath(ftsStart, pDest1, pDest2);
184 }
185 
188  const Surface& sDest) const {
189 
190  if (! sStart.isValid()){
191  if (sendLogWarning_){
192  edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: invalid input state"
193  <<std::endl;
194  }
195  return invalidState_;
196  }
197 
198  const Plane* pDest = dynamic_cast<const Plane*>(&sDest);
199  if (pDest != 0) return propagate(sStart, *pDest);
200 
201  const Cylinder* cDest = dynamic_cast<const Cylinder*>(&sDest);
202  if (cDest != 0) return propagate(sStart, *cDest);
203 
204  throw PropagationException("The surface is neither Cylinder nor Plane");
205 
206 }
207 
210  const Plane& pDest) const {
211 
212  if (! sStart.isValid()){
213  if (sendLogWarning_){
214  edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: invalid input state"
215  <<std::endl;
216  }
217  return invalidState_;
218  }
219  setIState(sStart);
220 
221  Point rPlane(pDest.position().x(), pDest.position().y(), pDest.position().z());
222  Vector nPlane(pDest.rotation().zx(), pDest.rotation().zy(), pDest.rotation().zz()); nPlane /= nPlane.mag();
223 
224  double pars[6] = { rPlane.x(), rPlane.y(), rPlane.z(),
225  nPlane.x(), nPlane.y(), nPlane.z() };
226 
227  propagate(PLANE_DT, pars);
228 
229  //(re)set it before leaving: dir =1 (-1) if path increased (decreased) and 0 if it didn't change
230  //need to implement this somewhere else as a separate function
231  double lDir = 0;
232  if (sStart.path() < svBuf_[cIndex_(nPoints_-1)].path()) lDir = 1.;
233  if (sStart.path() > svBuf_[cIndex_(nPoints_-1)].path()) lDir = -1.;
234  svBuf_[cIndex_(nPoints_-1)].dir = lDir;
235  return svBuf_[cIndex_(nPoints_-1)];
236 }
237 
240  const Cylinder& cDest) const {
241 
242  if (! sStart.isValid()){
243  if (sendLogWarning_){
244  edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: invalid input state"
245  <<std::endl;
246  }
247  return invalidState_;
248  }
249  setIState(sStart);
250 
251  double pars[6] = {0,0,0,0,0,0};
252  pars[RADIUS_P] = cDest.radius();
253 
254 
255  propagate(RADIUS_DT, pars);
256 
257  //(re)set it before leaving: dir =1 (-1) if path increased (decreased) and 0 if it didn't change
258  //need to implement this somewhere else as a separate function
259  double lDir = 0;
260  if (sStart.path() < svBuf_[cIndex_(nPoints_-1)].path()) lDir = 1.;
261  if (sStart.path() > svBuf_[cIndex_(nPoints_-1)].path()) lDir = -1.;
262  svBuf_[cIndex_(nPoints_-1)].dir = lDir;
263  return svBuf_[cIndex_(nPoints_-1)];
264 }
265 
268  const GlobalPoint& pDest) const {
269 
270  if (! sStart.isValid()){
271  if (sendLogWarning_){
272  edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: invalid input state"
273  <<std::endl;
274  }
275  return invalidState_;
276  }
277  setIState(sStart);
278 
279  double pars[6] = {pDest.x(), pDest.y(), pDest.z(), 0, 0, 0};
280 
281 
282  propagate(POINT_PCA_DT, pars);
283 
284  return svBuf_[cIndex_(nPoints_-1)];
285 }
286 
289  const GlobalPoint& pDest1, const GlobalPoint& pDest2) const {
290 
291  if ((pDest1-pDest2).mag() < 1e-10 || !sStart.isValid()){
292  if (sendLogWarning_){
293  if ((pDest1-pDest2).mag() < 1e-10)
294  edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: points are too close"
295  <<std::endl;
296  if (!sStart.isValid())
297  edm::LogWarning("SteppingHelixPropagator")<<"Can't propagate: invalid input state"
298  <<std::endl;
299  }
300  return invalidState_;
301  }
302  setIState(sStart);
303 
304  double pars[6] = {pDest1.x(), pDest1.y(), pDest1.z(),
305  pDest2.x(), pDest2.y(), pDest2.z()};
306 
307  propagate(LINE_PCA_DT, pars);
308 
309  return svBuf_[cIndex_(nPoints_-1)];
310 }
311 
313  nPoints_ = 0;
314  svBuf_[cIndex_(nPoints_)] = sStart; //do this anyways to have a fresh start
315  if (sStart.isComplete ) {
316  svBuf_[cIndex_(nPoints_)] = sStart;
317  nPoints_++;
318  } else {
319  loadState(svBuf_[cIndex_(nPoints_)], sStart.p3, sStart.r3, sStart.q,
320  propagationDirection(), sStart.covCurv);
321  nPoints_++;
322  }
324 }
325 
328  const double pars[6], double epsilon) const{
329 
330  static const std::string metname = "SteppingHelixPropagator";
331  StateInfo* svCurrent = &svBuf_[cIndex_(nPoints_-1)];
332 
333  //check if it's going to work at all
334  double tanDist = 0;
335  double dist = 0;
336  PropagationDirection refDirection = anyDirection;
337  Result result = refToDest(type, (*svCurrent), pars, dist, tanDist, refDirection);
338 
339  if (result != SteppingHelixStateInfo::OK || fabs(dist) > 1e6){
340  svCurrent->status_ = result;
341  if (fabs(dist) > 1e6) svCurrent->status_ = SteppingHelixStateInfo::INACC;
342  svCurrent->isValid_ = false;
343  svCurrent->field = field_;
344  if (sendLogWarning_){
345  edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Failed after first refToDest check with status "
347  <<std::endl;
348  }
349  return result;
350  }
351 
353  bool makeNextStep = true;
354  double dStep = defaultStep_;
355  PropagationDirection dir,oldDir;
356  dir = propagationDirection();
357  oldDir = dir;
358  int nOsc = 0;
359 
360  double distMag = 1e12;
361  double tanDistMag = 1e12;
362 
363  double distMat = 1e12;
364  double tanDistMat = 1e12;
365 
366  double tanDistNextCheck = -0.1;//just need a negative start val
367  double tanDistMagNextCheck = -0.1;
368  double tanDistMatNextCheck = -0.1;
369  double oldDStep = 0;
370  PropagationDirection oldRefDirection = propagationDirection();
371 
374 
375  bool isFirstStep = true;
376  bool expectNewMagVolume = false;
377 
378  int loopCount = 0;
379  while (makeNextStep){
380  dStep = defaultStep_;
381  svCurrent = &svBuf_[cIndex_(nPoints_-1)];
382  double curZ = svCurrent->r3.z();
383  double curR = svCurrent->r3.perp();
384  if ( fabs(curZ) < 440 && curR < 260) dStep = defaultStep_*2;
385 
386  //more such ifs might be scattered around
387  //even though dStep is large, it will still make stops at each volume boundary
389  dStep = 100.;
390  if (! useIsYokeFlag_ && fabs(curZ) < 667 && curR > 380 && curR < 850){
391  dStep = 5*(1+0.2*svCurrent->p3.mag());
392  }
393  }
394 
395  // refDirection = propagationDirection();
396 
397  tanDistNextCheck -= oldDStep;
398  tanDistMagNextCheck -= oldDStep;
399  tanDistMatNextCheck -= oldDStep;
400 
401  if (tanDistNextCheck < 0){
402  //use pre-computed values if it's the first step
403  if (! isFirstStep) refToDest(type, (*svCurrent), pars, dist, tanDist, refDirection);
404  // constrain allowed path for a tangential approach
405  if (fabs(tanDist) > 4.*(fabs(dist)) ) tanDist *= tanDist == 0 ? 0 :fabs(dist/tanDist*4.);
406 
407  tanDistNextCheck = fabs(tanDist)*0.5 - 0.5; //need a better guess (to-do)
408  //reasonable limit
409  if (tanDistNextCheck > defaultStep_*20. ) tanDistNextCheck = defaultStep_*20.;
410  oldRefDirection = refDirection;
411  } else {
412  tanDist = tanDist > 0. ? tanDist - oldDStep : tanDist + oldDStep;
413  refDirection = oldRefDirection;
414  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Skipped refToDest: guess tanDist = "<<tanDist
415  <<" next check at "<<tanDistNextCheck<<std::endl;
416  }
418  double fastSkipDist = fabs(dist) > fabs(tanDist) ? tanDist : dist;
419 
421  dir = refDirection;
422  } else {
423  dir = propagationDirection();
424  if (fabs(tanDist)<0.1 && refDirection != dir ){
425  //how did it get here? nOsc++;
426  dir = refDirection;
427  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;
428  }
429  }
430 
431  if (useMagVolumes_ && ! (fabs(dist) < fabs(epsilon))){//need to know the general direction
432  if (tanDistMagNextCheck < 0){
433  resultToMag = refToMagVolume((*svCurrent), dir, distMag, tanDistMag, fabs(fastSkipDist), expectNewMagVolume, fabs(tanDist));
434  // constrain allowed path for a tangential approach
435  if (fabs(tanDistMag) > 4.*(fabs(distMag)) ) tanDistMag *= tanDistMag == 0 ? 0 : fabs(distMag/tanDistMag*4.);
436 
437  tanDistMagNextCheck = fabs(tanDistMag)*0.5-0.5; //need a better guess (to-do)
438  //reasonable limit; "turn off" checking if bounds are further than the destination
439  if (tanDistMagNextCheck > defaultStep_*20.
440  || fabs(dist) < fabs(distMag)
441  || resultToMag ==SteppingHelixStateInfo::INACC)
442  tanDistMagNextCheck = defaultStep_*20 > fabs(fastSkipDist) ? fabs(fastSkipDist) : defaultStep_*20;
443  if (resultToMag != SteppingHelixStateInfo::INACC
444  && resultToMag != SteppingHelixStateInfo::OK) tanDistMagNextCheck = -1;
445  } else {
446  // resultToMag = SteppingHelixStateInfo::OK;
447  tanDistMag = tanDistMag > 0. ? tanDistMag - oldDStep : tanDistMag + oldDStep;
448  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Skipped refToMag: guess tanDistMag = "<<tanDistMag
449  <<" next check at "<<tanDistMagNextCheck;
450  }
451  }
452 
453  if (useMatVolumes_ && ! (fabs(dist) < fabs(epsilon))){//need to know the general direction
454  if (tanDistMatNextCheck < 0){
455  resultToMat = refToMatVolume((*svCurrent), dir, distMat, tanDistMat, fabs(fastSkipDist));
456  // constrain allowed path for a tangential approach
457  if (fabs(tanDistMat) > 4.*(fabs(distMat)) ) tanDistMat *= tanDistMat == 0 ? 0 : fabs(distMat/tanDistMat*4.);
458 
459  tanDistMatNextCheck = fabs(tanDistMat)*0.5-0.5; //need a better guess (to-do)
460  //reasonable limit; "turn off" checking if bounds are further than the destination
461  if (tanDistMatNextCheck > defaultStep_*20.
462  || fabs(dist) < fabs(distMat)
463  || resultToMat ==SteppingHelixStateInfo::INACC )
464  tanDistMatNextCheck = defaultStep_*20 > fabs(fastSkipDist) ? fabs(fastSkipDist) : defaultStep_*20;
465  if (resultToMat != SteppingHelixStateInfo::INACC
466  && resultToMat != SteppingHelixStateInfo::OK) tanDistMatNextCheck = -1;
467  } else {
468  // resultToMat = SteppingHelixStateInfo::OK;
469  tanDistMat = tanDistMat > 0. ? tanDistMat - oldDStep : tanDistMat + oldDStep;
470  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Skipped refToMat: guess tanDistMat = "<<tanDistMat
471  <<" next check at "<<tanDistMatNextCheck;
472  }
473  }
474 
475  double rDotP = svCurrent->r3.dot(svCurrent->p3);
476  if ((fabs(curZ) > 1.5e3 || curR >800.)
477  && ((dir == alongMomentum && rDotP > 0)
478  || (dir == oppositeToMomentum && rDotP < 0) )
479  ){
480  dStep = fabs(tanDist) -1e-12;
481  }
482  double tanDistMin = fabs(tanDist);
483  if (tanDistMin > fabs(tanDistMag)+0.05 &&
484  (resultToMag == SteppingHelixStateInfo::OK || resultToMag == SteppingHelixStateInfo::WRONG_VOLUME)){
485  tanDistMin = fabs(tanDistMag)+0.05; //try to step into the next volume
486  expectNewMagVolume = true;
487  } else expectNewMagVolume = false;
488 
489  if (tanDistMin > fabs(tanDistMat)+0.05 && resultToMat == SteppingHelixStateInfo::OK){
490  tanDistMin = fabs(tanDistMat)+0.05; //try to step into the next volume
491  if (expectNewMagVolume) expectNewMagVolume = false;
492  }
493 
494  if (tanDistMin*fabs(svCurrent->dEdx) > 0.5*svCurrent->p3.mag()){
495  tanDistMin = 0.5*svCurrent->p3.mag()/fabs(svCurrent->dEdx);
496  if (expectNewMagVolume) expectNewMagVolume = false;
497  }
498 
499 
500 
501  double tanDistMinLazy = fabs(tanDistMin);
502  if ((type == POINT_PCA_DT || type == LINE_PCA_DT)
503  && fabs(tanDist) < 2.*fabs(tanDistMin) ){
504  //being lazy here; the best is to take into account the curvature
505  tanDistMinLazy = fabs(tanDistMin)*0.5;
506  }
507 
508  if (fabs(tanDistMinLazy) < dStep){
509  dStep = fabs(tanDistMinLazy);
510  }
511 
512  //keep this path length for the next step
513  oldDStep = dStep;
514 
515  if (dStep > 1e-10 && ! (fabs(dist) < fabs(epsilon))){
516  StateInfo* svNext = &svBuf_[cIndex_(nPoints_)];
517  makeAtomStep((*svCurrent), (*svNext), dStep, dir, HEL_AS_F);
518 // if (useMatVolumes_ && expectNewMagVolume
519 // && svCurrent->magVol == svNext->magVol){
520 // double tmpDist=0;
521 // double tmpDistMag = 0;
522 // if (refToMagVolume((*svNext), dir, tmpDist, tmpDistMag, fabs(dist)) != SteppingHelixStateInfo::OK){
523 // //the point appears to be outside, but findVolume claims the opposite
524 // dStep += 0.05*fabs(tanDistMag/distMag); oldDStep = dStep; //do it again with a bigger step
525 // if (debug_) LogTrace(metname)
526 // <<"Failed to get into new mag volume: will try with new bigger step "<<dStep<<std::endl;
527 // makeAtomStep((*svCurrent), (*svNext), dStep, dir, HEL_AS_F);
528 // }
529 // }
530  nPoints_++; svCurrent = &svBuf_[cIndex_(nPoints_-1)];
531  if (oldDir != dir){
532  nOsc++;
533  tanDistNextCheck = -1;//check dist after osc
534  tanDistMagNextCheck = -1;
535  tanDistMatNextCheck = -1;
536  }
537  oldDir = dir;
538  }
539 
540  if (nOsc>1 && fabs(dStep)>epsilon){
541  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Ooops"<<std::endl;
542  }
543 
544  if (fabs(dist) < fabs(epsilon) ) result = SteppingHelixStateInfo::OK;
545 
546  if ((type == POINT_PCA_DT || type == LINE_PCA_DT )
547  && fabs(dStep) < fabs(epsilon) ){
548  //now check if it's not a branch point (peek ahead at 1 cm)
549  double nextDist = 0;
550  double nextTanDist = 0;
551  PropagationDirection nextRefDirection = anyDirection;
552  StateInfo* svNext = &svBuf_[cIndex_(nPoints_)];
553  makeAtomStep((*svCurrent), (*svNext), 1., dir, HEL_AS_F);
554  nPoints_++; svCurrent = &svBuf_[cIndex_(nPoints_-1)];
555  refToDest(type, (*svCurrent), pars, nextDist, nextTanDist, nextRefDirection);
556  if ( fabs(nextDist) > fabs(dist)){
557  nPoints_--; svCurrent = &svBuf_[cIndex_(nPoints_-1)];
559  if (debug_){
560  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Found real local minimum in PCA"<<std::endl;
561  }
562  } else {
563  //keep this trial point and continue
564  dStep = defaultStep_;
565  if (debug_){
566  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Found branch point in PCA"<<std::endl;
567  }
568  }
569  }
570 
571  if (nPoints_ > MAX_STEPS*1./defaultStep_ || loopCount > MAX_STEPS*100
572  || nOsc > 6 ) result = SteppingHelixStateInfo::FAULT;
573 
574  if (svCurrent->p3.mag() < 0.1 ) result = SteppingHelixStateInfo::RANGEOUT;
575 
576  curZ = svCurrent->r3.z();
577  curR = svCurrent->r3.perp();
578  if ( curR > 20000 || fabs(curZ) > 20000 ) result = SteppingHelixStateInfo::INACC;
579 
580  makeNextStep = result == SteppingHelixStateInfo::UNDEFINED;
581  svCurrent->status_ = result;
582  svCurrent->isValid_ = (result == SteppingHelixStateInfo::OK || makeNextStep );
583  svCurrent->field = field_;
584 
585  isFirstStep = false;
586  loopCount++;
587  }
588 
590  edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Propagation failed with status "
592  <<std::endl;
593  if (result == SteppingHelixStateInfo::RANGEOUT)
594  edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Momentum at last point is too low (<0.1) p_last = "
595  <<svCurrent->p3.mag()
596  <<std::endl;
597  if (result == SteppingHelixStateInfo::INACC)
598  edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Went too far: the last point is at "<<svCurrent->r3
599  <<std::endl;
600  if (result == SteppingHelixStateInfo::FAULT && nOsc > 6)
601  edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Infinite loop condidtion detected: going in cycles. Break after 6 cycles"
602  <<std::endl;
604  edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<" Tired to go farther. Made too many steps: more than "
606  <<std::endl;
607 
608  }
609 
610  if (debug_){
611  switch (type) {
612  case RADIUS_DT:
613  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to radius "<<pars[RADIUS_P]<<std::endl;
614  break;
615  case Z_DT:
616  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to z "<<pars[Z_P]<<std::endl;
617  break;
618  case PATHL_DT:
619  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to pathL "<<pars[PATHL_P]<<std::endl;
620  break;
621  case PLANE_DT:
622  {
623  Point rPlane(pars[0], pars[1], pars[2]);
624  Vector nPlane(pars[3], pars[4], pars[5]);
625  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to plane r0:"<<rPlane<<" n:"<<nPlane<<std::endl;
626  }
627  break;
628  case POINT_PCA_DT:
629  {
630  Point rDest(pars[0], pars[1], pars[2]);
631  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to PCA to point "<<rDest<<std::endl;
632  }
633  break;
634  case LINE_PCA_DT:
635  {
636  Point rDest1(pars[0], pars[1], pars[2]);
637  Point rDest2(pars[3], pars[4], pars[5]);
638  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to PCA to line "<<rDest1<<" - "<<rDest2<<std::endl;
639  }
640  break;
641  default:
642  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"going to NOT IMPLEMENTED"<<std::endl;
643  break;
644  }
645  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;
646  }
647 
648  return result;
649 }
650 
653  const SteppingHelixPropagator::Point& r3, int charge,
655  const AlgebraicSymMatrix55& covCurv) const{
656  static const std::string metname = "SteppingHelixPropagator";
657 
658  svCurrent.q = charge;
659  svCurrent.p3 = p3;
660  svCurrent.r3 = r3;
661  svCurrent.dir = dir == alongMomentum ? 1.: -1.;
662 
663  svCurrent.path_ = 0; // this could've held the initial path
664  svCurrent.radPath_ = 0;
665 
666  GlobalPoint gPointNegZ(svCurrent.r3.x(), svCurrent.r3.y(), -fabs(svCurrent.r3.z()));
667  GlobalPoint gPointNorZ(svCurrent.r3.x(), svCurrent.r3.y(), svCurrent.r3.z());
668 
669  float gpmag = gPointNegZ.mag2();
670  float pmag2 = p3.mag2();
671  if (gpmag > 1e20f ) {
672  LogTrace(metname)<<"Initial point is too far";
673  svCurrent.isValid_ = false;
674  return;
675  }
676  if (! (gpmag == gpmag) ) {
677  LogTrace(metname)<<"Initial point is a nan";
678  edm::LogWarning("SteppingHelixPropagatorNAN")<<"Initial point is a nan";
679  svCurrent.isValid_ = false;
680  return;
681  }
682  if (! (pmag2 == pmag2) ) {
683  LogTrace(metname)<<"Initial momentum is a nan";
684  edm::LogWarning("SteppingHelixPropagatorNAN")<<"Initial momentum is a nan";
685  svCurrent.isValid_ = false;
686  return;
687  }
688 
689  GlobalVector bf(0,0,0);
690  // = field_->inTesla(gPoint);
691  if (useMagVolumes_){
692  if (vbField_ ){
693  if (vbField_->isZSymmetric()){
694  svCurrent.magVol = vbField_->findVolume(gPointNegZ);
695  } else {
696  svCurrent.magVol = vbField_->findVolume(gPointNorZ);
697  }
698  if (useIsYokeFlag_){
699  double curRad = svCurrent.r3.perp();
700  if (curRad > 380 && curRad < 850 && fabs(svCurrent.r3.z()) < 667){
701  svCurrent.isYokeVol = isYokeVolume(svCurrent.magVol);
702  } else {
703  svCurrent.isYokeVol = false;
704  }
705  }
706  } else {
707  edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Failed to cast into VolumeBasedMagneticField: fall back to the default behavior"<<std::endl;
708  svCurrent.magVol = 0;
709  }
710  if (debug_){
711  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got volume at "<<svCurrent.magVol<<std::endl;
712  }
713  }
714 
715  if (useMagVolumes_ && svCurrent.magVol != 0 && ! useInTeslaFromMagField_){
716  if (vbField_->isZSymmetric()){
717  bf = svCurrent.magVol->inTesla(gPointNegZ);
718  if (debug_){
719  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Loaded bfield float: "<<bf
720  <<" at global float "<< gPointNegZ<<" double "<< svCurrent.r3<<std::endl;
721  LocalPoint lPoint(svCurrent.magVol->toLocal(gPointNegZ));
722  LocalVector lbf = svCurrent.magVol->fieldInTesla(lPoint);
723  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"\t cf in local locF: "<<lbf<<" at "<<lPoint<<std::endl;
724  }
725  } else {
726  bf = svCurrent.magVol->inTesla(gPointNorZ);
727  if (debug_){
728  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Loaded bfield float: "<<bf
729  <<" at global float "<< gPointNorZ<<" double "<< svCurrent.r3<<std::endl;
730  LocalPoint lPoint(svCurrent.magVol->toLocal(gPointNorZ));
731  LocalVector lbf = svCurrent.magVol->fieldInTesla(lPoint);
732  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"\t cf in local locF: "<<lbf<<" at "<<lPoint<<std::endl;
733  }
734  }
735  if (r3.z() > 0 && vbField_->isZSymmetric() ){
736  svCurrent.bf.set(-bf.x(), -bf.y(), bf.z());
737  } else {
738  svCurrent.bf.set(bf.x(), bf.y(), bf.z());
739  }
740  } else {
741  GlobalPoint gPoint(r3.x(), r3.y(), r3.z());
742  bf = field_->inTesla(gPoint);
743  svCurrent.bf.set(bf.x(), bf.y(), bf.z());
744  }
745  if (svCurrent.bf.mag2() < 1e-32) svCurrent.bf.set(0., 0., 1e-16);
746  if (debug_){
747  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific
748  <<"Loaded bfield double: "<<svCurrent.bf<<" from float: "<<bf
749  <<" at float "<< gPointNorZ<<" "<<gPointNegZ<<" double "<< svCurrent.r3<<std::endl;
750  }
751 
752 
753 
754  double dEdXPrime = 0;
755  double dEdx = 0;
756  double radX0 = 0;
757  MatBounds rzLims;
758  dEdx = getDeDx(svCurrent, dEdXPrime, radX0, rzLims);
759  svCurrent.dEdx = dEdx; svCurrent.dEdXPrime = dEdXPrime;
760  svCurrent.radX0 = radX0;
761  svCurrent.rzLims = rzLims;
762 
763  svCurrent.covCurv =covCurv;
764 
765  svCurrent.isComplete = true;
766  svCurrent.isValid_ = true;
767 
768  if (debug_){
769  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Loaded at path: "<<svCurrent.path()<<" radPath: "<<svCurrent.radPath()
770  <<" p3 "<<" pt: "<<svCurrent.p3.perp()<<" phi: "<<svCurrent.p3.phi()
771  <<" eta: "<<svCurrent.p3.eta()
772  <<" "<<svCurrent.p3
773  <<" r3: "<<svCurrent.r3
774  <<" bField: "<<svCurrent.bf.mag()
775  <<std::endl;
776  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Input Covariance in Curvilinear RF "<<covCurv<<std::endl;
777  }
778 }
779 
782  double dP, const SteppingHelixPropagator::Vector& tau,
783  const SteppingHelixPropagator::Vector& drVec, double dS, double dX0,
784  const AlgebraicMatrix55& dCovCurvTransform) const{
785  static const std::string metname = "SteppingHelixPropagator";
786  svNext.q = svPrevious.q;
787  svNext.dir = dS > 0.0 ? 1.: -1.;
788  svNext.p3 = tau; svNext.p3*=(svPrevious.p3.mag() - svNext.dir*fabs(dP));
789 
790  svNext.r3 = svPrevious.r3; svNext.r3 += drVec;
791 
792  svNext.path_ = svPrevious.path() + dS;
793  svNext.radPath_ = svPrevious.radPath() + dX0;
794 
795  GlobalPoint gPointNegZ(svNext.r3.x(), svNext.r3.y(), -fabs(svNext.r3.z()));
796  GlobalPoint gPointNorZ(svNext.r3.x(), svNext.r3.y(), svNext.r3.z());
797 
798  GlobalVector bf(0,0,0);
799 
800  if (useMagVolumes_){
801  if (vbField_ != 0){
802  if (vbField_->isZSymmetric()){
803  svNext.magVol = vbField_->findVolume(gPointNegZ);
804  } else {
805  svNext.magVol = vbField_->findVolume(gPointNorZ);
806  }
807  if (useIsYokeFlag_){
808  double curRad = svNext.r3.perp();
809  if (curRad > 380 && curRad < 850 && fabs(svNext.r3.z()) < 667){
810  svNext.isYokeVol = isYokeVolume(svNext.magVol);
811  } else {
812  svNext.isYokeVol = false;
813  }
814  }
815  } else {
816  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Failed to cast into VolumeBasedMagneticField"<<std::endl;
817  svNext.magVol = 0;
818  }
819  if (debug_){
820  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got volume at "<<svNext.magVol<<std::endl;
821  }
822  }
823 
824  if (useMagVolumes_ && svNext.magVol != 0 && ! useInTeslaFromMagField_){
825  if (vbField_->isZSymmetric()){
826  bf = svNext.magVol->inTesla(gPointNegZ);
827  } else {
828  bf = svNext.magVol->inTesla(gPointNorZ);
829  }
830  if (svNext.r3.z() > 0 && vbField_->isZSymmetric() ){
831  svNext.bf.set(-bf.x(), -bf.y(), bf.z());
832  } else {
833  svNext.bf.set(bf.x(), bf.y(), bf.z());
834  }
835  } else {
836  GlobalPoint gPoint(svNext.r3.x(), svNext.r3.y(), svNext.r3.z());
837  bf = field_->inTesla(gPoint);
838  svNext.bf.set(bf.x(), bf.y(), bf.z());
839  }
840  if (svNext.bf.mag2() < 1e-32) svNext.bf.set(0., 0., 1e-16);
841 
842 
843  double dEdXPrime = 0;
844  double dEdx = 0;
845  double radX0 = 0;
846  MatBounds rzLims;
847  dEdx = getDeDx(svNext, dEdXPrime, radX0, rzLims);
848  svNext.dEdx = dEdx; svNext.dEdXPrime = dEdXPrime;
849  svNext.radX0 = radX0;
850  svNext.rzLims = rzLims;
851 
852  //update Emat only if it's valid
853  svNext.hasErrorPropagated_ = svPrevious.hasErrorPropagated_;
854  if (svPrevious.hasErrorPropagated_){
855  {
856  AlgebraicMatrix55 tmp = dCovCurvTransform*svPrevious.covCurv;
857  ROOT::Math::AssignSym::Evaluate(svNext.covCurv, tmp*ROOT::Math::Transpose(dCovCurvTransform));
858 
859  svNext.covCurv += svPrevious.matDCovCurv;
860  }
861  } else {
862  //could skip dragging along the unprop. cov later .. now
863  // svNext.cov = svPrevious.cov;
864  }
865 
866  if (debug_){
867  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Now at path: "<<svNext.path()<<" radPath: "<<svNext.radPath()
868  <<" p3 "<<" pt: "<<svNext.p3.perp()<<" phi: "<<svNext.p3.phi()
869  <<" eta: "<<svNext.p3.eta()
870  <<" "<<svNext.p3
871  <<" r3: "<<svNext.r3
872  <<" dPhi: "<<acos(svNext.p3.unit().dot(svPrevious.p3.unit()))
873  <<" bField: "<<svNext.bf.mag()
874  <<" dedx: "<<svNext.dEdx
875  <<std::endl;
876  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"New Covariance "<<svNext.covCurv<<std::endl;
877  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Transf by dCovTransform "<<dCovCurvTransform<<std::endl;
878  }
879 }
880 
882  const SteppingHelixPropagator::Vector& tau) const{
883  Vector zRep(0., 0., 1.);
884  rep.lX = tau;
885  rep.lY = zRep.cross(tau); rep.lY *= 1./tau.perp();
886  rep.lZ = rep.lX.cross(rep.lY);
887 }
888 
891  double dS,
893  SteppingHelixPropagator::Fancy fancy) const{
894  static const std::string metname = "SteppingHelixPropagator";
895  if (debug_){
896  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Make atom step "<<svCurrent.path()<<" with step "<<dS<<" in direction "<<dir<<std::endl;
897  }
898 
899  double dP = 0;
900  double curP = svCurrent.p3.mag();
901  Vector tau = svCurrent.p3; tau *= 1./curP;
902  Vector tauNext(tau);
903  Vector drVec(0,0,0);
904 
905  dS = dir == alongMomentum ? fabs(dS) : -fabs(dS);
906 
907 
908  double radX0 = 1e24;
909 
910  switch (fancy){
911  case HEL_AS_F:
912  case HEL_ALL_F:{
913  double p0 = curP;// see above = svCurrent.p3.mag();
914  double p0Inv = 1./p0;
915  double b0 = svCurrent.bf.mag();
916 
917  //get to the mid-point first
918  double phi = (2.99792458e-3*svCurrent.q*b0*p0Inv)*dS/2.;
919  bool phiSmall = fabs(phi) < 1e-4;
920 
921  double cosPhi = 0;
922  double sinPhi = 0;
923 
924  double oneLessCosPhi=0;
925  double oneLessCosPhiOPhi=0;
926  double sinPhiOPhi=0;
927  double phiLessSinPhiOPhi=0;
928 
929  if (phiSmall){
930  double phi2 = phi*phi;
931  double phi3 = phi2*phi;
932  double phi4 = phi3*phi;
933  sinPhi = phi - phi3/6. + phi4*phi/120.;
934  cosPhi = 1. -phi2/2. + phi4/24.;
935  oneLessCosPhi = phi2/2. - phi4/24. + phi2*phi4/720.; // 0.5*phi*phi;//*(1.- phi*phi/12.);
936  oneLessCosPhiOPhi = 0.5*phi - phi3/24. + phi2*phi3/720.;//*(1.- phi*phi/12.);
937  sinPhiOPhi = 1. - phi*phi/6. + phi4/120.;
938  phiLessSinPhiOPhi = phi*phi/6. - phi4/120. + phi4*phi2/5040.;//*(1. - phi*phi/20.);
939  } else {
940  cosPhi = cos(phi);
941  sinPhi = sin(phi);
942  oneLessCosPhi = 1.-cosPhi;
943  oneLessCosPhiOPhi = oneLessCosPhi/phi;
944  sinPhiOPhi = sinPhi/phi;
945  phiLessSinPhiOPhi = 1 - sinPhiOPhi;
946  }
947 
948  Vector bHat = svCurrent.bf; bHat *= 1./b0; //bHat.mag();
949  // bool bAlongZ = fabs(bHat.z()) > 0.9999;
950 
951  Vector btVec(bHat.cross(tau)); // for balong z btVec.z()==0
952  double tauB = tau.dot(bHat);
953  Vector bbtVec(bHat*tauB - tau); // (-tau.x(), -tau.y(), 0)
954 
955  //don't need it here tauNext = tau + bbtVec*oneLessCosPhi - btVec*sinPhi;
956  drVec = bbtVec*phiLessSinPhiOPhi; drVec -= btVec*oneLessCosPhiOPhi; drVec += tau;
957  drVec *= dS/2.;
958 
959  double dEdx = svCurrent.dEdx;
960  double dEdXPrime = svCurrent.dEdXPrime;
961  radX0 = svCurrent.radX0;
962  dP = dEdx*dS;
963 
964  //improve with above values:
965  drVec += svCurrent.r3;
966  GlobalVector bfGV(0,0,0);
967  Vector bf(0,0,0); //(bfGV.x(), bfGV.y(), bfGV.z());
968  // = svCurrent.magVol->inTesla(GlobalPoint(drVec.x(), drVec.y(), -fabs(drVec.z())));
969  if (useMagVolumes_ && svCurrent.magVol != 0 && ! useInTeslaFromMagField_){
970  // this negative-z business will break at some point
971  if (vbField_->isZSymmetric()){
972  bfGV = svCurrent.magVol->inTesla(GlobalPoint(drVec.x(), drVec.y(), -fabs(drVec.z())));
973  } else {
974  bfGV = svCurrent.magVol->inTesla(GlobalPoint(drVec.x(), drVec.y(), drVec.z()));
975  }
976  if (drVec.z() > 0 && vbField_->isZSymmetric()){
977  bf.set(-bfGV.x(), -bfGV.y(), bfGV.z());
978  } else {
979  bf.set(bfGV.x(), bfGV.y(), bfGV.z());
980  }
981  } else {
982  bfGV = field_->inTesla(GlobalPoint(drVec.x(), drVec.y(), drVec.z()));
983  bf.set(bfGV.x(), bfGV.y(), bfGV.z());
984  }
985  double b0Init = b0;
986  b0 = bf.mag();
987  if (b0 < 1e-16) {
988  b0 = 1e-16;
989  bf.set(0., 0., 1e-16);
990  }
991  if (debug_){
992  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Improved b "<<b0
993  <<" at r3 "<<drVec<<std::endl;
994  }
995 
996  if (fabs((b0-b0Init)*dS) > 1){
997  //missed the mag volume boundary?
998  if (debug_){
999  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Large bf*dS change "<<fabs((b0-svCurrent.bf.mag())*dS)
1000  <<" --> recalc dedx"<<std::endl;
1001  }
1002  svNext.r3 = drVec;
1003  svNext.bf = bf;
1004  svNext.p3 = svCurrent.p3;
1005  svNext.isYokeVol = svCurrent.isYokeVol;
1006  MatBounds rzTmp;
1007  dEdx = getDeDx(svNext, dEdXPrime, radX0, rzTmp);
1008  dP = dEdx*dS;
1009  if (debug_){
1010  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"New dEdX= "<<dEdx
1011  <<" dP= "<<dP
1012  <<" for p0 "<<p0<<std::endl;
1013  }
1014  }
1015  //p0 is mid-way and b0 from mid-point
1016  p0 += dP/2.; if (p0 < 1e-2) p0 = 1e-2;
1017  p0Inv = 1./p0;
1018 
1019  phi = (2.99792458e-3*svCurrent.q*b0*p0Inv)*dS;
1020  phiSmall = fabs(phi) < 1e-4;
1021 
1022  if (phiSmall){
1023  double phi2 = phi*phi;
1024  double phi3 = phi2*phi;
1025  double phi4 = phi3*phi;
1026  sinPhi = phi - phi3/6. + phi4*phi/120.;
1027  cosPhi = 1. -phi2/2. + phi4/24.;
1028  oneLessCosPhi = phi2/2. - phi4/24. + phi2*phi4/720.; // 0.5*phi*phi;//*(1.- phi*phi/12.);
1029  oneLessCosPhiOPhi = 0.5*phi - phi3/24. + phi2*phi3/720.;//*(1.- phi*phi/12.);
1030  sinPhiOPhi = 1. - phi*phi/6. + phi4/120.;
1031  phiLessSinPhiOPhi = phi*phi/6. - phi4/120. + phi4*phi2/5040.;//*(1. - phi*phi/20.);
1032  }else {
1033  cosPhi = cos(phi);
1034  sinPhi = sin(phi);
1035  oneLessCosPhi = 1.-cosPhi;
1036  oneLessCosPhiOPhi = oneLessCosPhi/phi;
1037  sinPhiOPhi = sinPhi/phi;
1038  phiLessSinPhiOPhi = 1. - sinPhiOPhi;
1039  }
1040 
1041  bHat = bf; bHat *= 1./b0;// as above =1./bHat.mag();
1042  // bAlongZ = fabs(bHat.z()) > 0.9999;
1043  btVec = bHat.cross(tau); // for b||z (-tau.y(), tau.x() ,0)
1044  tauB = tau.dot(bHat);
1045  bbtVec = bHat*tauB - tau; //bHat.cross(btVec); for b||z: (-tau.x(), -tau.y(), 0)
1046 
1047  tauNext = bbtVec*oneLessCosPhi; tauNext -= btVec*sinPhi; tauNext += tau; //for b||z tauNext.z() == tau.z()
1048  double tauNextPerpInv = 1./tauNext.perp();
1049  drVec = bbtVec*phiLessSinPhiOPhi; drVec -= btVec*oneLessCosPhiOPhi; drVec += tau;
1050  drVec *= dS;
1051 
1052 
1053  if (svCurrent.hasErrorPropagated_){
1054  double theta02 = 0;
1055  double dX0 = fabs(dS)/radX0;
1056 
1057  if (applyRadX0Correction_){
1058  // this provides the integrand for theta^2
1059  // if summed up along the path, should result in
1060  // theta_total^2 = Int_0^x0{ f(x)dX} = (13.6/p0)^2*x0*(1+0.036*ln(x0+1))
1061  // x0+1 above is to make the result infrared safe.
1062  double x0 = fabs(svCurrent.radPath());
1063  double alphaX0 = 13.6e-3*p0Inv; alphaX0 *= alphaX0;
1064  double betaX0 = 0.038;
1065  double logx0p1 = log(x0+1);
1066  theta02 = dX0*alphaX0*(1+betaX0*logx0p1)*(1 + betaX0*logx0p1 + 2.*betaX0*x0/(x0+1) );
1067  } else {
1068  theta02 = 196e-6* p0Inv * p0Inv * dX0; //14.e-3/p0*sqrt(fabs(dS)/radX0); // .. drop log term (this is non-additive)
1069  }
1070 
1071  double epsilonP0 = 1.+ dP/(p0-0.5*dP);
1072  // double omegaP0 = -dP/(p0-0.5*dP) + dS*dEdXPrime;
1073  // double dsp = dS/(p0-0.5*dP); //use the initial p0 (not the mid-point) to keep the transport properly additive
1074 
1075  Vector tbtVec(bHat - tauB*tau); // for b||z tau.z()*(-tau.x(), -tau.y(), 1.-tau.z())
1076 
1078  {
1079  //Slightly modified copy of the curvilinear jacobian (don't use the original just because it's in float precision
1080  // and seems to have some assumptions about the field values
1081  // notation changes: p1--> tau, p2-->tauNext
1082  // theta --> phi
1083  // Vector p1 = tau;
1084  // Vector p2 = tauNext;
1085  Point xStart = svCurrent.r3;
1086  Vector dx = drVec;
1087  //GlobalVector h = MagneticField::inInverseGeV(xStart);
1088  // Martijn: field is now given as parameter.. GlobalVector h = globalParameters.magneticFieldInInverseGeV(xStart);
1089 
1090  //double qbp = fts.signedInverseMomentum();
1091  double qbp = svCurrent.q*p0Inv;
1092  // double absS = dS;
1093 
1094  // calculate transport matrix
1095  // Origin: TRPRFN
1096  double t11 = tau.x(); double t12 = tau.y(); double t13 = tau.z();
1097  double t21 = tauNext.x(); double t22 = tauNext.y(); double t23 = tauNext.z();
1098  double cosl0 = tau.perp();
1099  // double cosl1 = 1./tauNext.perp(); //not quite a cos .. it's a cosec--> change to cosecl1 below
1100  double cosecl1 = tauNextPerpInv;
1101  //AlgebraicMatrix a(5,5,1);
1102  // define average magnetic field and gradient
1103  // at initial point - inlike TRPRFN
1104  Vector hn = bHat;
1105  // double qp = -2.99792458e-3*b0;
1106  // double q = -h.mag()*qbp;
1107 
1108  double q = -phi/dS; //qp*qbp; // -phi/dS
1109  // double theta = -phi;
1110  double sint = -sinPhi; double cost = cosPhi;
1111  double hn1 = hn.x(); double hn2 = hn.y(); double hn3 = hn.z();
1112  double dx1 = dx.x(); double dx2 = dx.y(); double dx3 = dx.z();
1113  // double hnDt1 = hn1*t11 + hn2*t12 + hn3*t13;
1114 
1115  double gamma = hn1*t21 + hn2*t22 + hn3*t23;
1116  double an1 = hn2*t23 - hn3*t22;
1117  double an2 = hn3*t21 - hn1*t23;
1118  double an3 = hn1*t22 - hn2*t21;
1119  // double auInv = sqrt(1.- t13*t13); double au = auInv>0 ? 1./auInv : 1e24;
1120  double auInv = cosl0; double au = auInv>0 ? 1./auInv : 1e24;
1121  // double auInv = sqrt(t11*t11 + t12*t12); double au = auInv>0 ? 1./auInv : 1e24;
1122  double u11 = -au*t12; double u12 = au*t11;
1123  double v11 = -t13*u12; double v12 = t13*u11; double v13 = auInv;//t11*u12 - t12*u11;
1124  auInv = sqrt(1. - t23*t23); au = auInv>0 ? 1./auInv : 1e24;
1125  // auInv = sqrt(t21*t21 + t22*t22); au = auInv>0 ? 1./auInv : 1e24;
1126  double u21 = -au*t22; double u22 = au*t21;
1127  double v21 = -t23*u22; double v22 = t23*u21; double v23 = auInv;//t21*u22 - t22*u21;
1128  // now prepare the transport matrix
1129  // pp only needed in high-p case (WA)
1130  // double pp = 1./qbp;
1132  // moved up (where -h.mag() is needed()
1133  // double qp = q*pp;
1134  double anv = -(hn1*u21 + hn2*u22 );
1135  double anu = (hn1*v21 + hn2*v22 + hn3*v23);
1136  double omcost = oneLessCosPhi; double tmsint = -phi*phiLessSinPhiOPhi;
1137 
1138  double hu1 = - hn3*u12;
1139  double hu2 = hn3*u11;
1140  double hu3 = hn1*u12 - hn2*u11;
1141 
1142  double hv1 = hn2*v13 - hn3*v12;
1143  double hv2 = hn3*v11 - hn1*v13;
1144  double hv3 = hn1*v12 - hn2*v11;
1145 
1146  // 1/p - doesn't change since |tau| = |tauNext| ... not. It changes now
1147  dCCurvTransform_(0,0) = 1./(epsilonP0*epsilonP0)*(1. + dS*dEdXPrime);
1148 
1149  // lambda
1150 
1151  dCCurvTransform_(1,0) = phi*p0/svCurrent.q*cosecl1*
1152  (sinPhi*bbtVec.z() - cosPhi*btVec.z());
1153  //was dCCurvTransform_(1,0) = -qp*anv*(t21*dx1 + t22*dx2 + t23*dx3); //NOTE (SK) this was found to have an opposite sign
1154  //from independent re-calculation ... in fact the tauNext.dot.dR piece isnt reproduced
1155 
1156  dCCurvTransform_(1,1) = cost*(v11*v21 + v12*v22 + v13*v23) +
1157  sint*(hv1*v21 + hv2*v22 + hv3*v23) +
1158  omcost*(hn1*v11 + hn2*v12 + hn3*v13) * (hn1*v21 + hn2*v22 + hn3*v23) +
1159  anv*(-sint*(v11*t21 + v12*t22 + v13*t23) +
1160  omcost*(v11*an1 + v12*an2 + v13*an3) -
1161  tmsint*gamma*(hn1*v11 + hn2*v12 + hn3*v13) );
1162 
1163  dCCurvTransform_(1,2) = cost*(u11*v21 + u12*v22 ) +
1164  sint*(hu1*v21 + hu2*v22 + hu3*v23) +
1165  omcost*(hn1*u11 + hn2*u12 ) * (hn1*v21 + hn2*v22 + hn3*v23) +
1166  anv*(-sint*(u11*t21 + u12*t22 ) +
1167  omcost*(u11*an1 + u12*an2 ) -
1168  tmsint*gamma*(hn1*u11 + hn2*u12 ) );
1169  dCCurvTransform_(1,2) *= cosl0;
1170 
1171  // Commented out in part for reproducibility purposes: these terms are zero in cart->curv
1172  // dCCurvTransform_(1,3) = -q*anv*(u11*t21 + u12*t22 ); //don't show up in cartesian setup-->curv
1173  //why would lambdaNext depend explicitely on initial position ? any arbitrary init point can be chosen not
1174  // affecting the final state's momentum direction ... is this the field gradient in curvilinear coord?
1175  // dCCurvTransform_(1,4) = -q*anv*(v11*t21 + v12*t22 + v13*t23); //don't show up in cartesian setup-->curv
1176 
1177  // phi
1178 
1179  dCCurvTransform_(2,0) = - phi*p0/svCurrent.q*cosecl1*cosecl1*
1180  (oneLessCosPhi*bHat.z()*btVec.mag2() + sinPhi*btVec.z() + cosPhi*tbtVec.z()) ;
1181  //was dCCurvTransform_(2,0) = -qp*anu*(t21*dx1 + t22*dx2 + t23*dx3)*cosecl1;
1182 
1183  dCCurvTransform_(2,1) = cost*(v11*u21 + v12*u22 ) +
1184  sint*(hv1*u21 + hv2*u22 ) +
1185  omcost*(hn1*v11 + hn2*v12 + hn3*v13) *
1186  (hn1*u21 + hn2*u22 ) +
1187  anu*(-sint*(v11*t21 + v12*t22 + v13*t23) +
1188  omcost*(v11*an1 + v12*an2 + v13*an3) -
1189  tmsint*gamma*(hn1*v11 + hn2*v12 + hn3*v13) );
1190  dCCurvTransform_(2,1) *= cosecl1;
1191 
1192  dCCurvTransform_(2,2) = cost*(u11*u21 + u12*u22 ) +
1193  sint*(hu1*u21 + hu2*u22 ) +
1194  omcost*(hn1*u11 + hn2*u12 ) *
1195  (hn1*u21 + hn2*u22 ) +
1196  anu*(-sint*(u11*t21 + u12*t22 ) +
1197  omcost*(u11*an1 + u12*an2 ) -
1198  tmsint*gamma*(hn1*u11 + hn2*u12 ) );
1199  dCCurvTransform_(2,2) *= cosecl1*cosl0;
1200 
1201  // Commented out in part for reproducibility purposes: these terms are zero in cart->curv
1202  // dCCurvTransform_(2,3) = -q*anu*(u11*t21 + u12*t22 )*cosecl1;
1203  //why would lambdaNext depend explicitely on initial position ? any arbitrary init point can be chosen not
1204  // affecting the final state's momentum direction ... is this the field gradient in curvilinear coord?
1205  // dCCurvTransform_(2,4) = -q*anu*(v11*t21 + v12*t22 + v13*t23)*cosecl1;
1206 
1207  // yt
1208 
1209  double pp = 1./qbp;
1210  // (SK) these terms seem to consistently have a sign opp from private derivation
1211  dCCurvTransform_(3,0) = - pp*(u21*dx1 + u22*dx2 ); //NB: modified from the original: changed the sign
1212  dCCurvTransform_(4,0) = - pp*(v21*dx1 + v22*dx2 + v23*dx3);
1213 
1214 
1215  dCCurvTransform_(3,1) = (sint*(v11*u21 + v12*u22 ) +
1216  omcost*(hv1*u21 + hv2*u22 ) +
1217  tmsint*(hn1*u21 + hn2*u22 ) *
1218  (hn1*v11 + hn2*v12 + hn3*v13))/q;
1219 
1220  dCCurvTransform_(3,2) = (sint*(u11*u21 + u12*u22 ) +
1221  omcost*(hu1*u21 + hu2*u22 ) +
1222  tmsint*(hn1*u21 + hn2*u22 ) *
1223  (hn1*u11 + hn2*u12 ))*cosl0/q;
1224 
1225  dCCurvTransform_(3,3) = (u11*u21 + u12*u22 );
1226 
1227  dCCurvTransform_(3,4) = (v11*u21 + v12*u22 );
1228 
1229  // zt
1230 
1231  dCCurvTransform_(4,1) = (sint*(v11*v21 + v12*v22 + v13*v23) +
1232  omcost*(hv1*v21 + hv2*v22 + hv3*v23) +
1233  tmsint*(hn1*v21 + hn2*v22 + hn3*v23) *
1234  (hn1*v11 + hn2*v12 + hn3*v13))/q;
1235 
1236  dCCurvTransform_(4,2) = (sint*(u11*v21 + u12*v22 ) +
1237  omcost*(hu1*v21 + hu2*v22 + hu3*v23) +
1238  tmsint*(hn1*v21 + hn2*v22 + hn3*v23) *
1239  (hn1*u11 + hn2*u12 ))*cosl0/q;
1240 
1241  dCCurvTransform_(4,3) = (u11*v21 + u12*v22 );
1242 
1243  dCCurvTransform_(4,4) = (v11*v21 + v12*v22 + v13*v23);
1244  // end of TRPRFN
1245  }
1246 
1247  if (debug_){
1248  Basis rep; setRep(rep, tauNext);
1249  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"rep X: "<<rep.lX<<" "<<rep.lX.mag()
1250  <<"\t Y: "<<rep.lY<<" "<<rep.lY.mag()
1251  <<"\t Z: "<<rep.lZ<<" "<<rep.lZ.mag();
1252  }
1253  //mind the sign of dS and dP (dS*dP < 0 allways)
1254  //covariance should grow no matter which direction you propagate
1255  //==> take abs values.
1256  //reset not needed: fill all below svCurrent.matDCov *= 0.;
1257  double mulRR = theta02*dS*dS/3.;
1258  double mulRP = theta02*fabs(dS)*p0/2.;
1259  double mulPP = theta02*p0*p0;
1260  double losPP = dP*dP*1.6/fabs(dS)*(1.0 + p0*1e-3);
1261  //another guess .. makes sense for 1 cm steps 2./dS == 2 [cm] / dS [cm] at low pt
1262  //double it by 1TeV
1263  //not gaussian anyways
1264  // derived from the fact that sigma_p/eLoss ~ 0.08 after ~ 200 steps
1265 
1266  //curvilinear
1267  double sinThetaInv = tauNextPerpInv;
1268  double p0Mat = p0+ 0.5*dP; // FIXME change this to p0 after it's clear that there's agreement in everything else
1269  double p0Mat2 = p0Mat*p0Mat;
1270  // with 6x6 formulation
1271  svCurrent.matDCovCurv*=0;
1272 
1273  svCurrent.matDCovCurv(0,0) = losPP/(p0Mat2*p0Mat2);
1274  // svCurrent.matDCovCurv(0,1) = 0;
1275  // svCurrent.matDCovCurv(0,2) = 0;
1276  // svCurrent.matDCovCurv(0,3) = 0;
1277  // svCurrent.matDCovCurv(0,4) = 0;
1278 
1279  // svCurrent.matDCovCurv(1,0) = 0;
1280  svCurrent.matDCovCurv(1,1) = mulPP/p0Mat2;
1281  // svCurrent.matDCovCurv(1,2) = 0;
1282  // svCurrent.matDCovCurv(1,3) = 0;
1283  svCurrent.matDCovCurv(1,4) = mulRP/p0Mat;
1284 
1285  // svCurrent.matDCovCurv(2,0) = 0;
1286  // svCurrent.matDCovCurv(2,1) = 0;
1287  svCurrent.matDCovCurv(2,2) = mulPP/p0Mat2*(sinThetaInv*sinThetaInv);
1288  svCurrent.matDCovCurv(2,3) = mulRP/p0Mat*sinThetaInv;
1289  // svCurrent.matDCovCurv(2,4) = 0;
1290 
1291  // svCurrent.matDCovCurv(3,0) = 0;
1292  // svCurrent.matDCovCurv(3,1) = 0;
1293  svCurrent.matDCovCurv(3,2) = mulRP/p0Mat*sinThetaInv;
1294  // svCurrent.matDCovCurv(3,0) = 0;
1295  svCurrent.matDCovCurv(3,3) = mulRR;
1296  // svCurrent.matDCovCurv(3,4) = 0;
1297 
1298  // svCurrent.matDCovCurv(4,0) = 0;
1299  svCurrent.matDCovCurv(4,1) = mulRP/p0Mat;
1300  // svCurrent.matDCovCurv(4,2) = 0;
1301  // svCurrent.matDCovCurv(4,3) = 0;
1302  svCurrent.matDCovCurv(4,4) = mulRR;
1303  }
1304  break;
1305  }
1306  default:
1307  break;
1308  }
1309 
1310  double pMag = curP;//svCurrent.p3.mag();
1311 
1312  if (dir == oppositeToMomentum) dP = fabs(dP);
1313  else if( dP < 0) { //the case of negative dP
1314  dP = -dP > pMag ? -pMag+1e-5 : dP;
1315  }
1316 
1317  getNextState(svCurrent, svNext, dP, tauNext, drVec, dS, dS/radX0,
1319  return true;
1320 }
1321 
1323  double& dEdXPrime, double& radX0,
1324  MatBounds& rzLims) const{
1325  radX0 = 1.e24;
1326  dEdXPrime = 0.;
1327  rzLims = MatBounds();
1328  if (noMaterialMode_) return 0;
1329 
1330  double dEdx = 0.;
1331 
1332  double lR = sv.r3.perp();
1333  double lZ = fabs(sv.r3.z());
1334 
1335  //assume "Iron" .. seems to be quite the same for brass/iron/PbW04
1336  //good for Fe within 3% for 0.2 GeV to 10PeV
1337 
1338  double dEdX_HCal = 0.95; //extracted from sim
1339  double dEdX_ECal = 0.45;
1340  double dEdX_coil = 0.35; //extracted from sim .. closer to 40% in fact
1341  double dEdX_Fe = 1;
1342  double dEdX_MCh = 0.053; //chambers on average
1343  double dEdX_Trk = 0.0114;
1344  double dEdX_Air = 2E-4;
1345  double dEdX_Vac = 0.0;
1346 
1347  double radX0_HCal = 1.44/0.8; //guessing
1348  double radX0_ECal = 0.89/0.7;
1349  double radX0_coil = 4.; //
1350  double radX0_Fe = 1.76;
1351  double radX0_MCh = 1e3; //
1352  double radX0_Trk = 320.;
1353  double radX0_Air = 3.e4;
1354  double radX0_Vac = 3.e9; //"big" number for vacuum
1355 
1356 
1357  //not all the boundaries are set below: this will be a default
1358  if (! (lR < 380 && lZ < 785)){
1359  if (lZ > 785 ) rzLims = MatBounds(0, 1e4, 785, 1e4);
1360  if (lZ < 785 ) rzLims = MatBounds(380, 1e4, 0, 785);
1361  }
1362 
1363  //this def makes sense assuming we don't jump into endcap volume from the other z-side in one step
1364  //also, it is a positive shift only (at least for now): can't move ec further into the detector
1365  double ecShift = sv.r3.z() > 0 ? fabs(ecShiftPos_) : fabs(ecShiftNeg_);
1366 
1367  //this should roughly figure out where things are
1368  //(numbers taken from Fig1.1.2 TDR and from geom xmls)
1369  if (lR < 2.9){ //inside beampipe
1370  dEdx = dEdX_Vac; radX0 = radX0_Vac;
1371  rzLims = MatBounds(0, 2.9, 0, 1E4);
1372  }
1373  else if (lR < 129){
1374  if (lZ < 294){
1375  rzLims = MatBounds(2.9, 129, 0, 294); //tracker boundaries
1376  dEdx = dEdX_Trk; radX0 = radX0_Trk;
1377  //somewhat empirical formula that ~ matches the average if going from 0,0,0
1378  //assuming "uniform" tracker material
1379  //doesn't really track material layer to layer
1380  double lEtaDet = fabs(sv.r3.eta());
1381  double scaleRadX = lEtaDet > 1.5 ? 0.7724 : 1./cosh(0.5*lEtaDet);//sin(2.*atan(exp(-0.5*lEtaDet)));
1382  scaleRadX *= scaleRadX;
1383  if (lEtaDet > 2 && lZ > 20) scaleRadX *= (lEtaDet-1.);
1384  if (lEtaDet > 2.5 && lZ > 20) scaleRadX *= (lEtaDet-1.);
1385  radX0 *= scaleRadX;
1386  }
1387  //endcap part begins here
1388  else if ( lZ < 294 + ecShift ){
1389  //gap in front of EE here, piece inside 2.9<R<129
1390  rzLims = MatBounds(2.9, 129, 294, 294 + ecShift);
1391  dEdx = dEdX_Air; radX0 = radX0_Air;
1392  }
1393  else if (lZ < 372 + ecShift){
1394  rzLims = MatBounds(2.9, 129, 294 + ecShift, 372 + ecShift); //EE here, piece inside 2.9<R<129
1395  dEdx = dEdX_ECal; radX0 = radX0_ECal;
1396  }//EE averaged out over a larger space
1397  else if (lZ < 398 + ecShift){
1398  rzLims = MatBounds(2.9, 129, 372 + ecShift, 398 + ecShift); //whatever goes behind EE 2.9<R<129 is air up to Z=398
1399  dEdx = dEdX_HCal*0.05; radX0 = radX0_Air;
1400  }//betw EE and HE
1401  else if (lZ < 555 + ecShift){
1402  rzLims = MatBounds(2.9, 129, 398 + ecShift, 555 + ecShift); //HE piece 2.9<R<129;
1403  dEdx = dEdX_HCal*0.96; radX0 = radX0_HCal/0.96;
1404  } //HE calor abit less dense
1405  else {
1406  // rzLims = MatBounds(2.9, 129, 555, 785);
1407  // set the boundaries first: they serve as stop-points here
1408  // the material is set below
1409  if (lZ < 568 + ecShift) rzLims = MatBounds(2.9, 129, 555 + ecShift, 568 + ecShift); //a piece of HE support R<129, 555<Z<568
1410  else if (lZ < 625 + ecShift){
1411  if (lR < 85 + ecShift) rzLims = MatBounds(2.9, 85, 568 + ecShift, 625 + ecShift);
1412  else rzLims = MatBounds(85, 129, 568 + ecShift, 625 + ecShift);
1413  } else if (lZ < 785 + ecShift) rzLims = MatBounds(2.9, 129, 625 + ecShift, 785 + ecShift);
1414  else if (lZ < 1500 + ecShift) rzLims = MatBounds(2.9, 129, 785 + ecShift, 1500 + ecShift);
1415  else rzLims = MatBounds(2.9, 129, 1500 + ecShift, 1E4);
1416 
1417  //iron .. don't care about no material in front of HF (too forward)
1418  if (! (lZ > 568 + ecShift && lZ < 625 + ecShift && lR > 85 ) // HE support
1419  && ! (lZ > 785 + ecShift && lZ < 850 + ecShift && lR > 118)) {dEdx = dEdX_Fe; radX0 = radX0_Fe; }
1420  else { dEdx = dEdX_MCh; radX0 = radX0_MCh; } //ME at eta > 2.2
1421  }
1422  }
1423  else if (lR < 287){
1424  if (lZ < 372 + ecShift && lR < 177){ // 129<<R<177
1425  if (lZ < 304) rzLims = MatBounds(129, 177, 0, 304); //EB 129<<R<177 0<Z<304
1426  else if (lZ < 343){ // 129<<R<177 304<Z<343
1427  if (lR < 135 ) rzLims = MatBounds(129, 135, 304, 343);// tk piece 129<<R<135 304<Z<343
1428  else if (lR < 172 ){ //
1429  if (lZ < 311 ) rzLims = MatBounds(135, 172, 304, 311);
1430  else rzLims = MatBounds(135, 172, 311, 343);
1431  } else {
1432  if (lZ < 328) rzLims = MatBounds(172, 177, 304, 328);
1433  else rzLims = MatBounds(172, 177, 328, 343);
1434  }
1435  }
1436  else if ( lZ < 343 + ecShift){
1437  rzLims = MatBounds(129, 177, 343, 343 + ecShift); //gap
1438  }
1439  else {
1440  if (lR < 156 ) rzLims = MatBounds(129, 156, 343 + ecShift, 372 + ecShift);
1441  else if ( (lZ - 343 - ecShift) > (lR - 156)*1.38 )
1442  rzLims = MatBounds(156, 177, 127.72 + ecShift, 372 + ecShift, atan(1.38), 0);
1443  else rzLims = MatBounds(156, 177, 343 + ecShift, 127.72 + ecShift, 0, atan(1.38));
1444  }
1445 
1446  if (!(lR > 135 && lZ <343 + ecShift && lZ > 304 )
1447  && ! (lR > 156 && lZ < 372 + ecShift && lZ > 343 + ecShift && ((lZ-343. - ecShift )< (lR-156.)*1.38)))
1448  {
1449  //the crystals are the same length, but they are not 100% of material
1450  double cosThetaEquiv = 0.8/sqrt(1.+lZ*lZ/lR/lR) + 0.2;
1451  if (lZ > 343) cosThetaEquiv = 1.;
1452  dEdx = dEdX_ECal*cosThetaEquiv; radX0 = radX0_ECal/cosThetaEquiv;
1453  } //EB
1454  else {
1455  if ( (lZ > 304 && lZ < 328 && lR < 177 && lR > 135)
1456  && ! (lZ > 311 && lR < 172) ) {dEdx = dEdX_Fe; radX0 = radX0_Fe; } //Tk_Support
1457  else if ( lZ > 343 && lZ < 343 + ecShift) { dEdx = dEdX_Air; radX0 = radX0_Air; }
1458  else {dEdx = dEdX_ECal*0.2; radX0 = radX0_Air;} //cables go here <-- will be abit too dense for ecShift > 0
1459  }
1460  }
1461  else if (lZ < 554 + ecShift){ // 129<R<177 372<Z<554 AND 177<R<287 0<Z<554
1462  if (lR < 177){ // 129<R<177 372<Z<554
1463  if ( lZ > 372 + ecShift && lZ < 398 + ecShift )rzLims = MatBounds(129, 177, 372 + ecShift, 398 + ecShift); // EE 129<R<177 372<Z<398
1464  else if (lZ < 548 + ecShift) rzLims = MatBounds(129, 177, 398 + ecShift, 548 + ecShift); // HE 129<R<177 398<Z<548
1465  else rzLims = MatBounds(129, 177, 548 + ecShift, 554 + ecShift); // HE gap 129<R<177 548<Z<554
1466  }
1467  else if (lR < 193){ // 177<R<193 0<Z<554
1468  if ((lZ - 307) < (lR - 177.)*1.739) rzLims = MatBounds(177, 193, 0, -0.803, 0, atan(1.739));
1469  else if ( lZ < 389) rzLims = MatBounds(177, 193, -0.803, 389, atan(1.739), 0.);
1470  else if ( lZ < 389 + ecShift) rzLims = MatBounds(177, 193, 389, 389 + ecShift); // air insert
1471  else if ( lZ < 548 + ecShift ) rzLims = MatBounds(177, 193, 389 + ecShift, 548 + ecShift);// HE 177<R<193 389<Z<548
1472  else rzLims = MatBounds(177, 193, 548 + ecShift, 554 + ecShift); // HE gap 177<R<193 548<Z<554
1473  }
1474  else if (lR < 264){ // 193<R<264 0<Z<554
1475  double anApex = 375.7278 - 193./1.327; // 230.28695599096
1476  if ( (lZ - 375.7278) < (lR - 193.)/1.327) rzLims = MatBounds(193, 264, 0, anApex, 0, atan(1./1.327)); //HB
1477  else if ( (lZ - 392.7278 ) < (lR - 193.)/1.327)
1478  rzLims = MatBounds(193, 264, anApex, anApex+17., atan(1./1.327), atan(1./1.327)); // HB-HE gap
1479  else if ( (lZ - 392.7278 - ecShift ) < (lR - 193.)/1.327)
1480  rzLims = MatBounds(193, 264, anApex+17., anApex+17. + ecShift, atan(1./1.327), atan(1./1.327)); // HB-HE gap air insert
1481  // HE (372,193)-(517,193)-(517,264)-(417.8,264)
1482  else if ( lZ < 517 + ecShift ) rzLims = MatBounds(193, 264, anApex+17. + ecShift, 517 + ecShift, atan(1./1.327), 0);
1483  else if (lZ < 548 + ecShift){ // HE+gap 193<R<264 517<Z<548
1484  if (lR < 246 ) rzLims = MatBounds(193, 246, 517 + ecShift, 548 + ecShift); // HE 193<R<246 517<Z<548
1485  else rzLims = MatBounds(246, 264, 517 + ecShift, 548 + ecShift); // HE gap 246<R<264 517<Z<548
1486  }
1487  else rzLims = MatBounds(193, 264, 548 + ecShift, 554 + ecShift); // HE gap 193<R<246 548<Z<554
1488  }
1489  else if ( lR < 275){ // 264<R<275 0<Z<554
1490  if (lZ < 433) rzLims = MatBounds(264, 275, 0, 433); //HB
1491  else if (lZ < 554 ) rzLims = MatBounds(264, 275, 433, 554); // HE gap 264<R<275 433<Z<554
1492  else rzLims = MatBounds(264, 275, 554, 554 + ecShift); // HE gap 264<R<275 554<Z<554 air insert
1493  }
1494  else { // 275<R<287 0<Z<554
1495  if (lZ < 402) rzLims = MatBounds(275, 287, 0, 402);// HB 275<R<287 0<Z<402
1496  else if (lZ < 554) rzLims = MatBounds(275, 287, 402, 554);// //HE gap 275<R<287 402<Z<554
1497  else rzLims = MatBounds(275, 287, 554, 554 + ecShift); //HE gap 275<R<287 554<Z<554 air insert
1498  }
1499 
1500  if ((lZ < 433 || lR < 264) && (lZ < 402 || lR < 275) && (lZ < 517 + ecShift || lR < 246) //notches
1501  //I should've made HE and HF different .. now need to shorten HE to match
1502  && lZ < 548 + ecShift
1503  && ! (lZ < 389 + ecShift && lZ > 335 && lR < 193 ) //not a gap EB-EE 129<R<193
1504  && ! (lZ > 307 && lZ < 335 && lR < 193 && ((lZ - 307) > (lR - 177.)*1.739)) //not a gap
1505  && ! (lR < 177 && lZ < 398 + ecShift) //under the HE nose
1506  && ! (lR < 264 && lR > 193 && fabs(441.5+0.5*ecShift - lZ + (lR - 269.)/1.327) < (8.5 + ecShift*0.5)) ) //not a gap
1507  { dEdx = dEdX_HCal; radX0 = radX0_HCal; //hcal
1508  }
1509  else if( (lR < 193 && lZ > 389 && lZ < 389 + ecShift ) || (lR > 264 && lR < 287 && lZ > 554 && lZ < 554 + ecShift)
1510  ||(lR < 264 && lR > 193 && fabs(441.5+8.5+0.5*ecShift - lZ + (lR - 269.)/1.327) < ecShift*0.5) ) {
1511  dEdx = dEdX_Air; radX0 = radX0_Air;
1512  }
1513  else {dEdx = dEdX_HCal*0.05; radX0 = radX0_Air; }//endcap gap
1514  }
1515  //the rest is a tube of endcap volume 129<R<251 554<Z<largeValue
1516  else if (lZ < 564 + ecShift){ // 129<R<287 554<Z<564
1517  if (lR < 251) {
1518  rzLims = MatBounds(129, 251, 554 + ecShift, 564 + ecShift); // HE support 129<R<251 554<Z<564
1519  dEdx = dEdX_Fe; radX0 = radX0_Fe;
1520  }//HE support
1521  else {
1522  rzLims = MatBounds(251, 287, 554 + ecShift, 564 + ecShift); //HE/ME gap 251<R<287 554<Z<564
1523  dEdx = dEdX_MCh; radX0 = radX0_MCh;
1524  }
1525  }
1526  else if (lZ < 625 + ecShift){ // ME/1/1 129<R<287 564<Z<625
1527  rzLims = MatBounds(129, 287, 564 + ecShift, 625 + ecShift);
1528  dEdx = dEdX_MCh; radX0 = radX0_MCh;
1529  }
1530  else if (lZ < 785 + ecShift){ //129<R<287 625<Z<785
1531  if (lR < 275) rzLims = MatBounds(129, 275, 625 + ecShift, 785 + ecShift); //YE/1 129<R<275 625<Z<785
1532  else { // 275<R<287 625<Z<785
1533  if (lZ < 720 + ecShift) rzLims = MatBounds(275, 287, 625 + ecShift, 720 + ecShift); // ME/1/2 275<R<287 625<Z<720
1534  else rzLims = MatBounds(275, 287, 720 + ecShift, 785 + ecShift); // YE/1 275<R<287 720<Z<785
1535  }
1536  if (! (lR > 275 && lZ < 720 + ecShift)) { dEdx = dEdX_Fe; radX0 = radX0_Fe; }//iron
1537  else { dEdx = dEdX_MCh; radX0 = radX0_MCh; }
1538  }
1539  else if (lZ < 850 + ecShift){
1540  rzLims = MatBounds(129, 287, 785 + ecShift, 850 + ecShift);
1541  dEdx = dEdX_MCh; radX0 = radX0_MCh;
1542  }
1543  else if (lZ < 910 + ecShift){
1544  rzLims = MatBounds(129, 287, 850 + ecShift, 910 + ecShift);
1545  dEdx = dEdX_Fe; radX0 = radX0_Fe;
1546  }//iron
1547  else if (lZ < 975 + ecShift){
1548  rzLims = MatBounds(129, 287, 910 + ecShift, 975 + ecShift);
1549  dEdx = dEdX_MCh; radX0 = radX0_MCh;
1550  }
1551  else if (lZ < 1000 + ecShift){
1552  rzLims = MatBounds(129, 287, 975 + ecShift, 1000 + ecShift);
1553  dEdx = dEdX_Fe; radX0 = radX0_Fe;
1554  }//iron
1555  else if (lZ < 1063 + ecShift){
1556  rzLims = MatBounds(129, 287, 1000 + ecShift, 1063 + ecShift);
1557  dEdx = dEdX_MCh; radX0 = radX0_MCh;
1558  }
1559  else if ( lZ < 1073 + ecShift){
1560  rzLims = MatBounds(129, 287, 1063 + ecShift, 1073 + ecShift);
1561  dEdx = dEdX_Fe; radX0 = radX0_Fe;
1562  }
1563  else if (lZ < 1E4 ) {
1564  rzLims = MatBounds(129, 287, 1073 + ecShift, 1E4);
1565  dEdx = dEdX_Air; radX0 = radX0_Air;
1566  }
1567  else {
1568 
1569  dEdx = dEdX_Air; radX0 = radX0_Air;
1570  }
1571  }
1572  else if (lR <380 && lZ < 667){
1573  if (lZ < 630){
1574  if (lR < 315) rzLims = MatBounds(287, 315, 0, 630);
1575  else if (lR < 341 ) rzLims = MatBounds(315, 341, 0, 630); //b-field ~linear rapid fall here
1576  else rzLims = MatBounds(341, 380, 0, 630);
1577  } else rzLims = MatBounds(287, 380, 630, 667);
1578 
1579  if (lZ < 630) { dEdx = dEdX_coil; radX0 = radX0_coil; }//a guess for the solenoid average
1580  else {dEdx = dEdX_Air; radX0 = radX0_Air; }//endcap gap
1581  }
1582  else {
1583  if (lZ < 667) {
1584  if (lR < 850){
1585  bool isIron = false;
1586  //sanity check in addition to flags
1587  if (useIsYokeFlag_ && useMagVolumes_ && sv.magVol != 0){
1588  isIron = sv.isYokeVol;
1589  } else {
1590  double bMag = sv.bf.mag();
1591  isIron = (bMag > 0.75 && ! (lZ > 500 && lR <500 && bMag < 1.15)
1592  && ! (lZ < 450 && lR > 420 && bMag < 1.15 ) );
1593  }
1594  //tell the material stepper where mat bounds are
1595  rzLims = MatBounds(380, 850, 0, 667);
1596  if (isIron) { dEdx = dEdX_Fe; radX0 = radX0_Fe; }//iron
1597  else { dEdx = dEdX_MCh; radX0 = radX0_MCh; }
1598  } else {
1599  rzLims = MatBounds(850, 1E4, 0, 667);
1600  dEdx = dEdX_Air; radX0 = radX0_Air;
1601  }
1602  }
1603  else if (lR > 750 ){
1604  rzLims = MatBounds(750, 1E4, 667, 1E4);
1605  dEdx = dEdX_Air; radX0 = radX0_Air;
1606  }
1607  else if (lZ < 667 + ecShift){
1608  rzLims = MatBounds(287, 750, 667, 667 + ecShift);
1609  dEdx = dEdX_Air; radX0 = radX0_Air;
1610  }
1611  //the rest is endcap piece with 287<R<750 Z>667
1612  else if (lZ < 724 + ecShift){
1613  if (lR < 380 ) rzLims = MatBounds(287, 380, 667 + ecShift, 724 + ecShift);
1614  else rzLims = MatBounds(380, 750, 667 + ecShift, 724 + ecShift);
1615  dEdx = dEdX_MCh; radX0 = radX0_MCh;
1616  }
1617  else if (lZ < 785 + ecShift){
1618  if (lR < 380 ) rzLims = MatBounds(287, 380, 724 + ecShift, 785 + ecShift);
1619  else rzLims = MatBounds(380, 750, 724 + ecShift, 785 + ecShift);
1620  dEdx = dEdX_Fe; radX0 = radX0_Fe;
1621  }//iron
1622  else if (lZ < 850 + ecShift){
1623  rzLims = MatBounds(287, 750, 785 + ecShift, 850 + ecShift);
1624  dEdx = dEdX_MCh; radX0 = radX0_MCh;
1625  }
1626  else if (lZ < 910 + ecShift){
1627  rzLims = MatBounds(287, 750, 850 + ecShift, 910 + ecShift);
1628  dEdx = dEdX_Fe; radX0 = radX0_Fe;
1629  }//iron
1630  else if (lZ < 975 + ecShift){
1631  rzLims = MatBounds(287, 750, 910 + ecShift, 975 + ecShift);
1632  dEdx = dEdX_MCh; radX0 = radX0_MCh;
1633  }
1634  else if (lZ < 1000 + ecShift){
1635  rzLims = MatBounds(287, 750, 975 + ecShift, 1000 + ecShift);
1636  dEdx = dEdX_Fe; radX0 = radX0_Fe;
1637  }//iron
1638  else if (lZ < 1063 + ecShift){
1639  if (lR < 360){
1640  rzLims = MatBounds(287, 360, 1000 + ecShift, 1063 + ecShift);
1641  dEdx = dEdX_MCh; radX0 = radX0_MCh;
1642  }
1643  //put empty air where me4/2 should be (future)
1644  else {
1645  rzLims = MatBounds(360, 750, 1000 + ecShift, 1063 + ecShift);
1646  dEdx = dEdX_Air; radX0 = radX0_Air;
1647  }
1648  }
1649  else if ( lZ < 1073 + ecShift){
1650  rzLims = MatBounds(287, 750, 1063 + ecShift, 1073 + ecShift);
1651  //this plate does not exist: air
1652  dEdx = dEdX_Air; radX0 = radX0_Air;
1653  }
1654  else if (lZ < 1E4 ) {
1655  rzLims = MatBounds(287, 750, 1073 + ecShift, 1E4);
1656  dEdx = dEdX_Air; radX0 = radX0_Air;
1657  }
1658  else {dEdx = dEdX_Air; radX0 = radX0_Air; }//air
1659  }
1660 
1661  //dEdx so far is a relative number (relative to iron)
1662  //scale by what's expected for iron (the function was fit from pdg table)
1663  //0.065 (PDG) --> 0.044 to better match with MPV
1664  double p0 = sv.p3.mag();
1665  double logp0 = log(p0);
1666  double p0powN33 = 0;
1667  if (p0>3.) {
1668  // p0powN33 = exp(-0.33*logp0); //calculate for p>3GeV
1669  double xx=1./p0; xx=sqrt(sqrt(sqrt(sqrt(xx)))); xx=xx*xx*xx*xx*xx; // this is (p0)**(-5/16), close enough to -0.33
1670  p0powN33 = xx;
1671  }
1672  double dEdX_mat = -(11.4 + 0.96*fabs(logp0+log(2.8)) + 0.033*p0*(1.0 - p0powN33) )*1e-3;
1673  //in GeV/cm .. 0.8 to get closer to the median or MPV
1674 
1675  dEdXPrime = dEdx == 0 ? 0 : -dEdx*(0.96/p0 + 0.033 - 0.022*p0powN33)*1e-3; //== d(dEdX)/dp
1676  dEdx *= dEdX_mat;
1677 
1678  return dEdx;
1679 }
1680 
1681 
1683  int result = ind%MAX_POINTS;
1684  if (ind != 0 && result == 0){
1685  result = MAX_POINTS;
1686  }
1687  return result;
1688 }
1689 
1693  const double pars[6],
1694  double& dist, double& tanDist,
1695  PropagationDirection& refDirection,
1696  double fastSkipDist) const{
1697  static const std::string metname = "SteppingHelixPropagator";
1699 
1700  switch (dest){
1701  case RADIUS_DT:
1702  {
1703  double curR = sv.r3.perp();
1704  dist = pars[RADIUS_P] - curR;
1705  if (fabs(dist) > fastSkipDist){
1707  break;
1708  }
1709  double curP2 = sv.p3.mag2();
1710  double curPtPos2 = sv.p3.perp2(); if(curPtPos2< 1e-16) curPtPos2=1e-16;
1711 
1712  double cosDPhiPR = (sv.r3.x()*sv.p3.x()+sv.r3.y()*sv.p3.y());//only the sign is needed cos((sv.r3.deltaPhi(sv.p3)));
1713  refDirection = dist*cosDPhiPR > 0 ?
1715  tanDist = dist*sqrt(curP2/curPtPos2);
1716  result = SteppingHelixStateInfo::OK;
1717  }
1718  break;
1719  case Z_DT:
1720  {
1721  double curZ = sv.r3.z();
1722  dist = pars[Z_P] - curZ;
1723  if (fabs(dist) > fastSkipDist){
1725  break;
1726  }
1727  double curP = sv.p3.mag();
1728  refDirection = sv.p3.z()*dist > 0. ?
1730  tanDist = dist/sv.p3.z()*curP;
1731  result = SteppingHelixStateInfo::OK;
1732  }
1733  break;
1734  case PLANE_DT:
1735  {
1736  Point rPlane(pars[0], pars[1], pars[2]);
1737  Vector nPlane(pars[3], pars[4], pars[5]);
1738 
1739 
1740  // unfortunately this doesn't work: the numbers are too large
1741  // bool pVertical = fabs(pars[5])>0.9999;
1742  // double dRDotN = pVertical? (sv.r3.z() - rPlane.z())*nPlane.z() :(sv.r3 - rPlane).dot(nPlane);
1743  double dRDotN = (sv.r3.x()-rPlane.x())*nPlane.x() + (sv.r3.y()-rPlane.y())*nPlane.y() + (sv.r3.z()-rPlane.z())*nPlane.z();//(sv.r3 - rPlane).dot(nPlane);
1744 
1745  dist = fabs(dRDotN);
1746  if (dist > fastSkipDist){
1748  break;
1749  }
1750  double curP = sv.p3.mag();
1751  double p0 = curP;
1752  double p0Inv = 1./p0;
1753  Vector tau(sv.p3); tau *=p0Inv;
1754  double tN = tau.dot(nPlane);
1755  refDirection = tN*dRDotN < 0. ?
1757  double b0 = sv.bf.mag();
1758  if (fabs(tN)>1e-24){
1759  tanDist = -dRDotN/tN;
1760  } else {
1761  tN = 1e-24;
1762  if (fabs(dRDotN)>1e-24) tanDist = 1e6;
1763  else tanDist = 1;
1764  }
1765  if (fabs(tanDist) > 1e4) tanDist = 1e4;
1766  if (b0>1.5e-6){
1767  double b0Inv = 1./b0;
1768  double tNInv = 1./tN;
1769  Vector bHat(sv.bf); bHat *=b0Inv;
1770  double bHatN = bHat.dot(nPlane);
1771  double cosPB = bHat.dot(tau);
1772  double kVal = 2.99792458e-3*sv.q*p0Inv*b0;
1773  double aVal = tanDist*kVal;
1774  Vector lVec = bHat.cross(tau);
1775  double bVal = lVec.dot(nPlane)*tNInv;
1776  if (fabs(aVal*bVal)< 0.3){
1777  double cVal = 1. - bHatN*cosPB*tNInv; // - sv.bf.cross(lVec).dot(nPlane)*b0Inv*tNInv; //1- bHat_n*bHat_tau/tau_n;
1778  double aacVal = cVal*aVal*aVal;
1779  if (fabs(aacVal)<1){
1780  double tanDCorr = bVal/2. + (bVal*bVal/2. + cVal/6)*aVal;
1781  tanDCorr *= aVal;
1782  //+ (-bVal/24. + 0.625*bVal*bVal*bVal + 5./12.*bVal*cVal)*aVal*aVal*aVal
1783  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<tanDist<<" vs "
1784  <<tanDist*(1.+tanDCorr)<<" corr "<<tanDist*tanDCorr<<std::endl;
1785  tanDist *= (1.+tanDCorr);
1786  } else {
1787  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"AACVal "<< fabs(aacVal)
1788  <<" = "<<aVal<<"**2 * "<<cVal<<" too large:: will not converge"<<std::endl;
1789  }
1790  } else {
1791  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"ABVal "<< fabs(aVal*bVal)
1792  <<" = "<<aVal<<" * "<<bVal<<" too large:: will not converge"<<std::endl;
1793  }
1794  }
1795  result = SteppingHelixStateInfo::OK;
1796  }
1797  break;
1798  case CONE_DT:
1799  {
1800  Point cVertex(pars[0], pars[1], pars[2]);
1801  if (cVertex.perp2() < 1e-10){
1802  //assumes the cone axis/vertex is along z
1803  Vector relV3 = sv.r3 - cVertex;
1804  double relV3mag = relV3.mag();
1805  // double relV3Theta = relV3.theta();
1806  double theta(pars[3]);
1807  // double dTheta = theta-relV3Theta;
1808  double sinTheta = sin(theta);
1809  double cosTheta = cos(theta);
1810  double cosV3Theta = relV3.z()/relV3mag;
1811  if (cosV3Theta>1) cosV3Theta=1;
1812  if (cosV3Theta<-1) cosV3Theta=-1;
1813  double sinV3Theta = sqrt(1.-cosV3Theta*cosV3Theta);
1814 
1815  double sinDTheta = sinTheta*cosV3Theta - cosTheta*sinV3Theta;//sin(dTheta);
1816  double cosDTheta = cosTheta*cosV3Theta + sinTheta*sinV3Theta;//cos(dTheta);
1817  bool isInside = sinTheta > sinV3Theta && cosTheta*cosV3Theta > 0;
1818  dist = isInside || cosDTheta > 0 ? relV3mag*sinDTheta : relV3mag;
1819  if (fabs(dist) > fastSkipDist){
1821  break;
1822  }
1823 
1824  double relV3phi=relV3.phi();
1825  double normPhi = isInside ?
1826  Geom::pi() + relV3phi : relV3phi;
1827  double normTheta = theta > Geom::pi()/2. ?
1828  ( isInside ? 1.5*Geom::pi() - theta : theta - Geom::pi()/2. )
1829  : ( isInside ? Geom::pi()/2. - theta : theta + Geom::pi()/2 );
1830  //this is a normVector from the cone to the point
1831  Vector norm; norm.setRThetaPhi(fabs(dist), normTheta, normPhi);
1832  double curP = sv.p3.mag(); double cosp3theta = sv.p3.z()/curP;
1833  if (cosp3theta>1) cosp3theta=1;
1834  if (cosp3theta<-1) cosp3theta=-1;
1835  double sineConeP = sinTheta*cosp3theta - cosTheta*sqrt(1.-cosp3theta*cosp3theta);
1836 
1837  double sinSolid = norm.dot(sv.p3)/(fabs(dist)*curP);
1838  tanDist = fabs(sinSolid) > fabs(sineConeP) ? dist/fabs(sinSolid) : dist/fabs(sineConeP);
1839 
1840 
1841  refDirection = sinSolid > 0 ?
1843  if (debug_){
1844  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Cone pars: theta "<<theta
1845  <<" normTheta "<<norm.theta()
1846  <<" p3Theta "<<sv.p3.theta()
1847  <<" sinD: "<< sineConeP
1848  <<" sinSolid "<<sinSolid;
1849  }
1850  if (debug_){
1851  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"refToDest:toCone the point is "
1852  <<(isInside? "in" : "out")<<"side the cone"
1853  <<std::endl;
1854  }
1855  }
1856  result = SteppingHelixStateInfo::OK;
1857  }
1858  break;
1859  // case CYLINDER_DT:
1860  // break;
1861  case PATHL_DT:
1862  {
1863  double curS = fabs(sv.path());
1864  dist = pars[PATHL_P] - curS;
1865  if (fabs(dist) > fastSkipDist){
1867  break;
1868  }
1869  refDirection = pars[PATHL_P] > 0 ?
1871  tanDist = dist;
1872  result = SteppingHelixStateInfo::OK;
1873  }
1874  break;
1875  case POINT_PCA_DT:
1876  {
1877  Point pDest(pars[0], pars[1], pars[2]);
1878  double curP = sv.p3.mag();
1879  dist = (sv.r3 - pDest).mag()+ 1e-24;//add a small number to avoid 1/0
1880  tanDist = (sv.r3.dot(sv.p3) - pDest.dot(sv.p3))/curP;
1881  //account for bending in magnetic field (quite approximate)
1882  double b0 = sv.bf.mag();
1883  if (b0>1.5e-6){
1884  double p0 = curP;
1885  double kVal = 2.99792458e-3*sv.q/p0*b0;
1886  double aVal = fabs(dist*kVal);
1887  tanDist *= 1./(1.+ aVal);
1888  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"corrected by aVal "<<aVal<<" to "<<tanDist;
1889  }
1890  refDirection = tanDist < 0 ?
1892  result = SteppingHelixStateInfo::OK;
1893  }
1894  break;
1895  case LINE_PCA_DT:
1896  {
1897  Point rLine(pars[0], pars[1], pars[2]);
1898  Vector dLine(pars[3], pars[4], pars[5]);
1899  dLine = (dLine - rLine);
1900  dLine *= 1./dLine.mag(); if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"dLine "<<dLine;
1901 
1902  Vector dR = sv.r3 - rLine;
1903  double curP = sv.p3.mag();
1904  Vector dRPerp = dR - dLine*(dR.dot(dLine));
1905  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"dRperp "<<dRPerp;
1906 
1907  dist = dRPerp.mag() + 1e-24;//add a small number to avoid 1/0
1908  tanDist = dRPerp.dot(sv.p3)/curP;
1909  //angle wrt line
1910  double cosAlpha2 = dLine.dot(sv.p3)/curP; cosAlpha2 *= cosAlpha2;
1911  tanDist *= 1./sqrt(fabs(1.-cosAlpha2)+1e-96);
1912  //correct for dPhi in magnetic field: this isn't made quite right here
1913  //(the angle between the line and the trajectory plane is neglected .. conservative)
1914  double b0 = sv.bf.mag();
1915  if (b0>1.5e-6){
1916  double p0 = curP;
1917  double kVal = 2.99792458e-3*sv.q/p0*b0;
1918  double aVal = fabs(dist*kVal);
1919  tanDist *= 1./(1.+ aVal);
1920  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"corrected by aVal "<<aVal<<" to "<<tanDist;
1921  }
1922  refDirection = tanDist < 0 ?
1924  result = SteppingHelixStateInfo::OK;
1925  }
1926  break;
1927  default:
1928  {
1929  //some large number
1930  dist = 1e12;
1931  tanDist = 1e12;
1932  refDirection = anyDirection;
1934  }
1935  break;
1936  }
1937 
1938  double tanDistConstrained = tanDist;
1939  if (fabs(tanDist) > 4.*fabs(dist) ) tanDistConstrained *= tanDist == 0 ? 0 : fabs(dist/tanDist*4.);
1940 
1941  if (debug_){
1942  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"refToDest input: dest"<<dest<<" pars[]: ";
1943  for (int i = 0; i < 6; i++){
1944  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<", "<<i<<" "<<pars[i];
1945  }
1946  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<std::endl;
1947  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"refToDest output: "
1948  <<"\t dist" << dist
1949  <<"\t tanDist "<< tanDist
1950  <<"\t tanDistConstr "<< tanDistConstrained
1951  <<"\t refDirection "<< refDirection
1952  <<std::endl;
1953  }
1954  tanDist = tanDistConstrained;
1955 
1956  return result;
1957 }
1958 
1962  double& dist, double& tanDist,
1963  double fastSkipDist, bool expectNewMagVolume, double maxStep) const{
1964 
1965  static const std::string metname = "SteppingHelixPropagator";
1967  const MagVolume* cVol = sv.magVol;
1968 
1969  if (cVol == 0) return result;
1970  const std::vector<VolumeSide>& cVolFaces(cVol->faces());
1971 
1972  double distToFace[6] = {0,0,0,0,0,0};
1973  double tanDistToFace[6] = {0,0,0,0,0,0};
1974  PropagationDirection refDirectionToFace[6] = {anyDirection, anyDirection, anyDirection, anyDirection, anyDirection, anyDirection};
1975  Result resultToFace[6] = {result, result, result, result, result, result};
1976  int iFDest = -1;
1977  int iDistMin = -1;
1978 
1979  unsigned int iFDestSorted[6] = {0,0,0,0,0,0};
1980  unsigned int nDestSorted =0;
1981  unsigned int nearParallels = 0;
1982 
1983  double curP = sv.p3.mag();
1984 
1985  if (debug_){
1986  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Trying volume "<<DDSolidShapesName::name(cVol->shapeType())
1987  <<" with "<<cVolFaces.size()<<" faces"<<std::endl;
1988  }
1989 
1990  unsigned int nFaces = cVolFaces.size();
1991  for (unsigned int iFace = 0; iFace < nFaces; ++iFace){
1992  if (iFace > 5){
1993  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Too many faces"<<std::endl;
1994  }
1995  if (debug_){
1996  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Start with face "<<iFace<<std::endl;
1997  }
1998 // const Plane* cPlane = dynamic_cast<const Plane*>(&cVolFaces[iFace].surface());
1999 // const Cylinder* cCyl = dynamic_cast<const Cylinder*>(&cVolFaces[iFace].surface());
2000 // const Cone* cCone = dynamic_cast<const Cone*>(&cVolFaces[iFace].surface());
2001  const Surface* cPlane = 0; //only need to know the loc->glob transform
2002  const Cylinder* cCyl = 0;
2003  const Cone* cCone = 0;
2004  if (typeid(cVolFaces[iFace].surface()) == typeid(const Plane&)){
2005  cPlane = &cVolFaces[iFace].surface();
2006  } else if (typeid(cVolFaces[iFace].surface()) == typeid(const Cylinder&)){
2007  cCyl = dynamic_cast<const Cylinder*>(&cVolFaces[iFace].surface());
2008  } else if (typeid(cVolFaces[iFace].surface()) == typeid(const Cone&)){
2009  cCone = dynamic_cast<const Cone*>(&cVolFaces[iFace].surface());
2010  } else {
2011  edm::LogWarning(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Could not cast a volume side surface to a known type"<<std::endl;
2012  }
2013 
2014  if (debug_){
2015  if (cPlane!=0) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The face is a plane at "<<cPlane<<std::endl;
2016  if (cCyl!=0) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The face is a cylinder at "<<cCyl<<std::endl;
2017  }
2018 
2019  double pars[6] = {0,0,0,0,0,0};
2020  DestType dType = UNDEFINED_DT;
2021  if (cPlane != 0){
2022  Point rPlane(cPlane->position().x(),cPlane->position().y(),cPlane->position().z());
2023  // = cPlane->toGlobal(LocalVector(0,0,1.)); nPlane = nPlane.unit();
2024  Vector nPlane(cPlane->rotation().zx(), cPlane->rotation().zy(), cPlane->rotation().zz()); nPlane /= nPlane.mag();
2025 
2026  if (sv.r3.z() < 0 || !vbField_->isZSymmetric() ){
2027  pars[0] = rPlane.x(); pars[1] = rPlane.y(); pars[2] = rPlane.z();
2028  pars[3] = nPlane.x(); pars[4] = nPlane.y(); pars[5] = nPlane.z();
2029  } else {
2030  pars[0] = rPlane.x(); pars[1] = rPlane.y(); pars[2] = -rPlane.z();
2031  pars[3] = nPlane.x(); pars[4] = nPlane.y(); pars[5] = -nPlane.z();
2032  }
2033  dType = PLANE_DT;
2034  } else if (cCyl != 0){
2035  if (debug_){
2036  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Cylinder at "<<cCyl->position()
2037  <<" rorated by "<<cCyl->rotation()
2038  <<std::endl;
2039  }
2040  pars[RADIUS_P] = cCyl->radius();
2041  dType = RADIUS_DT;
2042  } else if (cCone != 0){
2043  if (debug_){
2044  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Cone at "<<cCone->position()
2045  <<" rorated by "<<cCone->rotation()
2046  <<" vertex at "<<cCone->vertex()
2047  <<" angle of "<<cCone->openingAngle()
2048  <<std::endl;
2049  }
2050  if (sv.r3.z() < 0 || !vbField_->isZSymmetric()){
2051  pars[0] = cCone->vertex().x(); pars[1] = cCone->vertex().y();
2052  pars[2] = cCone->vertex().z();
2053  pars[3] = cCone->openingAngle();
2054  } else {
2055  pars[0] = cCone->vertex().x(); pars[1] = cCone->vertex().y();
2056  pars[2] = -cCone->vertex().z();
2057  pars[3] = Geom::pi() - cCone->openingAngle();
2058  }
2059  dType = CONE_DT;
2060  } else {
2061  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Unknown surface"<<std::endl;
2062  resultToFace[iFace] = SteppingHelixStateInfo::UNDEFINED;
2063  continue;
2064  }
2065  resultToFace[iFace] =
2066  refToDest(dType, sv, pars,
2067  distToFace[iFace], tanDistToFace[iFace], refDirectionToFace[iFace], fastSkipDist);
2068 
2069 
2070  if (resultToFace[iFace] != SteppingHelixStateInfo::OK){
2071  if (resultToFace[iFace] == SteppingHelixStateInfo::INACC) result = SteppingHelixStateInfo::INACC;
2072  }
2073 
2074 
2075 
2076  //keep those in right direction for later use
2077  if (resultToFace[iFace] == SteppingHelixStateInfo::OK){
2078  double invDTFPosiF = 1./(1e-32+fabs(tanDistToFace[iFace]));
2079  double dSlope = fabs(distToFace[iFace])*invDTFPosiF;
2080  double maxStepL = maxStep> 100 ? 100 : maxStep; if (maxStepL < 10) maxStepL = 10.;
2081  bool isNearParallel = fabs(tanDistToFace[iFace]) + 100.*curP*dSlope < maxStepL //
2082  //a better choice is to use distance to next check of mag volume instead of 100cm; the last is for ~1.5arcLength(4T)+tandistance< maxStep
2083  && dSlope < 0.15 ; //
2084  if (refDirectionToFace[iFace] == dir || isNearParallel){
2085  if (isNearParallel) nearParallels++;
2086  iFDestSorted[nDestSorted] = iFace;
2087  nDestSorted++;
2088  }
2089  }
2090 
2091  //pick a shortest distance here (right dir only for now)
2092  if (refDirectionToFace[iFace] == dir){
2093  if (iDistMin == -1) iDistMin = iFace;
2094  else if (fabs(distToFace[iFace]) < fabs(distToFace[iDistMin])) iDistMin = iFace;
2095  }
2096  if (debug_)
2097  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<cVol<<" "<<iFace<<" "
2098  <<tanDistToFace[iFace]<<" "<<distToFace[iFace]<<" "<<refDirectionToFace[iFace]<<" "<<dir<<std::endl;
2099 
2100  }
2101 
2102  for (unsigned int i = 0;i<nDestSorted; ++i){
2103  unsigned int iMax = nDestSorted-i-1;
2104  for (unsigned int j=0;j<nDestSorted-i; ++j){
2105  if (fabs(tanDistToFace[iFDestSorted[j]]) > fabs(tanDistToFace[iFDestSorted[iMax]]) ){
2106  iMax = j;
2107  }
2108  }
2109  unsigned int iTmp = iFDestSorted[nDestSorted-i-1];
2110  iFDestSorted[nDestSorted-i-1] = iFDestSorted[iMax];
2111  iFDestSorted[iMax] = iTmp;
2112  }
2113 
2114  if (debug_){
2115  for (unsigned int i=0;i<nDestSorted;++i){
2116  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<cVol<<" "<<i<<" "<<iFDestSorted[i]<<" "<<tanDistToFace[iFDestSorted[i]]<<std::endl;
2117  }
2118  }
2119 
2120  //now go from the shortest to the largest distance hoping to get a point in the volume.
2121  //other than in case of a near-parallel travel this should stop after the first try
2122 
2123  for (unsigned int i=0; i<nDestSorted;++i){
2124  iFDest = iFDestSorted[i];
2125 
2126  double sign = dir == alongMomentum ? 1. : -1.;
2127  Point gPointEst(sv.r3);
2128  Vector lDelta(sv.p3); lDelta *= sign/curP*fabs(distToFace[iFDest]);
2129  gPointEst += lDelta;
2130  if (debug_){
2131  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Linear est point "<<gPointEst
2132  <<" for iFace "<<iFDest<<std::endl;
2133  }
2134  GlobalPoint gPointEstNegZ(gPointEst.x(), gPointEst.y(), -fabs(gPointEst.z()));
2135  GlobalPoint gPointEstNorZ(gPointEst.x(), gPointEst.y(), gPointEst.z() );
2136  if ( (vbField_->isZSymmetric() && cVol->inside(gPointEstNegZ))
2137  || ( !vbField_->isZSymmetric() && cVol->inside(gPointEstNorZ) ) ){
2138  if (debug_){
2139  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is inside the volume"<<std::endl;
2140  }
2141 
2142  result = SteppingHelixStateInfo::OK;
2143  dist = distToFace[iFDest];
2144  tanDist = tanDistToFace[iFDest];
2145  if (debug_){
2146  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got a point near closest boundary -- face "<<iFDest<<std::endl;
2147  }
2148  break;
2149  }
2150  }
2151 
2152  if (result != SteppingHelixStateInfo::OK && expectNewMagVolume){
2153  double sign = dir == alongMomentum ? 1. : -1.;
2154 
2155  //check if it's a wrong volume situation
2156  if (nDestSorted-nearParallels > 0 ) result = SteppingHelixStateInfo::WRONG_VOLUME;
2157  else {
2158  //get here if all faces in the corr direction were skipped
2159  Point gPointEst(sv.r3);
2160  double lDist = iDistMin == -1 ? fastSkipDist : fabs(distToFace[iDistMin]);
2161  if (lDist > fastSkipDist) lDist = fastSkipDist;
2162  Vector lDelta(sv.p3); lDelta *= sign/curP*lDist;
2163  gPointEst += lDelta;
2164  if (debug_){
2165  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Linear est point to shortest dist "<<gPointEst
2166  <<" for iFace "<<iDistMin<<" at distance "<<lDist*sign<<std::endl;
2167  }
2168  GlobalPoint gPointEstNegZ(gPointEst.x(), gPointEst.y(), -fabs(gPointEst.z()));
2169  GlobalPoint gPointEstNorZ(gPointEst.x(), gPointEst.y(), gPointEst.z() );
2170  if ( (vbField_->isZSymmetric() && cVol->inside(gPointEstNegZ))
2171  || ( !vbField_->isZSymmetric() && cVol->inside(gPointEstNorZ) ) ){
2172  if (debug_){
2173  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is inside the volume"<<std::endl;
2174  }
2175 
2176  }else {
2178  }
2179  }
2180 
2181  if (result == SteppingHelixStateInfo::WRONG_VOLUME){
2182  dist = sign*0.05;
2183  tanDist = dist*1.01;
2184  if( debug_){
2185  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Wrong volume located: return small dist, tandist"<<std::endl;
2186  }
2187  }
2188  }
2189 
2190  if (result == SteppingHelixStateInfo::INACC){
2191  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"All faces are too far"<<std::endl;
2192  } else if (result == SteppingHelixStateInfo::WRONG_VOLUME){
2193  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Appear to be in a wrong volume"<<std::endl;
2194  } else if (result != SteppingHelixStateInfo::OK){
2195  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Something else went wrong"<<std::endl;
2196  }
2197 
2198 
2199  return result;
2200 }
2201 
2202 
2206  double& dist, double& tanDist,
2207  double fastSkipDist) const{
2208 
2209  static const std::string metname = "SteppingHelixPropagator";
2211 
2212  double parLim[6] = {sv.rzLims.rMin, sv.rzLims.rMax,
2213  sv.rzLims.zMin, sv.rzLims.zMax,
2214  sv.rzLims.th1, sv.rzLims.th2 };
2215 
2216  double distToFace[4] = {0,0,0,0};
2217  double tanDistToFace[4] = {0,0,0,0};
2218  PropagationDirection refDirectionToFace[4] = {anyDirection, anyDirection, anyDirection, anyDirection};
2219  Result resultToFace[4] = {result, result, result, result};
2220  int iFDest = -1;
2221 
2222  double curP = sv.p3.mag();
2223 
2224  for (unsigned int iFace = 0; iFace < 4; iFace++){
2225  if (debug_){
2226  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Start with mat face "<<iFace<<std::endl;
2227  }
2228 
2229  double pars[6] = {0,0,0,0,0,0};
2230  DestType dType = UNDEFINED_DT;
2231  if (iFace > 1){
2232  if (fabs(parLim[iFace+2])< 1e-6){//plane
2233  if (sv.r3.z() < 0){
2234  pars[0] = 0; pars[1] = 0; pars[2] = -parLim[iFace];
2235  pars[3] = 0; pars[4] = 0; pars[5] = 1;
2236  } else {
2237  pars[0] = 0; pars[1] = 0; pars[2] = parLim[iFace];
2238  pars[3] = 0; pars[4] = 0; pars[5] = 1;
2239  }
2240  dType = PLANE_DT;
2241  } else {
2242  if (sv.r3.z() > 0){
2243  pars[0] = 0; pars[1] = 0;
2244  pars[2] = parLim[iFace];
2245  pars[3] = Geom::pi()/2. - parLim[iFace+2];
2246  } else {
2247  pars[0] = 0; pars[1] = 0;
2248  pars[2] = - parLim[iFace];
2249  pars[3] = Geom::pi()/2. + parLim[iFace+2];
2250  }
2251  dType = CONE_DT;
2252  }
2253  } else {
2254  pars[RADIUS_P] = parLim[iFace];
2255  dType = RADIUS_DT;
2256  }
2257 
2258  resultToFace[iFace] =
2259  refToDest(dType, sv, pars,
2260  distToFace[iFace], tanDistToFace[iFace], refDirectionToFace[iFace], fastSkipDist);
2261 
2262  if (resultToFace[iFace] != SteppingHelixStateInfo::OK){
2263  if (resultToFace[iFace] == SteppingHelixStateInfo::INACC) result = SteppingHelixStateInfo::INACC;
2264  continue;
2265  }
2266  if (refDirectionToFace[iFace] == dir || fabs(distToFace[iFace]) < 2e-2*fabs(tanDistToFace[iFace]) ){
2267  double sign = dir == alongMomentum ? 1. : -1.;
2268  Point gPointEst(sv.r3);
2269  Vector lDelta(sv.p3); lDelta *= sign*fabs(distToFace[iFace])/curP;
2270  gPointEst += lDelta;
2271  if (debug_){
2272  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Linear est point "<<gPointEst
2273  <<std::endl;
2274  }
2275  double lZ = fabs(gPointEst.z());
2276  double lR = gPointEst.perp();
2277  double tan4 = parLim[4] == 0 ? 0 : tan(parLim[4]);
2278  double tan5 = parLim[5] == 0 ? 0 : tan(parLim[5]);
2279  if ( (lZ - parLim[2]) > lR*tan4
2280  && (lZ - parLim[3]) < lR*tan5
2281  && lR > parLim[0] && lR < parLim[1]){
2282  if (debug_){
2283  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is inside the volume"<<std::endl;
2284  }
2285  //OK, guessed a point still inside the volume
2286  if (iFDest == -1){
2287  iFDest = iFace;
2288  } else {
2289  if (fabs(tanDistToFace[iFDest]) > fabs(tanDistToFace[iFace])){
2290  iFDest = iFace;
2291  }
2292  }
2293  } else {
2294  if (debug_){
2295  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"The point is NOT inside the volume"<<std::endl;
2296  }
2297  }
2298  }
2299 
2300  }
2301  if (iFDest != -1){
2302  result = SteppingHelixStateInfo::OK;
2303  dist = distToFace[iFDest];
2304  tanDist = tanDistToFace[iFDest];
2305  if (debug_){
2306  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Got a point near closest boundary -- face "<<iFDest<<std::endl;
2307  }
2308  } else {
2309  if (debug_){
2310  LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Failed to find a dest point inside the volume"<<std::endl;
2311  }
2312  }
2313 
2314  return result;
2315 }
2316 
2317 
2319  static const std::string metname = "SteppingHelixPropagator";
2320  if (vol == 0) return false;
2321  /*
2322  const MFGrid* mGrid = reinterpret_cast<const MFGrid*>(vol->provider());
2323  std::vector<int> dims(mGrid->dimensions());
2324 
2325  LocalVector lVCen(mGrid->nodeValue(dims[0]/2, dims[1]/2, dims[2]/2));
2326  LocalVector lVZLeft(mGrid->nodeValue(dims[0]/2, dims[1]/2, dims[2]/5));
2327  LocalVector lVZRight(mGrid->nodeValue(dims[0]/2, dims[1]/2, (dims[2]*4)/5));
2328 
2329  double mag2VCen = lVCen.mag2();
2330  double mag2VZLeft = lVZLeft.mag2();
2331  double mag2VZRight = lVZRight.mag2();
2332 
2333  bool result = false;
2334  if (mag2VCen > 0.6 && mag2VZLeft > 0.6 && mag2VZRight > 0.6){
2335  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Volume is magnetic, located at "<<vol->position()<<std::endl;
2336  result = true;
2337  } else {
2338  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Volume is not magnetic, located at "<<vol->position()<<std::endl;
2339  }
2340 
2341  */
2342  bool result = vol->isIron();
2343  if (result){
2344  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Volume is magnetic, located at "<<vol->position()<<std::endl;
2345  } else {
2346  if (debug_) LogTrace(metname)<<std::setprecision(17)<<std::setw(20)<<std::scientific<<"Volume is not magnetic, located at "<<vol->position()<<std::endl;
2347  }
2348 
2349  return result;
2350 }
type
Definition: HCALResponse.h:22
double z0() const
z coordinate
Definition: BeamSpot.h:69
int i
Definition: DBlmapReader.cc:9
T mag2() const
Definition: PV3DBase.h:65
Definition: Cone.h:20
DDSolidShape shapeType() const
Definition: MagVolume.h:31
virtual PropagationDirection propagationDirection() const
Definition: Propagator.h:143
tuple pp
Definition: createTree.py:15
Result refToDest(DestType dest, const SteppingHelixPropagator::StateInfo &sv, const double pars[6], double &dist, double &tanDist, PropagationDirection &refDirection, double fastSkipDist=1e12) const
(Internals) determine distance and direction from the current position to the plane ...
const std::string metname
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
void getFreeState(FreeTrajectoryState &fts) const
convert internal structure into the fts
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
ROOT::Math::SMatrixIdentity AlgebraicMatrixID
bool isYokeVolume(const MagVolume *vol) const
check if it&#39;s a yoke/iron based on this MagVol internals
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
SteppingHelixStateInfo::VolumeBounds MatBounds
Geom::Theta< T > theta() const
T y() const
Definition: PV3DBase.h:62
AlgebraicSymMatrix55 covCurv
const MagVolume * findVolume(const GlobalPoint &gp) const
std::pair< TrajectoryStateOnSurface, double > TsosPP
int cIndex_(int ind) const
(Internals) circular index for array of transients
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
double getDeDx(const SteppingHelixPropagator::StateInfo &sv, double &dEdXPrime, double &radX0, MatBounds &rzLims) const
estimate average (in fact smth. close to MPV and median) energy loss per unit path length ...
PropagationDirection
std::pair< FreeTrajectoryState, double > FtsPP
Result refToMatVolume(const SteppingHelixPropagator::StateInfo &sv, PropagationDirection dir, double &dist, double &tanDist, double fastSkipDist=1e12) const
double charge(const std::vector< uint8_t > &Ampls)
Definition: Plane.h:17
void loadState(SteppingHelixPropagator::StateInfo &svCurrent, const SteppingHelixPropagator::Vector &p3, const SteppingHelixPropagator::Point &r3, int charge, PropagationDirection dir, const AlgebraicSymMatrix55 &covCurv) const
virtual ::GlobalVector inTesla(const ::GlobalPoint &gp) const
Definition: MagVolume.h:46
T zx() const
double dydz() const
dydz slope
Definition: BeamSpot.h:85
GlobalPoint vertex() const
Global position of the cone vertex.
Definition: Cone.h:50
T zz() const
T mag() const
Definition: PV3DBase.h:66
virtual bool inside(const GlobalPoint &gp, double tolerance=0.) const =0
Scalar radius() const
Radius of the cylinder.
Definition: Cylinder.h:55
static const char * name(DDSolidShape s)
Definition: DDSolidShapes.h:21
TrajectoryStateOnSurface getStateOnSurface(const Surface &surf, bool returnTangentPlane=false) const
T sqrt(T t)
Definition: SSEVec.h:46
LocalPoint toLocal(const GlobalPoint &gp) const
const VolumeBasedMagneticField * vbField_
T z() const
Definition: PV3DBase.h:63
tuple result
Definition: query.py:137
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
T zy() const
int j
Definition: DBlmapReader.cc:9
AlgebraicSymMatrix55 matDCovCurv
virtual TrajectoryStateOnSurface propagate(const FreeTrajectoryState &ftsStart, const Plane &pDest) const
Propagate to Plane given a starting point.
Common base class.
#define LogTrace(id)
double dxdz() const
dxdz slope
Definition: BeamSpot.h:83
const MagneticField * field
bool isIron() const
Temporary hack to pass information on material. Will eventually be replaced!
Definition: MagVolume.h:51
bool makeAtomStep(SteppingHelixPropagator::StateInfo &svCurrent, SteppingHelixPropagator::StateInfo &svNext, double dS, PropagationDirection dir, SteppingHelixPropagator::Fancy fancy) const
main stepping function: compute next state vector after a step of length dS
Result refToMagVolume(const SteppingHelixPropagator::StateInfo &sv, PropagationDirection dir, double &dist, double &tanDist, double fastSkipDist=1e12, bool expectNewMagVolume=false, double maxStep=1e12) const
virtual std::pair< TrajectoryStateOnSurface, double > propagateWithPath(const FreeTrajectoryState &ftsStart, const Plane &pDest) const
virtual const std::vector< VolumeSide > & faces() const =0
Access to volume faces.
const AlgebraicSymMatrix55 unit55_
void setRep(SteppingHelixPropagator::Basis &rep, const SteppingHelixPropagator::Vector &tau) const
Set/compute basis vectors for local coordinates at current step (called by incrementState) ...
LocalVector fieldInTesla(const LocalPoint &lp) const
Definition: MagVolume.cc:11
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
StateInfo svBuf_[MAX_POINTS+1]
double pi()
Definition: Pi.h:31
double y0() const
y coordinate
Definition: BeamSpot.h:67
void setIState(const SteppingHelixStateInfo &sStart) const
(Internals) Init starting point
const RotationType & rotation() const
dbl *** dir
Definition: mlp_gen.cc:35
const double epsilon
void getNextState(const SteppingHelixPropagator::StateInfo &svPrevious, SteppingHelixPropagator::StateInfo &svNext, double dP, const SteppingHelixPropagator::Vector &tau, const SteppingHelixPropagator::Vector &drVec, double dS, double dX0, const AlgebraicMatrix55 &dCovCurv) const
T x() const
Definition: PV3DBase.h:61
const PositionType & position() const
const MagneticField * field_
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
static const std::string ResultName[MAX_RESULT]
Geom::Theta< float > openingAngle() const
Angle of the cone.
Definition: Cone.h:53
double p3[4]
Definition: TauolaWrapper.h:91
double x0() const
x coordinate
Definition: BeamSpot.h:65
Definition: DDAxes.h:10