OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
AmrMultiGrid.h
Go to the documentation of this file.
1 #ifndef AMR_MULTI_GRID_H
2 #define AMR_MULTI_GRID_H
3 
4 #include <vector>
5 #include <memory>
6 
7 #include <AMReX.H>
8 #include <AMReX_MultiFab.H>
9 #include <AMReX_Vector.H>
10 
11 #include "AmrMultiGridCore.h"
12 
14 #include "Amr/AmrBoxLib.h"
15 
16 #include "AmrMultiGridLevel.h"
17 
18 #include <fstream>
19 
20 #define AMR_MG_TIMER true
21 #define AMR_MG_WRITE false
22 
23 class AmrMultiGrid : public AmrPoissonSolver< AmrBoxLib > {
24 
25 public:
34 
36 
46 
47  typedef BottomSolver<
48  Teuchos::RCP<matrix_t>,
49  Teuchos::RCP<mv_t>,
52 
56 
58 
61 
62  typedef amrex::BoxArray boxarray_t;
63  typedef amrex::Box box_t;
64  typedef amrex::BaseFab<int> basefab_t;
65  typedef amrex::FArrayBox farraybox_t;
66 
68 
70 
72  enum Interpolater {
73  TRILINEAR = 0,
76  };
77 
79  enum BaseSolver {
80  // all Belos
84  CG,
89  // all Amesos2
90 #ifdef HAVE_AMESOS2_KLU2
91  , KLU2
92 #endif
93 #ifdef HAVE_AMESOS2_SUPERLU
94  , SUPERLU
95 #endif
96 #ifdef HAVE_AMESOS2_UMFPACK
97  , UMFPACK
98 #endif
99 #ifdef HAVE_AMESOS2_PARDISO_MKL
100  , PARDISO_MKL
101 #endif
102 #ifdef HAVE_AMESOS2_MUMPS
103  , MUMPS
104 #endif
105 #ifdef HAVE_AMESOS2_LAPACK
106  , LAPACK
107 #endif
108  // all MueLu
109  , SA
110  // add others ...
111  };
112 
114  enum Boundary {
118  };
119 
121  enum Norm {
122  L1,
123  L2,
125  };
126 
127 public:
128 
143  AmrMultiGrid(AmrBoxLib* itsAmrObject_p,
144  const std::string& bsolver,
145  const std::string& prec,
146  const bool& rebalance,
147  const std::string& reuse,
148  const std::string& bcx,
149  const std::string& bcy,
150  const std::string& bcz,
151  const std::string& smoother,
152  const std::size_t& nSweeps,
153  const std::string& interp,
154  const std::string& norm);
155 
169  unsigned short baseLevel,
170  unsigned short finestLevel,
171  bool prevAsGuess = true);
172 
177  void setNumberOfSweeps(const std::size_t& nSweeps);
178 
183  void setMaxNumberOfIterations(const std::size_t& maxiter);
184 
189  std::size_t getNumIters();
190 
197 
201  void setVerbose(bool verbose);
202 
203  double getXRangeMin(unsigned short level = 0);
204  double getXRangeMax(unsigned short level = 0);
205  double getYRangeMin(unsigned short level = 0);
206  double getYRangeMax(unsigned short level = 0);
207  double getZRangeMin(unsigned short level = 0);
208  double getZRangeMax(unsigned short level = 0);
209 
214  Inform &print(Inform &os) const;
215 
216 private:
217 
223  void initPhysicalBoundary_m(const Boundary* bc);
224 
231  void initLevels_m(const amrex::Vector<AmrField_u>& rho,
232  const amrex::Vector<AmrGeometry_t>& geom,
233  bool regrid);
234 
238  void clearMasks_m();
239 
244  void initGuess_m(bool reset);
245 
251 
256  bool isConverged_m(std::vector<scalar_t>& rhsNorms,
257  std::vector<scalar_t>& resNorms);
258 
266  void residual_m(const lo_t& level,
267  Teuchos::RCP<vector_t>& r,
268  const Teuchos::RCP<vector_t>& b,
269  const Teuchos::RCP<vector_t>& x);
270 
275  void relax_m(const lo_t& level);
276 
285  void residual_no_fine_m(const lo_t& level,
286  Teuchos::RCP<vector_t>& result,
287  const Teuchos::RCP<vector_t>& rhs,
288  const Teuchos::RCP<vector_t>& crs_rhs,
289  const Teuchos::RCP<vector_t>& b);
290 
291 #if AMR_MG_WRITE
292 
295  void writeResidualNorm_m();
296 #endif
297 
303  scalar_t evalNorm_m(const Teuchos::RCP<const vector_t>& x);
304 
310  void initResidual_m(std::vector<scalar_t>& rhsNorms,
311  std::vector<scalar_t>& resNorms);
312 
317 
322  void setup_m(const amrex::Vector<AmrField_u>& rho,
323  const amrex::Vector<AmrField_u>& phi,
324  const bool& matrices = true);
325 
329  void buildSingleLevel_m(const amrex::Vector<AmrField_u>& rho,
330  const amrex::Vector<AmrField_u>& phi,
331  const bool& matrices = true);
332 
336  void buildMultiLevel_m(const amrex::Vector<AmrField_u>& rho,
337  const amrex::Vector<AmrField_u>& phi,
338  const bool& matrices = true);
339 
345  void open_m(const lo_t& level, const bool& matrices);
346 
352  void close_m(const lo_t& level, const bool& matrices);
353 
364  void buildNoFinePoissonMatrix_m(const lo_t& level,
365  const go_t& gidx,
366  const AmrIntVect_t& iv,
367  const basefab_t& mfab,
368  const scalar_t* invdx2);
369 
381  void buildCompositePoissonMatrix_m(const lo_t& level,
382  const go_t& gidx,
383  const AmrIntVect_t& iv,
384  const basefab_t& mfab,
385  const basefab_t& rfab,
386  const basefab_t& cfab,
387  const scalar_t* invdx2);
388 
399  void buildRestrictionMatrix_m(const lo_t& level,
400  const go_t& gidx,
401  const AmrIntVect_t& iv,
402  D_DECL(const go_t& ii,
403  const go_t& jj,
404  const go_t& kk),
405  const basefab_t& rfab);
406 
418  void buildInterpolationMatrix_m(const lo_t& level,
419  const go_t& gidx,
420  const AmrIntVect_t& iv,
421  const basefab_t& cfab);
422 
435  void buildCrseBoundaryMatrix_m(const lo_t& level,
436  const go_t& gidx,
437  const AmrIntVect_t& iv,
438  const basefab_t& mfab,
439  const basefab_t& cfab,
440  const scalar_t* invdx2);
441 
455  void buildFineBoundaryMatrix_m(const lo_t& level,
456  const go_t& gidx,
457  const AmrIntVect_t& iv,
458  const basefab_t& mfab,
459  const basefab_t& rfab,
460  const basefab_t& cfab);
461 
467  void buildDensityVector_m(const lo_t& level,
468  const AmrField_t& rho);
469 
475  void buildPotentialVector_m(const lo_t& level,
476  const AmrField_t& phi);
477 
484  void buildGradientMatrix_m(const lo_t& level,
485  const go_t& gidx,
486  const AmrIntVect_t& iv,
487  const basefab_t& mfab,
488  const scalar_t* invdx);
489 
497  void amrex2trilinos_m(const lo_t& level,
498  const lo_t& comp,
499  const AmrField_t& mf,
500  Teuchos::RCP<vector_t>& mv);
501 
509  void trilinos2amrex_m(const lo_t& level,
510  const lo_t& comp,
511  AmrField_t& mf,
512  const Teuchos::RCP<vector_t>& mv);
513 
521  void map2vector_m(umap_t& map, indices_t& indices,
522  coefficients_t& values);
523 
530  void smooth_m(const lo_t& level,
531  Teuchos::RCP<vector_t>& e,
532  Teuchos::RCP<vector_t>& r);
533 
538  void restrict_m(const lo_t& level);
539 
544  void prolongate_m(const lo_t& level);
545 
550  void averageDown_m(const lo_t& level);
551 
556  void initInterpolater_m(const Interpolater& interp);
557 
562  void initCrseFineInterp_m(const Interpolater& interface);
563 
570  void initBaseSolver_m(const BaseSolver& solver,
571  const bool& rebalance,
572  const std::string& reuse);
573 
580  void initPrec_m(const Preconditioner& prec,
581  const bool& rebalance,
582  const std::string& reuse);
583 
588  Boundary convertToEnumBoundary_m(const std::string& bc);
589 
594  Interpolater convertToEnumInterpolater_m(const std::string& interp);
595 
600  BaseSolver convertToEnumBaseSolver_m(const std::string& bsolver);
601 
606  Preconditioner convertToEnumPreconditioner_m(const std::string& prec);
607 
612  Smoother convertToEnumSmoother_m(const std::string& smoother);
613 
618  Norm convertToEnumNorm_m(const std::string& norm);
619 
624  void writeSDDSHeader_m(std::ofstream& outfile);
625 
630  void writeSDDSData_m(const scalar_t& error);
631 
632 #if AMR_MG_TIMER
633 
636  void initTimer_m();
637 #endif
638 
639 private:
640  Teuchos::RCP<comm_t> comm_mp;
641  Teuchos::RCP<amr::node_t> node_mp;
642 
644  std::unique_ptr<AmrInterpolater<AmrMultiGridLevel_t> > interp_mp;
645 
647  std::unique_ptr<AmrInterpolater<AmrMultiGridLevel_t> > interface_mp;
648 
649  std::size_t nIter_m;
650  std::size_t bIter_m;
651  std::size_t maxiter_m;
652  std::size_t nSweeps_m;
654 
656  std::vector<std::unique_ptr<AmrMultiGridLevel_t > > mglevel_m;
657 
659  std::shared_ptr<bsolver_t> solver_mp;
660 
662  std::vector<std::shared_ptr<AmrSmoother> > smoother_m;
663 
665  std::shared_ptr<preconditioner_t> prec_mp;
666 
667  int lbase_m;
668  int lfine_m;
669  int nlevel_m;
670 
671  boundary_t bc_m[AMREX_SPACEDIM];
673 
675  std::string snorm_m;
676 
677  const scalar_t eps_m;
678 
679  bool verbose_m;
680  std::string fname_m;
681  std::ios_base::openmode flag_m;
682 
683 #if AMR_MG_TIMER
684  IpplTimings::TimerRef buildTimer_m;
685  IpplTimings::TimerRef restrictTimer_m;
686  IpplTimings::TimerRef smoothTimer_m;
687  IpplTimings::TimerRef interpTimer_m;
688  IpplTimings::TimerRef efieldTimer_m;
689  IpplTimings::TimerRef averageTimer_m;
690  IpplTimings::TimerRef bottomTimer_m;
691  IpplTimings::TimerRef dumpTimer_m;
692 #endif
693 };
694 
695 
696 inline Inform &operator<<(Inform &os, const AmrMultiGrid &fs) {
697  return fs.print(os);
698 }
699 
700 #endif
std::ostream & operator<<(std::ostream &os, const Attribute &attr)
Definition: Attribute.cpp:167
std::vector< scalar_t > coefficients_t
Type for matrix entries.
std::vector< std::unique_ptr< AmrMultiGridLevel_t > > mglevel_m
container for levels
Definition: AmrMultiGrid.h:656
Interface to Amesos2 solvers of the Trilinos package.
std::unique_ptr< AmrInterpolater< AmrMultiGridLevel_t > > interface_mp
interpolater for coarse-fine interface
Definition: AmrMultiGrid.h:647
AmrMultiGridLevel_t::boundary_t boundary_t
Definition: AmrMultiGrid.h:45
Bottom solver preconditioners.
constexpr double e
The value of .
Definition: Physics.h:40
void close_m(const lo_t &level, const bool &matrices)
std::unordered_map< go_t, scalar_t > umap_t
void buildNoFinePoissonMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const scalar_t *invdx2)
MueLuPreconditioner< AmrMultiGridLevel_t > MueLuPreconditioner_t
Definition: AmrMultiGrid.h:60
Norm convertToEnumNorm_m(const std::string &norm)
AmrPreconditioner< matrix_t, AmrMultiGridLevel_t > preconditioner_t
Definition: AmrMultiGrid.h:57
amr::AmrField_t AmrField_t
AmrMultiGrid(AmrBoxLib *itsAmrObject_p, const std::string &bsolver, const std::string &prec, const bool &rebalance, const std::string &reuse, const std::string &bcx, const std::string &bcy, const std::string &bcz, const std::string &smoother, const std::size_t &nSweeps, const std::string &interp, const std::string &norm)
void residual_m(const lo_t &level, Teuchos::RCP< vector_t > &r, const Teuchos::RCP< vector_t > &b, const Teuchos::RCP< vector_t > &x)
int lfine_m
fineste level
Definition: AmrMultiGrid.h:668
void buildDensityVector_m(const lo_t &level, const AmrField_t &rho)
amrex::FArrayBox farraybox_t
Definition: AmrMultiGrid.h:65
Teuchos::RCP< comm_t > comm_mp
communicator
Definition: AmrMultiGrid.h:640
void smooth_m(const lo_t &level, Teuchos::RCP< vector_t > &e, Teuchos::RCP< vector_t > &r)
double getXRangeMin(unsigned short level=0)
void buildRestrictionMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, D_DECL(const go_t &ii, const go_t &jj, const go_t &kk), const basefab_t &rfab)
int nBcPoints_m
maximum number of stencils points for BC
Definition: AmrMultiGrid.h:672
scalar_t getLevelResidualNorm(lo_t level)
int local_ordinal_t
bool verbose_m
If true, a SDDS file is written.
Definition: AmrMultiGrid.h:679
Smoother
All supported Ifpack2 smoothers.
Definition: AmrSmoother.h:27
void buildPotentialVector_m(const lo_t &level, const AmrField_t &phi)
void initCrseFineInterp_m(const Interpolater &interface)
void setVerbose(bool verbose)
FRONT * fs
Definition: hypervolume.cpp:59
void initGuess_m(bool reset)
Amesos2BottomSolver< AmrMultiGridLevel_t > Amesos2Solver_t
Definition: AmrMultiGrid.h:54
std::ios_base::openmode flag_m
std::ios::out or std::ios::app
Definition: AmrMultiGrid.h:681
void open_m(const lo_t &level, const bool &matrices)
MueLuBottomSolver< AmrMultiGridLevel_t > MueLuSolver_t
Definition: AmrMultiGrid.h:55
BaseSolver convertToEnumBaseSolver_m(const std::string &bsolver)
void initInterpolater_m(const Interpolater &interp)
std::size_t nIter_m
number of iterations till convergence
Definition: AmrMultiGrid.h:649
amrex::Box box_t
Definition: AmrMultiGrid.h:63
int nlevel_m
number of levelss
Definition: AmrMultiGrid.h:669
std::size_t nSweeps_m
number of smoothing iterations
Definition: AmrMultiGrid.h:652
AmrMultiGridLevel_t::AmrIntVect_t AmrIntVect_t
Definition: AmrMultiGrid.h:41
void solve(AmrScalarFieldContainer_t &rho, AmrScalarFieldContainer_t &phi, AmrVectorFieldContainer_t &efield, unsigned short baseLevel, unsigned short finestLevel, bool prevAsGuess=true)
amr::scalar_t scalar_t
Definition: AmrMultiGrid.h:33
double getYRangeMin(unsigned short level=0)
void setNumberOfSweeps(const std::size_t &nSweeps)
Abstract base class for all base level solvers.
Definition: BottomSolver.h:8
BottomSolver< Teuchos::RCP< matrix_t >, Teuchos::RCP< mv_t >, AmrMultiGridLevel_t > bsolver_t
Definition: AmrMultiGrid.h:51
int lbase_m
base level (currently only 0 supported)
Definition: AmrMultiGrid.h:667
void writeSDDSData_m(const scalar_t &error)
double getXRangeMax(unsigned short level=0)
Interpolater convertToEnumInterpolater_m(const std::string &interp)
AmrMultiGridLevel_t::AmrGeometry_t AmrGeometry_t
Definition: AmrMultiGrid.h:38
std::string snorm_m
norm for convergence criteria
Definition: AmrMultiGrid.h:675
void buildSingleLevel_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true)
AmrMultiGridLevel< matrix_t, vector_t > AmrMultiGridLevel_t
Definition: AmrMultiGrid.h:35
void buildCompositePoissonMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const basefab_t &rfab, const basefab_t &cfab, const scalar_t *invdx2)
double getZRangeMax(unsigned short level=0)
std::shared_ptr< bsolver_t > solver_mp
bottom solver
Definition: AmrMultiGrid.h:659
amr::dmap_t dmap_t
Definition: AmrMultiGrid.h:29
amr::AmrIntVect_t AmrIntVect_t
void buildCrseBoundaryMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const basefab_t &cfab, const scalar_t *invdx2)
AmrMultiGridLevel_t::coefficients_t coefficients_t
Definition: AmrMultiGrid.h:43
void initPrec_m(const Preconditioner &prec, const bool &rebalance, const std::string &reuse)
Teuchos::MpiComm< int > comm_t
double getYRangeMax(unsigned short level=0)
std::shared_ptr< AmrBoundary< AmrMultiGridLevel< MatrixType, VectorType > > > boundary_t
std::size_t getNumIters()
void map2vector_m(umap_t &map, indices_t &indices, coefficients_t &values)
amr::local_ordinal_t lo_t
Definition: AmrMultiGrid.h:31
BaseSolver
Supported bottom solvers.
Definition: AmrMultiGrid.h:79
std::vector< std::shared_ptr< AmrSmoother > > smoother_m
error smoother
Definition: AmrMultiGrid.h:662
amr::vector_t vector_t
Definition: AmrMultiGrid.h:27
AmrMultiGridLevel_t::AmrField_t AmrField_t
Definition: AmrMultiGrid.h:37
Norm norm_m
norm for convergence criteria (l1, l2, linf)
Definition: AmrMultiGrid.h:674
boundary_t bc_m[AMREX_SPACEDIM]
boundary conditions
Definition: AmrMultiGrid.h:671
Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > dmap_t
std::unique_ptr< AmrField_t > AmrField_u
void restrict_m(const lo_t &level)
double getZRangeMin(unsigned short level=0)
void averageDown_m(const lo_t &level)
void initResidual_m(std::vector< scalar_t > &rhsNorms, std::vector< scalar_t > &resNorms)
void prolongate_m(const lo_t &level)
Smoother convertToEnumSmoother_m(const std::string &smoother)
std::size_t maxiter_m
maximum number of iterations allowed
Definition: AmrMultiGrid.h:651
void clearMasks_m()
void buildInterpolationMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &cfab)
scalar_t iterate_m()
AmrMultiGridLevel_t::umap_t umap_t
Definition: AmrMultiGrid.h:44
void initPhysicalBoundary_m(const Boundary *bc)
Inform & print(Inform &os) const
Interface to Belos solvers of the Trilinos package.
void initBaseSolver_m(const BaseSolver &solver, const bool &rebalance, const std::string &reuse)
std::string fname_m
SDDS filename.
Definition: AmrMultiGrid.h:680
Norm
Supported convergence criteria.
Definition: AmrMultiGrid.h:121
Ifpack2Preconditioner< AmrMultiGridLevel_t > Ifpack2Preconditioner_t
Definition: AmrMultiGrid.h:59
Boundary convertToEnumBoundary_m(const std::string &bc)
Tpetra::CrsMatrix< scalar_t, local_ordinal_t, global_ordinal_t, node_t > matrix_t
Preconditioner convertToEnumPreconditioner_m(const std::string &prec)
AmrMultiGridLevel_t::AmrField_u AmrField_u
Definition: AmrMultiGrid.h:39
amr::comm_t comm_t
Definition: AmrMultiGrid.h:30
Tpetra::MultiVector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > multivector_t
void buildFineBoundaryMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const basefab_t &rfab, const basefab_t &cfab)
void setup_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true)
long global_ordinal_t
Smoother smootherType_m
type of smoother
Definition: AmrMultiGrid.h:653
amr::matrix_t matrix_t
Definition: AmrMultiGrid.h:26
void residual_no_fine_m(const lo_t &level, Teuchos::RCP< vector_t > &result, const Teuchos::RCP< vector_t > &rhs, const Teuchos::RCP< vector_t > &crs_rhs, const Teuchos::RCP< vector_t > &b)
amrex::BoxArray boxarray_t
Definition: AmrMultiGrid.h:62
void writeSDDSHeader_m(std::ofstream &outfile)
std::shared_ptr< AmrField_t > AmrField_s
Teuchos::RCP< amr::node_t > node_mp
kokkos node
Definition: AmrMultiGrid.h:641
Boundary
Supported physical boundaries.
Definition: AmrMultiGrid.h:114
void buildGradientMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const scalar_t *invdx)
scalar_t evalNorm_m(const Teuchos::RCP< const vector_t > &x)
void buildMultiLevel_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true)
amr::global_ordinal_t go_t
Definition: AmrMultiGrid.h:32
void trilinos2amrex_m(const lo_t &level, const lo_t &comp, AmrField_t &mf, const Teuchos::RCP< vector_t > &mv)
std::shared_ptr< preconditioner_t > prec_mp
preconditioner for bottom solver
Definition: AmrMultiGrid.h:665
BelosBottomSolver< AmrMultiGridLevel_t > BelosSolver_t
Definition: AmrMultiGrid.h:53
Timing::TimerRef TimerRef
Definition: IpplTimings.h:176
AmrMultiGridLevel_t::indices_t indices_t
Definition: AmrMultiGrid.h:42
void initLevels_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrGeometry_t > &geom, bool regrid)
amr::AmrScalarFieldContainer_t AmrScalarFieldContainer_t
Definition: PBunchDefs.h:60
void setMaxNumberOfIterations(const std::size_t &maxiter)
AmrMultiGridLevel_t::AmrField_s AmrField_s
Definition: AmrMultiGrid.h:40
Tpetra::Vector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > vector_t
std::size_t bIter_m
number of iterations of bottom solver
Definition: AmrMultiGrid.h:650
const scalar_t eps_m
rhs scale for convergence
Definition: AmrMultiGrid.h:677
Definition: Inform.h:41
void relax_m(const lo_t &level)
amr::Preconditioner Preconditioner
Definition: AmrMultiGrid.h:69
Interpolater
Supported interpolaters for prolongation operation.
Definition: AmrMultiGrid.h:72
void computeEfield_m(AmrVectorFieldContainer_t &efield)
bool isConverged_m(std::vector< scalar_t > &rhsNorms, std::vector< scalar_t > &resNorms)
std::vector< go_t > indices_t
Type for matrix indices.
amr::multivector_t mv_t
Definition: AmrMultiGrid.h:28
amrex::BaseFab< int > basefab_t
Definition: AmrMultiGrid.h:64
amr::AmrGeometry_t AmrGeometry_t
amr::AmrVectorFieldContainer_t AmrVectorFieldContainer_t
Definition: PBunchDefs.h:61
double scalar_t
AmrSmoother::Smoother Smoother
Definition: AmrMultiGrid.h:67
amr::AmrField_t AmrField_t
Definition: PBunchDefs.h:59
std::unique_ptr< AmrInterpolater< AmrMultiGridLevel_t > > interp_mp
interpolater without coarse-fine interface
Definition: AmrMultiGrid.h:644
void amrex2trilinos_m(const lo_t &level, const lo_t &comp, const AmrField_t &mf, Teuchos::RCP< vector_t > &mv)