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