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