IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
LandauDampingParameterList.cpp
Go to the documentation of this file.
1// Landau Damping Test that also has a parameterlist
2// Modifications are tagged with MOD START , MOD END
3// Usage:
4// srun ./LandauDampingParameterList
5// <nx> [<ny>...] <Np> <Nt> <stype [pre_args]>
6// <lbthres> --overallocate <ovfactor> --info 10
7// nx = No. cell-centered points in the x-direction
8// ny... = No. cell-centered points in the y-, z-, ...-direction
9// Np = Total no. of macro-particles in the simulation
10// Nt = Number of time steps
11// stype = Field solver type (FFT, CG and PCG supported) , PCG needs extra arguments
12// pre_args = when using the PCG, a preconditioner and its corresponding parameters need to
13// be defined. Supported pre_args are: [jacobi, chebyshev <degree>, newton <level>,
14// gauss-seidel <inner outer comm> , richardson <iter comm>]
15// jacobi = diagonal preconditioner
16// chebyshev = polynomial preconditioner (algo=chebyshev) of degree <degree>
17// newton = polynomial preconditioner (algo=newton) of degree 2^<level> - 1
18// gauss-seidel = symmetric 2-step gauss-seidel with <inner> inner iterations, <outer> outer
19// iterations richardson = jacobi-richardson preconditioner with <iter> iterations comm = 0
20// means that there is no halo-exchange during preconditioning, 1 means there is halo-exchange
21// lbthres = Load balancing threshold i.e., lbthres*100 is the maximum load imbalance
22// percentage which can be tolerated and beyond which
23// particle load balancing occurs. A value of 0.01 is good for many typical
24// simulations.
25// ovfactor = Over-allocation factor for the buffers used in the communication. Typical
26// values are 1.0, 2.0. Value 1.0 means no over-allocation.
27// Examples:
28// srun ./LandauDampingParameterList 64 64 64 10000 10 PCG chebyshev 31 0.01 --overallocate 2.0
29// --info 10 srun ./LandauDampingParameterList 64 64 64 10000 10 PCG gauss-seidel 5 2 0 0.01
30// --overallocate 2.0 --info 10
31
32#include <Kokkos_MathematicalConstants.hpp>
33#include <Kokkos_MathematicalFunctions.hpp>
34#include <Kokkos_Random.hpp>
35#include <chrono>
36#include <iostream>
37#include <random>
38#include <set>
39#include <string>
40#include <vector>
41
42#include "Utility/IpplTimings.h"
43
44#include "ChargedParticles.hpp"
45
46constexpr unsigned Dim = 3;
47
48template <typename T>
49struct Newton1D {
50 double tol = 1e-12;
51 int max_iter = 20;
52 double pi = Kokkos::numbers::pi_v<double>;
53
54 T k, alpha, u;
55
56 KOKKOS_INLINE_FUNCTION Newton1D() {}
57
58 KOKKOS_INLINE_FUNCTION Newton1D(const T& k_, const T& alpha_, const T& u_)
59 : k(k_)
60 , alpha(alpha_)
61 , u(u_) {}
62
63 KOKKOS_INLINE_FUNCTION ~Newton1D() {}
64
65 KOKKOS_INLINE_FUNCTION T f(T& x) {
66 T F;
67 F = x + (alpha * (Kokkos::sin(k * x) / k)) - u;
68 return F;
69 }
70
71 KOKKOS_INLINE_FUNCTION T fprime(T& x) {
72 T Fprime;
73 Fprime = 1 + (alpha * Kokkos::cos(k * x));
74 return Fprime;
75 }
76
77 KOKKOS_FUNCTION
78 void solve(T& x) {
79 int iterations = 0;
80 while (iterations < max_iter && Kokkos::fabs(f(x)) > tol) {
81 x = x - (f(x) / fprime(x));
82 iterations += 1;
83 }
84 }
85};
86
87template <typename T, class GeneratorPool, unsigned Dim>
88struct generate_random {
90 using value_type = typename T::value_type;
91 // Output View for the random numbers
92 view_type x, v;
93
94 // The GeneratorPool
95 GeneratorPool rand_pool;
96
98
99 T k, minU, maxU;
100
101 // Initialize all members
102 generate_random(view_type x_, view_type v_, GeneratorPool rand_pool_, value_type& alpha_, T& k_,
103 T& minU_, T& maxU_)
104 : x(x_)
105 , v(v_)
106 , rand_pool(rand_pool_)
107 , alpha(alpha_)
108 , k(k_)
109 , minU(minU_)
110 , maxU(maxU_) {}
111
112 KOKKOS_INLINE_FUNCTION void operator()(const size_t i) const {
113 // Get a random number state from the pool for the active thread
114 typename GeneratorPool::generator_type rand_gen = rand_pool.get_state();
115
116 value_type u;
117 for (unsigned d = 0; d < Dim; ++d) {
118 u = rand_gen.drand(minU[d], maxU[d]);
119 x(i)[d] = u / (1 + alpha);
120 Newton1D<value_type> solver(k[d], alpha, u);
121 solver.solve(x(i)[d]);
122 v(i)[d] = rand_gen.normal(0.0, 1.0);
123 }
124
125 // Give the state back, which will allow another thread to acquire it
126 rand_pool.free_state(rand_gen);
127 }
128};
129
130double CDF(const double& x, const double& alpha, const double& k) {
131 double cdf = x + (alpha / k) * std::sin(k * x);
132 return cdf;
133}
134
135KOKKOS_FUNCTION
136double PDF(const Vector_t<double, Dim>& xvec, const double& alpha, const Vector_t<double, Dim>& kw,
137 const unsigned Dim) {
138 double pdf = 1.0;
139
140 for (unsigned d = 0; d < Dim; ++d) {
141 pdf *= (1.0 + alpha * Kokkos::cos(kw[d] * xvec[d]));
142 }
143 return pdf;
144}
145
146const char* TestName = "LandauDamping";
147
148int main(int argc, char* argv[]) {
149 ippl::initialize(argc, argv);
150 {
152
153 Inform msg("LandauDamping");
154 Inform msg2all("LandauDamping", INFORM_ALL_NODES);
155
156 auto start = std::chrono::high_resolution_clock::now();
157
158 int arg = 1;
159
161 for (unsigned d = 0; d < Dim; d++) {
162 nr[d] = std::atoi(argv[arg++]);
163 }
164
165 static IpplTimings::TimerRef mainTimer = IpplTimings::getTimer("total");
166 static IpplTimings::TimerRef particleCreation = IpplTimings::getTimer("particlesCreation");
167 static IpplTimings::TimerRef dumpDataTimer = IpplTimings::getTimer("dumpData");
168 static IpplTimings::TimerRef PTimer = IpplTimings::getTimer("pushVelocity");
169 static IpplTimings::TimerRef RTimer = IpplTimings::getTimer("pushPosition");
170 static IpplTimings::TimerRef updateTimer = IpplTimings::getTimer("update");
171 static IpplTimings::TimerRef DummySolveTimer = IpplTimings::getTimer("solveWarmup");
172 static IpplTimings::TimerRef SolveTimer = IpplTimings::getTimer("solve");
173 static IpplTimings::TimerRef domainDecomposition = IpplTimings::getTimer("loadBalance");
174
175 IpplTimings::startTimer(mainTimer);
176
177 const size_type totalP = std::atoll(argv[arg++]);
178 const unsigned int nt = std::atoi(argv[arg++]);
179
180 msg << "Landau damping" << endl
181 << "nt " << nt << " Np= " << totalP << " grid = " << nr << endl;
182
183 using bunch_type = ChargedParticles<PLayout_t<double, Dim>, double, Dim>;
184
185 std::unique_ptr<bunch_type> P;
186
187 ippl::NDIndex<Dim> domain;
188 for (unsigned i = 0; i < Dim; i++) {
189 domain[i] = ippl::Index(nr[i]);
190 }
191
192 std::array<bool, Dim> isParallel;
193 isParallel.fill(true);
194
195 // create mesh and layout objects for this problem domain
196 Vector_t<double, Dim> kw = 0.5;
197 double alpha = 0.05;
198 Vector_t<double, Dim> rmin(0.0);
199 Vector_t<double, Dim> rmax = 2 * pi / kw;
200
201 Vector_t<double, Dim> hr = rmax / nr;
202 // Q = -\int\int f dx dv
203 double Q = std::reduce(rmax.begin(), rmax.end(), -1., std::multiplies<double>());
204 Vector_t<double, Dim> origin = rmin;
205 const double dt = std::min(.05, 0.5 * *std::min_element(hr.begin(), hr.end()));
206
207 const bool isAllPeriodic = true;
208 Mesh_t<Dim> mesh(domain, hr, origin);
209 FieldLayout_t<Dim> FL(MPI_COMM_WORLD, domain, isParallel, isAllPeriodic);
210 PLayout_t<double, Dim> PL(FL, mesh);
211
212 std::string solver = argv[arg++];
213 // MOD: Setup of Preconditioner Start
214 int gauss_seidel_inner_iterations;
215 int gauss_seidel_outer_iterations;
216 int newton_level;
217 int chebyshev_degree;
218 int richardson_iterations;
219 int communication;
220 double ssor_omega;
221 std::string preconditioner_type = "";
222
223 if (solver == "OPEN") {
224 throw IpplException("LandauDamping",
225 "Open boundaries solver incompatible with this simulation!");
226 }
227 ippl::ParameterList params;
228 if (solver == "PCG") {
229 preconditioner_type = argv[arg++];
230 if (preconditioner_type == "newton") {
231 newton_level = std::atoi(argv[arg++]);
232 } else if (preconditioner_type == "chebyshev") {
233 chebyshev_degree = std::atoi(argv[arg++]);
234 } else if (preconditioner_type == "richardson") {
235 richardson_iterations = std::atoi(argv[arg++]);
236 communication = std::atoi(argv[arg++]);
237 } else if (preconditioner_type == "gauss-seidel") {
238 gauss_seidel_inner_iterations = std::atoi(argv[arg++]);
239 gauss_seidel_outer_iterations = std::atoi(argv[arg++]);
240 communication = std::atoi(argv[arg++]);
241 } else if (preconditioner_type == "ssor") {
242 gauss_seidel_inner_iterations = std::atoi(argv[arg++]);
243 gauss_seidel_outer_iterations = std::atoi(argv[arg++]);
244 ssor_omega = std::stod(argv[arg++]);
245 }
246
247 params.add("preconditioner_type", preconditioner_type);
248 params.add("gauss_seidel_inner_iterations", gauss_seidel_inner_iterations);
249 params.add("gauss_seidel_outer_iterations", gauss_seidel_outer_iterations);
250 params.add("newton_level", newton_level);
251 params.add("chebyshev_degree", chebyshev_degree);
252 params.add("richardson_iterations", richardson_iterations);
253 params.add("communication", communication);
254 params.add("ssor_omega", ssor_omega);
255 }
256
257 P = std::make_unique<bunch_type>(PL, hr, rmin, rmax, isParallel, Q, solver);
258 // MOD END
259
260 P->nr_m = nr;
261
262 P->initializeFields(mesh, FL);
263
264 P->initSolver(params);
265 P->time_m = 0.0;
266 P->loadbalancethreshold_m = std::atof(argv[arg++]);
267
268 bool isFirstRepartition;
269
270 if ((P->loadbalancethreshold_m != 1.0) && (ippl::Comm->size() > 1)) {
271 msg << "Starting first repartition" << endl;
272 IpplTimings::startTimer(domainDecomposition);
273 isFirstRepartition = true;
274 const ippl::NDIndex<Dim>& lDom = FL.getLocalNDIndex();
275 const int nghost = P->rho_m.getNghost();
276 auto rhoview = P->rho_m.getView();
277
278 using index_array_type = typename ippl::RangePolicy<Dim>::index_array_type;
280 "Assign initial rho based on PDF", P->rho_m.getFieldRangePolicy(),
281 KOKKOS_LAMBDA(const index_array_type& args) {
282 // local to global index conversion
283 Vector_t<double, Dim> xvec = (args + lDom.first() - nghost + 0.5) * hr + origin;
284
285 // ippl::apply accesses the view at the given indices and obtains a
286 // reference; see src/Expression/IpplOperations.h
287 ippl::apply(rhoview, args) = PDF(xvec, alpha, kw, Dim);
288 });
289
290 Kokkos::fence();
291
292 P->initializeORB(FL, mesh);
293 P->repartition(FL, mesh, isFirstRepartition);
294 IpplTimings::stopTimer(domainDecomposition);
295 }
296
297 msg << "First domain decomposition done" << endl;
298 IpplTimings::startTimer(particleCreation);
299
300 typedef ippl::detail::RegionLayout<double, Dim, Mesh_t<Dim>>::uniform_type RegionLayout_t;
301 const RegionLayout_t& RLayout = PL.getRegionLayout();
302 const typename RegionLayout_t::host_mirror_type Regions = RLayout.gethLocalRegions();
303 Vector_t<double, Dim> Nr, Dr, minU, maxU;
304 int myRank = ippl::Comm->rank();
305 double factor = 1;
306 for (unsigned d = 0; d < Dim; ++d) {
307 Nr[d] = CDF(Regions(myRank)[d].max(), alpha, kw[d])
308 - CDF(Regions(myRank)[d].min(), alpha, kw[d]);
309 Dr[d] = CDF(rmax[d], alpha, kw[d]) - CDF(rmin[d], alpha, kw[d]);
310 minU[d] = CDF(Regions(myRank)[d].min(), alpha, kw[d]);
311 maxU[d] = CDF(Regions(myRank)[d].max(), alpha, kw[d]);
312 factor *= Nr[d] / Dr[d];
313 }
314
315 size_type nloc = (size_type)(factor * totalP);
316 size_type Total_particles = 0;
317
318 MPI_Allreduce(&nloc, &Total_particles, 1, MPI_UNSIGNED_LONG, MPI_SUM,
319 ippl::Comm->getCommunicator());
320
321 int rest = (int)(totalP - Total_particles);
322
323 if (ippl::Comm->rank() < rest) {
324 ++nloc;
325 }
326
327 P->create(nloc);
328 Kokkos::Random_XorShift64_Pool<> rand_pool64((size_type)(42 + 100 * ippl::Comm->rank()));
329 Kokkos::parallel_for(
330 nloc, generate_random<Vector_t<double, Dim>, Kokkos::Random_XorShift64_Pool<>, Dim>(
331 P->R.getView(), P->P.getView(), rand_pool64, alpha, kw, minU, maxU));
332
333 Kokkos::fence();
334 ippl::Comm->barrier();
335 IpplTimings::stopTimer(particleCreation);
336
337 P->q = P->Q_m / totalP;
338 msg << "particles created and initial conditions assigned " << endl;
339 isFirstRepartition = false;
340 // The update after the particle creation is not needed as the
341 // particles are generated locally
342
343 IpplTimings::startTimer(DummySolveTimer);
344 P->rho_m = 0.0;
345 P->runSolver();
346 IpplTimings::stopTimer(DummySolveTimer);
347
348 P->scatterCIC(totalP, 0, hr);
349
350 IpplTimings::startTimer(SolveTimer);
351 P->runSolver();
352 IpplTimings::stopTimer(SolveTimer);
353
354 P->gatherCIC();
355
356 IpplTimings::startTimer(dumpDataTimer);
357 P->dumpLandau();
358 P->gatherStatistics(totalP);
359 // P->dumpLocalDomains(FL, 0);
360 IpplTimings::stopTimer(dumpDataTimer);
361
362 // begin main timestep loop
363 msg << "Starting iterations ..." << endl;
364 for (unsigned int it = 0; it < nt; it++) {
365 // LeapFrog time stepping https://en.wikipedia.org/wiki/Leapfrog_integration
366 // Here, we assume a constant charge-to-mass ratio of -1 for
367 // all the particles hence eliminating the need to store mass as
368 // an attribute
369 // kick
370
372 P->P = P->P - 0.5 * dt * P->E;
374
375 // drift
377 P->R = P->R + dt * P->P;
379 // P->R.print();
380
381 // Since the particles have moved spatially update them to correct processors
382 IpplTimings::startTimer(updateTimer);
383 P->update();
384 IpplTimings::stopTimer(updateTimer);
385
386 // Domain Decomposition
387 if (P->balance(totalP, it + 1)) {
388 msg << "Starting repartition" << endl;
389 IpplTimings::startTimer(domainDecomposition);
390 P->repartition(FL, mesh, isFirstRepartition);
391 IpplTimings::stopTimer(domainDecomposition);
392 // IpplTimings::startTimer(dumpDataTimer);
393 // P->dumpLocalDomains(FL, it+1);
394 // IpplTimings::stopTimer(dumpDataTimer);
395 }
396
397 // scatter the charge onto the underlying grid
398 P->scatterCIC(totalP, it + 1, hr);
399
400 // Field solve
401 IpplTimings::startTimer(SolveTimer);
402 P->runSolver();
403 IpplTimings::stopTimer(SolveTimer);
404
405 // gather E field
406 P->gatherCIC();
407
408 // kick
410 P->P = P->P - 0.5 * dt * P->E;
412
413 P->time_m += dt;
414 IpplTimings::startTimer(dumpDataTimer);
415 P->dumpLandau();
416 P->gatherStatistics(totalP);
417 IpplTimings::stopTimer(dumpDataTimer);
418 msg << "Finished time step: " << it + 1 << " time: " << P->time_m << endl;
419
420 if (checkSignalHandler()) {
421 msg << "Aborting timestepping loop due to signal " << interruptSignalReceived
422 << endl;
423 break;
424 }
425 }
426
427 msg << "LandauDamping: End." << endl;
428 IpplTimings::stopTimer(mainTimer);
430 IpplTimings::print(std::string("timing.dat"));
431 auto end = std::chrono::high_resolution_clock::now();
432
433 std::chrono::duration<double> time_chrono =
434 std::chrono::duration_cast<std::chrono::duration<double>>(end - start);
435 std::cout << "Elapsed time: " << time_chrono.count() << std::endl;
436 }
438
439 return 0;
440}
constexpr unsigned Dim
const double pi
bool checkSignalHandler()
int interruptSignalReceived
void setSignalHandler()
int main(int argc, char *argv[])
KOKKOS_FUNCTION double PDF(const Vector_t< double, Dim > &xvec, const double &alpha, const Vector_t< double, Dim > &kw, const unsigned Dim)
double CDF(const double &x, const double &alpha, const double &k)
constexpr unsigned Dim
ippl::detail::size_type size_type
Definition datatypes.h:23
const char * TestName
ippl::FieldLayout< Dim > FieldLayout_t
Definition datatypes.h:21
ippl::Vector< T, Dim > Vector_t
Definition datatypes.h:38
typename ippl::ParticleSpatialLayout< T, Dim, Mesh_t< Dim > > PLayout_t
Definition datatypes.h:15
ippl::UniformCartesian< double, Dim > Mesh_t
Definition datatypes.h:12
Inform & endl(Inform &inf)
Definition Inform.cpp:42
#define INFORM_ALL_NODES
Definition Inform.h:38
void initialize(int &argc, char *argv[], MPI_Comm comm)
Definition Ippl.cpp:16
void finalize()
Definition Ippl.cpp:94
KOKKOS_INLINE_FUNCTION constexpr decltype(auto) apply(const View &view, const Coords &coords)
void parallel_for(const std::string &name, const ExecPolicy &policy, const FunctorType &functor)
std::unique_ptr< mpi::Communicator > Comm
Definition Ippl.h:22
const NDIndex_t & getLocalNDIndex() const
KOKKOS_INLINE_FUNCTION Vector< int, Dim > first() const
Definition NDIndex.hpp:170
const host_mirror_type gethLocalRegions() const
KOKKOS_INLINE_FUNCTION constexpr iterator begin()
Definition Vector.hpp:160
KOKKOS_INLINE_FUNCTION constexpr iterator end()
Definition Vector.hpp:165
Kokkos::View< typename NPtr< T, Dim >::type, Properties... > view_type
Definition ViewTypes.h:45
Timing::TimerRef TimerRef
static TimerRef getTimer(const char *nm)
static void stopTimer(TimerRef t)
static void print()
static void startTimer(TimerRef t)
::ippl::Vector< index_type, Dim > index_array_type
void add(const std::string &key, const T &value)
KOKKOS_FUNCTION void solve(T &x)
KOKKOS_INLINE_FUNCTION Newton1D(const T &k_, const T &alpha_, const T &u_)
KOKKOS_INLINE_FUNCTION T f(T &x)
KOKKOS_INLINE_FUNCTION ~Newton1D()
KOKKOS_INLINE_FUNCTION Newton1D()
KOKKOS_INLINE_FUNCTION T fprime(T &x)
KOKKOS_INLINE_FUNCTION void operator()(const size_t i) const
typename T::value_type value_type
generate_random(view_type x_, view_type v_, GeneratorPool rand_pool_, value_type &alpha_, T &k_, T &minU_, T &maxU_)
typename ippl::detail::ViewType< T, 1 >::view_type view_type