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