OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
ParallelTTracker.h
Go to the documentation of this file.
1 #ifndef OPAL_ParallelTTracker_HH
2 #define OPAL_ParallelTTracker_HH
3 
4 // ------------------------------------------------------------------------
5 // $RCSfile: ParallelTTracker.h,v $
6 // ------------------------------------------------------------------------
7 // $Revision: 1.1.2.1 $
8 // ------------------------------------------------------------------------
9 // Copyright: see Copyright.readme
10 // ------------------------------------------------------------------------
11 //
12 // Class: ParallelTTracker
13 //
14 // ------------------------------------------------------------------------
15 //
16 // $Date: 2004/11/12 20:10:11 $
17 // $Author: adelmann $
18 //
19 // ------------------------------------------------------------------------
20 
21 #include "Algorithms/Tracker.h"
22 #include "Steppers/BorisPusher.h"
23 #include "Structure/DataSink.h"
25 
26 #include "BasicActions/Option.h"
27 #include "Utilities/Options.h"
28 
29 #include "Physics/Physics.h"
30 
32 #include "Algorithms/IndexMap.h"
34 #include "AbsBeamline/BeamBeam.h"
37 #include "AbsBeamline/Corrector.h"
38 #include "AbsBeamline/Diagnostic.h"
39 #include "AbsBeamline/Degrader.h"
40 #include "AbsBeamline/Drift.h"
43 #include "AbsBeamline/Lambertson.h"
44 #include "AbsBeamline/Marker.h"
45 #include "AbsBeamline/Monitor.h"
46 #include "AbsBeamline/Multipole.h"
47 #include "AbsBeamline/MultipoleT.h"
48 #include "AbsBeamline/Probe.h"
49 #include "AbsBeamline/RBend.h"
50 #include "AbsBeamline/RBend3D.h"
51 #include "AbsBeamline/RFCavity.h"
54 #include "AbsBeamline/SBend.h"
55 #include "AbsBeamline/Separator.h"
56 #include "AbsBeamline/Septum.h"
57 #include "AbsBeamline/Solenoid.h"
60 
61 #include "Beamlines/Beamline.h"
62 #include "Elements/OpalBeamline.h"
63 #include "Solvers/WakeFunction.hh"
64 
65 #include <list>
66 #include <vector>
67 
69 
70 class ParallelTTracker: public Tracker {
71 
72 public:
74  // The beam line to be tracked is "bl".
75  // The particle reference data are taken from "data".
76  // The particle bunch tracked is initially empty.
77  // If [b]revBeam[/b] is true, the beam runs from s = C to s = 0.
78  // If [b]revTrack[/b] is true, we track against the beam.
79  explicit ParallelTTracker(const Beamline &bl,
80  const PartData &data,
81  bool revBeam,
82  bool revTrack);
83 
85  // The beam line to be tracked is "bl".
86  // The particle reference data are taken from "data".
87  // The particle bunch tracked is taken from [b]bunch[/b].
88  // If [b]revBeam[/b] is true, the beam runs from s = C to s = 0.
89  // If [b]revTrack[/b] is true, we track against the beam.
90  explicit ParallelTTracker(const Beamline &bl,
92  DataSink &ds,
93  const PartData &data,
94  bool revBeam,
95  bool revTrack,
96  const std::vector<unsigned long long> &maxSTEPS,
97  double zstart,
98  const std::vector<double> &zstop,
99  const std::vector<double> &dt);
100 
101 
102  virtual ~ParallelTTracker();
103 
104  virtual void visitAlignWrapper(const AlignWrapper &);
105 
107  virtual void visitBeamBeam(const BeamBeam &);
108 
110  virtual void visitBeamStripping(const BeamStripping &);
111 
113  virtual void visitCCollimator(const CCollimator &);
114 
115 
117  virtual void visitCorrector(const Corrector &);
118 
120  virtual void visitDegrader(const Degrader &);
121 
123  virtual void visitDiagnostic(const Diagnostic &);
124 
126  virtual void visitDrift(const Drift &);
127 
129  virtual void visitFlexibleCollimator(const FlexibleCollimator &);
130 
132  virtual void visitLambertson(const Lambertson &);
133 
135  virtual void visitMarker(const Marker &);
136 
138  virtual void visitMonitor(const Monitor &);
139 
141  virtual void visitMultipole(const Multipole &);
142 
144  virtual void visitMultipoleT(const MultipoleT &);
145 
147  virtual void visitProbe(const Probe &);
148 
150  virtual void visitRBend(const RBend &);
151 
153  virtual void visitRBend3D(const RBend3D &);
154 
156  virtual void visitRFCavity(const RFCavity &);
157 
159  virtual void visitTravelingWave(const TravelingWave &);
160 
162  virtual void visitRFQuadrupole(const RFQuadrupole &);
163 
165  virtual void visitSBend(const SBend &);
166 
168  virtual void visitSeparator(const Separator &);
169 
171  virtual void visitSeptum(const Septum &);
172 
174  virtual void visitSolenoid(const Solenoid &);
175 
177  virtual void visitSource(const Source &);
178 
180  virtual void visitParallelPlate(const ParallelPlate &);
181 
183  virtual void visitCyclotronValley(const CyclotronValley &);
184 
185 
187  // overwrite the execute-methode from DefaultVisitor
188  virtual void visitBeamline(const Beamline &);
189 
191  // overwrite the execute-methode from DefaultVisitor
192  virtual void execute();
193 
194 private:
195 
196  // Not implemented.
199  void operator=(const ParallelTTracker &);
200 
201  /******************** STATE VARIABLES ***********************************/
202 
204 
206 
208 
210 
212 
214 
215  double pathLength_m;
216 
218  double zstart_m;
219 
225 
227 
228  // This variable controls the minimal number of steps of emission (using bins)
229  // before we can merge the bins
231 
232  // The space charge solver crashes if we use less than ~10 particles.
233  // This variable controls the number of particles to be emitted before we use
234  // the space charge solver.
236 
237  // this variable controls the minimal number of steps until we repartition the particles
238  unsigned int repartFreq_m;
239 
240  unsigned int emissionSteps_m;
241 
243 
249 
250  std::set<ParticleMatterInteractionHandler*> activeParticleMatterInteractionHandlers_m;
252 
254 
255 #ifdef OPAL_DKS
256  DKSOPAL *dksbase;
257  void *r_ptr;
258  void *p_ptr;
259  void *Ef_ptr;
260  void *Bf_ptr;
261  void *dt_ptr;
262 
263  int stream1;
264  int stream2;
265 
266  int numDeviceElements;
267 
268  void setupDKS();
269  void allocateDeviceMemory();
270  void freeDeviceMemory();
271  void sendRPdt();
272  void sendEfBf();
273  void pushParticlesDKS(bool send = true);
274  void kickParticlesDKS();
275 
276 #endif
277 
278  /********************** END VARIABLES ***********************************/
279 
280  void kickParticles(const BorisPusher &pusher);
281  void pushParticles(const BorisPusher &pusher);
282  void updateReferenceParticle(const BorisPusher &pusher);
283 
284  void writePhaseSpace(const long long step, bool psDump, bool statDump);
285 
286  /********** BEGIN AUTOPHSING STUFF **********/
287  void updateRFElement(std::string elName, double maxPhi);
288  void printRFPhases();
289  void saveCavityPhases();
290  void restoreCavityPhases();
291  /************ END AUTOPHSING STUFF **********/
292 
293  void prepareSections();
294 
295  void timeIntegration1(BorisPusher & pusher);
296  void timeIntegration2(BorisPusher & pusher);
297  void selectDT(bool backTrack = false);
298  void changeDT(bool backTrack = false);
299  void emitParticles(long long step);
303  void computeSpaceChargeFields(unsigned long long step);
304  // void prepareOpalBeamlineSections();
305  void dumpStats(long long step, bool psDump, bool statDump);
306  void setOptionalVariables();
307  bool hasEndOfLineReached();
308  void handleRestartRun();
309  void prepareEmission();
310  void setTime();
311  void doBinaryRepartition();
312 
313  void transformBunch(const CoordinateSystemTrafo &trafo);
314 
315  void updateReference(const BorisPusher &pusher);
316  void updateRefToLabCSTrafo();
317  void applyFractionalStep(const BorisPusher &pusher, double tau);
318  void findStartPosition(const BorisPusher &pusher);
319  void autophaseCavities(const BorisPusher &pusher);
320 
322 };
323 
325  itsOpalBeamline_m.visit(wrap, *this, itsBunch_m);
326 }
327 
329  itsOpalBeamline_m.visit(bb, *this, itsBunch_m);
330 }
331 
333  itsOpalBeamline_m.visit(bstp, *this, itsBunch_m);
334 }
335 
337  itsOpalBeamline_m.visit(coll, *this, itsBunch_m);
338 }
339 
340 
341 inline void ParallelTTracker::visitCorrector(const Corrector &corr) {
342  itsOpalBeamline_m.visit(corr, *this, itsBunch_m);
343 }
344 
345 
347  itsOpalBeamline_m.visit(deg, *this, itsBunch_m);
348 }
349 
350 
352  itsOpalBeamline_m.visit(diag, *this, itsBunch_m);
353 }
354 
355 
356 inline void ParallelTTracker::visitDrift(const Drift &drift) {
357  itsOpalBeamline_m.visit(drift, *this, itsBunch_m);
358 }
359 
360 
362  itsOpalBeamline_m.visit(coll, *this, itsBunch_m);
363 }
364 
365 
367  itsOpalBeamline_m.visit(lamb, *this, itsBunch_m);
368 }
369 
370 
371 inline void ParallelTTracker::visitMarker(const Marker &marker) {
372  itsOpalBeamline_m.visit(marker, *this, itsBunch_m);
373 }
374 
375 
376 inline void ParallelTTracker::visitMonitor(const Monitor &mon) {
377  itsOpalBeamline_m.visit(mon, *this, itsBunch_m);
378 }
379 
380 
381 inline void ParallelTTracker::visitMultipole(const Multipole &mult) {
382  itsOpalBeamline_m.visit(mult, *this, itsBunch_m);
383 }
384 
386  itsOpalBeamline_m.visit(mult, *this, itsBunch_m);
387 }
388 
389 inline void ParallelTTracker::visitProbe(const Probe &prob) {
390  itsOpalBeamline_m.visit(prob, *this, itsBunch_m);
391 }
392 
393 
394 inline void ParallelTTracker::visitRBend(const RBend &bend) {
395  itsOpalBeamline_m.visit(bend, *this, itsBunch_m);
396 }
397 
398 inline void ParallelTTracker::visitRBend3D(const RBend3D &bend) {
399  itsOpalBeamline_m.visit(bend, *this, itsBunch_m);
400 }
401 
402 
404  itsOpalBeamline_m.visit(as, *this, itsBunch_m);
405 }
406 
408  itsOpalBeamline_m.visit(as, *this, itsBunch_m);
409 }
410 
411 
413  itsOpalBeamline_m.visit(rfq, *this, itsBunch_m);
414 }
415 
416 inline void ParallelTTracker::visitSBend(const SBend &bend) {
417  itsOpalBeamline_m.visit(bend, *this, itsBunch_m);
418 }
419 
420 
422  itsOpalBeamline_m.visit(sep, *this, itsBunch_m);
423 }
424 
425 
426 inline void ParallelTTracker::visitSeptum(const Septum &sept) {
427  itsOpalBeamline_m.visit(sept, *this, itsBunch_m);
428 }
429 
430 
431 inline void ParallelTTracker::visitSolenoid(const Solenoid &solenoid) {
432  itsOpalBeamline_m.visit(solenoid, *this, itsBunch_m);
433 }
434 
435 inline void ParallelTTracker::visitSource(const Source &source) {
436  itsOpalBeamline_m.visit(source, *this, itsBunch_m);
437 }
438 
440  itsOpalBeamline_m.visit(pplate, *this, itsBunch_m);
441 }
442 
444  itsOpalBeamline_m.visit(cv, *this, itsBunch_m);
445 }
446 
447 inline void ParallelTTracker::kickParticles(const BorisPusher &pusher) {
448  int localNum = itsBunch_m->getLocalNum();
449  for (int i = 0; i < localNum; ++i)
450  pusher.kick(itsBunch_m->R[i], itsBunch_m->P[i], itsBunch_m->Ef[i], itsBunch_m->Bf[i], itsBunch_m->dt[i]);
451 }
452 
453 inline void ParallelTTracker::pushParticles(const BorisPusher &pusher) {
455 
456  for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
457  pusher.push(itsBunch_m->R[i], itsBunch_m->P[i], itsBunch_m->dt[i]);
458  }
460 }
461 
462 #ifdef OPAL_DKS
463 inline void ParallelTTracker::setupDKS() {
464  dksbase = new DKSOPAL;
465  dksbase->setAPI("Cuda");
466  dksbase->setDevice("-gpu");
467  dksbase->initDevice();
468 
469  dksbase->createStream(stream1);
470  dksbase->createStream(stream2);
471 }
472 
473 
474 inline void ParallelTTracker::allocateDeviceMemory() {
475 
476  int ierr;
477  r_ptr = dksbase->allocateMemory<Vector_t>(itsBunch_m->getLocalNum(), ierr);
478  p_ptr = dksbase->allocateMemory<Vector_t>(itsBunch_m->getLocalNum(), ierr);
479  Ef_ptr = dksbase->allocateMemory<Vector_t>(itsBunch_m->getLocalNum(), ierr);
480  Bf_ptr = dksbase->allocateMemory<Vector_t>(itsBunch_m->getLocalNum(), ierr);
481  dt_ptr = dksbase->allocateMemory<double>(itsBunch_m->getLocalNum(), ierr);
482 
483  if (Ippl::getNodes() == 1) {
484  dksbase->registerHostMemory(&itsBunch_m->R[0], itsBunch_m->getLocalNum());
485  dksbase->registerHostMemory(&itsBunch_m->P[0], itsBunch_m->getLocalNum());
486  dksbase->registerHostMemory(&itsBunch_m->Ef[0], itsBunch_m->getLocalNum());
487  dksbase->registerHostMemory(&itsBunch_m->Bf[0], itsBunch_m->getLocalNum());
488  dksbase->registerHostMemory(&itsBunch_m->dt[0], itsBunch_m->getLocalNum());
489 
490  }
491 
492  numDeviceElements = itsBunch_m->getLocalNum();
493 }
494 
495 inline void ParallelTTracker::freeDeviceMemory() {
496  dksbase->freeMemory<Vector_t>(r_ptr, numDeviceElements);
497  dksbase->freeMemory<Vector_t>(p_ptr, numDeviceElements);
498  dksbase->freeMemory<Vector_t>(Ef_ptr, numDeviceElements);
499  dksbase->freeMemory<Vector_t>(Bf_ptr, numDeviceElements);
500  dksbase->freeMemory<double>(dt_ptr, numDeviceElements);
501 
502  if (Ippl::getNodes() == 1) {
503  dksbase->unregisterHostMemory(&itsBunch_m->R[0]);
504  dksbase->unregisterHostMemory(&itsBunch_m->P[0]);
505  dksbase->unregisterHostMemory(&itsBunch_m->Ef[0]);
506  dksbase->unregisterHostMemory(&itsBunch_m->Bf[0]);
507  dksbase->unregisterHostMemory(&itsBunch_m->dt[0]);
508  }
509 }
510 
511 inline void ParallelTTracker::sendRPdt() {
512  dksbase->writeDataAsync<Vector_t>(r_ptr, &itsBunch_m->R[0], itsBunch_m->getLocalNum(), stream1);
513  dksbase->writeDataAsync<Vector_t>(p_ptr, &itsBunch_m->P[0], itsBunch_m->getLocalNum(), stream1);
514  dksbase->writeDataAsync<double>(dt_ptr, &itsBunch_m->dt[0], itsBunch_m->getLocalNum(), stream1);
515 }
516 
517 inline void ParallelTTracker::sendEfBf() {
518  dksbase->writeDataAsync<Vector_t>(Ef_ptr, &itsBunch_m->Ef[0],
519  itsBunch_m->getLocalNum(), stream1);
520  dksbase->writeDataAsync<Vector_t>(Bf_ptr, &itsBunch_m->Bf[0],
521  itsBunch_m->getLocalNum(), stream1);
522 }
523 
524 inline void ParallelTTracker::pushParticlesDKS(bool send) {
525 
526  //send data to the GPU
527  if (send)
528  sendRPdt();
529  //execute particle push on the GPU
530  dksbase->callParallelTTrackerPush(r_ptr, p_ptr, dt_ptr, itsBunch_m->getLocalNum(),
531  Physics::c, stream1);
532  //get particles back to CPU
533  dksbase->readData<Vector_t>(r_ptr, &itsBunch_m->R[0], itsBunch_m->getLocalNum(), stream1);
534 }
535 
536 inline
537 void ParallelTTracker::kickParticlesDKS() {
538  //send data to the GPU
539  sendEfBf();
540  sendRPdt();
541 
542  //sync device
543  dksbase->syncDevice();
544 
545  //execute kick on the GPU
546  dksbase->callParallelTTrackerKick(r_ptr, p_ptr, Ef_ptr, Bf_ptr, dt_ptr,
548  itsBunch_m->getLocalNum(), Physics::c, stream2);
549 
550  dksbase->syncDevice();
551 
552  //get data back from GPU
553  dksbase->readDataAsync<Vector_t>(p_ptr, &itsBunch_m->P[0], itsBunch_m->getLocalNum(), stream2);
554 
555 }
556 #endif
557 
558 #endif // OPAL_ParallelTTracker_HH
static int getNodes()
Definition: IpplInfo.cpp:773
ParticleAttrib< Vector_t > P
Definition: Source.h:12
Interface for solenoids.
Definition: RBend3D.h:39
void switchOffUnitlessPositions(bool use_dt_per_particle=false)
elements
Definition: IndexMap.cpp:141
void pushParticles(const BorisPusher &pusher)
WakeFunction * wakeFunction_m
virtual void visitTravelingWave(const TravelingWave &)
Apply the algorithm to a RFCavity.
void findStartPosition(const BorisPusher &pusher)
ParticleAttrib< Vector_t > Ef
Definition: TSVMeta.h:24
Interface for septum magnet.
Definition: Septum.h:11
void transformBunch(const CoordinateSystemTrafo &trafo)
virtual void visitRBend(const RBend &)
Apply the algorithm to a RBend.
Interface for electrostatic separator.
Definition: Separator.h:33
unsigned long totalParticlesInSimulation_m
virtual void execute()
Apply the algorithm to the top-level beamline.
virtual void visitBeamline(const Beamline &)
Apply the algorithm to a beam line.
Interface for beam position monitors.
Definition: Monitor.h:41
Interface for RF Quadrupole.
Definition: RFQuadrupole.h:30
virtual void visitRBend3D(const RBend3D &)
Apply the algorithm to a RBend.
virtual void visitLambertson(const Lambertson &)
Apply the algorithm to a Lambertson.
virtual void visitBeamBeam(const BeamBeam &)
Apply the algorithm to a BeamBeam.
IpplTimings::TimerRef timeIntegrationTimer1_m
Define the position of a misaligned element.
Definition: AlignWrapper.h:39
Interface for RF cavity.
Definition: ParallelPlate.h:36
void updateReferenceParticle(const BorisPusher &pusher)
Particle reference data.
Definition: PartData.h:38
Interface for general corrector.
Definition: Corrector.h:35
Abstract collimator.
IpplTimings::TimerRef fieldEvaluationTimer_m
Definition: RBend.h:73
unsigned int emissionSteps_m
Interface for beam diagnostics.
Definition: Diagnostic.h:32
virtual void visitRFQuadrupole(const RFQuadrupole &)
Apply the algorithm to a RFQuadrupole.
Interface for a marker.
Definition: Marker.h:32
unsigned int repartFreq_m
Interface for drift space.
Definition: Drift.h:33
size_t numParticlesInSimulation_m
IpplTimings::TimerRef BinRepartTimer_m
void updateReference(const BorisPusher &pusher)
virtual void visitSeptum(const Septum &)
Apply the algorithm to a Septum.
virtual void visitCorrector(const Corrector &)
Apply the algorithm to a Corrector.
void kickParticles(const BorisPusher &pusher)
void computeParticleMatterInteraction(IndexMap::value_t elements, OrbitThreader &oth)
Interface for general multipole.
Definition: Multipole.h:46
double zstart_m
where to start
T deg(T x)
Convert radians to degrees.
Definition: matheval.hpp:82
void updateRFElement(std::string elName, double maxPhi)
OpalBeamline itsOpalBeamline_m
void computeWakefield(IndexMap::value_t &elements)
Interface for probe.
Definition: Probe.h:16
void computeExternalFields(OrbitThreader &oth)
void writePhaseSpace(const long long step, bool psDump, bool statDump)
virtual void visitDegrader(const Degrader &)
Apply the algorithm to a Degrader.
virtual void visitMonitor(const Monitor &)
Apply the algorithm to a Monitor.
void applyFractionalStep(const BorisPusher &pusher, double tau)
double getQ() const
The constant charge per particle.
Definition: PartData.h:107
void computeSpaceChargeFields(unsigned long long step)
virtual void visitDiagnostic(const Diagnostic &)
Apply the algorithm to a Diagnostic.
IpplTimings::TimerRef WakeFieldTimer_m
virtual void visitDrift(const Drift &)
Apply the algorithm to a Drift.
void changeDT(bool backTrack=false)
virtual void visitCyclotronValley(const CyclotronValley &)
Apply the algorithm to a CyclotronValley.
IpplTimings::TimerRef timeIntegrationTimer2_m
constexpr double c
The velocity of light in m/s.
Definition: Physics.h:52
Interface for cyclotron collimator.
Definition: CCollimator.h:13
virtual void visitMarker(const Marker &)
Apply the algorithm to a Marker.
void kick(const Vector_t &R, Vector_t &P, const Vector_t &Ef, const Vector_t &Bf, const double &dt) const
Definition: BorisPusher.h:48
Abstract beam-beam interaction.
Definition: BeamBeam.h:37
void visit(const T &, BeamlineVisitor &, PartBunchBase< double, 3 > *)
Definition: OpalBeamline.h:98
Interface for RF cavity.
Definition: TravelingWave.h:37
StepSizeConfig stepSizes_m
Definition: SBend.h:68
Interface for cyclotron valley.
void operator=(const ParallelTTracker &)
void switchToUnitlessPositions(bool use_dt_per_particle=false)
virtual void visitSeparator(const Separator &)
Apply the algorithm to a Separator.
virtual void visitMultipole(const Multipole &)
Apply the algorithm to a Multipole.
virtual void visitSolenoid(const Solenoid &)
Apply the algorithm to a Solenoid.
void emitParticles(long long step)
Interface for solenoids.
Definition: Solenoid.h:36
virtual void visitBeamStripping(const BeamStripping &)
Apply the algorithm to a BeamStripping.
An abstract sequence of beam line components.
Definition: Beamline.h:37
const PartData itsReference
The reference information.
size_t getLocalNum() const
void handleRestartRun()
void timeIntegration2(BorisPusher &pusher)
virtual void visitRFCavity(const RFCavity &)
Apply the algorithm to a RFCavity.
ParticleAttrib< double > dt
virtual void visitSource(const Source &)
Apply the algorithm to a Solenoid.
Interface for RF cavity.
Definition: RFCavity.h:37
double getM() const
The constant mass per particle.
Definition: PartData.h:112
Abstract collimator.
Definition: Degrader.h:37
virtual void visitSBend(const SBend &)
Apply the algorithm to a SBend.
virtual void visitCCollimator(const CCollimator &)
Apply the algorithm to a collimator.
virtual void visitProbe(const Probe &)
Apply the algorithm to a Probe.
virtual void visitParallelPlate(const ParallelPlate &)
Apply the algorithm to a ParallelPlate.
ParticleAttrib< Vector_t > Bf
ParticlePos_t & R
void timeIntegration1(BorisPusher &pusher)
virtual void visitFlexibleCollimator(const FlexibleCollimator &)
Apply the algorithm to a flexible collimator.
Timing::TimerRef TimerRef
Definition: IpplTimings.h:176
DataSink * itsDataSink_m
std::set< ParticleMatterInteractionHandler * > activeParticleMatterInteractionHandlers_m
void push(Vector_t &R, const Vector_t &P, const double &dt) const
Definition: BorisPusher.h:116
virtual void visitMultipoleT(const MultipoleT &)
Apply the algorithm to a MultipoleT.
PartBunchBase< double, 3 > * itsBunch_m
The bunch of particles to be tracked.
Definition: Tracker.h:174
void autophaseCavities(const BorisPusher &pusher)
virtual void visitAlignWrapper(const AlignWrapper &)
Apply the algorithm to an align wrapper.
void selectDT(bool backTrack=false)
Interface for a Lambertson septum.
Definition: Lambertson.h:33
std::set< std::shared_ptr< Component > > value_t
Definition: IndexMap.h:18
void dumpStats(long long step, bool psDump, bool statDump)
Track particles or bunches.
Definition: Tracker.h:84