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
37 #include "FastSimulation/Particle/interface/ParticleTable.h" // TODO: move this
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);
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  ParticleTable::Sentry ptable(&(*pdt));
185 
186  // Get the GenParticle collection
188  iEvent.getByToken(genParticlesToken_, genParticles);
189 
190  // Load the ParticleManager which returns the particles that have to be propagated
191  // Creates a fastsim::Particle out of a GenParticle/secondary
192  fastsim::ParticleManager particleManager(*genParticles->GetEvent()
193  ,*pdt
197  ,*simTracks_
198  ,*simVertices_);
199 
200  // Initialize the calorimeter geometry
202  {
203  if(watchCaloGeometry_.check(iSetup) || watchCaloTopology_.check(iSetup)){
205  iSetup.get<CaloGeometryRecord>().get(pG);
206  myCalorimetry->getCalorimeter()->setupGeometry(*pG);
207 
208  edm::ESHandle<CaloTopology> theCaloTopology;
209  iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
210  myCalorimetry->getCalorimeter()->setupTopology(*theCaloTopology);
211  myCalorimetry->getCalorimeter()->initialize(geometry_.getMagneticFieldZ(math::XYZTLorentzVector(0., 0., 0., 0.)));
212 
213  myCalorimetry->getHFShowerLibrary()->initHFShowerLibrary(iSetup);
214  }
215 
216  // Important: this also cleans the calorimetry information from the last event
217  myCalorimetry->initialize(_randomEngine.get());
218  }
219 
220  // The vector of SimTracks needed for the CalorimetryManager
221  std::vector<FSimTrack> myFSimTracks;
222 
223 
224  LogDebug(MESSAGECATEGORY) << "################################"
225  << "\n###############################";
226 
227  // loop over particles
228  for(std::unique_ptr<fastsim::Particle> particle = particleManager.nextParticle(*_randomEngine); particle != nullptr; particle=particleManager.nextParticle(*_randomEngine))
229  {
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  {
248  LogDebug(MESSAGECATEGORY) << " moved to next layer: " << *layer;
249  LogDebug(MESSAGECATEGORY) << " new state: " << *particle;
250 
251  // Hack to interface "old" calo to "new" tracking
252  // Particle reached calorimetry so stop further propagation
254  {
255  layer = nullptr;
256  // particle no longer is on a layer
257  particle->resetOnLayer();
258  break;
259  }
260 
261  // break after 25 ns: only happens for particles stuck in loops
262  if(particle->position().T() > 25)
263  {
264  layer = nullptr;
265  // particle no longer is on a layer
266  particle->resetOnLayer();
267  break;
268  }
269 
270  // perform interaction between layer and particle
271  // do only if there is actual material
272  if(layer->getThickness(particle->position(), particle->momentum()) > 1E-10){
273  int nSecondaries = 0;
274  // loop on interaction models
275  for(fastsim::InteractionModel * interactionModel : layer->getInteractionModels())
276  {
277  LogDebug(MESSAGECATEGORY) << " interact with " << *interactionModel;
278  std::vector<std::unique_ptr<fastsim::Particle> > secondaries;
279  interactionModel->interact(*particle,*layer,secondaries,*_randomEngine);
280  nSecondaries += secondaries.size();
281  particleManager.addSecondaries(particle->position(),particle->simTrackIndex(),secondaries,layer);
282  }
283 
284  // kinematic cuts: particle might e.g. lost all its energy
285  if(!particleFilter_.acceptsEn(*particle))
286  {
287  // Add endvertex if particle did not create any secondaries
288  if(nSecondaries==0) particleManager.addEndVertex(particle.get());
289  layer = nullptr;
290  break;
291  }
292  }
293 
294  LogDebug(MESSAGECATEGORY) << "--------------------------------"
295  << "\n-------------------------------";
296  }
297 
298  // do decays
299  if(!particle->isStable() && particle->remainingProperLifeTimeC() < 1E-10)
300  {
301  LogDebug(MESSAGECATEGORY) << "Decaying particle...";
302  std::vector<std::unique_ptr<fastsim::Particle> > secondaries;
303  decayer_.decay(*particle,secondaries, _randomEngine->theEngine());
304  LogDebug(MESSAGECATEGORY) << " decay has " << secondaries.size() << " products";
305  particleManager.addSecondaries(particle->position(), particle->simTrackIndex(),secondaries);
306  continue;
307  }
308 
309  LogDebug(MESSAGECATEGORY) << "################################"
310  << "\n###############################";
311  }
312 
313 
314  // -----------------------------
315  // Hack to interface "old" calorimetry with "new" propagation in tracker
316  // The CalorimetryManager has to know which particle could in principle hit which parts of the calorimeter
317  // I think it's a bit strange to propagate the particle even further (and even decay it) if it already hits
318  // some part of the calorimetry but this is how the code works...
319  // -----------------------------
320 
321  if(particle->position().Perp2() >= 128.*128. || std::abs(particle->position().Z()) >= 302.){
322 
323  LogDebug(MESSAGECATEGORY) << "\n moving particle to calorimetry: " << *particle;
324 
325  // create FSimTrack (this is the object the old propagation uses)
326  myFSimTracks.push_back(createFSimTrack(particle.get(), &particleManager));
327  // particle was decayed
328  if(!particle->isStable() && particle->remainingProperLifeTimeC() < 1E-10)
329  {
330  continue;
331  }
332 
333  LogDebug(MESSAGECATEGORY) << "################################"
334  << "\n###############################";
335  }
336 
337  // -----------------------------
338  // End Hack
339  // -----------------------------
340 
341 
342  LogDebug(MESSAGECATEGORY) << "################################"
343  << "\n###############################";
344 
345  }
346 
347  // store simTracks and simVertices
348  iEvent.put(std::move(simTracks_));
349  iEvent.put(std::move(simVertices_));
350  // store products of interaction models, i.e. simHits
351  for(auto & interactionModel : interactionModels_)
352  {
353  interactionModel->storeProducts(iEvent);
354  }
355 
356 
357  // -----------------------------
358  // Calorimetry Manager
359  // -----------------------------
361  {
362  for(auto myFSimTrack : myFSimTracks)
363  {
364  myCalorimetry->reconstructTrack(myFSimTrack, _randomEngine.get());
365  }
366  }
367 
368 
369  // -----------------------------
370  // Store Hits
371  // -----------------------------
372  std::unique_ptr<edm::PCaloHitContainer> p4(new edm::PCaloHitContainer);
373  std::unique_ptr<edm::PCaloHitContainer> p5(new edm::PCaloHitContainer);
374  std::unique_ptr<edm::PCaloHitContainer> p6(new edm::PCaloHitContainer);
375  std::unique_ptr<edm::PCaloHitContainer> p7(new edm::PCaloHitContainer);
376 
377  std::unique_ptr<edm::SimTrackContainer> m1(new edm::SimTrackContainer);
378 
380  {
381  myCalorimetry->loadFromEcalBarrel(*p4);
382  myCalorimetry->loadFromEcalEndcap(*p5);
383  myCalorimetry->loadFromPreshower(*p6);
384  myCalorimetry->loadFromHcal(*p7);
385  if(simulateMuons){
386  myCalorimetry->harvestMuonSimTracks(*m1);
387  }
388  }
389  iEvent.put(std::move(p4),"EcalHitsEB");
390  iEvent.put(std::move(p5),"EcalHitsEE");
391  iEvent.put(std::move(p6),"EcalHitsES");
392  iEvent.put(std::move(p7),"HcalHits");
393  iEvent.put(std::move(m1),"MuonSimTracks");
394 }
395 
396 void
398 {
399  _randomEngine.reset();
400 }
401 
402 FSimTrack
404 {
405  FSimTrack myFSimTrack(particle->pdgId(),
406  particleManager->getSimTrack(particle->simTrackIndex()).momentum(),
407  particle->simVertexIndex(),
408  particle->genParticleIndex(),
409  particle->simTrackIndex(),
410  particle->charge(),
411  particle->position(),
412  particle->momentum(),
413  particleManager->getSimVertex(particle->simVertexIndex()));
414 
415  // move the particle through the caloLayers
416  fastsim::LayerNavigator caloLayerNavigator(caloGeometry_);
417  const fastsim::SimplifiedGeometry * caloLayer = nullptr;
418 
419  // moveParticleToNextLayer(..) returns 0 in case that particle decays
420  // in this case particle is propagated up to its decay vertex
421  while(caloLayerNavigator.moveParticleToNextLayer(*particle,caloLayer))
422  {
423  LogDebug(MESSAGECATEGORY) << " moved to next caloLayer: " << *caloLayer;
424  LogDebug(MESSAGECATEGORY) << " new state: " << *particle;
425 
426  // break after 25 ns: only happens for particles stuck in loops
427  if(particle->position().T() > 50)
428  {
429  caloLayer = nullptr;
430  break;
431  }
432 
434  // Define ParticlePropagators (RawParticle) needed for CalorimetryManager and save them
436 
437  RawParticle PP(particle->pdgId(), particle->momentum());
438  PP.setVertex(particle->position());
439 
440  // no material
441  if(caloLayer->getThickness(particle->position(), particle->momentum()) < 1E-10)
442  {
443  // unfortunately needed for CalorimetryManager
445  if(!myFSimTrack.onEcal())
446  {
447  myFSimTrack.setEcal(PP, 0);
448  }
449  }
450  else if(caloLayer->getCaloType() == fastsim::SimplifiedGeometry::HCAL){
451  if(!myFSimTrack.onHcal())
452  {
453  myFSimTrack.setHcal(PP, 0);
454  }
455  }
456  else if(caloLayer->getCaloType() == fastsim::SimplifiedGeometry::VFCAL){
457  if(!myFSimTrack.onVFcal())
458  {
459  myFSimTrack.setVFcal(PP, 0);
460  }
461  }
462 
463  // not necessary to continue propagation
465  {
466  myFSimTrack.setGlobal();
467  caloLayer = nullptr;
468  break;
469  }
470 
471  continue;
472  }
473 
474  // Stupid variable used by the old propagator
475  // For details check BaseParticlePropagator.h
476  int success = 0;
477  if(caloLayer->isForward())
478  {
479  success = 2;
480  // particle moves inwards
481  if(particle->position().Z() * particle->momentum().Z() < 0)
482  {
483  success *= -1;
484  }
485  }
486  else
487  {
488  success = 1;
489  // particle moves inwards
490  if(particle->momentum().X() * particle->position().X() + particle->momentum().Y() * particle->position().Y() < 0)
491  {
492  success *= -1;
493  }
494  }
495 
496  // Save the hit
498  {
499  if(!myFSimTrack.onLayer1())
500  {
501  myFSimTrack.setLayer1(PP, success);
502  }
503  }
504 
506  {
507  if(!myFSimTrack.onLayer2())
508  {
509  myFSimTrack.setLayer2(PP, success);
510  }
511  }
512 
514  {
515  if(!myFSimTrack.onEcal())
516  {
517  myFSimTrack.setEcal(PP, success);
518  }
519  }
520 
522  {
523  if(!myFSimTrack.onHcal())
524  {
525  myFSimTrack.setHcal(PP, success);
526  }
527  }
528 
530  {
531  if(!myFSimTrack.onVFcal())
532  {
533  myFSimTrack.setVFcal(PP, success);
534  }
535  }
536 
537  // Particle reached end of detector
539  {
540  myFSimTrack.setGlobal();
541  caloLayer = nullptr;
542  break;
543  }
544 
545  LogDebug(MESSAGECATEGORY) << "--------------------------------"
546  << "\n-------------------------------";
547  }
548 
549  // do decays
550  // don't have to worry about daughters if particle already within the calorimetry
551  // since they will be rejected by the vertex cut of the ParticleFilter
552  if(!particle->isStable() && particle->remainingProperLifeTimeC() < 1E-10)
553  {
554  LogDebug(MESSAGECATEGORY) << "Decaying particle...";
555  std::vector<std::unique_ptr<fastsim::Particle> > secondaries;
556  decayer_.decay(*particle,secondaries, _randomEngine->theEngine());
557  LogDebug(MESSAGECATEGORY) << " decay has " << secondaries.size() << " products";
558  particleManager->addSecondaries(particle->position(), particle->simTrackIndex(),secondaries);
559  }
560 
561  return myFSimTrack;
562 }
563 
#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:137
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)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
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: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.
void setLayer2(const RawParticle &pp, int success)
Set the preshower layer2 variables.
Definition: FSimTrack.cc:74
virtual FSimTrack createFSimTrack(fastsim::Particle *particle, fastsim::ParticleManager *particleManager)
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:136
#define nullptr
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...
void getData(T &iHolder) const
Definition: EventSetup.h:98
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.
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:230
void setLayer1(const RawParticle &pp, int success)
Set the preshower layer1 variables.
Definition: FSimTrack.cc:67
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 getMagneticFieldZ(const math::XYZTLorentzVector &position) const
Initializes the tracker geometry.
Definition: Geometry.cc:161
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:33
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:57
Implementation of non-stable particle decays.
Definition: Decayer.h:37
HLT enums.
T get() const
Definition: EventSetup.h:63
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
fastsim::ParticleFilter particleFilter_
Decides which particles have to be propagated.
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
Definition: RawParticle.h:288
def move(src, dest)
Definition: eostools.py:510
T get(const Candidate &c)
Definition: component.h:55