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 
92  std::vector<std::unique_ptr<fastsim::InteractionModel> > interactionModels_;
93  std::map<std::string, fastsim::InteractionModel*> interactionModelMap_;
95 };
96 
97 const std::string FastSimProducer::MESSAGECATEGORY = "FastSimulation";
98 
100  : genParticlesToken_(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("src"))),
101  geometry_(iConfig.getParameter<edm::ParameterSet>("trackerDefinition")),
102  caloGeometry_(iConfig.getParameter<edm::ParameterSet>("caloDefinition")),
103  beamPipeRadius_(iConfig.getParameter<double>("beamPipeRadius")),
104  deltaRchargedMother_(iConfig.getParameter<double>("deltaRchargedMother")),
105  particleFilter_(iConfig.getParameter<edm::ParameterSet>("particleFilter")),
107  simulateCalorimetry(iConfig.getParameter<bool>("simulateCalorimetry")),
108  simulateMuons(iConfig.getParameter<bool>("simulateMuons")) {
109  //----------------
110  // define interaction models
111  //---------------
112 
113  const edm::ParameterSet& modelCfgs = iConfig.getParameter<edm::ParameterSet>("interactionModels");
114  for (const std::string& modelName : modelCfgs.getParameterNames()) {
115  const edm::ParameterSet& modelCfg = modelCfgs.getParameter<edm::ParameterSet>(modelName);
116  std::string modelClassName(modelCfg.getParameter<std::string>("className"));
117  // Use plugin-factory to create model
118  std::unique_ptr<fastsim::InteractionModel> interactionModel(
119  fastsim::InteractionModelFactory::get()->create(modelClassName, modelName, modelCfg));
120  if (!interactionModel.get()) {
121  throw cms::Exception("FastSimProducer") << "InteractionModel " << modelName << " could not be created";
122  }
123  // Add model to list
124  interactionModels_.push_back(std::move(interactionModel));
125  // and create the map
127  }
128 
129  //----------------
130  // calorimetry
131  //---------------
132 
133  if (simulateCalorimetry) {
134  myCalorimetry.reset(new CalorimetryManager(nullptr,
135  iConfig.getParameter<edm::ParameterSet>("Calorimetry"),
136  iConfig.getParameter<edm::ParameterSet>("MaterialEffectsForMuonsInECAL"),
137  iConfig.getParameter<edm::ParameterSet>("MaterialEffectsForMuonsInHCAL"),
138  iConfig.getParameter<edm::ParameterSet>("GFlash")));
139  }
140 
141  //----------------
142  // register products
143  //----------------
144 
145  // SimTracks and SimVertices
146  produces<edm::SimTrackContainer>();
147  produces<edm::SimVertexContainer>();
148  // products of interaction models, i.e. simHits
149  for (auto& interactionModel : interactionModels_) {
150  interactionModel->registerProducts(producesCollector());
151  }
152  produces<edm::PCaloHitContainer>("EcalHitsEB");
153  produces<edm::PCaloHitContainer>("EcalHitsEE");
154  produces<edm::PCaloHitContainer>("EcalHitsES");
155  produces<edm::PCaloHitContainer>("HcalHits");
156  produces<edm::SimTrackContainer>("MuonSimTracks");
157 }
158 
160  _randomEngine = std::make_unique<RandomEngineAndDistribution>(id);
161 }
162 
164  LogDebug(MESSAGECATEGORY) << " produce";
165 
168 
169  // Define containers for SimTracks, SimVertices
170  std::unique_ptr<edm::SimTrackContainer> simTracks_(new edm::SimTrackContainer);
171  std::unique_ptr<edm::SimVertexContainer> simVertices_(new edm::SimVertexContainer);
172 
173  // Get the particle data table (in case lifetime or charge of GenParticles not set)
175  iSetup.getData(pdt);
176 
177  // Get the GenParticle collection
179  iEvent.getByToken(genParticlesToken_, genParticles);
180 
181  // Load the ParticleManager which returns the particles that have to be propagated
182  // Creates a fastsim::Particle out of a GenParticle/secondary
183  fastsim::ParticleManager particleManager(*genParticles->GetEvent(),
184  *pdt,
188  *simTracks_,
189  *simVertices_);
190 
191  // Initialize the calorimeter geometry
192  if (simulateCalorimetry) {
193  if (watchCaloGeometry_.check(iSetup) || watchCaloTopology_.check(iSetup)) {
195  iSetup.get<CaloGeometryRecord>().get(pG);
196  myCalorimetry->getCalorimeter()->setupGeometry(*pG);
197 
198  edm::ESHandle<CaloTopology> theCaloTopology;
199  iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
200  myCalorimetry->getCalorimeter()->setupTopology(*theCaloTopology);
201  myCalorimetry->getCalorimeter()->initialize(geometry_.getMagneticFieldZ(math::XYZTLorentzVector(0., 0., 0., 0.)));
202 
203  myCalorimetry->getHFShowerLibrary()->initHFShowerLibrary(iSetup);
204  }
205 
206  // Important: this also cleans the calorimetry information from the last event
207  myCalorimetry->initialize(_randomEngine.get());
208  }
209 
210  // The vector of SimTracks needed for the CalorimetryManager
211  std::vector<FSimTrack> myFSimTracks;
212 
213  LogDebug(MESSAGECATEGORY) << "################################"
214  << "\n###############################";
215 
216  // loop over particles
217  for (std::unique_ptr<fastsim::Particle> particle = particleManager.nextParticle(*_randomEngine); particle != nullptr;
218  particle = particleManager.nextParticle(*_randomEngine)) {
219  LogDebug(MESSAGECATEGORY) << "\n moving NEXT particle: " << *particle;
220 
221  // -----------------------------
222  // This condition is necessary because of hack for calorimetry
223  // -> The CalorimetryManager should also be implemented based on this new FastSim classes (Particle.h) in a future project.
224  // A second loop (below) loops over all parts of the calorimetry in order to create a track of the old FastSim class FSimTrack.
225  // The condition below (R<128, z<302) makes sure that the particle geometrically is outside the tracker boundaries
226  // -----------------------------
227 
228  if (particle->position().Perp2() < 128. * 128. && std::abs(particle->position().Z()) < 302.) {
229  // move the particle through the layers
230  fastsim::LayerNavigator layerNavigator(geometry_);
231  const fastsim::SimplifiedGeometry* layer = nullptr;
232 
233  // moveParticleToNextLayer(..) returns 0 in case that particle decays
234  // in this case particle is propagated up to its decay vertex
235  while (layerNavigator.moveParticleToNextLayer(*particle, layer)) {
236  LogDebug(MESSAGECATEGORY) << " moved to next layer: " << *layer;
237  LogDebug(MESSAGECATEGORY) << " new state: " << *particle;
238 
239  // Hack to interface "old" calo to "new" tracking
240  // Particle reached calorimetry so stop further propagation
242  layer = nullptr;
243  // particle no longer is on a layer
244  particle->resetOnLayer();
245  break;
246  }
247 
248  // break after 25 ns: only happens for particles stuck in loops
249  if (particle->position().T() > 25) {
250  layer = nullptr;
251  // particle no longer is on a layer
252  particle->resetOnLayer();
253  break;
254  }
255 
256  // perform interaction between layer and particle
257  // do only if there is actual material
258  if (layer->getThickness(particle->position(), particle->momentum()) > 1E-10) {
259  int nSecondaries = 0;
260  // loop on interaction models
261  for (fastsim::InteractionModel* interactionModel : layer->getInteractionModels()) {
262  LogDebug(MESSAGECATEGORY) << " interact with " << *interactionModel;
263  std::vector<std::unique_ptr<fastsim::Particle> > secondaries;
264  interactionModel->interact(*particle, *layer, secondaries, *_randomEngine);
265  nSecondaries += secondaries.size();
266  particleManager.addSecondaries(particle->position(), particle->simTrackIndex(), secondaries, layer);
267  }
268 
269  // kinematic cuts: particle might e.g. lost all its energy
270  if (!particleFilter_.acceptsEn(*particle)) {
271  // Add endvertex if particle did not create any secondaries
272  if (nSecondaries == 0)
273  particleManager.addEndVertex(particle.get());
274  layer = nullptr;
275  break;
276  }
277  }
278 
279  LogDebug(MESSAGECATEGORY) << "--------------------------------"
280  << "\n-------------------------------";
281  }
282 
283  // do decays
284  if (!particle->isStable() && particle->remainingProperLifeTimeC() < 1E-10) {
285  LogDebug(MESSAGECATEGORY) << "Decaying particle...";
286  std::vector<std::unique_ptr<fastsim::Particle> > secondaries;
287  decayer_.decay(*particle, secondaries, _randomEngine->theEngine());
288  LogDebug(MESSAGECATEGORY) << " decay has " << secondaries.size() << " products";
289  particleManager.addSecondaries(particle->position(), particle->simTrackIndex(), secondaries);
290  continue;
291  }
292 
293  LogDebug(MESSAGECATEGORY) << "################################"
294  << "\n###############################";
295  }
296 
297  // -----------------------------
298  // Hack to interface "old" calorimetry with "new" propagation in tracker
299  // The CalorimetryManager has to know which particle could in principle hit which parts of the calorimeter
300  // I think it's a bit strange to propagate the particle even further (and even decay it) if it already hits
301  // some part of the calorimetry but this is how the code works...
302  // -----------------------------
303 
304  if (particle->position().Perp2() >= 128. * 128. || std::abs(particle->position().Z()) >= 302.) {
305  LogDebug(MESSAGECATEGORY) << "\n moving particle to calorimetry: " << *particle;
306 
307  // create FSimTrack (this is the object the old propagation uses)
308  myFSimTracks.push_back(createFSimTrack(particle.get(), &particleManager, *pdt));
309  // particle was decayed
310  if (!particle->isStable() && particle->remainingProperLifeTimeC() < 1E-10) {
311  continue;
312  }
313 
314  LogDebug(MESSAGECATEGORY) << "################################"
315  << "\n###############################";
316  }
317 
318  // -----------------------------
319  // End Hack
320  // -----------------------------
321 
322  LogDebug(MESSAGECATEGORY) << "################################"
323  << "\n###############################";
324  }
325 
326  // store simTracks and simVertices
327  iEvent.put(std::move(simTracks_));
328  iEvent.put(std::move(simVertices_));
329  // store products of interaction models, i.e. simHits
330  for (auto& interactionModel : interactionModels_) {
331  interactionModel->storeProducts(iEvent);
332  }
333 
334  // -----------------------------
335  // Calorimetry Manager
336  // -----------------------------
337  if (simulateCalorimetry) {
338  for (auto myFSimTrack : myFSimTracks) {
339  myCalorimetry->reconstructTrack(myFSimTrack, _randomEngine.get());
340  }
341  }
342 
343  // -----------------------------
344  // Store Hits
345  // -----------------------------
346  std::unique_ptr<edm::PCaloHitContainer> p4(new edm::PCaloHitContainer);
347  std::unique_ptr<edm::PCaloHitContainer> p5(new edm::PCaloHitContainer);
348  std::unique_ptr<edm::PCaloHitContainer> p6(new edm::PCaloHitContainer);
349  std::unique_ptr<edm::PCaloHitContainer> p7(new edm::PCaloHitContainer);
350 
351  std::unique_ptr<edm::SimTrackContainer> m1(new edm::SimTrackContainer);
352 
353  if (simulateCalorimetry) {
354  myCalorimetry->loadFromEcalBarrel(*p4);
355  myCalorimetry->loadFromEcalEndcap(*p5);
356  myCalorimetry->loadFromPreshower(*p6);
357  myCalorimetry->loadFromHcal(*p7);
358  if (simulateMuons) {
359  myCalorimetry->harvestMuonSimTracks(*m1);
360  }
361  }
362  iEvent.put(std::move(p4), "EcalHitsEB");
363  iEvent.put(std::move(p5), "EcalHitsEE");
364  iEvent.put(std::move(p6), "EcalHitsES");
365  iEvent.put(std::move(p7), "HcalHits");
366  iEvent.put(std::move(m1), "MuonSimTracks");
367 }
368 
370 
372  fastsim::ParticleManager* particleManager,
373  HepPDT::ParticleDataTable const& particleTable) {
374  FSimTrack myFSimTrack(particle->pdgId(),
375  particleManager->getSimTrack(particle->simTrackIndex()).momentum(),
376  particle->simVertexIndex(),
377  particle->genParticleIndex(),
378  particle->simTrackIndex(),
379  particle->charge(),
380  particle->position(),
381  particle->momentum(),
382  particleManager->getSimVertex(particle->simVertexIndex()));
383 
384  // move the particle through the caloLayers
385  fastsim::LayerNavigator caloLayerNavigator(caloGeometry_);
386  const fastsim::SimplifiedGeometry* caloLayer = nullptr;
387 
388  // moveParticleToNextLayer(..) returns 0 in case that particle decays
389  // in this case particle is propagated up to its decay vertex
390  while (caloLayerNavigator.moveParticleToNextLayer(*particle, caloLayer)) {
391  LogDebug(MESSAGECATEGORY) << " moved to next caloLayer: " << *caloLayer;
392  LogDebug(MESSAGECATEGORY) << " new state: " << *particle;
393 
394  // break after 25 ns: only happens for particles stuck in loops
395  if (particle->position().T() > 50) {
396  caloLayer = nullptr;
397  break;
398  }
399 
401  // Define ParticlePropagators (RawParticle) needed for CalorimetryManager and save them
403 
404  RawParticle PP = makeParticle(&particleTable, particle->pdgId(), particle->momentum(), particle->position());
405 
406  // no material
407  if (caloLayer->getThickness(particle->position(), particle->momentum()) < 1E-10) {
408  // unfortunately needed for CalorimetryManager
409  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::ECAL) {
410  if (!myFSimTrack.onEcal()) {
411  myFSimTrack.setEcal(PP, 0);
412  }
413  } else if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::HCAL) {
414  if (!myFSimTrack.onHcal()) {
415  myFSimTrack.setHcal(PP, 0);
416  }
417  } else if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::VFCAL) {
418  if (!myFSimTrack.onVFcal()) {
419  myFSimTrack.setVFcal(PP, 0);
420  }
421  }
422 
423  // not necessary to continue propagation
424  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::VFCAL) {
425  myFSimTrack.setGlobal();
426  caloLayer = nullptr;
427  break;
428  }
429 
430  continue;
431  }
432 
433  // Stupid variable used by the old propagator
434  // For details check BaseParticlePropagator.h
435  int success = 0;
436  if (caloLayer->isForward()) {
437  success = 2;
438  // particle moves inwards
439  if (particle->position().Z() * particle->momentum().Z() < 0) {
440  success *= -1;
441  }
442  } else {
443  success = 1;
444  // particle moves inwards
445  if (particle->momentum().X() * particle->position().X() + particle->momentum().Y() * particle->position().Y() <
446  0) {
447  success *= -1;
448  }
449  }
450 
451  // Save the hit
453  if (!myFSimTrack.onLayer1()) {
454  myFSimTrack.setLayer1(PP, success);
455  }
456  }
457 
459  if (!myFSimTrack.onLayer2()) {
460  myFSimTrack.setLayer2(PP, success);
461  }
462  }
463 
464  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::ECAL) {
465  if (!myFSimTrack.onEcal()) {
466  myFSimTrack.setEcal(PP, success);
467  }
468  }
469 
470  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::HCAL) {
471  if (!myFSimTrack.onHcal()) {
472  myFSimTrack.setHcal(PP, success);
473  }
474  }
475 
476  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::VFCAL) {
477  if (!myFSimTrack.onVFcal()) {
478  myFSimTrack.setVFcal(PP, success);
479  }
480  }
481 
482  // Particle reached end of detector
483  if (caloLayer->getCaloType() == fastsim::SimplifiedGeometry::VFCAL) {
484  myFSimTrack.setGlobal();
485  caloLayer = nullptr;
486  break;
487  }
488 
489  LogDebug(MESSAGECATEGORY) << "--------------------------------"
490  << "\n-------------------------------";
491  }
492 
493  // do decays
494  // don't have to worry about daughters if particle already within the calorimetry
495  // since they will be rejected by the vertex cut of the ParticleFilter
496  if (!particle->isStable() && particle->remainingProperLifeTimeC() < 1E-10) {
497  LogDebug(MESSAGECATEGORY) << "Decaying particle...";
498  std::vector<std::unique_ptr<fastsim::Particle> > secondaries;
499  decayer_.decay(*particle, secondaries, _randomEngine->theEngine());
500  LogDebug(MESSAGECATEGORY) << " decay has " << secondaries.size() << " products";
501  particleManager->addSecondaries(particle->position(), particle->simTrackIndex(), secondaries);
502  }
503 
504  return myFSimTrack;
505 }
506 
#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:49
std::vector< PCaloHit > PCaloHitContainer
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
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:140
int onLayer2() const
Definition: FSimTrack.h:106
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:525
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:150
int genParticleIndex() const
Return index of the particle in the genParticle vector.
Definition: Particle.h:165
edm::EDGetTokenT< edm::HepMCProduct > genParticlesToken_
Token to get the genParticles.
void setEcal(const RawParticle &pp, int success)
Set the ecal variables.
Definition: FSimTrack.cc:118
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:161
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:112
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:134
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:111
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:28
bool getData(T &iHolder) const
Definition: EventSetup.h:113
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:159
~FastSimProducer() override
int iEvent
Definition: GenABIO.cc:224
void setLayer1(const RawParticle &pp, int success)
Set the preshower layer1 variables.
Definition: FSimTrack.cc:106
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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:121
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:28
int simTrackIndex() const
Return index of the SimTrack.
Definition: Particle.h:153
void produce(edm::Event &, const edm::EventSetup &) override
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
edm::ESWatcher< CaloTopologyRecord > watchCaloTopology_
const std::vector< InteractionModel * > & getInteractionModels() const
Return the vector of all interaction models that are assigned with a layer.
double charge() const
Return charge of the particle.
Definition: Particle.h:137
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:32
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:34
HLT enums.
T get() const
Definition: EventSetup.h:73
int onLayer1() const
Definition: FSimTrack.h:101
bool isStable() const
Returns true if particle is considered stable.
Definition: Particle.h:171
int onHcal() const
Definition: FSimTrack.h:116
void setVFcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:130
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:143
(Kinematic) cuts on the particles that are propagated.
void setHcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:124
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:16
void endStream() override
std::vector< SimTrack > SimTrackContainer
double getMagneticFieldZ(const math::XYZTLorentzVector &position) const
Initializes the tracker geometry.
Definition: Geometry.cc:154
fastsim::ParticleFilter particleFilter_
Decides which particles have to be propagated.
def move(src, dest)
Definition: eostools.py:511