OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
OPAL
AmrMultiGrid.h
Go to the documentation of this file.
1 //
2 // Class AmrMultiGrid
3 // Main class of the AMR Poisson multigrid solver.
4 // It implements the multigrid solver described in https://doi.org/10.1016/j.cpc.2019.106912
5 //
6 // Copyright (c) 2017 - 2020, Matthias Frey, Paul Scherrer Institut, Villigen PSI, Switzerland
7 // All rights reserved
8 //
9 // Implemented as part of the PhD thesis
10 // "Precise Simulations of Multibunches in High Intensity Cyclotrons"
11 //
12 // This file is part of OPAL.
13 //
14 // OPAL is free software: you can redistribute it and/or modify
15 // it under the terms of the GNU General Public License as published by
16 // the Free Software Foundation, either version 3 of the License, or
17 // (at your option) any later version.
18 //
19 // You should have received a copy of the GNU General Public License
20 // along with OPAL. If not, see <https://www.gnu.org/licenses/>.
21 //
22 #ifndef AMR_MULTI_GRID_H
23 #define AMR_MULTI_GRID_H
24 
25 #include <vector>
26 #include <memory>
27 
28 #include <AMReX.H>
29 #include <AMReX_MultiFab.H>
30 #include <AMReX_Vector.H>
31 
32 #include "AmrMultiGridCore.h"
33 
35 #include "Amr/AmrBoxLib.h"
36 
37 #include "AmrMultiGridLevel.h"
38 
39 #include <fstream>
40 
41 #define AMR_MG_TIMER true
42 #define AMR_MG_WRITE false
43 
44 class AmrMultiGrid : public AmrPoissonSolver< AmrBoxLib > {
45 
46 public:
55 
57 
67 
68  typedef BottomSolver<
69  Teuchos::RCP<matrix_t>,
70  Teuchos::RCP<mv_t>,
73 
77 
79 
82 
83  typedef amrex::BoxArray boxarray_t;
84  typedef amrex::Box box_t;
85  typedef amrex::BaseFab<int> basefab_t;
86  typedef amrex::FArrayBox farraybox_t;
87 
89 
91 
93  enum Interpolater {
94  TRILINEAR = 0,
97  };
98 
100  enum BaseSolver {
101  // all Belos
105  CG,
110  // all Amesos2
111 #ifdef HAVE_AMESOS2_KLU2
112  , KLU2
113 #endif
114 #ifdef HAVE_AMESOS2_SUPERLU
115  , SUPERLU
116 #endif
117 #ifdef HAVE_AMESOS2_UMFPACK
118  , UMFPACK
119 #endif
120 #ifdef HAVE_AMESOS2_PARDISO_MKL
121  , PARDISO_MKL
122 #endif
123 #ifdef HAVE_AMESOS2_MUMPS
124  , MUMPS
125 #endif
126 #ifdef HAVE_AMESOS2_LAPACK
127  , LAPACK
128 #endif
129  // all MueLu
130  , SA
131  // add others ...
132  };
133 
135  enum Boundary {
138  PERIODIC
139  };
140 
142  enum Norm {
143  L1,
144  L2,
145  LINF
146  };
147 
148 public:
149 
164  AmrMultiGrid(AmrBoxLib* itsAmrObject_p,
165  const std::string& bsolver,
166  const std::string& prec,
167  const bool& rebalance,
168  const std::string& reuse,
169  const std::string& bcx,
170  const std::string& bcy,
171  const std::string& bcz,
172  const std::string& smoother,
173  const std::size_t& nSweeps,
174  const std::string& interp,
175  const std::string& norm);
176 
190  unsigned short baseLevel,
191  unsigned short finestLevel,
192  bool prevAsGuess = true);
193 
198  void setNumberOfSweeps(const std::size_t& nSweeps);
199 
204  void setMaxNumberOfIterations(const std::size_t& maxiter);
205 
210  std::size_t getNumIters();
211 
218 
222  void setVerbose(bool verbose);
223 
224 
228  void setTolerance(const scalar_t& eps);
229 
230  double getXRangeMin(unsigned short level = 0);
231  double getXRangeMax(unsigned short level = 0);
232  double getYRangeMin(unsigned short level = 0);
233  double getYRangeMax(unsigned short level = 0);
234  double getZRangeMin(unsigned short level = 0);
235  double getZRangeMax(unsigned short level = 0);
236 
241  Inform &print(Inform &os) const;
242 
243 private:
244 
250  void initPhysicalBoundary_m(const Boundary* bc);
251 
258  void initLevels_m(const amrex::Vector<AmrField_u>& rho,
259  const amrex::Vector<AmrGeometry_t>& geom,
260  bool regrid);
261 
265  void clearMasks_m();
266 
271  void initGuess_m(bool reset);
272 
278 
283  bool isConverged_m(std::vector<scalar_t>& rhsNorms,
284  std::vector<scalar_t>& resNorms);
285 
293  void residual_m(const lo_t& level,
294  Teuchos::RCP<vector_t>& r,
295  const Teuchos::RCP<vector_t>& b,
296  const Teuchos::RCP<vector_t>& x);
297 
302  void relax_m(const lo_t& level);
303 
312  void residual_no_fine_m(const lo_t& level,
313  Teuchos::RCP<vector_t>& result,
314  const Teuchos::RCP<vector_t>& rhs,
315  const Teuchos::RCP<vector_t>& crs_rhs,
316  const Teuchos::RCP<vector_t>& b);
317 
318 #if AMR_MG_WRITE
322  void writeResidualNorm_m();
323 #endif
324 
330  scalar_t evalNorm_m(const Teuchos::RCP<const vector_t>& x);
331 
337  void initResidual_m(std::vector<scalar_t>& rhsNorms,
338  std::vector<scalar_t>& resNorms);
339 
344 
349  void setup_m(const amrex::Vector<AmrField_u>& rho,
350  const amrex::Vector<AmrField_u>& phi,
351  const bool& matrices = true);
352 
356  void buildSingleLevel_m(const amrex::Vector<AmrField_u>& rho,
357  const amrex::Vector<AmrField_u>& phi,
358  const bool& matrices = true);
359 
363  void buildMultiLevel_m(const amrex::Vector<AmrField_u>& rho,
364  const amrex::Vector<AmrField_u>& phi,
365  const bool& matrices = true);
366 
372  void open_m(const lo_t& level, const bool& matrices);
373 
379  void close_m(const lo_t& level, const bool& matrices);
380 
391  void buildNoFinePoissonMatrix_m(const lo_t& level,
392  const go_t& gidx,
393  const AmrIntVect_t& iv,
394  const basefab_t& mfab,
395  const scalar_t* invdx2);
396 
408  void buildCompositePoissonMatrix_m(const lo_t& level,
409  const go_t& gidx,
410  const AmrIntVect_t& iv,
411  const basefab_t& mfab,
412  const basefab_t& rfab,
413  const scalar_t* invdx2);
414 
425  void buildRestrictionMatrix_m(const lo_t& level,
426  const go_t& gidx,
427  const AmrIntVect_t& iv,
428  D_DECL(const go_t& ii,
429  const go_t& jj,
430  const go_t& kk),
431  const basefab_t& rfab);
432 
444  void buildInterpolationMatrix_m(const lo_t& level,
445  const go_t& gidx,
446  const AmrIntVect_t& iv,
447  const basefab_t& cfab);
448 
461  void buildCrseBoundaryMatrix_m(const lo_t& level,
462  const go_t& gidx,
463  const AmrIntVect_t& iv,
464  const basefab_t& mfab,
465  const basefab_t& cfab,
466  const scalar_t* invdx2);
467 
481  void buildFineBoundaryMatrix_m(const lo_t& level,
482  const go_t& gidx,
483  const AmrIntVect_t& iv,
484  const basefab_t& mfab,
485  const basefab_t& rfab);
486 
492  void buildDensityVector_m(const lo_t& level,
493  const AmrField_t& rho);
494 
500  void buildPotentialVector_m(const lo_t& level,
501  const AmrField_t& phi);
502 
509  void buildGradientMatrix_m(const lo_t& level,
510  const go_t& gidx,
511  const AmrIntVect_t& iv,
512  const basefab_t& mfab,
513  const scalar_t* invdx);
514 
522  void amrex2trilinos_m(const lo_t& level,
523  const lo_t& comp,
524  const AmrField_t& mf,
525  Teuchos::RCP<vector_t>& mv);
526 
534  void trilinos2amrex_m(const lo_t& level,
535  const lo_t& comp,
536  AmrField_t& mf,
537  const Teuchos::RCP<vector_t>& mv);
538 
546  void map2vector_m(umap_t& map, indices_t& indices,
547  coefficients_t& values);
548 
555  void smooth_m(const lo_t& level,
556  Teuchos::RCP<vector_t>& e,
557  Teuchos::RCP<vector_t>& r);
558 
563  void restrict_m(const lo_t& level);
564 
569  void prolongate_m(const lo_t& level);
570 
575  void averageDown_m(const lo_t& level);
576 
581  void initInterpolater_m(const Interpolater& interp);
582 
587  void initCrseFineInterp_m(const Interpolater& interface);
588 
595  void initBaseSolver_m(const BaseSolver& solver,
596  const bool& rebalance,
597  const std::string& reuse);
598 
604  void initPrec_m(const Preconditioner& prec,
605  const std::string& reuse);
606 
611  Boundary convertToEnumBoundary_m(const std::string& bc);
612 
617  Interpolater convertToEnumInterpolater_m(const std::string& interp);
618 
623  BaseSolver convertToEnumBaseSolver_m(const std::string& bsolver);
624 
629  Preconditioner convertToEnumPreconditioner_m(const std::string& prec);
630 
635  Smoother convertToEnumSmoother_m(const std::string& smoother);
636 
641  Norm convertToEnumNorm_m(const std::string& norm);
642 
647  void writeSDDSHeader_m(std::ofstream& outfile);
648 
653  void writeSDDSData_m(const scalar_t& error);
654 
655 #if AMR_MG_TIMER
659  void initTimer_m();
660 #endif
661 
662 private:
663  Teuchos::RCP<comm_t> comm_mp;
664 
666  std::unique_ptr<AmrInterpolater<AmrMultiGridLevel_t> > interp_mp;
667 
669  std::unique_ptr<AmrInterpolater<AmrMultiGridLevel_t> > interface_mp;
670 
671  std::size_t nIter_m;
672  std::size_t bIter_m;
673  std::size_t maxiter_m;
674  std::size_t nSweeps_m;
676 
678  std::vector<std::unique_ptr<AmrMultiGridLevel_t > > mglevel_m;
679 
681  std::shared_ptr<bsolver_t> solver_mp;
682 
684  std::vector<std::shared_ptr<AmrSmoother> > smoother_m;
685 
687  std::shared_ptr<preconditioner_t> prec_mp;
688 
689  int lbase_m;
690  int lfine_m;
691  int nlevel_m;
692 
693  boundary_t bc_m[AMREX_SPACEDIM];
695 
697  std::string snorm_m;
698 
700 
701  bool verbose_m;
702  std::string fname_m;
703  std::ios_base::openmode flag_m;
704 
705 #if AMR_MG_TIMER
706  IpplTimings::TimerRef buildTimer_m;
707  IpplTimings::TimerRef restrictTimer_m;
708  IpplTimings::TimerRef smoothTimer_m;
709  IpplTimings::TimerRef interpTimer_m;
710  IpplTimings::TimerRef efieldTimer_m;
711  IpplTimings::TimerRef averageTimer_m;
712  IpplTimings::TimerRef bottomTimer_m;
713  IpplTimings::TimerRef dumpTimer_m;
714 #endif
715 };
716 
717 
718 inline Inform &operator<<(Inform &os, const AmrMultiGrid &fs) {
719  return fs.print(os);
720 }
721 
722 #endif
amr::AmrScalarFieldContainer_t AmrScalarFieldContainer_t
Definition: PBunchDefs.h:35
amr::AmrField_t AmrField_t
Definition: PBunchDefs.h:34
amr::AmrVectorFieldContainer_t AmrVectorFieldContainer_t
Definition: PBunchDefs.h:36
Inform & operator<<(Inform &os, const AmrMultiGrid &fs)
Definition: AmrMultiGrid.h:718
long global_ordinal_t
Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > dmap_t
int local_ordinal_t
Tpetra::Vector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > vector_t
double scalar_t
Tpetra::MultiVector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > multivector_t
Tpetra::CrsMatrix< scalar_t, local_ordinal_t, global_ordinal_t, node_t > matrix_t
Teuchos::MpiComm< int > comm_t
constexpr double e
The value of.
Definition: Physics.h:39
FRONT * fs
Definition: hypervolume.cpp:59
void trilinos2amrex_m(const lo_t &level, const lo_t &comp, AmrField_t &mf, const Teuchos::RCP< vector_t > &mv)
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)
amrex::BaseFab< int > basefab_t
Definition: AmrMultiGrid.h:85
void relax_m(const lo_t &level)
void initGuess_m(bool reset)
void close_m(const lo_t &level, const bool &matrices)
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)
double getZRangeMax(unsigned short level=0)
amr::Preconditioner Preconditioner
Definition: AmrMultiGrid.h:90
void restrict_m(const lo_t &level)
Boundary
Supported physical boundaries.
Definition: AmrMultiGrid.h:135
Boundary convertToEnumBoundary_m(const std::string &bc)
std::shared_ptr< bsolver_t > solver_mp
bottom solver
Definition: AmrMultiGrid.h:681
BaseSolver convertToEnumBaseSolver_m(const std::string &bsolver)
void setup_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true)
void clearMasks_m()
AmrMultiGridLevel_t::boundary_t boundary_t
Definition: AmrMultiGrid.h:66
double getXRangeMin(unsigned short level=0)
void initPhysicalBoundary_m(const Boundary *bc)
AmrMultiGridLevel_t::coefficients_t coefficients_t
Definition: AmrMultiGrid.h:64
AmrMultiGridLevel_t::AmrGeometry_t AmrGeometry_t
Definition: AmrMultiGrid.h:59
Teuchos::RCP< comm_t > comm_mp
communicator
Definition: AmrMultiGrid.h:663
std::string fname_m
SDDS filename.
Definition: AmrMultiGrid.h:702
BaseSolver
Supported bottom solvers.
Definition: AmrMultiGrid.h:100
AmrMultiGridLevel_t::umap_t umap_t
Definition: AmrMultiGrid.h:65
amrex::FArrayBox farraybox_t
Definition: AmrMultiGrid.h:86
std::string snorm_m
norm for convergence criteria
Definition: AmrMultiGrid.h:697
std::size_t nSweeps_m
number of smoothing iterations
Definition: AmrMultiGrid.h:674
void map2vector_m(umap_t &map, indices_t &indices, coefficients_t &values)
void buildFineBoundaryMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const basefab_t &rfab)
void writeSDDSHeader_m(std::ofstream &outfile)
std::vector< std::unique_ptr< AmrMultiGridLevel_t > > mglevel_m
container for levels
Definition: AmrMultiGrid.h:678
void amrex2trilinos_m(const lo_t &level, const lo_t &comp, const AmrField_t &mf, Teuchos::RCP< vector_t > &mv)
scalar_t getLevelResidualNorm(lo_t level)
scalar_t evalNorm_m(const Teuchos::RCP< const vector_t > &x)
void buildSingleLevel_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true)
int lbase_m
base level (currently only 0 supported)
Definition: AmrMultiGrid.h:689
void initPrec_m(const Preconditioner &prec, const std::string &reuse)
void setVerbose(bool verbose)
void initInterpolater_m(const Interpolater &interp)
amr::comm_t comm_t
Definition: AmrMultiGrid.h:51
std::unique_ptr< AmrInterpolater< AmrMultiGridLevel_t > > interp_mp
interpolater without coarse-fine interface
Definition: AmrMultiGrid.h:666
amr::scalar_t scalar_t
Definition: AmrMultiGrid.h:54
scalar_t iterate_m()
void buildGradientMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const scalar_t *invdx)
AmrMultiGridLevel_t::AmrField_u AmrField_u
Definition: AmrMultiGrid.h:60
BelosBottomSolver< AmrMultiGridLevel_t > BelosSolver_t
Definition: AmrMultiGrid.h:74
void initCrseFineInterp_m(const Interpolater &interface)
std::shared_ptr< preconditioner_t > prec_mp
preconditioner for bottom solver
Definition: AmrMultiGrid.h:687
amrex::Box box_t
Definition: AmrMultiGrid.h:84
std::unique_ptr< AmrInterpolater< AmrMultiGridLevel_t > > interface_mp
interpolater for coarse-fine interface
Definition: AmrMultiGrid.h:669
Norm convertToEnumNorm_m(const std::string &norm)
Interpolater convertToEnumInterpolater_m(const std::string &interp)
amr::global_ordinal_t go_t
Definition: AmrMultiGrid.h:53
void averageDown_m(const lo_t &level)
amr::matrix_t matrix_t
Definition: AmrMultiGrid.h:47
Smoother convertToEnumSmoother_m(const std::string &smoother)
amr::dmap_t dmap_t
Definition: AmrMultiGrid.h:50
double getZRangeMin(unsigned short level=0)
amr::vector_t vector_t
Definition: AmrMultiGrid.h:48
Amesos2BottomSolver< AmrMultiGridLevel_t > Amesos2Solver_t
Definition: AmrMultiGrid.h:75
void buildPotentialVector_m(const lo_t &level, const AmrField_t &phi)
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)
scalar_t eps_m
rhs scale for convergence
Definition: AmrMultiGrid.h:699
AmrMultiGridLevel_t::AmrField_t AmrField_t
Definition: AmrMultiGrid.h:58
std::size_t nIter_m
number of iterations till convergence
Definition: AmrMultiGrid.h:671
int lfine_m
fineste level
Definition: AmrMultiGrid.h:690
double getXRangeMax(unsigned short level=0)
Inform & print(Inform &os) const
AmrMultiGridLevel_t::AmrIntVect_t AmrIntVect_t
Definition: AmrMultiGrid.h:62
BottomSolver< Teuchos::RCP< matrix_t >, Teuchos::RCP< mv_t >, AmrMultiGridLevel_t > bsolver_t
Definition: AmrMultiGrid.h:72
bool verbose_m
If true, a SDDS file is written.
Definition: AmrMultiGrid.h:701
std::size_t bIter_m
number of iterations of bottom solver
Definition: AmrMultiGrid.h:672
int nlevel_m
number of levelss
Definition: AmrMultiGrid.h:691
amr::local_ordinal_t lo_t
Definition: AmrMultiGrid.h:52
std::ios_base::openmode flag_m
std::ios::out or std::ios::app
Definition: AmrMultiGrid.h:703
void setNumberOfSweeps(const std::size_t &nSweeps)
void initResidual_m(std::vector< scalar_t > &rhsNorms, std::vector< scalar_t > &resNorms)
AmrMultiGridLevel_t::AmrField_s AmrField_s
Definition: AmrMultiGrid.h:61
AmrSmoother::Smoother Smoother
Definition: AmrMultiGrid.h:88
void buildMultiLevel_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true)
void open_m(const lo_t &level, const bool &matrices)
void buildDensityVector_m(const lo_t &level, const AmrField_t &rho)
MueLuPreconditioner< AmrMultiGridLevel_t > MueLuPreconditioner_t
Definition: AmrMultiGrid.h:81
AmrMultiGridLevel< matrix_t, vector_t > AmrMultiGridLevel_t
Definition: AmrMultiGrid.h:56
std::size_t getNumIters()
std::size_t maxiter_m
maximum number of iterations allowed
Definition: AmrMultiGrid.h:673
void buildNoFinePoissonMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const scalar_t *invdx2)
double getYRangeMin(unsigned short level=0)
void writeSDDSData_m(const scalar_t &error)
int nBcPoints_m
maximum number of stencils points for BC
Definition: AmrMultiGrid.h:694
void setMaxNumberOfIterations(const std::size_t &maxiter)
void initBaseSolver_m(const BaseSolver &solver, const bool &rebalance, const std::string &reuse)
Smoother smootherType_m
type of smoother
Definition: AmrMultiGrid.h:675
Norm
Supported convergence criteria.
Definition: AmrMultiGrid.h:142
amrex::BoxArray boxarray_t
Definition: AmrMultiGrid.h:83
void initLevels_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrGeometry_t > &geom, bool regrid)
void smooth_m(const lo_t &level, Teuchos::RCP< vector_t > &e, Teuchos::RCP< vector_t > &r)
Ifpack2Preconditioner< AmrMultiGridLevel_t > Ifpack2Preconditioner_t
Definition: AmrMultiGrid.h:80
boundary_t bc_m[AMREX_SPACEDIM]
boundary conditions
Definition: AmrMultiGrid.h:693
void residual_m(const lo_t &level, Teuchos::RCP< vector_t > &r, const Teuchos::RCP< vector_t > &b, const Teuchos::RCP< vector_t > &x)
std::vector< std::shared_ptr< AmrSmoother > > smoother_m
error smoother
Definition: AmrMultiGrid.h:684
void buildInterpolationMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &cfab)
void computeEfield_m(AmrVectorFieldContainer_t &efield)
AmrMultiGridLevel_t::indices_t indices_t
Definition: AmrMultiGrid.h:63
Interpolater
Supported interpolaters for prolongation operation.
Definition: AmrMultiGrid.h:93
void prolongate_m(const lo_t &level)
bool isConverged_m(std::vector< scalar_t > &rhsNorms, std::vector< scalar_t > &resNorms)
Norm norm_m
norm for convergence criteria (l1, l2, linf)
Definition: AmrMultiGrid.h:696
double getYRangeMax(unsigned short level=0)
void setTolerance(const scalar_t &eps)
void solve(AmrScalarFieldContainer_t &rho, AmrScalarFieldContainer_t &phi, AmrVectorFieldContainer_t &efield, unsigned short baseLevel, unsigned short finestLevel, bool prevAsGuess=true)
AmrPreconditioner< matrix_t, AmrMultiGridLevel_t > preconditioner_t
Definition: AmrMultiGrid.h:78
MueLuBottomSolver< AmrMultiGridLevel_t > MueLuSolver_t
Definition: AmrMultiGrid.h:76
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)
Preconditioner convertToEnumPreconditioner_m(const std::string &prec)
void buildCompositePoissonMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const basefab_t &rfab, const scalar_t *invdx2)
amr::multivector_t mv_t
Definition: AmrMultiGrid.h:49
std::unordered_map< go_t, scalar_t > umap_t
std::unique_ptr< AmrField_t > AmrField_u
std::shared_ptr< AmrField_t > AmrField_s
amr::AmrField_t AmrField_t
std::vector< go_t > indices_t
Type for matrix indices.
amr::AmrIntVect_t AmrIntVect_t
amr::AmrGeometry_t AmrGeometry_t
std::shared_ptr< AmrBoundary< AmrMultiGridLevel< MatrixType, VectorType > > > boundary_t
std::vector< scalar_t > coefficients_t
Type for matrix entries.
Smoother
All supported Ifpack2 smoothers.
Definition: AmrSmoother.h:46
Definition: Inform.h:42
Timing::TimerRef TimerRef
Definition: IpplTimings.h:176