OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
ParallelCyclotronTracker.h
Go to the documentation of this file.
1 #ifndef OPAL_ParallelCyclotronTracker_HH
2 #define OPAL_ParallelCyclotronTracker_HH
3 
4 // ------------------------------------------------------------------------
5 // $RCSfile: ParallelCyclotronTracker.h,v $
6 // ------------------------------------------------------------------------
7 // $Revision: 1.1.2.1 $
8 // ------------------------------------------------------------------------
9 // Copyright: see Copyright.readme
10 // ------------------------------------------------------------------------
11 //
12 // Class: ParallelCyclotron
13 //
14 // ------------------------------------------------------------------------
15 //
16 // $Date: 2004/11/12 20:10:11 $
17 // $Author: adelmann $
18 //
19 // ------------------------------------------------------------------------
20 
21 #include "Algorithms/Tracker.h"
23 #include <vector>
24 #include <tuple>
25 #include <memory>
26 
27 #include "MultiBunchHandler.h"
28 #include "Steppers/Steppers.h"
29 
30 class DataSink;
31 class PluginElement;
32 class LossDataSink;
33 
34 template <class T, unsigned Dim>
35 class PartBunchBase;
36 
37 // Class ParallelCyclotronTracker
38 // ------------------------------------------------------------------------
39 
42  double sinAzimuth;
43  double cosAzimuth;
45 };
46 
48 
49 public:
50 
51  enum class MODE {
52  UNDEFINED = -1,
53  SINGLE = 0,
54  SEO = 1,
55  BUNCH = 2
56  };
57  typedef std::vector<double> dvector_t;
58  typedef std::vector<int> ivector_t;
59  typedef std::pair<double[8], Component *> element_pair;
60  typedef std::pair<ElementBase::ElementType, element_pair> type_pair;
61  typedef std::list<type_pair *> beamline_list;
62 
64  // The beam line to be tracked is "bl".
65  // The particle reference data are taken from "data".
66  // The particle bunch tracked is taken from [b]bunch[/b].
67  // If [b]revBeam[/b] is true, the beam runs from s = C to s = 0.
68  // If [b]revTrack[/b] is true, we track against the beam.
70  const PartData &data, bool revBeam,
71  bool revTrack, int maxSTEPS,
72  int timeIntegrator,
73  const int& numBunch,
74  const double& mbEta,
75  const double& mbPara,
76  const std::string& mbMode,
77  const std::string& mbBinning);
78 
79  virtual ~ParallelCyclotronTracker();
80 
82  virtual void visitRing(const Ring &ring);
83 
85  virtual void visitCyclotron(const Cyclotron &cycl);
86 
88  virtual void visitRFCavity(const RFCavity &);
89 
91  virtual void visitBeamBeam(const BeamBeam &);
92 
94  virtual void visitBeamStripping(const BeamStripping &);
95 
97  virtual void visitCCollimator(const CCollimator &);
98 
100  virtual void visitCorrector(const Corrector &);
101 
103  virtual void visitDegrader(const Degrader &);
104 
106  virtual void visitDiagnostic(const Diagnostic &);
107 
109  virtual void visitDrift(const Drift &);
110 
112  virtual void visitFlexibleCollimator(const FlexibleCollimator &);
113 
115  virtual void visitLambertson(const Lambertson &);
116 
118  virtual void visitMarker(const Marker &);
119 
121  virtual void visitMonitor(const Monitor &);
122 
124  virtual void visitMultipole(const Multipole &);
125 
127  virtual void visitMultipoleT (const MultipoleT &);
128 
130  virtual void visitMultipoleTStraight (const MultipoleTStraight &);
131 
134 
137 
139  virtual void visitOffset(const Offset &);
140 
142  virtual void visitProbe(const Probe &);
143 
145  virtual void visitRBend(const RBend &);
146 
148  virtual void visitRFQuadrupole(const RFQuadrupole &);
149 
151  virtual void visitSBend(const SBend &);
152 
154  virtual void visitSBend3D(const SBend3D &);
155 
157  virtual void visitScalingFFAMagnet(const ScalingFFAMagnet &bend);
158 
160  virtual void visitSeparator(const Separator &);
161 
163  virtual void visitSeptum(const Septum &);
164 
166  virtual void visitSolenoid(const Solenoid &);
167 
169  virtual void visitStripper(const Stripper &);
170 
172  virtual void visitParallelPlate(const ParallelPlate &);
173 
175  virtual void visitCyclotronValley(const CyclotronValley &);
176 
178  virtual void visitVariableRFCavity(const VariableRFCavity &cav);
179 
182  (const VariableRFCavityFringeField &cav);
183 
185  virtual void visitVerticalFFAMagnet(const VerticalFFAMagnet &bend);
186 
188  // overwrite the execute-methode from DefaultVisitor
189  virtual void execute();
190 
192  // overwrite the execute-methode from DefaultVisitor
193  virtual void visitBeamline(const Beamline &);
194 
196  inline void setLastDumpedStep(const int para) {lastDumpedStep_m = para ; }
197 
199  inline void setPr(double x) {referencePr = x;}
200  inline void setPt(double x) {referencePt = x;}
201  inline void setPz(double x) {referencePz = x;}
202  inline void setR(double x) {referenceR = x;}
203  inline void setTheta(double x) {referenceTheta = x;}
204  inline void setZ(double x) {referenceZ = x;}
205  inline void setBeGa(double x) {bega = x;}
206  inline void setPhi(double x) {referencePhi = x;}
207  inline void setPsi(double x) {referencePsi = x;}
208  inline void setPreviousH5Local(bool x) {previousH5Local = x;}
212 
213 private:
214 
215  // Not implemented.
218  void operator=(const ParallelCyclotronTracker &);
219 
221  std::list<Component *> myElements;
223  std::vector<PluginElement*> pluginElements_m;
224  std::vector<CavityCrossData> cavCrossDatas_m;
225 
227 
229 
231 
234 
236  static Vector_t const xaxis;
237  static Vector_t const yaxis;
238  static Vector_t const zaxis;
239 
241  double bega;
242  double referenceR;
244  double referenceZ = 0.0;
245 
246  double referencePr;
247  double referencePt;
248  double referencePz = 0.0;
250 
251  double referencePsi;
252  double referencePhi;
253 
254  bool spiral_flag = false;
255 
257 
259 
262 
263  // only used for multi-bunch mode
264  std::unique_ptr<MultiBunchHandler> mbHandler_m;
265 
267 
268  // for each bunch we have a path length
269  double pathLength_m;
270  // Multi time step tracker
271  void MtsTracker();
272 
273  void GenericTracker();
274  bool getFieldsAtPoint(const double &t, const size_t &Pindex, Vector_t &Efield, Vector_t &Bfield);
275 
276  /*
277  Local Variables both used by the integration methods
278  */
279 
280  long long step_m;
281  long long restartStep0_m;
282 
284 
285  // only for dumping
286  double azimuth_m;
288 
289  /* only for dumping to stat-file --> make azimuth monotonically increasing
290  * @param theta computed azimuth [deg]
291  * @param prevAzimuth previous angle [deg]
292  * @param azimuth to be updated [deg]
293  * @param bunchNr in restart mode only --> to compute initial azimuth
294  */
295  void dumpAngle(const double& theta,
296  double& prevAzimuth,
297  double& azimuth,
298  const short& bunchNr = 0);
299 
300  double computeRadius(const Vector_t& meanR) const;
301 
302  void computePathLengthUpdate(std::vector<double>& dl, const double& dt);
303 
304  // external field arrays for dumping
306 
307  const int myNode_m;
308  const size_t initialLocalNum_m;
309  const size_t initialTotalNum_m;
310 
312  std::vector<std::unique_ptr<std::ofstream> > outfTheta_m;
314  std::vector<double> azimuth_angle_m;
316  void openFiles(size_t numFiles, std::string fn);
317  void closeFiles();
320  std::ofstream outfTrackOrbit_m;
321 
322  void buildupFieldList(double BcParameter[], ElementBase::ElementType elementType, Component *elptr);
323 
324  // angle range [0~2PI) degree
325  double calculateAngle(double x, double y);
326  // angle range [-PI~PI) degree
327  double calculateAngle2(double x, double y);
328 
329  bool checkGapCross(Vector_t Rold, Vector_t Rnew, RFCavity * rfcavity, double &DistOld);
330  bool RFkick(RFCavity * rfcavity, const double t, const double dt, const int Pindex);
331 
332  bool getTunes(dvector_t &t, dvector_t &r, dvector_t &z, int lastTurn, double &nur, double &nuz);
333 
340 
341  /*
342  * @param bunchNr if <= -1 --> take all particles in simulation, if bunchNr > -1,
343  * take only particles of *this* bunch
344  */
345  Vector_t calcMeanR(short bunchNr = -1) const;
346 
347  Vector_t calcMeanP() const;
348 
349  void repartition(); // Do repartition between nodes if step_m is multiple of Options::repartFreq
350 
351  // Transform the x- and y-parts of a particle attribute (position, momentum, fields) from the
352  // global reference frame to the local reference frame.
353 
354  // phi is the angle of the bunch measured counter-clockwise from the positive x-axis.
355  void globalToLocal(ParticleAttrib<Vector_t> & vectorArray,
356  double phi, Vector_t const translationToGlobal = 0);
357 
358  // Transform the x- and y-parts of a particle attribute (position, momentum, fields) from the
359  // local reference frame to the global reference frame.
360  void localToGlobal(ParticleAttrib<Vector_t> & vectorArray,
361  double phi, Vector_t const translationToGlobal = 0);
362 
363  // Overloaded version of globalToLocal using a quaternion for 3D rotation
364  inline void globalToLocal(ParticleAttrib<Vector_t> & vectorArray,
365  Quaternion_t const quaternion,
366  Vector_t const meanR = Vector_t(0.0));
367 
368  // Overloaded version of localToGlobal using a quaternion for 3D rotation
369  inline void localToGlobal(ParticleAttrib<Vector_t> & vectorArray,
370  Quaternion_t const quaternion,
371  Vector_t const meanR = Vector_t(0.0));
372 
373  // Overloaded version of globalToLocal using phi and theta for pseudo 3D rotation
374  inline void globalToLocal(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 localToGlobal using phi and theta for pseudo 3D rotation
379  inline void localToGlobal(ParticleAttrib<Vector_t> & particleVectors,
380  double const phi, double const psi,
381  Vector_t const meanR = Vector_t(0.0));
382 
383  // Overloaded version of globalToLocal using phi and theta for pseudo 3D rotation, single vector
384  inline void globalToLocal(Vector_t & myVector,
385  double const phi, double const psi,
386  Vector_t const meanR = Vector_t(0.0));
387 
388  // Overloaded version of localToGlobal using phi and theta for pseudo 3D rotation, single vector
389  inline void localToGlobal(Vector_t & myVector,
390  double const phi, double const psi,
391  Vector_t const meanR = Vector_t(0.0));
392 
393  // Rotate the particles by an angle and axis defined in the quaternion (4-vector w,x,y,z)
394  inline void rotateWithQuaternion(ParticleAttrib<Vector_t> & vectorArray, Quaternion_t const quaternion);
395 
396  // Returns the quaternion of the rotation from vector u to vector v
397  inline void getQuaternionTwoVectors(Vector_t u, Vector_t v, Quaternion_t & quaternion);
398 
399  // Normalization of a quaternion
400  inline void normalizeQuaternion(Quaternion_t & quaternion);
401 
402  // Normalization of a quaternion
403  inline void normalizeVector(Vector_t & vector);
404 
405  // Some rotations
406  inline void rotateAroundZ(ParticleAttrib<Vector_t> & particleVectors, double const phi);
407  inline void rotateAroundX(ParticleAttrib<Vector_t> & particleVectors, double const psi);
408  inline void rotateAroundZ(Vector_t & myVector, double const phi);
409  inline void rotateAroundX(Vector_t & myVector, double const psi);
410 
411  // Push particles for time h.
412  // Apply effects of RF Gap Crossings.
413  // Unit assumptions: [itsBunch_m->R] = m, [itsBunch_m->P] = 1, [h] = s, [c] = m/s, [itsBunch_m->getT()] = s
414  bool push(double h);
415  // Kick particles for time h
416  // The fields itsBunch_m->Bf, itsBunch_m->Ef are used to calculate the forces
417  bool kick(double h);
418 
419  // Apply the trilogy half push - kick - half push,
420  // considering only external fields
421  void borisExternalFields(double h);
422 
423  // apply the plugin elements: probe, collimator, stripper, septum
424  bool applyPluginElements(const double dt);
425 
426  // destroy particles if they are marked as Bin=-1 in the plugin elements or out of global aperture
427  bool deleteParticle(bool = false);
428 
429  void initTrackOrbitFile();
430 
431  void singleParticleDump();
432 
433  void bunchDumpStatData();
434 
436 
438 
439  void initDistInGlobalFrame();
440 
441  void checkNumPart(std::string s);
442 
443  // we store a pointer explicitly to the Ring
445 
446  // to save geometry losses
447  std::unique_ptr<LossDataSink> lossDs_m;
448 
449  // If Ring is defined take the harmonic number from Ring; else use
450  // cyclotron
451  double getHarmonicNumber() const;
452 
453  typedef std::function<bool(const double&,
454  const size_t&,
455  Vector_t&,
457 
458  std::unique_ptr< Stepper<function_t> > itsStepper_mp;
459 
460  struct settings {
463  double deltaTheta;
465  } setup_m;
466 
468 
470 
472  bool isTurnDone();
473 
475  void update_m(double& t, const double& dt, const bool& finishedTurn);
476 
480  std::tuple<double, double, double> initializeTracking_m();
481 
482  void finalizeTracking_m(dvector_t& Ttime,
483  dvector_t& Tdeltr,
484  dvector_t& Tdeltz,
485  ivector_t& TturnNumber);
486 
487  void seoMode_m(double& t, const double dt, bool& finishedTurn,
488  dvector_t& Ttime, dvector_t& Tdeltr,
489  dvector_t& Tdeltz, ivector_t& TturnNumber);
490 
491  void singleMode_m(double& t, const double dt, bool& finishedTurn, double& oldReferenceTheta);
492 
493  void bunchMode_m(double& t, const double dt, bool& finishedTurn);
494 
495  void gapCrossKick_m(size_t i, double t, double dt,
496  const Vector_t& Rold, const Vector_t& Pold);
497 
498 
499  inline void dumpAzimuthAngles_m(const double& t,
500  const Vector_t& R,
501  const Vector_t& P,
502  const double& oldReferenceTheta,
503  const double& temp_meanTheta);
504 
505  inline void dumpThetaEachTurn_m(const double& t,
506  const Vector_t& R,
507  const Vector_t& P,
508  const double& temp_meanTheta,
509  bool& finishedTurn);
510 
512 
513  bool computeExternalFields_m(const size_t& i,
514  const double& t,
515  Vector_t& Efield,
516  Vector_t& Bfield);
517 
518  void injectBunch(bool& flagTransition);
519 
520  void saveInjectValues();
521 
522  bool isMultiBunch() const;
523 
524  bool hasMultiBunch() const;
525 
526  /*
527  * @param dt time step in ns
528  */
529  void updatePathLength(const double& dt);
530 
531  /*
532  * @param dt time step in ns
533  */
534  void updateTime(const double& dt);
535 
536  void updateAzimuthAndRadius();
537 
545  void initPathLength();
546 };
547 
554 inline
555 double ParallelCyclotronTracker::calculateAngle(double x, double y) {
556  double thetaXY = atan2(y, x);
557 
558  if (thetaXY < 0) return thetaXY + Physics::two_pi;
559  return thetaXY;
560 }
561 
568 inline
570 {
571  return atan2(y,x);
572 }
573 
574 
575 inline
577  return ( (mbHandler_m != nullptr) && itsBunch_m->weHaveBins() );
578 }
579 
580 
581 inline
583  return ( isMultiBunch() && mbHandler_m->getNumBunch() > 1 );
584 }
585 
586 #endif // OPAL_ParallelCyclotronTracker_HH
bool weHaveBins() const
void updatePathLength(const double &dt)
bool computeExternalFields_m(const size_t &i, const double &t, Vector_t &Efield, Vector_t &Bfield)
double calculateAngle(double x, double y)
virtual void visitOffset(const Offset &)
Apply the algorithm to a Offset.
void finalizeTracking_m(dvector_t &Ttime, dvector_t &Tdeltr, dvector_t &Tdeltz, ivector_t &TturnNumber)
void dumpAngle(const double &theta, double &prevAzimuth, double &azimuth, const short &bunchNr=0)
std::tuple< double, double, double > initializeTracking_m()
Definition: Offset.h:66
Definition: TSVMeta.h:24
std::vector< double > azimuth_angle_m
the different azimuthal angles for the outfTheta_m output files
Interface for septum magnet.
Definition: Septum.h:11
Interface for electrostatic separator.
Definition: Separator.h:33
virtual void visitMarker(const Marker &)
Apply the algorithm to a Marker.
Interface for a Cyclotron.
Definition: Cyclotron.h:91
std::list< Component * > myElements
virtual void visitCCollimator(const CCollimator &)
Apply the algorithm to a collimator.
Interface for beam position monitors.
Definition: Monitor.h:41
Interface for RF Quadrupole.
Definition: RFQuadrupole.h:30
constexpr double two_pi
The value of .
Definition: Physics.h:34
std::vector< PluginElement * > pluginElements_m
virtual void visitDegrader(const Degrader &)
Apply the algorithm to a Degrader.
virtual void visitVerticalFFAMagnet(const VerticalFFAMagnet &bend)
Apply the algorithm to a VerticalFFAMagnet.
Interface for RF cavity.
Definition: ParallelPlate.h:36
virtual void visitRFCavity(const RFCavity &)
Apply the algorithm to a RFCavity.
IpplTimings::TimerRef DumpTimer_m
Particle reference data.
Definition: PartData.h:38
Interface for general corrector.
Definition: Corrector.h:35
Abstract collimator.
Definition: RBend.h:73
Interface for beam diagnostics.
Definition: Diagnostic.h:32
double computeRadius(const Vector_t &meanR) const
std::list< type_pair * > beamline_list
Interface for a marker.
Definition: Marker.h:32
struct ParallelCyclotronTracker::settings setup_m
virtual void visitRFQuadrupole(const RFQuadrupole &)
Apply the algorithm to a RFQuadrupole.
IpplTimings::TimerRef IntegrationTimer_m
virtual void visitVariableRFCavity(const VariableRFCavity &cav)
Apply the algorithm to a VariabelRFCavity.
int maxSteps_m
The maximal number of steps the system is integrated.
void localToGlobal(ParticleAttrib< Vector_t > &vectorArray, double phi, Vector_t const translationToGlobal=0)
virtual void visitSeptum(const Septum &)
Apply the algorithm to a Septum.
void globalToLocal(ParticleAttrib< Vector_t > &vectorArray, double phi, Vector_t const translationToGlobal=0)
Interface for drift space.
Definition: Drift.h:33
virtual void visitDrift(const Drift &)
Apply the algorithm to a Drift.
void rotateAroundZ(ParticleAttrib< Vector_t > &particleVectors, double const phi)
void update_m(double &t, const double &dt, const bool &finishedTurn)
Update time and path length, write to output files.
void getQuaternionTwoVectors(Vector_t u, Vector_t v, Quaternion_t &quaternion)
virtual void visitSolenoid(const Solenoid &)
Apply the algorithm to a Solenoid.
Interface for general multipole.
Definition: Multipole.h:46
void gapCrossKick_m(size_t i, double t, double dt, const Vector_t &Rold, const Vector_t &Pold)
Vector_t calcMeanR(short bunchNr=-1) const
virtual void visitMultipoleTCurvedVarRadius(const MultipoleTCurvedVarRadius &)
Apply the algorithm to a MultipoleTCurvedVarRadius.
double bega
The reference variables.
void bunchMode_m(double &t, const double dt, bool &finishedTurn)
std::unique_ptr< Stepper< function_t > > itsStepper_mp
IpplTimings::TimerRef DelParticleTimer_m
Interface for probe.
Definition: Probe.h:16
virtual void visitProbe(const Probe &)
Apply the algorithm to a Probe.
virtual void execute()
Apply the algorithm to the top-level beamline.
void rotateAroundX(ParticleAttrib< Vector_t > &particleVectors, double const psi)
Class: DataSink.
Definition: OpalData.h:29
void setLastDumpedStep(const int para)
set last dumped step
virtual void visitFlexibleCollimator(const FlexibleCollimator &)
Apply the algorithm to a flexible collimator.
virtual void visitRBend(const RBend &)
Apply the algorithm to a RBend.
bool getFieldsAtPoint(const double &t, const size_t &Pindex, Vector_t &Efield, Vector_t &Bfield)
bool getTunes(dvector_t &t, dvector_t &r, dvector_t &z, int lastTurn, double &nur, double &nuz)
void injectBunch(bool &flagTransition)
virtual void visitMultipoleTCurvedConstRadius(const MultipoleTCurvedConstRadius &)
Apply the algorithm to a MultipoleTCurvedConstRadius.
virtual void visitMultipole(const Multipole &)
Apply the algorithm to a Multipole.
bool isTurnDone()
Check if turn done.
Interface for cyclotron collimator.
Definition: CCollimator.h:13
virtual void visitCyclotronValley(const CyclotronValley &)
Apply the algorithm to a CyclotronValley.it is empty for cyclotrontracker .
virtual void visitRing(const Ring &ring)
Apply the algorithm to an Ring.
virtual void visitScalingFFAMagnet(const ScalingFFAMagnet &bend)
Apply the algorithm to a ScalingFFAMagnet.
void computePathLengthUpdate(std::vector< double > &dl, const double &dt)
Abstract beam-beam interaction.
Definition: BeamBeam.h:37
Definition: SBend.h:68
IpplTimings::TimerRef TransformTimer_m
void dumpAzimuthAngles_m(const double &t, const Vector_t &R, const Vector_t &P, const double &oldReferenceTheta, const double &temp_meanTheta)
bool RFkick(RFCavity *rfcavity, const double t, const double dt, const int Pindex)
Interface for cyclotron valley.
PETE_TBTree< FnArcTan2, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > atan2(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
IpplTimings::TimerRef BinRepartTimer_m
Vektor< double, 3 > Vector_t
Definition: Vektor.h:6
std::unique_ptr< MultiBunchHandler > mbHandler_m
std::vector< std::unique_ptr< std::ofstream > > outfTheta_m
output coordinates at different azimuthal angles and one after every turn
virtual void visitCorrector(const Corrector &)
Apply the algorithm to a Corrector.
Interface for solenoids.
Definition: Solenoid.h:36
An abstract sequence of beam line components.
Definition: Beamline.h:37
std::pair< ElementBase::ElementType, element_pair > type_pair
virtual void visitVariableRFCavityFringeField(const VariableRFCavityFringeField &cav)
Apply the algorithm to a VariabelRFCavity.
void openFiles(size_t numFiles, std::string fn)
@ open / close output coordinate files
std::function< bool(const double &, const size_t &, Vector_t &, Vector_t &)> function_t
virtual void visitSBend3D(const SBend3D &)
Apply the algorithm to a SBend3D.
virtual void visitCyclotron(const Cyclotron &cycl)
Apply the algorithm to a Cyclotorn.
virtual void visitMultipoleT(const MultipoleT &)
Apply the algorithm to a MultipoleT.
void singleMode_m(double &t, const double dt, bool &finishedTurn, double &oldReferenceTheta)
IpplTimings::TimerRef PluginElemTimer_m
static Vector_t const xaxis
The positive axes unit vectors.
Interface for RF cavity.
Definition: RFCavity.h:37
virtual void visitBeamStripping(const BeamStripping &)
Apply the algorithm to a Beam Stripping.
Abstract collimator.
Definition: Degrader.h:37
virtual void visitBeamline(const Beamline &)
Apply the algorithm to a beam line.
std::vector< CavityCrossData > cavCrossDatas_m
void normalizeVector(Vector_t &vector)
virtual void visitBeamBeam(const BeamBeam &)
Apply the algorithm to a BeamBeam.
void setPr(double x)
Method for restart.
INTEGRATOR
Definition: Steppers.h:9
bool applyPluginElements(const double dt)
virtual void visitMonitor(const Monitor &)
Apply the algorithm to a Monitor.
Ring describes a ring type geometry for tracking.
Definition: Ring.h:64
virtual void visitSBend(const SBend &)
Apply the algorithm to a SBend.
virtual void visitSeparator(const Separator &)
Apply the algorithm to a Separator.
void rotateWithQuaternion(ParticleAttrib< Vector_t > &vectorArray, Quaternion_t const quaternion)
virtual void visitLambertson(const Lambertson &)
Apply the algorithm to a Lambertson.
Timing::TimerRef TimerRef
Definition: IpplTimings.h:176
Interface for a single beam element.
Definition: Component.h:51
void seoMode_m(double &t, const double dt, bool &finishedTurn, dvector_t &Ttime, dvector_t &Tdeltr, dvector_t &Tdeltz, ivector_t &TturnNumber)
void operator=(const ParallelCyclotronTracker &)
void dumpThetaEachTurn_m(const double &t, const Vector_t &R, const Vector_t &P, const double &temp_meanTheta, bool &finishedTurn)
bool checkGapCross(Vector_t Rold, Vector_t Rnew, RFCavity *rfcavity, double &DistOld)
void normalizeQuaternion(Quaternion_t &quaternion)
virtual void visitDiagnostic(const Diagnostic &)
Apply the algorithm to a Diagnostic.
std::pair< double[8], Component * > element_pair
virtual void visitMultipoleTStraight(const MultipoleTStraight &)
Apply the algorithm to a MultipoleTStraight.
void buildupFieldList(double BcParameter[], ElementBase::ElementType elementType, Component *elptr)
double calculateAngle2(double x, double y)
PartBunchBase< double, 3 > * itsBunch_m
The bunch of particles to be tracked.
Definition: Tracker.h:174
std::vector< double > dvector_t
virtual void visitStripper(const Stripper &)
Apply the algorithm to a charge stripper.
Interface for a Lambertson septum.
Definition: Lambertson.h:33
virtual void visitParallelPlate(const ParallelPlate &)
Apply the algorithm to a ParallelPlate, it is empty for cyclotrontracker .
std::unique_ptr< LossDataSink > lossDs_m
Track particles or bunches.
Definition: Tracker.h:84