OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
OPAL
ParallelTTracker.cpp
Go to the documentation of this file.
1 //
2 // Class ParallelTTracker
3 // OPAL-T tracker.
4 // The visitor class for tracking particles with time as independent
5 // variable.
6 //
7 // Copyright (c) 200x - 2014, Christof Kraus, Paul Scherrer Institut, Villigen PSI, Switzerland
8 // 2015 - 2016, Christof Metzger-Kraus, Helmholtz-Zentrum Berlin, Germany
9 // 2017 - 2020, Christof Metzger-Kraus
10 // All rights reserved
11 //
12 // This file is part of OPAL.
13 //
14 // OPAL is free software: you can redistribute it and/or modify
15 // it under the terms of the GNU General Public License as published by
16 // the Free Software Foundation, either version 3 of the License, or
17 // (at your option) any later version.
18 //
19 // You should have received a copy of the GNU General Public License
20 // along with OPAL. If not, see <https://www.gnu.org/licenses/>.
21 //
23 
24 #include <cfloat>
25 #include <cmath>
26 #include <fstream>
27 #include <limits>
28 #include <iomanip>
29 #include <iostream>
30 #include <sstream>
31 #include <string>
32 
33 #include "AbsBeamline/Monitor.h"
37 #include "BasicActions/Option.h"
38 #ifdef ENABLE_OPAL_FEL
40 #endif
41 #include "Beamlines/Beamline.h"
44 #include "Elements/OpalBeamline.h"
49 #include "Utilities/Options.h"
50 #include "Utilities/Timer.h"
51 #include "Utilities/Util.h"
53 
54 extern Inform* gmsg;
55 
56 class PartData;
57 
59  const PartData &reference,
60  bool revBeam,
61  bool revTrack):
62  Tracker(beamline, reference, revBeam, revTrack),
63  itsDataSink_m(NULL),
64  itsOpalBeamline_m(beamline.getOrigin3D(), beamline.getInitialDirection()),
65  globalEOL_m(false),
66  wakeStatus_m(false),
67  wakeFunction_m(NULL),
68  pathLength_m(0.0),
69  zstart_m(0.0),
70  dtCurrentTrack_m(0.0),
71  minStepforReBin_m(-1),
72  minBinEmitted_m(std::numeric_limits<size_t>::max()),
73  repartFreq_m(-1),
74  emissionSteps_m(std::numeric_limits<unsigned int>::max()),
75  numParticlesInSimulation_m(0),
76  timeIntegrationTimer1_m(IpplTimings::getTimer("TIntegration1")),
77  timeIntegrationTimer2_m(IpplTimings::getTimer("TIntegration2")),
78  fieldEvaluationTimer_m(IpplTimings::getTimer("External field eval")),
79  BinRepartTimer_m(IpplTimings::getTimer("Binaryrepart")),
80  WakeFieldTimer_m(IpplTimings::getTimer("WakeField")),
81  particleMatterStatus_m(false)
82 {}
83 
86  DataSink &ds,
87  const PartData &reference,
88  bool revBeam,
89  bool revTrack,
90  const std::vector<unsigned long long> &maxSteps,
91  double zstart,
92  const std::vector<double> &zstop,
93  const std::vector<double> &dt):
94  Tracker(beamline, bunch, reference, revBeam, revTrack),
95  itsDataSink_m(&ds),
96  itsOpalBeamline_m(beamline.getOrigin3D(), beamline.getInitialDirection()),
97  globalEOL_m(false),
98  wakeStatus_m(false),
99  wakeFunction_m(NULL),
100  pathLength_m(0.0),
101  zstart_m(zstart),
102  dtCurrentTrack_m(0.0),
103  minStepforReBin_m(-1),
104  minBinEmitted_m(std::numeric_limits<size_t>::max()),
105  repartFreq_m(-1),
106  emissionSteps_m(std::numeric_limits<unsigned int>::max()),
107  numParticlesInSimulation_m(0),
108  timeIntegrationTimer1_m(IpplTimings::getTimer("TIntegration1")),
109  timeIntegrationTimer2_m(IpplTimings::getTimer("TIntegration2")),
110  fieldEvaluationTimer_m(IpplTimings::getTimer("External field eval")),
111  BinRepartTimer_m(IpplTimings::getTimer("Binaryrepart")),
112  WakeFieldTimer_m(IpplTimings::getTimer("WakeField")),
113  particleMatterStatus_m(false)
114 {
115  for (unsigned int i = 0; i < zstop.size(); ++ i) {
116  stepSizes_m.push_back(dt[i], zstop[i], maxSteps[i]);
117  }
118 
121 }
122 
124 
126  const FlaggedBeamline* fbl = static_cast<const FlaggedBeamline*>(&bl);
127  if (fbl->getRelativeFlag()) {
128  OpalBeamline stash(fbl->getOrigin3D(), fbl->getInitialDirection());
129  stash.swap(itsOpalBeamline_m);
130  fbl->iterate(*this, false);
133  stash.merge(itsOpalBeamline_m);
134  stash.swap(itsOpalBeamline_m);
135  } else {
136  fbl->iterate(*this, false);
137  }
138 }
139 
140 void ParallelTTracker::updateRFElement(std::string elName, double maxPhase) {
143  cavities.insert(cavities.end(), travelingwaves.begin(), travelingwaves.end());
144 
145  for (FieldList::iterator fit = cavities.begin(); fit != cavities.end(); ++ fit) {
146  if ((*fit).getElement()->getName() == elName) {
147 
148  RFCavity *element = static_cast<RFCavity *>((*fit).getElement().get());
149 
150  element->setPhasem(maxPhase);
151  element->setAutophaseVeto();
152 
153  INFOMSG("Restored cavity phase from the h5 file. Name: " << element->getName() << ", phase: " << maxPhase << " rad" << endl);
154  return;
155  }
156  }
157 }
158 
161 }
162 
165 
166  if (OpalData::getInstance()->hasPriorTrack() ||
167  OpalData::getInstance()->inRestartRun()) {
170  for (; it < end; ++ it) {
171  updateRFElement((*it).first, (*it).second);
172  }
173  }
174 }
175 
177  Inform msg("ParallelTTracker ", *gmsg);
179 
180  BorisPusher pusher(itsReference);
181  const double globalTimeShift = itsBunch_m->weHaveEnergyBins()? OpalData::getInstance()->getGlobalPhaseShift(): 0.0;
183 
184  // the time step needs to be positive in the setup
187 
189 
190  if (OpalData::getInstance()->hasPriorTrack() || OpalData::getInstance()->inRestartRun()) {
192  }
193 
194  prepareSections();
195 
196  double minTimeStep = stepSizes_m.getMinTimeStep();
197 
199 
200  if ( OpalData::getInstance()->hasPriorTrack() ||
201  OpalData::getInstance()->inRestartRun()) {
202 
205 
207  } else {
208 
209  double momentum = euclidean_norm(itsBunch_m->get_pmean_Distribution());
211  itsBunch_m->toLabTrafo_m = beamlineToLab;
212 
213  itsBunch_m->RefPartR_m = beamlineToLab.transformTo(Vector_t(0.0));
214  itsBunch_m->RefPartP_m = beamlineToLab.rotateTo(momentum * Vector_t(0, 0, 1));
215 
216  if (itsBunch_m->getTotalNum() > 0) {
218  momentum = itsReference.getP() / itsBunch_m->getM();
219  itsBunch_m->RefPartP_m = beamlineToLab.rotateTo(momentum * Vector_t(0, 0, 1));
220  }
221 
222  if (zstart_m > pathLength_m) {
223  findStartPosition(pusher);
224  }
225 
227  }
228  }
229 
231  if (back_track) {
235  ++ stepSizes_m;
236  }
237  }
238 
239  Vector_t rmin(0.0), rmax(0.0);
240  if (itsBunch_m->getTotalNum() > 0) {
241  itsBunch_m->get_bounds(rmin, rmax);
242  }
243 
247  pathLength_m,
248  -rmin(2),
249  itsBunch_m->getT(),
250  (back_track? -minTimeStep: minTimeStep),
251  stepSizes_m,
253 
254  oth.execute();
255 
257 
259 
260  setTime();
261 
262  double time = itsBunch_m->getT() - globalTimeShift;
263  itsBunch_m->setT(time);
264 
265  *gmsg << level1 << *itsBunch_m << endl;
266 
267  unsigned long long step = itsBunch_m->getGlobalTrackStep();
268  OPALTimer::Timer myt1;
269  *gmsg << "* Track start at: " << myt1.time() << ", t= " << Util::getTimeString(time) << "; "
270  << "zstart at: " << Util::getLengthString(pathLength_m)
271  << endl;
272 
273  prepareEmission();
274 
275  *gmsg << "* Executing ParallelTTracker\n"
276  << "* Initial dt = " << Util::getTimeString(itsBunch_m->getdT()) << "\n"
277  << "* Max integration steps = " << stepSizes_m.getMaxSteps()
278  << ", next step = " << step << endl << endl;
279 
281 
282  globalEOL_m = false;
283  wakeStatus_m = false;
284  deletedParticles_m = false;
286  while (!stepSizes_m.reachedEnd()) {
287  unsigned long long trackSteps = stepSizes_m.getNumSteps() + step;
290 
291  for (; step < trackSteps; ++ step) {
292  Vector_t rmin(0.0), rmax(0.0);
293  if (itsBunch_m->getTotalNum() > 0) {
294  itsBunch_m->get_bounds(rmin, rmax);
295  }
296 
297  timeIntegration1(pusher);
298 
299  itsBunch_m->Ef = Vector_t(0.0);
300  itsBunch_m->Bf = Vector_t(0.0);
301 
303 
305  emitParticles(step);
307 
309 
310  timeIntegration2(pusher);
311 
313 
314  if (itsBunch_m->getT() > 0.0 ||
315  itsBunch_m->getdT() < 0.0) {
316  updateReference(pusher);
317  }
318 
319  if (deletedParticles_m) {
321  deletedParticles_m = false;
322  }
324 
325  if (hasEndOfLineReached()) break;
326 
327  bool const psDump = ((itsBunch_m->getGlobalTrackStep() % Options::psDumpFreq) + 1 == Options::psDumpFreq);
328  bool const statDump = ((itsBunch_m->getGlobalTrackStep() % Options::statDumpFreq) + 1 == Options::statDumpFreq);
329  dumpStats(step, psDump, statDump);
330 
332 
334  double driftPerTimeStep = std::abs(itsBunch_m->getdT()) * Physics::c * beta;
335  if (std::abs(stepSizes_m.getZStop() - pathLength_m) < 0.5 * driftPerTimeStep) {
336  break;
337  }
338  }
339 
340  if (globalEOL_m)
341  break;
342 
343  ++ stepSizes_m;
344  }
345 
347 
350  }
351 
352  bool const psDump = (((itsBunch_m->getGlobalTrackStep() - 1) % Options::psDumpFreq) + 1 != Options::psDumpFreq);
353  bool const statDump = (((itsBunch_m->getGlobalTrackStep() - 1) % Options::statDumpFreq) + 1 != Options::statDumpFreq);
354  writePhaseSpace((step + 1), psDump, statDump);
355 
356  msg << level2 << "Dump phase space of last step" << endl;
357 
359 
360  OPALTimer::Timer myt3;
361  *gmsg << endl << "* Done executing ParallelTTracker at "
362  << myt3.time() << endl << endl;
363 
365 
367 }
368 
370 
371  itsBeamline_m.accept(*this);
373 
377 }
378 
380 
382  pushParticles(pusher);
384 }
385 
386 
388  /*
389  transport and emit particles
390  that passed the cathode in the first
391  half-step or that would pass it in the
392  second half-step.
393 
394  to make IPPL and the field solver happy
395  make sure that at least 10 particles are emitted
396 
397  also remember that node 0 has
398  all the particles to be emitted
399 
400  this has to be done *after* the calculation of the
401  space charges! thereby we neglect space charge effects
402  in the very first step of a new-born particle.
403 
404  */
405 
407  kickParticles(pusher);
408  //switchElements();
409  pushParticles(pusher);
410 
411  const unsigned int localNum = itsBunch_m->getLocalNum();
412  for (unsigned int i = 0; i < localNum; ++ i) {
413  itsBunch_m->dt[i] = itsBunch_m->getdT();
414  }
415 
417 }
418 
419 void ParallelTTracker::selectDT(bool backTrack) {
420 
421  if (itsBunch_m->getIfBeamEmitting()) {
422  double dt = itsBunch_m->getEmissionDeltaT();
423  itsBunch_m->setdT(dt);
424  } else {
425  double dt = dtCurrentTrack_m;
426  itsBunch_m->setdT(dt);
427  }
428  if (backTrack) {
430  }
431 }
432 
433 void ParallelTTracker::changeDT(bool backTrack) {
434  selectDT(backTrack);
435  const unsigned int localNum = itsBunch_m->getLocalNum();
436  for (unsigned int i = 0; i < localNum; ++ i) {
437  itsBunch_m->dt[i] = itsBunch_m->getdT();
438  }
439 }
440 
441 void ParallelTTracker::emitParticles(long long step) {
442  if (!itsBunch_m->weHaveEnergyBins()) return;
443 
444  if (itsBunch_m->getIfBeamEmitting()) {
447 
448  transformBunch(refToGun);
449 
451 
452  Vector_t externalE = Vector_t(0.0);
453  Vector_t externalB = Vector_t(0.0);
455  gunToLab.rotateTo(Vector_t(0.0)),
456  itsBunch_m->getT() + 0.5 * itsBunch_m->getdT(),
457  externalE,
458  externalB);
459  itsBunch_m->emitParticles(externalE(2));
460 
464 
465  transformBunch(refToGun.inverted());
466  }
467 
468  if (step > minStepforReBin_m) {
469  itsBunch_m->Rebin();
471  }
472 }
473 
474 
475 void ParallelTTracker::computeSpaceChargeFields(unsigned long long step) {
476  if (numParticlesInSimulation_m <= minBinEmitted_m || !itsBunch_m->hasFieldSolver()) {
477  return;
478  }
479 
481  Quaternion alignment = getQuaternion(itsBunch_m->get_pmean(), Vector_t(0, 0, 1));
482  CoordinateSystemTrafo beamToReferenceCSTrafo(Vector_t(0, 0, pathLength_m), alignment.conjugate());
483  CoordinateSystemTrafo referenceToBeamCSTrafo = beamToReferenceCSTrafo.inverted();
484  const unsigned int localNum1 = itsBunch_m->getLocalNum();
485  for (unsigned int i = 0; i < localNum1; ++ i) {
486  itsBunch_m->R[i] = referenceToBeamCSTrafo.transformTo(itsBunch_m->R[i]);
487  }
488 
489  itsBunch_m->boundp();
490 
491  if (step % repartFreq_m + 1 == repartFreq_m) {
493  }
494 
495  if (itsBunch_m->weHaveEnergyBins()) {
499  for (int binNumber = 0; binNumber < itsBunch_m->getLastEmittedEnergyBin() &&
500  binNumber < itsBunch_m->getNumberOfEnergyBins(); ++binNumber) {
501 
502  itsBunch_m->setBinCharge(binNumber);
504  itsBunch_m->computeSelfFields(binNumber);
505  itsBunch_m->Q = Q_back;
506 
507  }
508 
509  } else {
512  }
513 
514  const unsigned int localNum2 = itsBunch_m->getLocalNum();
515  for (unsigned int i = 0; i < localNum2; ++ i) {
516  itsBunch_m->R[i] = beamToReferenceCSTrafo.transformTo(itsBunch_m->R[i]);
517  itsBunch_m->Ef[i] = beamToReferenceCSTrafo.rotateTo(itsBunch_m->Ef[i]);
518  itsBunch_m->Bf[i] = beamToReferenceCSTrafo.rotateTo(itsBunch_m->Bf[i]);
519  }
520 }
521 
522 
525  Inform msg("ParallelTTracker ", *gmsg);
526  const unsigned int localNum = itsBunch_m->getLocalNum();
527  bool locPartOutOfBounds = false, globPartOutOfBounds = false;
528  Vector_t rmin(0.0), rmax(0.0);
529  if (itsBunch_m->getTotalNum() > 0)
530  itsBunch_m->get_bounds(rmin, rmax);
532 
533  try {
534  elements = oth.query(pathLength_m + 0.5 * (rmax(2) + rmin(2)), rmax(2) - rmin(2));
535  } catch(IndexMap::OutOfBounds &e) {
536  globalEOL_m = true;
537  return;
538  }
539 
540  IndexMap::value_t::const_iterator it = elements.begin();
541  const IndexMap::value_t::const_iterator end = elements.end();
542 
543  for (; it != end; ++ it) {
544  CoordinateSystemTrafo refToLocalCSTrafo = (itsOpalBeamline_m.getMisalignment((*it)) *
546  CoordinateSystemTrafo localToRefCSTrafo = refToLocalCSTrafo.inverted();
547 
548  (*it)->setCurrentSCoordinate(pathLength_m + rmin(2));
549 
550  for (unsigned int i = 0; i < localNum; ++ i) {
551  if (itsBunch_m->Bin[i] < 0) continue;
552 
553  itsBunch_m->R[i] = refToLocalCSTrafo.transformTo(itsBunch_m->R[i]);
554  itsBunch_m->P[i] = refToLocalCSTrafo.rotateTo(itsBunch_m->P[i]);
555  double &dt = itsBunch_m->dt[i];
556 
557  Vector_t localE(0.0), localB(0.0);
558 
559  if ((*it)->apply(i, itsBunch_m->getT() + 0.5 * dt, localE, localB)) {
560  itsBunch_m->R[i] = localToRefCSTrafo.transformTo(itsBunch_m->R[i]);
561  itsBunch_m->P[i] = localToRefCSTrafo.rotateTo(itsBunch_m->P[i]);
562  itsBunch_m->Bin[i] = -1;
563  locPartOutOfBounds = true;
564 
565  continue;
566  }
567 
568  itsBunch_m->R[i] = localToRefCSTrafo.transformTo(itsBunch_m->R[i]);
569  itsBunch_m->P[i] = localToRefCSTrafo.rotateTo(itsBunch_m->P[i]);
570  itsBunch_m->Ef[i] += localToRefCSTrafo.rotateTo(localE);
571  itsBunch_m->Bf[i] += localToRefCSTrafo.rotateTo(localB);
572  }
573  }
574 
576 
578 #ifdef ENABLE_OPAL_FEL
579  computeUndulator(elements);
580 #endif
582 
583  reduce(locPartOutOfBounds, globPartOutOfBounds, OpOrAssign());
584 
585  size_t ne = 0;
586  if (globPartOutOfBounds) {
587  if (itsBunch_m->hasFieldSolver()) {
589  } else {
590  ne = itsBunch_m->destroyT();
591  }
593  deletedParticles_m = true;
594  }
595 
596  size_t totalNum = itsBunch_m->getTotalNum();
598  numParticlesInSimulation_m = totalNum;
599  }
600 
601  if (ne > 0) {
602  msg << level1 << "* Deleted " << ne << " particles, "
603  << "remaining " << numParticlesInSimulation_m << " particles" << endl;
604  }
605 }
606 
607 #ifdef ENABLE_OPAL_FEL
608 void ParallelTTracker::computeUndulator(IndexMap::value_t &elements) {
609  // Check if bunch has entered undulator field.
610  UndulatorRep* und;
611  IndexMap::value_t::const_iterator it = elements.begin();
612  for (; it != elements.end(); ++ it)
613  if ((*it)->getType() == ElementBase::UNDULATOR) {
614  und = dynamic_cast<UndulatorRep*>(it->get());
615  if (!und->getHasBeenSimulated())
616  break;
617  }
618  if (it == elements.end())
619  return;
620 
621  // Apply MITHRA full wave solver for undulator.
622  CoordinateSystemTrafo refToLocalCSTrafo = (itsOpalBeamline_m.getMisalignment((*it)) *
624 
625  und->apply(itsBunch_m, refToLocalCSTrafo);
626 
628 }
629 #endif
630 
632  bool hasWake = false;
633  WakeFunction *wfInstance;
634 
635  Inform msg("ParallelTTracker ", *gmsg);
636 
637  const unsigned int localNum = itsBunch_m->getLocalNum();
638 
639  IndexMap::value_t::const_iterator it = elements.begin();
640  const IndexMap::value_t::const_iterator end = elements.end();
641 
642  for (; it != end; ++ it) {
643  if ((*it)->hasWake() && !hasWake) {
645 
646  hasWake = true;
647 
648  if ((*it)->getWake()->getType() == "CSRWakeFunction" ||
649  (*it)->getWake()->getType() == "CSRIGFWakeFunction") {
650  if ((*it)->getType() == ElementBase::RBEND ||
651  (*it)->getType() == ElementBase::SBEND) {
652  wfInstance = (*it)->getWake();
653  wakeFunction_m = wfInstance;
654  } else {
655  wfInstance = wakeFunction_m;
656  }
657  } else {
658  wfInstance = (*it)->getWake();
659  }
660 
661  if (!wfInstance) {
662  throw OpalException("ParallelTTracker::computeWakefield",
663  "empty wake function");
664  }
665 
666  if (!wakeStatus_m) {
667  msg << level2 << "============== START WAKE CALCULATION =============" << endl;
668  wfInstance->initialize((*it).get());
669  wakeStatus_m = true;
670  }
671 
673 
674  Quaternion alignment = getQuaternion(itsBunch_m->get_pmean(), Vector_t(0, 0, 1));
675  CoordinateSystemTrafo referenceToBeamCSTrafo(Vector_t(0.0), alignment);
676  CoordinateSystemTrafo beamToReferenceCSTrafo = referenceToBeamCSTrafo.inverted();
677 
678  for (unsigned int i = 0; i < localNum; ++ i) {
679  itsBunch_m->R[i] = referenceToBeamCSTrafo.transformTo(itsBunch_m->R[i]);
680  itsBunch_m->P[i] = referenceToBeamCSTrafo.rotateTo(itsBunch_m->P[i]);
681  itsBunch_m->Ef[i] = referenceToBeamCSTrafo.rotateTo(itsBunch_m->Ef[i]);
682  }
683 
684  wfInstance->apply(itsBunch_m);
685 
686  for (unsigned int i = 0; i < localNum; ++ i) {
687  itsBunch_m->R[i] = beamToReferenceCSTrafo.transformTo(itsBunch_m->R[i]);
688  itsBunch_m->P[i] = beamToReferenceCSTrafo.rotateTo(itsBunch_m->P[i]);
689  itsBunch_m->Ef[i] = beamToReferenceCSTrafo.rotateTo(itsBunch_m->Ef[i]);
690  }
691 
693  }
694  }
695 
696  if (wakeStatus_m && !hasWake) {
697  msg << level2 << "=============== END WAKE CALCULATION ==============" << endl;
698  wakeStatus_m = false;
699  }
700 }
701 
703  Inform msg("ParallelTTracker ", *gmsg);
704  std::set<IndexMap::value_t::value_type> elementsWithParticleMatterInteraction;
705  std::set<ParticleMatterInteractionHandler*> particleMatterinteractionHandlers;
706  IndexMap::key_t currentRange{0.0, 0.0};
707 
708  while (elements.size() > 0) {
709  auto it = elements.begin();
710  if ((*it)->hasParticleMatterInteraction()) {
711  elementsWithParticleMatterInteraction.insert(*it);
712  particleMatterinteractionHandlers.insert((*it)->getParticleMatterInteraction());
713 
714  IndexMap::key_t range = oth.getRange(*it, pathLength_m);
715  currentRange.begin = std::min(currentRange.begin, range.begin);
716  currentRange.end = std::max(currentRange.end, range.end);
717 
718  IndexMap::value_t touching = oth.getTouchingElements(range);
719  elements.insert(touching.begin(), touching.end());
720  }
721 
722  elements.erase(it);
723  }
724 
725  if (elementsWithParticleMatterInteraction.size() > 0) {
726  std::set<ParticleMatterInteractionHandler*> oldSPHandlers;
727  std::vector<ParticleMatterInteractionHandler*> leftBehindSPHandlers, newSPHandlers;
729  oldSPHandlers.insert(it);
730  }
731 
732  leftBehindSPHandlers.resize(std::max(oldSPHandlers.size(),
733  particleMatterinteractionHandlers.size()));
734  auto last = std::set_difference(oldSPHandlers.begin(), oldSPHandlers.end(),
735  particleMatterinteractionHandlers.begin(), particleMatterinteractionHandlers.end(),
736  leftBehindSPHandlers.begin());
737  leftBehindSPHandlers.resize(last - leftBehindSPHandlers.begin());
738 
739  for (auto it: leftBehindSPHandlers) {
740  if (!it->stillActive()) {
742  }
743  }
744 
745  newSPHandlers.resize(std::max(oldSPHandlers.size(),
746  elementsWithParticleMatterInteraction.size()));
747  last = std::set_difference(particleMatterinteractionHandlers.begin(), particleMatterinteractionHandlers.end(),
748  oldSPHandlers.begin(), oldSPHandlers.end(),
749  newSPHandlers.begin());
750  newSPHandlers.resize(last - newSPHandlers.begin());
751 
752  for (auto it: newSPHandlers) {
754  }
755 
757  msg << level2 << "============== START PARTICLE MATTER INTERACTION CALCULATION =============" << endl;
758  particleMatterStatus_m = true;
759  }
760  }
761 
763  do {
765  //check if more than one degrader has particles
766  ParticleMatterInteractionHandler* onlyDegraderWithParticles = NULL;
767  int degradersWithParticlesCount = 0;
769  it->setFlagAllParticlesIn(false);
770  if (it->getParticlesInMat() > 0) {
771  onlyDegraderWithParticles = it;
772  ++ degradersWithParticlesCount;
773  }
774  }
775 
776  //if max particles per node is 2, and only one degrader has particles set
777  //AllParticlesIn for this degrader to true
778  unsigned int localNum = itsBunch_m->getLocalNum();
779  unsigned int totalNum = 0;
780  reduce(localNum, totalNum, OpAddAssign());
781  bool allParticlesInMat = (totalNum == 0 &&
782  degradersWithParticlesCount == 1);
783 
784  if (allParticlesInMat) {
785  onlyDegraderWithParticles->setFlagAllParticlesIn(true);
786  msg << "All particles in degrader" << endl;
787  }
788 
789  auto boundingSphere = itsBunch_m->getLocalBoundingSphere();
790  unsigned int rediffusedParticles = 0;
791  unsigned int numEnteredParticles = 0;
793  ElementBase* element = it->getElement();
794  CoordinateSystemTrafo refToLocalCSTrafo = (element->getMisalignment() *
796  CoordinateSystemTrafo localToRefCSTrafo = refToLocalCSTrafo.inverted();
797 
798  const unsigned int localNum = itsBunch_m->getLocalNum();
799  for (unsigned int i = 0; i < localNum; ++i) {
800  itsBunch_m->R[i] = refToLocalCSTrafo.transformTo(itsBunch_m->R[i]);
801  itsBunch_m->P[i] = refToLocalCSTrafo.rotateTo(itsBunch_m->P[i]);
802  }
803  boundingSphere.first = refToLocalCSTrafo.transformTo(boundingSphere.first);
804 
805  it->apply(itsBunch_m, boundingSphere);
806  it->print(msg);
807 
808  boundingSphere.first = localToRefCSTrafo.transformTo(boundingSphere.first);
809 
810  const unsigned int newLocalNum = itsBunch_m->getLocalNum();
811  for (unsigned int i = 0; i < newLocalNum; ++i) {
812  itsBunch_m->R[i] = localToRefCSTrafo.transformTo(itsBunch_m->R[i]);
813  itsBunch_m->P[i] = localToRefCSTrafo.rotateTo(itsBunch_m->P[i]);
814  }
815 
816  rediffusedParticles += it->getRediffused();
817  numEnteredParticles += it->getNumEntered();
818 
819  //if all particles were in material update time to time in degrader
820  if (it->getFlagAllParticlesIn()) {
821  double timeDifference = it->getTime() - itsBunch_m->getdT() - itsBunch_m->getT();
822  if (timeDifference > 0.0) {
823  const unsigned int numSteps = std::ceil(timeDifference / itsBunch_m->getdT());
824  const double origdT = itsBunch_m->getdT();
825  itsBunch_m->setdT(timeDifference / numSteps);
826  BorisPusher pusher(itsReference);
827  for (unsigned int i = 0; i < numSteps; ++ i) {
828  updateReference(pusher);
829  }
830  itsBunch_m->setdT(origdT);
831  }
832  itsBunch_m->setT(it->getTime() - itsBunch_m->getdT());
833  }
834  }
835 
836  //perform boundp only if there are particles in the bunch, or there are particles
837  //comming out of the degrader
838  if (numEnteredParticles > 0 || rediffusedParticles > 0) {
839  totalNum -= (numEnteredParticles + rediffusedParticles);
840  if (totalNum > minBinEmitted_m && itsBunch_m->hasFieldSolver()) {
841  itsBunch_m->boundp();
842  } else {
844  }
845  }
846 
847  //if bunch has no particles update time to time in degrader
848  if (itsBunch_m->getTotalNum() == 0)
850 
851  } while (itsBunch_m->getTotalNum() == 0);
852 
853 
855  msg << level2 << "============== END PARTICLE MATTER INTERACTION CALCULATION =============" << endl;
856  particleMatterStatus_m = false;
857  }
858  }
859 }
860 
862  size_t particles_or_bins = std::max(minBinEmitted_m, size_t(1000));
863  if (itsBunch_m->hasFieldSolver() && numParticlesInSimulation_m > particles_or_bins) {
864 
865  INFOMSG("*****************************************************************" << endl);
866  INFOMSG("do repartition because of repartFreq_m" << endl);
867  INFOMSG("*****************************************************************" << endl);
870  Ippl::Comm->barrier();
872  INFOMSG("*****************************************************************" << endl);
873  INFOMSG("do repartition done" << endl);
874  INFOMSG("*****************************************************************" << endl);
875  }
876 }
877 
878 void ParallelTTracker::dumpStats(long long step, bool psDump, bool statDump) {
879  OPALTimer::Timer myt2;
880  Inform msg("ParallelTTracker ", *gmsg);
881 
882  if (itsBunch_m->getGlobalTrackStep() % 1000 + 1 == 1000) {
883  msg << level1;
884  } else if (itsBunch_m->getGlobalTrackStep() % 100 + 1 == 100) {
885  msg << level2;
886  } else {
887  msg << level3;
888  }
889 
890  if (numParticlesInSimulation_m == 0) {
891  msg << myt2.time() << " "
892  << "Step " << std::setw(6) << itsBunch_m->getGlobalTrackStep() << "; "
893  << " -- no emission yet -- "
894  << "t= " << Util::getTimeString(itsBunch_m->getT())
895  << endl;
896  return;
897  }
898 
900  size_t totalParticles_f = numParticlesInSimulation_m;
901  if (totalParticles_f <= minBinEmitted_m) {
902  msg << myt2.time() << " "
903  << "Step " << std::setw(6) << itsBunch_m->getGlobalTrackStep() << "; "
904  << "only " << std::setw(4) << totalParticles_f << " particles emitted; "
905  << "t= " << Util::getTimeString(itsBunch_m->getT()) << ", "
907  << endl;
909  throw OpalException("ParallelTTracker::dumpStats()",
910  "there seems to be something wrong with the position of the bunch!");
911  } else {
912 
913  msg << myt2.time() << " "
914  << "Step " << std::setw(6) << itsBunch_m->getGlobalTrackStep() << " "
915  << "at " << Util::getLengthString(pathLength_m) << ", "
916  << "t= " << Util::getTimeString(itsBunch_m->getT()) << ", "
918  << endl;
919 
920  writePhaseSpace(step, psDump, statDump);
921  }
922 }
923 
924 
926  Inform msg("ParallelTTracker ", *gmsg);
927 
929  RealVariable *ar = dynamic_cast<RealVariable *>(OpalData::getInstance()->find("MINBINEMITTED"));
930  if (ar)
931  minBinEmitted_m = static_cast<size_t>(ar->getReal());
932  msg << level2 << "MINBINEMITTED " << minBinEmitted_m << endl;
933 
935  RealVariable *br = dynamic_cast<RealVariable *>(OpalData::getInstance()->find("MINSTEPFORREBIN"));
936  if (br)
937  minStepforReBin_m = static_cast<int>(br->getReal());
938  msg << level2 << "MINSTEPFORREBIN " << minStepforReBin_m << endl;
939 
940  // there is no point to do repartitioning with one node
941  if (Ippl::getNodes() == 1) {
943  } else {
945  RealVariable *rep = dynamic_cast<RealVariable *>(OpalData::getInstance()->find("REPARTFREQ"));
946  if (rep)
947  repartFreq_m = static_cast<int>(rep->getReal());
948  msg << level2 << "REPARTFREQ " << repartFreq_m << endl;
949  }
950 }
951 
952 
955  return globalEOL_m;
956 }
957 
959  const unsigned int localNum = itsBunch_m->getLocalNum();
960  for (unsigned int i = 0; i < localNum; ++i) {
961  itsBunch_m->dt[i] = itsBunch_m->getdT();
962  }
963 }
964 
966  Inform msg("ParallelTTracker ", *gmsg);
967 
968  if (!itsBunch_m->doEmission()) return;
969 
971  msg << level2 << "Do emission for " << itsBunch_m->getTEmission() << " [s] using "
972  << itsBunch_m->getNumberOfEnergyBins() << " energy bins " << endl
973  << "Change dT from " << itsBunch_m->getdT() << " [s] to "
974  << itsBunch_m->getEmissionDeltaT() << " [s] during emission " << endl;;
975 
976 }
977 
978 void ParallelTTracker::writePhaseSpace(const long long /*step*/, bool psDump, bool statDump) {
979  extern Inform *gmsg;
980  Inform msg("OPAL ", *gmsg);
981  Vector_t externalE, externalB;
982  Vector_t FDext[2]; // FDext = {BHead, EHead, BRef, ERef, BTail, ETail}.
983 
984  // Sample fields at (xmin, ymin, zmin), (xmax, ymax, zmax) and the centroid location. We
985  // are sampling the electric and magnetic fields at the back, front and
986  // center of the beam.
987  Vector_t rmin, rmax;
988  itsBunch_m->get_bounds(rmin, rmax);
989 
990  if (psDump || statDump) {
991  externalB = Vector_t(0.0);
992  externalE = Vector_t(0.0);
995  itsBunch_m->getT() - 0.5 * itsBunch_m->getdT(),
996  externalE,
997  externalB);
998  FDext[0] = itsBunch_m->toLabTrafo_m.rotateFrom(externalB);
999  FDext[1] = itsBunch_m->toLabTrafo_m.rotateFrom(externalE * 1e-6);
1000  }
1001 
1002  if (statDump) {
1003  std::vector<std::pair<std::string, unsigned int> > collimatorLosses;
1005  if (collimators.size() != 0) {
1006  for (FieldList::iterator it = collimators.begin(); it != collimators.end(); ++ it) {
1007  FlexibleCollimator* coll = static_cast<FlexibleCollimator*>(it->getElement().get());
1008  std::string name = coll->getName();
1009  unsigned int losses = coll->getLosses();
1010  collimatorLosses.push_back(std::make_pair(name, losses));
1011  }
1012  std::sort(collimatorLosses.begin(), collimatorLosses.end(),
1013  [](const std::pair<std::string, unsigned int>& a, const std::pair<std::string, unsigned int>& b) ->bool {
1014  return a.first < b.first;
1015  });
1016  std::vector<unsigned int> bareLosses(collimatorLosses.size(),0);
1017  for (size_t i = 0; i < collimatorLosses.size(); ++ i){
1018  bareLosses[i] = collimatorLosses[i].second;
1019  }
1020 
1021  reduce(&bareLosses[0], &bareLosses[0] + bareLosses.size(), &bareLosses[0], OpAddAssign());
1022 
1023  for (size_t i = 0; i < collimatorLosses.size(); ++ i){
1024  collimatorLosses[i].second = bareLosses[i];
1025  }
1026  }
1027  // Write statistical data.
1028  itsDataSink_m->dumpSDDS(itsBunch_m, FDext, collimatorLosses);
1029 
1030  msg << level3 << "* Wrote beam statistics." << endl;
1031  }
1032 
1033  if (psDump && (itsBunch_m->getTotalNum() > 0)) {
1034  // Write fields to .h5 file.
1035  const size_t localNum = itsBunch_m->getLocalNum();
1036  double distToLastStop = stepSizes_m.getFinalZStop() - pathLength_m;
1039  bool driftToCorrectPosition = std::abs(distToLastStop) < 0.5 * euclidean_norm(driftPerTimeStep);
1040  Ppos_t stashedR;
1041  Vector_t stashedRefPartR;
1042 
1043  if (driftToCorrectPosition) {
1044  const double tau = distToLastStop / euclidean_norm(driftPerTimeStep) * itsBunch_m->getdT();
1045  if (localNum > 0) {
1046  stashedR.create(localNum);
1047  stashedR = itsBunch_m->R;
1048  stashedRefPartR = itsBunch_m->RefPartR_m;
1049 
1050  for (size_t i = 0; i < localNum; ++ i) {
1051  itsBunch_m->R[i] += tau * (Physics::c * itsBunch_m->P[i] / Util::getGamma(itsBunch_m->P[i]) -
1052  driftPerTimeStep / itsBunch_m->getdT());
1053  }
1054  }
1055 
1056  driftPerTimeStep = itsBunch_m->toLabTrafo_m.rotateTo(driftPerTimeStep);
1057  itsBunch_m->RefPartR_m = itsBunch_m->RefPartR_m + tau * driftPerTimeStep / itsBunch_m->getdT();
1058  CoordinateSystemTrafo update(tau * driftPerTimeStep / itsBunch_m->getdT(),
1059  Quaternion(1.0, 0.0, 0.0, 0.0));
1061 
1063 
1065  }
1066  if (!statDump && !driftToCorrectPosition) itsBunch_m->calcBeamParameters();
1067 
1068  msg << *itsBunch_m << endl;
1069  itsDataSink_m->dumpH5(itsBunch_m, FDext);
1070 
1071  if (driftToCorrectPosition) {
1072  if (localNum > 0) {
1073  itsBunch_m->R = stashedR;
1074  stashedR.destroy(localNum, 0);
1075  }
1076 
1077  itsBunch_m->RefPartR_m = stashedRefPartR;
1079 
1081  }
1082 
1083  msg << level2 << "* Wrote beam phase space." << endl;
1084  }
1085 }
1086 
1088  updateReferenceParticle(pusher);
1090 }
1091 
1093  const double direction = back_track? -1: 1;
1094  const double dt = direction * std::min(itsBunch_m->getT(),
1095  direction * itsBunch_m->getdT());
1096  const double scaleFactor = Physics::c * dt;
1097  Vector_t Ef(0.0), Bf(0.0);
1098 
1099  itsBunch_m->RefPartR_m /= scaleFactor;
1101  itsBunch_m->RefPartR_m *= scaleFactor;
1102 
1104  IndexMap::value_t::const_iterator it = elements.begin();
1105  const IndexMap::value_t::const_iterator end = elements.end();
1106 
1107  for (; it != end; ++ it) {
1108  const CoordinateSystemTrafo &refToLocalCSTrafo = itsOpalBeamline_m.getCSTrafoLab2Local((*it));
1109 
1110  Vector_t localR = refToLocalCSTrafo.transformTo(itsBunch_m->RefPartR_m);
1111  Vector_t localP = refToLocalCSTrafo.rotateTo(itsBunch_m->RefPartP_m);
1112  Vector_t localE(0.0), localB(0.0);
1113 
1114  if ((*it)->applyToReferenceParticle(localR,
1115  localP,
1116  itsBunch_m->getT() - 0.5 * dt,
1117  localE,
1118  localB)) {
1119  *gmsg << level1 << "The reference particle hit an element" << endl;
1120  globalEOL_m = true;
1121  }
1122 
1123  Ef += refToLocalCSTrafo.rotateFrom(localE);
1124  Bf += refToLocalCSTrafo.rotateFrom(localB);
1125  }
1126 
1127  pusher.kick(itsBunch_m->RefPartR_m, itsBunch_m->RefPartP_m, Ef, Bf, dt);
1128 
1129  itsBunch_m->RefPartR_m /= scaleFactor;
1131  itsBunch_m->RefPartR_m *= scaleFactor;
1132 }
1133 
1135  const unsigned int localNum = itsBunch_m->getLocalNum();
1136  for (unsigned int i = 0; i < localNum; ++i) {
1137  itsBunch_m->R[i] = trafo.transformTo(itsBunch_m->R[i]);
1138  itsBunch_m->P[i] = trafo.rotateTo(itsBunch_m->P[i]);
1139  itsBunch_m->Ef[i] = trafo.rotateTo(itsBunch_m->Ef[i]);
1140  itsBunch_m->Bf[i] = trafo.rotateTo(itsBunch_m->Bf[i]);
1141  }
1142 }
1143 
1147 
1149 
1150  CoordinateSystemTrafo update(R, getQuaternion(P, Vector_t(0, 0, 1)));
1151 
1152  transformBunch(update);
1153 
1155 }
1156 
1157 void ParallelTTracker::applyFractionalStep(const BorisPusher &pusher, double tau) {
1158  double t = itsBunch_m->getT();
1159  t += tau;
1160  itsBunch_m->setT(t);
1161 
1162  // the push method below pushes for half a time step. Hence the ref particle
1163  // should be pushed for 2 * tau
1164  itsBunch_m->RefPartR_m /= (Physics::c * 2 * tau);
1166  itsBunch_m->RefPartR_m *= (Physics::c * 2 * tau);
1167 
1169 
1172  CoordinateSystemTrafo update(R, getQuaternion(P, Vector_t(0, 0, 1)));
1174 }
1175 
1177 
1178  StepSizeConfig stepSizesCopy(stepSizes_m);
1179  if (back_track) {
1180  stepSizesCopy.shiftZStopLeft(zstart_m);
1181  }
1182 
1183  double t = 0.0;
1184  itsBunch_m->setT(t);
1185 
1186  dtCurrentTrack_m = stepSizesCopy.getdT();
1187  selectDT();
1188 
1191  double gamma = 0.1 / itsBunch_m->getM() + 1.0;
1192  double beta = sqrt(1.0 - 1.0 / std::pow(gamma, 2));
1194  }
1195 
1196  while (true) {
1197  autophaseCavities(pusher);
1198 
1199  t += itsBunch_m->getdT();
1200  itsBunch_m->setT(t);
1201 
1202  Vector_t oldR = itsBunch_m->RefPartR_m;
1203  updateReferenceParticle(pusher);
1205 
1207 
1208  if (pathLength_m > stepSizesCopy.getZStop()) {
1209  ++ stepSizesCopy;
1210 
1211  if (stepSizesCopy.reachedEnd()) {
1212  -- stepSizesCopy;
1213  double tau = (stepSizesCopy.getZStop() - pathLength_m) / speed;
1214  applyFractionalStep(pusher, tau);
1215 
1216  break;
1217  }
1218 
1219  dtCurrentTrack_m = stepSizesCopy.getdT();
1220  selectDT();
1221  }
1222 
1223  if (std::abs(pathLength_m - zstart_m) <= 0.5 * itsBunch_m->getdT() * speed) {
1224  double tau = (zstart_m - pathLength_m) / speed;
1225  applyFractionalStep(pusher, tau);
1226 
1227  break;
1228  }
1229  }
1230 
1231  changeDT();
1232 }
1233 
1235 
1236  double t = itsBunch_m->getT();
1238  pusher.push(nextR, itsBunch_m->RefPartP_m, itsBunch_m->getdT());
1239  nextR *= Physics::c * itsBunch_m->getdT();
1240 
1241  auto elementSet = itsOpalBeamline_m.getElements(nextR);
1242  for (auto element: elementSet) {
1243  if (element->getType() == ElementBase::TRAVELINGWAVE) {
1244  const TravelingWave *TWelement = static_cast<const TravelingWave *>(element.get());
1245  if (!TWelement->getAutophaseVeto()) {
1246  CavityAutophaser ap(itsReference, element);
1249  t, itsBunch_m->getdT());
1250  }
1251 
1252  } else if (element->getType() == ElementBase::RFCAVITY) {
1253  const RFCavity *RFelement = static_cast<const RFCavity *>(element.get());
1254  if (!RFelement->getAutophaseVeto()) {
1255  CavityAutophaser ap(itsReference, element);
1258  t, itsBunch_m->getdT());
1259  }
1260  }
1261  }
1262 }
1263 
1265  unsigned int who;
1266  unsigned int whom;
1267  unsigned int howMany;
1268 };
1269 
1271  const int numNodes = Ippl::getNodes();
1272  if (itsBunch_m->hasFieldSolver() || numNodes == 1) return;
1273 
1274  long onAverage = itsBunch_m->getTotalNum() / Ippl::getNodes();
1275  if (itsBunch_m->getTotalNum() % Ippl::getNodes() > 0.5 * Ippl::getNodes())
1276  ++ onAverage;
1277 
1278  std::vector<long> localParticles(numNodes, 0);
1279  localParticles[Ippl::myNode()] = itsBunch_m->getLocalNum() - onAverage;
1280  allreduce(&(localParticles[0]),
1281  numNodes,
1282  std::plus<long>());
1283 
1284  std::vector<DistributionInfo> send;
1285  std::vector<DistributionInfo> receive;
1286 
1287  for (int i = 0; i < Ippl::getNodes(); ++ i) {
1288  if (localParticles[i] <= 0) continue;
1289 
1290  for (int j = 0; j < Ippl::getNodes(); ++ j) {
1291  if (j == i || localParticles[j] >= 0) continue;
1292 
1293  long numParts = std::min(localParticles[i], -localParticles[j]);
1294  localParticles[i] -= numParts;
1295  localParticles[j] += numParts;
1296 
1297  if (i == Ippl::myNode() || j == Ippl::myNode()) {
1298  DistributionInfo request;
1299  request.who = i;
1300  request.whom = j;
1301  request.howMany = numParts;
1302 
1303  if (i == Ippl::myNode()) {
1304  send.push_back(request);
1305  } else {
1306  receive.push_back(request);
1307  }
1308  }
1309 
1310  if (localParticles[i] == 0) break;
1311  }
1312  }
1313 
1314  std::vector<MPI_Request> requests;
1315  const long sizeSingleParticle = 9 * sizeof(double) + sizeof(short) + sizeof(int) + sizeof(PID_t::Return_t);
1316  long idx = itsBunch_m->getLocalNum() - 1;
1318 
1319  std::vector<char> send_msgbuf;
1320 
1321  if (send.size() > 0) {
1322  const char *buffer;
1323 
1324  unsigned int totalSend = 0, startIndex = 0;
1325  for (DistributionInfo &request: send) {
1326  totalSend += request.howMany;
1327  }
1328  send_msgbuf.reserve(totalSend * sizeSingleParticle);
1329 
1330  for (DistributionInfo &request: send) {
1331  size_t sizePrior = send_msgbuf.size();
1332  for (long i = 0; i < request.howMany; ++ i, -- idx) {
1333  buffer = reinterpret_cast<const char*>(&(itsBunch_m->R[idx](0)));
1334  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + 3 * sizeof(double));
1335  buffer = reinterpret_cast<const char*>(&(itsBunch_m->P[idx](0)));
1336  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + 3 * sizeof(double));
1337  buffer = reinterpret_cast<const char*>(&(itsBunch_m->Q[idx]));
1338  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(double));
1339  buffer = reinterpret_cast<const char*>(&(itsBunch_m->M[idx]));
1340  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(double));
1341  buffer = reinterpret_cast<const char*>(&(itsBunch_m->dt[idx]));
1342  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(double));
1343  buffer = reinterpret_cast<const char*>(&(itsBunch_m->POrigin[idx]));
1344  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(ParticleOrigin));
1345  buffer = reinterpret_cast<const char*>(&(itsBunch_m->TriID[idx]));
1346  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(int));
1347  buffer = reinterpret_cast<const char*>(&(itsBunch_m->ID[idx]));
1348  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(PID_t::Return_t));
1349  }
1350 
1351  size_t sendSizeThis = send_msgbuf.size() - sizePrior;
1352  MPI_Request req = Ippl::Comm->raw_isend(&(send_msgbuf[startIndex]),
1353  sendSizeThis,
1354  request.whom,
1355  tag);
1356 
1357  requests.push_back(req);
1358 
1359  startIndex += sendSizeThis;
1360  }
1361 
1362  itsBunch_m->destroy(totalSend, idx + 1, true);
1363  }
1364 
1365  for (unsigned int i = 0; i < receive.size(); ++ i) {
1366  int node = Communicate::COMM_ANY_NODE;
1367  char *recvbuf;
1368  const int bufsize = Ippl::Comm->raw_probe_receive(recvbuf, node, tag);
1369 
1370  int j = 0;
1371 
1372  while (j < bufsize) {
1373  ++ idx;
1374  itsBunch_m->create(1);
1375  {
1376  const double *buffer = reinterpret_cast<const double*>(recvbuf + j);
1377  itsBunch_m->R[idx] = Vector_t(buffer[0], buffer[1], buffer[2]);
1378  itsBunch_m->P[idx] = Vector_t(buffer[3], buffer[4], buffer[5]);
1379  itsBunch_m->Q[idx] = buffer[6];
1380  itsBunch_m->M[idx] = buffer[7];
1381  itsBunch_m->dt[idx] = buffer[8];
1382  }
1383  j += 9 * sizeof(double);
1384 
1385  {
1386  const ParticleOrigin *buffer = reinterpret_cast<const ParticleOrigin*>(recvbuf + j);
1387  itsBunch_m->POrigin[idx] = buffer[0];
1388  }
1389  j += sizeof(ParticleOrigin);
1390 
1391  {
1392  const int *buffer = reinterpret_cast<const int*>(recvbuf + j);
1393  itsBunch_m->TriID[idx] = buffer[0];
1394  }
1395  j += sizeof(int);
1396 
1397  {
1398  const PID_t::Return_t *buffer = reinterpret_cast<const PID_t::Return_t*>(recvbuf + j);
1399  itsBunch_m->ID[idx] = buffer[0];
1400  }
1401  j += sizeof(PID_t::Return_t);
1402  }
1403 
1404  delete[] recvbuf;
1405  }
1406 
1407  if (requests.size() > 0) {
1408  MPI_Waitall(requests.size(), &(requests[0]), MPI_STATUSES_IGNORE);
1409  }
1410 }
elements
Definition: IndexMap.cpp:163
Inform * gmsg
Definition: Main.cpp:62
Tps< T > pow(const Tps< T > &x, int y)
Integer power.
Definition: TpsMath.h:76
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
ParticleOrigin
Definition: PartBunchBase.h:46
PartBunchBase< T, Dim >::ConstIterator end(PartBunchBase< T, Dim > const &bunch)
Quaternion getQuaternion(Vector_t u, Vector_t ref)
Definition: Quaternion.cpp:34
std::list< ClassicField > FieldList
Definition: ClassicField.h:42
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:84
T::PETE_Expr_t::PETE_Return_t min(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:76
PETE_TBTree< OpNE, Index::PETE_Expr_t, PETE_Scalar< double > > ne(const Index &idx, double x)
Definition: IndexInlines.h:357
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
#define P_SPATIAL_TRANSFER_TAG
Definition: Tags.h:82
#define P_LAYOUT_CYCLE
Definition: Tags.h:86
std::complex< double > a
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
PETE_TBTree< FnCopysign, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > copysign(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
PETE_TUTree< FnCeil, typename T::PETE_Expr_t > ceil(const PETE_Expr< T > &l)
Definition: PETE.h:728
Inform & level1(Inform &inf)
Definition: Inform.cpp:45
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
Inform & level2(Inform &inf)
Definition: Inform.cpp:46
Inform & level3(Inform &inf)
Definition: Inform.cpp:47
#define INFOMSG(msg)
Definition: IpplInfo.h:348
const std::string name
std::string::const_iterator iterator_t
Definition: array.cpp:19
constexpr double e
The value of.
Definition: Physics.h:39
constexpr double c
The velocity of light in m/s.
Definition: Physics.h:51
std::string::iterator iterator
Definition: MSLang.h:16
T isinf(T x)
isinf function with adjusted return type
Definition: matheval.hpp:70
T isnan(T x)
isnan function with adjusted return type
Definition: matheval.hpp:66
int psDumpFreq
The frequency to dump the phase space, i.e.dump data when steppsDumpFreq==0.
Definition: Options.cpp:39
int minBinEmitted
The number of bins that have to be emitted before the bin are squashed into a single bin.
Definition: Options.cpp:51
int minStepForRebin
The number of steps into the simulation before the bins are squashed into a single bin.
Definition: Options.cpp:53
int repartFreq
The frequency to do particles repartition for better load balance between nodes.
Definition: Options.cpp:49
int statDumpFreq
The frequency to dump statistical values, e.e. dump data when stepstatDumpFreq==0.
Definition: Options.cpp:41
double getKineticEnergy(Vector_t p, double mass)
Definition: Util.h:37
std::string getEnergyString(double energyInMeV, unsigned int precision=3)
Definition: Util.h:128
std::string getTimeString(double time, unsigned int precision=3)
Definition: Util.h:55
double getGamma(Vector_t p)
Definition: Util.h:27
std::string getLengthString(double spos, unsigned int precision=3)
Definition: Util.h:78
ParticleAttrib< Vector_t > Ef
double get_meanKineticEnergy() const
ParticlePos_t & R
virtual void resetInterpolationCache(bool clearCache=false)
ParticleAttrib< int > Bin
Vector_t RefPartP_m
virtual void computeSelfFields()=0
size_t getLocalNum() const
void set_sPos(double s)
ParticleAttrib< double > M
size_t emitParticles(double eZ)
Emit particles in the given bin i.e. copy the particles from the bin structure into the particle cont...
size_t getTotalNum() const
bool weHaveEnergyBins()
size_t boundp_destroyT()
ParticleAttrib< Vector_t > P
void switchToUnitlessPositions(bool use_dt_per_particle=false)
double getTEmission()
double getdT() const
double getEmissionDeltaT()
ParticleAttrib< ParticleOrigin > POrigin
bool getIfBeamEmitting()
Vector_t RefPartR_m
ParticleAttrib< double > Q
int getLastEmittedEnergyBin()
void calcBeamParameters()
Vector_t get_pmean_Distribution() const
size_t getNumberOfEmissionSteps()
CoordinateSystemTrafo toLabTrafo_m
void setGlobalMeanR(Vector_t globalMeanR)
virtual void do_binaryRepart()
ParticleAttrib< int > TriID
void calcGammas()
Compute the gammas of all bins.
ParticleAttrib< double > dt
void setdT(double dt)
void get_bounds(Vector_t &rmin, Vector_t &rmax)
virtual void setBinCharge(int bin, double q)
Set the charge of one bin to the value of q and all other to zero.
int getNumberOfEnergyBins()
Vector_t get_centroid() const
Vector_t get_pmean() const
void destroy(size_t M, size_t I, bool doNow=false)
ParticleAttrib< Vector_t > Bf
virtual void boundp()
void create(size_t M)
long long getGlobalTrackStep() const
void setT(double t)
std::pair< Vector_t, double > getLocalBoundingSphere()
ParticleIndex_t & ID
double get_sPos() const
double getM() const
double getT() const
void switchOffUnitlessPositions(bool use_dt_per_particle=false)
void setOpenMode(OPENMODE openMode)
Definition: OpalData.cpp:348
std::vector< MaxPhasesT >::iterator getLastMaxPhases()
Definition: OpalData.cpp:401
double getGlobalPhaseShift()
units: (sec)
Definition: OpalData.cpp:451
std::vector< MaxPhasesT >::iterator getFirstMaxPhases()
Definition: OpalData.cpp:397
Object * find(const std::string &name)
Find entry.
Definition: OpalData.cpp:565
void setInPrepState(bool state)
Definition: OpalData.cpp:299
void setPriorTrack(const bool &value=true)
true if in follow-up track
Definition: OpalData.cpp:307
void setGlobalPhaseShift(double shift)
units: (sec)
Definition: OpalData.cpp:446
static OpalData * getInstance()
Definition: OpalData.cpp:195
double getPhaseAtMaxEnergy(const Vector_t &R, const Vector_t &P, double t, double dt)
std::set< std::shared_ptr< Component > > value_t
Definition: IndexMap.h:47
first_type begin
Definition: IndexMap.h:43
second_type end
Definition: IndexMap.h:44
IndexMap::value_t query(IndexMap::key_t::first_type step, IndexMap::key_t::second_type length)
IndexMap::value_t getTouchingElements(const IndexMap::key_t &range) const
IndexMap::key_t getRange(const IndexMap::value_t::value_type &element, double position) const
void selectDT(bool backTrack=false)
void computeExternalFields(OrbitThreader &oth)
IpplTimings::TimerRef timeIntegrationTimer2_m
IpplTimings::TimerRef WakeFieldTimer_m
OpalBeamline itsOpalBeamline_m
double zstart_m
where to start
void transformBunch(const CoordinateSystemTrafo &trafo)
void autophaseCavities(const BorisPusher &pusher)
void computeWakefield(IndexMap::value_t &elements)
void computeSpaceChargeFields(unsigned long long step)
virtual void visitBeamline(const Beamline &)
Apply the algorithm to a beam line.
void timeIntegration2(BorisPusher &pusher)
DataSink * itsDataSink_m
void updateReferenceParticle(const BorisPusher &pusher)
std::set< ParticleMatterInteractionHandler * > activeParticleMatterInteractionHandlers_m
WakeFunction * wakeFunction_m
size_t numParticlesInSimulation_m
virtual void execute()
Apply the algorithm to the top-level beamline.
IpplTimings::TimerRef timeIntegrationTimer1_m
void pushParticles(const BorisPusher &pusher)
void dumpStats(long long step, bool psDump, bool statDump)
void changeDT(bool backTrack=false)
void emitParticles(long long step)
void findStartPosition(const BorisPusher &pusher)
void computeParticleMatterInteraction(IndexMap::value_t elements, OrbitThreader &oth)
unsigned int repartFreq_m
void updateRFElement(std::string elName, double maxPhi)
void timeIntegration1(BorisPusher &pusher)
IpplTimings::TimerRef BinRepartTimer_m
void writePhaseSpace(const long long step, bool psDump, bool statDump)
StepSizeConfig stepSizes_m
IpplTimings::TimerRef fieldEvaluationTimer_m
void applyFractionalStep(const BorisPusher &pusher, double tau)
void kickParticles(const BorisPusher &pusher)
void updateReference(const BorisPusher &pusher)
unsigned int emissionSteps_m
double getMinTimeStep() const
StepSizeConfig & advanceToPos(double spos)
double getdT() const
void shiftZStopLeft(double back)
double getZStop() const
void push_back(double dt, double zstop, unsigned long numSteps)
bool reachedEnd() const
unsigned long getNumSteps() const
unsigned long long getMaxSteps() const
void sortAscendingZStop()
void resetIterator()
double getFinalZStop() const
virtual const std::string & getName() const
Get element name.
void getMisalignment(double &x, double &y, double &s) const
virtual void accept(BeamlineVisitor &visitor) const =0
Apply visitor.
CoordinateSystemTrafo getCSTrafoGlobal2Local() const
Definition: ElementBase.h:515
unsigned int getLosses() const
static void writeStatistics()
Definition: Monitor.cpp:206
Interface for RF cavity.
Definition: RFCavity.h:37
virtual bool getAutophaseVeto() const
Definition: RFCavity.h:453
virtual void setPhasem(double phase)
Definition: RFCavity.h:403
virtual void setAutophaseVeto(bool veto=true)
Definition: RFCavity.h:448
Interface for Traveling Wave.
Definition: TravelingWave.h:37
void apply(PartBunchBase< double, 3 > *itsBunch, CoordinateSystemTrafo const &refToLocalCSTrafo)
Definition: Undulator.cpp:86
bool getHasBeenSimulated() const
Definition: Undulator.cpp:346
const PartData itsReference
The reference information.
Vector_t rotateFrom(const Vector_t &r) const
Vector_t transformFrom(const Vector_t &r) const
Vector_t transformTo(const Vector_t &r) const
Vector_t rotateTo(const Vector_t &r) const
CoordinateSystemTrafo inverted() const
Particle reference data.
Definition: PartData.h:35
double getP() const
The constant reference momentum per particle.
Definition: PartData.h:114
double gamma
Definition: PartData.h:97
double beta
Definition: PartData.h:96
Quaternion conjugate() const
Definition: Quaternion.h:105
const Beamline & itsBeamline_m
Definition: Tracker.h:122
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
bool getRelativeFlag() const
Definition: TBeamline.h:451
Quaternion getInitialDirection() const
Definition: TBeamline.h:441
Vector_t getOrigin3D() const
Definition: TBeamline.h:431
virtual void iterate(BeamlineVisitor &, bool r2l) const
Apply visitor to all elements of the line.
Definition: TBeamline.h:203
virtual void initialize(const ElementBase *)
Definition: WakeFunction.h:35
virtual void apply(PartBunchBase< double, 3 > *bunch)=0
FieldList getElementByType(ElementBase::ElementType)
Vector_t rotateToLocalCS(const std::shared_ptr< Component > &comp, const Vector_t &r) const
Definition: OpalBeamline.h:178
CoordinateSystemTrafo getMisalignment(const std::shared_ptr< Component > &comp) const
Definition: OpalBeamline.h:200
void merge(OpalBeamline &rhs)
Vector_t transformToLocalCS(const std::shared_ptr< Component > &comp, const Vector_t &r) const
Definition: OpalBeamline.h:166
bool containsSource()
Definition: OpalBeamline.h:205
void prepareSections()
void compute3DLattice()
CoordinateSystemTrafo getCSTrafoLab2Local(const std::shared_ptr< Component > &comp) const
Definition: OpalBeamline.h:190
void switchElementsOff()
void save3DInput()
unsigned long getFieldAt(const unsigned int &, const Vector_t &, const long &, const double &, Vector_t &, Vector_t &)
void save3DLattice()
void activateElements()
std::set< std::shared_ptr< Component > > getElements(const Vector_t &x)
void swap(OpalBeamline &rhs)
void kick(const Vector_t &R, Vector_t &P, const Vector_t &Ef, const Vector_t &Bf, const double &dt) const
Definition: BorisPusher.h:65
void push(Vector_t &R, const Vector_t &P, const double &dt) const
Definition: BorisPusher.h:131
void dumpSDDS(PartBunchBase< double, 3 > *beam, Vector_t FDext[], const double &azimuth=-1) const
Definition: DataSink.cpp:109
void storeCavityInformation()
Write cavity information from H5 file.
Definition: DataSink.cpp:138
void dumpH5(PartBunchBase< double, 3 > *beam, Vector_t FDext[]) const
Definition: DataSink.cpp:90
The base class for all OPAL exceptions.
Definition: OpalException.h:28
std::string time() const
Return time.
Definition: Timer.cpp:42
virtual double getReal() const
Return value.
virtual void destroy(size_t M, size_t I, bool optDestroy=true)
virtual void create(size_t)
virtual MPI_Request raw_isend(void *, int, int, int)
Definition: Communicate.h:196
virtual int raw_probe_receive(char *&, int &, int &)
Definition: Communicate.h:208
void barrier(void)
int next_tag(int t, int s=1000)
Definition: TagMaker.h:39
Definition: Inform.h:42
static int getNodes()
Definition: IpplInfo.cpp:670
static int myNode()
Definition: IpplInfo.cpp:691
static Communicate * Comm
Definition: IpplInfo.h:84
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