OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
MLPoissonSolver.cpp
Go to the documentation of this file.
1 #include "MLPoissonSolver.h"
2 
4 
5 #include "Amr/AmrDefs.h"
6 
7 #include <AMReX_MultiFabUtil.H>
8 #include <AMReX_MLPoisson.H>
9 #include <AMReX_MLMG.H>
10 
12  : AmrPoissonSolver<AmrBoxLib>(itsAmrObject_p),
13  reltol_m(1.0e-10),
14  abstol_m(0.0)
15 { }
16 
20  unsigned short baseLevel,
21  unsigned short finestLevel,
22  bool prevAsGuess)
23 {
24  for (int i = 0; i <= finestLevel; ++i) {
25  phi[i]->setVal(0.0, 1);
26  for (int j = 0; j < AMREX_SPACEDIM; ++j) {
27  efield[i][j]->setVal(0.0, 1);
28  }
29  }
30 
31  this->mlmg_m(rho, phi, efield, baseLevel, finestLevel);
32 }
33 
34 
35 double MLPoissonSolver::getXRangeMin(unsigned short level) {
36  return itsAmrObject_mp->Geom(level).ProbLo(0);
37 }
38 
39 
40 double MLPoissonSolver::getXRangeMax(unsigned short level) {
41  return itsAmrObject_mp->Geom(level).ProbHi(0);
42 }
43 
44 
45 double MLPoissonSolver::getYRangeMin(unsigned short level) {
46  return itsAmrObject_mp->Geom(level).ProbLo(1);
47 }
48 
49 
50 double MLPoissonSolver::getYRangeMax(unsigned short level) {
51  return itsAmrObject_mp->Geom(level).ProbHi(1);
52 }
53 
54 
55 double MLPoissonSolver::getZRangeMin(unsigned short level) {
56  return itsAmrObject_mp->Geom(level).ProbLo(2);
57 }
58 
59 
60 double MLPoissonSolver::getZRangeMax(unsigned short level) {
61  return itsAmrObject_mp->Geom(level).ProbHi(2);
62 }
63 
64 
66  os << "* ************* A M R e X P o i s s o n S o l v e r ************************************ " << endl
67  << "* relative tolerance " << reltol_m << '\n'
68  << "* absolute tolerance " << abstol_m << '\n' << endl
69  << "* ******************************************************************** " << endl;
70  return os;
71 }
72 
73 
77  int baseLevel,
78  int finestLevel)
79 {
80  /* According to
81  * amrex/Tutorials/LinearSolvers/ABecLaplacian_C/MyTest.H
82  */
83  amrex::LPInfo info;
84  info.setAgglomeration(true);
85  info.setConsolidation(true);
86  info.setMaxCoarseningLevel(10);
87 
88  const GeomContainer_t& geom = itsAmrObject_mp->Geom();
89 
90  int nlevels = finestLevel - baseLevel + 1;
91  GeomContainer_t geom_p(nlevels);
92  AmrGridContainer_t ba(nlevels);
93  AmrProcMapContainer_t dm(nlevels);
94  AmrFieldContainer_pt phi_p(nlevels);
95  const_AmrFieldContainer_pt rho_p(nlevels);
96 
97  for (int ilev = 0; ilev < nlevels; ++ilev) {
98  geom_p[ilev] = geom[ilev + baseLevel];
99  ba[ilev] = rho[ilev + baseLevel]->boxArray();
100  dm[ilev] = rho[ilev + baseLevel]->DistributionMap();
101  rho_p[ilev] = rho[ilev + baseLevel].get();
102  phi_p[ilev] = phi[ilev + baseLevel].get();
103  }
104 
105  amrex::MLPoisson mlpoisson(geom_p, ba, dm, info);
106 
107  mlpoisson.setMaxOrder(3);
108 
109  mlpoisson.setDomainBC({AMREX_D_DECL(amrex::LinOpBCType::Dirichlet,
110  amrex::LinOpBCType::Dirichlet,
111  amrex::LinOpBCType::Dirichlet)},
112  {AMREX_D_DECL(amrex::LinOpBCType::Dirichlet,
113  amrex::LinOpBCType::Dirichlet,
114  amrex::LinOpBCType::Dirichlet)});
115 
116  for (int ilev = 0; ilev < nlevels; ++ilev) {
117  mlpoisson.setLevelBC(ilev, phi[ilev].get());
118  }
119 
120  amrex::MLMG mlmg(mlpoisson);
121  mlmg.setMaxIter(200);
122  mlmg.setMaxFmgIter(0);
123  mlmg.setVerbose(0);
124  mlmg.setCGVerbose(0);
125 
126 
127  mlmg.solve(phi_p,
128  rho_p,
129  reltol_m, abstol_m);
130 
131  amrex::Vector<std::array<AmrField_t*, AMREX_SPACEDIM> > grad(nlevels);
132 
133  for (int i = 0; i < nlevels; ++i) {
134  for (int j = 0; j < AMREX_SPACEDIM; ++j) {
135  grad[i][j] = efield[i][j].get();
136  }
137  }
138 
139  mlmg.getGradSolution(grad);
140 
141  for (int i = 0; i < nlevels; ++i) {
142  for (int j = 0; j < AMREX_SPACEDIM; ++j) {
143  efield[i][j]->mult(-1.0, 0, 1, 1);
144  }
145  }
146 }
double getZRangeMax(unsigned short level=0)
amrex::Vector< const AmrBoxLib::AmrField_t * > const_AmrFieldContainer_pt
constexpr double e
The value of .
Definition: Physics.h:40
AmrBoxLib::AmrGeomContainer_t GeomContainer_t
bool info
Info flag.
Definition: Options.cpp:8
double getXRangeMax(unsigned short level=0)
double getYRangeMin(unsigned short level=0)
AmrBoxLib::AmrScalarFieldContainer_t AmrScalarFieldContainer_t
MLPoissonSolver(AmrBoxLib *itsAmrObject_p)
double getXRangeMin(unsigned short level=0)
amrex::Vector< AmrBoxLib::AmrField_t * > AmrFieldContainer_pt
double reltol_m
Relative tolearance for solver.
Inform & print(Inform &os) const
AmrBoxLib::AmrGridContainer_t AmrGridContainer_t
AmrBoxLib::AmrProcMapContainer_t AmrProcMapContainer_t
AmrBoxLib::AmrVectorFieldContainer_t AmrVectorFieldContainer_t
double getZRangeMin(unsigned short level=0)
double getYRangeMax(unsigned short level=0)
void solve(AmrScalarFieldContainer_t &rho, AmrScalarFieldContainer_t &phi, AmrVectorFieldContainer_t &efield, unsigned short baseLevel, unsigned short finestLevel, bool prevAsGuess=true)
double abstol_m
Absolute tolerance for solver.
Definition: Inform.h:41
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
void mlmg_m(AmrScalarFieldContainer_t &rho, AmrScalarFieldContainer_t &phi, AmrVectorFieldContainer_t &efield, int baseLevel, int finestLevel)