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