OPAL (Object Oriented Parallel Accelerator Library)  2024.1
OPAL
ParallelTTracker.h
Go to the documentation of this file.
1 //
2 // Class ParallelTTracker
3 // OPAL-T tracker.
4 // The visitor class for tracking particles with time as independent
5 // variable.
6 //
7 // Copyright (c) 200x - 2014, Christof Kraus, Paul Scherrer Institut, Villigen PSI, Switzerland
8 // 2015 - 2016, Christof Metzger-Kraus, Helmholtz-Zentrum Berlin, Germany
9 // 2017 - 2020, Christof Metzger-Kraus
10 // All rights reserved
11 //
12 // This file is part of OPAL.
13 //
14 // OPAL is free software: you can redistribute it and/or modify
15 // it under the terms of the GNU General Public License as published by
16 // the Free Software Foundation, either version 3 of the License, or
17 // (at your option) any later version.
18 //
19 // You should have received a copy of the GNU General Public License
20 // along with OPAL. If not, see <https://www.gnu.org/licenses/>.
21 //
22 #ifndef OPAL_ParallelTTracker_HH
23 #define OPAL_ParallelTTracker_HH
24 
25 #include "Algorithms/Tracker.h"
26 #include "Steppers/BorisPusher.h"
27 #include "Structure/DataSink.h"
29 
30 #include "BasicActions/Option.h"
31 #include "Utilities/Options.h"
32 
33 #include "Physics/Physics.h"
34 
36 #include "Algorithms/IndexMap.h"
38 #include "AbsBeamline/Corrector.h"
39 #include "AbsBeamline/Degrader.h"
40 #include "AbsBeamline/Drift.h"
43 #include "AbsBeamline/Marker.h"
44 #include "AbsBeamline/Monitor.h"
45 #include "AbsBeamline/Multipole.h"
46 #include "AbsBeamline/MultipoleT.h"
47 #include "AbsBeamline/Probe.h"
48 #include "AbsBeamline/RBend.h"
49 #include "AbsBeamline/RBend3D.h"
50 #include "AbsBeamline/RFCavity.h"
51 #include "AbsBeamline/SBend.h"
52 #include "AbsBeamline/Septum.h"
53 #include "AbsBeamline/Solenoid.h"
55 #ifdef ENABLE_OPAL_FEL
56 #include "AbsBeamline/Undulator.h"
57 #endif
58 #include "AbsBeamline/Vacuum.h"
59 
60 #include "Beamlines/Beamline.h"
61 #include "Elements/OpalBeamline.h"
62 #include "Solvers/WakeFunction.h"
63 
64 #include <list>
65 #include <vector>
66 
68 
69 class ParallelTTracker: public Tracker {
70 
71 public:
73  // The beam line to be tracked is "bl".
74  // The particle reference data are taken from "data".
75  // The particle bunch tracked is initially empty.
76  // If [b]revBeam[/b] is true, the beam runs from s = C to s = 0.
77  // If [b]revTrack[/b] is true, we track against the beam.
78  explicit ParallelTTracker(const Beamline &bl,
79  const PartData &data,
80  bool revBeam,
81  bool revTrack);
82 
84  // The beam line to be tracked is "bl".
85  // The particle reference data are taken from "data".
86  // The particle bunch tracked is taken from [b]bunch[/b].
87  // If [b]revBeam[/b] is true, the beam runs from s = C to s = 0.
88  // If [b]revTrack[/b] is true, we track against the beam.
89  explicit ParallelTTracker(const Beamline &bl,
91  DataSink &ds,
92  const PartData &data,
93  bool revBeam,
94  bool revTrack,
95  const std::vector<unsigned long long> &maxSTEPS,
96  double zstart,
97  const std::vector<double> &zstop,
98  const std::vector<double> &dt);
99 
100 
101  virtual ~ParallelTTracker();
102 
104  // overwrite the execute-methode from DefaultVisitor
105  virtual void execute();
106 
108  // overwrite the execute-methode from DefaultVisitor
109  virtual void visitBeamline(const Beamline &);
110 
112  virtual void visitCCollimator(const CCollimator &);
113 
115  virtual void visitCorrector(const Corrector &);
116 
118  virtual void visitDegrader(const Degrader &);
119 
121  virtual void visitDrift(const Drift &);
122 
124  virtual void visitFlexibleCollimator(const FlexibleCollimator &);
125 
127  virtual void visitMarker(const Marker &);
128 
130  virtual void visitMonitor(const Monitor &);
131 
133  virtual void visitMultipole(const Multipole &);
134 
136  virtual void visitMultipoleT(const MultipoleT &);
137 
139  virtual void visitProbe(const Probe &);
140 
142  virtual void visitRBend(const RBend &);
143 
145  virtual void visitRBend3D(const RBend3D &);
146 
148  virtual void visitRFCavity(const RFCavity &);
149 
151  virtual void visitSBend(const SBend &);
152 
154  virtual void visitSeptum(const Septum &);
155 
157  virtual void visitSolenoid(const Solenoid &);
158 
160  virtual void visitSource(const Source &);
161 
163  virtual void visitTravelingWave(const TravelingWave &);
164 
165 #ifdef ENABLE_OPAL_FEL
166  virtual void visitUndulator(const Undulator &);
168 #endif
169 
171  virtual void visitVacuum(const Vacuum &);
172 
173 private:
174 
175  // Not implemented.
178  void operator=(const ParallelTTracker &);
179 
180  /******************** STATE VARIABLES ***********************************/
181 
183 
185 
187 
189 
191 
193 
194  double pathLength_m;
195 
197  double zstart_m;
198 
204 
206 
207  // This variable controls the minimal number of steps of emission (using bins)
208  // before we can merge the bins
210 
211  // The space charge solver crashes if we use less than ~10 particles.
212  // This variable controls the number of particles to be emitted before we use
213  // the space charge solver.
215 
216  // this variable controls the minimal number of steps until we repartition the particles
217  unsigned int repartFreq_m;
218 
219  unsigned int emissionSteps_m;
220 
222 
228 
229  std::set<ParticleMatterInteractionHandler*> activeParticleMatterInteractionHandlers_m;
231 
232  /********************** END VARIABLES ***********************************/
233 
234  void kickParticles(const BorisPusher &pusher);
235  void pushParticles(const BorisPusher &pusher);
236  void updateReferenceParticle(const BorisPusher &pusher);
237 
238  void writePhaseSpace(const long long step, bool psDump, bool statDump);
239 
240  /********** BEGIN AUTOPHSING STUFF **********/
241  void updateRFElement(std::string elName, double maxPhi);
242  void printRFPhases();
243  void saveCavityPhases();
244  void restoreCavityPhases();
245  /************ END AUTOPHSING STUFF **********/
246 
247  void prepareSections();
248 
249  void timeIntegration1(BorisPusher & pusher);
250  void timeIntegration2(BorisPusher & pusher);
251  void selectDT(bool backTrack = false);
252  void changeDT(bool backTrack = false);
253  void emitParticles(long long step);
257 #ifdef ENABLE_OPAL_FEL
258  void computeUndulator(IndexMap::value_t &elements);
259 #endif
260  void computeSpaceChargeFields(unsigned long long step);
261  // void prepareOpalBeamlineSections();
262  void dumpStats(long long step, bool psDump, bool statDump);
263  void setOptionalVariables();
264  bool hasEndOfLineReached(const BoundingBox& globalBoundingBox);
265  void handleRestartRun();
266  void prepareEmission();
267  void setTime();
268  void doBinaryRepartition();
269 
270  void transformBunch(const CoordinateSystemTrafo &trafo);
271 
272  void updateReference(const BorisPusher &pusher);
273  void updateRefToLabCSTrafo();
274  void applyFractionalStep(const BorisPusher &pusher, double tau);
275  void findStartPosition(const BorisPusher &pusher);
276  void autophaseCavities(const BorisPusher &pusher);
277 
279 };
280 
281 
283  itsOpalBeamline_m.visit(coll, *this, itsBunch_m);
284 }
285 
286 inline void ParallelTTracker::visitCorrector(const Corrector &corr) {
287  itsOpalBeamline_m.visit(corr, *this, itsBunch_m);
288 }
289 
291  itsOpalBeamline_m.visit(deg, *this, itsBunch_m);
292 }
293 
294 inline void ParallelTTracker::visitDrift(const Drift &drift) {
295  itsOpalBeamline_m.visit(drift, *this, itsBunch_m);
296 }
297 
299  itsOpalBeamline_m.visit(coll, *this, itsBunch_m);
300 }
301 
302 inline void ParallelTTracker::visitMarker(const Marker &marker) {
303  itsOpalBeamline_m.visit(marker, *this, itsBunch_m);
304 }
305 
306 inline void ParallelTTracker::visitMonitor(const Monitor &mon) {
307  itsOpalBeamline_m.visit(mon, *this, itsBunch_m);
308 }
309 
310 inline void ParallelTTracker::visitMultipole(const Multipole &mult) {
311  itsOpalBeamline_m.visit(mult, *this, itsBunch_m);
312 }
313 
315  itsOpalBeamline_m.visit(mult, *this, itsBunch_m);
316 }
317 
318 inline void ParallelTTracker::visitProbe(const Probe &prob) {
319  itsOpalBeamline_m.visit(prob, *this, itsBunch_m);
320 }
321 
322 inline void ParallelTTracker::visitRBend(const RBend &bend) {
323  itsOpalBeamline_m.visit(bend, *this, itsBunch_m);
324 }
325 
326 inline void ParallelTTracker::visitRBend3D(const RBend3D &bend) {
327  itsOpalBeamline_m.visit(bend, *this, itsBunch_m);
328 }
329 
331  itsOpalBeamline_m.visit(as, *this, itsBunch_m);
332 }
333 
334 inline void ParallelTTracker::visitSBend(const SBend &bend) {
335  itsOpalBeamline_m.visit(bend, *this, itsBunch_m);
336 }
337 
338 inline void ParallelTTracker::visitSeptum(const Septum &sept) {
339  itsOpalBeamline_m.visit(sept, *this, itsBunch_m);
340 }
341 
342 inline void ParallelTTracker::visitSolenoid(const Solenoid &solenoid) {
343  itsOpalBeamline_m.visit(solenoid, *this, itsBunch_m);
344 }
345 
346 inline void ParallelTTracker::visitSource(const Source &source) {
347  itsOpalBeamline_m.visit(source, *this, itsBunch_m);
348 }
349 
351  itsOpalBeamline_m.visit(as, *this, itsBunch_m);
352 }
353 
354 #ifdef ENABLE_OPAL_FEL
355 inline void ParallelTTracker::visitUndulator(const Undulator &u) {
357 }
358 #endif
359 
360 inline void ParallelTTracker::visitVacuum(const Vacuum &vac) {
361  itsOpalBeamline_m.visit(vac, *this, itsBunch_m);
362 }
363 
364 inline void ParallelTTracker::kickParticles(const BorisPusher &pusher) {
365  int localNum = itsBunch_m->getLocalNum();
366  for (int i = 0; i < localNum; ++i)
367  pusher.kick(itsBunch_m->R[i], itsBunch_m->P[i], itsBunch_m->Ef[i], itsBunch_m->Bf[i], itsBunch_m->dt[i]);
368 }
369 
370 inline void ParallelTTracker::pushParticles(const BorisPusher &pusher) {
372 
373  for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
374  pusher.push(itsBunch_m->R[i], itsBunch_m->P[i], itsBunch_m->dt[i]);
375  }
377 }
378 
379 #endif // OPAL_ParallelTTracker_HH
std::set< ParticleMatterInteractionHandler * > activeParticleMatterInteractionHandlers_m
void timeIntegration2(BorisPusher &pusher)
unsigned int repartFreq_m
virtual void visitTravelingWave(const TravelingWave &)
Apply the algorithm to a traveling wave.
virtual void visitDegrader(const Degrader &)
Apply the algorithm to a degrader.
void writePhaseSpace(const long long step, bool psDump, bool statDump)
void findStartPosition(const BorisPusher &pusher)
DataSink * itsDataSink_m
IpplTimings::TimerRef WakeFieldTimer_m
void computeWakefield(IndexMap::value_t &elements)
void handleRestartRun()
Interface for a marker.
Definition: Marker.h:32
IpplTimings::TimerRef BinRepartTimer_m
WakeFunction * wakeFunction_m
ParticleAttrib< Vector_t > P
virtual void visitFlexibleCollimator(const FlexibleCollimator &)
Apply the algorithm to a flexible collimator.
void switchOffUnitlessPositions(bool use_dt_per_particle=false)
void selectDT(bool backTrack=false)
void computeExternalFields(OrbitThreader &oth)
IpplTimings::TimerRef timeIntegrationTimer2_m
std::set< std::shared_ptr< Component > > value_t
Definition: IndexMap.h:47
virtual void visitMarker(const Marker &)
Apply the algorithm to a marker.
ParticleAttrib< Vector_t > Ef
void operator=(const ParallelTTracker &)
double zstart_m
where to start
Interface for solenoids.
Definition: Solenoid.h:36
virtual void visitVacuum(const Vacuum &)
Apply the algorithm to a vacuum space.
Interface for solenoids.
Definition: RBend3D.h:39
void dumpStats(long long step, bool psDump, bool statDump)
void push(Vector_t &R, const Vector_t &P, const double &dt) const
Definition: BorisPusher.h:131
An abstract sequence of beam line components.
Definition: Beamline.h:34
Interface for general multipole.
Definition: Multipole.h:47
Interface for drift space.
Definition: Drift.h:33
Definition: SBend.h:68
Interface for general corrector.
Definition: Corrector.h:35
Definition: Vacuum.h:61
virtual void visitRBend3D(const RBend3D &)
Apply the algorithm to a rectangular bend.
void kickParticles(const BorisPusher &pusher)
Timing::TimerRef TimerRef
Definition: IpplTimings.h:176
void changeDT(bool backTrack=false)
virtual void visitRFCavity(const RFCavity &)
Apply the algorithm to a RF cavity.
ParticleAttrib< Vector_t > Bf
virtual void visitCCollimator(const CCollimator &)
Apply the algorithm to a collimator.
virtual void visitSolenoid(const Solenoid &)
Apply the algorithm to a solenoid.
void switchToUnitlessPositions(bool use_dt_per_particle=false)
virtual void visitCorrector(const Corrector &)
Apply the algorithm to a closed orbit corrector.
virtual void visitSource(const Source &)
Apply the algorithm to a source.
void updateReferenceParticle(const BorisPusher &pusher)
size_t getLocalNum() const
virtual void visitMultipoleT(const MultipoleT &)
Apply the algorithm to an arbitrary multipole.
virtual void visitMultipole(const Multipole &)
Apply the algorithm to a multipole.
elements
Definition: IndexMap.cpp:163
virtual void visitDrift(const Drift &)
Apply the algorithm to a drift space.
virtual void visitSBend(const SBend &)
Apply the algorithm to a sector bend.
void emitParticles(long long step)
void pushParticles(const BorisPusher &pusher)
ParticleAttrib< double > dt
Definition: RBend.h:58
void computeParticleMatterInteraction(IndexMap::value_t elements, OrbitThreader &oth)
void applyFractionalStep(const BorisPusher &pusher, double tau)
Definition: Source.h:30
bool hasEndOfLineReached(const BoundingBox &globalBoundingBox)
IpplTimings::TimerRef timeIntegrationTimer1_m
ParticlePos_t & R
void autophaseCavities(const BorisPusher &pusher)
void updateRFElement(std::string elName, double maxPhi)
virtual void visitRBend(const RBend &)
Apply the algorithm to a rectangular bend.
void timeIntegration1(BorisPusher &pusher)
T deg(T x)
Convert radians to degrees.
Definition: matheval.hpp:75
virtual void execute()
Apply the algorithm to the top-level beamline.
PartBunchBase< double, 3 > * itsBunch_m
The bunch of particles to be tracked.
Definition: Tracker.h:151
OpalBeamline itsOpalBeamline_m
void visit(const T &, BeamlineVisitor &, PartBunchBase< double, 3 > *)
Definition: OpalBeamline.h:105
Definition: Probe.h:28
void transformBunch(const CoordinateSystemTrafo &trafo)
virtual void visitMonitor(const Monitor &)
Apply the algorithm to a beam position monitor.
StepSizeConfig stepSizes_m
Definition: Septum.h:23
unsigned int emissionSteps_m
void updateReference(const BorisPusher &pusher)
void kick(const Vector_t &R, Vector_t &P, const Vector_t &Ef, const Vector_t &Bf, const double &dt) const
Definition: BorisPusher.h:65
virtual void visitBeamline(const Beamline &)
Apply the algorithm to a beam line.
virtual void visitProbe(const Probe &)
Apply the algorithm to a probe.
virtual void visitSeptum(const Septum &)
Apply the algorithm to a septum.
size_t numParticlesInSimulation_m
IpplTimings::TimerRef fieldEvaluationTimer_m
void computeSpaceChargeFields(unsigned long long step)