OPALX (Object Oriented Parallel Accelerator Library for Exascale) MINIorX
OPALX
ParallelTracker.cpp
Go to the documentation of this file.
1//
2// Class ParallelTracker
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//
23
24#include <cfloat>
25#include <cmath>
26#include <fstream>
27#include <iomanip>
28#include <iostream>
29#include <limits>
30#include <sstream>
31#include <string>
32
33#include <boost/numeric/ublas/io.hpp>
34
38#include "BasicActions/Option.h"
39
40#include "Beamlines/Beamline.h"
43#include "Physics/Units.h"
44
48#include "Utilities/Options.h"
49#include "Utilities/Timer.h"
50#include "Utilities/Util.h"
52
53#include "AbsBeamline/Offset.h"
56
57extern Inform* gmsg;
58
59class PartData;
60
62 const Beamline& beamline, const PartData& reference, bool revBeam, bool revTrack)
63 : Tracker(beamline, reference, revBeam, revTrack),
64 itsDataSink_m(nullptr),
65 itsOpalBeamline_m(beamline.getOrigin3D(), beamline.getInitialDirection()),
66 opalRing_m(nullptr),
67 globalEOL_m(false),
68 wakeStatus_m(false),
69 wakeFunction_m(nullptr),
70 pathLength_m(0.0),
71 zstart_m(0.0),
72 dtCurrentTrack_m(0.0),
73 minStepforReBin_m(-1),
74 repartFreq_m(0),
75 emissionSteps_m(std::numeric_limits<unsigned int>::max()),
76 numParticlesInSimulation_m(0),
77 timeIntegrationTimer1_m(IpplTimings::getTimer("TIntegration1")),
78 timeIntegrationTimer2_m(IpplTimings::getTimer("TIntegration2")),
79 fieldEvaluationTimer_m(IpplTimings::getTimer("External field eval")),
80 PluginElemTimer_m(IpplTimings::getTimer("PluginElements")),
81 BinRepartTimer_m(IpplTimings::getTimer("Binaryrepart")),
82 OrbThreader_m(IpplTimings::getTimer("OrbThreader")) {
83}
84
86 const Beamline& beamline, PartBunch_t* bunch, DataSink& ds, const PartData& reference,
87 bool revBeam, bool revTrack, const std::vector<unsigned long long>& maxSteps, double zstart,
88 const std::vector<double>& zstop, const std::vector<double>& dt)
89 : Tracker(beamline, bunch, reference, revBeam, revTrack),
90 itsDataSink_m(&ds),
91 itsOpalBeamline_m(beamline.getOrigin3D(), beamline.getInitialDirection()),
92 opalRing_m(nullptr),
93 globalEOL_m(false),
94 wakeStatus_m(false),
95 wakeFunction_m(nullptr),
96 pathLength_m(0.0),
97 zstart_m(zstart),
98 dtCurrentTrack_m(0.0),
99 minStepforReBin_m(-1),
100 repartFreq_m(0),
101 emissionSteps_m(std::numeric_limits<unsigned int>::max()),
102 numParticlesInSimulation_m(0),
103 timeIntegrationTimer1_m(IpplTimings::getTimer("TIntegration1")),
104 timeIntegrationTimer2_m(IpplTimings::getTimer("TIntegration2")),
105 fieldEvaluationTimer_m(IpplTimings::getTimer("External field eval")),
106 BinRepartTimer_m(IpplTimings::getTimer("Binaryrepart")),
107 OrbThreader_m(IpplTimings::getTimer("OrbThreader")) {
108
109 for (unsigned int i = 0; i < zstop.size(); ++i) {
110 stepSizes_m.push_back(dt[i], zstop[i], maxSteps[i]);
111 }
112
115}
116
118}
119
121 *gmsg << "Adding ScalingFFAMagnet" << endl;
122 /*
123 if (opalRing_m != nullptr) {
124 ScalingFFAMagnet* newBend = bend.clone(); // setup the end field, if required
125 newBend->setupEndField();
126 opalRing_m->appendElement(*newBend);
127 } else {
128 throw OpalException("ParallelCyclotronTracker::visitScalingFFAMagnet",
129 "Need to define a RINGDEFINITION to use ScalingFFAMagnet element");
130 }
131 */
132}
133
135 double BcParameter[], ElementType elementType, Component* elptr) {
137
138 type_pair* localpair = new type_pair();
139 localpair->first = elementType;
140
141 for (int i = 0; i < 8; i++)
142 *(((localpair->second).first) + i) = *(BcParameter + i);
143
144 (localpair->second).second = elptr;
145
146 // always put cyclotron as the first element in the list.
147 if (elementType == ElementType::RING) {
148 sindex = FieldDimensions.begin();
149 } else {
150 sindex = FieldDimensions.end();
151 }
152 FieldDimensions.insert(sindex, localpair);
153}
154
155/*
156 * @param ring
157 */
158
160 *gmsg << "* ----------------------------- Ring ------------------------------------- *" << endl;
161
162 delete opalRing_m;
163
164 opalRing_m = dynamic_cast<Ring*>(ring.clone());
165
166 myElements.push_back(opalRing_m);
167
169
173
174 if (referenceTheta <= -180.0 || referenceTheta > 180.0) {
175 throw OpalException(
176 "Error in ParallelTracker::visitRing", "PHIINIT is out of [-180, 180)!");
177 }
178
179 referenceZ = 0.0;
180 referencePz = 0.0;
181
184
185 if (referencePtot < 0.0)
186 referencePt *= -1.0;
187
190
191 double BcParameter[8] = {}; // zero initialise array
192
194
195 // Finally print some diagnostic
196 *gmsg << "* Initial beam radius = " << referenceR << " [mm] " << endl;
197 *gmsg << "* Initial gamma = " << itsReference.getGamma() << endl;
198 *gmsg << "* Initial beta = " << itsReference.getBeta() << endl;
199 *gmsg << "* Total reference momentum = " << referencePtot << " [beta gamma]" << endl;
200 *gmsg << "* Reference azimuthal momentum = " << referencePt << " [beta gamma]" << endl;
201 *gmsg << "* Reference radial momentum = " << referencePr << " [beta gamma]" << endl;
202 *gmsg << "* " << opalRing_m->getSymmetry() << " fold field symmetry " << endl;
203 *gmsg << "* Harmonic number h = " << opalRing_m->getHarmonicNumber() << " " << endl;
204}
205
212 *gmsg << "Adding Vertical FFA Magnet" << endl;
213 if (opalRing_m != nullptr)
215 else
216 throw OpalException(
217 "ParallelCyclotronTracker::visitVerticalFFAMagnet",
218 "Need to define a RINGDEFINITION to use VerticalFFAMagnet element");
219}
220
222 const FlaggedBeamline* fbl = static_cast<const FlaggedBeamline*>(&bl);
223 if (fbl->getRelativeFlag()) {
224 OpalBeamline stash(fbl->getOrigin3D(), fbl->getInitialDirection());
225 stash.swap(itsOpalBeamline_m);
226 fbl->iterate(*this, false);
230 stash.swap(itsOpalBeamline_m);
231 } else {
232 fbl->iterate(*this, false);
233 }
234}
235
238 if (opalRing_m == nullptr)
239 throw OpalException(
240 "ParallelCylcotronTracker::visitOffset",
241 "Attempt to place an offset when Ring not defined");
242 Offset* offNonConst = const_cast<Offset*>(&off);
245}
246
247void ParallelTracker::updateRFElement(std::string elName, double maxPhase) {
250 cavities.insert(cavities.end(), travelingwaves.begin(), travelingwaves.end());
251
252 for (FieldList::iterator fit = cavities.begin(); fit != cavities.end(); ++fit) {
253 if ((*fit).getElement()->getName() == elName) {
254 RFCavity* element = static_cast<RFCavity*>((*fit).getElement().get());
255
256 element->setPhasem(maxPhase);
257 element->setAutophaseVeto();
258
259 *ippl::Info << "Restored cavity phase from the h5 file. Name: " << element->getName()
260 << ", phase: " << maxPhase << " rad" << endl;
261 return;
262 }
263 }
264}
265
268
269 bool flag = false;
270 for (PluginElement* element : pluginElements_m) {
271 bool tmp = element->check(itsBunch_m, turnnumber_m, itsBunch_m->getT(), dt);
272 flag |= tmp;
273
274 if (tmp) {
276 *gmsg << "* Total number of particles after PluginElement= "
277 << itsBunch_m->getTotalNum() << endl;
278 }
279 }
280
282 return flag;
283}
284
287}
288
291
292 if (OpalData::getInstance()->hasPriorTrack() || OpalData::getInstance()->inRestartRun()) {
295 for (; it < end; ++it) {
296 updateRFElement((*it).first, (*it).second);
297 }
298 }
299}
300
302 Inform msg("ParallelTracker ", *gmsg);
304 bool back_track = false;
305
307 const double globalTimeShift =
310
311 // the time step needs to be positive in the setup
312 itsBunch_m->setdT(std::abs(itsBunch_m->getdT()));
314
315 if (OpalData::getInstance()->hasPriorTrack() || OpalData::getInstance()->inRestartRun()) {
317 }
318
320
321 double minTimeStep = stepSizes_m.getMinTimeStep();
322
324
326
327 itsBunch_m->toLabTrafo_m = beamlineToLab;
328
329 itsBunch_m->RefPartR_m = beamlineToLab.transformTo(itsBunch_m->getParticleContainer()->getMeanR());
330 itsBunch_m->RefPartP_m = beamlineToLab.rotateTo(itsBunch_m->getParticleContainer()->getMeanP());
331
332 if (itsBunch_m->getTotalNum() > 0) {
333 if (zstart_m > pathLength_m) {
334 findStartPosition(pusher);
335 }
336
338 }
339
341
342 Vector_t<double, 3> rmin(0.0), rmax(0.0);
343 if (itsBunch_m->getTotalNum() > 0) {
344 itsBunch_m->get_bounds(rmin, rmax);
345 }
346
347 *gmsg << "ParallelTrack: momentum= " << itsBunch_m->getParticleContainer()->getMeanP()(2) << endl;
348 *gmsg << "itsBunch_m->RefPartR_m= " << itsBunch_m->RefPartR_m << endl;
349 *gmsg << "itsBunch_m->RefPartP_m= " << itsBunch_m->RefPartP_m << endl;
350
352
353 bool const psDump0 = 0;
354 bool const statDump0 = 0;
355
356 writePhaseSpace(0, psDump0, statDump0);
357 msg << level2 << "Dump initial phase space" << endl;
358
359
360 OrbitThreader oth(
362 itsBunch_m->getT(), (back_track ? -minTimeStep : minTimeStep), stepSizes_m,
364
365 oth.execute();
367
368 BoundingBox globalBoundingBox = oth.getBoundingBox();
369
371
372 setTime();
373
374 double time = itsBunch_m->getT() - globalTimeShift;
375 itsBunch_m->setT(time);
376
377 unsigned long long step = itsBunch_m->getGlobalTrackStep();
378 OPALTimer::Timer myt1;
379 *gmsg << "* Track start at: " << myt1.time() << ", t= " << Util::getTimeString(time) << "; "
380 << "zstart at: " << Util::getLengthString(pathLength_m) << endl;
381 *gmsg << "* Executing ParallelTracker\n"
382 << "* Initial dt = " << Util::getTimeString(itsBunch_m->getdT()) << "\n"
383 << "* Max integration steps = " << stepSizes_m.getMaxSteps() << ", next step = " << step
384 << endl
385 << endl;
386
388
389 globalEOL_m = false;
390 wakeStatus_m = false;
391 deletedParticles_m = false;
393
395
398
399 while (!stepSizes_m.reachedEnd()) {
400
401 unsigned long long trackSteps = stepSizes_m.getNumSteps() + step;
404
405 for (; step < trackSteps; ++step) {
406 Vector_t<double, 3> rmin(0.0), rmax(0.0);
407 if (itsBunch_m->getTotalNum() > 0) {
408 itsBunch_m->get_bounds(rmin, rmax);
409 {
410 *gmsg << "* Bunch bounds at step " << step + 1 << ": "
411 << " x[" << Util::getLengthString(rmin(0)) << ", "
412 << Util::getLengthString(rmax(0)) << "] "
413 << " y[" << Util::getLengthString(rmin(1)) << ", "
414 << Util::getLengthString(rmax(1)) << "] "
415 << " z[" << Util::getLengthString(rmin(2)) << ", "
416 << Util::getLengthString(rmax(2)) << "] " << endl;
417 }
418 {
419 // Print out mass and charge per particle
420 *gmsg << "* Particle mass at step " << step + 1 << ": " << itsBunch_m->getMassPerParticle() << endl;
421 *gmsg << "* Charge mass at step " << step + 1 << ": " << itsBunch_m->getChargePerParticle() << endl;
422 }
423 }
424 // ADA
425 timeIntegration1(pusher);
426
428
429 // \todo for a drift we can neglect that
430 // computeExternalFields(oth);
431
432 timeIntegration2(pusher);
433
434 /*
435 The following lines contain debugging output that was needed to fix units.
436 \todo I (Github @aliemen) will remove it once everything is correct.
437 At the moment, some units are better, but the self fields still seem off.
438 */
439 /*
440 {
441 // Calculate bunch size: rms in position
442 Vector_t<double, 3> meanR = itsBunch_m->getParticleContainer()->getMeanR();
443 Vector_t<double, 3> sumR2(0.0);
444 unsigned long long numParticles = itsBunch_m->getTotalNum();
445 auto Rview = itsBunch_m->getParticleContainer()->R.getView();
446 Kokkos::parallel_reduce(
447 "bunchSizeCalc", numParticles,
448 KOKKOS_LAMBDA(const size_t i, Vector_t<double,3>& localSumR2) {
449 Vector_t<double, 3> diff = Rview(i) - meanR;
450 localSumR2 += diff * diff;
451 }, sumR2);
452 ippl::Comm->allreduce(sumR2, 1, std::plus<Vector_t<double, 3>>());
453 Vector_t<double, 3> rmsSize = sqrt( sumR2 / static_cast<double>(numParticles) );
454 *gmsg << "* After step " << step + 1 << ": Bunch RMS size (x,y,z) = ("
455 << Util::getLengthString(rmsSize(0)) << ", "
456 << Util::getLengthString(rmsSize(1)) << ", "
457 << Util::getLengthString(rmsSize(2)) << ") " << endl;
458 }
459
460 {
461 // Output sum of particle charges
462 double localChargeSum = 0.0;
463 unsigned long long numParticles = itsBunch_m->getTotalNum();
464 auto chargeView = itsBunch_m->getParticleContainer()->Q.getView();
465 Kokkos::parallel_reduce(
466 "chargeSumCalc", numParticles,
467 KOKKOS_LAMBDA(const size_t i, double& localSum) {
468 localSum += chargeView(i);
469 }, localChargeSum);
470 double globalChargeSum = 0.0;
471 ippl::Comm->allreduce(&localChargeSum, &globalChargeSum, 1, std::plus<double>());
472 *gmsg << "* After step " << step + 1 << ": Total bunch charge = "
473 << globalChargeSum << endl; // Util::getChargeString(
474 }
475 {
476 // Extract charge of first particle for debugging
477 double firstParticleCharge = 0.0;
478 auto chargeView = itsBunch_m->getParticleContainer()->Q.getView();
479 Kokkos::parallel_reduce(
480 "firstParticleCharge", 1,
481 KOKKOS_LAMBDA(const size_t i, double& localCharge) {
482 localCharge = chargeView(0);
483 }, firstParticleCharge);
484 *gmsg << "* After step " << step + 1 << ": First particle charge = "
485 << firstParticleCharge << endl;
486 }
487 {
488 // Calculate mean electric field in x/y/z for debugging
489 Vector_t<double, 3> localSumE(0.0);
490 unsigned long long numParticles = itsBunch_m->getTotalNum();
491 auto Eview = itsBunch_m->getParticleContainer()->E.getView();
492 Kokkos::parallel_reduce(
493 "meanECalc", numParticles,
494 KOKKOS_LAMBDA(const size_t i, Vector_t<double,3>& localSum) {
495 localSum += Eview(i);
496 }, localSumE);
497 ippl::Comm->allreduce(localSumE, 1, std::plus<Vector_t<double, 3>>());
498 Vector_t<double, 3> meanE = localSumE / static_cast<double>(numParticles);
499 *gmsg << "* After step " << step + 1 << ": Mean electric field (x,y,z) = ("
500 << meanE(0) << ", "
501 << meanE(1) << ", "
502 << meanE(2) << ")" << endl;
503 }
504 {
505 // Compute mean z position for debugging
506 double localSumZ = 0.0;
507 unsigned long long numParticles = itsBunch_m->getTotalNum();
508 auto Rview = itsBunch_m->getParticleContainer()->R.getView();
509 Kokkos::parallel_reduce(
510 "meanZCalc", numParticles,
511 KOKKOS_LAMBDA(const size_t i, double& localSum) {
512 localSum += Rview(i)[2];
513 }, localSumZ);
514 double globalSumZ = 0.0;
515 ippl::Comm->allreduce(&localSumZ, &globalSumZ, 1, std::plus<double>());
516 double meanZ = globalSumZ / static_cast<double>(numParticles);
517 *gmsg << "* Mean z position at step " << step << ": " << Util::getLengthString(meanZ) << endl;
518 }
519 {
520 // Look at the first mass and the first charge in the bunch (deep copy)
521 double firstMass = 0.0;
522 double firstCharge = 0.0;
523 auto massView = itsBunch_m->getParticleContainer()->M.getView();
524 auto chargeView = itsBunch_m->getParticleContainer()->Q.getView();
525 Kokkos::parallel_reduce(
526 "firstMassCharge", 1,
527 KOKKOS_LAMBDA(const size_t i, double& localMass, double& localCharge) {
528 localMass = massView(0);
529 localCharge = chargeView(0);
530 }, firstMass, firstCharge);
531 *gmsg << "* After step " << step + 1 << ": First particle mass = "
532 << firstMass << ", charge = " << firstCharge << endl;
533 }
534 */
535
537 // \todo emitParticles(step);
538 //selectDT(back_track);
539
541
542 if (itsBunch_m->getT() > 0.0 || itsBunch_m->getdT() < 0.0) {
543 updateReference(pusher);
544 }
545
546 if (deletedParticles_m) {
547 // \todo doDelete
548 deletedParticles_m = false;
549 }
550
552
553 // if (hasEndOfLineReached(globalBoundingBox)) break;
554
555 bool const psDump =
558 bool const statDump =
561 dumpStats(step, psDump, statDump);
562
564
566 double beta = euclidean_norm(pdivg);
567 double driftPerTimeStep = std::abs(itsBunch_m->getdT()) * Physics::c * beta;
568
569 // Output driftPerTimeStep for debugging
570 *gmsg << "* Drift per time step: " << Util::getLengthString(driftPerTimeStep) << endl;
571
572 if (std::abs(stepSizes_m.getZStop() - pathLength_m) < 0.5 * driftPerTimeStep) {
573 break;
574 }
575 }
576
577 if (globalEOL_m)
578 break;
579 ++stepSizes_m;
580 }
582
584
585 bool const psDump =
587 bool const statDump =
590
591 writePhaseSpace((step + 1), psDump, statDump);
592
593 *gmsg << "* Dump phase space of last step" << endl;
594
596
597 OPALTimer::Timer myt3;
598 *gmsg << endl << "* Done executing ParallelTracker at " << myt3.time() << endl << endl;
599}
600 /*
601 OpalData::getInstance()->setPriorTrack();
602 */
603
605 itsBeamline_m.accept(*this);
610}
611
614 pushParticles(pusher);
616}
617
619 /*
620 transport and emit particles
621 that passed the cathode in the first
622 half-step or that would pass it in the
623 second half-step.
624
625 to make IPPL and the field solver happy
626 make sure that at least 10 particles are emitted
627
628 also remember that node 0 has
629 all the particles to be emitted
630
631 this has to be done *after* the calculation of the
632 space charges! thereby we neglect space charge effects
633 in the very first step of a new-born particle.
634
635 */
636
638 kickParticles(pusher);
639 // switchElements();
640 pushParticles(pusher);
641
642
643 auto dtview = itsBunch_m->getParticleContainer()->dt.getView();
644 double newdT = itsBunch_m->getdT();
645
647 "changeDT", ippl::getRangePolicy(dtview),
648 KOKKOS_LAMBDA(const size_t i) {
649 dtview(i) = newdT;
650 });
651
653}
654
655void ParallelTracker::selectDT(bool backTrack) {
656 double dt = dtCurrentTrack_m;
657 itsBunch_m->setdT(dt);
658
659 if (backTrack) {
660 itsBunch_m->setdT(-std::abs(itsBunch_m->getdT()));
661 }
662}
663
664void ParallelTracker::changeDT(bool backTrack) {
665
666 selectDT(backTrack);
667
668 auto dtview = itsBunch_m->getParticleContainer()->dt.getView();
669 double newdT = itsBunch_m->getdT();
670
672 "changeDT", ippl::getRangePolicy(dtview),
673 KOKKOS_LAMBDA(const size_t i) {
674 dtview(i) = newdT;
675 });
676
677}
678
679void ParallelTracker::emitParticles(long long step) {
680
681}
682
683void ParallelTracker::computeSpaceChargeFields(unsigned long long step) {
684
685 if (!itsBunch_m->hasFieldSolver()) {
686 *gmsg << "no solver avaidable " << endl;
687 return;
688 }
689
691
693
694 CoordinateSystemTrafo beamToReferenceCSTrafo(
695 Vector_t<double, 3>(0, 0, pathLength_m), alignment.conjugate());
696
697 CoordinateSystemTrafo referenceToBeamCSTrafo = beamToReferenceCSTrafo.inverted();
698
708 const matrix_t rot = referenceToBeamCSTrafo.getRotationMatrix();
709 const ippl::Vector<double, 3> org = referenceToBeamCSTrafo.getOrigin();
710
711
712 typedef Kokkos::View<double**> ViewMatrixType;
713 ViewMatrixType Rot("Rot", 3, 3);
714 ViewMatrixType::HostMirror h_Rot = Kokkos::create_mirror_view(Rot);
715
716 // Initialize M matrix on host.
717 for ( int i = 0; i < 3; ++i ) {
718 for ( int j = 0; j < 3; ++j ) {
719 h_Rot( i, j ) = rot(i, j);
720 }
721 }
722
723 Kokkos::deep_copy(Rot, h_Rot);
724
725 // Reset fields to zero
728
729 auto Rview = itsBunch_m->getParticleContainer()->R.getView();
730
732 "referenceToBeamCSTrafo", ippl::getRangePolicy(Rview),
733 KOKKOS_LAMBDA(const size_t k) {
734 ippl::Vector<double, 3> x = Rview(k); // x({Rview(i)[0],Rview(i)[1],Rview(i)[2]});
735 x = x - org;
736 for (size_t i = 0; i < 3; ++i ) {
737 Rview(k)[i] = 0.0;
738 for (size_t j = 0; j < 3; ++j ) {
739 Rview(k)[i] += Rot(i, j) * x(j);
740 }
741 }
742 });
743
745
746 if (step % repartFreq_m + 1 == repartFreq_m) {
748 }
749
751
753
767 auto Eview = itsBunch_m->getParticleContainer()->E.getView();
768 auto Bview = itsBunch_m->getParticleContainer()->B.getView();
770 "beamToReferenceCSTrafo", ippl::getRangePolicy(Rview),
771 KOKKOS_LAMBDA(const size_t k) {
772 ippl::Vector<double, 3> x = Rview(k); // ({Rview(k)[0],Rview(k)[1],Rview(k)[2]});
773 ippl::Vector<double, 3> e = Eview(k); // ({Eview(k)[0],Eview(k)[1],Eview(k)[2]});
774 ippl::Vector<double, 3> b = Bview(k); // ({Bview(k)[0],Bview(k)[1],Bview(k)[2]});
775
776 // beamToReferenceCSTrafo.rotateTo
777 for (size_t i = 0; i < 3; ++i ) {
778 Eview(k)[i] = 0.0;
779 Bview(k)[i] = 0.0;
780 for (size_t j = 0; j < 3; ++j ) {
781 Eview(k)[i] += Rot(j,i) * e(j);
782 Bview(k)[i] += Rot(j,i) * b(j);
783 }
784 }
785 // beamToReferenceCSTrafo.transformTo
786 x = x + org;
787 for (size_t i = 0; i < 3; ++i ) {
788 Rview(k)[i] = 0.0;
789 for (size_t j = 0; j < 3; ++j ) {
790 Rview(k)[i] += Rot(j,i) * x(j);
791 }
792 }
793 });
794
795}
796
799 Inform msg("ParallelTracker ", *gmsg);
800 const unsigned int localNum = itsBunch_m->getLocalNum();
801 bool locPartOutOfBounds = false, globPartOutOfBounds = false;
802 Vector_t<double, 3> rmin(0.0), rmax(0.0);
803 if (itsBunch_m->getTotalNum() > 0)
804 itsBunch_m->get_bounds(rmin, rmax);
806
807 try {
808 elements = oth.query(pathLength_m + 0.5 * (rmax(2) + rmin(2)), rmax(2) - rmin(2));
809 } catch (IndexMap::OutOfBounds& e) {
810 globalEOL_m = true;
811 return;
812 }
813
814 IndexMap::value_t::const_iterator it = elements.begin();
815 const IndexMap::value_t::const_iterator end = elements.end();
816
817 for (; it != end; ++it) {
818
819 CoordinateSystemTrafo refToLocalCSTrafo = (itsOpalBeamline_m.getMisalignment((*it))
822
823 CoordinateSystemTrafo localToRefCSTrafo = refToLocalCSTrafo.inverted();
824
825 (*it)->setCurrentSCoordinate(pathLength_m + rmin(2));
826
827 for (unsigned int i = 0; i < localNum; ++i) {
828 // \todo if (itsBunch_m->Bin[i] < 0)
829 // continue;
830
831 // \todo itsBunch_m->R[i] = refToLocalCSTrafo.transformTo(itsBunch_m->R[i]);
832 // \todo itsBunch_m->P[i] = refToLocalCSTrafo.rotateTo(itsBunch_m->P[i]);
833 double dt = 1.0; // \todo itsBunch_m->dt[i];
834
835 Vector_t<double, 3> localE(0.0), localB(0.0);
836
837 if ((*it)->apply(i, itsBunch_m->getT() + 0.5 * dt, localE, localB)) {
838 // itsBunch_m->R[i] = localToRefCSTrafo.transformTo(itsBunch_m->R[i]);
839 // itsBunch_m->P[i] = localToRefCSTrafo.rotateTo(itsBunch_m->P[i]);
840 // itsBunch_m->Bin[i] = -1;
841 locPartOutOfBounds = true;
842
843 continue;
844 }
845
846 // itsBunch_m->R[i] = localToRefCSTrafo.transformTo(itsBunch_m->R[i]);
847 // itsBunch_m->P[i] = localToRefCSTrafo.rotateTo(itsBunch_m->P[i]);
848 // itsBunch_m->Ef[i] += localToRefCSTrafo.rotateTo(localE);
849 // itsBunch_m->Bf[i] += localToRefCSTrafo.rotateTo(localB);
850 }
851 }
852
854
855 // \todo reduce(locPartOutOfBounds, globPartOutOfBounds, OpOrAssign());
856 ippl::Comm->reduce(locPartOutOfBounds, globPartOutOfBounds, 1, std::logical_or<bool>());
857
858 size_t ne = 0;
859 if (globPartOutOfBounds) {
860 if (itsBunch_m->hasFieldSolver()) {
862 }
863
864 // \todo
865 // else {
866 // ne = itsBunch_m->destroyT();
867 // }
869 deletedParticles_m = true;
870 }
871
872 size_t totalNum = itsBunch_m->getTotalNum();
874
875 if (ne > 0) {
876 msg << level1 << "* Deleted " << ne << " particles, "
877 << "remaining " << numParticlesInSimulation_m << " particles" << endl;
878 }
879}
880
882 if (itsBunch_m->hasFieldSolver()) {
883 *ippl::Info << "*****************************************************************" << endl;
884 *ippl::Info << "do repartition because of repartFreq_m" << endl;
885 *ippl::Info << "*****************************************************************" << endl;
887 ippl::Comm->barrier();
889 *ippl::Info << "*****************************************************************" << endl;
890 *ippl::Info << "do repartition done" << endl;
891 *ippl::Info << "*****************************************************************" << endl;
892 }
893}
894
895void ParallelTracker::dumpStats(long long step, bool psDump, bool statDump) {
896 OPALTimer::Timer myt2;
897
898 /*
899 if (itsBunch_m->getGlobalTrackStep() % 1000 + 1 == 1000) {
900 *gmsg << level1;
901 } else if (itsBunch_m->getGlobalTrackStep() % 100 + 1 == 100) {
902 *gmsg << level2;
903 } else {
904 *gmsg << level3;
905 }
906 */
907
909 *gmsg << "* " << myt2.time() << " "
910 << "Step " << std::setw(6) << itsBunch_m->getGlobalTrackStep() << "; "
911 << " -- no emission yet -- "
912 << "t= " << Util::getTimeString(itsBunch_m->getT()) << endl;
913 return;
914 }
915
916 // \todo itsBunch_m->calcEMean();
917 // size_t totalParticles_f = numParticlesInSimulation_m;
919 throw OpalException(
920 "ParallelTracker::dumpStats()",
921 "there seems to be something wrong with the position of the bunch!");
922 } else {
923 *gmsg << "* " << myt2.time() << " "
924 << "Step " << std::setw(6) << itsBunch_m->getGlobalTrackStep() << " "
925 << "at " << Util::getLengthString(pathLength_m) << ", "
926 << "t= " << Util::getTimeString(itsBunch_m->getT()) << ", "
928
929 writePhaseSpace(step, psDump, statDump);
930 }
931}
932
934
935 /*
936 minStepforReBin_m = Options::minStepForRebin;
937 RealVariable* br =
938 dynamic_cast<RealVariable*>(OpalData::getInstance()->find("MINSTEPFORREBIN"));
939 if (br)
940 minStepforReBin_m = static_cast<int>(br->getReal());
941 msg << level2 << "MINSTEPFORREBIN " << minStepforReBin_m << endl;
942 */
943
944 // there is no point to do repartitioning with one node
945 if (ippl::Comm->size() == 1) {
947 } else {
949 RealVariable* rep =
950 dynamic_cast<RealVariable*>(OpalData::getInstance()->find("REPARTFREQ"));
951 if (rep)
952 repartFreq_m = static_cast<int>(rep->getReal());
953 *gmsg << "* REPARTFREQ " << repartFreq_m << endl;
954 }
955}
956
957bool ParallelTracker::hasEndOfLineReached(const BoundingBox& globalBoundingBox) {
958 // \todo check in IPPL 1.0 it was OpBitwiseAndAssign()
959 ippl::Comm->reduce(globalEOL_m, globalEOL_m, 1, std::logical_and<bool>());
960 globalEOL_m = globalEOL_m || globalBoundingBox.isOutside(itsBunch_m->RefPartR_m);
961 return globalEOL_m;
962}
963
965
966 auto dtview = itsBunch_m->getParticleContainer()->dt.getView();
967 double newdT = itsBunch_m->getdT();
968
970 "changeDT", ippl::getRangePolicy(dtview),
971 KOKKOS_LAMBDA(const int i) {
972 dtview(i) = newdT;
973 });
974}
975
976void ParallelTracker::writePhaseSpace(const long long /*step*/, bool psDump, bool statDump) {
977 Vector_t<double, 3> externalE, externalB;
978 Vector_t<double, 3> FDext[2]; // FDext = {BHead, EHead, BRef, ERef, BTail, ETail}.
979
980 // Sample fields at (xmin, ymin, zmin), (xmax, ymax, zmax) and the centroid location. We
981 // are sampling the electric and magnetic fields at the back, front and
982 // center of the beam.
983 Vector_t<double, 3> rmin, rmax;
984 itsBunch_m->get_bounds(rmin, rmax);
985
986 if (psDump || statDump) {
987 externalB = Vector_t<double, 3>(0.0);
988 externalE = Vector_t<double, 3>(0.0);
991 itsBunch_m->getT() - 0.5 * itsBunch_m->getdT(), externalE, externalB);
992 FDext[0] = externalB; // \todo itsBunch_m->toLabTrafo_m.rotateFrom(externalB);
993 FDext[1] =
994 (externalE * Units::Vpm2MVpm); // \todo itsBunch_m->toLabTrafo_m.rotateFrom(externalE *
995 // Units::Vpm2MVpm);
996 }
997
998 if (statDump) {
999 itsDataSink_m->dumpSDDS(itsBunch_m, FDext, -1.0);
1000 *gmsg << "* Wrote beam statistics." << endl;
1001 }
1002
1003 if (psDump && (itsBunch_m->getTotalNum() > 0)) {
1004 // Write bunch to .h5 file.
1006
1007 /*
1008 // Write fields to .h5 file.
1009 const size_t localNum = itsBunch_m->getLocalNum();
1010 double distToLastStop = stepSizes_m.getFinalZStop() - pathLength_m;
1011 Vector_t<double, 3> beta = itsBunch_m->RefPartP_m / Util::getGamma(itsBunch_m->RefPartP_m);
1012 Vector_t<double, 3> driftPerTimeStep =
1013 itsBunch_m->getdT()
1014 * Physics::c; // \todo * itsBunch_m->toLabTrafo_m.rotateFrom(beta);
1015 bool driftToCorrectPosition =
1016 std::abs(distToLastStop) < 0.5 * euclidean_norm(driftPerTimeStep);
1017 // \todo Ppos_t stashedR;
1018 Vector_t<double, 3> stashedR;
1019 Vector_t<double, 3> stashedRefPartR;
1020
1021 if (driftToCorrectPosition) {
1022 const double tau =
1023 distToLastStop / euclidean_norm(driftPerTimeStep) * itsBunch_m->getdT();
1024
1025 if (localNum > 0) {
1026
1027 stashedR.create(localNum);
1028 stashedR = itsBunch_m->R;
1029 stashedRefPartR = itsBunch_m->RefPartR_m;
1030
1031 for (size_t i = 0; i < localNum; ++i) {
1032 itsBunch_m->R[i] +=
1033 tau
1034 * (Physics::c * itsBunch_m->P[i] / Util::getGamma(itsBunch_m->P[i])
1035 - driftPerTimeStep / itsBunch_m->getdT());
1036 }
1037
1038 }
1039
1040 driftPerTimeStep = itsBunch_m->toLabTrafo_m.rotateTo(driftPerTimeStep);
1041 itsBunch_m->RefPartR_m =
1042 itsBunch_m->RefPartR_m + tau * driftPerTimeStep / itsBunch_m->getdT();
1043 CoordinateSystemTrafo update(
1044 tau * driftPerTimeStep / itsBunch_m->getdT(), Quaternion(1.0, 0.0, 0.0, 0.0));
1045 itsBunch_m->toLabTrafo_m = itsBunch_m->toLabTrafo_m * update.inverted();
1046
1047 itsBunch_m->set_sPos(stepSizes_m.getFinalZStop());
1048
1049 itsBunch_m->calcBeamParameters();
1050 }
1051 if (!statDump && !driftToCorrectPosition)
1052 itsBunch_m->calcBeamParameters();
1053
1054 msg << *itsBunch_m << endl;
1055
1056 itsDataSink_m->dumpH5(itsBunch_m, FDext);
1057 */
1058
1059 /*
1060
1061 if (driftToCorrectPosition) {
1062 if (localNum > 0) {
1063 itsBunch_m->R = stashedR;
1064 stashedR.destroy(localNum, 0);
1065 }
1066
1067 itsBunch_m->RefPartR_m = stashedRefPartR;
1068 itsBunch_m->set_sPos(pathLength_m);
1069
1070 itsBunch_m->calcBeamParameters();
1071 }
1072
1073 */
1074
1075 *gmsg << level2 << "* Wrote beam phase space." << endl;
1076 }
1077}
1078
1082}
1083
1085 const double direction = back_track ? -1 : 1;
1086 const double dt = direction * std::min(itsBunch_m->getT(), direction * itsBunch_m->getdT());
1087 const double scaleFactor = Physics::c * dt;
1088 Vector_t<double, 3> Ef(0.0), Bf(0.0);
1089
1090 itsBunch_m->RefPartR_m /= scaleFactor;
1092 itsBunch_m->RefPartR_m *= scaleFactor;
1093
1095 IndexMap::value_t::const_iterator it = elements.begin();
1096 const IndexMap::value_t::const_iterator end = elements.end();
1097
1098 for (; it != end; ++it) {
1099 const CoordinateSystemTrafo& refToLocalCSTrafo =
1101
1102 Vector_t<double, 3> localR = refToLocalCSTrafo.transformTo(itsBunch_m->RefPartR_m);
1103 Vector_t<double, 3> localP = refToLocalCSTrafo.rotateTo(itsBunch_m->RefPartP_m);
1104 Vector_t<double, 3> localE(0.0), localB(0.0);
1105
1106 if ((*it)->applyToReferenceParticle(
1107 localR, localP, itsBunch_m->getT() - 0.5 * dt, localE, localB)) {
1108 *gmsg << level1 << "The reference particle hit an element" << endl;
1109 globalEOL_m = true;
1110 }
1111
1112 Ef += refToLocalCSTrafo.rotateFrom(localE);
1113 Bf += refToLocalCSTrafo.rotateFrom(localB);
1114 }
1115
1116 pusher.kick(itsBunch_m->RefPartR_m, itsBunch_m->RefPartP_m, Ef, Bf, dt);
1117
1118 itsBunch_m->RefPartR_m /= scaleFactor;
1120 itsBunch_m->RefPartR_m *= scaleFactor;
1121}
1122
1124 //const unsigned int localNum = itsBunch_m->getLocalNum();
1125 //for (unsigned int i = 0; i < localNum; ++i) {
1126 /* \todo host device ....
1127 itsBunch_m->R[i] = trafo.transformTo(itsBunch_m->R[i]);
1128 itsBunch_m->P[i] = trafo.rotateTo(itsBunch_m->P[i]);
1129 itsBunch_m->Ef[i] = trafo.rotateTo(itsBunch_m->Ef[i]);
1130 itsBunch_m->Bf[i] = trafo.rotateTo(itsBunch_m->Bf[i]);
1131 */
1132 //}
1133 *gmsg << "* ParallelTracker::transformBunch() is not implemented." << endl;
1134}
1135
1137 /*
1138 Inform m("updateRefToLabCSTrafo ");
1139 m << " initial RefPartR: " << itsBunch_m->RefPartR_m << endl;
1140 m << " RefPartR after transformFrom( " << itsBunch_m->RefPartR_m << endl;
1141 itsBunch_m->toLabTrafo_m.rotateFrom(itsBunch_m->RefPartP_m);
1142 Vector_t<double, 3> P = itsBunch_m->RefPartP_m;
1143 m << " CS " << itsBunch_m->toLabTrafo_m << endl;
1144 m << " rotMat= " << itsBunch_m->toLabTrafo_m.getRotationMatrix() << endl;
1145 m << " initial: path Lenght= " << pathLength_m << endl;
1146 m << " dt= " << itsBunch_m->getdT() << " R=" << R << endl;
1147 */
1148
1149 // std::cout << "Before updateRefToLabCSTrafo(): RefPartP_m: " << itsBunch_m->RefPartP_m << std::endl;
1150 //std::cout << "Rotation matrix = " << itsBunch_m->toLabTrafo_m.getRotationMatrix() << std::endl;
1151
1153
1155 //*gmsg << "* updateRefToLabCSTrafo(): P before transformFrom: " << itsBunch_m->RefPartP_m << endl;
1157 //*gmsg << "* updateRefToLabCSTrafo(): P after transformFrom: " << P << endl;
1158
1159 pathLength_m += std::copysign(1, itsBunch_m->getdT()) * euclidean_norm(R);
1160
1161 /*
1162 This will produce a NaN Quaternion if P is parallel to the z axis because the cross product
1163 of two parallel vectors is zero, leading to a zero norm when normalizing the axis of rotation.
1164 The fix is to check if P is parallel to the z axis and not rotate in that case.
1165 */
1166
1167 // First normalize P
1168 double normP = euclidean_norm(P);
1169 if (normP < 1e-12) {
1170 pathLength_m += 0.0;
1171 return;
1172 }
1173 P /= normP;
1174
1175 // After normalizing, we can check alignment with z-axis
1176 Vector_t<double, 3> target(0, 0, 1);
1177 double dot = P.dot(target);
1178 bool aligned = std::abs(std::abs(dot) - 1.0) < 1e-6;
1179 Quaternion Q = aligned ? Quaternion(0, 0, 0, 1.0) : getQuaternion(P, target);
1180
1181 if (aligned) {
1182 *gmsg << "* Warning: Reference particle momentum is aligned with z-axis; no quaternion rotation applied." << endl;
1183 }
1184
1185 CoordinateSystemTrafo update(R, Q);
1186
1188 transformBunch(update);
1189
1192}
1193
1194void ParallelTracker::applyFractionalStep(const BorisPusher& pusher, double tau) {
1195 double t = itsBunch_m->getT();
1196 t += tau;
1197 itsBunch_m->setT(t);
1198
1199 // the push method below pushes for half a time step. Hence the ref particle
1200 // should be pushed for 2 * tau
1201 itsBunch_m->RefPartR_m /= (Physics::c * 2 * tau);
1203 itsBunch_m->RefPartR_m *= (Physics::c * 2 * tau);
1204
1212}
1213
1215 StepSizeConfig stepSizesCopy(stepSizes_m);
1216 if (back_track) {
1217 stepSizesCopy.shiftZStopLeft(zstart_m);
1218 }
1219
1220 double t = 0.0;
1221 itsBunch_m->setT(t);
1222
1223 dtCurrentTrack_m = stepSizesCopy.getdT();
1224 selectDT();
1225
1226 while (true) {
1227 autophaseCavities(pusher);
1228
1229 t += itsBunch_m->getdT();
1230 itsBunch_m->setT(t);
1231
1234
1235 // \todo should not need the tmp
1238
1240 double speed = euclidean_norm(tmp);
1241
1242
1243 if (pathLength_m > stepSizesCopy.getZStop()) {
1244 ++stepSizesCopy;
1245
1246 if (stepSizesCopy.reachedEnd()) {
1247 --stepSizesCopy;
1248 double tau = (stepSizesCopy.getZStop() - pathLength_m) / speed;
1249 applyFractionalStep(pusher, tau);
1250
1251 break;
1252 }
1253
1254 dtCurrentTrack_m = stepSizesCopy.getdT();
1255 selectDT();
1256 }
1257
1258 if (std::abs(pathLength_m - zstart_m) <= 0.5 * itsBunch_m->getdT() * speed) {
1259 double tau = (zstart_m - pathLength_m) / speed;
1260 applyFractionalStep(pusher, tau);
1261
1262 break;
1263 }
1264 }
1265
1266 changeDT();
1267}
1268
1270 double t = itsBunch_m->getT();
1272 pusher.push(nextR, itsBunch_m->RefPartP_m, itsBunch_m->getdT());
1273 nextR *= Physics::c * itsBunch_m->getdT();
1274
1275 auto elementSet = itsOpalBeamline_m.getElements(nextR);
1276 for (auto element : elementSet) {
1277 if (element->getType() == ElementType::TRAVELINGWAVE) {
1278 const TravelingWave* TWelement = static_cast<const TravelingWave*>(element.get());
1279 if (!TWelement->getAutophaseVeto()) {
1280 CavityAutophaser ap(itsReference, element);
1284 itsBunch_m->getdT());
1285 }
1286
1287 } else if (element->getType() == ElementType::RFCAVITY) {
1288 const RFCavity* RFelement = static_cast<const RFCavity*>(element.get());
1289 if (!RFelement->getAutophaseVeto()) {
1290 CavityAutophaser ap(itsReference, element);
1294 itsBunch_m->getdT());
1295 }
1296 }
1297 }
1298}
1299
1301 unsigned int who;
1302 unsigned int whom;
1303 unsigned int howMany;
1304};
double dot(const Vector3D &lhs, const Vector3D &rhs)
Vector dot product.
Definition: Vector3D.cpp:118
ElementType
Definition: ElementBase.h:88
elements
Definition: IndexMap.cpp:162
Quaternion getQuaternion(ippl::Vector< double, 3 > u, ippl::Vector< double, 3 > ref)
Definition: Quaternion.cpp:35
boost::numeric::ublas::matrix< double > matrix_t
Definition: BoostMatrix.h:23
PartBunch< T, Dim >::ConstIterator end(PartBunch< T, Dim > const &bunch)
Inform * gmsg
Definition: changes.cpp:7
std::list< ClassicField > FieldList
Definition: ClassicField.h:43
std::string::const_iterator iterator_t
Definition: array.cpp:19
KOKKOS_INLINE_FUNCTION double euclidean_norm(const Vector_t< T, D > &v)
Definition: OPALTypes.h:35
constexpr KOKKOS_INLINE_FUNCTION auto second()
Definition: AbsorbingBC.h:14
Inform & level2(Inform &inf)
Definition: Inform.cpp:51
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
Inform & level1(Inform &inf)
Definition: Inform.cpp:48
STL namespace.
constexpr double e
The value of.
Definition: Physics.h:39
constexpr double c
The velocity of light in m/s.
Definition: Physics.h:45
constexpr double Vpm2MVpm
Definition: Units.h:125
constexpr double deg2rad
Definition: Units.h:143
std::string::iterator iterator
Definition: MSLang.h:15
T isinf(T x)
isinf function with adjusted return type
Definition: matheval.hpp:70
T isnan(T x)
isnan function with adjusted return type
Definition: matheval.hpp:66
int psDumpFreq
The frequency to dump the phase space, i.e.dump data when steppsDumpFreq==0.
Definition: Options.cpp:39
int repartFreq
The frequency to do particles repartition for better load balance between nodes.
Definition: Options.cpp:49
int statDumpFreq
The frequency to dump statistical values, e.e. dump data when stepstatDumpFreq==0.
Definition: Options.cpp:41
double getGamma(ippl::Vector< double, 3 > p)
Definition: Util.h:44
std::string getEnergyString(double energyInMeV, unsigned int precision=3)
Definition: Util.h:140
std::string getTimeString(double time, unsigned int precision=3)
Definition: Util.h:70
std::string getLengthString(double spos, unsigned int precision=3)
Definition: Util.h:93
std::unique_ptr< Inform > Info
Definition: Ippl.h:29
KOKKOS_INLINE_FUNCTION Vector< T, Dim > max(const Vector< T, Dim > &a, const Vector< T, Dim > &b)
Definition: Vector.hpp:226
RangePolicy< View::rank, typenameView::execution_space, PolicyArgs... >::policy_type getRangePolicy(const View &view, int shift=0)
KOKKOS_INLINE_FUNCTION Vector< T, Dim > min(const Vector< T, Dim > &a, const Vector< T, Dim > &b)
Definition: Vector.hpp:217
void parallel_for(const std::string &name, const ExecPolicy &policy, const FunctorType &functor)
std::unique_ptr< mpi::Communicator > Comm
Definition: Ippl.h:22
Interface for a single beam element.
Definition: Component.h:50
virtual const std::string & getName() const
Get element name.
virtual void accept(BeamlineVisitor &visitor) const =0
Apply visitor.
@Class Offset
Definition: Offset.h:66
void updateGeometry(Vector_t< double, 3 > startPosition, Vector_t< double, 3 > startDirection)
Definition: Offset.cpp:177
virtual bool getAutophaseVeto() const
Definition: RFCavity.h:380
virtual void setPhasem(double phase)
Definition: RFCavity.h:344
virtual void setAutophaseVeto(bool veto=true)
Definition: RFCavity.h:376
Ring describes a ring type geometry for tracking.
Definition: Ring.h:62
double getBeamRInit() const
Get the initial beam radius.
Definition: Ring.h:255
virtual void initialise(PartBunch_t *bunch, double &startField, double &endField) override
Initialise the Ring.
Definition: Ring.cpp:186
double getBeamPRInit() const
Get the initial beam radial momentum.
Definition: Ring.h:275
Vector_t< double, 3 > getNextPosition() const
Get the initial element's start position in cartesian coordinates.
Definition: Ring.cpp:234
double getSymmetry() const
Get the rotational symmetry of the ring (number of cells)
Definition: Ring.h:332
Vector_t< double, 3 > getNextNormal() const
Get the initial element's start normal in cartesian coordinates.
Definition: Ring.cpp:243
virtual ElementBase * clone() const override
Inherited copy constructor.
Definition: Ring.h:155
double getBeamPhiInit() const
Get the initial beam azimuthal angle.
Definition: Ring.h:265
double getHarmonicNumber()
Get the harmonic number for RF (number of bunches in the ring)
Definition: Ring.h:234
void appendElement(const Component &element)
Add element to the ring.
Definition: Ring.cpp:252
Sector bending magnet with an FFA-style field index and spiral end shape.
Bending magnet with an exponential dependence on field in the vertical plane.
std::vector< MaxPhasesT >::iterator getLastMaxPhases()
Definition: OpalData.cpp:401
double getGlobalPhaseShift()
units: (sec)
Definition: OpalData.cpp:450
std::vector< MaxPhasesT >::iterator getFirstMaxPhases()
Definition: OpalData.cpp:397
Object * find(const std::string &name)
Find entry.
Definition: OpalData.cpp:563
void setInPrepState(bool state)
Definition: OpalData.cpp:299
void setGlobalPhaseShift(double shift)
units: (sec)
Definition: OpalData.cpp:445
static OpalData * getInstance()
Definition: OpalData.cpp:195
void setOpenMode(OpenMode openMode)
Definition: OpalData.cpp:348
void setdT(double dt)
double getMassPerParticle() const
void performBunchSanityChecks() const
Definition: PartBunch.cpp:669
bool hasFieldSolver()
void computeSelfFields()
void incTrackSteps()
double getdT() const
void setT(double t)
size_t boundp_destroyT()
This is only temporary in order to get the collimator and pepperpot working.
void calcBeamParameters()
void boundp()
Vector_t< T, Dim > get_pmean() const
void do_binaryRepart()
double getT() const
void updateNumTotal()
Vector_t< T, Dim > RefPartR_m
Reference particle structures.
Definition: #PartBunch.hpp#:85
std::shared_ptr< ParticleContainer_t > getParticleContainer()
Definition: PartBunch.h:221
void incrementT()
double get_meanKineticEnergy()
void get_bounds(Vector_t< T, Dim > &rmin, Vector_t< T, Dim > &rmax) const
void set_sPos(double s)
get the spos of the primary beam
size_t getLocalNum() const
double getChargePerParticle() const
get the macro particle charge
Vector_t< T, Dim > get_centroid() const
long long getGlobalTrackStep() const
size_t getTotalNum() const
Vector_t< T, Dim > RefPartP_m
Definition: #PartBunch.hpp#:86
void setGlobalMeanR(Vector_t< T, Dim > globalMeanR)
bool weHaveEnergyBins()
CoordinateSystemTrafo toLabTrafo_m
Definition: #PartBunch.hpp#:88
const PartData itsReference
The reference information.
double getPhaseAtMaxEnergy(const Vector_t< double, 3 > &R, const Vector_t< double, 3 > &P, double t, double dt)
ippl::Vector< double, 3 > getOrigin() const
ippl::Vector< double, 3 > transformFrom(const ippl::Vector< double, 3 > &r) const
ippl::Vector< double, 3 > rotateFrom(const ippl::Vector< double, 3 > &r) const
ippl::Vector< double, 3 > rotateTo(const ippl::Vector< double, 3 > &r) const
matrix_t getRotationMatrix() const
ippl::Vector< double, 3 > transformTo(const ippl::Vector< double, 3 > &r) const
CoordinateSystemTrafo inverted() const
std::set< std::shared_ptr< Component > > value_t
Definition: IndexMap.h:47
IndexMap::value_t query(IndexMap::key_t::first_type step, IndexMap::key_t::second_type length)
BoundingBox getBoundingBox() const
virtual void visitScalingFFAMagnet(const ScalingFFAMagnet &bend)
Apply the algorithm to a scaling FFA.
void pushParticles(const BorisPusher &pusher)
void emitParticles(long long step)
void updateRFElement(std::string elName, double maxPhi)
void changeDT(bool backTrack=false)
IpplTimings::TimerRef BinRepartTimer_m
DataSink * itsDataSink_m
virtual void visitOffset(const Offset &)
Apply the algorithm to a offset (placement).
beamline_list FieldDimensions
void autophaseCavities(const BorisPusher &pusher)
bool applyPluginElements(const double dt)
IpplTimings::TimerRef OrbThreader_m
virtual void execute()
Apply the algorithm to the top-level beamline.
std::vector< PluginElement * > pluginElements_m
void kickParticles(const BorisPusher &pusher)
unsigned long long repartFreq_m
void transformBunch(const CoordinateSystemTrafo &trafo)
virtual void visitRing(const Ring &ring)
Apply the algorithm to a ring.
void writePhaseSpace(const long long step, bool psDump, bool statDump)
void timeIntegration2(BorisPusher &pusher)
virtual void visitBeamline(const Beamline &)
Apply the algorithm to a beam line.
double zstart_m
where to start
size_t numParticlesInSimulation_m
void applyFractionalStep(const BorisPusher &pusher, double tau)
std::list< Component * > myElements
void updateReference(const BorisPusher &pusher)
void updateReferenceParticle(const BorisPusher &pusher)
void buildupFieldList(double BcParameter[], ElementType elementType, Component *elptr)
IpplTimings::TimerRef timeIntegrationTimer2_m
void computeExternalFields(OrbitThreader &oth)
void computeSpaceChargeFields(unsigned long long step)
std::pair< ElementType, element_pair > type_pair
bool hasEndOfLineReached(const BoundingBox &globalBoundingBox)
void dumpStats(long long step, bool psDump, bool statDump)
IpplTimings::TimerRef timeIntegrationTimer1_m
IpplTimings::TimerRef PluginElemTimer_m
void selectDT(bool backTrack=false)
void findStartPosition(const BorisPusher &pusher)
IpplTimings::TimerRef fieldEvaluationTimer_m
void timeIntegration1(BorisPusher &pusher)
StepSizeConfig stepSizes_m
stores informations where to change the time step and where to stop the simulation,...
virtual ~ParallelTracker()
OpalBeamline itsOpalBeamline_m
virtual void visitVerticalFFAMagnet(const VerticalFFAMagnet &bend)
Apply the algorithm to a vertical FFA magnet.
Particle reference data.
Definition: PartData.h:37
double getBeta() const
The relativistic beta per particle.
Definition: PartData.h:126
double getGamma() const
The relativistic gamma per particle.
Definition: PartData.h:131
double beta
Definition: PartData.h:98
Quaternion conjugate() const
Definition: Quaternion.hpp:88
double getMinTimeStep() const
StepSizeConfig & advanceToPos(double spos)
double getdT() const
void printDirect(Inform &out) const
void shiftZStopLeft(double back)
double getZStop() const
void push_back(double dt, double zstop, unsigned long numSteps)
bool reachedEnd() const
unsigned long getNumSteps() const
unsigned long long getMaxSteps() const
void sortAscendingZStop()
const Beamline & itsBeamline_m
Definition: Tracker.h:119
PartBunch_t * itsBunch_m
The bunch of particles to be tracked.
Definition: Tracker.h:123
An abstract sequence of beam line components.
Definition: Beamline.h:34
bool getRelativeFlag() const
Definition: TBeamline.h:423
Quaternion getInitialDirection() const
Definition: TBeamline.h:413
Vector_t< double, 3 > getOrigin3D() const
Definition: TBeamline.h:403
virtual void iterate(BeamlineVisitor &, bool r2l) const
Apply visitor to all elements of the line.
Definition: TBeamline.h:197
FieldList getElementByType(ElementType)
Vector_t< double, 3 > rotateToLocalCS(const std::shared_ptr< Component > &comp, const Vector_t< double, 3 > &r) const
Definition: OpalBeamline.h:149
CoordinateSystemTrafo getMisalignment(const std::shared_ptr< Component > &comp) const
Definition: OpalBeamline.h:168
void merge(OpalBeamline &rhs)
unsigned long getFieldAt(const unsigned int &, const Vector_t< double, 3 > &, const long &, const double &, Vector_t< double, 3 > &, Vector_t< double, 3 > &)
void prepareSections()
void compute3DLattice()
CoordinateSystemTrafo getCSTrafoLab2Local(const std::shared_ptr< Component > &comp) const
Definition: OpalBeamline.h:159
Vector_t< double, 3 > transformToLocalCS(const std::shared_ptr< Component > &comp, const Vector_t< double, 3 > &r) const
Definition: OpalBeamline.h:139
std::set< std::shared_ptr< Component > > getElements(const Vector_t< double, 3 > &x)
void switchElementsOff()
void save3DInput()
void save3DLattice()
void activateElements()
void swap(OpalBeamline &rhs)
KOKKOS_INLINE_FUNCTION void kick(const Vector_t< double, 3 > &R, Vector_t< double, 3 > &P, const Vector_t< double, 3 > &Ef, const Vector_t< double, 3 > &Bf, const double &dt) const
Definition: BorisPusher.h:67
KOKKOS_INLINE_FUNCTION void push(Vector_t< double, 3 > &R, const Vector_t< double, 3 > &P, const double &dt) const
Definition: BorisPusher.h:129
bool isOutside(const Vector_t< double, 3 > &position) const
Definition: BoundingBox.h:60
void dumpH5(PartBunch_t *beam, Vector_t< double, 3 > FDext[]) const
Definition: DataSink.cpp:62
void storeCavityInformation()
Write cavity information from H5 file.
Definition: DataSink.cpp:107
void dumpSDDS(PartBunch_t *beam, Vector_t< double, 3 > FDext[], const double &azimuth=-1) const
Definition: DataSink.cpp:81
The base class for all OPAL exceptions.
Definition: OpalException.h:28
std::string time() const
Return time.
virtual double getReal() const
Return value.
KOKKOS_INLINE_FUNCTION T dot(const Vector< T, Dim > &rhs) const
Definition: Vector.hpp:182
Definition: Inform.h:40
static void stopTimer(TimerRef t)
Definition: IpplTimings.h:156
static void startTimer(TimerRef t)
Definition: IpplTimings.h:153