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