CMS 3D CMS Logo

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