OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
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 #include "AbsBeamline/Vacuum.h"
56 
57 #ifdef ENABLE_OPAL_FEL
58 #include "AbsBeamline/Undulator.h"
59 #endif
60 
61 #include "Beamlines/Beamline.h"
62 #include "Elements/OpalBeamline.h"
63 #include "Solvers/WakeFunction.h"
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 
105  // overwrite the execute-methode from DefaultVisitor
106  virtual void execute();
107 
109  // overwrite the execute-methode from DefaultVisitor
110  virtual void visitBeamline(const Beamline &);
111 
113  virtual void visitCCollimator(const CCollimator &);
114 
116  virtual void visitCorrector(const Corrector &);
117 
119  virtual void visitDegrader(const Degrader &);
120 
122  virtual void visitDrift(const Drift &);
123 
125  virtual void visitFlexibleCollimator(const FlexibleCollimator &);
126 
128  virtual void visitMarker(const Marker &);
129 
131  virtual void visitMonitor(const Monitor &);
132 
134  virtual void visitMultipole(const Multipole &);
135 
137  virtual void visitMultipoleT(const MultipoleT &);
138 
140  virtual void visitProbe(const Probe &);
141 
143  virtual void visitRBend(const RBend &);
144 
146  virtual void visitRBend3D(const RBend3D &);
147 
149  virtual void visitRFCavity(const RFCavity &);
150 
152  virtual void visitSBend(const SBend &);
153 
155  virtual void visitSeptum(const Septum &);
156 
158  virtual void visitSolenoid(const Solenoid &);
159 
161  virtual void visitSource(const Source &);
162 
164  virtual void visitTravelingWave(const TravelingWave &);
165 
166 #ifdef ENABLE_OPAL_FEL
168  virtual void visitUndulator(const Undulator &);
169 #endif
170 
172  virtual void visitVacuum(const Vacuum &);
173 
174 
175 private:
176 
177  // Not implemented.
181 
182  /******************** STATE VARIABLES ***********************************/
183 
185 
187 
189 
191 
193 
195 
196  double pathLength_m;
197 
199  double zstart_m;
200 
206 
208 
209  // This variable controls the minimal number of steps of emission (using bins)
210  // before we can merge the bins
212 
213  // The space charge solver crashes if we use less than ~10 particles.
214  // This variable controls the number of particles to be emitted before we use
215  // the space charge solver.
217 
218  // this variable controls the minimal number of steps until we repartition the particles
219  unsigned int repartFreq_m;
220 
221  unsigned int emissionSteps_m;
222 
224 
230 
231  std::set<ParticleMatterInteractionHandler*> activeParticleMatterInteractionHandlers_m;
233 
234  /********************** END VARIABLES ***********************************/
235 
236  void kickParticles(const BorisPusher &pusher);
237  void pushParticles(const BorisPusher &pusher);
238  void updateReferenceParticle(const BorisPusher &pusher);
239 
240  void writePhaseSpace(const long long step, bool psDump, bool statDump);
241 
242  /********** BEGIN AUTOPHSING STUFF **********/
243  void updateRFElement(std::string elName, double maxPhi);
245  void saveCavityPhases();
246  void restoreCavityPhases();
247  /************ END AUTOPHSING STUFF **********/
248 
249  void prepareSections();
250 
251  void timeIntegration1(BorisPusher & pusher);
252  void timeIntegration2(BorisPusher & pusher);
253  void selectDT(bool backTrack = false);
254  void changeDT(bool backTrack = false);
255  void emitParticles(long long step);
259 #ifdef ENABLE_OPAL_FEL
260  void computeUndulator(IndexMap::value_t &elements);
261 #endif
262  void computeSpaceChargeFields(unsigned long long step);
263  // void prepareOpalBeamlineSections();
264  void dumpStats(long long step, bool psDump, bool statDump);
265  void setOptionalVariables();
266  bool hasEndOfLineReached();
268  void prepareEmission();
269  void setTime();
270  void doBinaryRepartition();
271 
272  void transformBunch(const CoordinateSystemTrafo &trafo);
273 
274  void updateReference(const BorisPusher &pusher);
275  void updateRefToLabCSTrafo();
276  void applyFractionalStep(const BorisPusher &pusher, double tau);
277  void findStartPosition(const BorisPusher &pusher);
278  void autophaseCavities(const BorisPusher &pusher);
279 
281 };
282 
283 
285  itsOpalBeamline_m.visit(coll, *this, itsBunch_m);
286 }
287 
288 inline void ParallelTTracker::visitCorrector(const Corrector &corr) {
289  itsOpalBeamline_m.visit(corr, *this, itsBunch_m);
290 }
291 
294 }
295 
296 inline void ParallelTTracker::visitDrift(const Drift &drift) {
297  itsOpalBeamline_m.visit(drift, *this, itsBunch_m);
298 }
299 
301  itsOpalBeamline_m.visit(coll, *this, itsBunch_m);
302 }
303 
304 inline void ParallelTTracker::visitMarker(const Marker &marker) {
305  itsOpalBeamline_m.visit(marker, *this, itsBunch_m);
306 }
307 
308 inline void ParallelTTracker::visitMonitor(const Monitor &mon) {
309  itsOpalBeamline_m.visit(mon, *this, itsBunch_m);
310 }
311 
312 inline void ParallelTTracker::visitMultipole(const Multipole &mult) {
313  itsOpalBeamline_m.visit(mult, *this, itsBunch_m);
314 }
315 
317  itsOpalBeamline_m.visit(mult, *this, itsBunch_m);
318 }
319 
320 inline void ParallelTTracker::visitProbe(const Probe &prob) {
321  itsOpalBeamline_m.visit(prob, *this, itsBunch_m);
322 }
323 
324 inline void ParallelTTracker::visitRBend(const RBend &bend) {
325  itsOpalBeamline_m.visit(bend, *this, itsBunch_m);
326 }
327 
328 inline void ParallelTTracker::visitRBend3D(const RBend3D &bend) {
329  itsOpalBeamline_m.visit(bend, *this, itsBunch_m);
330 }
331 
333  itsOpalBeamline_m.visit(as, *this, itsBunch_m);
334 }
335 
336 inline void ParallelTTracker::visitSBend(const SBend &bend) {
337  itsOpalBeamline_m.visit(bend, *this, itsBunch_m);
338 }
339 
340 inline void ParallelTTracker::visitSeptum(const Septum &sept) {
341  itsOpalBeamline_m.visit(sept, *this, itsBunch_m);
342 }
343 
344 inline void ParallelTTracker::visitSolenoid(const Solenoid &solenoid) {
345  itsOpalBeamline_m.visit(solenoid, *this, itsBunch_m);
346 }
347 
348 inline void ParallelTTracker::visitSource(const Source &source) {
349  itsOpalBeamline_m.visit(source, *this, itsBunch_m);
350 }
351 
353  itsOpalBeamline_m.visit(as, *this, itsBunch_m);
354 }
355 
356 #ifdef ENABLE_OPAL_FEL
357 inline void ParallelTTracker::visitUndulator(const Undulator &u) {
359 }
360 #endif
361 
362 inline void ParallelTTracker::visitVacuum(const Vacuum &vac) {
363  itsOpalBeamline_m.visit(vac, *this, itsBunch_m);
364 }
365 
366 inline void ParallelTTracker::kickParticles(const BorisPusher &pusher) {
367  int localNum = itsBunch_m->getLocalNum();
368  for (int i = 0; i < localNum; ++i)
369  pusher.kick(itsBunch_m->R[i], itsBunch_m->P[i], itsBunch_m->Ef[i], itsBunch_m->Bf[i], itsBunch_m->dt[i]);
370 }
371 
372 inline void ParallelTTracker::pushParticles(const BorisPusher &pusher) {
374 
375  for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
376  pusher.push(itsBunch_m->R[i], itsBunch_m->P[i], itsBunch_m->dt[i]);
377  }
379 }
380 
381 #endif // OPAL_ParallelTTracker_HH
elements
Definition: IndexMap.cpp:163
T deg(T x)
Convert radians to degrees.
Definition: matheval.hpp:74
ParticleAttrib< Vector_t > Ef
ParticlePos_t & R
size_t getLocalNum() const
ParticleAttrib< Vector_t > P
void switchToUnitlessPositions(bool use_dt_per_particle=false)
ParticleAttrib< double > dt
ParticleAttrib< Vector_t > Bf
void switchOffUnitlessPositions(bool use_dt_per_particle=false)
std::set< std::shared_ptr< Component > > value_t
Definition: IndexMap.h:47
virtual void visitRFCavity(const RFCavity &)
Apply the algorithm to a RF cavity.
void selectDT(bool backTrack=false)
virtual void visitCCollimator(const CCollimator &)
Apply the algorithm to a collimator.
void computeExternalFields(OrbitThreader &oth)
virtual void visitRBend(const RBend &)
Apply the algorithm to a rectangular bend.
IpplTimings::TimerRef timeIntegrationTimer2_m
IpplTimings::TimerRef WakeFieldTimer_m
virtual void visitFlexibleCollimator(const FlexibleCollimator &)
Apply the algorithm to a flexible collimator.
OpalBeamline itsOpalBeamline_m
double zstart_m
where to start
void transformBunch(const CoordinateSystemTrafo &trafo)
void autophaseCavities(const BorisPusher &pusher)
void computeWakefield(IndexMap::value_t &elements)
void handleRestartRun()
void computeSpaceChargeFields(unsigned long long step)
virtual void visitTravelingWave(const TravelingWave &)
Apply the algorithm to a traveling wave.
virtual void visitBeamline(const Beamline &)
Apply the algorithm to a beam line.
void timeIntegration2(BorisPusher &pusher)
DataSink * itsDataSink_m
virtual void visitRBend3D(const RBend3D &)
Apply the algorithm to a rectangular bend.
void updateReferenceParticle(const BorisPusher &pusher)
std::set< ParticleMatterInteractionHandler * > activeParticleMatterInteractionHandlers_m
WakeFunction * wakeFunction_m
size_t numParticlesInSimulation_m
virtual void visitVacuum(const Vacuum &)
Apply the algorithm to a vacuum space.
virtual void execute()
Apply the algorithm to the top-level beamline.
virtual void visitSource(const Source &)
Apply the algorithm to a source.
IpplTimings::TimerRef timeIntegrationTimer1_m
void pushParticles(const BorisPusher &pusher)
void dumpStats(long long step, bool psDump, bool statDump)
void changeDT(bool backTrack=false)
virtual void visitMarker(const Marker &)
Apply the algorithm to a marker.
ParallelTTracker(const ParallelTTracker &)
void emitParticles(long long step)
void findStartPosition(const BorisPusher &pusher)
void computeParticleMatterInteraction(IndexMap::value_t elements, OrbitThreader &oth)
void operator=(const ParallelTTracker &)
unsigned int repartFreq_m
void updateRFElement(std::string elName, double maxPhi)
void timeIntegration1(BorisPusher &pusher)
IpplTimings::TimerRef BinRepartTimer_m
virtual void visitDegrader(const Degrader &)
Apply the algorithm to a degrader.
virtual void visitSBend(const SBend &)
Apply the algorithm to a sector bend.
void writePhaseSpace(const long long step, bool psDump, bool statDump)
virtual void visitSolenoid(const Solenoid &)
Apply the algorithm to a solenoid.
virtual void visitMultipoleT(const MultipoleT &)
Apply the algorithm to an arbitrary multipole.
StepSizeConfig stepSizes_m
virtual void visitDrift(const Drift &)
Apply the algorithm to a drift space.
IpplTimings::TimerRef fieldEvaluationTimer_m
virtual void visitMonitor(const Monitor &)
Apply the algorithm to a beam position monitor.
void applyFractionalStep(const BorisPusher &pusher, double tau)
virtual void visitSeptum(const Septum &)
Apply the algorithm to a septum.
void kickParticles(const BorisPusher &pusher)
void updateReference(const BorisPusher &pusher)
virtual void visitProbe(const Probe &)
Apply the algorithm to a probe.
unsigned int emissionSteps_m
virtual void visitCorrector(const Corrector &)
Apply the algorithm to a closed orbit corrector.
virtual void visitMultipole(const Multipole &)
Apply the algorithm to a multipole.
Interface for general corrector.
Definition: Corrector.h:35
Interface for drift space.
Definition: Drift.h:33
Interface for a marker.
Definition: Marker.h:32
Interface for general multipole.
Definition: Multipole.h:47
Definition: Probe.h:28
Definition: RBend.h:58
Interface for solenoids.
Definition: RBend3D.h:39
Interface for RF cavity.
Definition: RFCavity.h:37
Definition: SBend.h:68
Definition: Septum.h:23
Interface for solenoids.
Definition: Solenoid.h:36
Definition: Source.h:30
Interface for Traveling Wave.
Definition: TravelingWave.h:37
Definition: Vacuum.h:67
Particle reference data.
Definition: PartData.h:35
PartBunchBase< double, 3 > * itsBunch_m
The bunch of particles to be tracked.
Definition: Tracker.h:151
An abstract sequence of beam line components.
Definition: Beamline.h:34
void visit(const T &, BeamlineVisitor &, PartBunchBase< double, 3 > *)
Definition: OpalBeamline.h:105
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
void push(Vector_t &R, const Vector_t &P, const double &dt) const
Definition: BorisPusher.h:131
Timing::TimerRef TimerRef
Definition: IpplTimings.h:176