OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
ParallelCyclotronTracker.cpp
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//
29
34#include "AbsBeamline/Drift.h"
35#include "AbsBeamline/Marker.h"
36#include "AbsBeamline/Monitor.h"
43#include "AbsBeamline/Offset.h"
45#include "AbsBeamline/Probe.h"
46#include "AbsBeamline/RBend.h"
48#include "AbsBeamline/Ring.h"
49#include "AbsBeamline/SBend.h"
50#include "AbsBeamline/SBend3D.h"
52#include "AbsBeamline/Septum.h"
55#include "AbsBeamline/Vacuum.h"
59
62
63#include "Algorithms/Ctunes.h"
65
68
69#include "Beamlines/Beamline.h"
71
73
75
76#include "Physics/Physics.h"
77#include "Physics/Units.h"
78
80#include "Structure/DataSink.h"
82
84#include "Utilities/Options.h"
85
86#include <boost/filesystem.hpp>
87
88#include <cmath>
89#include <fstream>
90#include <iostream>
91#include <limits>
92#include <numeric>
93
95
99
100extern Inform *gmsg;
101
116 DataSink& ds,
117 const PartData& reference,
118 bool revBeam, bool revTrack,
119 int maxSTEPS,
120 Steppers::TimeIntegrator timeintegrator,
121 const int& numBunch,
122 const double& mbEta,
123 const double& mbPara,
124 const std::string& mbMode,
125 const std::string& mbBinning)
126 : Tracker(beamline, bunch, reference, revBeam, revTrack)
127 , bgf_m(nullptr)
128 , maxSteps_m(maxSTEPS)
129 , lastDumpedStep_m(0)
130 , myNode_m(Ippl::myNode())
131 , initialLocalNum_m(bunch->getLocalNum())
132 , initialTotalNum_m(bunch->getTotalNum())
133 , opalRing_m(nullptr)
134 , itsStepper_mp(nullptr)
135 , mode_m(TrackingMode::UNDEFINED)
136 , stepper_m(timeintegrator)
137{
138 itsBeamline = dynamic_cast<Beamline*>(beamline.clone());
139 itsDataSink = &ds;
140
141 if ( numBunch > 1 ) {
142 mbHandler_m = std::unique_ptr<MultiBunchHandler>(
143 new MultiBunchHandler(bunch, numBunch, mbEta,
144 mbPara, mbMode, mbBinning)
145 );
146 }
147
149 TransformTimer_m = IpplTimings::getTimer("Frametransform");
151 BinRepartTimer_m = IpplTimings::getTimer("Binaryrepart");
152 PluginElemTimer_m = IpplTimings::getTimer("PluginElements");
153 DelParticleTimer_m = IpplTimings::getTimer("DeleteParticles");
154
156}
157
163 if (bgf_m)
164 lossDs_m->save();
165 for (Component* component : myElements) {
166 delete(component);
167 }
168 for (auto fd : FieldDimensions) {
169 delete(fd);
170 }
171 delete itsBeamline;
172 // delete opalRing_m;
173}
174
175
177 if (!bgf_m) return;
178
179 Inform msg("bgf_main_collision_test ");
180
185 Vector_t intecoords = 0.0;
186
187 // This has to match the dT in the rk4 pusher
188 double dtime = itsBunch_m->getdT() * getHarmonicNumber();
189
190 int triId = 0;
191 for (size_t i = 0; i < itsBunch_m->getLocalNum(); i++) {
192 int res = bgf_m->partInside(itsBunch_m->R[i], itsBunch_m->P[i],
193 dtime, intecoords, triId);
194 if (res >= 0) {
195 lossDs_m->addParticle(OpalParticle(itsBunch_m->ID[i],
196 itsBunch_m->R[i], itsBunch_m->P[i],
197 itsBunch_m->getT(),
198 itsBunch_m->Q[i], itsBunch_m->M[i]),
199 std::make_pair(turnnumber_m, itsBunch_m->bunchNum[i]));
200 itsBunch_m->Bin[i] = -1;
201 Inform gmsgALL("OPAL", INFORM_ALL_NODES);
202 gmsgALL << level4 << "* Particle " << itsBunch_m->ID[i]
203 << " lost on boundary geometry" << endl;
204 }
205 }
206}
207
208// only used for dumping into stat file
209void ParallelCyclotronTracker::dumpAngle(const double& theta,
210 double& prevAzimuth,
211 double& azimuth,
212 const short& bunchNr) {
213 if ( prevAzimuth < 0.0 ) { // only at first occurrence
214 double plus = 0.0;
215 if ( OpalData::getInstance()->inRestartRun() ) {
216 plus = 360.0 * (turnnumber_m - bunchNr);
217 }
218 azimuth = theta + plus;
219 } else {
220 double dtheta = theta - prevAzimuth;
221 if ( dtheta < 0 ) {
222 dtheta += 360.0;
223 }
224 if ( dtheta > 180 ) { // rotating clockwise, reduce angle
225 dtheta -= 360;
226 }
227 azimuth += dtheta;
228 }
229 prevAzimuth = theta;
230}
231
232
234 return Units::m2mm * std::sqrt(meanR(0) * meanR(0) + meanR(1) * meanR(1));
235}
236
237
239 const double& dt)
240{
241 // the last element in dotP is the dot-product over all particles
242 std::vector<double> dotP(dl.size());
244
245 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
246 dotP[itsBunch_m->bunchNum[i]] += dot(itsBunch_m->P[i], itsBunch_m->P[i]);
247 }
248
249 allreduce(dotP.data(), dotP.size(), std::plus<double>());
250
251 // dot-product over all particles
252 double sum = std::accumulate(dotP.begin(), dotP.end() - 1, 0);
253 dotP.back() = sum / double(itsBunch_m->getTotalNum());
254
255 // bunch specific --> multi-bunches only
256 for (short b = 0; b < (short)dotP.size() - 1; ++b) {
257 dotP[b] /= double(itsBunch_m->getTotalNumPerBunch(b));
258 }
259
260 } else if ( itsBunch_m->getLocalNum() == 0 ) {
261 // here we are in DumpFrame::GLOBAL mode
262 dotP[0] = 0.0;
263 } else {
264 // here we are in DumpFrame::GLOBAL mode
265 dotP[0] = dot(itsBunch_m->P[0], itsBunch_m->P[0]);
266 }
267
268 for (size_t i = 0; i < dotP.size(); ++i) {
269 double const gamma = std::sqrt(1.0 + dotP[i]);
270 double const beta = std::sqrt(dotP[i]) / gamma;
271 dl[i] = c_mmtns * dt * Units::mm2m * beta;
272
273 }
274}
275
276
282void ParallelCyclotronTracker::openFiles(size_t numFiles, std::string SfileName) {
283
284 for (unsigned int i=0; i<numFiles; i++) {
285 std::string SfileName2 = SfileName;
286 if (i<=2) {
287 SfileName2 += std::string("-Angle" + std::to_string(int(i)) + ".dat");
288 }
289 else {
290 // for single Particle Mode, output after each turn, to define matched initial phase ellipse.
291 SfileName2 += std::string("-afterEachTurn.dat");
292 }
293
294 outfTheta_m.emplace_back(new std::ofstream(SfileName2.c_str()));
295 outfTheta_m.back()->precision(8);
296 outfTheta_m.back()->setf(std::ios::scientific, std::ios::floatfield);
297 *outfTheta_m.back() << "# r [mm] beta_r*gamma "
298 << "theta [deg] beta_theta*gamma "
299 << "z [mm] beta_z*gamma" << std::endl;
300 }
301}
302
309 for (auto & file : outfTheta_m) {
310 file->close();
311 }
312}
313
314
321 *gmsg << "* ----------------------------- Cyclotron -------------------------------- *" << endl;
322
323 cycl_m = dynamic_cast<Cyclotron*>(cycl.clone());
324 myElements.push_back(cycl_m);
325
326 // Is this a Spiral Inflector Simulation? If yes, we'll give the user some
327 // useful information
329 if (spiral_flag) {
330 *gmsg << endl << "* This is a Spiral Inflector Simulation! This means the following:" << endl;
331 *gmsg << "* 1.) It is up to the user to provide appropriate geometry, electric and magnetic fields!" << endl;
332 *gmsg << "* (Use BANDRF type cyclotron and use RFMAPFN to load both magnetic" << endl;
333 *gmsg << "* and electric fields, setting SUPERPOSE to an array of TRUE values.)" << endl;
334 *gmsg << "* 2.) For high currents it is strongly recommended to use the SAAMG fieldsolver," << endl;
335 *gmsg << "* FFT does not give the correct results (boundary conditions are missing)." << endl;
336 *gmsg << "* 3.) The whole geometry will be meshed and used for the fieldsolver." << endl;
337 *gmsg << "* There will be no transformations of the bunch into a local frame und consequently," << endl;
338 *gmsg << "* the problem will be treated non-relativistically!" << endl;
339 *gmsg << "* (This is not an issue for spiral inflectors as they are typically < 100 keV/amu.)" << endl;
340 *gmsg << endl << "* Note: For now, multi-bunch mode (MBM) needs to be de-activated for spiral inflector" << endl;
341 *gmsg << "* and space charge needs to be solved every time-step. numBunch_m and scSolveFreq are reset." << endl;
342 if (isMultiBunch()) {
343 mbHandler_m = nullptr;
344 }
345 }
346
347 // Fresh run (no restart):
349
350 // Get reference values from cyclotron element
351 // For now, these are still stored in mm. should be the only ones. -DW
357
358 if (referenceTheta <= -180.0 || referenceTheta > 180.0) {
359 throw OpalException("Error in ParallelCyclotronTracker::visitCyclotron",
360 "PHIINIT is out of [-180, 180)!");
361 }
362
364
365 // Calculate reference azimuthal (tangential) momentum from total-, z- and radial momentum:
366 float insqrt = referencePtot * referencePtot - \
367 referencePr * referencePr - referencePz * referencePz;
368
369 if (insqrt < 0) {
370 if (insqrt > -1.0e-10) {
371 referencePt = 0.0;
372 } else {
373 throw OpalException("Error in ParallelCyclotronTracker::visitCyclotron",
374 "Pt imaginary!");
375 }
376 } else {
377 referencePt = std::sqrt(insqrt);
378 }
379
380 if (referencePtot < 0.0) {
381 referencePt *= -1.0;
382 }
383 // End calculate referencePt
384
385 // Restart a run:
386 } else {
387
388 // If the user wants to save the restarted run in local frame,
389 // make sure the previous h5 file was local too
391 if (!previousH5Local) {
392 throw OpalException("Error in ParallelCyclotronTracker::visitCyclotron",
393 "You are trying a local restart from a global h5 file!");
394 }
395 // Else, if the user wants to save the restarted run in global frame,
396 // make sure the previous h5 file was global too
397 } else {
398 if (previousH5Local) {
399 throw OpalException("Error in ParallelCyclotronTracker::visitCyclotron",
400 "You are trying a global restart from a local h5 file!");
401 }
402 }
403
404 // Adjust some of the reference variables from the h5 file
408 if (referenceTheta <= -180.0 || referenceTheta > 180.0) {
409 throw OpalException("Error in ParallelCyclotronTracker::visitCyclotron",
410 "PHIINIT is out of [-180, 180)!");
411 }
412 }
413
416
417 *gmsg << endl;
418 *gmsg << "* Bunch global starting position:" << endl;
419 *gmsg << "* RINIT = " << referenceR << " [mm]" << endl;
420 *gmsg << "* PHIINIT = " << referenceTheta << " [deg]" << endl;
421 *gmsg << "* ZINIT = " << referenceZ << " [mm]" << endl;
422 *gmsg << endl;
423 *gmsg << "* Bunch global starting momenta:" << endl;
424 *gmsg << "* Initial gamma = " << itsReference.getGamma() << endl;
425 *gmsg << "* Initial beta = " << itsReference.getBeta() << endl;
426 *gmsg << "* Reference total momentum = " << referencePtot << " [beta gamma]" << endl;
427 *gmsg << "* Reference azimuthal momentum (Pt) = " << referencePt << " [beta gamma]" << endl;
428 *gmsg << "* Reference radial momentum (Pr) = " << referencePr << " [beta gamma]" << endl;
429 *gmsg << "* Reference axial momentum (Pz) = " << referencePz << " [beta gamma]" << endl;
430 *gmsg << endl;
431
432 double sym = cycl_m->getSymmetry();
433 *gmsg << "* " << sym << "-fold field symmetry " << endl;
434
435 // ckr: this just returned the default value as defined in Component.h
436 // double rff = cycl_m->getRfFrequ(0);
437 // *gmsg << "* Rf frequency= " << rff << " [MHz]" << endl;
438
439 std::string fmfn = cycl_m->getFieldMapFN();
440 *gmsg << "* Field map file = '" << fmfn << "'" << endl;
441
442 std::string type = cycl_m->getCyclotronType();
443 *gmsg << "* Type of cyclotron = " << type << " " << endl;
444
445 double rmin = cycl_m->getMinR();
446 double rmax = cycl_m->getMaxR();
447 *gmsg << "* Radial aperture = " << rmin << " ... " << rmax<<" [m] "<< endl;
448
449 double zmin = cycl_m->getMinZ();
450 double zmax = cycl_m->getMaxZ();
451 *gmsg << "* Vertical aperture = " << zmin << " ... " << zmax<<" [m]"<< endl;
452
453 double h = cycl_m->getCyclHarm();
454 *gmsg << "* Number of trimcoils = " << cycl_m->getNumberOfTrimcoils() << endl;
455 *gmsg << "* Harmonic number h = " << h << " " << endl;
456
457 std::vector<double> rffrequ = cycl_m->getRfFrequ();
458 *gmsg << "* RF frequency = " << Util::doubleVectorToString(rffrequ) << " [MHz]" << endl;
459
462 std::vector<double> rfphi = cycl_m->getRfPhi();
463 for (size_t i = 0; i < rfphi.size(); ++i) {
464 rfphi[i] = rfphi[i] * Units::rad2deg;
465 }
466 *gmsg << "* RF inital phase = " << Util::doubleVectorToString(rfphi) << " [deg]" << endl;
467
468 std::vector<double> escale = cycl_m->getEScale();
469 *gmsg << "* E-field scale factor = " << Util::doubleVectorToString(escale) << endl;
470
471 std::vector<bool> superpose = cycl_m->getSuperpose();
472 *gmsg << "* Superpose electric field maps -> " << Util::boolVectorToUpperString(superpose) << endl;
473 }
474
475 // Read in cyclotron field maps
477
478 double BcParameter[8] = {};
479 BcParameter[0] = Units::mm2m * cycl_m->getRmin();
480 BcParameter[1] = Units::mm2m * cycl_m->getRmax();
481
482 // Store inner radius and outer radius of cyclotron field map in the list
484}
485
492 *gmsg << "* ----------------------------- Collimator ------------------------------- *" << endl;
493
494 CCollimator* elptr = dynamic_cast<CCollimator*>(coll.clone());
495 myElements.push_back(elptr);
496
497 *gmsg << "* Name = " << elptr->getName() << endl;
498
499 double xstart = elptr->getXStart();
500 *gmsg << "* XStart = " << xstart << " [m]" << endl;
501
502 double xend = elptr->getXEnd();
503 *gmsg << "* XEnd = " << xend << " [m]" << endl;
504
505 double ystart = elptr->getYStart();
506 *gmsg << "* YStart = " << ystart << " [m]" << endl;
507
508 double yend = elptr->getYEnd();
509 *gmsg << "* YEnd = " << yend << " [m]" << endl;
510
511 double zstart = elptr->getZStart();
512 *gmsg << "* ZStart = " << zstart << " [m]" << endl;
513
514 double zend = elptr->getZEnd();
515 *gmsg << "* ZEnd = " << zend << " [m]" << endl;
516
517 double width = elptr->getWidth();
518 *gmsg << "* Width = " << width << " [m]" << endl;
519
520 elptr->initialise(itsBunch_m);
521
522 double BcParameter[8] = {};
523
524 BcParameter[0] = xstart;
525 BcParameter[1] = xend;
526 BcParameter[2] = ystart;
527 BcParameter[3] = yend;
528 BcParameter[4] = width;
529
530 buildupFieldList(BcParameter, ElementType::CCOLLIMATOR, elptr);
531}
532
539 *gmsg << "In Corrector; L= " << corr.getElementLength() << endl;
540 myElements.push_back(dynamic_cast<Corrector*>(corr.clone()));
541}
542
549 *gmsg << "In Degrader; L= " << deg.getElementLength() << endl;
550 myElements.push_back(dynamic_cast<Degrader*>(deg.clone()));
551
552}
553
560 *gmsg << "In drift L= " << drift.getElementLength() << endl;
561 myElements.push_back(dynamic_cast<Drift*>(drift.clone()));
562}
563
570
571}
572
579 if (opalRing_m == nullptr)
580 throw OpalException(
581 "ParallelCylcotronTracker::visitOffset",
582 "Attempt to place an offset when Ring not defined");
583 Offset* offNonConst = const_cast<Offset*>(&off);
587}
588
595 // *gmsg << "In Marker; L= " << marker.getElementLength() << endl;
596 myElements.push_back(dynamic_cast<Marker*>(marker.clone()));
597 // Do nothing.
598}
599
606 // *gmsg << "In Monitor; L= " << corr.getElementLength() << endl;
607 myElements.push_back(dynamic_cast<Monitor*>(corr.clone()));
608 // applyDrift(flip_s * corr.getElementLength());
609}
610
617 *gmsg << "In Multipole; L= " << mult.getElementLength() << " however the element is missing " << endl;
618 myElements.push_back(dynamic_cast<Multipole*>(mult.clone()));
619}
620
627 *gmsg << "Adding MultipoleT" << endl;
628 if (opalRing_m != nullptr) {
630 } else {
631 throw OpalException("ParallelCyclotronTracker::visitMultipoleT",
632 "Need to define a RINGDEFINITION to use MultipoleT element");
633 }
634 myElements.push_back(dynamic_cast<MultipoleT*>(multT.clone()));
635}
636
643 *gmsg << "Adding MultipoleTStraight" << endl;
644 if (opalRing_m != nullptr) {
645 opalRing_m->appendElement(multTstraight);
646 } else {
647 throw OpalException("ParallelCyclotronTracker::visitMultipoleTStraight",
648 "Need to define a RINGDEFINITION to use MultipoleTStraight element");
649 }
650 myElements.push_back(dynamic_cast<MultipoleTStraight*>(multTstraight.clone()));
651}
652
659 *gmsg << "Adding MultipoleTCurvedConstRadius" << endl;
660 if (opalRing_m != nullptr) {
661 opalRing_m->appendElement(multTccurv);
662 } else {
663 throw OpalException("ParallelCyclotronTracker::visitMultipoleTCurvedConstRadius",
664 "Need to define a RINGDEFINITION to use MultipoleTCurvedConstRadius element");
665 }
666 myElements.push_back(dynamic_cast<MultipoleTCurvedConstRadius*>(multTccurv.clone()));
667}
668
675 *gmsg << "Adding MultipoleTCurvedVarRadius" << endl;
676 if (opalRing_m != nullptr) {
677 opalRing_m->appendElement(multTvcurv);
678 } else {
679 throw OpalException("ParallelCyclotronTracker::visitMultipoleTCurvedVarRadius",
680 "Need to define a RINGDEFINITION to use MultipoleTCurvedVarRadius element");
681 }
682 myElements.push_back(dynamic_cast<MultipoleTCurvedVarRadius*>(multTvcurv.clone()));
683}
684
691 *gmsg << "* ----------------------------- Probe ------------------------------------ *" << endl;
692
693 Probe* elptr = dynamic_cast<Probe*>(prob.clone());
694 myElements.push_back(elptr);
695
696 *gmsg << "* Name = " << elptr->getName() << endl;
697
698 double xstart = elptr->getXStart();
699 *gmsg << "* XStart = " << xstart << " [m]" << endl;
700
701 double xend = elptr->getXEnd();
702 *gmsg << "* XEnd = " << xend << " [m]" << endl;
703
704 double ystart = elptr->getYStart();
705 *gmsg << "* YStart = " << ystart << " [m]" << endl;
706
707 double yend = elptr->getYEnd();
708 *gmsg << "* YEnd = " << yend << " [m]" << endl;
709
710 // initialise, do nothing
711 elptr->initialise(itsBunch_m);
712
713 double BcParameter[8] = {};
714
715 BcParameter[0] = xstart;
716 BcParameter[1] = xend;
717 BcParameter[2] = ystart;
718 BcParameter[3] = yend;
719 BcParameter[4] = 1 ; // width
720
721 buildupFieldList(BcParameter, ElementType::PROBE, elptr);
722}
723
730 *gmsg << "In RBend; L= " << bend.getElementLength() << " however the element is missing " << endl;
731 myElements.push_back(dynamic_cast<RBend*>(bend.clone()));
732}
733
740 *gmsg << "* ----------------------------- RFCavity --------------------------------- * " << endl;
741
742 RFCavity* elptr = dynamic_cast<RFCavity*>(as.clone());
743 myElements.push_back(elptr);
744
745 if ( elptr->getCavityType() != CavityType::SGSW ) {
746 throw OpalException("ParallelCyclotronTracker::visitRFCavity",
747 "\"" + elptr->getCavityTypeString() + "\" is not valid \"TYPE\" for RFCavity.\n"
748 "The ParallelCyclotronTracker can only play with cyclotron type RF system currently...");
749 }
750
751 *gmsg << "* Name = " << elptr->getName() << endl;
752
753 std::string fmfn = elptr->getFieldMapFN();
754 *gmsg << "* RF Field map file = '" << fmfn << "'" << endl;
755
756 double rmin = elptr->getRmin();
757 *gmsg << "* Minimal radius of cavity = " << rmin << " [mm]" << endl;
758
759 double rmax = elptr->getRmax();
760 *gmsg << "* Maximal radius of cavity = " << rmax << " [mm]" << endl;
761
762 double rff = elptr->getCycFrequency();
763 *gmsg << "* RF frequency (2*pi*f) = " << rff << " [rad/s]" << endl;
764
765 double angle = elptr->getAzimuth();
766 *gmsg << "* Cavity azimuth position = " << angle << " [deg] " << endl;
767
768 double gap = elptr->getGapWidth();
769 *gmsg << "* Cavity gap width = " << gap << " [mm] " << endl;
770
771 double pdis = elptr->getPerpenDistance();
772 *gmsg << "* Cavity Shift distance = " << pdis << " [mm] " << endl;
773
774 double phi0 = elptr->getPhi0();
775 *gmsg << "* Initial RF phase (t=0) = " << phi0 << " [deg] " << endl;
776
777 /*
778 Setup time dependence and in case of no
779 timedependence use a polynom with a_0 = 1 and a_k = 0, k = 1,2,3.
780 */
781
782 std::shared_ptr<AbstractTimeDependence> freq_atd = nullptr;
783 std::shared_ptr<AbstractTimeDependence> ampl_atd = nullptr;
784 std::shared_ptr<AbstractTimeDependence> phase_atd = nullptr;
785
786 dvector_t unityVec;
787 unityVec.push_back(1.);
788 unityVec.push_back(0.);
789 unityVec.push_back(0.);
790 unityVec.push_back(0.);
791
792 std::string frequencyModelName = elptr->getFrequencyModelName();
793 std::string amplitudeModelName = elptr->getAmplitudeModelName();
794 std::string phaseModelName = elptr->getPhaseModelName();
795
796 if (!frequencyModelName.empty()) {
797 freq_atd = AbstractTimeDependence::getTimeDependence(frequencyModelName);
798 *gmsg << "* Variable frequency RF Model name " << frequencyModelName << endl;
799 } else {
800 freq_atd = std::shared_ptr<AbstractTimeDependence>(new PolynomialTimeDependence(unityVec));
801 }
802
803 if (!amplitudeModelName.empty()) {
804 ampl_atd = AbstractTimeDependence::getTimeDependence(amplitudeModelName);
805 *gmsg << "* Variable amplitude RF Model name " << amplitudeModelName << endl;
806 } else {
807 ampl_atd = std::shared_ptr<AbstractTimeDependence>(new PolynomialTimeDependence(unityVec));
808 }
809
810 if (!phaseModelName.empty()) {
811 phase_atd = AbstractTimeDependence::getTimeDependence(phaseModelName);
812 *gmsg << "* Variable phase RF Model name " << phaseModelName << endl;
813 } else {
814 phase_atd = std::shared_ptr<AbstractTimeDependence>(new PolynomialTimeDependence(unityVec));
815 }
816 // read cavity voltage profile data from file.
817 elptr->initialise(itsBunch_m, freq_atd, ampl_atd, phase_atd);
818
819 double BcParameter[8] = {};
820
821 BcParameter[0] = Units::mm2m * rmin;
822 BcParameter[1] = Units::mm2m * rmax;
823 BcParameter[2] = Units::mm2m * pdis;
824 BcParameter[3] = angle;
825
826 buildupFieldList(BcParameter, ElementType::RFCAVITY, elptr);
827}
828
834 *gmsg << "* ----------------------------- Ring ------------------------------------- *" << endl;
835
836 delete opalRing_m;
837
838 opalRing_m = dynamic_cast<Ring*>(ring.clone());
839
840 myElements.push_back(opalRing_m);
841
843
847
848 if (referenceTheta <= -180.0 || referenceTheta > 180.0) {
849 throw OpalException("Error in ParallelCyclotronTracker::visitRing",
850 "PHIINIT is out of [-180, 180)!");
851 }
852
853 referenceZ = 0.0;
854 referencePz = 0.0;
855
858
859 if (referencePtot < 0.0)
860 referencePt *= -1.0;
861
864
865 double BcParameter[8] = {}; // zero initialise array
866
868
869 // Finally print some diagnostic
870 *gmsg << "* Initial beam radius = " << referenceR << " [mm] " << endl;
871 *gmsg << "* Initial gamma = " << itsReference.getGamma() << endl;
872 *gmsg << "* Initial beta = " << itsReference.getBeta() << endl;
873 *gmsg << "* Total reference momentum = " << referencePtot << " [beta gamma]" << endl;
874 *gmsg << "* Reference azimuthal momentum = " << referencePt << " [beta gamma]" << endl;
875 *gmsg << "* Reference radial momentum = " << referencePr << " [beta gamma]" << endl;
876 *gmsg << "* " << opalRing_m->getSymmetry() << " fold field symmetry " << endl;
877 *gmsg << "* Harmonic number h = " << opalRing_m->getHarmonicNumber() << " " << endl;
878}
879
886 *gmsg << "In SBend; L = " << bend.getElementLength() << " however the element is missing " << endl;
887 myElements.push_back(dynamic_cast<SBend*>(bend.clone()));
888}
889
891 *gmsg << "Adding SBend3D" << endl;
892 if (opalRing_m != nullptr)
894 else
895 throw OpalException("ParallelCyclotronTracker::visitSBend3D",
896 "Need to define a RINGDEFINITION to use SBend3D element");
897}
898
900 *gmsg << "Adding ScalingFFAMagnet" << endl;
901 if (opalRing_m != nullptr) {
902 ScalingFFAMagnet* newBend = bend.clone(); // setup the end field, if required
903 newBend->setupEndField();
904 opalRing_m->appendElement(*newBend);
905 } else {
906 throw OpalException("ParallelCyclotronTracker::visitScalingFFAMagnet",
907 "Need to define a RINGDEFINITION to use ScalingFFAMagnet element");
908 }
909}
910
917 *gmsg << "* ----------------------------- Septum ----------------------------------- *" << endl;
918
919 Septum* elptr = dynamic_cast<Septum*>(sept.clone());
920 myElements.push_back(elptr);
921
922 *gmsg << "* Name = " << elptr->getName() << endl;
923
924 double xstart = elptr->getXStart();
925 *gmsg << "* XStart = " << xstart << " [m]" << endl;
926
927 double xend = elptr->getXEnd();
928 *gmsg << "* XEnd = " << xend << " [m]" << endl;
929
930 double ystart = elptr->getYStart();
931 *gmsg << "* YStart = " << ystart << " [m]" << endl;
932
933 double yend = elptr->getYEnd();
934 *gmsg << "* YEnd = " << yend << " [m]" << endl;
935
936 double width = elptr->getWidth();
937 *gmsg << "* Width = " << width << " [m]" << endl;
938
939 // initialise, do nothing
940 elptr->initialise(itsBunch_m);
941
942 double BcParameter[8] = {};
943
944 BcParameter[0] = xstart;
945 BcParameter[1] = xend;
946 BcParameter[2] = ystart;
947 BcParameter[3] = yend;
948 BcParameter[4] = width;
949
950 buildupFieldList(BcParameter, ElementType::SEPTUM, elptr);
951}
952
959 myElements.push_back(dynamic_cast<Solenoid*>(solenoid.clone()));
960 Component* elptr = *(--myElements.end());
961 if (!elptr->hasAttribute("ELEMEDGE")) {
962 *gmsg << "Solenoid: no position of the element given!" << endl;
963 return;
964 }
965}
966
973 *gmsg << "* ----------------------------- Stripper --------------------------------- *" << endl;
974
975 Stripper* elptr = dynamic_cast<Stripper*>(stripper.clone());
976 myElements.push_back(elptr);
977
978 *gmsg << "* Name = " << elptr->getName() << endl;
979
980 double xstart = elptr->getXStart();
981 *gmsg << "* XStart = " << xstart << " [m]" << endl;
982
983 double xend = elptr->getXEnd();
984 *gmsg << "* XEnd = " << xend << " [m]" << endl;
985
986 double ystart = elptr->getYStart();
987 *gmsg << "* YStart = " << ystart << " [m]" << endl;
988
989 double yend = elptr->getYEnd();
990 *gmsg << "* YEnd = " << yend << " [m]" << endl;
991
992 double opcharge = elptr->getOPCharge();
993 *gmsg << "* Charge of outcoming particle = +e * " << opcharge << endl;
994
995 double opmass = elptr->getOPMass();
996 *gmsg << "* Mass of the outcoming particle = " << opmass << " [GeV/c^2]" << endl;
997
998 bool stop = elptr->getStop();
999 *gmsg << std::boolalpha
1000 << "* Particles stripped will be deleted after interaction -> "
1001 << stop << endl;
1002
1003 elptr->initialise(itsBunch_m);
1004
1005 double BcParameter[8] = {};
1006
1007 BcParameter[0] = xstart;
1008 BcParameter[1] = xend;
1009 BcParameter[2] = ystart;
1010 BcParameter[3] = yend;
1011 BcParameter[4] = 1; // width
1012 BcParameter[5] = opcharge;
1013 BcParameter[6] = opmass;
1014
1015 buildupFieldList(BcParameter, ElementType::STRIPPER, elptr);
1016}
1017
1024 *gmsg << "* ----------------------------- Vacuum ----------------------------------- *" << endl;
1025
1026 Vacuum* elptr = dynamic_cast<Vacuum*>(vac.clone());
1027 myElements.push_back(elptr);
1028
1029 double BcParameter[8] = {};
1030
1031 std::string gas = elptr->getResidualGasName();
1032 *gmsg << "* Residual gas = " << gas << endl;
1033
1034 double pressure = elptr->getPressure();
1035 if ( boost::filesystem::exists(elptr->getPressureMapFN()) ) {
1036 std::string pmfn = elptr->getPressureMapFN();
1037 *gmsg << "* Pressure field map file = '" << pmfn << "'" << endl;
1038 *gmsg << "* Default pressure = " << pressure << " [mbar]" << endl;
1039 } else {
1040 *gmsg << "* Pressure = " << pressure << " [mbar]" << endl;
1041 }
1042 double pscale = elptr->getPScale();
1043
1044 double temperature = elptr->getTemperature();
1045 *gmsg << "* Temperature = " << temperature << " [K]" << endl;
1046
1047 bool stop = elptr->getStop();
1048 *gmsg << std::boolalpha
1049 << "* Particles stripped will be deleted after interaction -> "
1050 << stop << endl;
1051
1052 elptr->initialise(itsBunch_m);
1053
1054 BcParameter[0] = pressure;
1055 BcParameter[1] = pscale;
1056 BcParameter[2] = temperature;
1057
1058 buildupFieldList(BcParameter, ElementType::VACUUM, elptr);
1059}
1060
1067 *gmsg << "Adding Variable RF Cavity" << endl;
1068 if (opalRing_m != nullptr)
1070 else
1071 throw OpalException("ParallelCyclotronTracker::visitVariableRFCavity",
1072 "Need to define a RINGDEFINITION to use VariableRFCavity element");
1073}
1074
1081 (const VariableRFCavityFringeField& cav) {
1082 *gmsg << "Adding Variable RF Cavity with Fringe Field" << endl;
1083 if (opalRing_m != nullptr)
1085 else
1086 throw OpalException("ParallelCyclotronTracker::visitVariableRFCavityFringeField",
1087 "Need to define a RINGDEFINITION to use VariableRFCavity element");
1088}
1089
1096 *gmsg << "Adding Vertical FFA Magnet" << endl;
1097 if (opalRing_m != nullptr)
1099 else
1100 throw OpalException("ParallelCyclotronTracker::visitVerticalFFAMagnet",
1101 "Need to define a RINGDEFINITION to use VerticalFFAMagnet element");
1102}
1103
1104
1112void ParallelCyclotronTracker::buildupFieldList(double BcParameter[], ElementType elementType, Component *elptr) {
1114
1115 type_pair *localpair = new type_pair();
1116 localpair->first = elementType;
1117
1118 for (int i = 0; i < 8; i++)
1119 *(((localpair->second).first) + i) = *(BcParameter + i);
1120
1121 (localpair->second).second = elptr;
1122
1123 // always put cyclotron as the first element in the list.
1124 if (elementType == ElementType::RING || elementType == ElementType::CYCLOTRON) {
1125 sindex = FieldDimensions.begin();
1126 } else {
1127 sindex = FieldDimensions.end();
1128 }
1129 FieldDimensions.insert(sindex, localpair);
1130}
1131
1138 const FlaggedBeamline* fbl = static_cast<const FlaggedBeamline*>(&bl);
1139 fbl->iterate(*this, false);//*dynamic_cast<BeamlineVisitor *>(this), false);
1140}
1141
1143 int nlp = itsBunch_m->getLocalNum();
1144 int minnlp = 0;
1145 int maxnlp = 111111;
1146 reduce(nlp, minnlp, OpMinAssign());
1147 reduce(nlp, maxnlp, OpMaxAssign());
1148 *gmsg << s << "min local particle number: "<< minnlp << endl;
1149 *gmsg << "* max local particle number: " << maxnlp << endl;
1150}
1151
1152
1154 /*
1155 Initialize common variables and structures
1156 for the integrators
1157 */
1158 if (OpalData::getInstance()->inRestartRun()) {
1160 }
1161
1162 step_m = 0;
1163 restartStep0_m = 0;
1164 turnnumber_m = 1;
1165 azimuth_m = -1.0;
1166 prevAzimuth_m = -1.0;
1167
1168 // Record how many bunches have already been injected. ONLY FOR MPM
1169 if (isMultiBunch())
1170 mbHandler_m->setNumBunch(itsBunch_m->getNumBunch());
1171
1172 itsBeamline->accept(*this);
1173 if (opalRing_m != nullptr)
1175
1176 // Display the selected elements
1177 *gmsg << "* ------------------------------------------------------------------------ *" << endl;
1178 *gmsg << "* The selected Beam line elements are :" << endl;
1179 for (auto fd : FieldDimensions) {
1180 ElementType type = fd->first;
1181 *gmsg << "* -> " << ElementBase::getTypeString(type) << endl;
1182 if (type == ElementType::RFCAVITY) {
1183 RFCavity* cav = static_cast<RFCavity*>((fd->second).second);
1184 CavityCrossData ccd = {cav, cav->getSinAzimuth(), cav->getCosAzimuth(), cav->getPerpenDistance() * Units::mm2m};
1185 cavCrossDatas_m.push_back(ccd);
1186 } else if ( type == ElementType::CCOLLIMATOR ||
1190 PluginElement* element = static_cast<PluginElement*>((fd->second).second);
1191 pluginElements_m.push_back(element);
1192 }
1193 }
1194 *gmsg << "* ------------------------------------------------------------------------ *" << endl << endl;
1195
1196 // Get BoundaryGeometry that is already initialized
1198 if (bgf_m)
1199 lossDs_m = std::unique_ptr<LossDataSink>(new LossDataSink(bgf_m->getOpalName(),!Options::asciidump));
1200
1201 // External field arrays for dumping
1202 for (int k = 0; k < 2; k++)
1203 FDext_m[k] = Vector_t(0.0, 0.0, 0.0);
1204
1205 extE_m = Vector_t(0.0, 0.0, 0.0);
1206 extB_m = Vector_t(0.0, 0.0, 0.0);
1207 DumpFields::writeFields((((*FieldDimensions.begin())->second).second));
1208 DumpEMFields::writeFields((((*FieldDimensions.begin())->second).second));
1209
1211 this,
1212 std::placeholders::_1,
1213 std::placeholders::_2,
1214 std::placeholders::_3,
1215 std::placeholders::_4);
1216
1217 switch ( stepper_m ) {
1219 *gmsg << "* 2nd order Leap-Frog integrator" << endl;
1220 itsStepper_mp.reset(new LF2<function_t>(func));
1221 break;
1222 }
1224 *gmsg << "* Multiple time stepping (MTS) integrator" << endl;
1225 break;
1226 }
1228 default: {
1229 *gmsg << "* 4th order Runge-Kutta integrator" << endl;
1230 itsStepper_mp.reset(new RK4<function_t>(func));
1231 break;
1232 }
1233 }
1234
1236 MtsTracker();
1237 } else {
1239 }
1240
1241 *gmsg << "* ------------------------------------------------------------------------ *" << endl;
1242 *gmsg << "* Finalizing i.e. write data and close files :" << endl;
1243 for (auto fd : FieldDimensions) {
1244 ((fd->second).second)->finalise();
1245 }
1246 *gmsg << "* ------------------------------------------------------------------------ *" << endl;
1247}
1248
1249
1251 /*
1252 * variable unit meaning
1253 * ------------------------------------------------
1254 * t [ns] time
1255 * dt [ns] time step
1256 * oldReferenceTheta [rad] azimuth angle
1257 * itsBunch_m->R [m] particle position
1258 *
1259 */
1260
1261 double t = 0, dt = 0, oldReferenceTheta = 0;
1262 std::tie(t, dt, oldReferenceTheta) = initializeTracking_m();
1263
1264 int const numSubsteps = std::max(Options::mtsSubsteps, 1);
1265 *gmsg << "* MTS: Number of substeps per step is " << numSubsteps << endl;
1266
1267 double const dt_inner = dt / double(numSubsteps);
1268 *gmsg << "* MTS: The inner time step is therefore " << dt_inner << " [ns]" << endl;
1269
1270// int SteptoLastInj = itsBunch_m->getSteptoLastInj();
1271
1272 bool flagTransition = false; // flag to determine when to transit from single-bunch to multi-bunches mode
1273
1274 *gmsg << "* ---------------------- Start tracking ---------------------------------- *" << endl;
1275
1276 if ( itsBunch_m->hasFieldSolver() )
1278
1279 for (; (step_m < maxSteps_m) && (itsBunch_m->getTotalNum()>0); step_m++) {
1280
1281 bool finishedTurn = false;
1282
1283 if (step_m % Options::sptDumpFreq == 0) {
1285 }
1286
1288
1289 // First half kick from space charge force
1290 if (itsBunch_m->hasFieldSolver()) {
1291 kick(0.5 * dt);
1292 }
1293
1294 // Substeps for external field integration
1295 for (int n = 0; n < numSubsteps; ++n)
1296 borisExternalFields(dt_inner);
1297
1298 // bunch injection
1299 // TODO: Where is correct location for this piece of code? Beginning/end of step? Before field solve?
1300 injectBunch(flagTransition);
1301
1302 if ( itsBunch_m->hasFieldSolver() ) {
1304 } else {
1305 // if field solver is not available , only update bunch, to transfer particles between nodes if needed,
1306 // reset parameters such as LocalNum, initialTotalNum_m.
1307 // INFOMSG("No space charge Effects are included!"<<endl;);
1308 if ((step_m % Options::repartFreq * 100) == 0) { //TODO: why * 100?
1309 Vector_t const meanP = calcMeanP();
1310 double const phi = calculateAngle(meanP(0), meanP(1)) - 0.5 * Physics::pi;
1311 Vector_t const meanR = calcMeanR();
1312 globalToLocal(itsBunch_m->R, phi, meanR);
1314 repartition();
1315 localToGlobal(itsBunch_m->R, phi, meanR);
1316 }
1317 }
1318
1319 // Second half kick from space charge force
1321 kick(0.5 * dt);
1322
1323 // recalculate bingamma and reset the BinID for each particles according to its current gamma
1324 if (isMultiBunch() && (step_m % Options::rebinFreq == 0)) {
1325 mbHandler_m->updateParticleBins(itsBunch_m);
1326 }
1327
1328 // dump some data after one push in single particle tracking
1329 if ( mode_m == TrackingMode::SINGLE ) {
1330
1331 unsigned int i = 0;
1332
1333 double temp_meanTheta = calculateAngle2(itsBunch_m->R[i](0),
1334 itsBunch_m->R[i](1)); // [-pi, pi]
1335
1337 temp_meanTheta, finishedTurn);
1338
1340 oldReferenceTheta, temp_meanTheta);
1341
1342 oldReferenceTheta = temp_meanTheta;
1343 } else if ( mode_m == TrackingMode::BUNCH ) {
1344 // both for single bunch and multi-bunch
1345 // avoid dump at the first step
1346 // finishedTurn has not been changed in first push
1347 if ( isTurnDone() ) {
1348 ++turnnumber_m;
1349 finishedTurn = true;
1350
1351 *gmsg << endl;
1352 *gmsg << "*** Finished turn " << turnnumber_m - 1
1353 << ", Total number of live particles: "
1354 << itsBunch_m->getTotalNum() << endl;
1355 }
1356
1357 // Recalculate bingamma and reset the BinID for each particles according to its current gamma
1358 if (isMultiBunch() && (step_m % Options::rebinFreq == 0)) {
1359 mbHandler_m->updateParticleBins(itsBunch_m);
1360 }
1361 }
1362
1363 // printing + updating bunch parameters + updating t
1364 update_m(t, dt, finishedTurn);
1365
1366 }
1367
1368 // Some post-integration stuff
1369 *gmsg << endl;
1370 *gmsg << "* ---------------------- DONE TRACKING PARTICLES ------------------------- *" << endl;
1371
1372
1373 //FIXME
1374 dvector_t Ttime = dvector_t();
1375 dvector_t Tdeltr = dvector_t();
1376 dvector_t Tdeltz = dvector_t();
1377 ivector_t TturnNumber = ivector_t();
1378
1379 finalizeTracking_m(Ttime, Tdeltr, Tdeltz, TturnNumber);
1380}
1381
1382
1384 /*
1385 * variable unit meaning
1386 * ------------------------------------------------
1387 * t [ns] time
1388 * dt [ns] time step
1389 * oldReferenceTheta [rad] azimuth angle
1390 * itsBunch_m->R [m] particle position
1391 *
1392 */
1393 // Generic Tracker that has three modes defined by timeIntegrator_m:
1394 // 0 --- RK-4 (default)
1395 // 1 --- LF-2
1396 // (2 --- MTS ... not yet implemented in generic tracker)
1397 // mbHandler_m->getNumBunch() determines the number of bunches in multibunch mode (MBM, 1 for OFF)
1398 // Total number of particles determines single particle mode (SPM, 1 particle) or
1399 // Static Equilibrium Orbit Mode (SEO, 2 particles)
1400
1401 double t = 0, dt = 0, oldReferenceTheta = 0;
1402 std::tie(t, dt, oldReferenceTheta) = initializeTracking_m();
1403
1404 // If initialTotalNum_m = 2, trigger SEO mode and prepare for transverse tuning calculation
1405 // Where is the IF here? -DW
1406 dvector_t Ttime, Tdeltr, Tdeltz;
1407 ivector_t TturnNumber;
1408
1409 // Apply the plugin elements: probe, collimator, stripper, septum once before first step
1410 bool flagNeedUpdate = applyPluginElements(dt);
1411
1412 // Destroy particles if they are marked as Bin = -1 in the plugin elements
1413 // or out of global aperture
1414 deleteParticle(flagNeedUpdate);
1415
1416 /********************************
1417 * Main integration loop *
1418 ********************************/
1419 *gmsg << endl;
1420 *gmsg << "* ---------------------- Start tracking ---------------------------------- *" << endl;
1421
1422 for (; (step_m < maxSteps_m) && (itsBunch_m->getTotalNum()>0); step_m++) {
1423
1424 bool finishedTurn = false;
1425
1426 switch (mode_m) {
1427 case TrackingMode::SEO: {
1428 // initialTotalNum_m == 2
1429 seoMode_m(t, dt, finishedTurn, Ttime, Tdeltr, Tdeltz, TturnNumber);
1430 break;
1431 }
1432 case TrackingMode::SINGLE: {
1433 // initialTotalNum_m == 1
1434 singleMode_m(t, dt, finishedTurn, oldReferenceTheta);
1435 break;
1436 }
1437 case TrackingMode::BUNCH: {
1438 // initialTotalNum_m > 2
1439 // Start Tracking for number of particles > 2 (i.e. not single and not SEO mode)
1440 bunchMode_m(t, dt, finishedTurn);
1441 break;
1442 }
1444 default:
1445 throw OpalException("ParallelCyclotronTracker::GenericTracker()",
1446 "No such tracking mode.");
1447 }
1448 }
1449 // Update bunch and some parameters and output some info
1450 update_m(t, dt, finishedTurn);
1451
1452 } // end for: the integration is DONE after maxSteps_m steps or if all particles are lost!
1453
1454 // Some post-integration stuff
1455 *gmsg << endl;
1456 *gmsg << "* ---------------------- DONE TRACKING PARTICLES ------------------------- *" << endl;
1457
1458 finalizeTracking_m(Ttime, Tdeltr, Tdeltz, TturnNumber);
1459}
1460
1461bool ParallelCyclotronTracker::getFieldsAtPoint(const double& t, const size_t& Pindex, Vector_t& Efield, Vector_t& Bfield) {
1462
1463 bool outOfBound = this->computeExternalFields_m(Pindex, t, Efield, Bfield);
1464
1465 // For runs without space charge effects, override this step to save time
1466 if (itsBunch_m->hasFieldSolver()) {
1467
1468 // Don't do for reference particle
1469 if (itsBunch_m->ID[Pindex] != 0) {
1470
1471 // add external Field and self space charge field
1472 Efield += itsBunch_m->Ef[Pindex];
1473 Bfield += itsBunch_m->Bf[Pindex];
1474 }
1475 }
1476
1477 return outOfBound;
1478}
1479
1480
1493 RFCavity * rfcavity, double& Dold)
1494{
1495 bool flag = false;
1496 double sinx = rfcavity->getSinAzimuth();
1497 double cosx = rfcavity->getCosAzimuth();
1498 // TODO: Presumably this is still in mm, so for now, change to m -DW
1499 double PerpenDistance = Units::mm2m * rfcavity->getPerpenDistance();
1500 double distNew = (Rnew[0] * sinx - Rnew[1] * cosx) - PerpenDistance;
1501 double distOld = (Rold[0] * sinx - Rold[1] * cosx) - PerpenDistance;
1502 if (distOld > 0.0 && distNew <= 0.0) flag = true;
1503 // This parameter is used correct cavity phase
1504 Dold = Units::m2mm * distOld;
1505 return flag;
1506}
1507
1508bool ParallelCyclotronTracker::RFkick(RFCavity * rfcavity, const double t, const double dt, const int Pindex){
1509 // For OPAL 2.0: As long as the RFCavity is in mm, we have to change R to mm here -DW
1510 double radius = std::sqrt(std::pow(Units::m2mm * itsBunch_m->R[Pindex](0), 2.0) + std::pow(Units::m2mm * itsBunch_m->R[Pindex](1), 2.0)
1511 - std::pow(rfcavity->getPerpenDistance() , 2.0));
1512 double rmin = rfcavity->getRmin();
1513 double rmax = rfcavity->getRmax();
1514 double nomalRadius = (radius - rmin) / (rmax - rmin);
1515 double tempP[3];
1516 if (nomalRadius <= 1.0 && nomalRadius >= 0.0) {
1517
1518 for (int j = 0; j < 3; j++)
1519 tempP[j] = itsBunch_m->P[Pindex](j); //[px,py,pz] units: dimensionless
1520
1521 // here evaluate voltage and conduct momenta kicking;
1522 rfcavity->getMomentaKick(nomalRadius, tempP, t, dt, itsBunch_m->ID[Pindex], itsBunch_m->getM(), itsBunch_m->getQ()); // t : ns
1523
1524 for (int j = 0; j < 3; j++)
1525 itsBunch_m->P[Pindex](j) = tempP[j];
1526 return true;
1527 }
1528 return false;
1529}
1530
1531
1532struct adder {
1533 adder() : sum(0) {}
1534 double sum;
1535 void operator()(double x) { sum += x; }
1536};
1537
1551 int lastTurn, double& /*nur*/, double& /*nuz*/) {
1552 TUNE_class *tune;
1553
1554 int Ndat = t.size();
1555
1556 /*
1557 remove mean
1558 */
1559 double rsum = for_each(r.begin(), r.end(), adder()).sum;
1560
1561 for (int i = 0; i < Ndat; i++)
1562 r[i] -= rsum;
1563
1564 double zsum = for_each(z.begin(), z.end(), adder()).sum;
1565
1566 for (int i = 0; i < Ndat; i++)
1567 z[i] -= zsum;
1568 double ti = *(t.begin());
1569 double tf = t[t.size()-1];
1570 double T = (tf - ti);
1571
1572 t.clear();
1573 for (int i = 0; i < Ndat; i++) {
1574 t.push_back(i);
1575 }
1576
1577 T = t[Ndat-1];
1578
1579 *gmsg << endl;
1580 *gmsg << "* ************************************* nuR ******************************************* *" << endl;
1581 *gmsg << endl << "* ===> " << Ndat << " data points Ti=" << ti << " Tf= " << tf << " -> T= " << T << endl;
1582
1583 int nhis_lomb = 10;
1584 int stat = 0;
1585 // book tune class
1586 tune = new TUNE_class();
1587 stat = tune->lombAnalysis(t, r, nhis_lomb, T / lastTurn);
1588 if (stat != 0)
1589 *gmsg << "* TUNE: Lomb analysis failed" << endl;
1590 *gmsg << "* ************************************************************************************* *" << endl;
1591
1592 delete tune;
1593 tune = nullptr;
1594 // FIXME: FixMe: need to come from the inputfile
1595 nhis_lomb = 100;
1596
1597 if (zsum != 0.0) {
1598
1599 *gmsg << endl;
1600 *gmsg << "* ************************************* nuZ ******************************************* *" << endl;
1601 *gmsg << endl << "* ===> " << Ndat << " data points Ti=" << ti << " Tf= " << tf << " -> T= " << T << endl;
1602
1603 // book tune class
1604 tune = new TUNE_class();
1605 stat = tune->lombAnalysis(t, z, nhis_lomb, T / lastTurn);
1606 if (stat != 0)
1607 *gmsg << "* TUNE: Lomb analysis failed" << endl;
1608 *gmsg << "* ************************************************************************************* *" << endl;
1609
1610 delete tune;
1611 tune = nullptr;
1612 }
1613 return true;
1614}
1615
1616
1618 if (opalRing_m != nullptr)
1619 return opalRing_m->getHarmonicNumber();
1620 Cyclotron* elcycl = dynamic_cast<Cyclotron*>(((*FieldDimensions.begin())->second).second);
1621 if (elcycl != nullptr)
1622 return elcycl->getCyclHarm();
1623 throw OpalException("ParallelCyclotronTracker::getHarmonicNumber()",
1624 std::string("The first item in the FieldDimensions list does not ")
1625 +std::string("seem to be a Ring or a Cyclotron element"));
1626}
1627
1628
1630 Vector_t meanR(0.0, 0.0, 0.0);
1631
1632 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
1633 // take all particles if bunchNr <= -1
1634 if ( bunchNr > -1 && itsBunch_m->bunchNum[i] != bunchNr)
1635 continue;
1636
1637 for (int d = 0; d < 3; ++d) {
1638 meanR(d) += itsBunch_m->R[i](d);
1639 }
1640 }
1641
1642 reduce(meanR, meanR, OpAddAssign());
1643
1644 size_t n = itsBunch_m->getTotalNum();
1645
1646 if ( bunchNr > -1 )
1647 n = itsBunch_m->getTotalNumPerBunch(bunchNr);
1648
1649 return meanR / Vector_t(n);
1650}
1651
1653 Vector_t meanP(0.0, 0.0, 0.0);
1654
1655 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
1656 for (int d = 0; d < 3; ++d) {
1657 meanP(d) += itsBunch_m->P[i](d);
1658 }
1659 }
1660
1661 reduce(meanP, meanP, OpAddAssign());
1662 return meanP / Vector_t(itsBunch_m->getTotalNum());
1663}
1664
1666 if ((step_m % Options::repartFreq) == 0) {
1671 }
1672}
1673
1675 double phi, Vector_t const translationToGlobal) {
1677 particleVectors -= translationToGlobal;
1678
1679 Tenzor<double, 3> const rotation( std::cos(phi), std::sin(phi), 0,
1680 -std::sin(phi), std::cos(phi), 0,
1681 0, 0, 1); // clockwise rotation
1682
1683 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
1684
1685 particleVectors[i] = dot(rotation, particleVectors[i]);
1686 }
1688}
1689
1691 double phi, Vector_t const translationToGlobal) {
1693 Tenzor<double, 3> const rotation(std::cos(phi), -std::sin(phi), 0,
1694 std::sin(phi), std::cos(phi), 0,
1695 0, 0, 1); // counter-clockwise rotation
1696
1697 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
1698
1699 particleVectors[i] = dot(rotation, particleVectors[i]);
1700 }
1701
1702 particleVectors += translationToGlobal;
1704}
1705
1706
1708 Quaternion_t const quaternion,
1709 Vector_t const meanR) {
1711
1712 // Translation from global to local
1713 particleVectors -= meanR;
1714
1715 // Rotation to align P_mean with x-axis
1716 rotateWithQuaternion(particleVectors, quaternion);
1718}
1719
1721 Quaternion_t const quaternion,
1722 Vector_t const meanR) {
1724 // Reverse the quaternion by multiplying the axis components (x,y,z) with -1
1725 Quaternion_t reverseQuaternion = quaternion * -1.0;
1726 reverseQuaternion(0) *= -1.0;
1727
1728 // Rotation back to original P_mean direction
1729 rotateWithQuaternion(particleVectors, reverseQuaternion);
1730
1731 // Translation from local to global
1732 particleVectors += meanR;
1734}
1735
1737 double const phi,
1738 double const psi,
1739 Vector_t const meanR) {
1740
1742 //double const tolerance = 1.0e-4; // TODO What is a good angle threshold? -DW
1743
1744 // Translation from global to local
1745 particleVectors -= meanR;
1746
1747 // Rotation to align P_mean with x-axis
1748 rotateAroundZ(particleVectors, phi);
1749
1750 // If theta is large enough (i.e. there is significant momentum in z direction)
1751 // rotate around x-axis next
1752 //if (fabs(psi) > tolerance)
1753 rotateAroundX(particleVectors, psi);
1755}
1756
1758 double const phi,
1759 double const psi,
1760 Vector_t const meanR) {
1761
1763 //double const tolerance = 1.0e-4; // TODO What is a good angle threshold? -DW
1764
1765 // Translation from global to local
1766 myVector -= meanR;
1767
1768 // Rotation to align P_mean with x-axis
1769 rotateAroundZ(myVector, phi);
1770
1771 // If theta is large enough (i.e. there is significant momentum in z direction)
1772 // rotate around x-axis next
1773 //if (fabs(psi) > tolerance)
1774 rotateAroundX(myVector, psi);
1776}
1777
1779 double const phi,
1780 double const psi,
1781 Vector_t const meanR) {
1782
1784 //double const tolerance = 1.0e-4; // TODO What is a good angle threshold? -DW
1785
1786 // If theta is large enough (i.e. there is significant momentum in z direction)
1787 // rotate back around x-axis next
1788 //if (fabs(psi) > tolerance)
1789 rotateAroundX(particleVectors, -psi);
1790
1791 // Rotation to align P_mean with x-axis
1792 rotateAroundZ(particleVectors, -phi);
1793
1794 // Translation from local to global
1795 particleVectors += meanR;
1797}
1798
1800 double const phi,
1801 double const psi,
1802 Vector_t const meanR) {
1803
1805 //double const tolerance = 1.0e-4; // TODO What is a good angle threshold? -DW
1806
1807 // If theta is large enough (i.e. there is significant momentum in z direction)
1808 // rotate back around x-axis next
1809 //if (fabs(psi) > tolerance)
1810 rotateAroundX(myVector, -psi);
1811
1812 // Rotation to align P_mean with x-axis
1813 rotateAroundZ(myVector, -phi);
1814
1815 // Translation from local to global
1816 myVector += meanR;
1818}
1819
1821
1822 Vector_t const quaternionVectorComponent = Vector_t(quaternion(1), quaternion(2), quaternion(3));
1823 double const quaternionScalarComponent = quaternion(0);
1824
1825 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
1826
1827 particleVectors[i] = 2.0f * dot(quaternionVectorComponent, particleVectors[i]) * quaternionVectorComponent +
1828 (quaternionScalarComponent * quaternionScalarComponent -
1829 dot(quaternionVectorComponent, quaternionVectorComponent)) * particleVectors[i] + 2.0f *
1830 quaternionScalarComponent * cross(quaternionVectorComponent, particleVectors[i]);
1831 }
1832}
1833
1835
1836 double tolerance = 1.0e-10;
1837 double length2 = dot(quaternion, quaternion);
1838
1839 if (std::abs(length2) > tolerance && std::abs(length2 - 1.0f) > tolerance) {
1840
1841 double length = std::sqrt(length2);
1842 quaternion /= length;
1843 }
1844}
1845
1847
1848 double tolerance = 1.0e-10;
1849 double length2 = dot(vector, vector);
1850
1851 if (std::abs(length2) > tolerance && std::abs(length2 - 1.0f) > tolerance) {
1852
1853 double length = std::sqrt(length2);
1854 vector /= length;
1855 }
1856}
1857
1858inline void ParallelCyclotronTracker::rotateAroundZ(ParticleAttrib<Vector_t>& particleVectors, double const phi) {
1859 // Clockwise rotation of particles 'particleVectors' by 'phi' around Z axis
1860
1861 Tenzor<double, 3> const rotation( std::cos(phi), std::sin(phi), 0,
1862 -std::sin(phi), std::cos(phi), 0,
1863 0, 0, 1);
1864
1865 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
1866 particleVectors[i] = dot(rotation, particleVectors[i]);
1867 }
1868}
1869
1870inline void ParallelCyclotronTracker::rotateAroundZ(Vector_t& myVector, double const phi) {
1871 // Clockwise rotation of single vector 'myVector' by 'phi' around Z axis
1872
1873 Tenzor<double, 3> const rotation( std::cos(phi), std::sin(phi), 0,
1874 -std::sin(phi), std::cos(phi), 0,
1875 0, 0, 1);
1876
1877 myVector = dot(rotation, myVector);
1878}
1879
1880inline void ParallelCyclotronTracker::rotateAroundX(ParticleAttrib<Vector_t>& particleVectors, double const psi) {
1881 // Clockwise rotation of particles 'particleVectors' by 'psi' around X axis
1882
1883 Tenzor<double, 3> const rotation(1, 0, 0,
1884 0, std::cos(psi), std::sin(psi),
1885 0, -std::sin(psi), std::cos(psi));
1886
1887 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
1888
1889 particleVectors[i] = dot(rotation, particleVectors[i]);
1890 }
1891}
1892
1893inline void ParallelCyclotronTracker::rotateAroundX(Vector_t& myVector, double const psi) {
1894 // Clockwise rotation of single vector 'myVector' by 'psi' around X axis
1895
1896 Tenzor<double, 3> const rotation(1, 0, 0,
1897 0, std::cos(psi), std::sin(psi),
1898 0, -std::sin(psi), std::cos(psi));
1899
1900 myVector = dot(rotation, myVector);
1901}
1902
1904 // four vector (w,x,y,z) of the quaternion of P_mean with the positive x-axis
1905
1906 normalizeVector(u);
1907 normalizeVector(v);
1908
1909 double k_cos_theta = dot(u, v);
1910 double k = std::sqrt(dot(u, u) * dot(v, v));
1911 double tolerance1 = 1.0e-5;
1912 double tolerance2 = 1.0e-8;
1913 Vector_t resultVectorComponent;
1914
1915 if (std::abs(k_cos_theta / k + 1.0) < tolerance1) {
1916 // u and v are almost exactly antiparallel so we need to do
1917 // 180 degree rotation around any vector orthogonal to u
1918
1919 resultVectorComponent = cross(u, xaxis);
1920
1921 // If by chance u is parallel to xaxis, use zaxis instead
1922 if (dot(resultVectorComponent, resultVectorComponent) < tolerance2) {
1923 resultVectorComponent = cross(u, zaxis);
1924 }
1925
1926 double halfAngle = 0.5 * Physics::pi;
1927 double sinHalfAngle = std::sin(halfAngle);
1928
1929 resultVectorComponent *= sinHalfAngle;
1930
1931 k = 0.0;
1932 k_cos_theta = std::cos(halfAngle);
1933
1934 } else {
1935 resultVectorComponent = cross(u, v);
1936 }
1937
1938 quaternion(0) = k_cos_theta + k;
1939 quaternion(1) = resultVectorComponent(0);
1940 quaternion(2) = resultVectorComponent(1);
1941 quaternion(3) = resultVectorComponent(2);
1942
1943 normalizeQuaternion(quaternion);
1944}
1945
1946
1948 /* h [ns]
1949 * dt1 [ns]
1950 * dt2 [ns]
1951 */
1953
1954 h *= Units::ns2s;
1955
1956 bool flagNeedUpdate = false;
1957 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
1958 Vector_t const oldR = itsBunch_m->R[i];
1959 double const gamma = std::sqrt(1.0 + dot(itsBunch_m->P[i], itsBunch_m->P[i]));
1960 double const c_gamma = Physics::c / gamma;
1961 Vector_t const v = itsBunch_m->P[i] * c_gamma;
1962 itsBunch_m->R[i] += h * v;
1963 for (const auto & ccd : cavCrossDatas_m) {
1964 double const distNew = (itsBunch_m->R[i][0] * ccd.sinAzimuth - itsBunch_m->R[i][1] * ccd.cosAzimuth) - ccd.perpenDistance;
1965 bool tagCrossing = false;
1966 double distOld;
1967 if (distNew <= 0.0) {
1968 distOld = (oldR[0] * ccd.sinAzimuth - oldR[1] * ccd.cosAzimuth) - ccd.perpenDistance;
1969 if (distOld > 0.0) tagCrossing = true;
1970 }
1971 if (tagCrossing) {
1972 double const dt1 = distOld / std::sqrt(dot(v, v));
1973 double const dt2 = h - dt1;
1974
1975 // Retrack particle from the old postion to cavity gap point
1976 itsBunch_m->R[i] = oldR + dt1 * v;
1977
1978 // Momentum kick
1979 //itsBunch_m->R[i] *= 1000.0; // RFkick uses [itsBunch_m->R] == mm
1980 RFkick(ccd.cavity, itsBunch_m->getT() * Units::s2ns, dt1, i);
1981 //itsBunch_m->R[i] *= 0.001;
1982
1983 itsBunch_m->R[i] += dt2 * itsBunch_m->P[i] * c_gamma;
1984 }
1985 }
1986 flagNeedUpdate |= (itsBunch_m->Bin[i] < 0);
1987 }
1988
1990 return flagNeedUpdate;
1991}
1992
1993
1996
1997 bool flagNeedUpdate = false;
1998 BorisPusher pusher;
1999 double const q = itsBunch_m->Q[0] / Physics::q_e; // For now all particles have the same charge
2000 double const M = itsBunch_m->M[0] * Units::GeV2eV; // For now all particles have the same rest energy
2001
2002 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
2003
2004 pusher.kick(itsBunch_m->R[i], itsBunch_m->P[i],
2005 itsBunch_m->Ef[i], itsBunch_m->Bf[i],
2006 h * Units::ns2s, M, q);
2007 flagNeedUpdate |= (itsBunch_m->Bin[i] < 0);
2008 }
2010 return flagNeedUpdate;
2011}
2012
2013
2015 // h in [ns]
2016
2017 // push particles for first half step
2018 bool flagNeedUpdate = push(0.5 * h);
2019
2020 // Evaluate external fields
2022 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); ++i) {
2023
2024 itsBunch_m->Ef[i] = Vector_t(0.0, 0.0, 0.0);
2025 itsBunch_m->Bf[i] = Vector_t(0.0, 0.0, 0.0);
2026
2028 itsBunch_m->Ef[i], itsBunch_m->Bf[i]);
2029 }
2031
2032 // Kick particles for full step
2033 flagNeedUpdate |= kick(h);
2034
2035 // push particles for second half step
2036 flagNeedUpdate |= push(0.5 * h);
2037
2038 // apply the plugin elements: probe, collimator, stripper, septum
2039 flagNeedUpdate |= applyPluginElements(h);
2040 // destroy particles if they are marked as Bin=-1 in the plugin elements or out of global apeture
2041 deleteParticle(flagNeedUpdate);
2042}
2043
2044
2047
2048 for (beamline_list::iterator sindex = ++(FieldDimensions.begin()); sindex != FieldDimensions.end(); ++sindex) {
2049 if (((*sindex)->first) == ElementType::VACUUM) {
2050 Vacuum* vac = static_cast<Vacuum*>(((*sindex)->second).second);
2052 }
2053 }
2054
2055 bool flag = false;
2056 for (PluginElement* element : pluginElements_m) {
2057 bool tmp = element->check(itsBunch_m,
2059 itsBunch_m->getT(),
2060 dt);
2061 flag |= tmp;
2062
2063 if ( tmp ) {
2065 *gmsg << "* Total number of particles after PluginElement= "
2066 << itsBunch_m->getTotalNum() << endl;
2067 }
2068 }
2069
2071 return flag;
2072}
2073
2076 // Update immediately if any particles are lost during this step
2077
2078 if ((step_m + 1) % Options::delPartFreq != 0) {
2080 return false;
2081 }
2082
2083 allreduce(flagNeedUpdate, 1, std::logical_or<bool>());
2084
2085 if (flagNeedUpdate) {
2086 short bunchCount = itsBunch_m->getNumBunch();
2087 std::vector<size_t> locLostParticleNum(bunchCount, 0);
2088
2089 const int leb = itsBunch_m->getLastemittedBin();
2090 std::unique_ptr<size_t[]> localBinCount;
2091
2092 if ( isMultiBunch() ) {
2093 localBinCount = std::unique_ptr<size_t[]>(new size_t[leb]);
2094 for (int i = 0; i < leb; ++i)
2095 localBinCount[i] = 0;
2096 }
2097
2098 for (unsigned int i = 0; i < itsBunch_m->getLocalNum(); i++) {
2099 if (itsBunch_m->Bin[i] < 0) {
2100 ++locLostParticleNum[itsBunch_m->bunchNum[i]];
2101 itsBunch_m->destroy(1, i);
2102 } else if ( isMultiBunch() ) {
2103 /* we need to count the local number of particles
2104 * per energy bin
2105 */
2106 ++localBinCount[itsBunch_m->Bin[i]];
2107 }
2108 }
2109
2110 if ( isMultiBunch() ) {
2111 // set the local bin count
2112 for (int i = 0; i < leb; ++i) {
2113 itsBunch_m->setLocalBinCount(localBinCount[i], i);
2114 }
2115 }
2116
2117 std::vector<size_t> localnum(bunchCount + 1);
2118 for (size_t i = 0; i < localnum.size() - 1; ++i) {
2119 localnum[i] = itsBunch_m->getLocalNumPerBunch(i) - locLostParticleNum[i];
2120 itsBunch_m->setLocalNumPerBunch(localnum[i], i);
2121 }
2122
2123 /* We need to destroy the particles now
2124 * before we compute the means. We also
2125 * have to update the total number of particles
2126 * otherwise the statistics are wrong.
2127 */
2129
2130 /* total number of particles of individual bunches
2131 * last index of vector contains total number over all
2132 * bunches, used as a check
2133 */
2134 std::vector<size_t> totalnum(bunchCount + 1);
2135 localnum[bunchCount] = itsBunch_m->getLocalNum();
2136
2137 allreduce(localnum.data(), totalnum.data(), localnum.size(), std::plus<size_t>());
2138 itsBunch_m->setTotalNum(totalnum[bunchCount]);
2139
2140 for (short i = 0; i < bunchCount; ++i) {
2141 itsBunch_m->setTotalNumPerBunch(totalnum[i], i);
2142 }
2143
2144 size_t sum = std::accumulate(totalnum.begin(),
2145 totalnum.end() - 1, 0);
2146
2147 if ( sum != totalnum[bunchCount] ) {
2148 throw OpalException("ParallelCyclotronTracker::deleteParticle()",
2149 "Total number of particles " + std::to_string(totalnum[bunchCount]) +
2150 " != " + std::to_string(sum) + " (sum over all bunches)");
2151 }
2152
2153 size_t globLostParticleNum = 0;
2154 size_t locNumLost = std::accumulate(locLostParticleNum.begin(),
2155 locLostParticleNum.end(), 0);
2156
2157 reduce(locNumLost, globLostParticleNum, 1, std::plus<size_t>());
2158
2159 if ( globLostParticleNum > 0 ) {
2160 *gmsg << level3 << "At step " << step_m
2161 << ", lost " << globLostParticleNum << " particles" << endl;
2162 }
2163
2164 if (totalnum[bunchCount] == 0) {
2166 return flagNeedUpdate;
2167 }
2168
2169 Vector_t const meanR = calcMeanR();
2170 Vector_t const meanP = calcMeanP();
2171
2172 // Bunch (local) azimuth at meanR w.r.t. y-axis
2173 double const phi = calculateAngle(meanP(0), meanP(1)) - 0.5 * Physics::pi;
2174
2175 // Bunch (local) elevation at meanR w.r.t. xy plane
2176 double const psi = 0.5 * Physics::pi - std::acos(meanP(2) / std::sqrt(dot(meanP, meanP)));
2177
2178 // For statistics data, the bunch is transformed into a local coordinate system
2179 // with meanP in direction of y-axis -DW
2180 globalToLocal(itsBunch_m->R, phi, psi, meanR);
2181 globalToLocal(itsBunch_m->P, phi, psi, Vector_t(0.0)); // P should be rotated, but not shifted
2182
2183 //itsBunch_m->R *= Vector_t(0.001); // mm --> m
2184
2185 // Now destroy particles and update pertinent parameters in local frame
2186 // Note that update() will be called within boundp() -DW
2187 itsBunch_m->boundp();
2188 //itsBunch_m->update();
2189
2191
2192 //itsBunch_m->R *= Vector_t(1000.0); // m --> mm
2193
2194 localToGlobal(itsBunch_m->R, phi, psi, meanR);
2195 localToGlobal(itsBunch_m->P, phi, psi, Vector_t(0.0));
2196
2197 // If particles were deleted, recalculate bingamma and reset BinID for remaining particles
2198 if ( isMultiBunch() )
2199 mbHandler_m->updateParticleBins(itsBunch_m);
2200 }
2201
2203 return flagNeedUpdate;
2204}
2205
2207
2208 std::string f = OpalData::getInstance()->getInputBasename() + std::string("-trackOrbit.dat");
2209
2210 outfTrackOrbit_m.setf(std::ios::scientific, std::ios::floatfield);
2211 outfTrackOrbit_m.precision(8);
2212
2213 if (myNode_m == 0) {
2214
2215 if (OpalData::getInstance()->inRestartRun()) {
2216 outfTrackOrbit_m.open(f.c_str(), std::ios::app);
2217 outfTrackOrbit_m << "# Restart at integration step " << itsBunch_m->getLocalTrackStep() << std::endl;
2218 } else {
2219 outfTrackOrbit_m.open(f.c_str());
2220 outfTrackOrbit_m << "# The six-dimensional phase space data in the global Cartesian coordinates" << std::endl;
2221 outfTrackOrbit_m << "# Part. ID x [m] beta_x*gamma y [m] beta_y*gamma z [m] beta_z*gamma" << std::endl;
2222 }
2223 }
2224}
2225
2227
2228 if (!OpalData::getInstance()->inRestartRun()) {
2229 // Start a new run (no restart)
2230
2231 double const initialReferenceTheta = referenceTheta * Units::deg2rad;
2232
2233 // TODO: Replace with TracerParticle
2234 // Force the initial phase space values of the particle with ID = 0 to zero,
2235 // to set it as a reference particle.
2236 if (initialTotalNum_m > 2) {
2237 for (size_t i = 0; i < initialLocalNum_m; ++i) {
2238 if (itsBunch_m->ID[i] == 0) {
2239 itsBunch_m->R[i] = Vector_t(0.0);
2240 itsBunch_m->P[i] = Vector_t(0.0);
2241 }
2242 }
2243 }
2244
2245 // Initialize global R
2246 //itsBunch_m->R *= Vector_t(1000.0); // m --> mm
2247
2251
2252 localToGlobal(itsBunch_m->R, initialReferenceTheta, initMeanR);
2253
2254 // Initialize global P (Cartesian, but input P_ref is in Pr, Ptheta, Pz,
2255 // so translation has to be done before the rotation this once)
2256 // Cave: In the local frame, the positive y-axis is the direction of movement -DW
2257 for (size_t i = 0; i < initialLocalNum_m; ++i) {
2258 itsBunch_m->P[i](0) += referencePr;
2259 itsBunch_m->P[i](1) += referencePt;
2260 itsBunch_m->P[i](2) += referencePz;
2261 }
2262
2263 // Out of the three coordinates of meanR (R, Theta, Z) only the angle
2264 // changes the momentum vector...
2265 localToGlobal(itsBunch_m->P, initialReferenceTheta);
2266
2268 if (distType == DistributionType::FROMFILE) {
2270 }
2271
2272 // Initialize the bin number of the first bunch to 0
2273 for (size_t i = 0; i < initialLocalNum_m; ++i) {
2274 itsBunch_m->Bin[i] = 0;
2275 }
2276
2277 // Set time step per particle
2278 setTimeStep();
2279
2280 // Backup initial distribution if multi bunch mode
2281 if ((initialTotalNum_m > 2) && isMultiBunch() && mbHandler_m->isForceMode()) {
2282 mbHandler_m->saveBunch(itsBunch_m);
2283 }
2284
2285 // Else: Restart from the distribution in the h5 file
2286 } else {
2287
2288 // Do a local frame restart (we have already checked that the old h5 file was saved in local
2289 // frame as well).
2291
2292 *gmsg << "* Restart in the local frame" << endl;
2293 //itsBunch_m->R *= Vector_t(1000.0); // m --> mm
2294
2298
2299 // Do the transformations
2302
2303 // Initialize the bin number of the first bunch to 0
2304 for (size_t i = 0; i < initialLocalNum_m; ++i) {
2305 itsBunch_m->Bin[i] = 0;
2306 }
2307
2308 // Or do a global frame restart (no transformations necessary)
2309 } else {
2310 *gmsg << "* Restart in the global frame" << endl;
2311
2313 //itsBunch_m->R *= Vector_t(1000.0); // m --> mm
2314 }
2315 }
2316
2317 // set the number of particles per bunch
2319
2320 // ------------ Get some Values ---------------------------------------------------------- //
2321 Vector_t const meanR = calcMeanR();
2322 Vector_t const meanP = calcMeanP();
2323
2324 // Bunch (local) azimuth at meanR w.r.t. y-axis
2325 double const phi = calculateAngle(meanP(0), meanP(1)) - 0.5 * Physics::pi;
2326
2327 // Bunch (local) elevation at meanR w.r.t. xy plane
2328 double const psi = 0.5 * Physics::pi - std::acos(meanP(2) / std::sqrt(dot(meanP, meanP)));
2329
2330 double radius = std::sqrt(meanR[0] * meanR[0] + meanR[1] * meanR[1]); // [m]
2331
2332 if ( isMultiBunch() ) {
2333 mbHandler_m->setRadiusTurns(radius);
2334 }
2335
2336 // Do boundp and repartition in the local frame at beginning of this run
2337 globalToLocal(itsBunch_m->R, phi, psi, meanR);
2338 globalToLocal(itsBunch_m->P, phi, psi); // P should be rotated, but not shifted
2339
2340 //itsBunch_m->R *= Vector_t(0.001); // mm --> m
2341
2342 itsBunch_m->boundp();
2343
2344 checkNumPart(std::string("\n* Before repartition: "));
2345 repartition();
2346 checkNumPart(std::string("\n* After repartition: "));
2347
2349
2350 *gmsg << endl << "* *********************** Bunch information in local frame: ************************";
2351 *gmsg << *itsBunch_m << endl;
2352
2353 //itsBunch_m->R *= Vector_t(1000.0); // m --> mm
2354
2355 localToGlobal(itsBunch_m->R, phi, psi, meanR);
2356 localToGlobal(itsBunch_m->P, phi, psi);
2357
2358 // Save initial distribution if not a restart
2359 if (!OpalData::getInstance()->inRestartRun()) {
2360 step_m -= 1;
2361
2364
2365 step_m += 1;
2366 }
2367
2368 // Print out the Bunch information at beginning of the run. Because the bunch information
2369 // displays in units of m we have to change back and forth one more time.
2370 //itsBunch_m->R *= Vector_t(0.001); // mm --> m
2371
2373
2374 // multi-bunch simulation only
2376
2377 *gmsg << endl << "* *********************** Bunch information in global frame: ***********************";
2378 *gmsg << *itsBunch_m << endl;
2379
2380 //itsBunch_m->R *= Vector_t(1000.0); // m --> mm
2381}
2382
2384 for (size_t i = 0; i < initialLocalNum_m; ++i) {
2385 itsBunch_m->dt[i] = itsBunch_m->getdT();
2386 }
2387}
2388
2390
2391 double pTotalMean = 0.0;
2392 for (size_t i = 0; i < initialLocalNum_m; ++i) {
2393 pTotalMean += euclidean_norm(itsBunch_m->P[i]);
2394 }
2395
2396 allreduce(pTotalMean, 1, std::plus<double>());
2397
2398 pTotalMean /= initialTotalNum_m;
2399
2400 if (std::abs(pTotalMean - referencePtot) / pTotalMean > 1e-2) {
2401 throw OpalException("ParallelCyclotronTracker::checkFileMomentum",
2402 "The total momentum of the particle distribution\n"
2403 "in the global reference frame: " +
2404 std::to_string(pTotalMean) + ",\n"
2405 "is different from the momentum given\n"
2406 "in the \"BEAM\" command: " +
2407 std::to_string(referencePtot) + ".\n"
2408 "In Opal-cycl the initial distribution\n"
2409 "is specified in the local reference frame.\n"
2410 "When using a \"FROMFILE\" type distribution, the momentum \n"
2411 "must be the same as the specified in the \"BEAM\" command,\n"
2412 "which is in global reference frame.");
2413 }
2414}
2415
2416
2417//TODO: This can be completely rewritten with TracerParticle -DW
2420
2421 if (Ippl::getNodes() > 1 ) {
2422
2423 double x;
2424 int id;
2425 dvector_t tmpr;
2426 ivector_t tmpi;
2427
2429
2430 // for all nodes, find the location of particle with ID = 0 & 1 in bunch containers
2431 int found[2] = {-1, -1};
2432 int counter = 0;
2433
2434 for (size_t i = 0; i < itsBunch_m->getLocalNum(); ++i) {
2435 if (itsBunch_m->ID[i] == 0) {
2436 found[counter] = i;
2437 counter++;
2438 }
2439 if (itsBunch_m->ID[i] == 1) {
2440 found[counter] = i;
2441 counter++;
2442 }
2443 }
2444
2445 if (myNode_m == 0) {
2446 int notReceived = Ippl::getNodes() - 1;
2447 int numberOfPart = 0;
2448 // receive other nodes
2449 while(notReceived > 0) {
2450
2451 int node = COMM_ANY_NODE;
2452 Message *rmsg = Ippl::Comm->receive_block(node, tag);
2453
2454 if (rmsg == nullptr)
2455 ERRORMSG("Could not receive from client nodes in main." << endl);
2456
2457 --notReceived;
2458
2459 rmsg->get(&numberOfPart);
2460
2461 for (int i = 0; i < numberOfPart; ++i) {
2462 rmsg->get(&id);
2463 tmpi.push_back(id);
2464 for (int ii = 0; ii < 6; ii++) {
2465 rmsg->get(&x);
2466 tmpr.push_back(x);
2467 }
2468 }
2469 delete rmsg;
2470 }
2471 // own node
2472 for (int i = 0; i < counter; ++i) {
2473
2474 tmpi.push_back(itsBunch_m->ID[found[i]]);
2475
2476 for (int j = 0; j < 3; ++j) {
2477 tmpr.push_back(itsBunch_m->R[found[i]](j));
2478 tmpr.push_back(itsBunch_m->P[found[i]](j));
2479 }
2480 }
2481 // store
2482 dvector_t::iterator itParameter = tmpr.begin();
2483
2484 for (auto tmpid : tmpi) {
2485
2486 outfTrackOrbit_m << "ID" << tmpid;
2487
2488 if (tmpid == 0) { // for stat file
2489 itsBunch_m->RefPartR_m[0] = *itParameter;
2490 itsBunch_m->RefPartR_m[1] = *(itParameter + 2);
2491 itsBunch_m->RefPartR_m[2] = *(itParameter + 4);
2492 itsBunch_m->RefPartP_m[0] = *(itParameter + 1);
2493 itsBunch_m->RefPartP_m[1] = *(itParameter + 3);
2494 itsBunch_m->RefPartP_m[2] = *(itParameter + 5);
2495 }
2496 for (int ii = 0; ii < 6; ii++) {
2497 outfTrackOrbit_m << " " << *itParameter;
2498 ++itParameter;
2499 }
2500
2502 }
2503 } else {
2504 // for other nodes
2505 Message *smsg = new Message();
2506 smsg->put(counter);
2507
2508 for (int i = 0; i < counter; i++) {
2509
2510 smsg->put(itsBunch_m->ID[found[i]]);
2511
2512 for (int j = 0; j < 3; j++) {
2513 smsg->put(itsBunch_m->R[found[i]](j));
2514 smsg->put(itsBunch_m->P[found[i]](j));
2515 }
2516 }
2517
2518 if (!Ippl::Comm->send(smsg, 0, tag)) {
2519 ERRORMSG("Ippl::Comm->send(smsg, 0, tag) failed " << endl);
2520 }
2521 }
2522
2524
2525 } else {
2526
2527 for (size_t i = 0; i < itsBunch_m->getLocalNum(); i++) {
2528 if (itsBunch_m->ID[i] == 0 || itsBunch_m->ID[i] == 1) {
2529
2530 outfTrackOrbit_m << "ID" << itsBunch_m->ID[i] << " ";
2531 outfTrackOrbit_m << itsBunch_m->R[i](0) << " " << itsBunch_m->P[i](0) << " ";
2532 outfTrackOrbit_m << itsBunch_m->R[i](1) << " " << itsBunch_m->P[i](1) << " ";
2533 outfTrackOrbit_m << itsBunch_m->R[i](2) << " " << itsBunch_m->P[i](2) << std::endl;
2534
2535 if (itsBunch_m->ID[i] == 0) { // for stat file
2538 }
2539 }
2540 }
2541 }
2542
2544}
2545
2547
2549
2550 // dump stat file per bunch in case of multi-bunch mode
2551 if (isMultiBunch()) {
2552 double phi = 0.0, psi = 0.0;
2553 Vector_t meanR = calcMeanR();
2554
2555 // Bunch (global) angle w.r.t. x-axis (cylinder coordinates)
2556 double theta = calculateAngle(meanR(0), meanR(1)) * Units::rad2deg;
2557
2558 // fix azimuth_m --> make monotonically increasing
2560
2562
2564 Vector_t meanP = calcMeanP();
2565
2566 // Bunch (local) azimuth at meanR w.r.t. y-axis
2567 phi = calculateAngle(meanP(0), meanP(1)) - 0.5 * Physics::pi;
2568
2569 // Bunch (local) elevation at meanR w.r.t. xy plane
2570 psi = 0.5 * Physics::pi - std::acos(meanP(2) / std::sqrt(dot(meanP, meanP)));
2571
2572 // Rotate so Pmean is in positive y direction. No shift, so that normalized emittance and
2573 // unnormalized emittance as well as centroids are calculated correctly in
2574 // PartBunch::calcBeamParameters()
2575 globalToLocal(itsBunch_m->R, phi, psi, meanR);
2576 globalToLocal(itsBunch_m->P, phi, psi);
2577 }
2578
2580
2582 localToGlobal(itsBunch_m->R, phi, psi, meanR);
2583 localToGlobal(itsBunch_m->P, phi, psi);
2584 }
2585
2587 return;
2588 }
2589
2590 // --------------------------------- Get some Values ---------------------------------------- //
2591 double const temp_t = itsBunch_m->getT() * Units::s2ns;
2592 Vector_t meanR;
2593 Vector_t meanP;
2595 meanR = calcMeanR();
2596 meanP = calcMeanP();
2597 } else if (itsBunch_m->getLocalNum() > 0) {
2598 meanR = itsBunch_m->R[0];
2599 meanP = itsBunch_m->P[0];
2600 }
2601 double phi = 0;
2602 double psi = 0;
2603
2604 // Bunch (global) angle w.r.t. x-axis (cylinder coordinates)
2605 double azimuth = calculateAngle(meanR(0), meanR(1)) * Units::rad2deg;
2606
2607 // fix azimuth_m --> make monotonically increasing
2609
2610 // -------------- Calculate the external fields at the center of the bunch ----------------- //
2611 beamline_list::iterator DumpSindex = FieldDimensions.begin();
2612
2613 extE_m = Vector_t(0.0, 0.0, 0.0);
2614 extB_m = Vector_t(0.0, 0.0, 0.0);
2615
2616 (((*DumpSindex)->second).second)->apply(meanR, meanP, temp_t, extE_m, extB_m);
2617
2618 // If we are saving in local frame, bunch and fields at the bunch center have to be rotated
2619 // TODO: Make decision if we maybe want to always save statistics data in local frame? -DW
2621 // -------------------- ----------- Do Transformations ---------------------------------- //
2622 // Bunch (local) azimuth at meanR w.r.t. y-axis
2623 phi = calculateAngle(meanP(0), meanP(1)) - 0.5 * Physics::pi;
2624
2625 // Bunch (local) elevation at meanR w.r.t. xy plane
2626 psi = 0.5 * Physics::pi - std::acos(meanP(2) / std::sqrt(dot(meanP, meanP)));
2627
2628 // Rotate so Pmean is in positive y direction. No shift, so that normalized emittance and
2629 // unnormalized emittance as well as centroids are calculated correctly in
2630 // PartBunch::calcBeamParameters()
2631 globalToLocal(extB_m, phi, psi);
2632 globalToLocal(extE_m, phi, psi);
2633 globalToLocal(itsBunch_m->R, phi, psi);
2634 globalToLocal(itsBunch_m->P, phi, psi);
2635 }
2636
2637 //itsBunch_m->R *= Vector_t(0.001); // mm -> m
2638
2639 FDext_m[0] = extB_m * Units::kG2T;
2640 FDext_m[1] = extE_m; // kV/mm? -DW
2641
2642 // Save the stat file
2644
2645 //itsBunch_m->R *= Vector_t(1000.0); // m -> mm
2646
2647 // If we are in local mode, transform back after saving
2649 localToGlobal(itsBunch_m->R, phi, psi);
2650 localToGlobal(itsBunch_m->P, phi, psi);
2651 }
2652
2654}
2655
2656
2658 // --------------------------------- Particle dumping --------------------------------------- //
2659 // Note: Don't dump when
2660 // 1. after one turn
2661 // in order to sychronize the dump step for multi-bunch and single bunch for comparison
2662 // with each other during post-process phase.
2663 // ------------------------------------------------------------------------------------------ //
2665
2666 // --------------------------------- Get some Values ---------------------------------------- //
2667 double const temp_t = itsBunch_m->getT() * Units::s2ns;
2668
2669 Vector_t meanR;
2670 Vector_t meanP;
2671
2672 // in case of multi-bunch mode take always bunch mean (although it takes all bunches)
2674 meanR = calcMeanR();
2675 meanP = calcMeanP();
2676 } else if (itsBunch_m->getLocalNum() > 0) {
2677 meanR = itsBunch_m->R[0];
2678 meanP = itsBunch_m->P[0];
2679 }
2680
2681 double const betagamma_temp = std::sqrt(dot(meanP, meanP));
2682
2683 double const E = itsBunch_m->get_meanKineticEnergy();
2684
2685 // Bunch (global) angle w.r.t. x-axis (cylinder coordinates)
2686 double const theta = std::atan2(meanR(1), meanR(0));
2687
2688 // Bunch (local) azimuth at meanR w.r.t. y-axis
2689 double const phi = calculateAngle(meanP(0), meanP(1)) - 0.5 * Physics::pi;
2690
2691 // Bunch (local) elevation at meanR w.r.t. xy plane
2692 double const psi = 0.5 * Physics::pi - std::acos(meanP(2) / std::sqrt(dot(meanP, meanP)));
2693
2694 // ---------------- Re-calculate reference values in format of input values ----------------- //
2695 // Position:
2696 // New OPAL 2.0: Init in m (back to mm just for output) -DW
2697 referenceR = computeRadius(meanR); // includes m --> mm conversion
2699 referenceZ = Units::m2mm * meanR(2);
2700
2701 // Momentum in Theta-hat, R-hat, Z-hat coordinates at position meanR:
2702 referencePtot = betagamma_temp;
2703 referencePz = meanP(2);
2704 referencePr = meanP(0) * std::cos(theta) + meanP(1) * std::sin(theta);
2707
2708 // ----- Calculate the external fields at the center of the bunch (Cave: Global Frame) ----- //
2709 beamline_list::iterator DumpSindex = FieldDimensions.begin();
2710
2711 extE_m = Vector_t(0.0, 0.0, 0.0);
2712 extB_m = Vector_t(0.0, 0.0, 0.0);
2713
2714 (((*DumpSindex)->second).second)->apply(meanR, meanP, temp_t, extE_m, extB_m);
2715 FDext_m[0] = extB_m * Units::kG2T;
2716 FDext_m[1] = extE_m;
2717
2718 // -------------- If flag psDumpFrame is global, dump bunch in global frame ------------- //
2720 bool dumpLocalFrame = true;
2721 std::string dumpString = "local";
2723 dumpLocalFrame = false;
2724 dumpString = "global";
2725 }
2726
2727 if (dumpLocalFrame == true) {
2728 // ---------------- If flag psDumpFrame is local, dump bunch in local frame ---------------- //
2729
2730 // The bunch is transformed into a local coordinate system with meanP in direction of y-axis
2731 globalToLocal(itsBunch_m->R, phi, psi, meanR);
2732 //globalToLocal(itsBunch_m->R, phi, psi, meanR * Vector_t(0.001));
2733 globalToLocal(itsBunch_m->P, phi, psi); // P should only be rotated
2734
2735 globalToLocal(extB_m, phi, psi);
2736 globalToLocal(extE_m, phi, psi);
2737 }
2738
2739 FDext_m[0] = extB_m * Units::kG2T;
2740 FDext_m[1] = extE_m;
2741
2742 lastDumpedStep_m = itsDataSink->dumpH5(itsBunch_m, // Local and in m
2743 FDext_m, E,
2747 referenceR,
2749 referenceZ,
2750 phi / Units::deg2rad, // P_mean azimuth
2751 // at ref. R/Th/Z
2752 psi / Units::deg2rad, // P_mean elevation
2753 // at ref. R/Th/Z
2754 dumpLocalFrame); // Flag localFrame
2755
2756 if (dumpLocalFrame == true) {
2757 // Return to global frame
2758 localToGlobal(itsBunch_m->R, phi, psi, meanR);
2759 //localToGlobal(itsBunch_m->R, phi, psi, meanR * Vector_t(0.001));
2760 localToGlobal(itsBunch_m->P, phi, psi);
2761 }
2762
2763 // Tell user in which mode we are dumping
2764 // New: no longer dumping for num part < 3, omit phase space dump number info
2765 if (lastDumpedStep_m == -1) {
2766 *gmsg << endl << "* Integration step " << step_m + 1
2767 << " (no phase space dump for <= 2 particles)" << endl;
2768 } else {
2769 *gmsg << endl << "* Phase space dump " << lastDumpedStep_m
2770 << " (" << dumpString << " frame) at integration step " << step_m + 1 << endl;
2771 }
2772 }
2773
2774 // Print dump information on screen
2775 *gmsg << "* T = " << temp_t << " ns"
2776 << ", Live Particles: " << itsBunch_m->getTotalNum() << endl;
2777 *gmsg << "* E = " << E << " MeV"
2778 << ", beta * gamma = " << betagamma_temp << endl;
2779 *gmsg << "* Bunch position: R = " << referenceR << " mm"
2780 << ", Theta = " << referenceTheta << " Deg"
2781 << ", Z = " << referenceZ << " mm" << endl;
2782 *gmsg << "* Local Azimuth = " << phi / Units::deg2rad << " Deg"
2783 << ", Local Elevation = " << psi / Units::deg2rad << " Deg" << endl;
2784
2786}
2787
2789 return (step_m > 10) && (((step_m + 1) %setup_m.stepsPerTurn) == 0);
2790}
2791
2792void ParallelCyclotronTracker::update_m(double& t, const double& dt,
2793 const bool& finishedTurn)
2794{
2795 // Reference parameters
2796 t += dt;
2797
2798 updateTime(dt);
2799
2801 if (!(step_m + 1 % 1000)) {
2802 *gmsg << "Step " << step_m + 1 << endl;
2803 }
2804
2805 updatePathLength(dt);
2806
2807 // Here is global frame, don't do itsBunch_m->boundp();
2808
2809 if (itsBunch_m->getTotalNum()>0) {
2810 // Only dump last step if we have particles left.
2811 // Check separately for phase space (ps) and statistics (stat) data dump frequency
2812 if ( mode_m != TrackingMode::SEO && ( ((step_m + 1) % Options::psDumpFreq == 0) ||
2813 (Options::psDumpEachTurn && finishedTurn)))
2814 {
2815 // Write phase space data to h5 file
2817 }
2818
2819 if ( mode_m != TrackingMode::SEO && ( ((step_m + 1) % Options::statDumpFreq == 0) ||
2820 (Options::psDumpEachTurn && finishedTurn)))
2821 {
2822 // Write statistics data to stat file
2824 }
2825 }
2826
2827 if (Options::psDumpEachTurn && finishedTurn) {
2828 for (PluginElement* element : pluginElements_m) {
2829 element->save();
2830 }
2831 }
2832}
2833
2834
2835std::tuple<double, double, double> ParallelCyclotronTracker::initializeTracking_m() {
2836 // Read in some control parameters
2839
2840 // Define 3 special azimuthal angles where dump particle's six parameters
2841 // at each turn into 3 ASCII files. only important in single particle tracking
2842 azimuth_angle_m.resize(3);
2843 azimuth_angle_m[0] = 0.0;
2844 azimuth_angle_m[1] = 22.5 * Units::deg2rad;
2845 azimuth_angle_m[2] = 45.0 * Units::deg2rad;
2846
2847 double harm = getHarmonicNumber();
2848 double dt = itsBunch_m->getdT() * Units::s2ns * harm;
2849 double t = itsBunch_m->getT() * Units::s2ns;
2850
2851 double oldReferenceTheta = referenceTheta * Units::deg2rad; // init here, reset each step
2852 setup_m.deltaTheta = Physics::pi / (setup_m.stepsPerTurn); // half of the average angle per step
2853
2854 //int stepToLastInj = itsBunch_m->getSteptoLastInj(); // TODO: Do we need this? -DW
2855
2856 // Record how many bunches have already been injected. ONLY FOR MBM
2857 if (isMultiBunch()) {
2858 mbHandler_m->setNumBunch(itsBunch_m->getNumBunch());
2859 }
2860
2862
2863 // Get data from h5 file for restart run and reset current step
2864 // to last step of previous simulation
2865 if (OpalData::getInstance()->inRestartRun()) {
2866
2870
2871 *gmsg << "* Restart at integration step " << restartStep0_m
2872 << " at turn " << turnnumber_m - 1 << endl;
2873
2875 }
2876
2877 setup_m.stepsNextCheck = step_m + setup_m.stepsPerTurn; // Steps to next check for transition
2878
2880
2881 if (isMultiBunch()) {
2882 mbHandler_m->updateParticleBins(itsBunch_m);
2883 }
2884
2885 // --- Output to user --- //
2886 *gmsg << "* Beginning of this run is at t = " << t << " [ns]" << endl;
2887 *gmsg << "* The time step is set to dt = " << dt << " [ns]" << endl;
2888
2889 if ( isMultiBunch() ) {
2890 *gmsg << "* MBM: Time interval between neighbour bunches is set to "
2891 << setup_m.stepsPerTurn * dt << "[ns]" << endl;
2892 *gmsg << "* MBM: The particles energy bin reset frequency is set to "
2894 }
2895
2896 *gmsg << "* Single particle trajectory dump frequency is set to " << Options::sptDumpFreq << endl;
2897 *gmsg << "* The frequency to solve space charge fields is set to " << setup_m.scSolveFreq << endl;
2898 *gmsg << "* The repartition frequency is set to " << Options::repartFreq << endl;
2899
2900 switch ( mode_m ) {
2901 case TrackingMode::SEO: {
2902 *gmsg << endl;
2903 *gmsg << "* ------------------------- STATIC EQUILIBRIUM ORBIT MODE ----------------------------- *" << endl
2904 << "* Instruction: When the total particle number is equal to 2, SEO mode is triggered *" << endl
2905 << "* automatically. This mode does NOT include any RF cavities. The initial distribution *" << endl
2906 << "* file must be specified. In the file the first line is for reference particle and the *" << endl
2907 << "* second line is for off-center particle. The tune is calculated by FFT routines based *" << endl
2908 << "* on these two particles. *" << endl
2909 << "* ---------------- NOTE: SEO MODE ONLY WORKS SERIALLY ON SINGLE NODE ------------------ *" << endl;
2910
2911 if (Ippl::getNodes() != 1)
2912 throw OpalException("Error in ParallelCyclotronTracker::initializeTracking_m",
2913 "SEO MODE ONLY WORKS SERIALLY ON SINGLE NODE!");
2914 break;
2915 }
2916 case TrackingMode::SINGLE: {
2917 *gmsg << endl;
2918 *gmsg << "* ------------------------------ SINGLE PARTICLE MODE --------------------------------- *" << endl
2919 << "* Instruction: When the total particle number is equal to 1, single particle mode is *" << endl
2920 << "* triggered automatically. The initial distribution file must be specified which should *" << endl
2921 << "* contain only one line for the single particle *" << endl
2922 << "* ---------NOTE: SINGLE PARTICLE MODE ONLY WORKS SERIALLY ON A SINGLE NODE ------------ *" << endl;
2923
2924 if (Ippl::getNodes() != 1)
2925 throw OpalException("Error in ParallelCyclotronTracker::initializeTracking_m",
2926 "SINGLE PARTICLE MODE ONLY WORKS SERIALLY ON A SINGLE NODE!");
2927
2928 // For single particle mode open output files
2930 break;
2931 }
2932 case TrackingMode::BUNCH: {
2933 break;
2934 }
2936 default: {
2937 throw OpalException("ParallelCyclotronTracker::initializeTracking_m()",
2938 "No such tracking mode.");
2939 }
2940 }
2941
2942 return std::make_tuple(t, dt, oldReferenceTheta);
2943}
2944
2945
2947 dvector_t& Tdeltr,
2948 dvector_t& Tdeltz, ivector_t& TturnNumber) {
2949
2950 for (size_t ii = 0; ii < (itsBunch_m->getLocalNum()); ii++) {
2951 if (itsBunch_m->ID[ii] == 0) {
2952 double FinalMomentum2 = std::pow(itsBunch_m->P[ii](0), 2.0) + std::pow(itsBunch_m->P[ii](1), 2.0) + std::pow(itsBunch_m->P[ii](2), 2.0);
2953 double FinalEnergy = (std::sqrt(1.0 + FinalMomentum2) - 1.0) * itsBunch_m->getM() * Units::eV2MeV;
2954 *gmsg << "* Final energy of reference particle = " << FinalEnergy << " [MeV]" << endl;
2955 *gmsg << "* Total phase space dump number(includes the initial distribution) = " << lastDumpedStep_m + 1 << endl;
2956 *gmsg << "* One can restart simulation from the last dump step (--restart " << lastDumpedStep_m << ")" << endl;
2957 }
2958 }
2959
2961
2962 switch ( mode_m ) {
2963 case TrackingMode::SEO: {
2964 // Calculate tunes after tracking.
2965 *gmsg << endl;
2966 *gmsg << "* **************** The result for tune calculation (NO space charge) ******************* *" << endl
2967 << "* Number of tracked turns: " << TturnNumber.back() << endl;
2968 double nur, nuz;
2969 getTunes(Ttime, Tdeltr, Tdeltz, TturnNumber.back(), nur, nuz);
2970 break;
2971 }
2973 closeFiles();
2974 // fall through
2975 case TrackingMode::BUNCH: // we do nothing
2977 default: {
2978 // not for multibunch
2979 if (!isMultiBunch()) {
2980 *gmsg << "*" << endl;
2981 *gmsg << "* Finished during turn " << turnnumber_m << " (" << turnnumber_m - 1 << " turns completed)" << endl;
2982 *gmsg << "* Cave: Turn number is not correct for restart mode"<< endl;
2983 }
2984 }
2985 }
2986
2988
2989 if (myNode_m == 0) {
2990 outfTrackOrbit_m.close();
2991 }
2992
2993 *gmsg << endl << "* *********************** Bunch information in global frame: ***********************";
2994
2995 if (itsBunch_m->getTotalNum() > 0) {
2996 // Print out the Bunch information at end of the run.
2998 *gmsg << *itsBunch_m << endl;
2999 } else {
3000 *gmsg << endl << "* No Particles left in bunch!" << endl;
3001 *gmsg << "* **********************************************************************************" << endl;
3002 }
3003}
3004
3006 if ( initialTotalNum_m == 1 ) {
3008 } else if ( initialTotalNum_m == 2 ) {
3010 } else if ( initialTotalNum_m > 2 ) {
3012 } else {
3014 }
3015}
3016
3017void ParallelCyclotronTracker::seoMode_m(double& t, const double dt, bool& /*finishedTurn*/,
3018 dvector_t& Ttime, dvector_t& Tdeltr,
3019 dvector_t& Tdeltz, ivector_t& TturnNumber) {
3020
3021 // 2 particles: Trigger SEO mode
3022 // (Switch off cavity and calculate betatron oscillation tuning)
3023 double r_tuning[2], z_tuning[2];
3024
3026 for (size_t i = 0; i < (itsBunch_m->getLocalNum()); i++) {
3027
3028 if ((step_m % Options::sptDumpFreq == 0)) {
3029 outfTrackOrbit_m << "ID" << (itsBunch_m->ID[i]);
3030 outfTrackOrbit_m << " " << itsBunch_m->R[i](0)
3031 << " " << itsBunch_m->P[i](0)
3032 << " " << itsBunch_m->R[i](1)
3033 << " " << itsBunch_m->P[i](1)
3034 << " " << itsBunch_m->R[i](2)
3035 << " " << itsBunch_m->P[i](2)
3036 << std::endl;
3037 }
3038
3039 double OldTheta = calculateAngle(itsBunch_m->R[i](0), itsBunch_m->R[i](1));
3040 r_tuning[i] = itsBunch_m->R[i](0) * std::cos(OldTheta) +
3041 itsBunch_m->R[i](1) * std::sin(OldTheta);
3042
3043 z_tuning[i] = itsBunch_m->R[i](2);
3044
3045 // Integrate for one step in the lab Cartesian frame (absolute value).
3046 itsStepper_mp->advance(itsBunch_m, i, t, dt);
3047
3048 if ( (i == 0) && isTurnDone() ) {
3049 turnnumber_m++;
3050 }
3051
3052 } // end for: finish one step tracking for all particles for initialTotalNum_m == 2 mode
3054
3055 // store dx and dz for future tune calculation if higher precision needed, reduce freqSample.
3056 if (step_m % Options::sptDumpFreq == 0) {
3057 Ttime.push_back(t * Units::ns2s);
3058 Tdeltz.push_back(z_tuning[1]);
3059 Tdeltr.push_back(r_tuning[1] - r_tuning[0]);
3060 TturnNumber.push_back(turnnumber_m);
3061 }
3062}
3063
3064
3065void ParallelCyclotronTracker::singleMode_m(double& t, const double dt,
3066 bool& finishedTurn, double& oldReferenceTheta) {
3067 // 1 particle: Trigger single particle mode
3068
3069 // ********************************************************************************** //
3070 // * This was moved here b/c collision should be tested before the actual * //
3071 // * timestep (bgf_main_collision_test() predicts the next step automatically) * //
3072
3073 // apply the plugin elements: probe, collimator, stripper, septum
3074 bool flagNeedUpdate = applyPluginElements(dt);
3075
3076 // check if we lose particles at the boundary
3078 // ********************************************************************************** //
3079
3080 if (itsBunch_m->getLocalNum() == 0) return; // might happen if particle is in collimator
3081
3083
3084 unsigned int i = 0; // we only have a single particle
3085
3086 if ( step_m % Options::sptDumpFreq == 0 ) {
3087 outfTrackOrbit_m << "ID" <<itsBunch_m->ID[i]
3088 << " " << itsBunch_m->R[i](0)
3089 << " " << itsBunch_m->P[i](0)
3090 << " " << itsBunch_m->R[i](1)
3091 << " " << itsBunch_m->P[i](1)
3092 << " " << itsBunch_m->R[i](2)
3093 << " " << itsBunch_m->P[i](2)
3094 << std::endl;
3095 }
3096
3097 double temp_meanTheta = calculateAngle2(itsBunch_m->R[i](0),
3098 itsBunch_m->R[i](1)); // [-pi, pi]
3099
3101 temp_meanTheta, finishedTurn);
3102
3104 oldReferenceTheta, temp_meanTheta);
3105
3106 oldReferenceTheta = temp_meanTheta;
3107
3108 // used for gap crossing checking
3109 Vector_t Rold = itsBunch_m->R[i]; // [x,y,z] (mm)
3110 Vector_t Pold = itsBunch_m->P[i]; // [px,py,pz] (beta*gamma)
3111
3112 // integrate for one step in the lab Cartesian frame (absolute value).
3113 itsStepper_mp->advance(itsBunch_m, i, t, dt);
3114
3115 // If gap crossing happens, do momenta kicking (not if gap crossing just happened)
3116 if (itsBunch_m->cavityGapCrossed[i] == true)
3117 itsBunch_m->cavityGapCrossed[i] = false;
3118 else
3119 gapCrossKick_m(i, t, dt, Rold, Pold);
3121
3122 // Destroy particles if they are marked as Bin = -1 in the plugin elements
3123 // or out of the global aperture
3124 flagNeedUpdate |= (itsBunch_m->Bin[i] < 0);
3125 deleteParticle(flagNeedUpdate);
3126}
3127
3128
3129void ParallelCyclotronTracker::bunchMode_m(double& t, const double dt, bool& finishedTurn) {
3130
3131 // Flag for transition from single-bunch to multi-bunches mode
3132 static bool flagTransition = false;
3133
3134 // single particle dumping
3135 if (step_m % Options::sptDumpFreq == 0)
3137
3138 // Find out if we need to do bunch injection
3139 injectBunch(flagTransition);
3140
3141// oldReferenceTheta = calculateAngle2(meanR(0), meanR(1));
3142
3143 // Calculate SC field before each time step and keep constant during integration.
3144 // Space Charge effects are included only when total macropaticles number is NOT LESS THAN 1000.
3145 if (itsBunch_m->hasFieldSolver()) {
3146
3147 if (step_m % setup_m.scSolveFreq == 0) {
3149 } else {
3150 // If we are not solving for the space charge fields at this time step
3151 // we will apply the fields from the previous step and have to rotate them
3152 // accordingly. For this we find the quaternion between the previous mean momentum (PreviousMeanP)
3153 // and the current mean momentum (meanP) and rotate the fields with this quaternion.
3154
3155 Vector_t const meanP = calcMeanP();
3156
3157 Quaternion_t quaternionToNewMeanP;
3158
3159 getQuaternionTwoVectors(PreviousMeanP, meanP, quaternionToNewMeanP);
3160
3161 // Reset PreviousMeanP. Cave: This HAS to be after the quaternion is calculated!
3163
3164 // Rotate the fields
3165 globalToLocal(itsBunch_m->Ef, quaternionToNewMeanP);
3166 globalToLocal(itsBunch_m->Bf, quaternionToNewMeanP);
3167 }
3168 }
3169
3170 // *** This was moved here b/c collision should be tested before the **********************
3171 // *** actual timestep (bgf_main_collision_test() predicts the next step automatically) -DW
3172 // Apply the plugin elements: probe, collimator, stripper, septum
3173 bool flagNeedUpdate = applyPluginElements(dt);
3174
3175 // check if we lose particles at the boundary
3177
3179
3180 for (size_t i = 0; i < itsBunch_m->getLocalNum(); i++) {
3181
3182 // used for gap crossing checking
3183 Vector_t Rold = itsBunch_m->R[i]; // [x,y,z] (mm)
3184 Vector_t Pold = itsBunch_m->P[i]; // [px,py,pz] (beta*gamma)
3185
3186 // Integrate for one step in the lab Cartesian frame (absolute value).
3187 itsStepper_mp->advance(itsBunch_m, i, t, dt);
3188
3189 // If gap crossing happens, do momenta kicking (not if gap crossing just happened)
3190 if (itsBunch_m->cavityGapCrossed[i] == true) {
3191 itsBunch_m->cavityGapCrossed[i] = false;
3192 } else {
3193 gapCrossKick_m(i, t, dt, Rold, Pold);
3194 }
3195 flagNeedUpdate |= (itsBunch_m->Bin[i] < 0);
3196 }
3197
3199
3200 // Destroy particles if they are marked as Bin = -1 in the plugin elements
3201 // or out of global aperture
3202 deleteParticle(flagNeedUpdate);
3203
3204 // Recalculate bingamma and reset the BinID for each particles according to its current gamma
3205 if (isMultiBunch() && step_m % Options::rebinFreq == 0) {
3206 mbHandler_m->updateParticleBins(itsBunch_m);
3207 }
3208
3209 // Some status output for user after each turn
3210 if ( isTurnDone() ) {
3211 turnnumber_m++;
3212 finishedTurn = true;
3213
3214 *gmsg << endl;
3215 *gmsg << "*** Finished turn " << turnnumber_m - 1
3216 << ", Total number of live particles: "
3217 << itsBunch_m->getTotalNum() << endl;
3218 }
3219
3221}
3222
3223
3225 double dt,
3226 const Vector_t& Rold,
3227 const Vector_t& Pold) {
3228
3229 for (beamline_list::iterator sindex = ++(FieldDimensions.begin());
3230 sindex != FieldDimensions.end(); ++sindex)
3231 {
3232 bool tag_crossing = false;
3233 double DistOld = 0.0; //mm
3234 RFCavity * rfcav;
3235
3236 if (((*sindex)->first) == ElementType::RFCAVITY) {
3237 // here check gap cross in the list, if do , set tag_crossing to TRUE
3238 rfcav = static_cast<RFCavity *>(((*sindex)->second).second);
3239 tag_crossing = checkGapCross(Rold, itsBunch_m->R[i], rfcav, DistOld);
3240 }
3241
3242 if ( tag_crossing ) {
3243 itsBunch_m->cavityGapCrossed[i] = true;
3244
3245 double oldMomentum2 = dot(Pold, Pold);
3246 double oldBetgam = std::sqrt(oldMomentum2);
3247 double oldGamma = std::sqrt(1.0 + oldMomentum2);
3248 double oldBeta = oldBetgam / oldGamma;
3249 double dt1 = DistOld / (Physics::c * oldBeta * Units::m2mm / Units::s2ns); // ns, c in [mm/ns]
3250 double dt2 = dt - dt1;
3251
3252 // retrack particle from the old postion to cavity gap point
3253 // restore the old coordinates and momenta
3254 itsBunch_m->R[i] = Rold;
3255 itsBunch_m->P[i] = Pold;
3256
3257 if (dt / dt1 < 1.0e9) {
3258 itsStepper_mp->advance(itsBunch_m, i, t, dt1);
3259 }
3260 // Momentum kick
3261 RFkick(rfcav, t, dt1, i);
3262
3263 /* Retrack particle from cavity gap point for
3264 * the left time to finish the entire timestep
3265 */
3266 if (dt / dt2 < 1.0e9) {
3267 itsStepper_mp->advance(itsBunch_m, i, t, dt2);
3268 }
3269 }
3270 }
3271}
3272
3273
3275 const Vector_t& R,
3276 const Vector_t& P,
3277 const double& oldReferenceTheta,
3278 const double& temp_meanTheta) {
3279
3280 for (unsigned int i=0; i<=2; i++) {
3281 if ((oldReferenceTheta < azimuth_angle_m[i] - setup_m.deltaTheta) &&
3282 ( temp_meanTheta >= azimuth_angle_m[i] - setup_m.deltaTheta))
3283 {
3284 *(outfTheta_m[i]) << "#Turn number = " << turnnumber_m
3285 << ", Time = " << t << " [ns]" << std::endl
3286 << " " << std::hypot(R(0), R(1))
3287 << " " << P(0) * std::cos(temp_meanTheta) + P(1) * std::sin(temp_meanTheta)
3288 << " " << temp_meanTheta * Units::rad2deg
3289 << " " << - P(0) * std::sin(temp_meanTheta) + P(1) * std::cos(temp_meanTheta)
3290 << " " << R(2)
3291 << " " << P(2) << std::endl;
3292 }
3293 }
3294}
3295
3297 const Vector_t& R,
3298 const Vector_t& P,
3299 const double& temp_meanTheta,
3300 bool& finishedTurn) {
3301 if ( isTurnDone() ) {
3302 ++turnnumber_m;
3303 finishedTurn = true;
3304
3305 *gmsg << "* SPT: Finished turn " << turnnumber_m - 1 << endl;
3306
3307 *(outfTheta_m[3]) << "#Turn number = " << turnnumber_m
3308 << ", Time = " << t << " [ns]" << std::endl
3309 << " " << std::sqrt(R(0) * R(0) + R(1) * R(1))
3310 << " " << P(0) * std::cos(temp_meanTheta) +
3311 P(1) * std::sin(temp_meanTheta)
3312 << " " << temp_meanTheta / Units::deg2rad
3313 << " " << - P(0) * std::sin(temp_meanTheta) +
3314 P(1) * std::cos(temp_meanTheta)
3315 << " " << R(2)
3316 << " " << P(2) << std::endl;
3317 }
3318}
3319
3320
3322 // Firstly reset E and B to zero before fill new space charge field data for each track step
3323 itsBunch_m->Bf = Vector_t(0.0);
3324 itsBunch_m->Ef = Vector_t(0.0);
3325
3327 // --- Single bunch mode with spiral inflector --- //
3328
3329 // If we are doing a spiral inflector simulation and are using the SAAMG solver
3330 // we don't rotate or translate the bunch and gamma is 1.0 (non-relativistic).
3331
3332 //itsBunch_m->R *= Vector_t(0.001); // mm --> m
3333
3334 itsBunch_m->setGlobalMeanR(Vector_t(0.0, 0.0, 0.0));
3336
3338
3339 //itsBunch_m->R *= Vector_t(1000.0); // m --> mm
3340
3341 } else {
3342
3343 Vector_t const meanR = calcMeanR(); // (m)
3344
3346
3347 // Since calcMeanP takes into account all particles of all bins (TODO: Check this! -DW)
3348 // Using the quaternion method with PreviousMeanP and yaxis should give the correct result
3349
3350 Quaternion_t quaternionToYAxis;
3351
3352 getQuaternionTwoVectors(PreviousMeanP, yaxis, quaternionToYAxis);
3353
3354 globalToLocal(itsBunch_m->R, quaternionToYAxis, meanR);
3355
3356 //itsBunch_m->R *= Vector_t(0.001); // mm --> m
3357
3358 if ((step_m + 1) % Options::boundpDestroyFreq == 0) {
3360 } else {
3361 itsBunch_m->boundp();
3362 }
3363
3364 if (hasMultiBunch()) {
3365 // --- Multibunch mode --- //
3366
3367 // Calculate gamma for each energy bin
3369
3370 repartition();
3371
3372 // Calculate space charge field for each energy bin
3373 for (int b = 0; b < itsBunch_m->getLastemittedBin(); b++) {
3374
3376 //itsBunch_m->setGlobalMeanR(0.001 * meanR);
3377 itsBunch_m->setGlobalMeanR(meanR);
3378 itsBunch_m->setGlobalToLocalQuaternion(quaternionToYAxis);
3380 }
3381
3383
3384 } else {
3385 // --- Single bunch mode --- //
3386
3387 double temp_meangamma = std::sqrt(1.0 + dot(PreviousMeanP, PreviousMeanP));
3388
3389 repartition();
3390
3391 //itsBunch_m->setGlobalMeanR(0.001 * meanR);
3392 itsBunch_m->setGlobalMeanR(meanR);
3393 itsBunch_m->setGlobalToLocalQuaternion(quaternionToYAxis);
3394
3395 itsBunch_m->computeSelfFields_cycl(temp_meangamma);
3396 }
3397
3398 //scale coordinates back
3399 //itsBunch_m->R *= Vector_t(1000.0); // m --> mm
3400
3401 // Transform coordinates back to global
3402 localToGlobal(itsBunch_m->R, quaternionToYAxis, meanR);
3403
3404 // Transform self field back to global frame (rotate only)
3405 localToGlobal(itsBunch_m->Ef, quaternionToYAxis);
3406 localToGlobal(itsBunch_m->Bf, quaternionToYAxis);
3407 }
3408}
3409
3410
3411bool ParallelCyclotronTracker::computeExternalFields_m(const size_t& i, const double& t,
3412 Vector_t& Efield, Vector_t& Bfield) {
3413
3415
3416 // Flag whether a particle is out of field
3417 bool outOfBound = (((*sindex)->second).second)->apply(i, t, Efield, Bfield);
3418
3419 Bfield *= Units::kG2T;
3420 Efield *= Units::kV2V / Units::mm2m;
3421
3422 return outOfBound;
3423}
3424
3426 Vector_t& Efield, Vector_t& Bfield) {
3427
3429 // Flag whether a particle is out of field
3430 bool outOfBound = (((*sindex)->second).second)->apply(R, P, t, Efield, Bfield);
3431
3432 Bfield *= Units::kG2T;
3433 Efield *= Units::kV2V / Units::mm2m;
3434
3435 return outOfBound;
3436}
3437
3438
3439
3440void ParallelCyclotronTracker::injectBunch(bool& flagTransition) {
3442 return;
3443 }
3444
3445 const short result = mbHandler_m->injectBunch(itsBunch_m, itsReference,
3446 flagTransition);
3447
3448 switch ( result ) {
3449 case 0: {
3450 // nothing happened
3451 break;
3452 }
3453 case 1: {
3454 // bunch got saved
3457 if ( flagTransition ) {
3458 *gmsg << "* MBM: Saving beam distribution at turn " << turnnumber_m << endl;
3459 *gmsg << "* MBM: After one revolution, Multi-Bunch Mode will be invoked" << endl;
3460 }
3461 break;
3462 }
3463 case 2: {
3464 // bunch got injected
3466 break;
3467 }
3468 default: {
3469 throw OpalException("ParallelCyclotronTracker::injectBunch()",
3470 "Unknown return value " + std::to_string(result));
3471 }
3472 }
3473}
3474
3475
3477 if ( !isMultiBunch() ) {
3478 return;
3479 }
3480
3481 Vector_t meanR = calcMeanR();
3482
3483 // Bunch (global) angle w.r.t. x-axis (cylinder coordinates)
3484 double theta = calculateAngle(meanR(0), meanR(1)) * Units::rad2deg;
3485
3486 // fix azimuth_m --> make monotonically increasing
3488
3489 const double radius = computeRadius(meanR);
3490
3491 MultiBunchHandler::injection_t& inj = mbHandler_m->getInjectionValues();
3492
3493 inj.time = itsBunch_m->getT() * Units::s2ns;
3495 inj.azimuth = azimuth_m;
3496 inj.radius = radius;
3497}
3498
3499
3501 /* the last element includes all particles,
3502 * all other elements are bunch number specific
3503 */
3504 std::vector<double> lpaths(1);
3505
3506 if ( isMultiBunch() ) {
3507 lpaths.resize(mbHandler_m->getNumBunch() + 1);
3508 }
3509
3510 computePathLengthUpdate(lpaths, dt);
3511
3512 pathLength_m += lpaths.back();
3514
3515 if ( isMultiBunch() ) {
3516 mbHandler_m->updatePathLength(lpaths);
3517 }
3518}
3519
3520
3522 // t is in seconds
3523 double t = itsBunch_m->getT();
3524
3525 itsBunch_m->setT(t + dt * Units::ns2s);
3526
3527 if ( isMultiBunch() ) {
3528 mbHandler_m->updateTime(dt);
3529 }
3530}
3531
3532
3534 if (!isMultiBunch()) {
3535 return;
3536 }
3537
3538 for (short b = 0; b < mbHandler_m->getNumBunch(); ++b) {
3539 Vector_t meanR = calcMeanR(b);
3540 MultiBunchHandler::beaminfo_t& binfo = mbHandler_m->getBunchInfo(b);
3541
3542 binfo.radius = computeRadius(meanR);
3543 double azimuth = calculateAngle(meanR(0), meanR(1)) * Units::rad2deg;
3544 dumpAngle(azimuth, binfo.prevAzimuth, binfo.azimuth, b);
3545 }
3546}
3547
3548
3550 if ( isMultiBunch() ) {
3551 // we need to reset the path length of each bunch
3553 }
3554}
Vector3D cross(const Vector3D &lhs, const Vector3D &rhs)
Vector cross product.
Definition: Vector3D.cpp:111
double dot(const Vector3D &lhs, const Vector3D &rhs)
Vector dot product.
Definition: Vector3D.cpp:118
ElementType
Definition: ElementBase.h:88
Quaternion Quaternion_t
Definition: Quaternion.h:44
Tps< T > cos(const Tps< T > &x)
Cosine.
Definition: TpsMath.h:129
Tps< T > pow(const Tps< T > &x, int y)
Integer power.
Definition: TpsMath.h:76
Tps< T > sin(const Tps< T > &x)
Sine.
Definition: TpsMath.h:111
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition: TpsMath.h:91
T euclidean_norm(const Vector< T > &)
Euclidean norm.
Definition: Vector.h:243
constexpr double c_mmtns
Inform * gmsg
Definition: Main.cpp:61
DistributionType
Definition: Distribution.h:50
void allreduce(const T *input, T *output, int count, Op op)
Definition: GlobalComm.hpp:510
bool reduce(Communicate &, InputIterator, InputIterator, OutputIterator, const ReduceOp &, bool *IncludeVal=0)
Definition: GlobalComm.hpp:55
const int COMM_ANY_NODE
Definition: Communicate.h:40
#define IPPL_APP_CYCLE
Definition: Tags.h:103
#define IPPL_APP_TAG4
Definition: Tags.h:97
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:84
bool for_each(const BareFieldIterator< T, D > &p, SameFieldID s, C)
Definition: AssignDefs.h:30
T::PETE_Expr_t::PETE_Return_t sum(const PETE_Expr< T > &expr)
Definition: PETE.h:1111
PETE_TUTree< FnArcCos, typename T::PETE_Expr_t > acos(const PETE_Expr< T > &l)
Definition: PETE.h:725
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
PETE_TBTree< FnArcTan2, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > atan2(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
Inform & level4(Inform &inf)
Definition: Inform.cpp:48
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
Inform & level3(Inform &inf)
Definition: Inform.cpp:47
#define INFORM_ALL_NODES
Definition: Inform.h:39
#define ERRORMSG(msg)
Definition: IpplInfo.h:350
constexpr double q_e
The elementary charge in As.
Definition: Physics.h:69
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 pi
The value of.
Definition: Physics.h:30
constexpr double mm2m
Definition: Units.h:29
constexpr double m2mm
Definition: Units.h:26
constexpr double ns2s
Definition: Units.h:47
constexpr double eV2MeV
Definition: Units.h:77
constexpr double kV2V
Definition: Units.h:62
constexpr double GeV2eV
Definition: Units.h:68
constexpr double kG2T
Definition: Units.h:59
constexpr double rad2deg
Definition: Units.h:146
constexpr double deg2rad
Definition: Units.h:143
constexpr double s2ns
Definition: Units.h:44
std::string::iterator iterator
Definition: MSLang.h:16
T deg(T x)
Convert radians to degrees.
Definition: matheval.hpp:74
int mtsSubsteps
Definition: Options.cpp:59
int psDumpFreq
The frequency to dump the phase space, i.e.dump data when steppsDumpFreq==0.
Definition: Options.cpp:39
int boundpDestroyFreq
Definition: Options.cpp:87
bool psDumpEachTurn
phase space dump flag for OPAL-cycl
Definition: Options.cpp:43
bool asciidump
Definition: Options.cpp:85
int sptDumpFreq
The frequency to dump single particle trajectory of particles with ID = 0 & 1.
Definition: Options.cpp:47
unsigned int delPartFreq
The frequency to delete particles (currently: OPAL-cycl only)
Definition: Options.cpp:109
int scSolveFreq
The frequency to solve space charge fields.
Definition: Options.cpp:57
int repartFreq
The frequency to do particles repartition for better load balance between nodes.
Definition: Options.cpp:49
int rebinFreq
The frequency to reset energy bin ID for all particles.
Definition: Options.cpp:55
DumpFrame psDumpFrame
flag to decide in which coordinate frame the phase space will be dumped for OPAL-cycl
Definition: Options.cpp:45
int statDumpFreq
The frequency to dump statistical values, e.e. dump data when stepstatDumpFreq==0.
Definition: Options.cpp:41
std::string doubleVectorToString(const std::vector< double > &v)
Definition: Util.cpp:175
std::string boolVectorToUpperString(const std::vector< bool > &b)
Definition: Util.cpp:160
TimeIntegrator
Definition: Steppers.h:25
@ UNDEFINED
Definition: Steppers.h:26
boost::function< boost::tuple< double, bool >(arguments_t)> type
Definition: function.hpp:21
const std::string & getOpalName() const
Return object name.
Definition: Object.cpp:310
ParticleAttrib< Vector_t > Ef
double get_meanKineticEnergy() const
ParticlePos_t & R
ParticleAttrib< int > Bin
Vector_t RefPartP_m
void boundp_destroyCycl()
double getQ() const
Access to reference data.
double getChargePerParticle() const
get the macro particle charge
size_t getLocalNum() const
void set_sPos(double s)
void setLocalTrackStep(long long n)
step in a TRACK command
void setLocalBinCount(size_t num, int bin)
ParticleAttrib< double > M
size_t getTotalNum() const
ParticleAttrib< Vector_t > P
double getdT() const
FieldSolverType getFieldSolverType() const
Return the fieldsolver type if we have a fieldsolver.
void setLocalNumPerBunch(size_t numpart, short n)
size_t getLocalNumPerBunch(short n) const
void setGlobalToLocalQuaternion(Quaternion_t globalToLocalQuaternion)
virtual void computeSelfFields_cycl(double gamma)=0
Vector_t RefPartR_m
ParticleAttrib< double > Q
void calcBeamParameters()
size_t getTotalNumPerBunch(short n) const
ParticleAttrib< short > cavityGapCrossed
void setTotalNum(size_t n)
long long getLocalTrackStep() const
void setGlobalMeanR(Vector_t globalMeanR)
short getNumBunch() const
virtual void do_binaryRepart()
ParticleAttrib< double > dt
int getStepsPerTurn() const
void countTotalNumPerBunch()
DistributionType getDistType() const
void calcGammas_cycl()
virtual void setBinCharge(int bin, double q)
Set the charge of one bin to the value of q and all other to zero.
void performDestroy(bool updateLocalNum=false)
ParticleAttrib< short > bunchNum
void setTotalNumPerBunch(size_t numpart, short n)
void destroy(size_t M, size_t I, bool doNow=false)
ParticleAttrib< Vector_t > Bf
virtual void boundp()
void setT(double t)
ParticleIndex_t & ID
double get_sPos() const
double getM() const
double getT() const
std::string getInputBasename()
get input file name without extension
Definition: OpalData.cpp:669
static OpalData * getInstance()
Definition: OpalData.cpp:196
BoundaryGeometry * getGlobalGeometry()
Definition: OpalData.cpp:461
void setOpenMode(OpenMode openMode)
Definition: OpalData.cpp:349
bool inRestartRun()
true if we do a restart run
Definition: OpalData.cpp:312
int lombAnalysis(double *x, double *y, int Ndat, int nhis)
Definition: Ctunes.cpp:146
void operator()(double x)
void computePathLengthUpdate(std::vector< double > &dl, const double &dt)
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.
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)
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
void normalizeVector(Vector_t &vector)
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)
virtual void visitMultipoleTStraight(const MultipoleTStraight &)
Apply the algorithm to an arbitrary straight multipole.
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)
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)
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
std::vector< double > dvector_t
IpplTimings::TimerRef TransformTimer_m
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 visitMultipoleTCurvedVarRadius(const MultipoleTCurvedVarRadius &)
Apply the algorithm to an arbitrary curved multipole of variable radius.
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.
IpplTimings::TimerRef DumpTimer_m
virtual void visitMultipoleTCurvedConstRadius(const MultipoleTCurvedConstRadius &)
Apply the algorithm to an arbitrary curved multipole of constant radius.
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 visitMonitor(const Monitor &)
Apply the algorithm to a beam position monitor.
static void writeFields(Component *field)
static void writeFields(Component *field)
Definition: DumpFields.cpp:129
double getWidth()
Definition: CCollimator.h:100
double getZStart()
Member variable access.
Definition: CCollimator.h:90
double getZEnd()
Definition: CCollimator.h:95
Interface for a single beam element.
Definition: Component.h:50
Interface for general corrector.
Definition: Corrector.h:35
virtual std::vector< double > getEScale() const
Definition: Cyclotron.cpp:284
virtual double getCyclHarm() const
Definition: Cyclotron.cpp:297
void setBFieldType()
Definition: Cyclotron.cpp:369
virtual bool getSpiralFlag() const
Definition: Cyclotron.cpp:175
virtual double getPRinit() const
Definition: Cyclotron.cpp:135
virtual double getPZinit() const
Definition: Cyclotron.cpp:159
virtual std::string getFieldMapFN() const
Definition: Cyclotron.cpp:183
virtual double getRmin() const
Definition: Cyclotron.cpp:301
virtual double getMaxR() const
Definition: Cyclotron.cpp:325
virtual std::vector< double > getRfFrequ() const
Definition: Cyclotron.cpp:226
virtual double getRmax() const
Definition: Cyclotron.cpp:305
virtual double getMaxZ() const
Definition: Cyclotron.cpp:345
virtual double getBScale() const
Definition: Cyclotron.cpp:276
virtual double getZinit() const
Definition: Cyclotron.cpp:151
virtual double getPHIinit() const
Definition: Cyclotron.cpp:143
virtual double getSymmetry() const
Definition: Cyclotron.cpp:252
virtual double getMinZ() const
Definition: Cyclotron.cpp:335
virtual std::vector< double > getRfPhi() const
Definition: Cyclotron.cpp:213
BFieldType getBFieldType() const
Definition: Cyclotron.cpp:389
virtual void initialise(PartBunchBase< double, 3 > *bunch, double &startField, double &endField)
Definition: Cyclotron.cpp:941
virtual std::vector< bool > getSuperpose() const
Definition: Cyclotron.cpp:239
unsigned int getNumberOfTrimcoils() const
Definition: Cyclotron.cpp:293
virtual double getRinit() const
Definition: Cyclotron.cpp:127
virtual double getMinR() const
Definition: Cyclotron.cpp:321
const std::string & getCyclotronType() const
Definition: Cyclotron.cpp:260
Interface for drift space.
Definition: Drift.h:33
virtual const std::string & getName() const
Get element name.
virtual double getElementLength() const
Get design length.
Definition: ElementBase.h:414
virtual void accept(BeamlineVisitor &visitor) const =0
Apply visitor.
virtual ElementBase * clone() const =0
Return clone.
virtual bool hasAttribute(const std::string &aKey) const
Test for existence of an attribute.
std::string getTypeString() const
Definition: ElementBase.h:578
Interface for a marker.
Definition: Marker.h:32
Interface for general multipole.
Definition: Multipole.h:47
ElementBase * clone() const override
Definition: MultipoleT.cpp:84
virtual ElementBase * clone() const override
virtual ElementBase * clone() const override
virtual ElementBase * clone() const override
Definition: Offset.h:66
void updateGeometry(Vector_t startPosition, Vector_t startDirection)
Definition: Offset.cpp:183
double getYStart() const
virtual void initialise(PartBunchBase< double, 3 > *bunch, double &startField, double &endField) override
Pure virtual implementation of Component.
double getXEnd() const
double getXStart() const
Member variable access.
double getYEnd() const
Definition: Probe.h:28
Definition: RBend.h:58
std::string getPhaseModelName()
Definition: RFCavity.h:465
virtual double getRmax() const
Definition: RFCavity.cpp:294
void getMomentaKick(const double normalRadius, double momentum[], const double t, const double dtCorrt, const int PID, const double restMass, const int chargenumber)
used in OPAL-cycl
Definition: RFCavity.cpp:363
virtual double getAzimuth() const
Definition: RFCavity.cpp:298
std::string getAmplitudeModelName()
Definition: RFCavity.h:450
virtual double getCosAzimuth() const
Definition: RFCavity.cpp:306
virtual void initialise(PartBunchBase< double, 3 > *bunch, double &startField, double &endField) override
Definition: RFCavity.cpp:167
virtual double getSinAzimuth() const
Definition: RFCavity.cpp:302
virtual std::string getFieldMapFN() const
Definition: RFCavity.cpp:335
virtual double getCycFrequency() const
Definition: RFCavity.cpp:349
virtual double getGapWidth() const
Definition: RFCavity.cpp:314
std::string getFrequencyModelName()
Definition: RFCavity.h:480
virtual double getPerpenDistance() const
Definition: RFCavity.cpp:310
CavityType getCavityType() const
Definition: RFCavity.h:415
virtual double getPhi0() const
Definition: RFCavity.cpp:318
std::string getCavityTypeString() const
Definition: RFCavity.cpp:331
virtual double getRmin() const
Definition: RFCavity.cpp:290
Ring describes a ring type geometry for tracking.
Definition: Ring.h:64
double getBeamRInit() const
Definition: Ring.h:235
Vector_t getNextPosition() const
Definition: Ring.cpp:211
Vector_t getNextNormal() const
Definition: Ring.cpp:220
double getBeamPRInit() const
Definition: Ring.h:247
double getSymmetry() const
Definition: Ring.h:286
virtual ElementBase * clone() const override
Definition: Ring.h:155
virtual void initialise(PartBunchBase< double, 3 > *bunch, double &startField, double &endField) override
Definition: Ring.cpp:162
double getBeamPhiInit() const
Definition: Ring.h:241
double getHarmonicNumber()
Definition: Ring.h:222
void lockRing()
Definition: Ring.cpp:289
void appendElement(const Component &element)
Definition: Ring.cpp:229
Definition: SBend.h:68
ScalingFFAMagnet * clone() const override
Definition: Septum.h:23
virtual void initialise(PartBunchBase< double, 3 > *bunch, double &startField, double &endField) override
Pure virtual implementation of Component.
Definition: Septum.cpp:47
double getWidth() const
Definition: Septum.cpp:57
Interface for solenoids.
Definition: Solenoid.h:36
double getOPMass() const
Definition: Stripper.cpp:85
double getOPCharge() const
Definition: Stripper.cpp:81
bool getStop() const
Definition: Stripper.cpp:93
Definition: Vacuum.h:61
double getPressure() const
Definition: Vacuum.cpp:128
std::string getPressureMapFN() const
Definition: Vacuum.cpp:141
double getPScale() const
Definition: Vacuum.cpp:149
virtual void initialise(PartBunchBase< double, 3 > *bunch, double &startField, double &endField) override
Definition: Vacuum.cpp:233
virtual bool checkVacuum(PartBunchBase< double, 3 > *bunch, Cyclotron *cycl)
Definition: Vacuum.cpp:201
double getTemperature() const
Definition: Vacuum.cpp:162
std::string getResidualGasName()
Definition: Vacuum.cpp:115
bool getStop() const
Definition: Vacuum.cpp:175
static std::shared_ptr< AbstractTimeDependence > getTimeDependence(std::string name)
const PartData itsReference
The reference information.
Particle reference data.
Definition: PartData.h:35
double getBeta() const
The relativistic beta per particle.
Definition: PartData.h:124
double getGamma() const
The relativistic gamma per particle.
Definition: PartData.h:129
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
virtual void iterate(BeamlineVisitor &, bool r2l) const
Apply visitor to all elements of the line.
Definition: TBeamline.h:205
void kick(const Vector_t &R, Vector_t &P, const Vector_t &Ef, const Vector_t &Bf, const double &dt) const
Definition: BorisPusher.h:65
Leap-Frog 2nd order.
Definition: LF2.h:26
4-th order Runnge-Kutta stepper
Definition: RK4.h:27
int partInside(const Vector_t &r, const Vector_t &v, const double dt, Vector_t &intecoords, int &triId)
void setMultiBunchInitialPathLengh(MultiBunchHandler *mbhandler_p)
Definition: DataSink.cpp:218
void dumpSDDS(PartBunchBase< double, 3 > *beam, Vector_t FDext[], const double &azimuth=-1) const
Definition: DataSink.cpp:105
void dumpH5(PartBunchBase< double, 3 > *beam, Vector_t FDext[]) const
Definition: DataSink.cpp:86
void writeMultiBunchStatistics(PartBunchBase< double, 3 > *beam, MultiBunchHandler *mbhandler)
Definition: DataSink.cpp:197
The base class for all OPAL exceptions.
Definition: OpalException.h:28
Message * receive_block(int &node, int &tag)
void barrier(void)
Message & put(const T &val)
Definition: Message.h:406
Message & get(const T &cval)
Definition: Message.h:476
int next_tag(int t, int s=1000)
Definition: TagMaker.h:39
Definition: Inform.h:42
static int getNodes()
Definition: IpplInfo.cpp:670
static Communicate * Comm
Definition: IpplInfo.h:84
static TimerRef getTimer(const char *nm)
Definition: IpplTimings.h:182
static void stopTimer(TimerRef t)
Definition: IpplTimings.h:192
static void startTimer(TimerRef t)
Definition: IpplTimings.h:187
Vektor< double, 3 > Vector_t
Definition: Vektor.h:6