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