CMS 3D CMS Logo

FastSimProducer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 #include <string>
4 
5 // framework
17 
18 // data formats
25 
26 // fastsim
38 
39 // Hack for calorimetry
47 
49 // Author: L. Vanelderen, S. Kurz
50 // Date: 29 May 2017
52 
54 
65 public:
66  explicit FastSimProducer(const edm::ParameterSet&);
67  ~FastSimProducer() override { ; }
68 
69 private:
70  void beginStream(edm::StreamID id) override;
71  void produce(edm::Event&, const edm::EventSetup&) override;
72  void endStream() override;
74  fastsim::ParticleManager* particleManager,
75  HepPDT::ParticleDataTable const& particleTable);
76 
80  double beamPipeRadius_;
83  std::unique_ptr<RandomEngineAndDistribution> _randomEngine;
84 
88  std::unique_ptr<CalorimetryManager> myCalorimetry; // unfortunately, default constructor cannot be called
90 
93 
95  std::vector<std::unique_ptr<fastsim::InteractionModel> > interactionModels_;
96  std::map<std::string, fastsim::InteractionModel*> interactionModelMap_;
98 };
99 
100 const std::string FastSimProducer::MESSAGECATEGORY = "FastSimulation";
101 
103  : genParticlesToken_(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("src"))),
104  geometry_(iConfig.getParameter<edm::ParameterSet>("trackerDefinition")),
105  caloGeometry_(iConfig.getParameter<edm::ParameterSet>("caloDefinition")),
106  beamPipeRadius_(iConfig.getParameter<double>("beamPipeRadius")),
107  deltaRchargedMother_(iConfig.getParameter<double>("deltaRchargedMother")),
108  particleFilter_(iConfig.getParameter<edm::ParameterSet>("particleFilter")),
110  simulateCalorimetry(iConfig.getParameter<bool>("simulateCalorimetry")),
111  simulateMuons(iConfig.getParameter<bool>("simulateMuons")),
112  fixLongLivedBug_(iConfig.getParameter<bool>("fixLongLivedBug")),
113  useFastSimsDecayer_(iConfig.getParameter<bool>("useFastSimsDecayer")) {
114  //----------------
115  // define interaction models
116  //---------------
117 
118  const edm::ParameterSet& modelCfgs = iConfig.getParameter<edm::ParameterSet>("interactionModels");
119  for (const std::string& modelName : modelCfgs.getParameterNames()) {
120  const edm::ParameterSet& modelCfg = modelCfgs.getParameter<edm::ParameterSet>(modelName);
121  std::string modelClassName(modelCfg.getParameter<std::string>("className"));
122  // Use plugin-factory to create model
123  std::unique_ptr<fastsim::InteractionModel> interactionModel(
124  fastsim::InteractionModelFactory::get()->create(modelClassName, modelName, modelCfg));
125  if (!interactionModel.get()) {
126  throw cms::Exception("FastSimProducer") << "InteractionModel " << modelName << " could not be created";
127  }
128  // Add model to list
129  interactionModels_.push_back(std::move(interactionModel));
130  // and create the map
132  }
133 
134  //----------------
135  // calorimetry
136  //---------------
137 
138  if (simulateCalorimetry) {
139  myCalorimetry.reset(new CalorimetryManager(nullptr,
140  iConfig.getParameter<edm::ParameterSet>("Calorimetry"),
141  iConfig.getParameter<edm::ParameterSet>("MaterialEffectsForMuonsInECAL"),
142  iConfig.getParameter<edm::ParameterSet>("MaterialEffectsForMuonsInHCAL"),
143  iConfig.getParameter<edm::ParameterSet>("GFlash")));
144  }
145 
146  //----------------
147  // register products
148  //----------------
149 
150  // SimTracks and SimVertices
151  produces<edm::SimTrackContainer>();
152  produces<edm::SimVertexContainer>();
153  // products of interaction models, i.e. simHits
154  for (auto& interactionModel : interactionModels_) {
155  interactionModel->registerProducts(*this);
156  }
157  produces<edm::PCaloHitContainer>("EcalHitsEB");
158  produces<edm::PCaloHitContainer>("EcalHitsEE");
159  produces<edm::PCaloHitContainer>("EcalHitsES");
160  produces<edm::PCaloHitContainer>("HcalHits");
161  produces<edm::SimTrackContainer>("MuonSimTracks");
162 }
163 
165  _randomEngine = std::make_unique<RandomEngineAndDistribution>(id);
166 }
167 
169  LogDebug(MESSAGECATEGORY) << " produce";
170 
173 
174  // Define containers for SimTracks, SimVertices
175  std::unique_ptr<edm::SimTrackContainer> simTracks_(new edm::SimTrackContainer);
176  std::unique_ptr<edm::SimVertexContainer> simVertices_(new edm::SimVertexContainer);
177 
178  // Get the particle data table (in case lifetime or charge of GenParticles not set)
180  iSetup.getData(pdt);
181 
182  // Get the GenParticle collection
184  iEvent.getByToken(genParticlesToken_, genParticles);
185 
186  // Load the ParticleManager which returns the particles that have to be propagated
187  // Creates a fastsim::Particle out of a GenParticle/secondary
188  fastsim::ParticleManager particleManager(*genParticles->GetEvent(),
189  *pdt,
193  *simTracks_,
194  *simVertices_,
197 
198  // Initialize the calorimeter geometry
199  if (simulateCalorimetry) {
200  if (watchCaloGeometry_.check(iSetup) || watchCaloTopology_.check(iSetup)) {
202  iSetup.get<CaloGeometryRecord>().get(pG);
203  myCalorimetry->getCalorimeter()->setupGeometry(*pG);
204 
205  edm::ESHandle<CaloTopology> theCaloTopology;
206  iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
207  myCalorimetry->getCalorimeter()->setupTopology(*theCaloTopology);
208  myCalorimetry->getCalorimeter()->initialize(geometry_.getMagneticFieldZ(math::XYZTLorentzVector(0., 0., 0., 0.)));
209 
210  myCalorimetry->getHFShowerLibrary()->initHFShowerLibrary(iSetup);
211  }
212 
213  // Important: this also cleans the calorimetry information from the last event
214  myCalorimetry->initialize(_randomEngine.get());
215  }
216 
217  // The vector of SimTracks needed for the CalorimetryManager
218  std::vector<FSimTrack> myFSimTracks;
219 
220  LogDebug(MESSAGECATEGORY) << "################################"
221  << "\n###############################";
222 
223  // loop over particles
224  for (std::unique_ptr<fastsim::Particle> particle = particleManager.nextParticle(*_randomEngine); particle != nullptr;
225  particle = particleManager.nextParticle(*_randomEngine)) {
226  LogDebug(MESSAGECATEGORY) << "\n moving NEXT particle: " << *particle;
227 
228  // -----------------------------
229  // This condition is necessary because of hack for calorimetry
230  // -> The CalorimetryManager should also be implemented based on this new FastSim classes (Particle.h) in a future project.
231  // A second loop (below) loops over all parts of the calorimetry in order to create a track of the old FastSim class FSimTrack.
232  // The condition below (R<128, z<302) makes sure that the particle geometrically is outside the tracker boundaries
233  // -----------------------------
234 
235  if (particle->position().Perp2() < 128. * 128. && std::abs(particle->position().Z()) < 302.) {
236  // move the particle through the layers
237  fastsim::LayerNavigator layerNavigator(geometry_);
238  const fastsim::SimplifiedGeometry* layer = nullptr;
239 
240  // moveParticleToNextLayer(..) returns 0 in case that particle decays
241  // in this case particle is propagated up to its decay vertex
242  while (layerNavigator.moveParticleToNextLayer(*particle, layer)) {
243  LogDebug(MESSAGECATEGORY) << " moved to next layer: " << *layer;
244  LogDebug(MESSAGECATEGORY) << " new state: " << *particle;
245 
246  // Hack to interface "old" calo to "new" tracking
247  // Particle reached calorimetry so stop further propagation
249  layer = nullptr;
250  // particle no longer is on a layer
251  particle->resetOnLayer();
252  break;
253  }
254 
255  // break after 25 ns: only happens for particles stuck in loops
256  if (particle->position().T() > 25) {
257  layer = nullptr;
258  // particle no longer is on a layer
259  particle->resetOnLayer();
260  break;
261  }
262 
263  // perform interaction between layer and particle
264  // do only if there is actual material
265  if (layer->getThickness(particle->position(), particle->momentum()) > 1E-10) {
266  int nSecondaries = 0;
267  // loop on interaction models
268  for (fastsim::InteractionModel* interactionModel : layer->getInteractionModels()) {
269  LogDebug(MESSAGECATEGORY) << " interact with " << *interactionModel;
270  std::vector<std::unique_ptr<fastsim::Particle> > secondaries;
271  interactionModel->interact(*particle, *layer, secondaries, *_randomEngine);
272  nSecondaries += secondaries.size();
273  particleManager.addSecondaries(particle->position(), particle->simTrackIndex(), secondaries, layer);
274  }
275 
276  // kinematic cuts: particle might e.g. lost all its energy
277  if (!particleFilter_.acceptsEn(*particle)) {
278  // Add endvertex if particle did not create any secondaries
279  if (nSecondaries == 0)
280  particleManager.addEndVertex(particle.get());
281  layer = nullptr;
282  break;
283  }
284  }
285 
286  LogDebug(MESSAGECATEGORY) << "--------------------------------"
287  << "\n-------------------------------";
288  }
289 
290  // do decays
291  if (!particle->isStable() && particle->remainingProperLifeTimeC() < 1E-10) {
292  LogDebug(MESSAGECATEGORY) << "Decaying particle...";
293  std::vector<std::unique_ptr<fastsim::Particle> > secondaries;
294  if (useFastSimsDecayer_) decayer_.decay(*particle, secondaries, _randomEngine->theEngine());
295  LogDebug(MESSAGECATEGORY) << " decay has " << secondaries.size() << " products";
296  particleManager.addSecondaries(particle->position(), particle->simTrackIndex(), secondaries);
297  continue;
298  }
299 
300  LogDebug(MESSAGECATEGORY) << "################################"
301  << "\n###############################";
302  }
303 
304  // -----------------------------
305  // Hack to interface "old" calorimetry with "new" propagation in tracker
306  // The CalorimetryManager has to know which particle could in principle hit which parts of the calorimeter
307  // I think it's a bit strange to propagate the particle even further (and even decay it) if it already hits
308  // some part of the calorimetry but this is how the code works...
309  // -----------------------------
310 
311  if (particle->position().Perp2() >= 128. * 128. || std::abs(particle->position().Z()) >= 302.) {
312  LogDebug(MESSAGECATEGORY) << "\n moving particle to calorimetry: " << *particle;
313 
314  // create FSimTrack (this is the object the old propagation uses)
315  myFSimTracks.push_back(createFSimTrack(particle.get(), &particleManager, *pdt));
316  // particle was decayed
317  if (!particle->isStable() && particle->remainingProperLifeTimeC() < 1E-10) {
318  continue;
319  }
320 
321  LogDebug(MESSAGECATEGORY) << "################################"
322  << "\n###############################";
323  }
324 
325  // -----------------------------
326  // End Hack
327  // -----------------------------
328 
329  LogDebug(MESSAGECATEGORY) << "################################"
330  << "\n###############################";
331  }
332 
333  // store simTracks and simVertices
334  iEvent.put(std::move(simTracks_));
335  iEvent.put(std::move(simVertices_));
336  // store products of interaction models, i.e. simHits
337  for (auto& interactionModel : interactionModels_) {
338  interactionModel->storeProducts(iEvent);
339  }
340 
341  // -----------------------------
342  // Calorimetry Manager
343  // -----------------------------
344  if (simulateCalorimetry) {
345  for (auto myFSimTrack : myFSimTracks) {
346  myCalorimetry->reconstructTrack(myFSimTrack, _randomEngine.get());
347  }
348  }
349 
350  // -----------------------------
351  // Store Hits
352  // -----------------------------
353  std::unique_ptr<edm::PCaloHitContainer> p4(new edm::PCaloHitContainer);
354  std::unique_ptr<edm::PCaloHitContainer> p5(new edm::PCaloHitContainer);
355  std::unique_ptr<edm::PCaloHitContainer> p6(new edm::PCaloHitContainer);
356  std::unique_ptr<edm::PCaloHitContainer> p7(new edm::PCaloHitContainer);
357 
358  std::unique_ptr<edm::SimTrackContainer> m1(new edm::SimTrackContainer);
359 
360  if (simulateCalorimetry) {
361  myCalorimetry->loadFromEcalBarrel(*p4);
362  myCalorimetry->loadFromEcalEndcap(*p5);
363  myCalorimetry->loadFromPreshower(*p6);
364  myCalorimetry->loadFromHcal(*p7);
365  if (simulateMuons) {
366  myCalorimetry->harvestMuonSimTracks(*m1);
367  }
368  }
369  iEvent.put(std::move(p4), "EcalHitsEB");
370  iEvent.put(std::move(p5), "EcalHitsEE");
371  iEvent.put(std::move(p6), "EcalHitsES");
372  iEvent.put(std::move(p7), "HcalHits");
373  iEvent.put(std::move(m1), "MuonSimTracks");
374 }
375 
377 
379  fastsim::ParticleManager* particleManager,
380  HepPDT::ParticleDataTable const& particleTable) {
381  FSimTrack myFSimTrack(particle->pdgId(),
382  particleManager->getSimTrack(particle->simTrackIndex()).momentum(),
383  particle->simVertexIndex(),
384  particle->genParticleIndex(),
385  particle->simTrackIndex(),
386  particle->charge(),
387  particle->position(),
388  particle->momentum(),
389  particleManager->getSimVertex(particle->simVertexIndex()));
390 
391  // move the particle through the caloLayers
392  fastsim::LayerNavigator caloLayerNavigator(caloGeometry_);
393  const fastsim::SimplifiedGeometry* caloLayer = nullptr;
394 
395  // moveParticleToNextLayer(..) returns 0 in case that particle decays
396  // in this case particle is propagated up to its decay vertex
397  while (caloLayerNavigator.moveParticleToNextLayer(*particle, caloLayer)) {
398  LogDebug(MESSAGECATEGORY) << " moved to next caloLayer: " << *caloLayer;
399  LogDebug(MESSAGECATEGORY) << " new state: " << *particle;
400 
401  // break after 25 ns: only happens for particles stuck in loops
402  if (particle->position().T() > 50) {
403  caloLayer = nullptr;
404  break;
405  }
406 
408  // Define ParticlePropagators (RawParticle) needed for CalorimetryManager and save them
410 
411  RawParticle PP = makeParticle(&particleTable, particle->pdgId(), particle->momentum(), particle->position());
412 
413  // no material
414  if (caloLayer->getThickness(particle->position(), particle->momentum()) < 1E-10) {
415  // unfortunately needed for CalorimetryManager
416  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::ECAL) {
417  if (!myFSimTrack.onEcal()) {
418  myFSimTrack.setEcal(PP, 0);
419  }
420  } else if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::HCAL) {
421  if (!myFSimTrack.onHcal()) {
422  myFSimTrack.setHcal(PP, 0);
423  }
424  } else if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::VFCAL) {
425  if (!myFSimTrack.onVFcal()) {
426  myFSimTrack.setVFcal(PP, 0);
427  }
428  }
429 
430  // not necessary to continue propagation
431  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::VFCAL) {
432  myFSimTrack.setGlobal();
433  caloLayer = nullptr;
434  break;
435  }
436 
437  continue;
438  }
439 
440  // Stupid variable used by the old propagator
441  // For details check BaseParticlePropagator.h
442  int success = 0;
443  if (caloLayer->isForward()) {
444  success = 2;
445  // particle moves inwards
446  if (particle->position().Z() * particle->momentum().Z() < 0) {
447  success *= -1;
448  }
449  } else {
450  success = 1;
451  // particle moves inwards
452  if (particle->momentum().X() * particle->position().X() + particle->momentum().Y() * particle->position().Y() <
453  0) {
454  success *= -1;
455  }
456  }
457 
458  // Save the hit
460  if (!myFSimTrack.onLayer1()) {
461  myFSimTrack.setLayer1(PP, success);
462  }
463  }
464 
466  if (!myFSimTrack.onLayer2()) {
467  myFSimTrack.setLayer2(PP, success);
468  }
469  }
470 
471  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::ECAL) {
472  if (!myFSimTrack.onEcal()) {
473  myFSimTrack.setEcal(PP, success);
474  }
475  }
476 
477  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::HCAL) {
478  if (!myFSimTrack.onHcal()) {
479  myFSimTrack.setHcal(PP, success);
480  }
481  }
482 
483  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::VFCAL) {
484  if (!myFSimTrack.onVFcal()) {
485  myFSimTrack.setVFcal(PP, success);
486  }
487  }
488 
489  // Particle reached end of detector
490  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::VFCAL) {
491  myFSimTrack.setGlobal();
492  caloLayer = nullptr;
493  break;
494  }
495 
496  LogDebug(MESSAGECATEGORY) << "--------------------------------"
497  << "\n-------------------------------";
498  }
499 
500  // do decays
501  // don't have to worry about daughters if particle already within the calorimetry
502  // since they will be rejected by the vertex cut of the ParticleFilter
503  if (!particle->isStable() && particle->remainingProperLifeTimeC() < 1E-10) {
504  LogDebug(MESSAGECATEGORY) << "Decaying particle...";
505  std::vector<std::unique_ptr<fastsim::Particle> > secondaries;
506  if (useFastSimsDecayer_) decayer_.decay(*particle, secondaries, _randomEngine->theEngine());
507  LogDebug(MESSAGECATEGORY) << " decay has " << secondaries.size() << " products";
508  particleManager->addSecondaries(particle->position(), particle->simTrackIndex(), secondaries);
509  }
510 
511  return myFSimTrack;
512 }
513 
#define LogDebug(id)
void beginStream(edm::StreamID id) override
std::vector< std::unique_ptr< fastsim::InteractionModel > > interactionModels_
All defined interaction models.
edm::ESWatcher< CaloGeometryRecord > watchCaloGeometry_
T getParameter(std::string const &) const
void update(const edm::EventSetup &iSetup, const std::map< std::string, InteractionModel * > &interactionModelMap)
Initializes the tracker geometry.
Definition: Geometry.cc:44
std::vector< PCaloHit > PCaloHitContainer
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
fastsim::Geometry caloGeometry_
Hack to interface "old" calo to "new" tracking.
std::unique_ptr< RandomEngineAndDistribution > _randomEngine
The random engine.
The core class of the new SimplifiedGeometryPropagator.
Implementation of a generic detector layer (base class for forward/barrel layers).
const math::XYZTLorentzVector & position() const
Return position of the particle.
Definition: Particle.h:142
int onLayer2() const
Definition: FSimTrack.h:101
Handles/tracks (possible) intersections of particle&#39;s trajectory and tracker layers.
Manages GenParticles and Secondaries from interactions.
double deltaRchargedMother_
Cut on deltaR for ClosestChargedDaughter algorithm (FastSim tracking)
def create(alignables, pedeDump, additionalData, outputFile, config)
HepPDT::ParticleDataTable ParticleDataTable
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
fastsim::Geometry geometry_
The definition of the tracker according to python config.
#define nullptr
double remainingProperLifeTimeC() const
Return the particle&#39;s remaining proper lifetime[in ct].
Definition: Particle.h:152
int genParticleIndex() const
Return index of the particle in the genParticle vector.
Definition: Particle.h:167
edm::EDGetTokenT< edm::HepMCProduct > genParticlesToken_
Token to get the genParticles.
void setEcal(const RawParticle &pp, int success)
Set the ecal variables.
Definition: FSimTrack.cc:81
fastsim::Decayer decayer_
Handles decays of non-stable particles using pythia.
void setGlobal()
particle did not decay before more detectors (useful for newProducer)
Definition: FSimTrack.h:156
virtual bool isForward() const =0
Returns false/true depending if the object is a (non-abstract) barrel/forward layer.
virtual FSimTrack createFSimTrack(fastsim::Particle *particle, fastsim::ParticleManager *particleManager, HepPDT::ParticleDataTable const &particleTable)
void setLayer2(const RawParticle &pp, int success)
Set the preshower layer2 variables.
Definition: FSimTrack.cc:74
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:136
bool moveParticleToNextLayer(Particle &particle, const SimplifiedGeometry *&layer)
Move particle along its trajectory to the next intersection with any of the tracker layers...
const SimVertex getSimVertex(unsigned i)
Returns the position of a given SimVertex. Needed for interfacing the code with the old calorimetry...
int onEcal() const
Definition: FSimTrack.h:106
const SimTrack getSimTrack(unsigned i)
Returns a given SimTrack. Needed for interfacing the code with the old calorimetry.
RawParticle makeParticle(HepPDT::ParticleDataTable const *, int id, const math::XYZTLorentzVector &p)
Definition: makeParticle.cc:29
bool getData(T &iHolder) const
Definition: EventSetup.h:111
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Base class for any interaction model between a particle and a tracker layer.
int simVertexIndex() const
Return index of the origin vertex.
Definition: Particle.h:161
~FastSimProducer() override
int iEvent
Definition: GenABIO.cc:224
void setLayer1(const RawParticle &pp, int success)
Set the preshower layer1 variables.
Definition: FSimTrack.cc:67
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const std::vector< InteractionModel * > & getInteractionModels() const
Return the vector of all interaction models that are assigned with a layer.
CaloType getCaloType() const
Hack to interface "old" Calorimetry with "new" Tracker.
double p4[4]
Definition: TauolaWrapper.h:92
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int onVFcal() const
Definition: FSimTrack.h:116
void addSecondaries(const math::XYZTLorentzVector &vertexPosition, int motherSimTrackId, std::vector< std::unique_ptr< Particle > > &secondaries, const SimplifiedGeometry *layer=nullptr)
Adds secondaries that are produced by any of the interactions (or particle decay) to the buffer...
double beamPipeRadius_
The radius of the beampipe.
std::vector< std::string > getParameterNames() const
void decay(const Particle &particle, std::vector< std::unique_ptr< Particle > > &secondaries, CLHEP::HepRandomEngine &engine) const
Decay particle using pythia.
Definition: Decayer.cc:29
int simTrackIndex() const
Return index of the SimTrack.
Definition: Particle.h:155
void produce(edm::Event &, const edm::EventSetup &) override
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
edm::ESWatcher< CaloTopologyRecord > watchCaloTopology_
double charge() const
Return charge of the particle.
Definition: Particle.h:139
std::vector< SimVertex > SimVertexContainer
FastSimProducer(const edm::ParameterSet &)
virtual const double getThickness(const math::XYZTLorentzVector &position) const =0
Return thickness of the layer at a given position.
Definition the tracker geometry (vectors of forward/barrel layers).
Definition: Geometry.h:34
std::map< std::string, fastsim::InteractionModel * > interactionModelMap_
Each interaction model has a unique name.
static const std::string MESSAGECATEGORY
Category of debugging messages ("FastSimulation")
std::unique_ptr< CalorimetryManager > myCalorimetry
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:52
Implementation of non-stable particle decays.
Definition: Decayer.h:37
HLT enums.
T get() const
Definition: EventSetup.h:71
int onLayer1() const
Definition: FSimTrack.h:96
bool isStable() const
Returns true if particle is considered stable.
Definition: Particle.h:175
int onHcal() const
Definition: FSimTrack.h:111
void setVFcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:95
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:145
(Kinematic) cuts on the particles that are propagated.
void setHcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:88
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:19
void endStream() override
std::vector< SimTrack > SimTrackContainer
double getMagneticFieldZ(const math::XYZTLorentzVector &position) const
Initializes the tracker geometry.
Definition: Geometry.cc:161
fastsim::ParticleFilter particleFilter_
Decides which particles have to be propagated.
def move(src, dest)
Definition: eostools.py:511
T get(const Candidate &c)
Definition: component.h:55