CMS 3D CMS Logo

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