CMS 3D CMS Logo

TrajectoryManager.cc
Go to the documentation of this file.
1 //Framework Headers
3 
4 //CMSSW Headers
9 
10 // Tracker reco geometry headers
19 //#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
20 
21 //FAMOS Headers
32 
33 //#include "FastSimulation/Utilities/interface/Histos.h"
34 //#include "FastSimulation/Utilities/interface/FamosLooses.h"
35 // Numbering scheme
36 
37 //#define FAMOS_DEBUG
38 #ifdef FAMOS_DEBUG
41 #endif
42 
43 #include <list>
44 
46 
48  const edm::ParameterSet& matEff,
50  const edm::ParameterSet& decays) :
51  mySimEvent(aSimEvent),
52  _theGeometry(nullptr),
53  _theFieldMap(nullptr),
54  theMaterialEffects(nullptr),
55  myDecayEngine(nullptr),
56  theGeomTracker(nullptr),
57  theGeomSearchTracker(nullptr),
58  theLayerMap(56, static_cast<const DetLayer*>(nullptr)), // reserve space for layers here
59  theNegLayerOffset(27),
60  // myHistos(0),
61  use_hardcoded(true)
62 
63 {
64  //std::cout << "TrajectoryManager.cc 1 use_hardcoded = " << use_hardcoded << std::endl;
65  use_hardcoded = matEff.getParameter<bool>("use_hardcoded_geometry");
66 
67  // Initialize Bthe stable particle decay engine
68  if ( decays.getParameter<bool>("ActivateDecays")) {
70  distCut = decays.getParameter<double>("DistCut");
71  }
72  // Initialize the Material Effects updator, if needed
73  if ( matEff.getParameter<bool>("PairProduction") ||
74  matEff.getParameter<bool>("Bremsstrahlung") ||
75  matEff.getParameter<bool>("MuonBremsstrahlung") ||
76  matEff.getParameter<bool>("EnergyLoss") ||
77  matEff.getParameter<bool>("MultipleScattering") ||
78  matEff.getParameter<bool>("NuclearInteraction")
79  )
80  theMaterialEffects = new MaterialEffects(matEff);
81 
82  // Save SimHits according to Optiom
83  // Only the hits from first half loop is saved
84  firstLoop = simHits.getUntrackedParameter<bool>("firstLoop",true);
85  // Only if pT>pTmin are the hits saved
86  pTmin = simHits.getUntrackedParameter<double>("pTmin",0.5);
87 
88  /*
89  // Get the Famos Histos pointer
90  myHistos = Histos::instance();
91 
92  // Initialize a few histograms
93 
94  myHistos->book("h302",1210,-121.,121.,1210,-121.,121.);
95  myHistos->book("h300",1210,-121.,121.,1210,-121.,121.);
96  myHistos->book("h301",1200,-300.,300.,1210,-121.,121.);
97  myHistos->book("h303",1200,-300.,300.,1210,-121.,121.);
98  */
99 }
100 
101 void
103  const TrackerInteractionGeometry* interactionGeometry,
104  const MagneticFieldMap* aFieldMap)
105 {
106 
107  // Initialize the reco tracker geometry
108  theGeomSearchTracker = geomSearchTracker;
109 
110  // Initialize the simplified tracker geometry
111  _theGeometry = interactionGeometry;
112 
114 
115  // Initialize the magnetic field
116  _theFieldMap = aFieldMap;
117 
118 }
119 
120 void
122 
123  theGeomTracker = geomTracker;
124 
125 }
126 
129  return _theGeometry;
130 }
131 
133 
134  if ( myDecayEngine ) delete myDecayEngine;
136 
137  //Write the histograms
138  /*
139  myHistos->put("histos.root");
140  if ( myHistos ) delete myHistos;
141  */
142 }
143 
144 void
146 {
147 
148  // Clear the hits of the previous event
149  // thePSimHits->clear();
150  thePSimHits.clear();
151 
152  // The new event
153  XYZTLorentzVector myBeamPipe = XYZTLorentzVector(0.,2.5, 9999999.,0.);
154 
155  std::list<TrackerLayer>::const_iterator cyliter;
156 
157  // bool debug = mySimEvent->id().event() == 8;
158 
159  // Loop over the particles (watch out: increasing upper limit!)
160  for( int fsimi=0; fsimi < (int) mySimEvent->nTracks(); ++fsimi) {
161  // If the particle has decayed inside the beampipe, or decays
162  // immediately, there is nothing to do
163  //if ( debug ) std::cout << mySimEvent->track(fsimi) << std::endl;
164  //if ( debug ) std::cout << "Not yet at end vertex ? " << mySimEvent->track(fsimi).notYetToEndVertex(myBeamPipe) << std::endl;
165  if( !mySimEvent->track(fsimi).notYetToEndVertex(myBeamPipe) ) continue;
166  mySimEvent->track(fsimi).setPropagate();
167 
168  // Get the geometry elements
169  cyliter = _theGeometry->cylinderBegin();
170  // Prepare the propagation
172  //The real work starts here
173  int success = 1;
174  int sign = +1;
175  int loop = 0;
176  int cyl = 0;
177 
178  // Find the initial cylinder to propagate to.
179  for ( ; cyliter != _theGeometry->cylinderEnd() ; ++cyliter ) {
180 
181  PP.setPropagationConditions(*cyliter);
182  if ( PP.inside() && !PP.onSurface() ) break;
183  ++cyl;
184 
185  }
186 
187  // The particle has a pseudo-rapidity (position or momentum direction)
188  // in excess of 3.0. Just simply go to the last tracker layer
189  // without bothering with all the details of the propagation and
190  // material effects.
191  // 08/02/06 - pv: increase protection from 0.99 (eta=2.9932) to 0.9998 (eta=4.9517)
192  // to simulate material effects at large eta
193  // if above 0.99: propagate to the last tracker cylinder where the material is concentrated!
194  double ppcos2T = PP.particle().cos2Theta();
195  double ppcos2V = PP.particle().cos2ThetaV();
196 
197  if(use_hardcoded){
198  if ( ( ppcos2T > 0.99 && ppcos2T < 0.9998 ) && ( cyl == 0 || ( ppcos2V > 0.99 && ppcos2V < 0.9998 ) ) ){
199  if ( cyliter != _theGeometry->cylinderEnd() ) {
200  cyliter = _theGeometry->cylinderEnd();
201  --cyliter;
202  }
203  // if above 0.9998: don't propagate at all (only to the calorimeters directly)
204  } else if ( ppcos2T > 0.9998 && ( cyl == 0 || ppcos2V > 0.9998 ) ) {
205  cyliter = _theGeometry->cylinderEnd();
206  }
207  }
208  else {
209  if ( ppcos2T > 0.9998 && ( cyl == 0 || ppcos2V > 0.9998 ) ) {
210  cyliter = _theGeometry->cylinderEnd();
211  }
212  }
213 
214  // Loop over the cylinders
215  while ( cyliter != _theGeometry->cylinderEnd() &&
216  loop<100 && // No more than 100 loops
217  mySimEvent->track(fsimi).notYetToEndVertex(PP.particle().vertex())) { // The particle decayed
218 
219  // Skip layers with no material (kept just for historical reasons)
220  if ( cyliter->surface().mediumProperties().radLen() < 1E-10 ) {
221  ++cyliter; ++cyl;
222  continue;
223  }
224 
225  // Pathological cases:
226  // To prevent from interacting twice in a row with the same layer
227  // bool escapeBarrel = (PP.getSuccess() == -1 && success == 1);
228  bool escapeBarrel = PP.getSuccess() == -1;
229  bool escapeEndcap = (PP.getSuccess() == -2 && success == 1);
230  // To break the loop
231  bool fullPropagation =
232  (PP.getSuccess() <= 0 && success==0) || escapeEndcap;
233 
234  if ( escapeBarrel ) {
235  ++cyliter; ++cyl;
236  while (cyliter != _theGeometry->cylinderEnd() && cyliter->forward() ) {
237  sign=1; ++cyliter; ++cyl;
238  }
239 
240  if ( cyliter == _theGeometry->cylinderEnd() ) {
241  --cyliter; --cyl; fullPropagation=true;
242  }
243 
244  }
245 
246  // Define the propagation conditions
247  PP.setPropagationConditions(*cyliter,!fullPropagation);
248  if ( escapeEndcap ) PP.increaseRCyl(0.0005);
249 
250  // Remember last propagation outcome
251  success = PP.getSuccess();
252 
253  // Propagation was not successful :
254  // Change the sign of the cylinder increment and count the loops
255  if ( !PP.propagateToBoundSurface(*cyliter) ||
256  PP.getSuccess()<=0) {
257  sign = -sign;
258  ++loop;
259  }
260 
261  // The particle may have decayed on its way... in which the daughters
262  // have to be added to the event record
263  if ( PP.hasDecayed() || (!mySimEvent->track(fsimi).nDaughters() && pdg::cTau(PP.particle().pid(),mySimEvent->theTable())<1E-3 ) ) {
264  updateWithDaughters(PP, fsimi, random);
265  break;
266  }
267 
268  // Exit by the endcaps or innermost cylinder :
269  // Positive cylinder increment
270  if ( PP.getSuccess()==2 || cyliter==_theGeometry->cylinderBegin() )
271  sign = +1;
272 
273  // Successful propagation to a cylinder, with some Material :
274  if( PP.getSuccess() > 0 && PP.onFiducial() ) {
275 
276  bool saveHit =
277  ( (loop==0 && sign>0) || !firstLoop ) && // Save only first half loop
278  PP.particle().charge()!=0. && // Consider only charged particles
279  cyliter->sensitive() && // Consider only sensitive layers
280  PP.particle().Perp2()>pTmin*pTmin; // Consider only pT > pTmin
281 
282  // Material effects are simulated there
283  if ( theMaterialEffects )
284  theMaterialEffects->interact(*mySimEvent,*cyliter,PP,fsimi, random);
285 
286  // There is a PP.setXYZT=(0,0,0,0) if bremss fails
287  saveHit &= PP.particle().E()>1E-6;
288 
289  if ( saveHit ) {
290  // Consider only active layers
291  if ( cyliter->sensitive() ) {
292  // Add information to the FSimTrack (not yet available)
293  // myTrack.addSimHit(PP.particle(),layer);
294 
295  // Return one or two (for overlap regions) PSimHits in the full
296  // tracker geometry
297  if ( theGeomTracker )
298  createPSimHits(*cyliter, PP, thePSimHits[fsimi], fsimi,mySimEvent->track(fsimi).type(), tTopo);
299 
300  /*
301  myHistos->fill("h302",PP.particle().X() ,PP.particle().Y());
302  if ( sin(PP.particle().vertex().Phi()) > 0. )
303  myHistos->fill("h303",PP.particle().Z(),PP.particle().R());
304  else
305  myHistos->fill("h303",PP.Z(),-PP.particle().R());
306  */
307 
308  }
309  }
310 
311  // Fill Histos (~poor man event display)
312  /*
313  myHistos->fill("h300",PP.particle().x(),PP.particle().y());
314  if ( sin(PP.particle().vertex().phi()) > 0. )
315  myHistos->fill("h301",PP.particle().z(),sqrt(PP.particle().vertex().Perp2()));
316  else
317  myHistos->fill("h301",PP.particle().z(),-sqrt(PP.particle().vertex().Perp2()));
318  */
319 
320  //The particle may have lost its energy in the material
321  if ( mySimEvent->track(fsimi).notYetToEndVertex(PP.particle().vertex()) &&
322  !mySimEvent->filter().acceptParticle(PP.particle()) )
323  mySimEvent->addSimVertex(PP.particle().vertex(),fsimi, FSimVertexType::END_VERTEX);
324 
325  }
326 
327  // Stop here if the particle has reached an end
328  if ( mySimEvent->track(fsimi).notYetToEndVertex(PP.particle().vertex()) ) {
329 
330  // Otherwise increment the cylinder iterator
331  // do {
332  if (sign==1) {++cyliter;++cyl;}
333  else {--cyliter;--cyl;}
334 
335  // Check if the last surface has been reached
336  if( cyliter==_theGeometry->cylinderEnd()) {
337 
338  // Try to propagate to the ECAL in half a loop
339  // Note: Layer1 = ECAL Barrel entrance, or Preshower
340  // entrance, or ECAL Endcap entrance (in the corner)
341  PP.propagateToEcal();
342  // PP.propagateToPreshowerLayer1();
343 
344  // If it is not possible, try go back to the last cylinder
345  if(PP.getSuccess()==0) {
346  --cyliter; --cyl; sign = -sign;
347  PP.setPropagationConditions(*cyliter);
348  PP.propagateToBoundSurface(*cyliter);
349 
350  // If there is definitely no way, leave it here.
351  if(PP.getSuccess()<0) {++cyliter; ++cyl;}
352 
353  }
354 
355  // Check if the particle has decayed on the way to ECAL
356  if ( PP.hasDecayed() )
357  updateWithDaughters(PP, fsimi, random);
358 
359  }
360  }
361 
362  }
363 
364  // Propagate all particles without a end vertex to the Preshower,
365  // theECAL and the HCAL.
366  if ( mySimEvent->track(fsimi).notYetToEndVertex(PP.particle().vertex()) )
367  propagateToCalorimeters(PP, fsimi, random);
368 
369  }
370 
371  // Save the information from Nuclear Interaction Generation
373 
374 }
375 
376 void
378 
379  FSimTrack& myTrack = mySimEvent->track(fsimi);
380 
381  // Set the position and momentum at the end of the tracker volume
382  myTrack.setTkPosition(PP.particle().vertex().Vect());
383  myTrack.setTkMomentum(PP.particle().momentum());
384 
385  // Propagate to Preshower Layer 1
386  PP.propagateToPreshowerLayer1(false);
387  if ( PP.hasDecayed() ) {
388  updateWithDaughters(PP, fsimi, random);
389  return;
390  }
391  if ( myTrack.notYetToEndVertex(PP.particle().vertex()) && PP.getSuccess() > 0 )
392  myTrack.setLayer1(PP.particle(),PP.getSuccess());
393 
394  // Propagate to Preshower Layer 2
395  PP.propagateToPreshowerLayer2(false);
396  if ( PP.hasDecayed() ) {
397  updateWithDaughters(PP, fsimi, random);
398  return;
399  }
400  if ( myTrack.notYetToEndVertex(PP.particle().vertex()) && PP.getSuccess() > 0 )
401  myTrack.setLayer2(PP.particle(),PP.getSuccess());
402 
403  // Propagate to Ecal Endcap
404  PP.propagateToEcalEntrance(false);
405  if ( PP.hasDecayed() ) {
406  updateWithDaughters(PP, fsimi, random);
407  return;
408  }
409  if ( myTrack.notYetToEndVertex(PP.particle().vertex()) )
410  myTrack.setEcal(PP.particle(),PP.getSuccess());
411 
412  // Propagate to HCAL entrance
413  PP.propagateToHcalEntrance(false);
414  if ( PP.hasDecayed() ) {
415  updateWithDaughters(PP,fsimi, random);
416  return;
417  }
418  if ( myTrack.notYetToEndVertex(PP.particle().vertex()) )
419  myTrack.setHcal(PP.particle(),PP.getSuccess());
420 
421  // Propagate to VFCAL entrance
422  PP.propagateToVFcalEntrance(false);
423  if ( PP.hasDecayed() ) {
424  updateWithDaughters(PP,fsimi, random);
425  return;
426  }
427  if ( myTrack.notYetToEndVertex(PP.particle().vertex()) )
428  myTrack.setVFcal(PP.particle(),PP.getSuccess());
429 
430 }
431 
432 bool
434 
435  std::list<TrackerLayer>::const_iterator cyliter;
436  bool done = false;
437 
438  // Get the geometry elements
439  cyliter = _theGeometry->cylinderBegin();
440 
441  // Find the layer to propagate to.
442  for ( ; cyliter != _theGeometry->cylinderEnd() ; ++cyliter ) {
443 
444  if ( layer != cyliter->layerNumber() ) continue;
445 
446  PP.setPropagationConditions(*cyliter);
447 
448  done =
449  PP.propagateToBoundSurface(*cyliter) &&
450  PP.getSuccess() > 0 &&
451  PP.onFiducial();
452 
453  break;
454 
455  }
456 
457  return done;
458 
459 }
460 
461 void
463 
464  // The particle was already decayed in the GenEvent, but still the particle was
465  // allowed to propagate (for magnetic field bending, for material effects, etc...)
466  // Just modify the momentum of the daughters in that case
467  unsigned nDaugh = mySimEvent->track(fsimi).nDaughters();
468  if ( nDaugh ) {
469 
470  // Move the vertex
471  unsigned vertexId = mySimEvent->track(fsimi).endVertex().id();
472  mySimEvent->vertex(vertexId).setPosition(PP.particle().vertex());
473 
474  // Before-propagation and after-propagation momentum and vertex position
475  XYZTLorentzVector momentumBefore = mySimEvent->track(fsimi).momentum();
476  const XYZTLorentzVector& momentumAfter = PP.particle().momentum();
477  double magBefore = std::sqrt(momentumBefore.Vect().mag2());
478  double magAfter = std::sqrt(momentumAfter.Vect().mag2());
479  // Rotation to be applied
480  XYZVector axis = momentumBefore.Vect().Cross(momentumAfter.Vect());
481  double angle = std::acos(momentumBefore.Vect().Dot(momentumAfter.Vect())/(magAfter*magBefore));
482  Rotation r(axis,angle);
483  // Rescaling to be applied
484  double rescale = magAfter/magBefore;
485 
486  // Move, rescale and rotate daugthers, grand-daughters, etc.
487  moveAllDaughters(fsimi,r,rescale);
488 
489  // The particle is not decayed in the GenEvent, decay it with PYTHIA
490  } else {
491 
492  // Decays are not activated : do nothing
493  if ( !myDecayEngine ) return;
494 
495  // Invoke PYDECY (Pythia6) or Pythia8 to decay the particle and get the daughters
497 
498  // Update the FSimEvent with an end vertex and with the daughters
499  if ( !daughters.empty() ) {
500  double distMin = 1E99;
501  int theClosestChargedDaughterId = -1;
502  DaughterParticleIterator daughter = daughters.begin();
503 
504  int ivertex = mySimEvent->addSimVertex(daughter->vertex(),fsimi,
506 
507  if ( ivertex != -1 ) {
508  for ( ; daughter != daughters.end(); ++daughter) {
509  int theDaughterId = mySimEvent->addSimTrack(&(*daughter), ivertex);
510  // Find the closest charged daughter (if charged mother)
511  if ( PP.particle().charge() * daughter->charge() > 1E-10 ) {
512  double dist = (daughter->Vect().Unit().Cross(PP.particle().Vect().Unit())).R();
513  if ( dist < distCut && dist < distMin ) {
514  distMin = dist;
515  theClosestChargedDaughterId = theDaughterId;
516  }
517  }
518  }
519  }
520  // Attach mother and closest daughter sp as to cheat tracking ;-)
521  if ( theClosestChargedDaughterId >=0 )
522  mySimEvent->track(fsimi).setClosestDaughterId(theClosestChargedDaughterId);
523  }
524 
525  }
526 
527 }
528 
529 
530 void
532 
533  //
534  for ( unsigned idaugh=0; idaugh < (unsigned)(mySimEvent->track(fsimi).nDaughters()); ++idaugh) {
535  // Initial momentum of the daughter
536  XYZTLorentzVector daughMomentum (mySimEvent->track(fsimi).daughter(idaugh).momentum());
537  // Rotate and rescale
538  XYZVector newMomentum (r * daughMomentum.Vect());
539  newMomentum *= rescale;
540  double newEnergy = std::sqrt(newMomentum.mag2() + daughMomentum.mag2());
541  // Set the new momentum
542  mySimEvent->track(fsimi).setMomentum(XYZTLorentzVector(newMomentum.X(),newMomentum.Y(),newMomentum.Z(),newEnergy));
543  // Watch out : recursive call to get all grand-daughters
544  int fsimDaug = mySimEvent->track(fsimi).daughter(idaugh).id();
545  moveAllDaughters(fsimDaug,r,rescale);
546  }
547 }
548 
549 void
551  const ParticlePropagator& PP,
552  std::map<double,PSimHit>& theHitMap,
553  int trackID, int partID, const TrackerTopology *tTopo) {
554 
555  // Propagate the particle coordinates to the closest tracker detector(s)
556  // in this layer and create the PSimHit(s)
557 
558  // const MagneticField& mf = MagneticFieldMap::instance()->magneticField();
559  // This solution is actually much faster !
561  AnalyticalPropagator alongProp(&mf, anyDirection);
563 
564  // std::cout << "PP.X() = " << PP.X() << std::endl;
565  // std::cout << "PP.Y() = " << PP.Y() << std::endl;
566  // std::cout << "PP.Z() = " << PP.Z() << std::endl;
567 
569  const DetLayer* tkLayer = detLayer(layer,PP.particle().Z());
570 
571  TrajectoryStateOnSurface trajState = makeTrajectoryState( tkLayer, PP, &mf);
572  float thickness = theMaterialEffects ? theMaterialEffects->thickness() : 0.;
573  float eloss = theMaterialEffects ? theMaterialEffects->energyLoss() : 0.;
574 
575  // Find, in the corresponding layers, the detectors compatible
576  // with the current track
577  std::vector<DetWithState> compat
578  = tkLayer->compatibleDets( trajState, alongProp, est);
579 
580  // And create the corresponding PSimHits
581  std::map<double,PSimHit> theTrackHits;
582  for (std::vector<DetWithState>::const_iterator i=compat.begin(); i!=compat.end(); i++) {
583  // Correct Eloss for last 3 rings of TEC (thick sensors, 0.05 cm)
584  // Disgusting fudge factor !
585  makePSimHits( i->first, i->second, theHitMap, trackID, eloss, thickness, partID,tTopo);
586  }
587 
588 }
589 
592  const ParticlePropagator& pp,
593  const MagneticField* field) const
594 {
595  GlobalPoint pos( pp.particle().X(), pp.particle().Y(), pp.particle().Z());
596  GlobalVector mom( pp.particle().Px(), pp.particle().Py(), pp.particle().Pz());
597  auto plane = layer->surface().tangentPlane(pos);
599  (GlobalTrajectoryParameters( pos, mom, TrackCharge( pp.particle().charge()), field), *plane);
600 }
601 
602 void
604  const TrajectoryStateOnSurface& ts,
605  std::map<double,PSimHit>& theHitMap,
606  int tkID, float el, float thick, int pID,
607  const TrackerTopology *tTopo)
608 {
609 
610  std::vector< const GeomDet*> comp = det->components();
611  if (!comp.empty()) {
612  for (std::vector< const GeomDet*>::const_iterator i = comp.begin();
613  i != comp.end(); i++) {
614  auto du = (*i);
615  if (du->isLeaf()) // not even needed (or it should iterate if really not leaf)
616  theHitMap.insert(theHitMap.end(),makeSinglePSimHit( *du, ts, tkID, el, thick, pID,tTopo));
617  }
618  }
619  else {
620  auto du = (det);
621  theHitMap.insert(theHitMap.end(),makeSinglePSimHit( *du, ts, tkID, el, thick, pID,tTopo));
622  }
623 
624 }
625 
626 std::pair<double,PSimHit>
628  const TrajectoryStateOnSurface& ts,
629  int tkID, float el, float thick, int pID,
630  const TrackerTopology *tTopo) const
631 {
632 
633  const float onSurfaceTolarance = 0.01; // 10 microns
634 
635  LocalPoint lpos;
636  LocalVector lmom;
637  if ( fabs( det.toLocal(ts.globalPosition()).z()) < onSurfaceTolarance) {
638  lpos = ts.localPosition();
639  lmom = ts.localMomentum();
640  }
641  else {
644  ts.transverseCurvature(),
645  anyDirection);
646  std::pair<bool,double> path = crossing.pathLength(det.surface());
647  if (!path.first) {
648  // edm::LogWarning("FastTracking") << "TrajectoryManager ERROR: crossing with det failed, skipping PSimHit";
649  return std::pair<double,PSimHit>(0.,PSimHit());
650  }
651  lpos = det.toLocal( GlobalPoint( crossing.position(path.second)));
652  lmom = det.toLocal( GlobalVector( crossing.direction(path.second)));
653  lmom = lmom.unit() * ts.localMomentum().mag();
654  }
655 
656  // The module (half) thickness
657  const BoundPlane& theDetPlane = det.surface();
658  float halfThick = 0.5*theDetPlane.bounds().thickness();
659  // The Energy loss rescaled to the module thickness
660  float eloss = el;
661  if ( thick > 0. ) {
662  // Total thickness is in radiation lengths, 1 radlen = 9.36 cm
663  // Sensitive module thickness is about 30 microns larger than
664  // the module thickness itself
665  eloss *= (2.* halfThick - 0.003) / (9.36 * thick);
666  }
667  // The entry and exit points, and the time of flight
668  float pZ = lmom.z();
669  LocalPoint entry = lpos + (-halfThick/pZ) * lmom;
670  LocalPoint exit = lpos + halfThick/pZ * lmom;
671  float tof = ts.globalPosition().mag() / 30. ; // in nanoseconds, FIXME: very approximate
672 
673  // If a hadron suffered a nuclear interaction, just assign the hits of the closest
674  // daughter to the mother's track. The same applies to a charged particle decay into
675  // another charged particle.
676  int localTkID = tkID;
677  if ( !mySimEvent->track(tkID).noMother() && mySimEvent->track(tkID).mother().closestDaughterId() == tkID )
678  localTkID = mySimEvent->track(tkID).mother().id();
679 
680  // FIXME: fix the track ID and the particle ID
681  PSimHit hit( entry, exit, lmom.mag(), tof, eloss, pID,
682  det.geographicalId().rawId(), localTkID,
683  lmom.theta(),
684  lmom.phi());
685 
686  // Check that the PSimHit is physically on the module!
687  unsigned subdet = DetId(hit.detUnitId()).subdetId();
688  double boundX = theDetPlane.bounds().width()/2.;
689  double boundY = theDetPlane.bounds().length()/2.;
690 
691  // Special treatment for TID and TEC trapeziodal modules
692  if ( subdet == 4 || subdet == 6 )
693  boundX *= 1. - hit.localPosition().y()/theDetPlane.position().perp();
694 
695 #ifdef FAMOS_DEBUG
696  unsigned detid = DetId(hit.detUnitId()).rawId();
697  unsigned stereo = 0;
698  unsigned theLayer = 0;
699  unsigned theRing = 0;
700  switch (subdet) {
701  case 1:
702  {
703 
704  theLayer = tTopo->pxbLayer(detid);
705  std::cout << "\tPixel Barrel Layer " << theLayer << std::endl;
706  stereo = 1;
707  break;
708  }
709  case 2:
710  {
711 
712  theLayer = tTopo->pxfDisk(detid);
713  std::cout << "\tPixel Forward Disk " << theLayer << std::endl;
714  stereo = 1;
715  break;
716  }
717  case 3:
718  {
719 
720  theLayer = tTopo->tibLayer(detid);
721  std::cout << "\tTIB Layer " << theLayer << std::endl;
722  stereo = module.stereo();
723  break;
724  }
725  case 4:
726  {
727 
728  theLayer = tTopo->tidWheel(detid);
729  theRing = tTopo->tidRing(detid);
730  unsigned int theSide = module.side();
731  if ( theSide == 1 )
732  std::cout << "\tTID Petal Back " << std::endl;
733  else
734  std::cout << "\tTID Petal Front" << std::endl;
735  std::cout << "\tTID Layer " << theLayer << std::endl;
736  std::cout << "\tTID Ring " << theRing << std::endl;
737  stereo = module.stereo();
738  break;
739  }
740  case 5:
741  {
742 
743  theLayer = tTopo->tobLayer(detid);
744  stereo = tTopo->tobStereo(detid);
745  std::cout << "\tTOB Layer " << theLayer << std::endl;
746  break;
747  }
748  case 6:
749  {
750 
751  theLayer = tTopo->tecWheel(detid);
752  theRing = tTopo->tecRing(detid);
753  unsigned int theSide = module.petal()[0];
754  if ( theSide == 1 )
755  std::cout << "\tTEC Petal Back " << std::endl;
756  else
757  std::cout << "\tTEC Petal Front" << std::endl;
758  std::cout << "\tTEC Layer " << theLayer << std::endl;
759  std::cout << "\tTEC Ring " << theRing << std::endl;
760  stereo = module.stereo();
761  break;
762  }
763  default:
764  {
765  stereo = 0;
766  break;
767  }
768  }
769 
770  std::cout << "Thickness = " << 2.*halfThick-0.003 << "; " << thick * 9.36 << std::endl
771  << "Length = " << det.surface().bounds().length() << std::endl
772  << "Width = " << det.surface().bounds().width() << std::endl;
773 
774  std::cout << "Hit position = "
775  << hit.localPosition().x() << " "
776  << hit.localPosition().y() << " "
777  << hit.localPosition().z() << std::endl;
778 #endif
779 
780  // Check if the hit is on the physical volume of the module
781  // (It happens that it is not, in the case of double sided modules,
782  // because the envelope of the gluedDet is larger than each of
783  // the mono and the stereo modules)
784 
785  double dist = 0.;
786  GlobalPoint IP (mySimEvent->track(localTkID).vertex().position().x(),
787  mySimEvent->track(localTkID).vertex().position().y(),
788  mySimEvent->track(localTkID).vertex().position().z());
789 
790  dist = ( fabs(hit.localPosition().x()) > boundX ||
791  fabs(hit.localPosition().y()) > boundY ) ?
792  // Will be used later as a flag to reject the PSimHit!
793  -( det.surface().toGlobal(hit.localPosition()) - IP ).mag2()
794  :
795  // These hits are kept!
796  ( det.surface().toGlobal(hit.localPosition()) - IP ).mag2();
797 
798  // Fill Histos (~poor man event display)
799  /*
800  GlobalPoint gpos( det.toGlobal(hit.localPosition()));
801  // std::cout << "gpos.x() = " << gpos.x() << std::endl;
802  // std::cout << "gpos.y() = " << gpos.y() << std::endl;
803 
804  myHistos->fill("h300",gpos.x(),gpos.y());
805  if ( sin(gpos.phi()) > 0. )
806  myHistos->fill("h301",gpos.z(),gpos.perp());
807  else
808  myHistos->fill("h301",gpos.z(),-gpos.perp());
809  */
810  return std::pair<double,PSimHit>(dist,hit);
811 
812 }
813 
814 void
816 {
817 
818  // These are the BoundSurface&, the BoundDisk* and the BoundCylinder* for that layer
819  // const BoundSurface& theSurface = layer.surface();
820  // BoundDisk* theDisk = layer.disk(); // non zero for endcaps
821  // BoundCylinder* theCylinder = layer.cylinder(); // non zero for barrel
822  // int theLayer = layer.layerNumber(); // 1->3 PixB, 4->5 PixD,
823  // // 6->9 TIB, 10->12 TID,
824  // // 13->18 TOB, 19->27 TEC
825 
828 
829  const std::vector< const BarrelDetLayer*>& barrelLayers =
831  LogDebug("FastTracking") << "Barrel DetLayer dump: ";
832  for (auto bl=barrelLayers.begin();
833  bl != barrelLayers.end(); ++bl) {
834  LogDebug("FastTracking")<< "radius " << (**bl).specificSurface().radius();
835  }
836 
837  const std::vector< const ForwardDetLayer*>& posForwardLayers =
839  LogDebug("FastTracking") << "Positive Forward DetLayer dump: ";
840  for (auto fl=posForwardLayers.begin();
841  fl != posForwardLayers.end(); ++fl) {
842  LogDebug("FastTracking") << "Z pos "
843  << (**fl).surface().position().z()
844  << " radii "
845  << (**fl).specificSurface().innerRadius()
846  << ", "
847  << (**fl).specificSurface().outerRadius();
848  }
849 
850  const float rTolerance = 1.5;
851  const float zTolerance = 3.;
852 
853  LogDebug("FastTracking")<< "Dump of TrackerInteractionGeometry cylinders:";
854  for( std::list<TrackerLayer>::const_iterator i=_theGeometry->cylinderBegin();
855  i!=_theGeometry->cylinderEnd(); ++i) {
856  const BoundCylinder* cyl = i->cylinder();
857  const BoundDisk* disk = i->disk();
858 
859  LogDebug("FastTracking") << "Famos Layer no " << i->layerNumber()
860  << " is sensitive? " << i->sensitive()
861  << " pos " << i->surface().position();
862  if (!i->sensitive()) continue;
863 
864  if (cyl != nullptr) {
865  LogDebug("FastTracking") << " cylinder radius " << cyl->radius();
866  bool found = false;
867  for (auto
868  bl=barrelLayers.begin(); bl != barrelLayers.end(); ++bl) {
869 
870  if (fabs( cyl->radius() - (**bl).specificSurface().radius()) < rTolerance) {
871  theLayerMap[i->layerNumber()] = *bl;
872  found = true;
873  LogDebug("FastTracking")<< "Corresponding DetLayer found with radius "
874  << (**bl).specificSurface().radius();
875  break;
876  }
877  }
878  if (!found) {
879  edm::LogWarning("FastTracking") << " Trajectory manager FAILED to find a corresponding DetLayer!";
880  }
881  }
882  else {
883  LogDebug("FastTracking") << " disk radii " << disk->innerRadius()
884  << ", " << disk->outerRadius();
885  bool found = false;
886  for (auto fl=posForwardLayers.begin();
887  fl != posForwardLayers.end(); ++fl) {
888 
889  if (fabs( disk->position().z() - (**fl).surface().position().z()) < zTolerance) {
890  theLayerMap[i->layerNumber()] = *fl;
891  found = true;
892  LogDebug("FastTracking") << "Corresponding DetLayer found with Z pos "
893  << (**fl).surface().position().z()
894  << " and radii "
895  << (**fl).specificSurface().innerRadius()
896  << ", "
897  << (**fl).specificSurface().outerRadius();
898  break;
899  }
900  }
901  if (!found) {
902  edm::LogWarning("FastTracking") << "FAILED to find a corresponding DetLayer!";
903  }
904  }
905  }
906 
907  // Put the negative layers in the same map but with an offset
908  const std::vector< const ForwardDetLayer*>& negForwardLayers = theGeomSearchTracker->negForwardLayers();
909  for (auto nl=negForwardLayers.begin();
910  nl != negForwardLayers.end(); ++nl) {
911  for (int i=0; i<=theNegLayerOffset; i++) {
912  if (theLayerMap[i] == nullptr) continue;
913  if ( fabs( (**nl).surface().position().z() +theLayerMap[i]-> surface().position().z()) < zTolerance) {
915  break;
916  }
917  }
918  }
919 
920 }
921 
922 const DetLayer*
923 TrajectoryManager::detLayer( const TrackerLayer& layer, float zpos) const
924 {
925  if (zpos > 0 || !layer.forward() ) return theLayerMap[layer.layerNumber()];
926  else return theLayerMap[layer.layerNumber()+theNegLayerOffset];
927 }
928 
929 void
931 {
932 
933  std::map<unsigned,std::map<double,PSimHit> >::const_iterator itrack = thePSimHits.begin();
934  std::map<unsigned,std::map<double,PSimHit> >::const_iterator itrackEnd = thePSimHits.end();
935  for ( ; itrack != itrackEnd; ++itrack ) {
936  std::map<double,PSimHit>::const_iterator it = (itrack->second).begin();
937  std::map<double,PSimHit>::const_iterator itEnd = (itrack->second).end();
938  for( ; it!= itEnd; ++it) {
939  /*
940  DetId theDetUnitId((it->second).detUnitId());
941  const GeomDet* theDet = theGeomTracker->idToDet(theDetUnitId);
942  std::cout << "Track/z/r after : "
943  << (it->second).trackId() << " "
944  << theDet->surface().toGlobal((it->second).localPosition()).z() << " "
945  << theDet->surface().toGlobal((it->second).localPosition()).perp() << std::endl;
946  */
947  // Keep only those hits that are on the physical volume of a module
948  // (The other hits have been assigned a negative <double> value.
949  if ( it->first > 0. ) c.push_back(it->second);
950  }
951  }
952 
953 }
#define LogDebug(id)
void reconstruct(const TrackerTopology *tTopo, RandomEngineAndDistribution const *)
Does the real job.
bool hasDecayed() const
Has the particle decayed while propagated ?
int id() const
the index in FBaseSimEvent
Definition: FSimVertex.h:44
void initializeLayerMap()
Initialize correspondence map between Famos interaction geometry and tracker reco geometry...
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
T getParameter(std::string const &) const
int addSimVertex(const XYZTLorentzVector &decayVertex, int im=-1, FSimVertexType::VertexType type=FSimVertexType::ANY)
Add a new vertex to the Event and to the various lists.
T getUntrackedParameter(std::string const &, T const &) const
const DetLayer * detLayer(const TrackerLayer &layer, float zpos) const
Returns the DetLayer pointer corresponding to the FAMOS layer.
virtual float length() const =0
void initializeTrackerGeometry(const TrackerGeometry *geomTracker)
Initialize the full Tracker Geometry.
void propagateToCalorimeters(ParticlePropagator &PP, int fsimi, RandomEngineAndDistribution const *)
Propagate the particle through the calorimeters.
XYZVector Vect() const
the momentum threevector
Definition: RawParticle.h:342
PythiaDecays * myDecayEngine
std::pair< double, PSimHit > makeSinglePSimHit(const GeomDetUnit &det, const TrajectoryStateOnSurface &ts, int tkID, float el, float thick, int pID, const TrackerTopology *tTopo) const
and there
bool acceptParticle(const RawParticle &p) const
virtual ConstReferenceCountingPointer< TangentPlane > tangentPlane(const GlobalPoint &) const =0
int addSimTrack(const RawParticle *p, int iv, int ig=-1, const HepMC::GenVertex *ev=0)
Add a new track to the Event and to the various lists.
bool propagateToPreshowerLayer1(bool first=true)
void makePSimHits(const GeomDet *det, const TrajectoryStateOnSurface &ts, std::map< double, PSimHit > &theHitMap, int tkID, float el, float thick, int pID, const TrackerTopology *tTopo)
and there
unsigned int tibLayer(const DetId &id) const
double cTau(int pdgID, const HepPDT::ParticleDataTable *pdt)
unsigned int tidRing(const DetId &id) const
const FSimVertex & endVertex() const
end vertex
std::vector< ForwardDetLayer const * > const & posForwardLayers() const
void createPSimHits(const TrackerLayer &layer, const ParticlePropagator &P_before, std::map< double, PSimHit > &theHitMap, int trackID, int partID, const TrackerTopology *tTopo)
Create a vector of PSimHits.
const MagneticFieldMap * _theFieldMap
unsigned int layerNumber() const
Returns the layer number.
Definition: TrackerLayer.h:82
bool forward() const
Is the layer forward ?
Definition: TrackerLayer.h:70
void setPosition(const math::XYZTLorentzVector &newPosition)
Reset the position (to be used with care)
Definition: FSimVertex.h:52
const TrackerGeometry * theGeomTracker
void initializeRecoGeometry(const GeometricSearchTracker *geomSearchTracker, const TrackerInteractionGeometry *interactionGeometry, const MagneticFieldMap *aFieldMap)
Initialize the Reconstruction Geometry.
unsigned int pxfDisk(const DetId &id) const
const XYZTLorentzVector & momentum() const
Temporary (until move of SimTrack to Mathcore) - No! Actually very useful.
Definition: FSimTrack.h:204
TrajectoryStateOnSurface makeTrajectoryState(const DetLayer *layer, const ParticlePropagator &pp, const MagneticField *field) const
Teddy, you must put comments there.
unsigned int tecRing(const DetId &id) const
ring id
const HepPDT::ParticleDataTable * theTable() const
Get the pointer to the particle data table.
Definition: FBaseSimEvent.h:57
#define nullptr
uint32_t tobStereo(const DetId &id) const
void moveAllDaughters(int fsimi, const Rotation &r, double rescale)
Move, rescale and rotate all daughters after propagation, material effects and decay of the mother...
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
const FSimTrack & daughter(int i) const
Ith daughter.
const GeometricSearchTracker * theGeomSearchTracker
void setEcal(const RawParticle &pp, int success)
Set the ecal variables.
Definition: FSimTrack.cc:81
std::list< TrackerLayer >::const_iterator cylinderEnd() const
Returns the last pointer in the cylinder list.
int getSuccess() const
Has propagation been performed and was barrel or endcap reached ?
const DaughterParticleList & particleDaughters(ParticlePropagator &particle, CLHEP::HepRandomEngine *)
Definition: PythiaDecays.cc:35
void setPropagate()
The particle has been propgated through the tracker.
Definition: FSimTrack.cc:61
GlobalPoint globalPosition() const
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:69
bool propagateToBoundSurface(const TrackerLayer &)
int nDaughters() const
Number of daughters.
const Bounds & bounds() const
Definition: Surface.h:120
DaughterParticleList::const_iterator DaughterParticleIterator
Definition: PythiaDecays.h:26
unsigned int tidWheel(const DetId &id) const
void setClosestDaughterId(int id)
Set the index of the closest charged daughter.
Definition: FSimTrack.h:198
void setLayer2(const RawParticle &pp, int success)
Set the preshower layer2 variables.
Definition: FSimTrack.cc:74
TRandom random
Definition: MVATrainer.cc:138
RawParticle const & particle() const
The particle being propagated.
ROOT::Math::AxisAngle Rotation
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
bool notYetToEndVertex(const XYZTLorentzVector &pos) const
Compare the end vertex position with another position.
Definition: FSimTrack.cc:45
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
FSimVertex & vertex(int id) const
Return vertex with given Id.
LocalVector localMomentum() const
virtual float width() const =0
int TrackCharge
Definition: TrackCharge.h:4
bool propagateToVFcalEntrance(bool first=true)
void setLayer1(const RawParticle &pp, int success)
Set the preshower layer1 variables.
Definition: FSimTrack.cc:67
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:340
~TrajectoryManager()
Default Destructor.
T mag() const
Definition: PV3DBase.h:67
const TrackerInteractionGeometry * _theGeometry
const math::XYZTLorentzVector & position() const
Temporary (until CMSSW moves to Mathcore) - No ! Actually very useful.
Definition: FSimVertex.h:49
MaterialEffects * theMaterialEffects
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
std::vector< ForwardDetLayer const * > const & negForwardLayers() const
double Y() const
y of vertex
Definition: RawParticle.h:306
double Py() const
y of the momentum
Definition: RawParticle.h:319
int closestDaughterId() const
Get the index of the closest charged daughter.
Definition: FSimTrack.h:201
double Z() const
z of vertex
Definition: RawParticle.h:307
void setPropagationConditions(const TrackerLayer &, bool firstLoop=true)
virtual const BoundSurface & surface() const =0
The surface of the GeometricSearchDet.
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
void save()
Save nuclear interaction information.
CLHEP::HepRandomEngine & theEngine() const
#define end
Definition: vmac.h:39
double Pz() const
z of the momentum
Definition: RawParticle.h:322
std::pair< const GeomDet *, TrajectoryStateOnSurface > DetWithState
double charge() const
get the MEASURED charge
Definition: RawParticle.h:313
TrajectoryManager()
Default Constructor.
bool propagateToEcalEntrance(bool first=true)
unsigned int nTracks() const
Number of tracks.
Definition: FSimEvent.cc:36
std::map< unsigned, std::map< double, PSimHit > > thePSimHits
unsigned int pxbLayer(const DetId &id) const
virtual std::vector< const GeomDet * > components() const
Returns direct components, if any.
Definition: GeomDet.h:88
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:339
double energyLoss() const
Return the energy loss by ionization in the current layer.
void interact(FSimEvent &simEvent, const TrackerLayer &layer, ParticlePropagator &PP, unsigned i, RandomEngineAndDistribution const *)
Vector3DBase unit() const
Definition: Vector3DBase.h:57
Definition: DetId.h:18
double thickness() const
Return the thickness of the current layer.
void loadSimHits(edm::PSimHitContainer &c) const
std::list< TrackerLayer >::const_iterator cylinderBegin() const
Returns the first pointer in the cylinder list.
void setTkPosition(const math::XYZVectorD &pos)
Definition: SimTrack.h:41
const TrackerInteractionGeometry * theGeometry()
Returns the pointer to geometry.
double X() const
x of vertex
Definition: RawParticle.h:305
double getMagneticField() const
Get the magnetic field.
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:22
std::vector< const DetLayer * > theLayerMap
bool propagateToHcalEntrance(bool first=true)
#define begin
Definition: vmac.h:32
GlobalVector globalMomentum() const
const KineParticleFilter & filter() const
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:91
double Px() const
x of the momentum
Definition: RawParticle.h:316
void setTkMomentum(const math::XYZTLorentzVectorD &mom)
Definition: SimTrack.h:43
void setVFcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:95
std::pair< const GeomDet *, TrajectoryStateOnSurface > DetWithState
bool noMother() const
no mother particle
bool propagateToLayer(ParticlePropagator &PP, unsigned layer)
std::vector< PSimHit > PSimHitContainer
void setHcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:88
const FSimVertex vertex() const
Origin vertex.
bool onFiducial() const
Is the vertex on some material ?
void updateWithDaughters(ParticlePropagator &PP, int fsimi, RandomEngineAndDistribution const *)
Decay the particle and update the SimEvent with daughters.
math::XYZVector XYZVector
Definition: RawParticle.h:28
void setMomentum(const math::XYZTLorentzVector &newMomentum)
Reset the momentum (to be used with care)
Definition: FSimTrack.h:207
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:56
unsigned int tecWheel(const DetId &id) const
Definition: vlib.h:208
bool propagateToPreshowerLayer2(bool first=true)
const FSimTrack & mother() const
mother
unsigned int tobLayer(const DetId &id) const
Global3DVector GlobalVector
Definition: GlobalVector.h:10
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:27
std::vector< BarrelDetLayer const * > const & barrelLayers() const
std::vector< RawParticle > DaughterParticleList
Definition: PythiaDecays.h:25
FSimTrack & track(int id) const
Return track with given Id.
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11