OPALX (Object Oriented Parallel Accelerator Library for Exascal) MINIorX
OPALX
FlatTop.cpp
Go to the documentation of this file.
1#include "Distribution.h"
2#include "SamplingBase.hpp"
3#include "FlatTop.h"
4#include <memory>
5#include <cmath>
6
10using GeneratorPool = typename Kokkos::Random_XorShift64_Pool<>;
11using Dist_t = ippl::random::NormalDistribution<double, 3>;
12
13FlatTop::FlatTop(std::shared_ptr<ParticleContainer_t> &pc,
14 std::shared_ptr<FieldContainer_t> &fc,
15 std::shared_ptr<Distribution_t> &opalDist)
16 : SamplingBase(pc, fc, opalDist), rand_pool_m(determineRandInit()) {
17 setParameters(opalDist);
18}
19
20void FlatTop::setWithDomainDecomp(bool withDomainDecomp) {
21 withDomainDecomp_m = withDomainDecomp;
22}
23
25 extern Inform* gmsg;
26 size_t randInit;
27 if (Options::seed == -1) {
28 randInit = 1234567;
29 *gmsg << "* Seed = " << randInit << " on all ranks" << endl;
30 } else {
31 randInit = static_cast<size_t>(Options::seed + 100 * ippl::Comm->rank());
32 }
33 return randInit;
34}
35
36void FlatTop::setParameters(const std::shared_ptr<Distribution_t> &opalDist) {
37 emitting_m = opalDist->emitting_m;
38 // time span of fall is [0, riseTime, riseTime+flattopTime, fallTime+flattopTime+riseTime ]
39 sigmaTFall_m = opalDist_m->getSigmaTFall();
40 sigmaTRise_m = opalDist_m->getSigmaTRise();
41 cutoffR_m = opalDist_m->getCutoffR();
42
43 fallTime_m = sigmaTFall_m * cutoffR_m[2]; // fall is [0, fallTime]
44 flattopTime_m = opalDist->getTPulseLengthFWHM()
45 - std::sqrt(2.0 * std::log(2.0)) * (sigmaTRise_m + sigmaTFall_m);
46 if (flattopTime_m < 0.0) {
47 flattopTime_m = 0.0;
48 }
50
52 opalDist_m->setTEmission(emissionTime_m);
53
54 // These expression are take from the old OPAL
55 // I think normalizedFlankArea is int_0^{cutoff} exp(-(x/sigma)^2/2 ) / sigma
56 // Instead of int_0^{cutoff} exp(-(x/sigma)^2/2 ) / sqrt(2*pi) / sigma, which is strange!
57 // So the distribution of tails are exp(-(x/sigma)^2/2 ) and not Gaussian!
58 normalizedFlankArea_m = 0.5 * std::sqrt(Physics::two_pi) * std::erf(cutoffR_m[2] / std::sqrt(2.0));
60
61 // make sure only z direction is decomposed
62 fc_m->setDecomp({false, false, true});
63}
64
65void FlatTop::generateUniformDisk(size_type nlocal, size_t nNew) {
66
67 GeneratorPool rand_pool = rand_pool_m;
71
72 view_type Rview = pc_m->R.getView();
73 view_type Pview = pc_m->P.getView();
74
75 double pi = Physics::pi;
76 Vector_t<double, 3> sigmaR = opalDist_m->getSigmaR();
77 // Sample (Rx,Ry) on a unit ring, then scale with sigmaRx and sigmaRy, set Px=Py=0
78 Kokkos::parallel_for(
79 "unitDisk", Kokkos::RangePolicy<>(nlocal, nlocal+nNew), KOKKOS_LAMBDA(const size_t j) {
80 auto generator = rand_pool.get_state();
81 double r = Kokkos::sqrt( generator.drand(0., 1.) );
82 double theta = 2.0 * pi * generator.drand(0., 1.);
83 rand_pool.free_state(generator);
84
85 Rview(j)[0] = r * Kokkos::cos(theta) * sigmaR[0];
86 Rview(j)[1] = r * Kokkos::sin(theta) * sigmaR[1];
87 Rview(j)[2] = 0.0;
88 Pview(j)[0] = 0.0;
89 Pview(j)[1] = 0.0;
90 Pview(j)[2] = 0.0;
91 });
92 Kokkos::fence();
93}
94
98
99void FlatTop::generateParticles(size_t& numberOfParticles, Vector_t<double, 3> nr) {
100
101 setNr(nr);
102
103 // initial allocation is similar for both emitting and non-emitting cases
104 allocateParticles(numberOfParticles);
105
106 if(emitting_m){
107 // set nlocal to 0 for the very first time step, before sampling particles
108 //pc_m->setLocalNum(0);
109
110 Kokkos::View<bool*> tmp_invalid("tmp_invalid", 0);
111 // \todo might be abuse of semantics: maybe think about new pc_m->setTotalNum or pc_m->updateTotal function instead?
112 pc_m->destroy(tmp_invalid, pc_m->getLocalNum());
113 }
114}
115
116double FlatTop::FlatTopProfile(double t){
117 double t0;
118 if(t<riseTime_m){
119 t0 = riseTime_m;
120 return exp( -pow((t-t0)/sigmaTRise_m,2) /2. );
121 // In the old opal, tails seem to be exp(-x^2/sigma^2/2) rather than Gaussian with normalizing factor.
122 }
123 else if( t>riseTime_m && t<riseTime_m + flattopTime_m){
124 return 1.;
125 }
128 return exp( -pow((t-t0)/sigmaTFall_m,2)/2. );
129 // In the old opal, tails seem to be exp(-x^2/sigma^2/2) rather than Gaussian with normalizing factor.
130 }
131 else
132 return 0.;
133}
134
135size_t FlatTop::computeNlocalUniformly(size_t nglobal){
136 MPI_Comm comm = MPI_COMM_WORLD;
137 int nranks;
138 int rank;
139 MPI_Comm_size(comm, &nranks);
140 MPI_Comm_rank(comm, &rank);
141
142 size_type nlocal = floor(nglobal/nranks);
143
144 size_t remaining = nglobal - nlocal*nranks;
145
146 if(remaining>0 && rank==0){
147 nlocal += remaining;
148 }
149
150 return nlocal;
151}
152
153double FlatTop::integrateTrapezoidal(double x1, double x2, double y1, double y2){
154 return 0.5 * (y1+y2) * fabs(x2-x1);
155}
156
157void FlatTop::initDomainDecomp(double BoxIncr) {
158 auto *mesh = &fc_m->getMesh();
159 auto *FL = &fc_m->getFL();
160 Vector_t<double, 3> sigmaR = opalDist_m->getSigmaR();
161 ippl::Vector<double, 3> o;
162 ippl::Vector<double, 3> e;
163 double tol = 1e-15; // enlarge grid by tol to avoid missing particles on boundaries
164 o[0] = -sigmaR[0] - tol;
165 e[0] = sigmaR[0] + tol;
166 o[1] = -sigmaR[1] - tol;
167 e[1] = sigmaR[1] + tol;
168 o[2] = 0.0 - tol;
169 e[2] = Physics::c * emissionTime_m + tol;
170
171 ippl::Vector<double, 3> l = e - o;
172 hr_m = (1.0+BoxIncr/100.)*(l / nr_m);
173 mesh->setMeshSpacing(hr_m);
174 mesh->setOrigin(o-0.5*hr_m*BoxIncr/100.);
175 pc_m->getLayout().updateLayout(*FL, *mesh);
176}
177
178double FlatTop::countEnteringParticlesPerRank(double t0, double tf){
179 size_type nlocalNew=0;
180 double tArea = 0.0;
181 tArea = integrateTrapezoidal(t0, tf, FlatTopProfile(t0), FlatTopProfile(tf));
182 size_type totalNew = floor(totalN_m * tArea / distArea_m);
183 nlocalNew = 0;
184
185 if(totalNew>0){
187 // the same number of particles per rank
188 nlocalNew = computeNlocalUniformly(totalNew);
189 }
190 else{
191 // select number of particles per rank using estimated domain decomposition at final emission time
192 // find min/max of particle positions for [t,t+dt]
193 Vector_t<double, 3> prmin, prmax;
194 Vector_t<double, 3> sigmaR = opalDist_m->getSigmaR();
195 prmin[0] = -sigmaR[0];
196 prmax[0] = sigmaR[0];
197 prmin[1] = -sigmaR[1];
198 prmax[1] = sigmaR[1];
199 prmin[2] = std::max(0.0, Physics::c * t0);
200 prmax[2] = Physics::c*tf;
201
202 double dx = prmax[0] - prmin[0];
203 double dy = prmax[1] - prmin[1];
204 double dz = prmax[2] - prmin[2];
205
206 if (dx <= 0 || dy <= 0 || dz <= 0) {
207 throw std::runtime_error("Invalid global particle volume: prmax must be greater than prmin.");
208 }
209
210 double globalpvolume = dx * dy * dz;
211
212 // find min/max of subdomains for the current rank
213 auto regions = pc_m->getLayout().getRegionLayout().gethLocalRegions();
214 int rank = ippl::Comm->rank();
215
216 if (rank < 0 || static_cast<size_t>(rank) >= regions.size()) {
217 throw std::runtime_error("Invalid rank index in gethLocalRegions()");
218 }
219
220 Vector_t<double, 3> locrmin, locrmax;
221 for (unsigned d = 0; d < Dim; ++d) {
222 locrmax[d] = regions(rank)[d].max();
223 locrmin[d] = regions(rank)[d].min();
224 }
225
226 if (prmax[0] >= locrmin[0] && prmin[0] <= locrmax[0] &&
227 prmax[1] >= locrmin[1] && prmin[1] <= locrmax[1] &&
228 prmax[2] >= locrmin[2] && prmin[2] <= locrmax[2]) {
229
230 double x1 = std::max(prmin[0], locrmin[0]);
231 double x2 = std::min(prmax[0], locrmax[0]);
232 double y1 = std::max(prmin[1], locrmin[1]);
233 double y2 = std::min(prmax[1], locrmax[1]);
234 double z1 = std::max(prmin[2], locrmin[2]);
235 double z2 = std::min(prmax[2], locrmax[2]);
236
237 if (x2 >= x1 && y2 >= y1 && z2 >= z1) {
238 double locpvolume = (x2 - x1) * (y2 - y1) * (z2 - z1);
239 if (globalpvolume > 0) {
240 nlocalNew = static_cast<int>(totalNew * locpvolume / globalpvolume);
241 }
242 else{
243 nlocalNew = 0;
244 }
245 } else {
246 nlocalNew = 0;
247 }
248 }
249 }
250 }
251 return nlocalNew;
252}
253
254void FlatTop::allocateParticles(size_t numberOfParticles){
255 totalN_m = numberOfParticles;
256
257 size_type nlocal;
258
260
261 pc_m->create(nlocal);
262}
263
264void FlatTop::emitParticles(double t, double dt) {
265 extern Inform* gmsg;
266
267 // count number of new particles to be emitted
269
270 // current number of particles per rank
271 size_type nlocal = pc_m->getLocalNum();
272
273 // extend particle container to accomodate new particles
274 pc_m->create(nNew);
275
276 // generate new particles on uniform disc
277 *gmsg << "* generate particles on a disc" << endl;
278 generateUniformDisk(nlocal, nNew);
279
280 *gmsg << "* new particles emmitted" << endl;
281}
282
284 size_type nNew;
285 int rank, numRanks;
286 double t = 0.0;
287 double c = Physics::c;
288
289 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
290 MPI_Comm_size(MPI_COMM_WORLD, &numRanks);
291
292 std::string filename = "timeNpart_" + std::to_string(rank) + ".txt";
293 std::ofstream file(filename);
294
295 // Check if the file opened successfully
296 if (!file.is_open()) {
297 throw std::runtime_error("Failed to open file: " + filename);
298 }
299
300 for(size_type i=0; i<nsteps; i++){
301 nNew = countEnteringParticlesPerRank(t, t + dt);
302
303 // current number of particles per rank
304 size_type nlocal = pc_m->getLocalNum();
305
306 // extend particle container to accomodate new particles
307 pc_m->create(nNew);
308
309 // generate new particles on uniform disc
310 generateUniformDisk(nlocal, nNew);
311
312 // write to a file
313 auto rViewDevice = pc_m->R.getView();
314 auto rView = Kokkos::create_mirror_view(rViewDevice);
315 Kokkos::deep_copy(rView,rViewDevice);
316
317 for(size_type j=nlocal; j<nlocal+nNew; j++){
318 file << t << " " << (emissionTime_m-t)*c << " " << rView(j)[0] << " " << rView(j)[1] << "\n";
319 }
320 file.flush(); // Ensure data is written to disk
321
322 t = t + dt;
323 }
324 file.close();
325 ippl::Comm->barrier();
326}
327
328void FlatTop::testEmitParticles(size_type nsteps, double dt) {
329 double t = 0.0;
330
331 for(size_type i=0; i<nsteps; i++){
332 emitParticles(t, dt);
333
334 t = t + dt;
335 }
336}
ParticleContainer< double, 3 > ParticleContainer_t
Definition Component.h:30
Inform * gmsg
Definition changes.cpp:7
ippl::Vector< T, Dim > Vector_t
typename Kokkos::Random_XorShift64_Pool<> GeneratorPool
typename ippl::detail::ViewType< ippl::Vector< double, Dim >, 1 >::view_type view_type
Distribution Distribution_t
Definition FlatTop.cpp:9
ippl::random::NormalDistribution< double, 3 > Dist_t
Definition FlatTop.cpp:11
FieldContainer< double, 3 > FieldContainer_t
Defines the FlatTop class used for sampling emitting particles.
const int nr
const double pi
Definition datatypes.h:4
constexpr unsigned Dim
Definition datatypes.h:6
constexpr double two_pi
The value of.
Definition Physics.h:33
constexpr double c
The velocity of light in m/s.
Definition Physics.h:45
constexpr double pi
The value of.
Definition Physics.h:30
int seed
The current random seed.
Definition Options.cpp:37
void setNr(Vector_t< double, 3 > nr)
Sets the number of grid points per direction.
Definition FlatTop.cpp:95
void testEmitParticles(size_type nsteps, double dt) override
Tests particle emission over a given number of steps.
Definition FlatTop.cpp:328
double distArea_m
Total area of the flattop distribution.
Definition FlatTop.h:65
double FlatTopProfile(double t)
Computes the flat-top profile value at a given time.
Definition FlatTop.cpp:116
bool emitting_m
Flag for particle emission status.
Definition FlatTop.h:71
Vector_t< double, 3 > nr_m
Number of grid points per direction.
Definition FlatTop.h:75
static size_t determineRandInit()
Determines the random seed initialization.
Definition FlatTop.cpp:24
void allocateParticles(size_t numberOfParticles)
Allocates memory for a given number of particles.
Definition FlatTop.cpp:254
double emissionTime_m
Total emission time.
Definition FlatTop.h:74
double sigmaTRise_m
Standard deviation for rise time profile.
Definition FlatTop.h:67
size_type totalN_m
Total number of particles.
Definition FlatTop.h:72
GeneratorPool rand_pool_m
Random number generator pool.
Definition FlatTop.h:62
double sigmaTFall_m
Standard deviation for fall time profile.
Definition FlatTop.h:66
bool withDomainDecomp_m
Flag for domain decomposition.
Definition FlatTop.h:73
void initDomainDecomp(double BoxIncr) override
Initializes the domain decomposition.
Definition FlatTop.cpp:157
Vector_t< double, 3 > hr_m
Grid spacing.
Definition FlatTop.h:76
double integrateTrapezoidal(double x1, double x2, double y1, double y2)
Integrates using the trapezoidal rule.
Definition FlatTop.cpp:153
double fallTime_m
Time duration for the fall phase.
Definition FlatTop.h:69
size_t computeNlocalUniformly(size_t nglobal)
Computes the local number of particles uniformly distributed among ranks.
Definition FlatTop.cpp:135
double normalizedFlankArea_m
Normalized area of the distribution flanks.
Definition FlatTop.h:64
void setWithDomainDecomp(bool withDomainDecomp) override
Sets whether to use domain decomposition.
Definition FlatTop.cpp:20
void emitParticles(double t, double dt) override
Emits new particles within a given time interval.
Definition FlatTop.cpp:264
double riseTime_m
Time duration for the rise phase.
Definition FlatTop.h:70
void generateParticles(size_t &numberOfParticles, Vector_t< double, 3 > nr) override
Generates particles with a given number and grid configuration.
Definition FlatTop.cpp:99
double countEnteringParticlesPerRank(double t0, double tf)
Counts the number of particles entering per rank in a given time interval.
Definition FlatTop.cpp:178
void generateUniformDisk(size_type nlocal, size_t nNew)
Generates particles (x,y) uniformly on a disk distribution.
Definition FlatTop.cpp:65
void testNumEmitParticles(size_type nsteps, double dt) override
Tests the number of emitted particles over a given number of steps.
Definition FlatTop.cpp:283
Vector_t< double, 3 > cutoffR_m
Cutoff radius.
Definition FlatTop.h:68
double flattopTime_m
Time duration of when the time profile is flat.
Definition FlatTop.h:63
FlatTop(std::shared_ptr< ParticleContainer_t > &pc, std::shared_ptr< FieldContainer_t > &fc, std::shared_ptr< Distribution_t > &opalDist)
Constructor for FlatTop.
Definition FlatTop.cpp:13
ippl::detail::size_type size_type
Definition FlatTop.h:61
void setParameters(const std::shared_ptr< Distribution_t > &opalDist)
Sets distribution parameters.
Definition FlatTop.cpp:36
SamplingBase(std::shared_ptr< ParticleContainer_t > &pc, std::shared_ptr< FieldContainer_t > &fc, std::shared_ptr< Distribution_t > &dist)
std::shared_ptr< FieldContainer_t > fc_m
std::shared_ptr< Distribution_t > opalDist_m
std::shared_ptr< ParticleContainer_t > pc_m