OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
ParallelTTracker.cpp
Go to the documentation of this file.
1 // ------------------------------------------------------------------------
2 // $RCSfile: ParallelTTracker.cpp,v $
3 // ------------------------------------------------------------------------
4 // $Revision: 1.1.2.1 $
5 // ------------------------------------------------------------------------
6 // Copyright: see Copyright.readme
7 // ------------------------------------------------------------------------
8 //
9 // Class: ParallelTTracker
10 // The visitor class for tracking particles with time as independent
11 // variable.
12 //
13 // ------------------------------------------------------------------------
14 //
15 // $Date: 2004/11/12 20:10:11 $
16 // $Author: adelmann $
17 //
18 // ------------------------------------------------------------------------
19 
21 
22 #include <cfloat>
23 #include <iostream>
24 #include <fstream>
25 #include <iomanip>
26 #include <sstream>
27 #include <string>
28 #include <limits>
29 #include <cmath>
30 
33 #include "Beamlines/Beamline.h"
35 #include "Lines/Sequence.h"
36 
38 
40 
41 #include "BasicActions/Option.h"
42 #include "Utilities/Options.h"
43 #include "Utilities/Util.h"
44 
47 #include "Utilities/Timer.h"
51 #include "AbsBeamline/Monitor.h"
52 
53 class PartData;
54 
56  const PartData &reference,
57  bool revBeam,
58  bool revTrack):
59  Tracker(beamline, reference, revBeam, revTrack),
60  itsDataSink_m(NULL),
61  itsOpalBeamline_m(beamline.getOrigin3D(), beamline.getInitialDirection()),
62  globalEOL_m(false),
63  wakeStatus_m(false),
64  wakeFunction_m(NULL),
65  pathLength_m(0.0),
66  zstart_m(0.0),
67  dtCurrentTrack_m(0.0),
68  minStepforReBin_m(-1),
69  minBinEmitted_m(std::numeric_limits<size_t>::max()),
70  repartFreq_m(-1),
71  emissionSteps_m(std::numeric_limits<unsigned int>::max()),
72  numParticlesInSimulation_m(0),
73  timeIntegrationTimer1_m(IpplTimings::getTimer("TIntegration1")),
74  timeIntegrationTimer2_m(IpplTimings::getTimer("TIntegration2")),
75  fieldEvaluationTimer_m(IpplTimings::getTimer("External field eval")),
76  BinRepartTimer_m(IpplTimings::getTimer("Binaryrepart")),
77  WakeFieldTimer_m(IpplTimings::getTimer("WakeField")),
78  particleMatterStatus_m(false),
79  totalParticlesInSimulation_m(0)
80 {
81 
82 #ifdef OPAL_DKS
84  setupDKS();
85 #endif
86 }
87 
90  DataSink &ds,
91  const PartData &reference,
92  bool revBeam,
93  bool revTrack,
94  const std::vector<unsigned long long> &maxSteps,
95  double zstart,
96  const std::vector<double> &zstop,
97  const std::vector<double> &dt):
98  Tracker(beamline, bunch, reference, revBeam, revTrack),
99  itsDataSink_m(&ds),
100  itsOpalBeamline_m(beamline.getOrigin3D(), beamline.getInitialDirection()),
101  globalEOL_m(false),
102  wakeStatus_m(false),
103  wakeFunction_m(NULL),
104  pathLength_m(0.0),
105  zstart_m(zstart),
106  dtCurrentTrack_m(0.0),
107  minStepforReBin_m(-1),
108  minBinEmitted_m(std::numeric_limits<size_t>::max()),
109  repartFreq_m(-1),
110  emissionSteps_m(std::numeric_limits<unsigned int>::max()),
111  numParticlesInSimulation_m(0),
112  timeIntegrationTimer1_m(IpplTimings::getTimer("TIntegration1")),
113  timeIntegrationTimer2_m(IpplTimings::getTimer("TIntegration2")),
114  fieldEvaluationTimer_m(IpplTimings::getTimer("External field eval")),
115  BinRepartTimer_m(IpplTimings::getTimer("Binaryrepart")),
116  WakeFieldTimer_m(IpplTimings::getTimer("WakeField")),
117  particleMatterStatus_m(false),
118  totalParticlesInSimulation_m(0)
119 {
120  for (unsigned int i = 0; i < zstop.size(); ++ i) {
121  stepSizes_m.push_back(dt[i], zstop[i], maxSteps[i]);
122  }
123 
126 
127 #ifdef OPAL_DKS
129  setupDKS();
130 #endif
131 }
132 
134 #ifdef OPAL_DKS
136  delete dksbase;
137 #endif
138 }
139 
141  const FlaggedBeamline* fbl = static_cast<const FlaggedBeamline*>(&bl);
142  if (fbl->getRelativeFlag()) {
143  OpalBeamline stash(fbl->getOrigin3D(), fbl->getInitialDirection());
144  stash.swap(itsOpalBeamline_m);
145  fbl->iterate(*this, false);
148  stash.merge(itsOpalBeamline_m);
149  stash.swap(itsOpalBeamline_m);
150  } else {
151  fbl->iterate(*this, false);
152  }
153 }
154 
155 void ParallelTTracker::updateRFElement(std::string elName, double maxPhase) {
158  cavities.insert(cavities.end(), travelingwaves.begin(), travelingwaves.end());
159 
160  for (FieldList::iterator fit = cavities.begin(); fit != cavities.end(); ++ fit) {
161  if ((*fit).getElement()->getName() == elName) {
162 
163  RFCavity *element = static_cast<RFCavity *>((*fit).getElement().get());
164 
165  element->setPhasem(maxPhase);
166  element->setAutophaseVeto();
167 
168  INFOMSG("Restored cavity phase from the h5 file. Name: " << element->getName() << ", phase: " << maxPhase << " rad" << endl);
169  return;
170  }
171  }
172 }
173 
176 }
177 
180 
181  if (OpalData::getInstance()->hasPriorTrack() ||
182  OpalData::getInstance()->inRestartRun()) {
183  iterator_t it = OpalData::getInstance()->getFirstMaxPhases();
184  iterator_t end = OpalData::getInstance()->getLastMaxPhases();
185  for (; it < end; ++ it) {
186  updateRFElement((*it).first, (*it).second);
187  }
188  }
189 }
190 
192  Inform msg("ParallelTTracker ", *gmsg);
194 
195  BorisPusher pusher(itsReference);
196  const double globalTimeShift = itsBunch_m->weHaveEnergyBins()? OpalData::getInstance()->getGlobalPhaseShift(): 0.0;
198 
199  // the time step needs to be positive in the setup
202 
204 
205  if (OpalData::getInstance()->hasPriorTrack() || OpalData::getInstance()->inRestartRun()) {
207  }
208 
209  prepareSections();
210 
211  double minTimeStep = stepSizes_m.getMinTimeStep();
212  unsigned long long totalNumSteps = stepSizes_m.getNumStepsFinestResolution();
213 
215 
216  if ( OpalData::getInstance()->hasPriorTrack() ||
217  OpalData::getInstance()->inRestartRun()) {
218 
221 
223  } else {
224 
225  double momentum = euclidean_norm(itsBunch_m->get_pmean_Distribution());
227  itsBunch_m->toLabTrafo_m = beamlineToLab;
228 
229  itsBunch_m->RefPartR_m = beamlineToLab.transformTo(Vector_t(0.0));
230  itsBunch_m->RefPartP_m = beamlineToLab.rotateTo(momentum * Vector_t(0, 0, 1));
231 
232  if (itsBunch_m->getTotalNum() > 0) {
234  momentum = itsReference.getP() / itsBunch_m->getM();
235  itsBunch_m->RefPartP_m = beamlineToLab.rotateTo(momentum * Vector_t(0, 0, 1));
236  }
237 
238  if (zstart_m > pathLength_m) {
239  findStartPosition(pusher);
240  }
241 
243  }
244  }
245 
247  if (back_track) {
251  ++ stepSizes_m;
252  }
253  }
254 
255  Vector_t rmin(0.0), rmax(0.0);
256  if (itsBunch_m->getTotalNum() > 0) {
257  itsBunch_m->get_bounds(rmin, rmax);
258  }
259 
263  pathLength_m,
264  -rmin(2),
265  itsBunch_m->getT(),
266  (back_track? -minTimeStep: minTimeStep),
267  totalNumSteps,
268  stepSizes_m.getFinalZStop() + 2 * rmax(2),
270 
271  oth.execute();
272 
274 
277 
278  setTime();
279 
280  double time = itsBunch_m->getT() - globalTimeShift;
281  itsBunch_m->setT(time);
282 
283  *gmsg << level1 << *itsBunch_m << endl;
284 
285  unsigned long long step = itsBunch_m->getGlobalTrackStep();
286  OPALTimer::Timer myt1;
287  *gmsg << "Track start at: " << myt1.time() << ", t= " << Util::getTimeString(time) << "; "
288  << "zstart at: " << Util::getLengthString(pathLength_m)
289  << endl;
290 
291  prepareEmission();
292 
293  *gmsg << level1
294  << "Executing ParallelTTracker, initial dt= " << Util::getTimeString(itsBunch_m->getdT()) << ";\n"
295  << "max integration steps " << stepSizes_m.getMaxSteps() << ", next step= " << step << endl;
296 
298 
299 #ifdef OPAL_DKS
301  allocateDeviceMemory();
302 #endif
303 
304  globalEOL_m = false;
305  wakeStatus_m = false;
306  deletedParticles_m = false;
308  while (!stepSizes_m.reachedEnd()) {
309  unsigned long long trackSteps = stepSizes_m.getNumSteps() + step;
312 
313  for (; step < trackSteps; ++ step) {
314  Vector_t rmin(0.0), rmax(0.0);
315  if (itsBunch_m->getTotalNum() > 0) {
316  itsBunch_m->get_bounds(rmin, rmax);
317  }
318 
319  timeIntegration1(pusher);
320 
321  itsBunch_m->Ef = Vector_t(0.0);
322  itsBunch_m->Bf = Vector_t(0.0);
323 
325 
327  emitParticles(step);
329 
331 
332  timeIntegration2(pusher);
333 
335 
336  if (itsBunch_m->getT() > 0.0 ||
337  itsBunch_m->getdT() < 0.0) {
338  updateReference(pusher);
339  }
340 
341  if (deletedParticles_m) {
343  deletedParticles_m = false;
344  }
346 
347  if (hasEndOfLineReached()) break;
348 
349  bool const psDump = ((itsBunch_m->getGlobalTrackStep() % Options::psDumpFreq) + 1 == Options::psDumpFreq);
350  bool const statDump = ((itsBunch_m->getGlobalTrackStep() % Options::statDumpFreq) + 1 == Options::statDumpFreq);
351  dumpStats(step, psDump, statDump);
352 
354 
356  double driftPerTimeStep = std::abs(itsBunch_m->getdT()) * Physics::c * beta;
357  if (std::abs(stepSizes_m.getZStop() - pathLength_m) < 0.5 * driftPerTimeStep) {
358  break;
359  }
360  }
361 
362  if (globalEOL_m)
363  break;
364 
365  ++ stepSizes_m;
366  }
367 
369 
372  }
373 
374  bool const psDump = (((itsBunch_m->getGlobalTrackStep() - 1) % Options::psDumpFreq) + 1 != Options::psDumpFreq);
375  bool const statDump = (((itsBunch_m->getGlobalTrackStep() - 1) % Options::statDumpFreq) + 1 != Options::statDumpFreq);
376  writePhaseSpace((step + 1), psDump, statDump);
377 
378  msg << level2 << "Dump phase space of last step" << endl;
379 
381 
382  OPALTimer::Timer myt3;
383  *gmsg << "done executing ParallelTTracker at " << myt3.time() << endl;
384 
385 #ifdef OPAL_DKS
387  freeDeviceMemory();
388 #endif
389 
391 
393 }
394 
396 
397  itsBeamline_m.accept(*this);
399 
403 }
404 
406 
408 #ifdef OPAL_DKS
410  pushParticlesDKS();
411  else
412  pushParticles(pusher);
413 #else
414  pushParticles(pusher);
415 #endif
416 
418 }
419 
420 
422  /*
423  transport and emit particles
424  that passed the cathode in the first
425  half-step or that would pass it in the
426  second half-step.
427 
428  to make IPPL and the field solver happy
429  make sure that at least 10 particles are emitted
430 
431  also remember that node 0 has
432  all the particles to be emitted
433 
434  this has to be done *after* the calculation of the
435  space charges! thereby we neglect space charge effects
436  in the very first step of a new-born particle.
437 
438  */
439 
441 #ifdef OPAL_DKS
442  if (IpplInfo::DKSEnabled) {
443  kickParticlesDKS();
444  pushParticlesDKS(false);
445  } else {
446  kickParticles(pusher);
447  pushParticles(pusher);
448  }
449 #else
450  kickParticles(pusher);
451 
452  //switchElements();
453  pushParticles(pusher);
454 #endif
455 
456  const unsigned int localNum = itsBunch_m->getLocalNum();
457  for (unsigned int i = 0; i < localNum; ++ i) {
458  itsBunch_m->dt[i] = itsBunch_m->getdT();
459  }
460 
462 }
463 
464 void ParallelTTracker::selectDT(bool backTrack) {
465 
466  if (itsBunch_m->getIfBeamEmitting()) {
467  double dt = itsBunch_m->getEmissionDeltaT();
468  itsBunch_m->setdT(dt);
469  } else {
470  double dt = dtCurrentTrack_m;
471  itsBunch_m->setdT(dt);
472  }
473  if (backTrack) {
475  }
476 }
477 
478 void ParallelTTracker::changeDT(bool backTrack) {
479  selectDT(backTrack);
480  const unsigned int localNum = itsBunch_m->getLocalNum();
481  for (unsigned int i = 0; i < localNum; ++ i) {
482  itsBunch_m->dt[i] = itsBunch_m->getdT();
483  }
484 }
485 
486 void ParallelTTracker::emitParticles(long long step) {
487  if (!itsBunch_m->weHaveEnergyBins()) return;
488 
489  if (itsBunch_m->getIfBeamEmitting()) {
492 
493  transformBunch(refToGun);
494 
496 
497  Vector_t externalE = Vector_t(0.0);
498  Vector_t externalB = Vector_t(0.0);
500  gunToLab.rotateTo(Vector_t(0.0)),
501  itsBunch_m->getT() + 0.5 * itsBunch_m->getdT(),
502  externalE,
503  externalB);
504  itsBunch_m->emitParticles(externalE(2));
505 
509 
510  transformBunch(refToGun.inverted());
511  }
512 
513  if (step > minStepforReBin_m) {
514  itsBunch_m->Rebin();
516  }
517 }
518 
519 
520 void ParallelTTracker::computeSpaceChargeFields(unsigned long long step) {
522 
523  if (!itsBunch_m->hasFieldSolver()) return;
524 
526  Quaternion alignment = getQuaternion(itsBunch_m->get_pmean(), Vector_t(0, 0, 1));
527  CoordinateSystemTrafo beamToReferenceCSTrafo(Vector_t(0, 0, pathLength_m), alignment.conjugate());
528  CoordinateSystemTrafo referenceToBeamCSTrafo = beamToReferenceCSTrafo.inverted();
529 
530  const unsigned int localNum1 = itsBunch_m->getLocalNum();
531  for (unsigned int i = 0; i < localNum1; ++ i) {
532  itsBunch_m->R[i] = referenceToBeamCSTrafo.transformTo(itsBunch_m->R[i]);
533  }
534 
535  itsBunch_m->boundp();
536 
537  if (step % repartFreq_m + 1 == repartFreq_m) {
539  }
540 
541  if (itsBunch_m->weHaveEnergyBins()) {
545  for (int binNumber = 0; binNumber < itsBunch_m->getLastEmittedEnergyBin() &&
546  binNumber < itsBunch_m->getNumberOfEnergyBins(); ++binNumber) {
547 
548  itsBunch_m->setBinCharge(binNumber);
550  itsBunch_m->computeSelfFields(binNumber);
551  itsBunch_m->Q = Q_back;
552 
553  }
554 
555  } else {
558  }
559 
560  const unsigned int localNum2 = itsBunch_m->getLocalNum();
561  for (unsigned int i = 0; i < localNum2; ++ i) {
562  itsBunch_m->R[i] = beamToReferenceCSTrafo.transformTo(itsBunch_m->R[i]);
563  itsBunch_m->Ef[i] = beamToReferenceCSTrafo.rotateTo(itsBunch_m->Ef[i]);
564  itsBunch_m->Bf[i] = beamToReferenceCSTrafo.rotateTo(itsBunch_m->Bf[i]);
565  }
566 }
567 
568 
571  Inform msg("ParallelTTracker ", *gmsg);
572  const unsigned int localNum = itsBunch_m->getLocalNum();
573  bool locPartOutOfBounds = false, globPartOutOfBounds = false;
574  Vector_t rmin(0.0), rmax(0.0);
575  if (itsBunch_m->getTotalNum() > 0)
576  itsBunch_m->get_bounds(rmin, rmax);
578 
579  try {
580  elements = oth.query(pathLength_m + 0.5 * (rmax(2) + rmin(2)), rmax(2) - rmin(2));
581  } catch(IndexMap::OutOfBounds &e) {
582  globalEOL_m = true;
583  return;
584  }
585 
586  IndexMap::value_t::const_iterator it = elements.begin();
587  const IndexMap::value_t::const_iterator end = elements.end();
588 
589  for (; it != end; ++ it) {
590  CoordinateSystemTrafo refToLocalCSTrafo = (itsOpalBeamline_m.getMisalignment((*it)) *
592  CoordinateSystemTrafo localToRefCSTrafo = refToLocalCSTrafo.inverted();
593 
594  (*it)->setCurrentSCoordinate(pathLength_m + rmin(2));
595 
596  for (unsigned int i = 0; i < localNum; ++ i) {
597  if (itsBunch_m->Bin[i] < 0) continue;
598 
599  itsBunch_m->R[i] = refToLocalCSTrafo.transformTo(itsBunch_m->R[i]);
600  itsBunch_m->P[i] = refToLocalCSTrafo.rotateTo(itsBunch_m->P[i]);
601  double &dt = itsBunch_m->dt[i];
602 
603  Vector_t localE(0.0), localB(0.0);
604 
605  if ((*it)->apply(i, itsBunch_m->getT() + 0.5 * dt, localE, localB)) {
606  itsBunch_m->R[i] = localToRefCSTrafo.transformTo(itsBunch_m->R[i]);
607  itsBunch_m->P[i] = localToRefCSTrafo.rotateTo(itsBunch_m->P[i]);
608  itsBunch_m->Bin[i] = -1;
609  locPartOutOfBounds = true;
610 
611  continue;
612  }
613 
614  itsBunch_m->R[i] = localToRefCSTrafo.transformTo(itsBunch_m->R[i]);
615  itsBunch_m->P[i] = localToRefCSTrafo.rotateTo(itsBunch_m->P[i]);
616  itsBunch_m->Ef[i] += localToRefCSTrafo.rotateTo(localE);
617  itsBunch_m->Bf[i] += localToRefCSTrafo.rotateTo(localB);
618  }
619  }
620 
622 
623  computeWakefield(elements);
624  computeParticleMatterInteraction(elements, oth);
625 
626  reduce(locPartOutOfBounds, globPartOutOfBounds, OpOrAssign());
627 
628  size_t ne = 0;
629  if (globPartOutOfBounds) {
630  if (itsBunch_m->hasFieldSolver()) {
631  ne = itsBunch_m->boundp_destroyT();
632  } else {
633  ne = itsBunch_m->destroyT();
634  }
637  deletedParticles_m = true;
638  }
639 
640  size_t totalNum = itsBunch_m->getTotalNum();
642  numParticlesInSimulation_m = totalNum;
643  }
644 
645  if (ne > 0) {
646  msg << level1 << "* Deleted " << ne << " particles, "
647  << "remaining " << totalParticlesInSimulation_m << " particles" << endl;
648  }
649 }
650 
652  bool hasWake = false;
653  WakeFunction *wfInstance;
654 
655  Inform msg("ParallelTTracker ", *gmsg);
656 
657  const unsigned int localNum = itsBunch_m->getLocalNum();
658 
659  IndexMap::value_t::const_iterator it = elements.begin();
660  const IndexMap::value_t::const_iterator end = elements.end();
661 
662  for (; it != end; ++ it) {
663  if ((*it)->hasWake() && !hasWake) {
665 
666  hasWake = true;
667 
668  if ((*it)->getWake()->getType() == "CSRWakeFunction" ||
669  (*it)->getWake()->getType() == "CSRIGFWakeFunction") {
670  if ((*it)->getType() == ElementBase::RBEND ||
671  (*it)->getType() == ElementBase::SBEND) {
672  wfInstance = (*it)->getWake();
673  wakeFunction_m = wfInstance;
674  } else {
675  wfInstance = wakeFunction_m;
676  }
677  } else {
678  wfInstance = (*it)->getWake();
679  }
680 
681  if (!wfInstance) {
682  throw OpalException("ParallelTTracker::computeWakefield",
683  "empty wake function");
684  }
685 
686  if (!wakeStatus_m) {
687  msg << level2 << "============== START WAKE CALCULATION =============" << endl;
688  wfInstance->initialize((*it).get());
689  wakeStatus_m = true;
690  }
691 
693 
694  Quaternion alignment = getQuaternion(itsBunch_m->get_pmean(), Vector_t(0, 0, 1));
695  CoordinateSystemTrafo referenceToBeamCSTrafo(Vector_t(0.0), alignment);
696  CoordinateSystemTrafo beamToReferenceCSTrafo = referenceToBeamCSTrafo.inverted();
697 
698  for (unsigned int i = 0; i < localNum; ++ i) {
699  itsBunch_m->R[i] = referenceToBeamCSTrafo.transformTo(itsBunch_m->R[i]);
700  itsBunch_m->P[i] = referenceToBeamCSTrafo.rotateTo(itsBunch_m->P[i]);
701  itsBunch_m->Ef[i] = referenceToBeamCSTrafo.rotateTo(itsBunch_m->Ef[i]);
702  }
703 
704  wfInstance->apply(itsBunch_m);
705 
706  for (unsigned int i = 0; i < localNum; ++ i) {
707  itsBunch_m->R[i] = beamToReferenceCSTrafo.transformTo(itsBunch_m->R[i]);
708  itsBunch_m->P[i] = beamToReferenceCSTrafo.rotateTo(itsBunch_m->P[i]);
709  itsBunch_m->Ef[i] = beamToReferenceCSTrafo.rotateTo(itsBunch_m->Ef[i]);
710  }
711 
713  }
714  }
715 
716  if (wakeStatus_m && !hasWake) {
717  msg << level2 << "=============== END WAKE CALCULATION ==============" << endl;
718  wakeStatus_m = false;
719  }
720 }
721 
723  Inform msg("ParallelTTracker ", *gmsg);
724  std::set<IndexMap::value_t::value_type> elementsWithParticleMatterInteraction;
725  std::set<ParticleMatterInteractionHandler*> particleMatterinteractionHandlers;
726  std::pair<double, double> currentRange(0.0, 0.0);
727 
728  while (elements.size() > 0) {
729  auto it = elements.begin();
730  if ((*it)->hasParticleMatterInteraction()) {
731  elementsWithParticleMatterInteraction.insert(*it);
732  particleMatterinteractionHandlers.insert((*it)->getParticleMatterInteraction());
733 
734  std::pair<double, double> range = oth.getRange(*it, pathLength_m);
735  currentRange.first = std::min(currentRange.first, range.first);
736  currentRange.second = std::max(currentRange.second, range.second);
737 
738  IndexMap::value_t touching = oth.getTouchingElements(range);
739  elements.insert(touching.begin(), touching.end());
740  }
741 
742  elements.erase(it);
743  }
744 
745  if (elementsWithParticleMatterInteraction.size() > 0) {
746  std::set<ParticleMatterInteractionHandler*> oldSPHandlers;
747  std::vector<ParticleMatterInteractionHandler*> leftBehindSPHandlers, newSPHandlers;
749  oldSPHandlers.insert(it);
750  }
751 
752  leftBehindSPHandlers.resize(std::max(oldSPHandlers.size(),
753  particleMatterinteractionHandlers.size()));
754  auto last = std::set_difference(oldSPHandlers.begin(), oldSPHandlers.end(),
755  particleMatterinteractionHandlers.begin(), particleMatterinteractionHandlers.end(),
756  leftBehindSPHandlers.begin());
757  leftBehindSPHandlers.resize(last - leftBehindSPHandlers.begin());
758 
759  for (auto it: leftBehindSPHandlers) {
760  if (!it->stillActive()) {
761  activeParticleMatterInteractionHandlers_m.erase(it);
762  }
763  }
764 
765  newSPHandlers.resize(std::max(oldSPHandlers.size(),
766  elementsWithParticleMatterInteraction.size()));
767  last = std::set_difference(particleMatterinteractionHandlers.begin(), particleMatterinteractionHandlers.end(),
768  oldSPHandlers.begin(), oldSPHandlers.end(),
769  newSPHandlers.begin());
770  newSPHandlers.resize(last - newSPHandlers.begin());
771 
772  for (auto it: newSPHandlers) {
773  activeParticleMatterInteractionHandlers_m.insert(it);
774  }
775 
777  msg << level2 << "============== START PARTICLE MATTER INTERACTION CALCULATION =============" << endl;
778  particleMatterStatus_m = true;
779  }
780  }
781 
783  do {
785  //check if more than one degrader has particles
786  ParticleMatterInteractionHandler* onlyDegraderWithParticles = NULL;
787  int degradersWithParticlesCount = 0;
789  it->setFlagAllParticlesIn(false);
790  if (it->getParticlesInMat() > 0) {
791  onlyDegraderWithParticles = it;
792  ++ degradersWithParticlesCount;
793  }
794  }
795 
796  //if max particles per node is 2, and only one degrader has particles set
797  //AllParticlesIn for this degrader to true
798  unsigned int localNum = itsBunch_m->getLocalNum();
799  unsigned int totalNum = 0;
800  reduce(localNum, totalNum, OpAddAssign());
801  bool allParticlesInMat = (totalNum == 0 &&
802  degradersWithParticlesCount == 1);
803 
804  if (allParticlesInMat) {
805  onlyDegraderWithParticles->setFlagAllParticlesIn(true);
806  msg << "All particles in degrader" << endl;
807  }
808 
809  auto boundingSphere = itsBunch_m->getLocalBoundingSphere();
810  unsigned int rediffusedParticles = 0;
811  unsigned int numEnteredParticles = 0;
812  for (auto it: activeParticleMatterInteractionHandlers_m) {
813  ElementBase* element = it->getElement();
814  CoordinateSystemTrafo refToLocalCSTrafo = (element->getMisalignment() *
816  CoordinateSystemTrafo localToRefCSTrafo = refToLocalCSTrafo.inverted();
817 
818  const unsigned int localNum = itsBunch_m->getLocalNum();
819  for (unsigned int i = 0; i < localNum; ++i) {
820  itsBunch_m->R[i] = refToLocalCSTrafo.transformTo(itsBunch_m->R[i]);
821  itsBunch_m->P[i] = refToLocalCSTrafo.rotateTo(itsBunch_m->P[i]);
822  }
823  boundingSphere.first = refToLocalCSTrafo.transformTo(boundingSphere.first);
824 
825  it->apply(itsBunch_m, boundingSphere, totalParticlesInSimulation_m);
826  it->print(msg);
827 
828  boundingSphere.first = localToRefCSTrafo.transformTo(boundingSphere.first);
829 
830  const unsigned int newLocalNum = itsBunch_m->getLocalNum();
831  for (unsigned int i = 0; i < newLocalNum; ++i) {
832  itsBunch_m->R[i] = localToRefCSTrafo.transformTo(itsBunch_m->R[i]);
833  itsBunch_m->P[i] = localToRefCSTrafo.rotateTo(itsBunch_m->P[i]);
834  }
835 
836  rediffusedParticles += it->getRediffused();
837  numEnteredParticles += it->getNumEntered();
838 
839  //if all particles were in material update time to time in degrader
840  if (it->getFlagAllParticlesIn()) {
841  double timeDifference = it->getTime() - itsBunch_m->getdT() - itsBunch_m->getT();
842  if (timeDifference > 0.0) {
843  const unsigned int numSteps = std::ceil(timeDifference / itsBunch_m->getdT());
844  const double origdT = itsBunch_m->getdT();
845  itsBunch_m->setdT(timeDifference / numSteps);
846  BorisPusher pusher(itsReference);
847  for (unsigned int i = 0; i < numSteps; ++ i) {
848  updateReference(pusher);
849  }
850  itsBunch_m->setdT(origdT);
851  }
852  itsBunch_m->setT(it->getTime() - itsBunch_m->getdT());
853  }
854  }
855 
856  //perform boundp only if there are particles in the bunch, or there are particles
857  //comming out of the degrader
858  if (numEnteredParticles > 0 || rediffusedParticles > 0) {
859  totalNum -= (numEnteredParticles + rediffusedParticles);
860  if (totalNum > minBinEmitted_m && itsBunch_m->hasFieldSolver()) {
861  itsBunch_m->boundp();
862  } else {
864  }
865  }
866 
867  //if bunch has no particles update time to time in degrader
868  if (itsBunch_m->getTotalNum() == 0)
870 
871  } while (itsBunch_m->getTotalNum() == 0);
872 
873 
875  msg << level2 << "============== END PARTICLE MATTER INTERACTION CALCULATION =============" << endl;
876  particleMatterStatus_m = false;
877  }
878  }
879 }
880 
882  size_t particles_or_bins = std::max(minBinEmitted_m, size_t(1000));
883  if (itsBunch_m->hasFieldSolver() && numParticlesInSimulation_m > particles_or_bins) {
884 
885  INFOMSG("*****************************************************************" << endl);
886  INFOMSG("do repartition because of repartFreq_m" << endl);
887  INFOMSG("*****************************************************************" << endl);
890  Ippl::Comm->barrier();
892  INFOMSG("*****************************************************************" << endl);
893  INFOMSG("do repartition done" << endl);
894  INFOMSG("*****************************************************************" << endl);
895  }
896 }
897 
898 void ParallelTTracker::dumpStats(long long step, bool psDump, bool statDump) {
899  OPALTimer::Timer myt2;
900  Inform msg("ParallelTTracker ", *gmsg);
901 
902  if (itsBunch_m->getGlobalTrackStep() % 1000 + 1 == 1000) {
903  msg << level1;
904  } else if (itsBunch_m->getGlobalTrackStep() % 100 + 1 == 100) {
905  msg << level2;
906  } else {
907  msg << level3;
908  }
909 
910  if (numParticlesInSimulation_m == 0) {
911  msg << myt2.time() << " "
912  << "Step " << std::setw(6) << itsBunch_m->getGlobalTrackStep() << "; "
913  << " -- no emission yet -- "
914  << "t= " << Util::getTimeString(itsBunch_m->getT())
915  << endl;
916  return;
917  }
918 
920  size_t totalParticles_f = numParticlesInSimulation_m;
921  if (totalParticles_f <= minBinEmitted_m) {
922  msg << myt2.time() << " "
923  << "Step " << std::setw(6) << itsBunch_m->getGlobalTrackStep() << "; "
924  << "only " << std::setw(4) << totalParticles_f << " particles emitted; "
925  << "t= " << Util::getTimeString(itsBunch_m->getT()) << ", "
927  << endl;
929  throw OpalException("ParallelTTracker::dumpStats()",
930  "there seems to be something wrong with the position of the bunch!");
931  } else {
932 
933  msg << myt2.time() << " "
934  << "Step " << std::setw(6) << itsBunch_m->getGlobalTrackStep() << " "
935  << "at " << Util::getLengthString(pathLength_m) << ", "
936  << "t= " << Util::getTimeString(itsBunch_m->getT()) << ", "
938  << endl;
939 
940  writePhaseSpace(step, psDump, statDump);
941  }
942 }
943 
944 
946  Inform msg("ParallelTTracker ", *gmsg);
947 
949  RealVariable *ar = dynamic_cast<RealVariable *>(OpalData::getInstance()->find("MINBINEMITTED"));
950  if (ar)
951  minBinEmitted_m = static_cast<size_t>(ar->getReal());
952  msg << level2 << "MINBINEMITTED " << minBinEmitted_m << endl;
953 
955  RealVariable *br = dynamic_cast<RealVariable *>(OpalData::getInstance()->find("MINSTEPFORREBIN"));
956  if (br)
957  minStepforReBin_m = static_cast<int>(br->getReal());
958  msg << level2 << "MINSTEPFORREBIN " << minStepforReBin_m << endl;
959 
960  // there is no point to do repartitioning with one node
961  if (Ippl::getNodes() == 1) {
963  } else {
965  RealVariable *rep = dynamic_cast<RealVariable *>(OpalData::getInstance()->find("REPARTFREQ"));
966  if (rep)
967  repartFreq_m = static_cast<int>(rep->getReal());
968  msg << level2 << "REPARTFREQ " << repartFreq_m << endl;
969  }
970 }
971 
972 
975  return globalEOL_m;
976 }
977 
979  const unsigned int localNum = itsBunch_m->getLocalNum();
980  for (unsigned int i = 0; i < localNum; ++i) {
981  itsBunch_m->dt[i] = itsBunch_m->getdT();
982  }
983 }
984 
986  Inform msg("ParallelTTracker ", *gmsg);
987 
988  if (!itsBunch_m->doEmission()) return;
989 
991  msg << level2 << "Do emission for " << itsBunch_m->getTEmission() << " [s] using "
992  << itsBunch_m->getNumberOfEnergyBins() << " energy bins " << endl
993  << "Change dT from " << itsBunch_m->getdT() << " [s] to "
994  << itsBunch_m->getEmissionDeltaT() << " [s] during emission " << endl;;
995 
996 }
997 
998 void ParallelTTracker::writePhaseSpace(const long long step, bool psDump, bool statDump) {
999  extern Inform *gmsg;
1000  Inform msg("OPAL ", *gmsg);
1001  Vector_t externalE, externalB;
1002  Vector_t FDext[2]; // FDext = {BHead, EHead, BRef, ERef, BTail, ETail}.
1003 
1004  // Sample fields at (xmin, ymin, zmin), (xmax, ymax, zmax) and the centroid location. We
1005  // are sampling the electric and magnetic fields at the back, front and
1006  // center of the beam.
1007  Vector_t rmin, rmax;
1008  itsBunch_m->get_bounds(rmin, rmax);
1009 
1010  if (psDump || statDump) {
1011  externalB = Vector_t(0.0);
1012  externalE = Vector_t(0.0);
1015  itsBunch_m->getT() - 0.5 * itsBunch_m->getdT(),
1016  externalE,
1017  externalB);
1018  FDext[0] = itsBunch_m->toLabTrafo_m.rotateFrom(externalB);
1019  FDext[1] = itsBunch_m->toLabTrafo_m.rotateFrom(externalE * 1e-6);
1020  }
1021 
1022  if (statDump) {
1023  std::vector<std::pair<std::string, unsigned int> > collimatorLosses;
1025  if (collimators.size() != 0) {
1026  for (FieldList::iterator it = collimators.begin(); it != collimators.end(); ++ it) {
1027  FlexibleCollimator* coll = static_cast<FlexibleCollimator*>(it->getElement().get());
1028  std::string name = coll->getName();
1029  unsigned int losses = coll->getLosses();
1030  collimatorLosses.push_back(std::make_pair(name, losses));
1031  }
1032  std::sort(collimatorLosses.begin(), collimatorLosses.end(),
1033  [](const std::pair<std::string, unsigned int>& a, const std::pair<std::string, unsigned int>& b) ->bool {
1034  return a.first < b.first;
1035  });
1036  std::vector<unsigned int> bareLosses(collimatorLosses.size(),0);
1037  for (size_t i = 0; i < collimatorLosses.size(); ++ i){
1038  bareLosses[i] = collimatorLosses[i].second;
1039  }
1040 
1041  reduce(&bareLosses[0], &bareLosses[0] + bareLosses.size(), &bareLosses[0], OpAddAssign());
1042 
1043  for (size_t i = 0; i < collimatorLosses.size(); ++ i){
1044  collimatorLosses[i].second = bareLosses[i];
1045  }
1046  }
1047  // Write statistical data.
1048  itsDataSink_m->dumpSDDS(itsBunch_m, FDext, collimatorLosses);
1049 
1050  msg << level3 << "* Wrote beam statistics." << endl;
1051  }
1052 
1053  if (psDump && (itsBunch_m->getTotalNum() > 0)) {
1054  // Write fields to .h5 file.
1055  const size_t localNum = itsBunch_m->getLocalNum();
1056  double distToLastStop = stepSizes_m.getFinalZStop() - pathLength_m;
1058  Vector_t driftPerTimeStep = itsBunch_m->getdT() * Physics::c * itsBunch_m->toLabTrafo_m.rotateFrom(beta);
1059  bool driftToCorrectPosition = std::abs(distToLastStop) < 0.5 * euclidean_norm(driftPerTimeStep);
1060  Ppos_t stashedR;
1061  Vector_t stashedRefPartR;
1062 
1063  if (driftToCorrectPosition) {
1064  const double tau = distToLastStop / euclidean_norm(driftPerTimeStep) * itsBunch_m->getdT();
1065  if (localNum > 0) {
1066  stashedR.create(localNum);
1067  stashedR = itsBunch_m->R;
1068  stashedRefPartR = itsBunch_m->RefPartR_m;
1069 
1070  for (size_t i = 0; i < localNum; ++ i) {
1071  itsBunch_m->R[i] += tau * (Physics::c * itsBunch_m->P[i] / Util::getGamma(itsBunch_m->P[i]) -
1072  driftPerTimeStep / itsBunch_m->getdT());
1073  }
1074  }
1075 
1076  driftPerTimeStep = itsBunch_m->toLabTrafo_m.rotateTo(driftPerTimeStep);
1077  itsBunch_m->RefPartR_m = itsBunch_m->RefPartR_m + tau * driftPerTimeStep / itsBunch_m->getdT();
1078  CoordinateSystemTrafo update(tau * driftPerTimeStep / itsBunch_m->getdT(),
1079  Quaternion(1.0, 0.0, 0.0, 0.0));
1081 
1083 
1085  }
1086  if (!statDump && !driftToCorrectPosition) itsBunch_m->calcBeamParameters();
1087 
1088  msg << *itsBunch_m << endl;
1089  itsDataSink_m->dumpH5(itsBunch_m, FDext);
1090 
1091  if (driftToCorrectPosition) {
1092  if (localNum > 0) {
1093  itsBunch_m->R = stashedR;
1094  stashedR.destroy(localNum, 0);
1095  }
1096 
1097  itsBunch_m->RefPartR_m = stashedRefPartR;
1099 
1101  }
1102 
1103  msg << level2 << "* Wrote beam phase space." << endl;
1104  }
1105 }
1106 
1108  updateReferenceParticle(pusher);
1110 }
1111 
1113  const double direction = back_track? -1: 1;
1114  const double dt = direction * std::min(itsBunch_m->getT(),
1115  direction * itsBunch_m->getdT());
1116  const double scaleFactor = Physics::c * dt;
1117  Vector_t Ef(0.0), Bf(0.0);
1118 
1119  itsBunch_m->RefPartR_m /= scaleFactor;
1121  itsBunch_m->RefPartR_m *= scaleFactor;
1122 
1124  IndexMap::value_t::const_iterator it = elements.begin();
1125  const IndexMap::value_t::const_iterator end = elements.end();
1126 
1127  for (; it != end; ++ it) {
1128  const CoordinateSystemTrafo &refToLocalCSTrafo = itsOpalBeamline_m.getCSTrafoLab2Local((*it));
1129 
1130  Vector_t localR = refToLocalCSTrafo.transformTo(itsBunch_m->RefPartR_m);
1131  Vector_t localP = refToLocalCSTrafo.rotateTo(itsBunch_m->RefPartP_m);
1132  Vector_t localE(0.0), localB(0.0);
1133 
1134  if ((*it)->applyToReferenceParticle(localR,
1135  localP,
1136  itsBunch_m->getT() - 0.5 * dt,
1137  localE,
1138  localB)) {
1139  *gmsg << level1 << "The reference particle hit an element" << endl;
1140  globalEOL_m = true;
1141  }
1142 
1143  Ef += refToLocalCSTrafo.rotateFrom(localE);
1144  Bf += refToLocalCSTrafo.rotateFrom(localB);
1145  }
1146 
1147  pusher.kick(itsBunch_m->RefPartR_m, itsBunch_m->RefPartP_m, Ef, Bf, dt);
1148 
1149  itsBunch_m->RefPartR_m /= scaleFactor;
1151  itsBunch_m->RefPartR_m *= scaleFactor;
1152 }
1153 
1155  const unsigned int localNum = itsBunch_m->getLocalNum();
1156  for (unsigned int i = 0; i < localNum; ++i) {
1157  itsBunch_m->R[i] = trafo.transformTo(itsBunch_m->R[i]);
1158  itsBunch_m->P[i] = trafo.rotateTo(itsBunch_m->P[i]);
1159  itsBunch_m->Ef[i] = trafo.rotateTo(itsBunch_m->Ef[i]);
1160  itsBunch_m->Bf[i] = trafo.rotateTo(itsBunch_m->Bf[i]);
1161  }
1162 }
1163 
1167 
1169 
1170  CoordinateSystemTrafo update(R, getQuaternion(P, Vector_t(0, 0, 1)));
1171 
1172  transformBunch(update);
1173 
1175 }
1176 
1177 void ParallelTTracker::applyFractionalStep(const BorisPusher &pusher, double tau) {
1178  double t = itsBunch_m->getT();
1179  t += tau;
1180  itsBunch_m->setT(t);
1181 
1182  // the push method below pushes for half a time step. Hence the ref particle
1183  // should be pushed for 2 * tau
1184  itsBunch_m->RefPartR_m /= (Physics::c * 2 * tau);
1186  itsBunch_m->RefPartR_m *= (Physics::c * 2 * tau);
1187 
1189 
1192  CoordinateSystemTrafo update(R, getQuaternion(P, Vector_t(0, 0, 1)));
1194 }
1195 
1197 
1198  StepSizeConfig stepSizesCopy(stepSizes_m);
1199  if (back_track) {
1200  stepSizesCopy.shiftZStopLeft(zstart_m);
1201  }
1202 
1203  double t = 0.0;
1204  itsBunch_m->setT(t);
1205 
1206  dtCurrentTrack_m = stepSizesCopy.getdT();
1207  selectDT();
1208 
1211  double gamma = 0.1 / itsBunch_m->getM() + 1.0;
1212  double beta = sqrt(1.0 - 1.0 / std::pow(gamma, 2));
1213  itsBunch_m->RefPartP_m = itsBunch_m->toLabTrafo_m.rotateTo(beta * gamma * Vector_t(0, 0, 1));
1214  }
1215 
1216  while (true) {
1217  autophaseCavities(pusher);
1218 
1219  t += itsBunch_m->getdT();
1220  itsBunch_m->setT(t);
1221 
1222  Vector_t oldR = itsBunch_m->RefPartR_m;
1223  updateReferenceParticle(pusher);
1225 
1227 
1228  if (pathLength_m > stepSizesCopy.getZStop()) {
1229  ++ stepSizesCopy;
1230 
1231  if (stepSizesCopy.reachedEnd()) {
1232  -- stepSizesCopy;
1233  double tau = (stepSizesCopy.getZStop() - pathLength_m) / speed;
1234  applyFractionalStep(pusher, tau);
1235 
1236  break;
1237  }
1238 
1239  dtCurrentTrack_m = stepSizesCopy.getdT();
1240  selectDT();
1241  }
1242 
1243  if (std::abs(pathLength_m - zstart_m) <= 0.5 * itsBunch_m->getdT() * speed) {
1244  double tau = (zstart_m - pathLength_m) / speed;
1245  applyFractionalStep(pusher, tau);
1246 
1247  break;
1248  }
1249  }
1250 
1251  changeDT();
1252 }
1253 
1255 
1256  double t = itsBunch_m->getT();
1258  pusher.push(nextR, itsBunch_m->RefPartP_m, itsBunch_m->getdT());
1259  nextR *= Physics::c * itsBunch_m->getdT();
1260 
1261  auto elementSet = itsOpalBeamline_m.getElements(nextR);
1262  for (auto element: elementSet) {
1263  if (element->getType() == ElementBase::TRAVELINGWAVE) {
1264  const TravelingWave *TWelement = static_cast<const TravelingWave *>(element.get());
1265  if (!TWelement->getAutophaseVeto()) {
1266  CavityAutophaser ap(itsReference, element);
1269  t, itsBunch_m->getdT());
1270  }
1271 
1272  } else if (element->getType() == ElementBase::RFCAVITY) {
1273  const RFCavity *RFelement = static_cast<const RFCavity *>(element.get());
1274  if (!RFelement->getAutophaseVeto()) {
1275  CavityAutophaser ap(itsReference, element);
1278  t, itsBunch_m->getdT());
1279  }
1280  }
1281  }
1282 }
1283 
1285  unsigned int who;
1286  unsigned int whom;
1287  unsigned int howMany;
1288 };
1289 
1291  const int numNodes = Ippl::getNodes();
1292  if (itsBunch_m->hasFieldSolver() || numNodes == 1) return;
1293 
1294  long onAverage = itsBunch_m->getTotalNum() / Ippl::getNodes();
1295  if (itsBunch_m->getTotalNum() % Ippl::getNodes() > 0.5 * Ippl::getNodes())
1296  ++ onAverage;
1297 
1298  std::vector<long> localParticles(numNodes, 0);
1299  localParticles[Ippl::myNode()] = itsBunch_m->getLocalNum() - onAverage;
1300  allreduce(&(localParticles[0]),
1301  numNodes,
1302  std::plus<long>());
1303 
1304  std::vector<DistributionInfo> send;
1305  std::vector<DistributionInfo> receive;
1306 
1307  for (int i = 0; i < Ippl::getNodes(); ++ i) {
1308  if (localParticles[i] <= 0) continue;
1309 
1310  for (int j = 0; j < Ippl::getNodes(); ++ j) {
1311  if (j == i || localParticles[j] >= 0) continue;
1312 
1313  long numParts = std::min(localParticles[i], -localParticles[j]);
1314  localParticles[i] -= numParts;
1315  localParticles[j] += numParts;
1316 
1317  if (i == Ippl::myNode() || j == Ippl::myNode()) {
1318  DistributionInfo request;
1319  request.who = i;
1320  request.whom = j;
1321  request.howMany = numParts;
1322 
1323  if (i == Ippl::myNode()) {
1324  send.push_back(request);
1325  } else {
1326  receive.push_back(request);
1327  }
1328  }
1329 
1330  if (localParticles[i] == 0) break;
1331  }
1332  }
1333 
1334  std::vector<MPI_Request> requests;
1335  const long sizeSingleParticle = 9 * sizeof(double) + sizeof(short) + sizeof(int) + sizeof(PID_t::Return_t);
1336  long idx = itsBunch_m->getLocalNum() - 1;
1338 
1339  std::vector<char> send_msgbuf;
1340 
1341  if (send.size() > 0) {
1342  const char *buffer;
1343 
1344  unsigned int totalSend = 0, startIndex = 0;
1345  for (DistributionInfo &request: send) {
1346  totalSend += request.howMany;
1347  }
1348  send_msgbuf.reserve(totalSend * sizeSingleParticle);
1349 
1350  for (DistributionInfo &request: send) {
1351  size_t sizePrior = send_msgbuf.size();
1352  for (long i = 0; i < request.howMany; ++ i, -- idx) {
1353  buffer = reinterpret_cast<const char*>(&(itsBunch_m->R[idx](0)));
1354  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + 3 * sizeof(double));
1355  buffer = reinterpret_cast<const char*>(&(itsBunch_m->P[idx](0)));
1356  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + 3 * sizeof(double));
1357  buffer = reinterpret_cast<const char*>(&(itsBunch_m->Q[idx]));
1358  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(double));
1359  buffer = reinterpret_cast<const char*>(&(itsBunch_m->M[idx]));
1360  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(double));
1361  buffer = reinterpret_cast<const char*>(&(itsBunch_m->dt[idx]));
1362  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(double));
1363  buffer = reinterpret_cast<const char*>(&(itsBunch_m->PType[idx]));
1364  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(short));
1365  buffer = reinterpret_cast<const char*>(&(itsBunch_m->TriID[idx]));
1366  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(int));
1367  buffer = reinterpret_cast<const char*>(&(itsBunch_m->ID[idx]));
1368  send_msgbuf.insert(send_msgbuf.end(), buffer, buffer + sizeof(PID_t::Return_t));
1369  }
1370 
1371  size_t sendSizeThis = send_msgbuf.size() - sizePrior;
1372  MPI_Request req = Ippl::Comm->raw_isend(&(send_msgbuf[startIndex]),
1373  sendSizeThis,
1374  request.whom,
1375  tag);
1376 
1377  requests.push_back(req);
1378 
1379  startIndex += sendSizeThis;
1380  }
1381 
1382  itsBunch_m->destroy(totalSend, idx + 1, true);
1383  }
1384 
1385  for (unsigned int i = 0; i < receive.size(); ++ i) {
1386  int node = Communicate::COMM_ANY_NODE;
1387  char *recvbuf;
1388  const int bufsize = Ippl::Comm->raw_probe_receive(recvbuf, node, tag);
1389 
1390  int j = 0;
1391 
1392  while (j < bufsize) {
1393  ++ idx;
1394  itsBunch_m->create(1);
1395  {
1396  const double *buffer = reinterpret_cast<const double*>(recvbuf + j);
1397  itsBunch_m->R[idx] = Vector_t(buffer[0], buffer[1], buffer[2]);
1398  itsBunch_m->P[idx] = Vector_t(buffer[3], buffer[4], buffer[5]);
1399  itsBunch_m->Q[idx] = buffer[6];
1400  itsBunch_m->M[idx] = buffer[7];
1401  itsBunch_m->dt[idx] = buffer[8];
1402  }
1403  j += 9 * sizeof(double);
1404 
1405  {
1406  const short *buffer = reinterpret_cast<const short*>(recvbuf + j);
1407  itsBunch_m->PType[idx] = buffer[0];
1408  }
1409  j += sizeof(short);
1410 
1411  {
1412  const int *buffer = reinterpret_cast<const int*>(recvbuf + j);
1413  itsBunch_m->TriID[idx] = buffer[0];
1414  }
1415  j += sizeof(int);
1416 
1417  {
1418  const PID_t::Return_t *buffer = reinterpret_cast<const PID_t::Return_t*>(recvbuf + j);
1419  itsBunch_m->ID[idx] = buffer[0];
1420  }
1421  j += sizeof(PID_t::Return_t);
1422  }
1423 
1424  delete[] recvbuf;
1425  }
1426 
1427  if (requests.size() > 0) {
1428  MPI_Waitall(requests.size(), &(requests[0]), MPI_STATUSES_IGNORE);
1429  }
1430 }
1431 
1432 // vi: set et ts=4 sw=4 sts=4:
1433 // Local Variables:
1434 // mode:c++
1435 // c-basic-offset: 4
1436 // indent-tabs-mode:nil
1437 // End:
static int getNodes()
Definition: IpplInfo.cpp:773
ParticleAttrib< Vector_t > P
void switchOffUnitlessPositions(bool use_dt_per_particle=false)
elements
Definition: IndexMap.cpp:141
double getMinTimeStep() const
void getMisalignment(double &x, double &y, double &s) const
void pushParticles(const BorisPusher &pusher)
double getT() const
WakeFunction * wakeFunction_m
std::vector< MaxPhasesT >::iterator getLastMaxPhases()
Definition: OpalData.cpp:447
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
void setPriorTrack(const bool &value=true)
true if in follow-up track
Definition: OpalData.cpp:332
void findStartPosition(const BorisPusher &pusher)
int getLastEmittedEnergyBin()
constexpr double e
The value of .
Definition: Physics.h:40
void destroy(size_t M, size_t I, bool doNow=false)
ParticleAttrib< Vector_t > Ef
Vector_t RefPartR_m
Interface for basic beam line object.
Definition: ElementBase.h:128
void setOpenMode(OPENMODE openMode)
Definition: OpalData.cpp:377
virtual void create(size_t)
bool reachedEnd() const
void prepareSections()
bool getRelativeFlag() const
Definition: TBeamline.h:470
void create(size_t M)
void transformBunch(const CoordinateSystemTrafo &trafo)
Vector_t rotateTo(const Vector_t &r) const
const Beamline & itsBeamline_m
Definition: Tracker.h:145
unsigned long totalParticlesInSimulation_m
virtual void execute()
Apply the algorithm to the top-level beamline.
virtual void visitBeamline(const Beamline &)
Apply the algorithm to a beam line.
void setGlobalMeanR(Vector_t globalMeanR)
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:123
The base class for all OPAL exceptions.
Definition: OpalException.h:28
virtual double getReal() const
Return value.
PETE_TBTree< FnCopysign, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > copysign(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
ParticleAttrib< double > Q
void compute3DLattice()
std::list< ClassicField > FieldList
Definition: ClassicField.h:47
PETE_TUTree< FnCeil, typename T::PETE_Expr_t > ceil(const PETE_Expr< T > &l)
Definition: PETE.h:811
IpplTimings::TimerRef timeIntegrationTimer1_m
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 updateReferenceParticle(const BorisPusher &pusher)
int getNumberOfEnergyBins()
Particle reference data.
Definition: PartData.h:38
double getdT() const
The REAL VARIABLE definition.
Definition: RealVariable.h:26
Inform * gmsg
Definition: Main.cpp:21
virtual void do_binaryRepart()
double getM() const
void sortAscendingZStop()
double getPhaseAtMaxEnergy(const Vector_t &R, const Vector_t &P, double t, double dt)
static bool DKSEnabled
Definition: IpplInfo.h:285
Abstract collimator.
IpplTimings::TimerRef fieldEvaluationTimer_m
unsigned int emissionSteps_m
int minStepForRebin
The number of steps into the simulation before the bins are squashed into a single bin...
Definition: Options.cpp:66
void barrier(void)
double get_meanKineticEnergy() const
static int myNode()
Definition: IpplInfo.cpp:794
virtual const std::string & getName() const
Get element name.
Definition: ElementBase.cpp:95
Vector_t transformTo(const Vector_t &r) const
unsigned int repartFreq_m
Quaternion getQuaternion(Vector_t u, Vector_t ref)
Definition: Quaternion.cpp:34
virtual bool getAutophaseVeto() const
Definition: RFCavity.h:458
Vector_t get_pmean_Distribution() const
void setdT(double dt)
PETE_TBTree< OpNE, Index::PETE_Expr_t, PETE_Scalar< double > > ne(const Index &idx, double x)
Definition: IndexInlines.h:357
unsigned long long getNumStepsFinestResolution() const
size_t numParticlesInSimulation_m
IpplTimings::TimerRef BinRepartTimer_m
unsigned long getFieldAt(const unsigned int &, const Vector_t &, const long &, const double &, Vector_t &, Vector_t &)
void updateReference(const BorisPusher &pusher)
std::string getLengthString(double spos, unsigned int precision=3)
Definition: Util.h:63
unsigned long long getMaxSteps() const
Vector_t getOrigin3D() const
Definition: TBeamline.h:450
void kickParticles(const BorisPusher &pusher)
size_t getTotalNum() const
void computeParticleMatterInteraction(IndexMap::value_t elements, OrbitThreader &oth)
double getFinalZStop() const
virtual void resetInterpolationCache(bool clearCache=false)
int next_tag(int t, int s=1000)
Definition: TagMaker.h:43
IndexMap::value_t getTouchingElements(const std::pair< double, double > &range)
bool containsSource()
Definition: OpalBeamline.h:224
int statDumpFreq
The frequency to dump statistical values, e.e. dump data when stepstatDumpFreq==0.
Definition: Options.cpp:54
double zstart_m
where to start
Vector_t transformToLocalCS(const std::shared_ptr< Component > &comp, const Vector_t &r) const
Definition: OpalBeamline.h:185
Inform & level2(Inform &inf)
Definition: Inform.cpp:46
void updateRFElement(std::string elName, double maxPhi)
OpalBeamline itsOpalBeamline_m
virtual void apply(PartBunchBase< double, 3 > *bunch)=0
void set_sPos(double s)
void computeWakefield(IndexMap::value_t &elements)
bool reduce(Communicate &, InputIterator, InputIterator, OutputIterator, const ReduceOp &, bool *IncludeVal=0)
Definition: GlobalComm.hpp:55
static OpalData * getInstance()
Definition: OpalData.cpp:209
void computeExternalFields(OrbitThreader &oth)
void writePhaseSpace(const long long step, bool psDump, bool statDump)
void activateElements()
void calcBeamParameters()
Vector_t get_pmean() const
Template class for beam lines.
Definition: TBeamline.h:40
double getGamma(Vector_t p)
Definition: Util.h:24
ParticleAttrib< short > PType
ParticleIndex_t & ID
double getdT() const
double getGlobalPhaseShift()
units: (sec)
Definition: OpalData.cpp:497
#define P_LAYOUT_CYCLE
Definition: Tags.h:86
void applyFractionalStep(const BorisPusher &pusher, double tau)
void computeSpaceChargeFields(unsigned long long step)
virtual void accept(BeamlineVisitor &visitor) const =0
Apply visitor.
IpplTimings::TimerRef WakeFieldTimer_m
void shiftZStopLeft(double back)
bool getIfBeamEmitting()
int repartFreq
The frequency to do particles repartition for better load balance between nodes.
Definition: Options.cpp:60
std::string getTimeString(double time, unsigned int precision=3)
Definition: Util.h:40
void changeDT(bool backTrack=false)
ParticleAttrib< int > TriID
IpplTimings::TimerRef timeIntegrationTimer2_m
void allreduce(const T *input, T *output, int count, Op op)
Definition: GlobalComm.hpp:510
constexpr double c
The velocity of light in m/s.
Definition: Physics.h:52
static void startTimer(TimerRef t)
Definition: IpplTimings.h:187
void switchElementsOff(const double &, ElementBase::ElementType eltype=ElementBase::ANY)
void kick(const Vector_t &R, Vector_t &P, const Vector_t &Ef, const Vector_t &Bf, const double &dt) const
Definition: BorisPusher.h:48
CoordinateSystemTrafo toLabTrafo_m
Interface for RF cavity.
Definition: TravelingWave.h:37
StepSizeConfig stepSizes_m
#define INFOMSG(msg)
Definition: IpplInfo.h:397
double getTEmission()
virtual void iterate(BeamlineVisitor &, bool r2l) const
Apply visitor to all elements of the line.
Definition: TBeamline.h:215
IndexMap::value_t query(IndexMap::key_t::first_type step, IndexMap::key_t::second_type length)
Definition: OrbitThreader.h:94
Vector_t RefPartP_m
void dumpH5(PartBunchBase< double, 3 > *beam, Vector_t FDext[]) const
Definition: DataSink.cpp:69
FieldList getElementByType(ElementBase::ElementType)
void switchToUnitlessPositions(bool use_dt_per_particle=false)
Tps< T > pow(const Tps< T > &x, int y)
Integer power.
Definition: TpsMath.h:76
virtual void computeSelfFields()=0
Vektor< double, 3 > Vector_t
Definition: Vektor.h:6
void swap(OpalBeamline &rhs)
void emitParticles(long long step)
An abstract sequence of beam line components.
Definition: Beamline.h:37
const PartData itsReference
The reference information.
size_t getLocalNum() const
Vector_t get_centroid() const
virtual void destroy(size_t M, size_t I, bool optDestroy=true)
void timeIntegration2(BorisPusher &pusher)
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition: TpsMath.h:91
std::vector< MaxPhasesT >::iterator getFirstMaxPhases()
Definition: OpalData.cpp:443
double beta
Definition: PartData.h:99
int psDumpFreq
The frequency to dump the phase space, i.e.dump data when steppsDumpFreq==0.
Definition: Options.cpp:48
ParticleAttrib< double > dt
double getP() const
The constant reference momentum per particle.
Definition: PartData.h:117
size_t boundp_destroyT()
void push_back(double dt, double zstop, unsigned long numSteps)
virtual void setAutophaseVeto(bool veto=true)
Definition: RFCavity.h:453
virtual void setPhasem(double phase)
Definition: RFCavity.h:408
std::pair< Vector_t, double > getLocalBoundingSphere()
static void writeStatistics()
Definition: Monitor.cpp:209
std::string time() const
Return time.
Definition: Timer.cpp:42
Vector_t rotateToLocalCS(const std::shared_ptr< Component > &comp, const Vector_t &r) const
Definition: OpalBeamline.h:197
std::set< std::shared_ptr< Component > > getElements(const Vector_t &x)
double getZStop() const
Interface for RF cavity.
Definition: RFCavity.h:37
Object * find(const std::string &name)
Find entry.
Definition: OpalData.cpp:618
void setGlobalPhaseShift(double shift)
units: (sec)
Definition: OpalData.cpp:492
void resetIterator()
int minBinEmitted
The number of bins that have to be emitted before the bin are squashed into a single bin...
Definition: Options.cpp:63
std::pair< double, double > getRange(const IndexMap::value_t::value_type &element, double position) const
unsigned long getNumSteps() const
bool weHaveEnergyBins()
CoordinateSystemTrafo getCSTrafoLab2Local(const std::shared_ptr< Component > &comp) const
Definition: OpalBeamline.h:209
T euclidean_norm(const Vector< T > &)
Euclidean norm.
Definition: Vector.h:243
Vector_t rotateFrom(const Vector_t &r) const
void setT(double t)
const std::string name
virtual void initialize(const ElementBase *ref)
Definition: WakeFunction.hh:19
CoordinateSystemTrafo getMisalignment(const std::shared_ptr< Component > &comp) const
Definition: OpalBeamline.h:219
size_t getNumberOfEmissionSteps()
T isnan(T x)
isnan function with adjusted return type
Definition: matheval.hpp:74
std::string::const_iterator iterator_t
Definition: array.cpp:8
ParticleAttrib< int > Bin
CoordinateSystemTrafo inverted() const
void dumpSDDS(PartBunchBase< double, 3 > *beam, Vector_t FDext[], const double &azimuth=-1) const
Definition: DataSink.cpp:99
virtual void boundp()
void calcGammas()
Compute the gammas of all bins.
double getEnergy(Vector_t p, double mass)
Definition: Util.h:29
virtual void setBinCharge(int bin, double q)
Set the charge of one bin to the value of q and all other to zero.
Quaternion getInitialDirection() const
Definition: TBeamline.h:460
ParticleAttrib< Vector_t > Bf
ParticleAttrib< double > M
ParticlePos_t & R
void save3DLattice()
Quaternion conjugate() const
Definition: Quaternion.h:104
Inform & level3(Inform &inf)
Definition: Inform.cpp:47
std::string::iterator iterator
Definition: MSLang.h:16
void timeIntegration1(BorisPusher &pusher)
DataSink * itsDataSink_m
double getEmissionDeltaT()
void get_bounds(Vector_t &rmin, Vector_t &rmax)
StepSizeConfig & advanceToPos(double spos)
static void stopTimer(TimerRef t)
Definition: IpplTimings.h:192
Definition: Inform.h:41
static Communicate * Comm
Definition: IpplInfo.h:93
#define P_SPATIAL_TRANSFER_TAG
Definition: Tags.h:82
Inform & level1(Inform &inf)
Definition: Inform.cpp:45
std::string getEnergyString(double energyInMeV, unsigned int precision=3)
Definition: Util.h:113
std::set< ParticleMatterInteractionHandler * > activeParticleMatterInteractionHandlers_m
T::PETE_Expr_t::PETE_Return_t min(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:95
void push(Vector_t &R, const Vector_t &P, const double &dt) const
Definition: BorisPusher.h:116
PartBunchBase< double, 3 > * itsBunch_m
The bunch of particles to be tracked.
Definition: Tracker.h:174
void autophaseCavities(const BorisPusher &pusher)
void selectDT(bool backTrack=false)
unsigned int getLosses() const
void setInPrepState(bool state)
Definition: OpalData.cpp:324
double gamma
Definition: PartData.h:100
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
void storeCavityInformation()
Write cavity information from H5 file.
Definition: DataSink.cpp:131
std::set< std::shared_ptr< Component > > value_t
Definition: IndexMap.h:18
T isinf(T x)
isinf function with adjusted return type
Definition: matheval.hpp:78
void dumpStats(long long step, bool psDump, bool statDump)
size_t emitParticles(double eZ)
Emit particles in the given bin i.e. copy the particles from the bin structure into the particle cont...
long long getGlobalTrackStep() const
Vector_t transformFrom(const Vector_t &r) const
void save3DInput()
Track particles or bunches.
Definition: Tracker.h:84
CoordinateSystemTrafo getCSTrafoGlobal2Local() const
Definition: ElementBase.h:603