OPALX (Object Oriented Parallel Accelerator Library for Exascal) master (dc2a29eed580)
OPALX
Loading...
Searching...
No Matches
ParallelTracker.h
Go to the documentation of this file.
1
22#ifndef OPALX_ParallelTracker_HH
23#define OPALX_ParallelTracker_HH
24
26#include "Algorithms/Tracker.h"
29#include "Structure/DataSink.h"
30
31#include "BasicActions/Option.h"
32#include "Utilities/Options.h"
33
34#include "Physics/Physics.h"
35
36#include "Algorithms/IndexMap.h"
38
40#include "AbsBeamline/Drift.h"
42#include "AbsBeamline/Laser.h"
43#include "AbsBeamline/Marker.h"
44#include "AbsBeamline/Monitor.h"
47#include "AbsBeamline/RBend.h"
49#include "AbsBeamline/SBend.h"
53#include "Beamlines/Beamline.h"
56
57#include <list>
58#include <memory>
59#include <tuple>
60#include <vector>
61
62class ParticleMatterInteractionHandler;
63class PluginElement;
64
71class ParallelTracker : public Tracker {
72private:
76 double zstart_m;
77
80
82 unsigned long long repartFreq_m;
83 std::vector<std::vector<std::shared_ptr<SamplingBase>>>
85
86 // --- Timers ---
87 IpplTimings::TimerRef timeIntegrationTimer1_m;
88 IpplTimings::TimerRef timeIntegrationTimer2_m;
89 IpplTimings::TimerRef fieldEvaluationTimer_m;
90 IpplTimings::TimerRef WakeFieldTimer_m;
91 IpplTimings::TimerRef PluginElemTimer_m;
92 IpplTimings::TimerRef BinRepartTimer_m;
93 IpplTimings::TimerRef OrbThreader_m;
94
95public:
102 explicit ParallelTracker(const Beamline& bl, bool revBeam);
103
116 explicit ParallelTracker(
117 const Beamline& bl, PartBunch_t& bunch, DataSink* ds, bool revBeam,
118 const std::vector<unsigned long long>& maxSTEPS, double zstart,
119 const std::vector<double>& zstop, const std::vector<double>& dt,
120 const std::vector<std::vector<std::shared_ptr<SamplingBase>>>& emittingSamplers = {});
121
123 virtual ~ParallelTracker();
124
127 virtual void visitBeamline(const Beamline&);
128
130 virtual void visitComponent(const Component&);
131
134
136 virtual void visitDrift(const Drift&);
137
139 virtual void visitLaser(const Laser&);
140
142 virtual void visitMonitor(const Monitor&);
143
145 virtual void visitMarker(const Marker&);
146
148 virtual void visitMultipole(const Multipole&);
149
151 virtual void visitMultipoleT(const MultipoleT&);
152
154 virtual void visitRBend(const RBend&);
155
157 virtual void visitRFCavity(const RFCavity&);
158
160 virtual void visitSBend(const SBend&);
161
163 virtual void visitTravelingWave(const TravelingWave&);
164
166 virtual void visitSolenoid(const Solenoid&);
167
169 virtual void execute();
170
177
184
186 void timeIntegration1(BorisPusher& pusher);
187
189 void timeIntegration2(BorisPusher& pusher);
190
195 void evolveSpinTBMT();
196
201 void computeSpaceChargeFields(unsigned long long step);
202
206
210 void emitFromEmissionSources(double t, double dt);
211
213 size_t applyGlobalProcesses(double dt);
214
216 void resetFields();
217
219 void changeDT();
220
222 void setTime();
223
224private:
226 void updateReference(const BorisPusher& pusher);
227
229 void updateReferenceParticles(const BorisPusher& pusher);
230
233
240 void writePhaseSpace(const long long step, bool psDump, bool statDump);
241
248 void dumpStats(long long step, bool psDump, bool statDump);
249
251 void prepareSections();
252
254 void selectDT();
255
258
263 bool hasEndOfLineReached(const BoundingBox& globalBoundingBox);
264
266 void doBinaryRepartition();
267
269 size_t deleteInvalidParticles(bool activeOnly, Inform& m, const std::string& reason);
270
271public:
274
275private:
277 void activateEmittingContainers(double t);
278
285
290 void printInitialContainerRefs(Inform& m) const;
291
293 void findStartPositions(const BorisPusher& pusher);
294
296 void autophaseCavities(const BorisPusher& pusher);
297
303 void updateRFElement(std::string elName, double maxPhi);
304
307
309 void saveCavityPhases();
310
312 void restoreCavityPhases();
313};
314
318
319inline void ParallelTracker::visitDrift(const Drift& drift) {
320 itsOpalBeamline_m.visit(drift, *this, *itsBunch_m);
321}
322
323inline void ParallelTracker::visitMonitor(const Monitor& monitor) {
324 itsOpalBeamline_m.visit(monitor, *this, *itsBunch_m);
325}
326
327inline void ParallelTracker::visitMarker(const Marker& marker) {
328 itsOpalBeamline_m.visit(marker, *this, *itsBunch_m);
329}
330
332 itsOpalBeamline_m.visit(mult, *this, *itsBunch_m);
333}
334
336 itsOpalBeamline_m.visit(mult, *this, *itsBunch_m);
337}
338
339inline void ParallelTracker::visitRBend(const RBend& bend) {
340 itsOpalBeamline_m.visit(bend, *this, *itsBunch_m);
341}
342
345}
346
347inline void ParallelTracker::visitSBend(const SBend& bend) {
348 itsOpalBeamline_m.visit(bend, *this, *itsBunch_m);
349}
350
354
357}
358
359#endif // OPALX_ParallelTracker_HH
ippl::Vector< T, Dim > Vector_t
An abstract sequence of beam line components.
Definition Beamline.h:34
Component applying a constant accelerating electric field (Ex,Ey,Ez).
Interface for drift space.
Definition Drift.h:31
Passive OPALX laser element.
Definition Laser.h:25
Interface for a marker.
Definition Marker.h:31
Interface for general multipole.
Definition Multipole.h:30
void visit(const T &, BeamlineVisitor &, PartBunch_t &)
Implements the main time-based simulation loop for parallel tracking.
void activateEmittingContainers(double t)
Force-activate containers whose emitting samplers have not yet finished.
void saveCavityPhases()
Persist cavity phases to the data sink.
void updateRFElement(std::string elName, double maxPhi)
Set stored RF phase on the named cavity or traveling-wave element.
virtual void visitRBend(const RBend &)
Apply the algorithm to a rectangular bend.
double dtCurrentTrack_m
Global for the current track segment.
IpplTimings::TimerRef BinRepartTimer_m
void printInitialContainerRefs(Inform &m) const
Log reference state for each container at track start.
void restoreCavityPhases()
Restore cavity phases from a prior track or restart.
void evolveSpinTBMT()
Thomas-BMT spin precession across all active containers that store Pol. Must be called after external...
virtual void visitSolenoid(const Solenoid &)
Apply the algorithm to a solenoid.
size_t deleteInvalidParticles(bool activeOnly, Inform &m, const std::string &reason)
Delete particles marked invalid by the central per-container mask.
void prepareSections()
Accept beamline visitor, prepare sections, compute and save 3D lattice.
DataSink * itsDataSink_m
Borrowed beam statistics and phase-space output sink.
void resetFields()
Zero E and B on all active particle containers.
virtual void visitMonitor(const Monitor &)
Apply the algorithm to a monitor.
void doBinaryRepartition()
Trigger binary repartition for the field solver if configured.
void updateRefToLabCSTrafo()
Refresh each container's reference-to-lab transform from current state.
void kickParticles(const BorisPusher &pusher, PartBunch_t::ParticleContainer_t &pc)
Boris half-kick using E, B and per-particle dt on one container.
void autophaseCavities(const BorisPusher &pusher)
Autophase TRAVELINGWAVE and RFCAVITY elements along the reference orbit.
void printRFPhases()
Print RF phases (debug/diagnostic hook).
IpplTimings::TimerRef WakeFieldTimer_m
IpplTimings::TimerRef OrbThreader_m
void emitFromEmissionSources(double t, double dt)
Emit macroparticles from configured samplers per container.
virtual void visitMultipole(const Multipole &)
Apply the algorithm to a multipole.
void setTime()
Reset per-particle dt views to the current global bunch dt.
virtual void visitDrift(const Drift &)
Apply the algorithm to a drift.
virtual void execute()
Run the main tracking loop until all step-size segments complete.
std::vector< std::vector< std::shared_ptr< SamplingBase > > > emittingSamplers_m
Per-container emitters.
virtual void visitMultipoleT(const MultipoleT &)
Apply the algorithm to a multipole (templated type).
void setOptionalVariables()
Load REPARTFREQ and related options from input.
unsigned long long repartFreq_m
Space-charge repartition period (steps); 0 disables it.
void updateReferenceParticles(const BorisPusher &pusher)
Advance reference positions/momenta through the beamline for one step.
virtual void visitSBend(const SBend &)
Apply the algorithm to a sector bend.
void writePhaseSpace(const long long step, bool psDump, bool statDump)
Write phase space and/or statistics when flags request it.
void timeIntegration2(BorisPusher &pusher)
Second half: kick then push all active containers.
size_t applyGlobalProcesses(double dt)
Apply global processes and return the global number of particles marked invalid.
virtual void visitBeamline(const Beamline &)
Visit the full beamline (iterates elements into OpalBeamline). Overrides DefaultVisitor.
virtual void visitLaser(const Laser &)
Reject laser tracking until dedicated laser tracking is implemented.
void selectDT()
Set global bunch dt to dtCurrentTrack_m.
double zstart_m
Path-length start position for the track (m).
void changeDT()
Set bunch dt from StepSizeConfig and copy to all container dt views.
size_t markBackwardParticlesAtSourcePlane()
Mark particles moving backward behind an active source/cathode plane.
void updateReference(const BorisPusher &pusher)
Update reference trajectories and lab/reference coordinate transforms.
IpplTimings::TimerRef timeIntegrationTimer2_m
void computeExternalFields(OrbitThreader &oth)
Apply external fields from elements intersecting each active container.
virtual void visitTravelingWave(const TravelingWave &)
Apply the algorithm to a traveling wave cavity.
void computeSpaceChargeFields(unsigned long long step)
Self-fields in beam frame (primary container); optional binary repartition.
bool hasEndOfLineReached(const BoundingBox &globalBoundingBox)
Whether tracking should stop at end-of-line (global reduction).
void findStartPositions(const BorisPusher &pusher)
Integrate references in time until path length reaches zstart_m.
bool globalEOL_m
End-of-line flag (e.g. orbit threader out of bounds).
void dumpStats(long long step, bool psDump, bool statDump)
Log per-container stats and trigger dumps according to dump flags.
virtual void visitComponent(const Component &)
Visit a generic component using the base tracker behavior.
IpplTimings::TimerRef timeIntegrationTimer1_m
IpplTimings::TimerRef PluginElemTimer_m
virtual void visitMarker(const Marker &)
Apply the algorithm to a marker.
void computeInitialBounds(Vector_t< double, 3 > &rmin, Vector_t< double, 3 > &rmax)
Union of per-container spatial bounds over MPI.
virtual void visitRFCavity(const RFCavity &)
Apply the algorithm to an RF cavity.
virtual void visitConstantEFieldCavity(const ConstantEFieldCavity &)
Apply the algorithm to a constant E-field cavity.
IpplTimings::TimerRef fieldEvaluationTimer_m
void pushParticles(const BorisPusher &pusher, PartBunch_t::ParticleContainer_t &pc)
Boris position push (unitless positions) on one container.
void timeIntegration1(BorisPusher &pusher)
First half of the leapfrog step: push all active containers.
StepSizeConfig stepSizes_m
virtual ~ParallelTracker()
Destructor; releases tracker resources.
OpalBeamline itsOpalBeamline_m
Cloned field elements and coordinate transforms.
Container for all per-particle (and per-simulation) fields tracked during OPALX tracking.
Abstract rectangular bend with straight body and curved reference path.
Definition RBend.h:16
Interface for standing wave cavities.
Definition RFCavity.h:34
Abstract sector bend with planar-arc body geometry.
Definition SBend.h:15
Abstract class for a solenoid magnet.
Definition Solenoid.h:33
PartBunch_t * itsBunch_m
The bunch of particles to be tracked. Borrowed; lifetime is managed by TrackRun.
Definition Tracker.h:119
Interface for traveling wave cavities.