OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
FieldSolver.cpp
Go to the documentation of this file.
1//
2// Class FieldSolver
3// The class for the OPAL FIELDSOLVER command.
4// A FieldSolver definition is used by most physics commands to define the
5// particle charge and the reference momentum, together with some other data.
6//
7// Copyright (c) 200x - 2022, Paul Scherrer Institut, Villigen PSI, Switzerland
8//
9// All rights reserved
10//
11// This file is part of OPAL.
12//
13// OPAL is free software: you can redistribute it and/or modify
14// it under the terms of the GNU General Public License as published by
15// the Free Software Foundation, either version 3 of the License, or
16// (at your option) any later version.
17//
18// You should have received a copy of the GNU General Public License
19// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
20//
22
30#include "Physics/Physics.h"
33#ifdef HAVE_SAAMG_SOLVER
35#endif
39
40#ifdef ENABLE_AMR
41 #include "Amr/AmrBoxLib.h"
42#ifdef AMREX_ENABLE_FBASELIB
44#endif
46 #include "Amr/AmrDefs.h"
48#endif
49
50#ifdef HAVE_AMR_MG_SOLVER
52#endif
53
54#include <map>
55
56using namespace Expressions;
57
58//TODO: o add a FIELD for DISCRETIZATION, MAXITERS, TOL...
59
60// The attributes of class FieldSolver.
61namespace {
62 namespace deprecated {
63 enum {
64 BCFFTT,
65 SIZE
66 };
67 }
68 enum {
69 FSTYPE = deprecated::SIZE, // The field solver name
70 // FOR FFT BASED SOLVER
71 MX, // mesh sixe in x
72 MY, // mesh sixe in y
73 MT, // mesh sixe in z
74 PARFFTX, // parallelized grid in x
75 PARFFTY, // parallelized grid in y
76 PARFFTT, // parallelized grid in z
77 BCFFTX, // boundary condition in x [FFT + AMR_MG only]
78 BCFFTY, // boundary condition in y [FFT + AMR_MG only]
79 BCFFTZ, // boundary condition in z [FFT + AMR_MG only]
80 GREENSF, // holds greensfunction to be used [FFT + P3M only]
81 BBOXINCR, // how much the boundingbox is increased
82 GEOMETRY, // geometry of boundary [SAAMG only]
83 ITSOLVER, // iterative solver [SAAMG + AMR_MG]
84 INTERPL, // interpolation used for boundary points [SAAMG only]
85 TOL, // tolerance of the SAAMG preconditioned solver [SAAMG only]
86 MAXITERS, // max number of iterations [SAAMG only]
87 PRECMODE, // preconditioner mode [SAAMG only]
88 RC, // cutoff radius for PP interactions
89 ALPHA, // Green’s function splitting parameter
90#ifdef ENABLE_AMR
91 AMR_MAXLEVEL, // AMR, maximum refinement level
92 AMR_REFX, // AMR, refinement ratio in x
93 AMR_REFY, // AMR, refinement ratio in y
94 AMR_REFZ, // AMR, refinement ration in z
95 AMR_MAXGRIDX, // AMR, maximum grid size in x (default: 16)
96 AMR_MAXGRIDY, // AMR, maximum grid size in y (default: 16)
97 AMR_MAXGRIDZ, // AMR, maximum grid size in z (default: 16)
98 AMR_BFX, // AMR, blocking factor in x (maxgrid needs to be a multiple, default: 8)
99 AMR_BFY, // AMR, blocking factor in y (maxgrid needs to be a multiple, default: 8)
100 AMR_BFZ, // AMR, blocking factor in z (maxgrid needs to be a multiple, default: 8)
101 AMR_TAGGING,
102 AMR_DENSITY,
103 AMR_MAX_NUM_PART,
104 AMR_MIN_NUM_PART,
105 AMR_SCALING,
106 AMR_DOMAIN_RATIO,
107#endif
108#ifdef HAVE_AMR_MG_SOLVER
109 // AMR_MG = Adaptive-Mesh-Refinement Multi-Grid
110 AMR_MG_SMOOTHER, // AMR, smoother for level solution
111 AMR_MG_NSWEEPS, // AMR, number of smoothing sweeps
112 AMR_MG_PREC, // AMR, preconditioner for bottom solver
113 AMR_MG_INTERP, // AMR, interpolater for solution from level l to l+1
114 AMR_MG_NORM, // AMR, norm convergence criteria
115 AMR_MG_VERBOSE, // AMR, enable solver info writing (SDDS file)
116 AMR_MG_REBALANCE, // AMR, rebalance smoothed aggregation (SA) preconditioner
117 AMR_MG_REUSE, // AMR, reuse type of SA (NONE, RP, RAP, SYMBOLIC or FULL)
118 AMR_MG_TOL, // AMR, tolerance of solver
119#endif
120 // FOR XXX BASED SOLVER
121 SIZE
122 };
123}
124
125
127 Definition(SIZE, "FIELDSOLVER",
128 "The \"FIELDSOLVER\" statement defines data for a the field solver") {
129
130 itsAttr[FSTYPE] = Attributes::makePredefinedString("FSTYPE", "Name of the attached field solver.",
131 {"FFT", "FFTPERIODIC", "SAAMG", "FMG", "ML", "AMR_MG", "NONE","P3M"});
132
133 itsAttr[MX] = Attributes::makeReal("MX", "Meshsize in x");
134 itsAttr[MY] = Attributes::makeReal("MY", "Meshsize in y");
135 itsAttr[MT] = Attributes::makeReal("MT", "Meshsize in z(t)");
136
137 itsAttr[PARFFTX] = Attributes::makeBool("PARFFTX",
138 "True, dimension 0 i.e x is parallelized",
139 false);
140
141 itsAttr[PARFFTY] = Attributes::makeBool("PARFFTY",
142 "True, dimension 1 i.e y is parallelized",
143 false);
144
145 itsAttr[PARFFTT] = Attributes::makeBool("PARFFTT",
146 "True, dimension 2 i.e z(t) is parallelized",
147 true);
148
149 //FFT ONLY:
150 itsAttr[BCFFTX] = Attributes::makePredefinedString("BCFFTX",
151 "Boundary conditions in x.",
152 {"OPEN", "DIRICHLET", "PERIODIC"},
153 "OPEN");
154
155 itsAttr[BCFFTY] = Attributes::makePredefinedString("BCFFTY",
156 "Boundary conditions in y.",
157 {"OPEN", "DIRICHLET", "PERIODIC"},
158 "OPEN");
159
160 itsAttr[BCFFTZ] = Attributes::makePredefinedString("BCFFTZ",
161 "Boundary conditions in z(t).",
162 {"OPEN", "DIRICHLET", "PERIODIC"},
163 "OPEN");
164
165 itsAttr[deprecated::BCFFTT] = Attributes::makePredefinedString("BCFFTT",
166 "Boundary conditions in z(t).",
167 {"OPEN", "DIRICHLET", "PERIODIC"},
168 "OPEN");
169
170 //GREENSF is also used in P3M
171 itsAttr[GREENSF] = Attributes::makePredefinedString("GREENSF",
172 "Which Greensfunction to be used.",
173 {"STANDARD", "INTEGRATED"},
174 "INTEGRATED");
175
176 itsAttr[BBOXINCR] = Attributes::makeReal("BBOXINCR",
177 "Increase of bounding box in % ",
178 2.0);
179
180 // P3M only:
181 itsAttr[RC] = Attributes::makeReal("RC",
182 "cutoff radius for PP interactions",
183 0.0);
184
186 "Standard Ewald Green’s function splitting parameter",
187 1e8);
188 //SAAMG and in case of FFT with dirichlet BC in x and y
189 itsAttr[GEOMETRY] = Attributes::makeUpperCaseString("GEOMETRY",
190 "GEOMETRY to be used as domain boundary",
191 "");
192
193 itsAttr[ITSOLVER] = Attributes::makePredefinedString("ITSOLVER",
194 "Type of iterative solver.",
195 {"CG",
196 "BICGSTAB",
197 "GMRES",
198 "MINRES",
199 "PCPG",
200 "STOCHASTIC_CG",
201 "RECYCLING_CG",
202 "RECYCLING_GMRES",
203 "KLU2",
204 "SUPERLU",
205 "UMFPACK",
206 "PARDISO_MKL",
207 "MUMPS",
208 "LAPACK",
209 "SA"},
210 "CG");
211
212 itsAttr[INTERPL] = Attributes::makePredefinedString("INTERPL",
213 "interpolation used for boundary points.",
214 {"CONSTANT", "LINEAR", "QUADRATIC"},
215 "LINEAR");
216
217 itsAttr[TOL] = Attributes::makeReal("TOL",
218 "Tolerance for iterative solver",
219 1e-8);
220
221 itsAttr[MAXITERS] = Attributes::makeReal("MAXITERS",
222 "Maximum number of iterations of iterative solver",
223 100);
224
225 itsAttr[PRECMODE] = Attributes::makePredefinedString("PRECMODE",
226 "Preconditioner Mode.",
227 {"STD", "HIERARCHY", "REUSE"},
228 "HIERARCHY");
229
230 // AMR
231#ifdef ENABLE_AMR
232 itsAttr[AMR_MAXLEVEL] = Attributes::makeReal("AMR_MAXLEVEL",
233 "Maximum number of levels in AMR",
234 0);
235
236 itsAttr[AMR_REFX] = Attributes::makeReal("AMR_REFX",
237 "Refinement ration in x-direction in AMR",
238 2);
239
240 itsAttr[AMR_REFY] = Attributes::makeReal("AMR_REFY",
241 "Refinement ration in y-direction in AMR",
242 2);
243
244 itsAttr[AMR_REFZ] = Attributes::makeReal("AMR_REFZ",
245 "Refinement ration in z-direction in AMR",
246 2);
247
248 itsAttr[AMR_MAXGRIDX] = Attributes::makeReal("AMR_MAXGRIDX",
249 "Maximum grid size in x for AMR",
250 16);
251
252 itsAttr[AMR_MAXGRIDY] = Attributes::makeReal("AMR_MAXGRIDY",
253 "Maximum grid size in y for AMR",
254 16);
255
256 itsAttr[AMR_MAXGRIDZ] = Attributes::makeReal("AMR_MAXGRIDZ",
257 "Maximum grid size in z for AMR",
258 16);
259
260 itsAttr[AMR_BFX] = Attributes::makeReal("AMR_BFX",
261 "Blocking factor in x for AMR (AMR_MAXGRIDX needs to be a multiple)",
262 8);
263
264 itsAttr[AMR_BFY] = Attributes::makeReal("AMR_BFY",
265 "Blocking factor in y for AMR (AMR_MAXGRIDY needs to be a multiple)",
266 8);
267
268 itsAttr[AMR_BFZ] = Attributes::makeReal("AMR_BFZ",
269 "Blocking factor in y for AMR (AMR_MAXGRIDZ needs to be a multiple)",
270 8);
271
272 itsAttr[AMR_TAGGING] = Attributes::makeUpperCaseString("AMR_TAGGING",
273 "Refinement criteria [CHARGE_DENSITY | POTENTIAL | EFIELD | MIN_NUM_PARTICLES | MAX_NUM_PARTICLES]",
274 "CHARGE_DENSITY");
275
276 itsAttr[AMR_DENSITY] = Attributes::makeReal("AMR_DENSITY",
277 "Tagging value for charge density refinement [C / cell volume]",
278 1.0e-14);
279
280 itsAttr[AMR_MAX_NUM_PART] = Attributes::makeReal("AMR_MAX_NUM_PART",
281 "Tagging value for max. #particles",
282 1);
283
284 itsAttr[AMR_MIN_NUM_PART] = Attributes::makeReal("AMR_MIN_NUM_PART",
285 "Tagging value for min. #particles",
286 1);
287
288 itsAttr[AMR_SCALING] = Attributes::makeReal("AMR_SCALING",
289 "Scaling value for maximum value tagging "
290 "(only POTENTIAL / CHARGE_DENSITY / "
291 "MOMENTA)", 0.75);
292
293 itsAttr[AMR_DOMAIN_RATIO] = Attributes::makeRealArray("AMR_DOMAIN_RATIO",
294 "Box ratio of AMR computation domain. Default: [-1, 1]^3");
295
296 // default
297 Attributes::setRealArray(itsAttr[AMR_DOMAIN_RATIO], {1.0, 1.0, 1.0});
298#endif
299
300#ifdef HAVE_AMR_MG_SOLVER
301 itsAttr[AMR_MG_SMOOTHER] = Attributes::makePredefinedString("AMR_MG_SMOOTHER",
302 "Smoothing of level solution.",
303 {"GS", "SGS", "JACOBI"}, "GS");
304
305 itsAttr[AMR_MG_NSWEEPS] = Attributes::makeReal("AMR_MG_NSWEEPS",
306 "Number of relaxation steps",
307 8);
308
309 itsAttr[AMR_MG_PREC] = Attributes::makePredefinedString("AMR_MG_PREC",
310 "Preconditioner of bottom solver.",
311 {"NONE", "ILUT", "CHEBYSHEV", "RILUK", "JACOBI", "BLOCK_JACOBI", "GS", "BLOCK_GS", "SA"},
312 "NONE");
313
314 itsAttr[AMR_MG_INTERP] = Attributes::makePredefinedString("AMR_MG_INTERP",
315 "Interpolater between levels.",
316 {"TRILINEAR", "LAGRANGE", "PC"},
317 "PC");
318
319 itsAttr[AMR_MG_NORM] = Attributes::makePredefinedString("AMR_MG_NORM",
320 "Norm for convergence criteria.",
321 {"L1_NORM", "L2_NORM", "LINF_NORM"},
322 "LINF_NORM");
323
324 itsAttr[AMR_MG_VERBOSE] = Attributes::makeBool("AMR_MG_VERBOSE",
325 "Write solver info in SDDS format (*.solver)",
326 false);
327
328 itsAttr[AMR_MG_REBALANCE] = Attributes::makeBool("AMR_MG_REBALANCE",
329 "Rebalancing of Smoothed Aggregation "
330 "Preconditioner",
331 false);
332
333 itsAttr[AMR_MG_REUSE] = Attributes::makePredefinedString("AMR_MG_REUSE",
334 "Reuse type of Smoothed Aggregation.",
335 {"NONE", "RP", "RAP", "SYMBOLIC", "FULL"},
336 "RAP");
337
338 itsAttr[AMR_MG_TOL] = Attributes::makeReal("AMR_MG_TOL",
339 "AMR MG solver tolerance (default: 1.0e-10)",
340 1.0e-10);
341#endif
342
343 mesh_m = 0;
344 FL_m = 0;
345 PL_m.reset(nullptr);
346 solver_m = 0;
347
349}
350
351
352FieldSolver::FieldSolver(const std::string& name, FieldSolver* parent):
353 Definition(name, parent)
354{
355 mesh_m = 0;
356 FL_m = 0;
357 PL_m.reset(nullptr);
358 solver_m = 0;
359}
360
361
363 if (mesh_m) {
364 delete mesh_m;
365 mesh_m = 0;
366 }
367 if (FL_m) {
368 delete FL_m;
369 FL_m = 0;
370 }
371 if (solver_m) {
372 delete solver_m;
373 solver_m = 0;
374 }
375}
376
377FieldSolver* FieldSolver::clone(const std::string& name) {
378 return new FieldSolver(name, this);
379}
380
385
386FieldSolver* FieldSolver::find(const std::string& name) {
387 FieldSolver* fs = dynamic_cast<FieldSolver*>(OpalData::getInstance()->find(name));
388
389 if (fs == 0) {
390 throw OpalException("FieldSolver::find()", "FieldSolver \"" + name + "\" not found.");
391 }
392 return fs;
393}
394
395std::string FieldSolver::getType() {
396 return Attributes::getString(itsAttr[FSTYPE]);
397}
398
399double FieldSolver::getMX() const {
400 return Attributes::getReal(itsAttr[MX]);
401}
402
403double FieldSolver::getMY() const {
404 return Attributes::getReal(itsAttr[MY]);
405}
406
407double FieldSolver::getMT() const {
408 return Attributes::getReal(itsAttr[MT]);
409}
410
411void FieldSolver::setMX(double value) {
412 Attributes::setReal(itsAttr[MX], value);
413}
414
415void FieldSolver::setMY(double value) {
416 Attributes::setReal(itsAttr[MY], value);
417}
418
419void FieldSolver::setMT(double value) {
420 Attributes::setReal(itsAttr[MT], value);
421}
422
424
425}
426
428
429 e_dim_tag decomp[3] = {SERIAL, SERIAL, SERIAL};
430
431 NDIndex<3> domain;
432 domain[0] = Index((int)getMX() + 1);
433 domain[1] = Index((int)getMY() + 1);
434 domain[2] = Index((int)getMT() + 1);
435
436 if (Attributes::getBool(itsAttr[PARFFTX]))
437 decomp[0] = PARALLEL;
438 if (Attributes::getBool(itsAttr[PARFFTY]))
439 decomp[1] = PARALLEL;
440 if (Attributes::getBool(itsAttr[PARFFTT]))
441 decomp[2] = PARALLEL;
442
443 if (Attributes::getString(itsAttr[FSTYPE]) == "FFTPERIODIC") {
444 decomp[0] = decomp[1] = SERIAL;
445 decomp[2] = PARALLEL;
446 }
447 // create prototype mesh and layout objects for this problem domain
448
449 if ( !isAmrSolverType() ) {
450 mesh_m = new Mesh_t(domain);
451 FL_m = new FieldLayout_t(*mesh_m, decomp);
452 PL_m.reset(new Layout_t(*FL_m, *mesh_m));
453 // OpalData::getInstance()->setMesh(mesh_m);
454 // OpalData::getInstance()->setFieldLayout(FL_m);
455 // OpalData::getInstance()->setLayout(PL_m);
456 }
457}
458
460 if (itsAttr[deprecated::BCFFTT])
461 return (Attributes::getString(itsAttr[deprecated::BCFFTT]) == "PERIODIC");
462
463 return (Attributes::getString(itsAttr[BCFFTZ]) == "PERIODIC");
464}
465
470 throw OpalException("FieldSolver::isAmrSolverType",
471 "The attribute \"FSTYPE\" has not correct type for AMR solver");
472 } else {
473 return Options::amr;
474 }
475}
476
478 static const std::map<std::string, FieldSolverType> stringFSType_s = {
479 {"NONE", FieldSolverType::NONE},
480 {"FFT", FieldSolverType::FFT},
481 {"FFTBOX", FieldSolverType::FFTBOX},
482 {"SAAMG", FieldSolverType::SAAMG},
483 {"P3M", FieldSolverType::P3M},
484 {"FMG", FieldSolverType::FMG},
485 {"ML", FieldSolverType::ML},
486 {"AMR_MG", FieldSolverType::AMRMG},
487 {"HYPRE", FieldSolverType::HYPRE},
488 {"HPGMG", FieldSolverType::HPGMG}
489 };
490 fsName_m = getType();
491 if (fsName_m.empty()) {
492 throw OpalException("FieldSolver::setFieldSolverType",
493 "The attribute \"FSTYPE\" isn't set for \"FIELDSOLVER\"!");
494 } else {
495 fsType_m = stringFSType_s.at(fsName_m);
496 }
497}
498
500 itsBunch_m = b;
501
502 std::string greens = Attributes::getString(itsAttr[GREENSF]);
503 std::string bcx = Attributes::getString(itsAttr[BCFFTX]);
504 std::string bcy = Attributes::getString(itsAttr[BCFFTY]);
505 std::string bcz = Attributes::getString(itsAttr[deprecated::BCFFTT]);
506 if (bcz.empty()) {
507 bcz = Attributes::getString(itsAttr[BCFFTZ]);
508 }
509
510 if ( isAmrSolverType() ) {
511 Inform m("FieldSolver::initAmrSolver");
512
513#ifdef ENABLE_AMR
515
517#endif
518 } else if (fsType_m == FieldSolverType::FFT) {
519 bool sinTrafo = ((bcx == "DIRICHLET") && (bcy == "DIRICHLET") && (bcz == "DIRICHLET"));
520 if (sinTrafo) {
521 std::cout << "FFTBOX ACTIVE" << std::endl;
522 //we go over all geometries and add the Geometry Elements to the geometry list
523 std::string geoms = Attributes::getString(itsAttr[GEOMETRY]);
524 std::string tmp;
525 //split and add all to list
526 std::vector<BoundaryGeometry*> geometries;
527 for (unsigned int i = 0; i <= geoms.length(); i++) {
528 if (i == geoms.length() || geoms[i] == ',') {
530 if (geom != 0)
531 geometries.push_back(geom);
532 tmp.clear();
533 } else {
534 tmp += geoms[i];
535 }
536 }
537 BoundaryGeometry* ttmp = geometries[0];
538 solver_m = new FFTBoxPoissonSolver(mesh_m, FL_m, greens, ttmp->getA());
539 itsBunch_m->set_meshEnlargement(Attributes::getReal(itsAttr[BBOXINCR]) / 100.0);
540 fsName_m = "FFTBOX";
542 } else {
543 solver_m = new FFTPoissonSolver(mesh_m, FL_m, greens, bcz);
544 itsBunch_m->set_meshEnlargement(Attributes::getReal(itsAttr[BBOXINCR]) / 100.0);
545 }
546 } else if (fsType_m == FieldSolverType::P3M) {
548 FL_m,
551 greens);
552
553 itsBunch_m->set_meshEnlargement(Attributes::getReal(itsAttr[BBOXINCR]) / 100.0);
554
555 } else if (fsType_m == FieldSolverType::SAAMG) {
556#ifdef HAVE_SAAMG_SOLVER
557 //we go over all geometries and add the Geometry Elements to the geometry list
558 std::string geoms = Attributes::getString(itsAttr[GEOMETRY]);
559 std::string tmp;
560 //split and add all to list
561 std::vector<BoundaryGeometry*> geometries;
562 for (unsigned int i = 0; i <= geoms.length(); i++) {
563 if (i == geoms.length() || geoms[i] == ',') {
565 if (geom != 0) {
566 geometries.push_back(geom);
567 }
568 tmp.clear();
569 } else
570 tmp += geoms[i];
571 }
572 solver_m = new MGPoissonSolver(dynamic_cast<PartBunch*>(itsBunch_m), mesh_m, FL_m,
573 geometries,
577 Attributes::getReal(itsAttr[MAXITERS]),
578 Attributes::getString(itsAttr[PRECMODE]));
579 itsBunch_m->set_meshEnlargement(Attributes::getReal(itsAttr[BBOXINCR]) / 100.0);
580#else
581 throw OpalException("FieldSolver::initSolver",
582 "SAAMG Solver not enabled! Please build OPAL with -DENABLE_SAAMG_SOLVER=1");
583#endif
584 } else {
585 solver_m = 0;
586 INFOMSG("No solver attached" << endl);
587 }
588}
589
591 return (solver_m != 0);
592}
593
595 os << "* ************* F I E L D S O L V E R ********************************************** " << endl;
596 os << "* FIELDSOLVER " << getOpalName() << '\n'
597 << "* TYPE " << fsName_m << '\n'
598 << "* N-PROCESSORS " << Ippl::getNodes() << '\n'
599 << "* MX " << Attributes::getReal(itsAttr[MX]) << '\n'
600 << "* MY " << Attributes::getReal(itsAttr[MY]) << '\n'
601 << "* MT " << Attributes::getReal(itsAttr[MT]) << '\n'
602 << "* BBOXINCR " << Attributes::getReal(itsAttr[BBOXINCR]) << endl;
604 os << "* RC " << Attributes::getReal(itsAttr[RC]) << '\n'
605 << "* ALPHA " << Attributes::getReal(itsAttr[ALPHA]) << '\n'
606 << "* GREENSF " << Attributes::getString(itsAttr[GREENSF]) << endl;
607 } else if (fsType_m == FieldSolverType::FFT) {
608 os << "* GREENSF " << Attributes::getString(itsAttr[GREENSF]) << endl;
609 } else if (fsType_m == FieldSolverType::SAAMG) {
610 os << "* GEOMETRY " << Attributes::getString(itsAttr[GEOMETRY]) << '\n'
611 << "* ITSOLVER " << Attributes::getString(itsAttr[ITSOLVER]) << '\n'
612 << "* INTERPL " << Attributes::getString(itsAttr[INTERPL]) << '\n'
613 << "* TOL " << Attributes::getReal(itsAttr[TOL]) << '\n'
614 << "* MAXITERS " << Attributes::getReal(itsAttr[MAXITERS]) << '\n'
615 << "* PRECMODE " << Attributes::getString(itsAttr[PRECMODE]) << endl;
616 } else if (Options::amr) {
617#ifdef ENABLE_AMR
618 os << "* AMR_MAXLEVEL " << Attributes::getReal(itsAttr[AMR_MAXLEVEL]) << '\n'
619 << "* AMR_REFX " << Attributes::getReal(itsAttr[AMR_REFX]) << '\n'
620 << "* AMR_REFY " << Attributes::getReal(itsAttr[AMR_REFY]) << '\n'
621 << "* AMR_REFZ " << Attributes::getReal(itsAttr[AMR_REFZ]) << '\n'
622 << "* AMR_MAXGRIDX " << Attributes::getReal(itsAttr[AMR_MAXGRIDX]) << '\n'
623 << "* AMR_MAXGRIDY " << Attributes::getReal(itsAttr[AMR_MAXGRIDY]) << '\n'
624 << "* AMR_MAXGRIDZ " << Attributes::getReal(itsAttr[AMR_MAXGRIDZ]) << '\n'
625 << "* AMR_BFX " << Attributes::getReal(itsAttr[AMR_BFX]) << '\n'
626 << "* AMR_BFY " << Attributes::getReal(itsAttr[AMR_BFY]) << '\n'
627 << "* AMR_BFZ " << Attributes::getReal(itsAttr[AMR_BFZ]) << '\n'
628 << "* AMR_TAGGING " << this->getTagging_m() <<'\n'
629 << "* AMR_DENSITY " << Attributes::getReal(itsAttr[AMR_DENSITY]) << '\n'
630 << "* AMR_MAX_NUM_PART " << Attributes::getReal(itsAttr[AMR_MAX_NUM_PART]) << '\n'
631 << "* AMR_MIN_NUM_PART " << Attributes::getReal(itsAttr[AMR_MIN_NUM_PART]) << '\n'
632 << "* AMR_DENSITY " << Attributes::getReal(itsAttr[AMR_DENSITY]) << '\n'
633 << "* AMR_SCALING " << Attributes::getReal(itsAttr[AMR_SCALING]) << endl;
634 auto length = Attributes::getRealArray(itsAttr[AMR_DOMAIN_RATIO]);
635 os << "* AMR_DOMAIN_RATIO ( ";
636 for (auto& l : length) {
637 os << l << " ";
638 }
639 os << ")" << endl;
640#endif
641 }
642
643#ifdef HAVE_AMR_MG_SOLVER
645 os << "* ITSOLVER (AMR_MG) "
646 << Attributes::getString(itsAttr[ITSOLVER]) << '\n'
647 << "* AMR_MG_PREC "
648 << Attributes::getString(itsAttr[AMR_MG_PREC]) << '\n'
649 << "* AMR_MG_REBALANCE "
650 << Attributes::getBool(itsAttr[AMR_MG_REBALANCE]) << '\n'
651 << "* AMR_MG_REUSE "
652 << Attributes::getString(itsAttr[AMR_MG_REUSE]) << '\n'
653 << "* AMR_MG_SMOOTHER "
654 << Attributes::getString(itsAttr[AMR_MG_SMOOTHER]) << '\n'
655 << "* AMR_MG_NSWEEPS "
656 << Attributes::getReal(itsAttr[AMR_MG_NSWEEPS]) << '\n'
657 << "* AMR_MG_INTERP "
658 << Attributes::getString(itsAttr[AMR_MG_INTERP]) << '\n'
659 << "* AMR_MG_NORM "
660 << Attributes::getString(itsAttr[AMR_MG_NORM]) << '\n'
661 << "* AMR_MG_TOL "
662 << Attributes::getReal(itsAttr[AMR_MG_TOL]) << '\n'
663 << "* AMR_MG_VERBOSE "
664 << Attributes::getBool(itsAttr[AMR_MG_VERBOSE]) << '\n'
665 << "* BCFFTX "
666 << Attributes::getString(itsAttr[BCFFTX]) << '\n'
667 << "* BCFFTY "
668 << Attributes::getString(itsAttr[BCFFTY]) << '\n';
669 if (itsAttr[deprecated::BCFFTT]) {
670 os << "* BCFFTT (deprec.) "
671 << Attributes::getString(itsAttr[deprecated::BCFFTT]) << endl;
672 } else {
673 os << "* BCFFTZ "
674 << Attributes::getString(itsAttr[BCFFTZ]) << endl;
675 }
676 }
677#endif
678
679 if (Attributes::getBool(itsAttr[PARFFTX])) {
680 os << "* XDIM parallel " << endl;
681 } else {
682 os << "* XDIM serial " << endl;
683 }
684
685 if (Attributes::getBool(itsAttr[PARFFTY])) {
686 os << "* YDIM parallel " << endl;
687 } else {
688 os << "* YDIM serial " << endl;
689 }
690
691 if (Attributes::getBool(itsAttr[PARFFTT])) {
692 os << "* Z(T)DIM parallel " << endl;
693 }else {
694 os << "* Z(T)DIM serial " << endl;
695 }
696
697 if ( !isAmrSolverType() ) {
698 INFOMSG(level3 << *mesh_m << endl);
699 INFOMSG(level3 << *PL_m << endl);
700 }
701
702 if (solver_m)
703 os << *solver_m << endl;
704 os << "* ********************************************************************************** " << endl;
705 return os;
706}
707
708#ifdef ENABLE_AMR
709std::string FieldSolver::getTagging_m() const {
710 std::string tagging = Attributes::getString(itsAttr[AMR_TAGGING]);
711
712 // This conversion was introduced to allow numbers for tagging (useful for the sampler).
713 if ( Util::isAllDigits(tagging) ) {
714 tagging = AmrObject::getTaggingString(std::stoi(tagging));
715 }
716 return tagging;
717}
718
719
721
722 auto domain = Attributes::getRealArray(itsAttr[AMR_DOMAIN_RATIO]);
723 dynamic_cast<AmrPartBunch*>(itsBunch_m)->setAmrDomainRatio(
724 Attributes::getRealArray(itsAttr[AMR_DOMAIN_RATIO])
725 );
726 itsBunch_m->set_meshEnlargement(Attributes::getReal(itsAttr[BBOXINCR]) * 0.01);
727
728 // setup initial info for creating the object
730 info.grid[0] = (int)this->getMX();
731 info.grid[1] = (int)this->getMY();
732 info.grid[2] = (int)this->getMT();
733 info.maxgrid[0] = Attributes::getReal(itsAttr[AMR_MAXGRIDX]);
734 info.maxgrid[1] = Attributes::getReal(itsAttr[AMR_MAXGRIDY]);
735 info.maxgrid[2] = Attributes::getReal(itsAttr[AMR_MAXGRIDZ]);
736 info.bf[0] = Attributes::getReal(itsAttr[AMR_BFX]);
737 info.bf[1] = Attributes::getReal(itsAttr[AMR_BFY]);
738 info.bf[2] = Attributes::getReal(itsAttr[AMR_BFZ]);
739 info.maxlevel = Attributes::getReal(itsAttr[AMR_MAXLEVEL]);
740 info.refratio[0] = Attributes::getReal(itsAttr[AMR_REFX]);
741 info.refratio[1] = Attributes::getReal(itsAttr[AMR_REFY]);
742 info.refratio[2] = Attributes::getReal(itsAttr[AMR_REFZ]);
743
745
746 itsAmrObject_mp->setTagging( this->getTagging_m() );
747
748 itsAmrObject_mp->setScalingFactor( Attributes::getReal(itsAttr[AMR_SCALING]) );
749
750 itsAmrObject_mp->setChargeDensity( Attributes::getReal(itsAttr[AMR_DENSITY]) );
751
752 itsAmrObject_mp->setMaxNumParticles(
753 Attributes::getReal(itsAttr[AMR_MAX_NUM_PART])
754 );
755
756 itsAmrObject_mp->setMinNumParticles(
757 Attributes::getReal(itsAttr[AMR_MIN_NUM_PART])
758 );
759}
760
761
763
765 if ( dynamic_cast<AmrBoxLib*>( itsAmrObject_mp.get() ) == 0 )
766 throw OpalException("FieldSolver::initAmrSolver_m()",
767 "ML solver requires AMReX.");
768 solver_m = new MLPoissonSolver(static_cast<AmrBoxLib*>(itsAmrObject_mp.get()));
769#ifdef AMREX_ENABLE_FBASELIB
770 } else if (fsType_m == FieldSolverType::FMG) {
771
772 if ( dynamic_cast<AmrBoxLib*>( itsAmrObject_mp.get() ) == 0 )
773 throw OpalException("FieldSolver::initAmrSolver_m()",
774 "FMultiGrid solver requires AMReX.");
775
776 solver_m = new FMGPoissonSolver(static_cast<AmrBoxLib*>(itsAmrObject_mp.get()));
777#endif
778 } else if (fsType_m == FieldSolverType::HYPRE) {
779 throw OpalException("FieldSolver::initAmrSolver_m()",
780 "HYPRE solver not yet implemented.");
781 } else if (fsType_m == FieldSolverType::HPGMG) {
782 throw OpalException("FieldSolver::initAmrSolver_m()",
783 "HPGMG solver not yet implemented.");
784 } else if (fsType_m == FieldSolverType::AMRMG) {
785#ifdef HAVE_AMR_MG_SOLVER
786 if ( dynamic_cast<AmrBoxLib*>( itsAmrObject_mp.get() ) == 0 )
787 throw OpalException("FieldSolver::initAmrSolver_m()",
788 "FMultiGrid solver requires AMReX.");
789
790 std::string bcz = Attributes::getString(itsAttr[deprecated::BCFFTT]);
791 if (bcz.empty()) {
792 bcz = Attributes::getString(itsAttr[BCFFTZ]);
793 }
794 solver_m = new AmrMultiGrid(static_cast<AmrBoxLib*>(itsAmrObject_mp.get()),
796 Attributes::getString(itsAttr[AMR_MG_PREC]),
797 Attributes::getBool(itsAttr[AMR_MG_REBALANCE]),
798 Attributes::getString(itsAttr[AMR_MG_REUSE]),
801 bcz,
802 Attributes::getString(itsAttr[AMR_MG_SMOOTHER]),
803 Attributes::getReal(itsAttr[AMR_MG_NSWEEPS]),
804 Attributes::getString(itsAttr[AMR_MG_INTERP]),
805 Attributes::getString(itsAttr[AMR_MG_NORM]));
806
807 dynamic_cast<AmrMultiGrid*>(solver_m)->setVerbose(
808 Attributes::getBool(itsAttr[AMR_MG_VERBOSE]));
809
810 dynamic_cast<AmrMultiGrid*>(solver_m)->setTolerance(
811 Attributes::getReal(itsAttr[AMR_MG_TOL]));
812#else
813 throw OpalException("FieldSolver::initAmrSolver_m()",
814 "Multigrid solver not enabled! "
815 "Please build OPAL with -DENABLE_AMR_MG_SOLVER=1");
816#endif
817 } else
818 throw OpalException("FieldSolver::initAmrSolver_m()",
819 "Unknown solver " + getType() + ".");
820}
821#endif
UniformCartesian< 3, double > Mesh_t
Definition PBunchDefs.h:22
ParticleSpatialLayout< double, 3, Mesh_t > Layout_t
Definition PBunchDefs.h:24
CenteredFieldLayout< 3, Mesh_t, Center_t > FieldLayout_t
Definition PBunchDefs.h:28
@ SIZE
Definition IndexMap.cpp:174
e_dim_tag
Definition FieldLayout.h:55
@ PARALLEL
Definition FieldLayout.h:55
@ SERIAL
Definition FieldLayout.h:55
Inform & endl(Inform &inf)
Definition Inform.cpp:42
Inform & level3(Inform &inf)
Definition Inform.cpp:47
#define INFOMSG(msg)
Definition IpplInfo.h:348
const std::string name
Representation objects and parsers for attribute expressions.
Attribute makeBool(const std::string &name, const std::string &help)
Make logical attribute.
double getReal(const Attribute &attr)
Return real value.
void setRealArray(Attribute &attr, const std::vector< double > &value)
Set array value.
Attribute makeUpperCaseString(const std::string &name, const std::string &help)
Make uppercase string attribute.
Attribute makePredefinedString(const std::string &name, const std::string &help, const std::initializer_list< std::string > &predefinedStrings)
Make predefined string attribute.
Attribute makeReal(const std::string &name, const std::string &help)
Make real attribute.
bool getBool(const Attribute &attr)
Return logical value.
void setReal(Attribute &attr, double val)
Set real value.
Attribute makeRealArray(const std::string &name, const std::string &help)
Create real array attribute.
std::vector< double > getRealArray(const Attribute &attr)
Get array value.
std::string getString(const Attribute &attr)
Get string value.
bool amr
Enable AMR if true.
Definition Options.cpp:99
bool isAllDigits(const std::string &str)
Definition Util.cpp:218
Definition(int size, const char *name, const char *help)
Constructor for exemplars.
void registerOwnership(const AttributeHandler::OwnerType &itsClass) const
Definition Object.cpp:191
const std::string & getOpalName() const
Return object name.
Definition Object.cpp:310
std::vector< Attribute > itsAttr
The object attributes.
Definition Object.h:216
Object * find(const std::string &name)
Find entry.
Definition OpalData.cpp:571
static OpalData * getInstance()
Definition OpalData.cpp:196
BoundaryGeometry * getGlobalGeometry()
Definition OpalData.cpp:461
static std::unique_ptr< AmrBoxLib > create(const AmrInfo &info, AmrPartBunch *bunch_p)
Definition AmrBoxLib.cpp:74
int grid[3]
Number of grid points in x-, y- and z-direction.
Definition AmrObject.h:54
static std::string getTaggingString(int number)
int maxgrid[3]
Maximum grid size in x-, y- and z-direction.
Definition AmrObject.h:55
static BoundaryGeometry * find(const std::string &name)
virtual void update()
Update the field solver data.
std::unique_ptr< Layout_t > PL_m
The particle layout.
FieldSolverType fsType_m
void setMT(double)
Store emittance for mode 3.
void initCartesianFields()
double getMY() const
Return meshsize.
std::string getTagging_m() const
Mesh_t * mesh_m
The cartesian mesh.
std::string fsName_m
PartBunchBase< double, 3 > * itsBunch_m
all the particles are here ...
std::string getType()
PoissonSolver * solver_m
the actual solver, should be a base object
Inform & printInfo(Inform &os) const
bool isAmrSolverType() const
static FieldSolver * find(const std::string &name)
Find named FieldSolver.
FieldSolver()
Exemplar constructor.
void setFieldSolverType()
virtual void execute()
Execute (init) the field solver data.
double getMX() const
Return meshsize.
FieldLayout_t * FL_m
The field layout f.
bool hasPeriodicZ()
void initAmrObject_m()
std::unique_ptr< AmrObject > itsAmrObject_mp
void setMY(double)
Store emittance for mode 2.
double getMT() const
Return meshsize.
virtual ~FieldSolver()
bool hasValidSolver()
void setMX(double)
Store emittance for mode 1.
virtual FieldSolver * clone(const std::string &name)
Make clone.
void initSolver(PartBunchBase< double, 3 > *b)
void initAmrSolver_m()
The base class for all OPAL exceptions.
Definition Index.h:237
static int getNodes()
Definition IpplInfo.cpp:670