OPALX (Object Oriented Parallel Accelerator Library for Exascal) MINIorX
OPALX
#PartBunch.hpp#
Go to the documentation of this file.
1//
2// Class PartBunch.hpp
3// Particle Bunch.
4// A representation of a particle bunch as a vector of particles.
5//
6// Copyright (c) 2008 - 2020, Paul Scherrer Institut, Villigen PSI, Switzerland
7// All rights reserved
8//
9// This file is part of OPAL.
10//
11// OPAL is free software: you can redistribute it and/or modify
12// it under the terms of the GNU General Public License as published by
13// the Free Software Foundation, either version 3 of the License, or
14// (at your option) any later version.
15//
16// You should have received a copy of the GNU General Public License
17// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
18//
19#ifndef OPAL_PartBunch_HPP
20#define OPAL_PartBunch_HPP
21
22class Distribution;
23
25#include "Algorithms/PartBins.h"
26#include "Algorithms/PartData.h"
27#include "Ippl.h"
28#include "OPALtypes.h"
29#include "Particle/ParticleAttrib.h"
30#include "Particle/ParticleLayout.h"
31#include "Utilities/Util.h"
32
33/*
34
35 Will change
36
37*/
38#include "PoissonSolvers/FFTOpenPoissonSolver.h"
39#include "PoissonSolvers/FFTPeriodicPoissonSolver.h"
40#include "PoissonSolvers/P3MSolver.h"
41
42#include "PoissonSolvers/PoissonCG.h"
43
44using ippl::detail::ConditionalType, ippl::detail::VariantFromConditionalTypes;
45
46template <typename T = double, unsigned Dim = 3>
47using FFTSolver_t = ConditionalType<
48 Dim == 2 || Dim == 3, ippl::FFTPeriodicPoissonSolver<VField_t<T, Dim>, Field_t<Dim>>>;
49
50template <typename T = double, unsigned Dim = 3>
51using P3MSolver_t = ConditionalType<Dim == 3, ippl::P3MSolver<VField_t<T, Dim>, Field_t<Dim>>>;
52
53template <typename T = double, unsigned Dim = 3>
55 ConditionalType<Dim == 3, ippl::FFTOpenPoissonSolver<VField_t<T, Dim>, Field_t<Dim>>>;
56
57template <typename T = double, unsigned Dim = 3>
58using Solver_t =
59 VariantFromConditionalTypes<FFTSolver_t<T, Dim>, P3MSolver_t<T, Dim>, OpenSolver_t<T, Dim>>;
60
61template <class PLayout, typename T, unsigned Dim = 3>
62class PartBunch : public ippl::ParticleBase<PLayout> {
63 using Base = ippl::ParticleBase<PLayout>;
64
65public:
67
68 ParticleAttrib<double> Q; // holds the charge for particle
69 ParticleAttrib<double> M; // holds the mass for particle
70 ParticleAttrib<int> bunchNum; // holds the mass for particle
72 ParticleAttrib<Vector_t<T, Dim>> E; // E field vector f
73 ParticleAttrib<Vector_t<T, Dim>> Ef; // e field vector for gun simulations
74 ParticleAttrib<Vector_t<T, Dim>> Bf; // b field vector for gun simulations
75 ParticleAttrib<int> Bin; // energy bin in which the particle is in, if zero -> rebin
76 ParticleAttrib<double> dt; // holds the dt timestep for particle
77 ParticleAttrib<double> Phi; // the electric potential
78 ParticleAttrib<Vector_t<T, Dim>> Eftmp; // e field vector for gun simulations
79 // ParticleAttrib<ParticleType> PType; // particle names
80
84
87
89
90 // ParticleOrigin refPOrigin_m;
91 // ParticleType refPType_m;
92
96 // Vector_t globalMeanR_m = Vector_t(0.0, 0.0, 0.0);
97 // Quaternion_t globalToLocalQuaternion_m = Quaternion_t(1.0, 0.0, 0.0, 0.0);
99 // Quaternion_t globalToLocalQuaternion_m;
100
104
106
109
111 std::unique_ptr<double[]> bingamma_m;
112
113 // FIXME: this should go into the Bin class!
114 // holds number of emitted particles of the bin
115 // jjyang: opal-cycl use *nBin_m of pbin_m
116 std::unique_ptr<size_t[]> binemitted_m;
117
120
123
126
129
131 std::vector<size_t> bunchTotalNum_m;
132 std::vector<size_t> bunchLocalNum_m;
133
138
141 double dh_m;
142
144
146 double qi_m;
148
149 // unit state of PartBunch
152
154 double centroid_m[2 * Dim];
155
157 // matrix_t moments_m;
158
160 double dt_m;
162 double t_m;
164 double spos_m;
165
168
171
173
174 typedef ippl::BConds<Field<T, Dim>, Dim> bc_type;
176
177 /*
178 Data structure for particle load balance information
179 */
180
182
183 std::unique_ptr<size_t[]> globalPartPerNode_m;
184
186
187 std::array<bool, Dim> isParallel_m;
188
192
193 // flag to tell if we are a DC-beam
196
197 double Q_m;
198
199public:
201 PLayout& pl, Vector_t<double, Dim> hr, Vector_t<double, Dim> rmin,
202 Vector_t<double, Dim> rmax, std::array<bool, Dim> decomp, double Qtot)
203 : ippl::ParticleBase<PLayout>(pl), hr_m(hr), rmin_m(rmin), rmax_m(rmax), Q_m(Qtot) {
204 // register the particle attributes
205 this->addAttribute(Q);
206 this->addAttribute(M);
207 this->addAttribute(P);
208 this->addAttribute(E);
209 this->addAttribute(Phi);
210 this->addAttribute(Ef);
211 this->addAttribute(Eftmp);
212 this->addAttribute(Bf);
213 this->addAttribute(Bin);
214 this->addAttribute(dt);
215 // this->addAttribute(PType);
216 this->addAttribute(bunchNum);
217
219 for (unsigned int i = 0; i < Dim; i++) {
220 isParallel_m[i] = decomp[i];
221 }
222 }
223
224 // PartBunch() = delete;
225
226 PartBunch& operator=(const PartBunch&) = delete;
227
229 }
230
232 updateLayout(fl, mesh);
233 }
234
235 /*
236
237 Physics
238
239 */
240
241 double getCouplingConstant() const {
242 return couplingConstant_m;
243 }
244
245 void setCouplingConstant(double c) {
247 }
248
251 // momentsComputer_m.compute(*this);
252 }
253
254 void calcBeamParametersInitial(); // Calculate initial beam parameters before emission.
255
256 // set the charge per simulation particle
257 void setCharge(double q) {
258 qi_m = q;
259 if (this->getTotalNum() != 0)
260 Q = qi_m;
261 else
262 *ippl::Warn
263 << "Could not set total charge in PartBunch::setCharge based on this->getTotalNum"
264 << endl;
265 }
266
267 // set the charge per simulation particle when total particle number equals 0
268 void setChargeZeroPart(double q) {
269 qi_m = q;
270 }
271
272 // set the mass per simulation particle
273 void setMass(double mass) {
274 massPerParticle_m = mass;
275 if (this->getTotalNum() != 0)
276 M = mass;
277 }
278
279 void setMassZeroPart(double mass) {
280 massPerParticle_m = mass;
281 }
282
284 double getCharge() const {
285 return 0.0; // ADA sum(this->Q);
286 }
287
289 double getChargePerParticle() const {
290 return qi_m;
291 }
292
293 double getMassPerParticle() const {
294 return massPerParticle_m;
295 }
296
298 double getQ() const {
299 return reference->getQ();
300 }
301 double getM() const {
302 return reference->getM();
303 }
304 double getP() const {
305 return reference->getP();
306 }
307
308 double getE() const {
309 return reference->getE();
310 }
311
312 // ParticleOrigin getPOrigin() const;
313 // ParticleType getPType() const;
314
315 double getInitialBeta() const {
316 return reference->getBeta();
317 }
318
319 double getInitialGamma() const {
320 return reference->getGamma();
321 }
322
325
326 void resetQ(double q) {
327 const_cast<PartData*>(reference)->setQ(q);
328 }
329 void resetM(double m) {
330 const_cast<PartData*>(reference)->setM(m);
331 }
332
333 // void setPOrigin(ParticleOrigin);
334 void setPType(const std::string& type) {
335 // refPType_m = ParticleProperties::getParticleType(type);
336 }
337
339 double getdE() const {
340 return 0.0;
341 }
342 double getGamma(int i) const {
343 return 0.0;
344 }
345
346 double getBeta(int i) const {
347 return 0.0;
348 }
349
350 void actT() {
351 }
352
353 const PartData* getReference() const {
354 return reference;
355 }
356
358 return 0.0;
359 }
360
361 /*
362
363 Domain Decomposition and Repartitioning
364
365 */
366
368 }
369
371 ippl::NDIndex<3> domain = getFieldLayout().getDomain();
372 for (unsigned int i = 0; i < Dim; i++)
373 grid[i] = domain[i].length();
374 }
375
376 void updateFields(const Vector_t<T, Dim>& /*hr*/, const Vector_t<T, Dim>& origin) {
377 /* \todo
378 this->getMesh().set_meshSpacing(&(hr_m[0]));
379 this->getMesh().set_origin(origin);
380 rho_m.initialize(this->getMesh(), getFieldLayout(), 1);
381 eg_m.initialize(this->getMesh(), getFieldLayout(), 1);
382 */
383 }
384
386 orb.initialize(fl, mesh, EFDMag_m);
387 }
388
390 // Repartition the domains
391 bool fromAnalyticDensity = false;
392 bool res = orb.binaryRepartition(this->R, fl, fromAnalyticDensity);
393
394 if (res != true) {
395 std::cout << "Could not repartition!" << std::endl;
396 return;
397 }
398 // Update
399 this->updateLayout(fl, mesh);
400 }
401
402 bool balance(unsigned int totalP) { //, int timestep = 1) {
403 int local = 0;
404 std::vector<int> res(ippl::Comm->size());
405 double threshold = 1.0;
406 double equalPart = (double)totalP / ippl::Comm->size();
407 double dev = std::abs((double)this->getLocalNum() - equalPart) / totalP;
408 if (dev > threshold) {
409 local = 1;
410 }
411 MPI_Allgather(&local, 1, MPI_INT, res.data(), 1, MPI_INT, ippl::Comm->getCommunicator());
412
413 for (unsigned int i = 0; i < res.size(); i++) {
414 if (res[i] == 1) {
415 return true;
416 }
417 }
418 return false;
419 }
420
422 for (int i = 0; i < ippl::Comm->size(); i++)
423 globalPartPerNode_m[i] = 0;
424
425 std::size_t localnum = this->getLocalNum();
426 // ADA gather(&localnum, &globalPartPerNode_m[0], 1);
427 }
428
429 size_t getLoadBalance(int p) const {
430 return globalPartPerNode_m[p];
431 }
432
433 /*
434
435 Mesh and Field Layout related functions
436
437 */
438
440
441 const Mesh_t<Dim>& getMesh() const {
442 return &getFieldLayout().get_Mesh();
443 }
444
446 // PLayout_t* layout = static_cast<Layout_t*>(&getLayout());
447 // return dynamic_cast<FieldLayout_t<Dim>&>(layout->getLayout().getFieldLayout());
448 return this->getFieldLayout();
449 }
450
451 // void setFieldLayout(FieldLayout_t* fLayout);
452
454 // Update local fields
455 static IpplTimings::TimerRef tupdateLayout = IpplTimings::getTimer("updateLayout");
456 IpplTimings::startTimer(tupdateLayout);
457 this->EFD_m.updateLayout(fl);
458 this->EFDMag_m.updateLayout(fl);
459
460 // Update layout with new FieldLayout
461 PLayout& layout = this->getLayout();
462 layout.updateLayout(fl, mesh);
463 IpplTimings::stopTimer(tupdateLayout);
464 static IpplTimings::TimerRef tupdatePLayout = IpplTimings::getTimer("updatePB");
465 IpplTimings::startTimer(tupdatePLayout);
466 this->update();
467 IpplTimings::stopTimer(tupdatePLayout);
468 }
469
470 bool isGridFixed() const {
471 return fixed_grid;
472 }
473
474 void boundp() {
475 /*
476 Assume rmin_m < 0.0
477 */
478
479 // if (!R.isDirty() && stateOfLastBoundP_ == unit_state_) return;
480 // if (!(this->R.isDirty() || this->ID.isDirty()) && stateOfLastBoundP_m == unit_state_m)
481 // return; //-DW
482
484
485 if (!isGridFixed()) {
486 const int dimIdx = (dcBeam_m ? 2 : 3);
487
493
494 // \todo this->updateDomainLength(nr_m);
497
498 double volume = 1.0;
499 for (int i = 0; i < dimIdx; i++) {
500 double length = std::abs(rmax_m[i] - rmin_m[i]);
501 if (length < 1e-10) {
502 rmax_m[i] += 1e-10;
503 rmin_m[i] -= 1e-10;
504 } else {
505 rmax_m[i] += dh_m * length;
506 rmin_m[i] -= dh_m * length;
507 }
508 hr_m[i] = (rmax_m[i] - rmin_m[i]) / (nr_m[i] - 1);
509 }
510 if (dcBeam_m) {
511 rmax_m[2] = rmin_m[2] + periodLength_m;
512 hr_m[2] = periodLength_m / (nr_m[2] - 1);
513 }
514 for (int i = 0; i < dimIdx; ++i) {
515 volume *= std::abs(rmax_m[i] - rmin_m[i]);
516 }
517
518 // ada if (volume < 1e-21 && this->this->getTotalNum() > 1 && std::abs(sum(Q)) >
519 // 0.0) {
520 // WARNMSG(level1 << "!!! Extremely high particle density detected !!!" << endl);
521 // }
522 // *ippl::Info << "It is a full boundp hz= " << hr_m << " rmax= " << rmax_m << "
523 // rmin= " << rmin_m
524 // << endl;
525
526 // ADA if (hr_m[0] * hr_m[1] * hr_m[2] <= 0) {
527 // throw GeneralClassicException("boundp() ", "h<0, can not build a mesh");
528 // }
529
530 Vector_t<T, Dim> origin =
531 rmin_m - Vector_t<T, Dim>(hr_m[0] / 2.0, hr_m[1] / 2.0, hr_m[2] / 2.0);
532 this->updateFields(hr_m, origin);
533 }
534 this->update();
535 // this->R.resetDirtyFlag();
536 }
537
540 // \todo this->updateDomainLength(nr_m);
541
542 std::vector<size_t> tmpbinemitted;
543
544 boundp();
545
546 size_t ne = 0;
547 const size_t localNum = this->getLocalNum();
548
549 double rzmean = 0.0; // momentsComputer_m.getMeanPosition()(2);
550 double rzrms = 0.0; // momentsComputer_m.getStandardDeviationPosition()(2);
551 const bool haveEnergyBins = weHaveEnergyBins();
552 if (haveEnergyBins) {
553 tmpbinemitted.resize(getNumberOfEnergyBins(), 0.0);
554 }
555 for (unsigned int i = 0; i < localNum; i++) {
556 if (this->Bin(i) < 0) { // ADA || (Options::remotePartDel > 0
557 // && std::abs(R(i)(2) - rzmean) < Options::remotePartDel * rzrms)) {
558 ne++;
559 // \todo this->destroy(1, i);
560 } else if (haveEnergyBins) {
561 tmpbinemitted[this->Bin(i)]++;
562 }
563 }
564
565 boundp();
566
569
570 // ippl::Comm->reduce(ne, ne, 1, std::plus<size_t>());
571 return ne;
572 }
573
574 size_t destroyT() {
575 std::unique_ptr<size_t[]> tmpbinemitted;
576
577 const size_t localNum = this->getLocalNum();
578 const size_t totalNum = this->getTotalNum();
579 size_t ne = 0;
580
581 if (weHaveEnergyBins()) {
582 tmpbinemitted = std::unique_ptr<size_t[]>(new size_t[getNumberOfEnergyBins()]);
583 for (int i = 0; i < getNumberOfEnergyBins(); i++) {
584 tmpbinemitted[i] = 0;
585 }
586 for (unsigned int i = 0; i < localNum; i++) {
587 if (this->Bin(i) < 0) {
588 // ADA destroy(1, i);
589 ++ne;
590 } else
591 tmpbinemitted[this->Bin(i)]++;
592 }
593 } else {
594 Inform dmsg("destroy: ", INFORM_ALL_NODES);
595 for (size_t i = 0; i < localNum; i++) {
596 if ((this->Bin(i) < 0)) {
597 ne++;
598 // ADA destroy(1, i);
599 }
600 }
601 }
602 /* ADA
603 if (ne > 0) {
604 performDestroy(true);
605 }
606 */
607
610
611 size_t newTotalNum = this->getLocalNum();
612
613 //
614 ippl::Comm->reduce(newTotalNum, newTotalNum, 1, std::plus<size_t>());
615
618
619 return totalNum - newTotalNum;
620 }
621
622 void set_meshEnlargement(double dh) {
623 dh_m = dh;
624 }
625
626 /*
627
628 Boundary Conditions
629
630 */
631
633 }
634
636 }
637
638 void setupBCs() {
640 }
641
643 this->setParticleBC(ippl::BC::PERIODIC);
644 }
645
646 /*
647 Misc Stuff
648 */
649
650 void runTests() { // Field Solver
651 Vector_t<T, Dim> ll(-0.005);
652 Vector_t<T, Dim> ur(0.005);
653
655
656 ippl::NDIndex<Dim> domain = getFieldLayout().getDomain();
657 for (unsigned int i = 0; i < Dim; i++)
658 nr_m[i] = domain[i].length();
659
660 for (int i = 0; i < 3; i++)
661 hr_m[i] = (ur[i] - ll[i]) / nr_m[i];
662
663 getFieldLayout().get_Mesh().set_meshSpacing(&(hr_m[0]));
664 getFieldLayout().get_Mesh().set_origin(ll);
665
666 rho_m.initialize(getFieldLayout().getMesh(), getFieldLayout(), 1);
667 eg_m.initialize(getFieldLayout().getMesh(), getFieldLayout(), 1);
668
669 // this->fs_m->solver_m->test(this);
670 }
671
672 void resetInterpolationCache(bool clearCache = false) {
673 // if (clearCache)
674 // interpolationCache_m.destroy(interpolationCache_m.size(), 0, true);
675 //
676 }
677
678 void swap(unsigned int i, unsigned int j) {
679 if (i >= this->getLocalNum() || j >= this->getLocalNum() || i == j)
680 return;
681
682 std::swap(this->R(i), this->R(j));
683 std::swap(this->P(i), this->P(j));
684 std::swap(this->Q(i), this->Q(j));
685 std::swap(this->M(i), this->M(j));
686 std::swap(this->Phi(i), this->Phi(j));
687 /*
688 std::swap(Ef(i), Ef[j]);
689 std::swap(Eftmp(i), Eftmp[j]);
690 std::swap(Bf(i), Bf[j]);
691 std::swap(Bin(i), Bin[j]);
692 std::swap(dt(i), dt[j]);
693 std::swap(PType(i), PType[j]);
694 std::swap(POrigin(i), POrigin[j]);
695 std::swap(TriID(i), TriID[j]);
696 std::swap(cavityGapCrossed(i), cavityGapCrossed[j]);
697 std::swap(bunchNum(i), bunchNum[j]);
698 */
699 }
700
701 double getRho(int x, int y, int z) {
702 // return rho_m[x][y][z].get();
703 return 0.0;
704 }
705
706 void gatherStatistics(unsigned int totalP) {
707 ippl::Comm->barrier();
708 std::cout << "Rank " << ippl::Comm->rank() << " has "
709 << (double)this->getLocalNum() / totalP * 100.0
710 << " percent of the total particles " << std::endl;
711 ippl::Comm->barrier();
712 }
713
714 // FIXME: unify methods, use convention that all particles have own dt
715 void switchToUnitlessPositions(bool use_dt_per_particle = false) {
716 bool hasToReset = false;
717 // if (!this->R.isDirty())
718 // hasToReset = true;
719
720 for (size_t i = 0; i < this->getLocalNum(); i++) {
721 double dt = getdT();
722 if (use_dt_per_particle)
723 dt = this->dt(i);
724
725 // this->R(i) /= Vector_t(Physics::c * dt);
726 this->R(i) /= Vector_t<T, Dim>(dt);
727 }
728
730
731 // if (hasToReset)
732 // this->R.resetDirtyFlag();
733 }
734
735 // FIXME: unify methods, use convention that all particles have own dt
736 void switchOffUnitlessPositions(bool use_dt_per_particle = false) {
737 }
738
741 std::size_t localnum = 0;
742 const Vector_t<T, Dim> meanR = get_rmean();
743
744 for (unsigned long k = 0; k < this->getLocalNum(); ++k)
745 if (std::abs(this->R(k)(0) - meanR(0)) > x(0)
746 || std::abs(this->R(k)(1) - meanR(1)) > x(1)
747 || std::abs(this->R(k)(2) - meanR(2)) > x(2)) {
748 ++localnum;
749 }
750
751 // ADA gather(&localnum, &globalPartPerNode_m[0], 1);
752
753 size_t npOutside = std::accumulate(
754 globalPartPerNode_m.get(), globalPartPerNode_m.get() + ippl::Comm->size(), 0,
755 std::plus<size_t>());
756
757 return npOutside;
758 }
759
770 unsigned int nBins, std::vector<double>& lineDensity, std::pair<double, double>& meshInfo) {
771 Vector_t<T, Dim> rmin, rmax;
772 get_bounds(rmin, rmax);
773
774 if (nBins < 2) {
775 Vector<int, 3> /*NDIndex<3>*/ grid;
776 // \todo this->updateDomainLength(grid);
777 nBins = grid[2];
778 }
779
780 double length = rmax(2) - rmin(2);
781 double zmin = rmin(2) - dh_m * length, zmax = rmax(2) + dh_m * length;
782 double hz = (zmax - zmin) / (nBins - 2);
783 double perMeter = 1.0 / hz; //(zmax - zmin);
784 zmin -= hz;
785
786 lineDensity.resize(nBins, 0.0);
787 std::fill(lineDensity.begin(), lineDensity.end(), 0.0);
788
789 const unsigned int lN = this->getLocalNum();
790 for (unsigned int i = 0; i < lN; ++i) {
791 const double z = this->R(i)(2) - 0.5 * hz;
792 unsigned int idx = (z - zmin) / hz;
793 double tau = (z - zmin) / hz - idx;
794
795 // ADA lineDensity[idx] += Q(i) * (1.0 - tau) * perMeter;
796 // lineDensity[idx + 1] += Q(i) * tau * perMeter;
797 }
798
799 // ADA reduce(&(lineDensity[0]), &(lineDensity[0]) + nBins, &(lineDensity[0]),
800 // OpAddAssign());
801
802 meshInfo.first = zmin;
803 meshInfo.second = hz;
804 }
805
806 void setBeamFrequency(double v) {
807 periodLength_m = 1.; // ADA Physics::c / f;
808 }
809
810 /*
811 Compatibility function
812
813 */
814
816 const Vector_t<T, Dim> maxE = max(eg_m);
817 // const double maxL = max(dot(eg_m,eg_m));
818 const Vector_t<T, Dim> minE = min(eg_m);
819 // *ippl::Info << "MaxE= " << maxE << " MinE= " << minE << endl;
820 return VectorPair_t(maxE, minE);
821 }
822
823 /*
824
825 Field Computations
826
827 */
828
830 }
831
833 void computeSelfFields(int b) {
834 }
835
836 // void setSolver(FieldSolver* fs);
837
839 // if (this->fs_m)
840 // return this->fs_m->hasValidSolver();
841 // else
842 return false;
843 }
844
845 bool getFieldSolverType() const {
846 // if (fs_m) {
847 // return true; // fs_m->getFieldSolverType();
848 // } else {
849 return false; // FieldSolverType::NONE;
850 //}
851 }
852
853 /*
854 I / O
855
856
857
858 Inform& print(Inform& os) {
859 return print(os);
860 }
861 */
862 /*
863 Inform& operator<<(Inform& os, PartBunch<PLayout, T, Dim>& p) {
864 return p.print(os);
865 }
866
867 */
868
869 Inform& print(Inform& os) {
870 // if (this->getLocalNum() != 0) { // to suppress Nans
871 Inform::FmtFlags_t ff = os.flags();
872
873 double pathLength = get_sPos();
874
875 os << std::scientific;
876 os << level1 << "\n";
877 os << "* ************** B U N C H "
878 "********************************************************* \n";
879 os << "* NP = " << this->getLocalNum() << "\n";
880
881 // os << "* Qtot = " << std::setw(17) << Util::getChargeString(std::abs(sum(Q)))
882 // << " "
883 // os << "Qi = " << std::setw(17) << Util::getChargeString(std::abs(qi_m)) << "\n";
884 // os << "* Ekin = " << std::setw(17)
885 // << Util::getEnergyString(get_meanKineticEnergy()) << " "
886 // << "dEkin = " << std::setw(17) << Util::getEnergyString(getdE()) << "\n";
887 // os << "* rmax = " << Util::getLengthString(rmax_m, 5) << "\n";
888 // os << "* rmin = " << Util::getLengthString(rmin_m, 5) << "\n";
889 /*
890 if (this->getTotalNum() >= 2) { // to suppress Nans
891 os << "* rms beam size = " << Util::getLengthString(get_rrms(), 5) << "\n";
892 os << "* rms momenta = " << std::setw(12) << std::setprecision(5) << get_prms()
893 << " [beta gamma]\n";
894 os << "* mean position = " << Util::getLengthString(get_rmean(), 5) << "\n";
895 os << "* mean momenta = " << std::setw(12) << std::setprecision(5) << get_pmean)
896 << " [beta gamma]\n";
897 os << "* rms emittance = " << std::setw(12) << std::setprecision(5) << get_emit()
898 << " (not normalized)\n";
899 os << "* rms correlation = " << std::setw(12) << std::setprecision(5) << get_rprms()
900 << "\n";
901 }
902
903 os << "* hr = " << Util::getLengthString(get_hr(), 5) << "\n";
904 os << "* dh = " << std::setw(13) << std::setprecision(5) << dh_m * 100
905 << " [%]\n";
906 os << "* t = " << std::setw(17) << Util::getTimeString(getT()) << " "
907 << "dT = " << std::setw(17) << Util::getTimeString(getdT()) << "\n";
908 os << "* spos = " << std::setw(17) << Util::getLengthString(pathLength) << "\n";
909 */
910 os << "* "
911 "********************************************************************************** "
912 << endl;
913 os.flags(ff);
914 // }
915 return os;
916 }
917
918 // save particles in case of one core
919 std::unique_ptr<Inform> pmsg_m;
920 std::unique_ptr<std::ofstream> f_stream;
921
922 void dumpData(int iteration) {
923 auto view = P.getView();
924
925 double Energy = 0.0;
926
927 Kokkos::parallel_reduce(
928 "Particle Energy", view.extent(0),
929 KOKKOS_LAMBDA(const int i, double& valL) {
930 double myVal = dot(view(i), view(i)).apply();
931 valL += myVal;
932 },
933 Kokkos::Sum<double>(Energy));
934
935 Energy *= 0.5;
936 double gEnergy = 0.0;
937
938 MPI_Reduce(&Energy, &gEnergy, 1, MPI_DOUBLE, MPI_SUM, 0, ippl::Comm->getCommunicator());
939
940 Inform csvout(NULL, "data/energy.csv", Inform::APPEND);
941 csvout.precision(10);
942 csvout.setf(std::ios::scientific, std::ios::floatfield);
943
944 csvout << iteration << " " << gEnergy << endl;
945
946 ippl::Comm->barrier();
947 }
948
949 /*
950
951 Field Particle Interpolation
952
953 */
954
955 void gatherCIC() {
956 // static IpplTimings::TimerRef gatherTimer = IpplTimings::getTimer("gather");
957 // IpplTimings::startTimer(gatherTimer);
958 gather(this->E, EFD_m, this->R);
959 // IpplTimings::stopTimer(gatherTimer);
960 }
961
962 void scatterCIC(unsigned int totalP, int iteration) {
963 // static IpplTimings::TimerRef scatterTimer = IpplTimings::getTimer("scatter");
964 // IpplTimings::startTimer(scatterTimer);
965 Inform m("scatter ");
966 EFDMag_m = 0.0;
967 scatter(Q, EFDMag_m, this->R);
968 // IpplTimings::stopTimer(scatterTimer);
969
970 static IpplTimings::TimerRef sumTimer = IpplTimings::getTimer("CheckCharge");
971 IpplTimings::startTimer(sumTimer);
972 double Q_grid = EFDMag_m.sum();
973
974 unsigned int Total_particles = 0;
975 unsigned int local_particles = this->getLocalNum();
976
977 MPI_Reduce(
978 &local_particles, &Total_particles, 1, MPI_UNSIGNED, MPI_SUM, 0,
979 ippl::Comm->getCommunicator());
980
981 double rel_error = std::fabs((Q_m - Q_grid) / Q_m);
982 m << "Rel. error in charge conservation = " << rel_error << endl;
983
984 if (ippl::Comm->rank() == 0) {
985 if (Total_particles != totalP || rel_error > 1e-10) {
986 std::cout << "Total particles in the sim. " << totalP << " "
987 << "after update: " << Total_particles << std::endl;
988 std::cout << "Total particles not matched in iteration: " << iteration << std::endl;
989 std::cout << "Q grid: " << Q_grid << "Q particles: " << Q_m << std::endl;
990 std::cout << "Rel. error in charge conservation: " << rel_error << std::endl;
991 exit(1);
992 }
993 }
994
995 IpplTimings::stopTimer(sumTimer);
996 }
997
998 /*
999
1000 Energy Bin Stuff
1001
1002 */
1003
1005 return true;
1006 }
1007
1009 return 1;
1010 }
1011
1013 return 1;
1014 }
1015
1017 return 1;
1018 }
1019
1020 void Rebin() {
1021 }
1022
1023 void setEnergyBins(int numberOfEnergyBins) {
1024 }
1025
1027 return false;
1028 }
1029
1030 void setTEmission(double t) {
1031 }
1032
1033 double getTEmission() {
1034 return 1.0;
1035 }
1036
1037 bool weHaveBins() const {
1038 return false;
1039 }
1040
1041 void setPBins(PartBins* pbin) {
1042 }
1043
1048
1049 size_t emitParticles(double eZ) {
1050 return 0;
1051 }
1052
1054 size_type total_particles = 0;
1055 size_type local_particles = this->getLocalNum();
1056
1057 MPI_Reduce(
1058 &local_particles, &total_particles, 1, MPI_UNSIGNED_LONG, MPI_SUM, 0,
1059 ippl::Comm->getCommunicator());
1060 // need to be stored
1061 }
1062
1063 void rebin() {
1064 this->Bin = 0;
1065 // ADA pbin_m->resetBins();
1066 // delete pbin_m; we did not allocate it!
1067 pbin_m = nullptr;
1068 }
1069
1071 if (pbin_m != nullptr)
1072 return pbin_m->getLastemittedBin();
1073 else
1074 return 0;
1075 }
1076
1077 void setLocalBinCount(size_t num, int bin) {
1078 if (pbin_m != nullptr) {
1079 // ADA pbin_m->setPartNum(bin, num);
1080 }
1081 }
1082
1084 void calcGammas() {
1085 }
1086
1088 double getBinGamma(int bin) {
1089 return 1.0;
1090 }
1091
1092 bool hasBinning() const {
1093 return (pbin_m != nullptr);
1094 }
1095
1097 void setBinCharge(int bin, double q) {
1098 // this->Q = where(eq(this->Bin, bin), q, 0.0);
1099 }
1100
1102 void setBinCharge(int bin) {
1103 // this->Q = where(eq(this->Bin, bin), this->qi_m, 0.0);
1104 }
1105
1107 double calcMeanPhi() {
1108 }
1109
1112 bool resetPartBinID2(const double eta) {
1113 return true;
1114 }
1115
1117 return true;
1118 }
1119
1120 /*
1121 Part Bunch Access Functions
1122
1123 */
1124
1125 double getPx(int i) {
1126 return 0.0;
1127 }
1128 double getPy(int i) {
1129 return 0.0;
1130 }
1131 double getPz(int i) {
1132 return 0.0;
1133 }
1134
1135 double getPx0(int i) {
1136 return 0.0;
1137 }
1138 double getPy0(int i) {
1139 return 0.0;
1140 }
1141
1142 double getX(int i) {
1143 return 0.0;
1144 }
1145 double getY(int i) {
1146 return 0.0;
1147 }
1148 double getZ(int i) {
1149 return 0.0;
1150 }
1151
1152 double getX0(int i) {
1153 return 0.0;
1154 }
1155 double getY0(int i) {
1156 return 0.0;
1157 }
1158
1159 void setZ(int i, double zcoo) {
1160 }
1161
1163 this->getLocalBounds(rmin, rmax);
1164
1165 double min[2 * Dim];
1166
1167 for (unsigned int i = 0; i < Dim; ++i) {
1168 min[2 * i] = rmin[i];
1169 min[2 * i + 1] = -rmax[i];
1170 }
1171
1172 // ippl::Comm->allreduce(min, 2 * Dim, std::less<double>());
1173
1174 for (unsigned int i = 0; i < Dim; ++i) {
1175 rmin[i] = min[2 * i];
1176 rmax[i] = -min[2 * i + 1];
1177 }
1178 }
1179
1181 const size_t localNum = this->getLocalNum();
1182 if (localNum == 0) {
1183 double maxValue = 1e8;
1184 rmin = Vector_t<T, Dim>({maxValue, maxValue, maxValue});
1185 rmax = Vector_t<T, Dim>({-maxValue, -maxValue, -maxValue});
1186 return;
1187 }
1188
1189 rmin = this->R(0);
1190 rmax = this->R(0);
1191 for (size_t i = 1; i < localNum; ++i) {
1192 for (unsigned short d = 0; d < 3u; ++d) {
1193 if (rmin(d) > this->R(i)(d))
1194 rmin(d) = this->R(i)(d);
1195 else if (rmax(d) < this->R(i)(d))
1196 rmax(d) = this->R(i)(d);
1197 }
1198 }
1199 }
1200
1201 std::pair<Vector_t<T, Dim>, double> getBoundingSphere() {
1202 Vector_t<T, Dim> rmin, rmax;
1203 get_bounds(rmin, rmax);
1204
1205 std::pair<Vector_t<T, Dim>, double> sphere;
1206 sphere.first = 0.5 * (rmin + rmax);
1207 sphere.second = std::sqrt(dot(rmax - sphere.first, rmax - sphere.first));
1208
1209 return sphere;
1210 }
1211
1212 std::pair<Vector_t<T, Dim>, double> getLocalBoundingSphere() {
1213 Vector_t<T, Dim> rmin, rmax;
1214 getLocalBounds(rmin, rmax);
1215
1216 std::pair<Vector_t<T, Dim>, double> sphere;
1217 sphere.first = 0.5 * (rmin + rmax);
1218 sphere.second = std::sqrt(dot(rmax - sphere.first, rmax - sphere.first));
1219
1220 return sphere;
1221 }
1222
1224 bounds(this->P, min, max);
1225 }
1226
1227 /*
1228 Misc Bunch related quantities
1229
1230 */
1231
1232 // get 2nd order momentum matrix
1233 // matrix_t getSigmaMatrix() const;
1234
1235 void setdT(double dt) {
1236 dt_m = dt;
1237 }
1238
1239 double getdT() const {
1240 return dt_m;
1241 }
1242
1243 void setT(double t) {
1244 t_m = t;
1245 }
1246
1247 void incrementT() {
1248 t_m += dt_m;
1249 }
1250
1251 double getT() const {
1252 return t_m;
1253 }
1254
1261
1262 void set_sPos(double s) {
1263 spos_m = s;
1264 }
1265
1266 double get_gamma() const {
1267 return 1.00;
1268 }
1269
1271 return 0.0;
1272 }
1273
1275 return rmin_m;
1276 }
1278 return rmax_m;
1279 }
1280
1282 return Vector_t<T, Dim>(0.0);
1283 }
1285 return Vector_t<T, Dim>(0.0);
1286 }
1288 return Vector_t<T, Dim>(0.0);
1289 }
1291 return Vector_t<T, Dim>(0.0);
1292 }
1294 return Vector_t<T, Dim>(0.0);
1295 }
1297 return Vector_t<T, Dim>(0.0);
1298 }
1303 return Vector_t<T, Dim>(0.0);
1304 }
1306 return Vector_t<T, Dim>(0.0);
1307 }
1309 return Vector_t<T, Dim>(0.0);
1310 }
1312 return Vector_t<T, Dim>(0.0);
1313 }
1315 return Vector_t<T, Dim>(0.0);
1316 }
1318 return Vector_t<T, Dim>(0.0);
1319 }
1321 return Vector_t<T, Dim>(0.0);
1322 }
1336 return Vector_t<T, Dim>(0.0);
1337 }
1338
1339 double get_Dx() const {
1340 return 0.0;
1341 }
1342 double get_Dy() const {
1343 return 0.0;
1344 }
1345 double get_DDx() const {
1346 return 0.0;
1347 }
1348 double get_DDy() const {
1349 return 0.0;
1350 }
1351
1352 /*
1353 Some quantities related to integrations/tracking
1354 */
1355
1356 void setStepsPerTurn(int n) {
1357 stepsPerTurn_m = n;
1358 }
1359
1360 int getStepsPerTurn() const {
1361 return stepsPerTurn_m;
1362 }
1363
1365 void setGlobalTrackStep(long long n) {
1367 }
1368
1369 long long getGlobalTrackStep() const {
1370 return globalTrackStep_m;
1371 }
1372
1374 void setLocalTrackStep(long long n) {
1375 localTrackStep_m = n;
1376 }
1377
1381 }
1382
1383 long long getLocalTrackStep() const {
1384 return localTrackStep_m;
1385 }
1386
1387 void setNumBunch(short n) {
1388 numBunch_m = n;
1389 bunchTotalNum_m.resize(n);
1390 bunchLocalNum_m.resize(n);
1391 }
1392
1393 short getNumBunch() const {
1394 return numBunch_m;
1395 }
1396
1398 globalMeanR_m = globalMeanR;
1399 }
1400
1404
1405 // void setGlobalToLocalQuaternion(Quaternion_t globalToLocalQuaternion);
1406 // Quaternion_t getGlobalToLocalQuaternion();
1407
1408 void setSteptoLastInj(int n) {
1409 SteptoLastInj_m = n;
1410 }
1411
1412 int getSteptoLastInj() const {
1413 return SteptoLastInj_m;
1414 }
1415
1416 double calculateAngle(double x, double y) {
1417 return 0.0;
1418 }
1419
1420 double get_sPos() const {
1421 return spos_m;
1422 }
1423};
1424
1426
1427#endif
double dot(const Vector3D &lhs, const Vector3D &rhs)
Vector dot product.
Definition Vector3D.cpp:118
VariantFromConditionalTypes< FFTSolver_t< T, Dim >, P3MSolver_t< T, Dim >, OpenSolver_t< T, Dim > > Solver_t
ConditionalType< Dim==2||Dim==3, ippl::FFTPeriodicPoissonSolver< VField_t< T, Dim >, Field_t< Dim > > > FFTSolver_t
PartBunch< PLayout_t< double, 3 >, double, 3 > PartBunch_t
ConditionalType< Dim==3, ippl::FFTOpenPoissonSolver< VField_t< T, Dim >, Field_t< Dim > > > OpenSolver_t
ConditionalType< Dim==3, ippl::P3MSolver< VField_t< T, Dim >, Field_t< Dim > > > P3MSolver_t
ippl::FieldLayout< Dim > FieldLayout_t
Definition PBunchDefs.h:27
ippl::Field< double, Dim, ViewArgs... > Field_t
Definition PBunchDefs.h:30
ippl::Field< ippl::Vector< T, Dim >, Dim, ViewArgs... > VField_t
Definition PBunchDefs.h:33
ippl::UniformCartesian< double, 3 > Mesh_t
Definition PBunchDefs.h:19
ippl::Vector< T, Dim > Vector_t
ippl::Field< T, Dim, Mesh_t< Dim >, Centering_t< Dim >, ViewArgs... > Field
ippl::detail::size_type size_type
ippl::ParticleAttrib< T > ParticleAttrib
constexpr unsigned Dim
Definition datatypes.h:6
ippl::OrthogonalRecursiveBisection< Field< double, Dim >, T > ORB
int getSteptoLastInj() const
double getX(int i)
double get_Dy() const
void updateDomainLength(Vector_t< int, 3 > &grid)
void gatherLoadBalanceStatistics()
void setdT(double dt)
double get_sPos() const
void initializeORB(FieldLayout_t< Dim > &fl, Mesh_t< Dim > &mesh)
Vector_t< T, Dim > get_hr() const
double getMassPerParticle() const
void setBeamFrequency(double v)
void setMass(double mass)
Vector_t< T, Dim > get_prms() const
void getLocalBounds(Vector_t< T, Dim > &rmin, Vector_t< T, Dim > &rmax) const
void updateFields(const Vector_t< T, Dim > &, const Vector_t< T, Dim > &origin)
bool hasFieldSolver()
bool weHaveBins() const
void computeSelfFields()
ParticleLayout< double, 3 > & getLayout()
Definition .PartBunch.h:115
size_t getNumberOfEmissionSteps()
void repartition(FieldLayout_t< Dim > &fl, Mesh_t< Dim > &mesh)
void incTrackSteps()
void setBinCharge(int bin, double q)
Set the charge of one bin to the value of q and all other to zero.
double getPx0(int i)
size_t destroyT()
double getdT() const
const Mesh_t< Dim > & getMesh() const
void setT(double t)
void setChargeZeroPart(double q)
Vector_t< T, Dim > get_99_99Percentile() const
std::pair< Vector_t< double, 3 >, Vector_t< double, 3 > > VectorPair_t
ParticleAttrib< Vector_t< double, Dim > > Ef
size_t calcNumPartsOutside(Vector_t< T, Dim > x)
returns the number of particles outside of a box defined by x
Vector_t< T, Dim > get_normalizedEps_99Percentile() const
size_t boundp_destroyT()
PartBunch & operator=(const PartBunch &)=delete
void setLocalTrackStep(long long n)
step in a TRACK command
int getLastemittedBin()
int getLastEmittedEnergyBin()
bool resetPartBinBunch()
double get_Dx() const
size_t getLoadBalance(int p) const
void setPType(const std::string &type)
void calcBeamParameters()
Vector_t< T, Dim > get_maxExtent() const
void boundp()
double getdE() const
Vector_t< T, Dim > get_pmean() const
double getQ() const
Access to reference data.
void swap(unsigned int i, unsigned int j)
std::pair< Vector_t< T, Dim >, double > getLocalBoundingSphere()
bool resetPartBinID2(const double eta)
void setBinCharge(int bin)
Set the charge of all other the ones in bin to zero.
void resetInterpolationCache(bool clearCache=false)
double getE() const
void do_binaryRepart()
double getT() const
void resetQ(double q)
Set reference data.
void setPBins(PartBins *pbin)
std::pair< Vector_t< T, Dim >, double > getBoundingSphere()
void switchToUnitlessPositions(bool use_dt_per_particle=false)
double getCharge() const
get the total charge per simulation particle
double getY(int i)
double getInitialBeta() const
void updateNumTotal()
void get_PBounds(Vector_t< T, Dim > &min, Vector_t< T, Dim > &max) const
Vector_t< T, Dim > get_rrms() const
double getTEmission()
void setBCForDCBeam()
Vector_t< T, Dim > get_halo() const
Vector_t< T, Dim > get_origin() const
short getNumBunch() const
double getRho(int x, int y, int z)
void calcBeamParametersInitial()
double getGamma(int i) const
Vector_t< T, Dim > get_68Percentile() const
bool balance(unsigned int totalP)
void setTEmission(double t)
bool hasBinning() const
void incrementT()
void setLocalBinCount(size_t num, int bin)
void setBCAllPeriodic()
double getX0(int i)
double get_meanKineticEnergy()
int getNumberOfEnergyBins()
long long getLocalTrackStep() const
double getP() const
void gatherStatistics(unsigned int totalP)
void get_bounds(Vector_t< T, Dim > &rmin, Vector_t< T, Dim > &rmax) const
void set_sPos(double s)
double calculateAngle(double x, double y)
size_t getLocalNum() const
void setZ(int i, double zcoo)
double calcMeanPhi()
calculate average angle of longitudinal direction of bins
void resizeMesh()
double get_DDy() const
bool isGridFixed() const
double getM() const
Vector_t< double, Dim > R(size_t i)
Definition PartBunch.h:276
Vector_t< T, Dim > get_99Percentile() const
void resetM(double m)
Vector_t< T, Dim > get_normalizedEps_99_99Percentile() const
void calcLineDensity(unsigned int nBins, std::vector< double > &lineDensity, std::pair< double, double > &meshInfo)
calculates the 1d line density (not normalized) and append it to a file.
void setupBCs()
Vector_t< T, Dim > get_rmean() const
double getChargePerParticle() const
get the macro particle charge
double getPy(int i)
ParticleAttrib< Vector_t< double, Dim > > P
Vector_t< T, Dim > getGlobalMeanR()
void setSteptoLastInj(int n)
double get_DDx() const
ParticleAttrib< Vector_t< double, Dim > > Eftmp
double getPz(int i)
void setStepsPerTurn(int n)
Vector_t< T, Dim > get_pmean_Distribution() const
Vector_t< T, Dim > get_emit() const
ippl::BConds< Field< T, Dim >, Dim > bc_type
const PartData * getReference() const
Vector_t< T, Dim > get_95Percentile() const
VectorPair_t getEExtrema()
double get_gamma() const
std::unique_ptr< std::ofstream > f_stream
void setEnergyBins(int numberOfEnergyBins)
Vector_t< T, Dim > get_centroid() const
double getY0(int i)
void setGlobalTrackStep(long long n)
step in multiple TRACK commands
PartBunch(PLayout &pl, Vector_t< double, Dim > hr, Vector_t< double, Dim > rmin, Vector_t< double, Dim > rmax, std::array< bool, Dim > decomp, double Qtot)
bool getIfBeamEmitting()
double getPy0(int i)
long long getGlobalTrackStep() const
void scatterCIC(unsigned int totalP, int iteration)
void initialize(FieldLayout_t< Dim > &fl, Mesh_t< Dim > &mesh)
ParticleAttrib< Vector_t< double, Dim > > E
double getInitialGamma() const
void updateLayout(FieldLayout_t< Dim > &fl, Mesh_t< Dim > &mesh)
Inform & print(Inform &os)
double getBeta(int i) const
void setBCAllOpen()
double getZ(int i)
void setCharge(double q)
FieldLayout_t< Dim > & getFieldLayout()
Vector_t< T, Dim > get_normalizedEps_95Percentile() const
void switchOffUnitlessPositions(bool use_dt_per_particle=false)
void setMassZeroPart(double mass)
void gatherCIC()
Vector_t< T, Dim > get_normalizedEps_68Percentile() const
double getCouplingConstant() const
size_t getTotalNum() const
void calcGammas()
Compute the gammas of all bins.
void set_meshEnlargement(double dh)
void computeSelfFields(int b)
double getEmissionDeltaT()
Vector_t< T, Dim > get_rprms() const
double getPx(int i)
void setCouplingConstant(double c)
void setGlobalMeanR(Vector_t< T, Dim > globalMeanR)
bool getFieldSolverType() const
Return the fieldsolver type if we have a fieldsolver.
void dumpData(int iteration)
bool weHaveEnergyBins()
ippl::ParticleBase< PLayout > Base
void setNumBunch(short n)
int getStepsPerTurn() const
double getBinGamma(int bin)
Get gamma of one bin.
void runTests()
Vector_t< T, Dim > get_norm_emit() const
ParticleAttrib< Vector_t< double, Dim > > Bf
size_t emitParticles(double eZ)
Emit particles in the given bin i.e. copy the particles from the bin structure into the particle cont...
Particle reference data.
Definition PartData.h:37