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