OPAL (Object Oriented Parallel Accelerator Library)
2021.1.99
OPAL
|
#include <AmrMultiGrid.h>
Public Member Functions | |
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 | solve (AmrScalarFieldContainer_t &rho, AmrScalarFieldContainer_t &phi, AmrVectorFieldContainer_t &efield, unsigned short baseLevel, unsigned short finestLevel, bool prevAsGuess=true) |
void | setNumberOfSweeps (const std::size_t &nSweeps) |
void | setMaxNumberOfIterations (const std::size_t &maxiter) |
std::size_t | getNumIters () |
scalar_t | getLevelResidualNorm (lo_t level) |
void | setVerbose (bool verbose) |
void | setTolerance (const scalar_t &eps) |
double | getXRangeMin (unsigned short level=0) |
double | getXRangeMax (unsigned short level=0) |
double | getYRangeMin (unsigned short level=0) |
double | getYRangeMax (unsigned short level=0) |
double | getZRangeMin (unsigned short level=0) |
double | getZRangeMax (unsigned short level=0) |
Inform & | print (Inform &os) const |
Public Member Functions inherited from AmrPoissonSolver< AmrBoxLib > | |
AmrPoissonSolver (AmrBoxLib *itsAmrObject_p) | |
virtual | ~AmrPoissonSolver () |
void | computePotential (Field_t &, Vector_t) |
void | computePotential (Field_t &, Vector_t, double) |
void | test (PartBunchBase< double, 3 > *) |
void | hasToRegrid () |
Public Member Functions inherited from PoissonSolver | |
virtual | ~PoissonSolver () |
virtual void | resizeMesh (Vector_t &, Vector_t &, const Vector_t &, const Vector_t &, double) |
Private Member Functions | |
void | initPhysicalBoundary_m (const Boundary *bc) |
void | initLevels_m (const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrGeometry_t > &geom, bool regrid) |
void | clearMasks_m () |
void | initGuess_m (bool reset) |
scalar_t | iterate_m () |
bool | isConverged_m (std::vector< scalar_t > &rhsNorms, std::vector< scalar_t > &resNorms) |
void | residual_m (const lo_t &level, Teuchos::RCP< vector_t > &r, const Teuchos::RCP< vector_t > &b, const Teuchos::RCP< vector_t > &x) |
void | relax_m (const lo_t &level) |
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 | evalNorm_m (const Teuchos::RCP< const vector_t > &x) |
void | initResidual_m (std::vector< scalar_t > &rhsNorms, std::vector< scalar_t > &resNorms) |
void | computeEfield_m (AmrVectorFieldContainer_t &efield) |
void | setup_m (const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true) |
void | buildSingleLevel_m (const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true) |
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 | close_m (const lo_t &level, const bool &matrices) |
void | buildNoFinePoissonMatrix_m (const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const scalar_t *invdx2) |
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) |
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) |
void | buildInterpolationMatrix_m (const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &cfab) |
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) |
void | buildFineBoundaryMatrix_m (const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const basefab_t &rfab) |
void | buildDensityVector_m (const lo_t &level, const AmrField_t &rho) |
void | buildPotentialVector_m (const lo_t &level, const AmrField_t &phi) |
void | buildGradientMatrix_m (const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const scalar_t *invdx) |
void | amrex2trilinos_m (const lo_t &level, const lo_t &comp, const AmrField_t &mf, Teuchos::RCP< vector_t > &mv) |
void | trilinos2amrex_m (const lo_t &level, const lo_t &comp, AmrField_t &mf, const Teuchos::RCP< vector_t > &mv) |
void | map2vector_m (umap_t &map, indices_t &indices, coefficients_t &values) |
void | smooth_m (const lo_t &level, Teuchos::RCP< vector_t > &e, Teuchos::RCP< vector_t > &r) |
void | restrict_m (const lo_t &level) |
void | prolongate_m (const lo_t &level) |
void | averageDown_m (const lo_t &level) |
void | initInterpolater_m (const Interpolater &interp) |
void | initCrseFineInterp_m (const Interpolater &interface) |
void | initBaseSolver_m (const BaseSolver &solver, const bool &rebalance, const std::string &reuse) |
void | initPrec_m (const Preconditioner &prec, const std::string &reuse) |
Boundary | convertToEnumBoundary_m (const std::string &bc) |
Interpolater | convertToEnumInterpolater_m (const std::string &interp) |
BaseSolver | convertToEnumBaseSolver_m (const std::string &bsolver) |
Preconditioner | convertToEnumPreconditioner_m (const std::string &prec) |
Smoother | convertToEnumSmoother_m (const std::string &smoother) |
Norm | convertToEnumNorm_m (const std::string &norm) |
void | writeSDDSHeader_m (std::ofstream &outfile) |
void | writeSDDSData_m (const scalar_t &error) |
Private Attributes | |
Teuchos::RCP< comm_t > | comm_mp |
communicator More... | |
std::unique_ptr< AmrInterpolater< AmrMultiGridLevel_t > > | interp_mp |
interpolater without coarse-fine interface More... | |
std::unique_ptr< AmrInterpolater< AmrMultiGridLevel_t > > | interface_mp |
interpolater for coarse-fine interface More... | |
std::size_t | nIter_m |
number of iterations till convergence More... | |
std::size_t | bIter_m |
number of iterations of bottom solver More... | |
std::size_t | maxiter_m |
maximum number of iterations allowed More... | |
std::size_t | nSweeps_m |
number of smoothing iterations More... | |
Smoother | smootherType_m |
type of smoother More... | |
std::vector< std::unique_ptr< AmrMultiGridLevel_t > > | mglevel_m |
container for levels More... | |
std::shared_ptr< bsolver_t > | solver_mp |
bottom solver More... | |
std::vector< std::shared_ptr< AmrSmoother > > | smoother_m |
error smoother More... | |
std::shared_ptr< preconditioner_t > | prec_mp |
preconditioner for bottom solver More... | |
int | lbase_m |
base level (currently only 0 supported) More... | |
int | lfine_m |
fineste level More... | |
int | nlevel_m |
number of levelss More... | |
boundary_t | bc_m [AMREX_SPACEDIM] |
boundary conditions More... | |
int | nBcPoints_m |
maximum number of stencils points for BC More... | |
Norm | norm_m |
norm for convergence criteria (l1, l2, linf) More... | |
std::string | snorm_m |
norm for convergence criteria More... | |
scalar_t | eps_m |
rhs scale for convergence More... | |
bool | verbose_m |
If true, a SDDS file is written. More... | |
std::string | fname_m |
SDDS filename. More... | |
std::ios_base::openmode | flag_m |
std::ios::out or std::ios::app More... | |
Additional Inherited Members | |
Protected Types inherited from PoissonSolver | |
typedef Field< int, 3, Mesh_t, Center_t > | IField_t |
typedef Field< std::complex< double >, 3, Mesh_t, Center_t > | CxField_t |
Protected Attributes inherited from AmrPoissonSolver< AmrBoxLib > | |
AmrBoxLib * | itsAmrObject_mp |
bool | regrid_m |
is set to true by itsAmrObject_mp and reset to false by solver More... | |
Definition at line 44 of file AmrMultiGrid.h.
Definition at line 75 of file AmrMultiGrid.h.
Definition at line 61 of file AmrMultiGrid.h.
Definition at line 58 of file AmrMultiGrid.h.
Definition at line 60 of file AmrMultiGrid.h.
Definition at line 59 of file AmrMultiGrid.h.
Definition at line 62 of file AmrMultiGrid.h.
Definition at line 56 of file AmrMultiGrid.h.
typedef amrex::BaseFab<int> AmrMultiGrid::basefab_t |
Definition at line 85 of file AmrMultiGrid.h.
Definition at line 74 of file AmrMultiGrid.h.
Definition at line 66 of file AmrMultiGrid.h.
typedef amrex::Box AmrMultiGrid::box_t |
Definition at line 84 of file AmrMultiGrid.h.
typedef amrex::BoxArray AmrMultiGrid::boxarray_t |
Definition at line 83 of file AmrMultiGrid.h.
typedef BottomSolver< Teuchos::RCP<matrix_t>, Teuchos::RCP<mv_t>, AmrMultiGridLevel_t > AmrMultiGrid::bsolver_t |
Definition at line 72 of file AmrMultiGrid.h.
Definition at line 64 of file AmrMultiGrid.h.
typedef amr::comm_t AmrMultiGrid::comm_t |
Definition at line 51 of file AmrMultiGrid.h.
typedef amr::dmap_t AmrMultiGrid::dmap_t |
Definition at line 50 of file AmrMultiGrid.h.
typedef amrex::FArrayBox AmrMultiGrid::farraybox_t |
Definition at line 86 of file AmrMultiGrid.h.
Definition at line 53 of file AmrMultiGrid.h.
Definition at line 80 of file AmrMultiGrid.h.
Definition at line 63 of file AmrMultiGrid.h.
Definition at line 52 of file AmrMultiGrid.h.
typedef amr::matrix_t AmrMultiGrid::matrix_t |
Definition at line 47 of file AmrMultiGrid.h.
Definition at line 81 of file AmrMultiGrid.h.
Definition at line 76 of file AmrMultiGrid.h.
typedef amr::multivector_t AmrMultiGrid::mv_t |
Definition at line 49 of file AmrMultiGrid.h.
Definition at line 90 of file AmrMultiGrid.h.
Definition at line 78 of file AmrMultiGrid.h.
typedef amr::scalar_t AmrMultiGrid::scalar_t |
Definition at line 54 of file AmrMultiGrid.h.
Definition at line 88 of file AmrMultiGrid.h.
Definition at line 65 of file AmrMultiGrid.h.
typedef amr::vector_t AmrMultiGrid::vector_t |
Definition at line 48 of file AmrMultiGrid.h.
Supported bottom solvers.
Enumerator | |
---|---|
BICGSTAB | |
MINRES | |
PCPG | |
CG | |
GMRES | |
STOCHASTIC_CG | |
RECYCLING_CG | |
RECYCLING_GMRES | |
SA |
Definition at line 100 of file AmrMultiGrid.h.
Supported physical boundaries.
Enumerator | |
---|---|
DIRICHLET | |
OPEN | |
PERIODIC |
Definition at line 135 of file AmrMultiGrid.h.
Supported interpolaters for prolongation operation.
Enumerator | |
---|---|
TRILINEAR | |
LAGRANGE | |
PIECEWISE_CONST |
Definition at line 93 of file AmrMultiGrid.h.
enum AmrMultiGrid::Norm |
AmrMultiGrid::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 | ||
) |
Instantiation used in Structure/FieldSolver.cpp
itsAmrObject_p | has information about refinemen ratios, etc. |
bsolver | bottom solver |
prec | preconditioner for bottom solver |
rebalance | of preconditioner (SA only) |
bcx | boundary condition in x |
bcy | boundary condition in y |
bcz | boundary condition in z |
smoother | for level solution |
nSweeps | when smoothing |
interp | interpolater between levels |
norm | for convergence criteria |
Definition at line 42 of file AmrMultiGrid.cpp.
References convertToEnumBaseSolver_m(), convertToEnumBoundary_m(), convertToEnumInterpolater_m(), convertToEnumNorm_m(), convertToEnumPreconditioner_m(), convertToEnumSmoother_m(), endl(), flag_m, fname_m, INFOMSG, initBaseSolver_m(), initCrseFineInterp_m(), initInterpolater_m(), initPhysicalBoundary_m(), initPrec_m(), norm_m, and smootherType_m.
|
private |
Data transfer from AMReX to Trilinos.
mf | is the multifab of a level |
comp | component to copy |
mv | is the vector to be filled |
level | where to perform |
Definition at line 1656 of file AmrMultiGrid.cpp.
References mglevel_m.
Referenced by buildDensityVector_m(), and buildPotentialVector_m().
|
private |
Average data from fine level to coarse
level | finest level is omitted |
Definition at line 1844 of file AmrMultiGrid.cpp.
References lfine_m, mglevel_m, IpplTimings::startTimer(), and IpplTimings::stopTimer().
Referenced by solve().
|
private |
Build the Poisson matrix for a level that got refined (it does not take the covered cells (covered by fine cells) into account). The finest level does not build such a matrix. It takes care of physical boundaries (i.e. mesh boundary). Internal boundaries (i.e. boundaries due to crse-fine interfaces) are treated by the boundary matrix.
iv | is the current cell |
mfab | is the mask (internal cell, boundary cell, ...) of that level |
rfab | is the mask between levels |
level | for which we build the special Poisson matrix |
Definition at line 1111 of file AmrMultiGrid.cpp.
References interface_mp, lbase_m, and mglevel_m.
Referenced by buildMultiLevel_m().
|
private |
The boundary values at the crse-fine-interface need to be taken into account properly. This matrix is used to update the fine boundary values from the coarse values, i.e.
\[ x^{(l)} = B_{crse}\cdot x^{(l-1)} \]
Dirichlet boundary condition
iv | is the current cell |
mfab | is the mask (internal cell, boundary cell, ...) of that level |
cells | all fine cells that are at the crse-fine interface |
level | the base level is omitted |
Definition at line 1358 of file AmrMultiGrid.cpp.
References interface_mp, lbase_m, and mglevel_m.
Referenced by buildMultiLevel_m().
|
inlineprivate |
Copy data from AMReX to Trilinos
rho | is the charge density |
level | for which to copy |
Definition at line 1564 of file AmrMultiGrid.cpp.
References amrex2trilinos_m(), and mglevel_m.
Referenced by buildMultiLevel_m(), and buildSingleLevel_m().
|
private |
The boundary values at the crse-fine-interface need to be taken into account properly. This matrix is used to update the coarse boundary values from fine values, i.e.
\[ x^{(l)} = B_{fine}\cdot x^{(l+1)} \]
Dirichlet boundary condition. Flux matching.
level | the finest level is omitted |
gidx | the global index |
iv | is the current cell |
mfab | is the mask (internal cell, boundary cell, ...) of that level |
rfab | is the mask between levels |
Definition at line 1423 of file AmrMultiGrid.cpp.
References begin(), end(), interface_mp, lfine_m, and mglevel_m.
Referenced by buildMultiLevel_m().
|
private |
Gradient matrix is used to compute the electric field
iv | is the current cell |
mfab | is the mask (internal cell, boundary cell, ...) of that level |
level | for which to compute |
Definition at line 1578 of file AmrMultiGrid.cpp.
References lbase_m, map2vector_m(), and mglevel_m.
Referenced by buildMultiLevel_m(), and buildSingleLevel_m().
|
private |
Interpolate data from coarse cells to appropriate refined cells. The interpolation scheme is allowed only to have local cells in the stencil, i.e. 2D --> 4, 3D --> 8.
\[ x^{(l)} = I\cdot x^{(l-1)} \]
iv | is the current cell |
level | for which to build the interpolation matrix. The finest level does not build such a matrix. |
Definition at line 1321 of file AmrMultiGrid.cpp.
References interp_mp, lbase_m, and mglevel_m.
Referenced by buildMultiLevel_m().
|
private |
Build all matrices and vectors needed for multi-level computation
Definition at line 778 of file AmrMultiGrid.cpp.
References buildCompositePoissonMatrix_m(), buildCrseBoundaryMatrix_m(), buildDensityVector_m(), buildFineBoundaryMatrix_m(), buildGradientMatrix_m(), buildInterpolationMatrix_m(), buildNoFinePoissonMatrix_m(), buildRestrictionMatrix_m(), close_m(), lbase_m, mglevel_m, nlevel_m, nSweeps_m, open_m(), smoother_m, smootherType_m, and solver_mp.
Referenced by setup_m().
|
private |
Build the Poisson matrix for a level assuming no finer level (i.e. the whole fine mesh is taken into account). It takes care of physical boundaries (i.e. mesh boundary). Internal boundaries (i.e. boundaries due to crse-fine interfaces) are treated by the boundary matrix.
iv | is the current cell |
mfab | is the mask (internal cell, boundary cell, ...) of that level |
level | for which we build the Poisson matrix |
Definition at line 1030 of file AmrMultiGrid.cpp.
References interface_mp, lbase_m, and mglevel_m.
Referenced by buildMultiLevel_m(), and buildSingleLevel_m().
|
inlineprivate |
Copy data from AMReX to Trilinos
phi | is the potential |
level | for which to copy |
Definition at line 1571 of file AmrMultiGrid.cpp.
References amrex2trilinos_m(), and mglevel_m.
|
private |
Build a matrix that averages down the data of the fine cells down to the corresponding coarse cells. The base level does not build such a matrix.
\[ x^{(l)} = R\cdot x^{(l+1)} \]
iv | is the current cell |
rfab | is the mask between levels |
level | for which to build restriction matrix |
Definition at line 1269 of file AmrMultiGrid.cpp.
References lfine_m, mglevel_m, and serialize().
Referenced by buildMultiLevel_m().
|
private |
Build all matrices and vectors needed for single-level computation
Definition at line 717 of file AmrMultiGrid.cpp.
References buildDensityVector_m(), buildGradientMatrix_m(), buildNoFinePoissonMatrix_m(), close_m(), lbase_m, mglevel_m, open_m(), and solver_mp.
Referenced by setup_m().
|
private |
Clear masks (required to build matrices) no longer needed.
Definition at line 317 of file AmrMultiGrid.cpp.
References mglevel_m, and nlevel_m.
Referenced by setup_m().
|
private |
Call fill complete
level | for which we filled matrix |
matrices | if we set matrices as well. |
Definition at line 992 of file AmrMultiGrid.cpp.
References lbase_m, lfine_m, mglevel_m, and solver_mp.
Referenced by buildMultiLevel_m(), and buildSingleLevel_m().
|
private |
efield | to compute |
Definition at line 664 of file AmrMultiGrid.cpp.
References lbase_m, mglevel_m, nlevel_m, IpplTimings::startTimer(), IpplTimings::stopTimer(), and trilinos2amrex_m().
Referenced by solve().
|
private |
Converts string to enum BaseSolver
bsolver | bottom solver |
Definition at line 2041 of file AmrMultiGrid.cpp.
References amr::SA.
Referenced by AmrMultiGrid().
|
private |
Convertstring to enum Boundary
bc | boundary condition |
Definition at line 2008 of file AmrMultiGrid.cpp.
Referenced by AmrMultiGrid().
|
private |
Converts string to enum Interpolater
interp | interpolater |
Definition at line 2024 of file AmrMultiGrid.cpp.
Referenced by AmrMultiGrid().
|
private |
Converts string to enum Norm
norm | either L1, L2, LInf |
Definition at line 2107 of file AmrMultiGrid.cpp.
References Hypervolume::n.
Referenced by AmrMultiGrid().
|
private |
Converts string to enum Preconditioner
prec | preconditioner |
Definition at line 2082 of file AmrMultiGrid.cpp.
References Ifpack2Preconditioner< Level >::fillMap(), MueLuPreconditioner< Level >::fillMap(), and EmissionModelT::NONE.
Referenced by AmrMultiGrid().
|
private |
Converts string to enum Smoother
smoother | of level solution |
Definition at line 2101 of file AmrMultiGrid.cpp.
References AmrSmoother::convertToEnumSmoother().
Referenced by AmrMultiGrid().
|
private |
Vector norm computation.
x | is the vector for which we compute the norm |
Definition at line 595 of file AmrMultiGrid.cpp.
References norm_m.
Referenced by getLevelResidualNorm(), and initResidual_m().
AmrMultiGrid::scalar_t AmrMultiGrid::getLevelResidualNorm | ( | lo_t | level | ) |
Obtain the residual norm of a level
level | for which error is requested |
Definition at line 184 of file AmrMultiGrid.cpp.
References evalNorm_m(), and mglevel_m.
Referenced by iterate_m().
std::size_t AmrMultiGrid::getNumIters | ( | ) |
Obtain some convergence info
Definition at line 179 of file AmrMultiGrid.cpp.
References nIter_m.
|
virtual |
Implements PoissonSolver.
Definition at line 2241 of file AmrMultiGrid.cpp.
References AmrPoissonSolver< AmrBoxLib >::itsAmrObject_mp.
|
virtual |
Implements PoissonSolver.
Definition at line 2236 of file AmrMultiGrid.cpp.
References AmrPoissonSolver< AmrBoxLib >::itsAmrObject_mp.
|
virtual |
Implements PoissonSolver.
Definition at line 2251 of file AmrMultiGrid.cpp.
References AmrPoissonSolver< AmrBoxLib >::itsAmrObject_mp.
|
virtual |
Implements PoissonSolver.
Definition at line 2246 of file AmrMultiGrid.cpp.
References AmrPoissonSolver< AmrBoxLib >::itsAmrObject_mp.
|
virtual |
Implements PoissonSolver.
Definition at line 2261 of file AmrMultiGrid.cpp.
References AmrPoissonSolver< AmrBoxLib >::itsAmrObject_mp.
|
virtual |
Implements PoissonSolver.
Definition at line 2256 of file AmrMultiGrid.cpp.
References AmrPoissonSolver< AmrBoxLib >::itsAmrObject_mp.
|
private |
Instantiate a bottom solver
solver | type |
rebalance | solver (SA only) |
reuse | types of SA hierarchy |
Definition at line 1904 of file AmrMultiGrid.cpp.
References MueLuBottomSolver< Level >::convertToMueLuReuseOption(), prec_mp, amr::SA, and solver_mp.
Referenced by AmrMultiGrid().
|
private |
Instantiate interface interpolater
interface | handler |
Definition at line 1885 of file AmrMultiGrid.cpp.
References interface_mp.
Referenced by AmrMultiGrid().
|
private |
Reset potential to zero (currently)
reset | solution to initial guess (zero) |
Definition at line 326 of file AmrMultiGrid.cpp.
References mglevel_m, and nlevel_m.
Referenced by solve().
|
private |
Instantiate interpolater
interp | interpolater type |
Definition at line 1867 of file AmrMultiGrid.cpp.
References interp_mp.
Referenced by AmrMultiGrid().
|
private |
Instantiate all levels and set boundary conditions
rho | is the charge density |
geom | is the geometry |
regrid | was performed |
Definition at line 227 of file AmrMultiGrid.cpp.
References bc_m, comm_mp, AmrBoxLib::getMeshScaling(), AmrPoissonSolver< AmrBoxLib >::itsAmrObject_mp, lbase_m, mglevel_m, and nlevel_m.
Referenced by solve().
|
private |
Instantiate boundary object
bc | boundary conditions @precondition length must be equal to AMREX_SPACEDIM |
Definition at line 199 of file AmrMultiGrid.cpp.
References bc_m, and nBcPoints_m.
Referenced by AmrMultiGrid().
|
private |
Instantiate a preconditioner for the bottom solver
precond | type |
reuse | types of SA hierarchy |
Definition at line 1978 of file AmrMultiGrid.cpp.
References amr::BLOCK_GS, amr::BLOCK_JACOBI, amr::CHEBYSHEV, MueLuPreconditioner< Level >::convertToMueLuReuseOption(), amr::GS, amr::ILUT, amr::JACOBI, EmissionModelT::NONE, prec_mp, amr::RILUK, and amr::SA.
Referenced by AmrMultiGrid().
|
private |
Initial convergence criteria values.
rhsNorms | per level of right-hand side (is filled) |
resNorms | per level of residual (is filled) |
Definition at line 623 of file AmrMultiGrid.cpp.
References comm_mp, endl(), evalNorm_m(), mglevel_m, nlevel_m, and residual_m().
Referenced by iterate_m().
|
private |
Compute norms / level and check convergence
Definition at line 391 of file AmrMultiGrid.cpp.
Referenced by iterate_m().
|
private |
Actual solve.
Definition at line 336 of file AmrMultiGrid.cpp.
References bIter_m, for_each(), getLevelResidualNorm(), initResidual_m(), isConverged_m(), lfine_m, maxiter_m, mglevel_m, nIter_m, nlevel_m, relax_m(), residual_m(), and solver_mp.
Referenced by solve().
|
inlineprivate |
Some indices might occur several times due to boundary conditions, etc. We avoid this by filling a map and then copying the data to a vector for filling the matrices. The map gets cleared inside the function.
indices | in matrix |
values | are the coefficients |
Definition at line 1725 of file AmrMultiGrid.cpp.
References for_each().
Referenced by buildGradientMatrix_m().
|
private |
Set matrix and vector pointer
level | for which we fill matrix + vector |
matrices | if we need to set matrices as well or only vectors. |
Definition at line 876 of file AmrMultiGrid.cpp.
References interface_mp, interp_mp, lbase_m, lfine_m, mglevel_m, nBcPoints_m, and solver_mp.
Referenced by buildMultiLevel_m(), and buildSingleLevel_m().
Print information abour tolerances.
os | output stream where to write to |
Definition at line 2266 of file AmrMultiGrid.cpp.
References endl().
|
private |
Update error of fine level based on error of coarse level
level | to update |
Definition at line 1827 of file AmrMultiGrid.cpp.
References mglevel_m, IpplTimings::startTimer(), and IpplTimings::stopTimer().
Referenced by relax_m().
|
private |
Recursive call.
level | to relax |
Definition at line 463 of file AmrMultiGrid.cpp.
References lbase_m, lfine_m, mglevel_m, prolongate_m(), residual_no_fine_m(), restrict_m(), smooth_m(), solver_mp, IpplTimings::startTimer(), and IpplTimings::stopTimer().
Referenced by iterate_m().
|
private |
Compute composite residual of a level
r | is the residual to compute |
b | is the right-hand side |
x | is the left-hand side |
level | to solve for |
Definition at line 400 of file AmrMultiGrid.cpp.
References lfine_m, and mglevel_m.
Referenced by initResidual_m(), and iterate_m().
|
private |
Compute the residual of a level without considering refined level.
result | is computed |
rhs | is the right-hand side |
crs_rhs | is the coarse right-hand side for internal boundary |
b | is the left-hand side |
level | to solver for |
Definition at line 550 of file AmrMultiGrid.cpp.
References mglevel_m.
Referenced by relax_m(), and restrict_m().
|
private |
Restrict coarse level residual based on fine level residual
level | to restrict |
Definition at line 1770 of file AmrMultiGrid.cpp.
References mglevel_m, residual_no_fine_m(), IpplTimings::startTimer(), and IpplTimings::stopTimer().
Referenced by relax_m().
void AmrMultiGrid::setMaxNumberOfIterations | ( | const std::size_t & | maxiter | ) |
Specify the maximum number of iterations
maxiter | \( [0, \infty[ \) |
Definition at line 170 of file AmrMultiGrid.cpp.
References maxiter_m.
void AmrMultiGrid::setNumberOfSweeps | ( | const std::size_t & | nSweeps | ) |
Specify the number of smoothing steps
nSweeps | for each smoothing step |
Definition at line 165 of file AmrMultiGrid.cpp.
References nSweeps_m.
void AmrMultiGrid::setTolerance | ( | const scalar_t & | eps | ) |
|
private |
Build all matrices and vectors, i.e. AMReX to Trilinos
matrices | if we need to build matrices as well or only vectors. |
Definition at line 688 of file AmrMultiGrid.cpp.
References buildMultiLevel_m(), buildSingleLevel_m(), clearMasks_m(), lbase_m, lfine_m, mglevel_m, solver_mp, and IpplTimings::startTimer().
Referenced by solve().
void AmrMultiGrid::setVerbose | ( | bool | verbose | ) |
Enable solver info dumping into SDDS file
Definition at line 189 of file AmrMultiGrid.cpp.
References verbose_m.
|
private |
Perform one smoothing step
e | error to update (left-hand side) |
r | residual (right-hand side) |
level | on which to relax |
Definition at line 1753 of file AmrMultiGrid.cpp.
References Physics::e, smoother_m, IpplTimings::startTimer(), and IpplTimings::stopTimer().
Referenced by relax_m().
|
virtual |
Used in OPAL
rho | right-hand side charge density on grid [C / m] |
phi | electrostatic potential (unknown) [V] |
efield | electric field [V / m] |
baseLevel | for solve |
finestLevel | for solve |
prevAsGuess | use of previous solution as initial guess |
Reimplemented from PoissonSolver.
Definition at line 110 of file AmrMultiGrid.cpp.
References averageDown_m(), computeEfield_m(), initGuess_m(), initLevels_m(), iterate_m(), AmrPoissonSolver< AmrBoxLib >::itsAmrObject_mp, lbase_m, lfine_m, mglevel_m, nlevel_m, AmrPoissonSolver< AmrBoxLib >::regrid_m, setup_m(), trilinos2amrex_m(), verbose_m, and writeSDDSData_m().
|
private |
Data transfer from Trilinos to AMReX.
level | to copy |
comp | component to copy |
mf | is the multifab to be filled |
mv | is the corresponding Trilinos vector |
Definition at line 1690 of file AmrMultiGrid.cpp.
References mglevel_m.
Referenced by computeEfield_m(), and solve().
|
private |
SDDS data write (done by root core)
error | to write |
Definition at line 2192 of file AmrMultiGrid.cpp.
References bIter_m, comm_mp, flag_m, fname_m, AmrBoxLib::getT(), AmrPoissonSolver< AmrBoxLib >::itsAmrObject_mp, nIter_m, AmrPoissonSolver< AmrBoxLib >::regrid_m, IpplTimings::startTimer(), IpplTimings::stopTimer(), and writeSDDSHeader_m().
Referenced by solve().
|
private |
SDDS header is written by root core
outfile | output stream |
Definition at line 2123 of file AmrMultiGrid.cpp.
References comm_mp, OPALTimer::Timer::date(), endl(), Util::getGitRevision(), OpalData::getInputFn(), OpalData::getInstance(), OpalData::isInOPALCyclMode(), OpalData::isInOPALTMode(), OPAL_PROJECT_NAME, OPAL_PROJECT_VERSION, snorm_m, and OPALTimer::Timer::time().
Referenced by writeSDDSData_m().
|
private |
boundary conditions
Definition at line 693 of file AmrMultiGrid.h.
Referenced by initLevels_m(), and initPhysicalBoundary_m().
|
private |
number of iterations of bottom solver
Definition at line 672 of file AmrMultiGrid.h.
Referenced by iterate_m(), and writeSDDSData_m().
|
private |
communicator
Definition at line 663 of file AmrMultiGrid.h.
Referenced by initLevels_m(), initResidual_m(), writeSDDSData_m(), and writeSDDSHeader_m().
|
private |
rhs scale for convergence
Definition at line 699 of file AmrMultiGrid.h.
Referenced by setTolerance().
|
private |
std::ios::out or std::ios::app
Definition at line 703 of file AmrMultiGrid.h.
Referenced by AmrMultiGrid(), and writeSDDSData_m().
|
private |
SDDS filename.
Definition at line 702 of file AmrMultiGrid.h.
Referenced by AmrMultiGrid(), and writeSDDSData_m().
|
private |
interpolater for coarse-fine interface
Definition at line 669 of file AmrMultiGrid.h.
Referenced by buildCompositePoissonMatrix_m(), buildCrseBoundaryMatrix_m(), buildFineBoundaryMatrix_m(), buildNoFinePoissonMatrix_m(), initCrseFineInterp_m(), and open_m().
|
private |
interpolater without coarse-fine interface
Definition at line 666 of file AmrMultiGrid.h.
Referenced by buildInterpolationMatrix_m(), initInterpolater_m(), and open_m().
|
private |
base level (currently only 0 supported)
Definition at line 689 of file AmrMultiGrid.h.
Referenced by buildCompositePoissonMatrix_m(), buildCrseBoundaryMatrix_m(), buildGradientMatrix_m(), buildInterpolationMatrix_m(), buildMultiLevel_m(), buildNoFinePoissonMatrix_m(), buildSingleLevel_m(), close_m(), computeEfield_m(), initLevels_m(), open_m(), relax_m(), setup_m(), and solve().
|
private |
fineste level
Definition at line 690 of file AmrMultiGrid.h.
Referenced by averageDown_m(), buildFineBoundaryMatrix_m(), buildRestrictionMatrix_m(), close_m(), iterate_m(), open_m(), relax_m(), residual_m(), setup_m(), and solve().
|
private |
maximum number of iterations allowed
Definition at line 673 of file AmrMultiGrid.h.
Referenced by iterate_m(), and setMaxNumberOfIterations().
|
private |
container for levels
Definition at line 678 of file AmrMultiGrid.h.
Referenced by amrex2trilinos_m(), averageDown_m(), buildCompositePoissonMatrix_m(), buildCrseBoundaryMatrix_m(), buildDensityVector_m(), buildFineBoundaryMatrix_m(), buildGradientMatrix_m(), buildInterpolationMatrix_m(), buildMultiLevel_m(), buildNoFinePoissonMatrix_m(), buildPotentialVector_m(), buildRestrictionMatrix_m(), buildSingleLevel_m(), clearMasks_m(), close_m(), computeEfield_m(), getLevelResidualNorm(), initGuess_m(), initLevels_m(), initResidual_m(), iterate_m(), open_m(), prolongate_m(), relax_m(), residual_m(), residual_no_fine_m(), restrict_m(), setup_m(), solve(), and trilinos2amrex_m().
|
private |
maximum number of stencils points for BC
Definition at line 694 of file AmrMultiGrid.h.
Referenced by initPhysicalBoundary_m(), and open_m().
|
private |
number of iterations till convergence
Definition at line 671 of file AmrMultiGrid.h.
Referenced by getNumIters(), iterate_m(), and writeSDDSData_m().
|
private |
number of levelss
Definition at line 691 of file AmrMultiGrid.h.
Referenced by buildMultiLevel_m(), clearMasks_m(), computeEfield_m(), initGuess_m(), initLevels_m(), initResidual_m(), iterate_m(), and solve().
|
private |
norm for convergence criteria (l1, l2, linf)
Definition at line 696 of file AmrMultiGrid.h.
Referenced by AmrMultiGrid(), and evalNorm_m().
|
private |
number of smoothing iterations
Definition at line 674 of file AmrMultiGrid.h.
Referenced by buildMultiLevel_m(), and setNumberOfSweeps().
|
private |
preconditioner for bottom solver
Definition at line 687 of file AmrMultiGrid.h.
Referenced by initBaseSolver_m(), and initPrec_m().
|
private |
error smoother
Definition at line 684 of file AmrMultiGrid.h.
Referenced by buildMultiLevel_m(), and smooth_m().
|
private |
type of smoother
Definition at line 675 of file AmrMultiGrid.h.
Referenced by AmrMultiGrid(), and buildMultiLevel_m().
|
private |
norm for convergence criteria
Definition at line 697 of file AmrMultiGrid.h.
Referenced by writeSDDSHeader_m().
|
private |
bottom solver
Definition at line 681 of file AmrMultiGrid.h.
Referenced by buildMultiLevel_m(), buildSingleLevel_m(), close_m(), initBaseSolver_m(), iterate_m(), open_m(), relax_m(), and setup_m().
|
private |
If true, a SDDS file is written.
Definition at line 701 of file AmrMultiGrid.h.
Referenced by setVerbose(), and solve().