OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
ParallelCyclotronTracker.h
Go to the documentation of this file.
1//
2// Class ParallelCyclotronTracker
3// Tracker for OPAL-Cycl
4//
5// Copyright (c) 2007 - 2014, Jianjun Yang, Andreas Adelmann and Matthias Toggweiler,
6// Paul Scherrer Institut, Villigen PSI, Switzerland
7// Copyright (c) 2014, Daniel Winklehner, MIT, Cambridge, MA, USA
8// Copyright (c) 2012 - 2022, Paul Scherrer Institut, Villigen PSI, Switzerland
9// All rights reserved
10//
11// Implemented as part of the PhD thesis
12// "Beam dynamics in high intensity cyclotrons including neighboring bunch effects"
13// and the paper
14// "Beam dynamics in high intensity cyclotrons including neighboring bunch effects:
15// Model, implementation, and application"
16// (https://journals.aps.org/prab/pdf/10.1103/PhysRevSTAB.13.064201)
17//
18// This file is part of OPAL.
19//
20// OPAL is free software: you can redistribute it and/or modify
21// it under the terms of the GNU General Public License as published by
22// the Free Software Foundation, either version 3 of the License, or
23// (at your option) any later version.
24//
25// You should have received a copy of the GNU General Public License
26// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
27//
28#ifndef OPAL_ParallelCyclotronTracker_HH
29#define OPAL_ParallelCyclotronTracker_HH
30
34#include "Algorithms/Tracker.h"
35#include "Steppers/Steppers.h"
36
37#include <memory>
38#include <tuple>
39#include <vector>
40
41class DataSink;
42class PluginElement;
43class LossDataSink;
44
45template <class T, unsigned Dim>
46class PartBunchBase;
47
54
56
57public:
58 typedef std::vector<double> dvector_t;
59 typedef std::vector<int> ivector_t;
60 typedef std::pair<double[8], Component*> element_pair;
61 typedef std::pair<ElementType, element_pair> type_pair;
62 typedef std::list<type_pair*> beamline_list;
63
65 // The beam line to be tracked is "bl".
66 // The particle reference data are taken from "data".
67 // The particle bunch tracked is taken from [b]bunch[/b].
68 // If [b]revBeam[/b] is true, the beam runs from s = C to s = 0.
69 // If [b]revTrack[/b] is true, we track against the beam.
71 const PartData& data, bool revBeam,
72 bool revTrack, int maxSTEPS,
73 Steppers::TimeIntegrator timeintegrator,
74 const int& numBunch,
75 const double& mbEta,
76 const double& mbPara,
77 const std::string& mbMode,
78 const std::string& mbBinning);
79
81
83 // overwrite the execute-methode from DefaultVisitor
84 virtual void execute();
85
87 // R - position
88 // P - momentum; note the field is returned in the lab frame
89 // t - time
90 // Efield - filled with values for the electric field
91 // Bfield - filled with values for the magnetic field
92 // Returns a boolean value, true if the particle was out of the nominal
93 // field map boundary, else false.
95 const Vector_t& P,
96 const double& t,
97 Vector_t& Efield,
98 Vector_t& Bfield);
99
101 // overwrite the execute-methode from DefaultVisitor
102 virtual void visitBeamline(const Beamline&);
103
105 virtual void visitCCollimator(const CCollimator&);
106
108 virtual void visitCorrector(const Corrector&);
109
111 virtual void visitCyclotron(const Cyclotron& cycl);
112
114 virtual void visitDegrader(const Degrader&);
115
117 virtual void visitDrift(const Drift&);
118
120 virtual void visitFlexibleCollimator(const FlexibleCollimator&);
121
123 virtual void visitMarker(const Marker&);
124
126 virtual void visitMonitor(const Monitor&);
127
129 virtual void visitMultipole(const Multipole&);
130
132 virtual void visitMultipoleT(const MultipoleT&);
133
135 virtual void visitOffset(const Offset&);
136
138 virtual void visitOutputPlane(const OutputPlane&);
139
141 virtual void visitProbe(const Probe&);
142
144 virtual void visitRBend(const RBend&);
145
147 virtual void visitRFCavity(const RFCavity&);
148
150 virtual void visitRing(const Ring& ring);
151
153 virtual void visitSBend(const SBend&);
154
156 virtual void visitSBend3D(const SBend3D&);
157
159 virtual void visitScalingFFAMagnet(const ScalingFFAMagnet& bend);
160
162 virtual void visitSeptum(const Septum&);
163
165 virtual void visitSolenoid(const Solenoid&);
166
168 virtual void visitStripper(const Stripper&);
169
171 virtual void visitVacuum(const Vacuum&);
172
174 virtual void visitVariableRFCavity(const VariableRFCavity& cav);
175
178
180 virtual void visitVerticalFFAMagnet(const VerticalFFAMagnet& bend);
181
183 inline void setLastDumpedStep(const int para) { lastDumpedStep_m = para; }
184
186 inline void setPr(double x) {referencePr = x;}
187 inline void setPt(double x) {referencePt = x;}
188 inline void setPz(double x) {referencePz = x;}
189 inline void setR(double x) {referenceR = x;}
190 inline void setTheta(double x) {referenceTheta = x;}
191 inline void setZ(double x) {referenceZ = x;}
192 inline void setBeGa(double x) {bega = x;}
193 inline void setPhi(double x) {referencePhi = x;}
194 inline void setPsi(double x) {referencePsi = x;}
195 inline void setPreviousH5Local(bool x) {previousH5Local = x;}
199
200 Ring* getRing() const {return opalRing_m;}
201private:
202 enum class TrackingMode: unsigned short {
207 };
208
209 // Not implemented.
213
215 std::list<Component*> myElements;
217 std::vector<PluginElement*> pluginElements_m;
218 std::vector<CavityCrossData> cavCrossDatas_m;
219
221
223
225
228
230 static Vector_t const xaxis;
231 static Vector_t const yaxis;
232 static Vector_t const zaxis;
233
235 double bega;
239 double referenceZ = 0.0;
240
243 double referencePz = 0.0;
245
248
249 bool spiral_flag = false;
250
252
254
257
258 // only used for multi-bunch mode
259 std::unique_ptr<MultiBunchHandler> mbHandler_m;
260
262
263 // for each bunch we have a path length
265 // Multi time step tracker
266 void MtsTracker();
267
268 void GenericTracker();
269 bool getFieldsAtPoint(const double& t, const size_t& Pindex, Vector_t& Efield, Vector_t& Bfield);
270
271 /*
272 Local Variables both used by the integration methods
273 */
274
275 long long step_m;
276 long long restartStep0_m;
277
279
280 // only for dumping
281 double azimuth_m;
283
284 /* only for dumping to stat-file --> make azimuth monotonically increasing
285 * @param theta computed azimuth [deg]
286 * @param prevAzimuth previous angle [deg]
287 * @param azimuth to be updated [deg]
288 * @param bunchNr in restart mode only --> to compute initial azimuth
289 */
290 void dumpAngle(const double& theta,
291 double& prevAzimuth,
292 double& azimuth,
293 const short& bunchNr = 0);
294
295 double computeRadius(const Vector_t& meanR) const;
296
297 void computePathLengthUpdate(std::vector<double>& dl, const double& dt);
298
299 // external field arrays for dumping
301
302 const int myNode_m;
303 const size_t initialLocalNum_m;
304 const size_t initialTotalNum_m;
305
307 std::vector<std::unique_ptr<std::ofstream> > outfTheta_m;
309 std::vector<double> azimuth_angle_m;
311 void openFiles(size_t numFiles, std::string fn);
312 void closeFiles();
315 std::ofstream outfTrackOrbit_m;
316
317 void buildupFieldList(double BcParameter[], ElementType elementType, Component* elptr);
318
319 // angle range [0~2PI) degree
320 double calculateAngle(double x, double y);
321 // angle range [-PI~PI) degree
322 double calculateAngle2(double x, double y);
323
324 bool checkGapCross(Vector_t Rold, Vector_t Rnew, RFCavity* rfcavity, double& DistOld);
325 bool RFkick(RFCavity* rfcavity, const double t, const double dt, const int Pindex);
326
327 bool getTunes(dvector_t& t, dvector_t& r, dvector_t& z, int lastTurn, double& nur, double& nuz);
328
335
336 /*
337 * @param bunchNr if <= -1 --> take all particles in simulation, if bunchNr > -1,
338 * take only particles of *this* bunch
339 */
340 Vector_t calcMeanR(short bunchNr = -1) const;
341
342 Vector_t calcMeanP() const;
343
344 void repartition(); // Do repartition between nodes if step_m is multiple of Options::repartFreq
345
346 // Transform the x- and y-parts of a particle attribute (position, momentum, fields) from the
347 // global reference frame to the local reference frame.
348
349 // phi is the angle of the bunch measured counter-clockwise from the positive x-axis.
350 void globalToLocal(ParticleAttrib<Vector_t>& vectorArray,
351 double phi, Vector_t const translationToGlobal = 0);
352
353 // Transform the x- and y-parts of a particle attribute (position, momentum, fields) from the
354 // local reference frame to the global reference frame.
355 void localToGlobal(ParticleAttrib<Vector_t>& vectorArray,
356 double phi, Vector_t const translationToGlobal = 0);
357
358 // Overloaded version of globalToLocal using a quaternion for 3D rotation
359 inline void globalToLocal(ParticleAttrib<Vector_t>& vectorArray,
360 Quaternion_t const quaternion,
361 Vector_t const meanR = Vector_t(0.0));
362
363 // Overloaded version of localToGlobal using a quaternion for 3D rotation
364 inline void localToGlobal(ParticleAttrib<Vector_t>& vectorArray,
365 Quaternion_t const quaternion,
366 Vector_t const meanR = Vector_t(0.0));
367
368 // Overloaded version of globalToLocal using phi and theta for pseudo 3D rotation
369 inline void globalToLocal(ParticleAttrib<Vector_t>& particleVectors,
370 double const phi, double const psi,
371 Vector_t const meanR = Vector_t(0.0));
372
373 // Overloaded version of localToGlobal using phi and theta for pseudo 3D rotation
374 inline void localToGlobal(ParticleAttrib<Vector_t>& particleVectors,
375 double const phi, double const psi,
376 Vector_t const meanR = Vector_t(0.0));
377
378 // Overloaded version of globalToLocal using phi and theta for pseudo 3D rotation, single vector
379 inline void globalToLocal(Vector_t& myVector,
380 double const phi, double const psi,
381 Vector_t const meanR = Vector_t(0.0));
382
383 // Overloaded version of localToGlobal using phi and theta for pseudo 3D rotation, single vector
384 inline void localToGlobal(Vector_t& myVector,
385 double const phi, double const psi,
386 Vector_t const meanR = Vector_t(0.0));
387
388 // Rotate the particles by an angle and axis defined in the quaternion (4-vector w,x,y,z)
389 inline void rotateWithQuaternion(ParticleAttrib<Vector_t>& vectorArray, Quaternion_t const quaternion);
390
391 // Returns the quaternion of the rotation from vector u to vector v
392 inline void getQuaternionTwoVectors(Vector_t u, Vector_t v, Quaternion_t& quaternion);
393
394 // Normalization of a quaternion
395 inline void normalizeQuaternion(Quaternion_t& quaternion);
396
397 // Normalization of a quaternion
398 inline void normalizeVector(Vector_t& vector);
399
400 // Some rotations
401 inline void rotateAroundZ(ParticleAttrib<Vector_t>& particleVectors, double const phi);
402 inline void rotateAroundX(ParticleAttrib<Vector_t>& particleVectors, double const psi);
403 inline void rotateAroundZ(Vector_t& myVector, double const phi);
404 inline void rotateAroundX(Vector_t& myVector, double const psi);
405
406 // Push particles for time h.
407 // Apply effects of RF Gap Crossings.
408 // Unit assumptions: [itsBunch_m->R] = m, [itsBunch_m->P] = 1, [h] = s, [c] = m/s, [itsBunch_m->getT()] = s
409 bool push(double h);
410 // Kick particles for time h
411 // The fields itsBunch_m->Bf, itsBunch_m->Ef are used to calculate the forces
412 bool kick(double h);
413
414 // Apply the trilogy half push - kick - half push,
415 // considering only external fields
416 void borisExternalFields(double h);
417
418 // apply the plugin elements: probe, collimator, stripper, septum
419 bool applyPluginElements(const double dt);
420
421 // destroy particles if they are marked as Bin=-1 in the plugin elements or out of global aperture
422 bool deleteParticle(bool = false);
423
424 void initTrackOrbitFile();
425
426 void singleParticleDump();
427
428 void bunchDumpStatData();
429
431
433
435
436 void setTimeStep();
437
438 void checkFileMomentum();
439
440 void checkNumPart(std::string s);
441
442 // we store a pointer explicitly to the Ring
444
445 // to save geometry losses
446 std::unique_ptr<LossDataSink> lossDs_m;
447
448 // If Ring is defined take the harmonic number from Ring; else use
449 // cyclotron
450 double getHarmonicNumber() const;
451
452 typedef std::function<bool(const double&,
453 const size_t&,
454 Vector_t&,
456
457 std::unique_ptr< Stepper<function_t> > itsStepper_mp;
458
465
467
469
470 // for change of coordinates
472
474 bool isTurnDone();
475
477 void update_m(double& t, const double& dt, const bool& finishedTurn);
478
482 std::tuple<double, double, double> initializeTracking_m();
483
484 void finalizeTracking_m(dvector_t& Ttime,
485 dvector_t& Tdeltr,
486 dvector_t& Tdeltz,
487 ivector_t& TturnNumber);
488
489 void setTrackingMode();
490
491 void seoMode_m(double& t, const double dt, bool& finishedTurn,
492 dvector_t& Ttime, dvector_t& Tdeltr,
493 dvector_t& Tdeltz, ivector_t& TturnNumber);
494
495 void singleMode_m(double& t, const double dt, bool& finishedTurn, double& oldReferenceTheta);
496
497 void bunchMode_m(double& t, const double dt, bool& finishedTurn);
498
499 void gapCrossKick_m(size_t i, double t, double dt,
500 const Vector_t& Rold, const Vector_t& Pold);
501
502
503 inline void dumpAzimuthAngles_m(const double& t,
504 const Vector_t& R,
505 const Vector_t& P,
506 const double& oldReferenceTheta,
507 const double& temp_meanTheta);
508
509 inline void dumpThetaEachTurn_m(const double& t,
510 const Vector_t& R,
511 const Vector_t& P,
512 const double& temp_meanTheta,
513 bool& finishedTurn);
514
516
517 bool computeExternalFields_m(const size_t& i,
518 const double& t,
519 Vector_t& Efield,
520 Vector_t& Bfield);
521
522 void injectBunch(bool& flagTransition);
523
524 void saveInjectValues();
525
526 bool isMultiBunch() const;
527
528 bool hasMultiBunch() const;
529
530 /*
531 * @param dt time step in s
532 */
533 void updatePathLength(const double& dt);
534
535 /*
536 * @param dt time step in s
537 */
538 void updateTime(const double& dt);
539
541
549 void initPathLength();
550};
551
558inline
560 double thetaXY = std::atan2(y, x);
561
562 if (thetaXY < 0) return thetaXY + Physics::two_pi;
563 return thetaXY;
564}
565
572inline
574 return std::atan2(y,x);
575}
576
577
578inline
580 return ( (mbHandler_m != nullptr) && itsBunch_m->weHaveBins() );
581}
582
583
584inline
586 return ( isMultiBunch() && mbHandler_m->getNumBunch() > 1 );
587}
588
589#endif // OPAL_ParallelCyclotronTracker_HH
ElementType
Definition ElementBase.h:88
Quaternion Quaternion_t
Definition Quaternion.h:42
boost::numeric::ublas::matrix< double > matrix_t
Definition BoostMatrix.h:23
constexpr double two_pi
The value of.
Definition Physics.h:33
TimeIntegrator
Definition Steppers.h:25
void computePathLengthUpdate(std::vector< double > &dl, const double &dt)
void setLastDumpedStep(const int para)
set last dumped step
std::vector< PluginElement * > pluginElements_m
void dumpAzimuthAngles_m(const double &t, const Vector_t &R, const Vector_t &P, const double &oldReferenceTheta, const double &temp_meanTheta)
IpplTimings::TimerRef IntegrationTimer_m
IpplTimings::TimerRef PluginElemTimer_m
virtual void visitDegrader(const Degrader &)
Apply the algorithm to a degrader.
virtual void visitRBend(const RBend &)
Apply the algorithm to a rectangular bend.
static Vector_t const xaxis
The positive axes unit vectors.
virtual void visitProbe(const Probe &)
Apply the algorithm to a probe.
virtual void visitVacuum(const Vacuum &)
Apply the algorithm to a vacuum space.
void setPr(double x)
Method for restart.
void operator=(const ParallelCyclotronTracker &)
std::function< bool(const double &, const size_t &, Vector_t &, Vector_t &)> function_t
bool computeExternalFields_m(const Vector_t &R, const Vector_t &P, const double &t, Vector_t &Efield, Vector_t &Bfield)
Calculate the field map at an arbitrary point.
std::vector< double > azimuth_angle_m
the different azimuthal angles for the outfTheta_m output files
std::unique_ptr< MultiBunchHandler > mbHandler_m
void injectBunch(bool &flagTransition)
std::pair< double[8], Component * > element_pair
void rotateAroundZ(ParticleAttrib< Vector_t > &particleVectors, double const phi)
virtual void visitStripper(const Stripper &)
Apply the algorithm to a particle stripper.
void normalizeQuaternion(Quaternion_t &quaternion)
virtual void visitCCollimator(const CCollimator &)
Apply the algorithm to a collimator.
void singleMode_m(double &t, const double dt, bool &finishedTurn, double &oldReferenceTheta)
void update_m(double &t, const double &dt, const bool &finishedTurn)
Update time and path length, write to output files.
std::vector< std::unique_ptr< std::ofstream > > outfTheta_m
output coordinates at different azimuthal angles and one after every turn
void rotateAroundX(ParticleAttrib< Vector_t > &particleVectors, double const psi)
double calculateAngle(double x, double y)
bool isTurnDone()
Check if turn done.
virtual void visitBeamline(const Beamline &)
Apply the algorithm to a beam line.
bool checkGapCross(Vector_t Rold, Vector_t Rnew, RFCavity *rfcavity, double &DistOld)
virtual void visitOffset(const Offset &)
Apply the algorithm to a offset (placement).
IpplTimings::TimerRef DelParticleTimer_m
bool RFkick(RFCavity *rfcavity, const double t, const double dt, const int Pindex)
virtual void visitScalingFFAMagnet(const ScalingFFAMagnet &bend)
Apply the algorithm to a scaling FFA magnet.
double computeRadius(const Vector_t &meanR) const
void updatePathLength(const double &dt)
void openFiles(size_t numFiles, std::string fn)
@ open / close output coordinate files
virtual void visitSBend3D(const SBend3D &)
Apply the algorithm to a sector bend with 3D field map.
double bega
The reference variables.
void dumpAngle(const double &theta, double &prevAzimuth, double &azimuth, const short &bunchNr=0)
bool applyPluginElements(const double dt)
std::list< type_pair * > beamline_list
virtual void visitCorrector(const Corrector &)
Apply the algorithm to a closed orbit corrector.
std::tuple< double, double, double > initializeTracking_m()
std::pair< ElementType, element_pair > type_pair
virtual void visitRing(const Ring &ring)
Apply the algorithm to a ring.
IpplTimings::TimerRef BinRepartTimer_m
void globalToLocal(ParticleAttrib< Vector_t > &vectorArray, double phi, Vector_t const translationToGlobal=0)
ParallelCyclotronTracker(const Beamline &bl, PartBunchBase< double, 3 > *bunch, DataSink &ds, const PartData &data, bool revBeam, bool revTrack, int maxSTEPS, Steppers::TimeIntegrator timeintegrator, const int &numBunch, const double &mbEta, const double &mbPara, const std::string &mbMode, const std::string &mbBinning)
Constructor.
bool getFieldsAtPoint(const double &t, const size_t &Pindex, Vector_t &Efield, Vector_t &Bfield)
std::vector< CavityCrossData > cavCrossDatas_m
virtual void visitDrift(const Drift &)
Apply the algorithm to a drift space.
std::unique_ptr< Stepper< function_t > > itsStepper_mp
struct ParallelCyclotronTracker::settings setup_m
void finalizeTracking_m(dvector_t &Ttime, dvector_t &Tdeltr, dvector_t &Tdeltz, ivector_t &TturnNumber)
virtual void visitVariableRFCavityFringeField(const VariableRFCavityFringeField &cav)
Apply the algorithm to a variable RF cavity with Fringe Field.
void gapCrossKick_m(size_t i, double t, double dt, const Vector_t &Rold, const Vector_t &Pold)
void getQuaternionTwoVectors(Vector_t u, Vector_t v, Quaternion_t &quaternion)
virtual void visitSBend(const SBend &)
Apply the algorithm to a sector bend.
virtual void visitMarker(const Marker &)
Apply the algorithm to a marker.
virtual void execute()
Apply the algorithm to the top-level beamline.
void dumpThetaEachTurn_m(const double &t, const Vector_t &R, const Vector_t &P, const double &temp_meanTheta, bool &finishedTurn)
void bunchMode_m(double &t, const double dt, bool &finishedTurn)
std::list< Component * > myElements
std::unique_ptr< LossDataSink > lossDs_m
virtual void visitMultipoleT(const MultipoleT &)
Apply the algorithm to an arbitrary multipole.
virtual void visitVariableRFCavity(const VariableRFCavity &cav)
Apply the algorithm to a variabel RF cavity.
void buildupFieldList(double BcParameter[], ElementType elementType, Component *elptr)
void rotateWithQuaternion(ParticleAttrib< Vector_t > &vectorArray, Quaternion_t const quaternion)
virtual void visitMultipole(const Multipole &)
Apply the algorithm to a multipole.
virtual void visitSolenoid(const Solenoid &)
Apply the algorithm to a solenoid.
Steppers::TimeIntegrator stepper_m
IpplTimings::TimerRef TransformTimer_m
ParallelCyclotronTracker(const ParallelCyclotronTracker &)
Vector_t calcMeanR(short bunchNr=-1) const
double calculateAngle2(double x, double y)
virtual void visitCyclotron(const Cyclotron &cycl)
Apply the algorithm to a cyclotron.
void localToGlobal(ParticleAttrib< Vector_t > &vectorArray, double phi, Vector_t const translationToGlobal=0)
virtual void visitRFCavity(const RFCavity &)
Apply the algorithm to a RF cavity.
bool getTunes(dvector_t &t, dvector_t &r, dvector_t &z, int lastTurn, double &nur, double &nuz)
virtual void visitFlexibleCollimator(const FlexibleCollimator &)
Apply the algorithm to a flexible collimator.
virtual void visitVerticalFFAMagnet(const VerticalFFAMagnet &bend)
Apply the algorithm to a vertical FFA magnet.
virtual void visitSeptum(const Septum &)
Apply the algorithm to a septum.
void seoMode_m(double &t, const double dt, bool &finishedTurn, dvector_t &Ttime, dvector_t &Tdeltr, dvector_t &Tdeltz, ivector_t &TturnNumber)
int maxSteps_m
The maximal number of steps the system is integrated.
virtual void visitOutputPlane(const OutputPlane &)
Apply the algorithm to a outputplane.
virtual void visitMonitor(const Monitor &)
Apply the algorithm to a beam position monitor.
Interface for a single beam element.
Definition Component.h:50
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
Ring describes a ring type geometry for tracking.
Definition Ring.h:53
Definition SBend.h:68
Interface for solenoids.
Definition Solenoid.h:36
Tracker(const Beamline &, const PartData &, bool backBeam, bool backTrack)
Constructor.
Definition Tracker.cpp:82
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
Timing::TimerRef TimerRef
Vektor< double, 3 > Vector_t