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