#include <source/solvers/poisson/CellPoissonHypreSolver.h>
Public Member Functions | |
CellPoissonHypreSolver (const std::string &object_name, tbox::Pointer< tbox::Database > database=NULL) | |
Constructor. | |
~CellPoissonHypreSolver () | |
void | initializeSolverState (tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy, int ln=0) |
Initialize to a given hierarchy. | |
void | deallocateSolverState () |
Reset to an uninitialized state. | |
void | setMatrixCoefficients (const PoissonSpecifications &spec) |
Set the matrix coefficients. | |
void | setSolnIdDepth (const int depth) |
Set default depth of the solution data involved in the solve. | |
void | setRhsIdDepth (const int depth) |
Set default depth of the rhs data involved in the solve. | |
void | setStoppingCriteria (const int max_iterations=10, const double relative_residual_tol=1.0e-6) |
Set the stopping criteria (max iterations and residual tolerance) for the linear solver. | |
int | solveSystem (const int u, const int f, bool homogeneous_bc=false) |
Solve the linear system Au=f. | |
int | getNumberOfIterations () const |
Return the number of iterations taken by the solver to converge. | |
void | setNumPreRelaxSteps (const int steps) |
Set the number of pre-relax steps used by the Hypre solve. | |
void | setNumPostRelaxSteps (const int steps) |
Set the number of post-relax steps used by the Hypre solve. | |
double | getRelativeResidualNorm () const |
Return the final residual norm returned by the Hypre solve. | |
void | setUseSMG (bool use_smg) |
Set whether to use Hypre's PFMG algorithm instead of the SMG algorithm. | |
void | setBoundaries (const std::string &boundary_type, const int fluxes=-1, const int flags=-1, int *bdry_types=NULL) |
Specify boundary condition directly, without using a RobinBcCoefStrategy<DIM> object. | |
void | setPhysicalBcCoefObject (const RobinBcCoefStrategy< DIM > *physical_bc_coef_strategy, const tbox::Pointer< hier::Variable< DIM > > variable=NULL) |
Specify boundary condition through the use of a Robin boundary condition object. | |
void | setPrintSolverInfo (const bool print) |
Set the flag for printing solver information. |
Class CellPoissonHypreSolver<DIM> uses the HYPRE preconditioner library to solve linear equations of the form , where C is a cell-centered array, D is a face-centered array, and u and f are cell-centered arrays (see PoissonSpecifications). The discretization is the standard second order finite difference stencil.
Robin boundary conditions are used through the interface class RobinBcCoefStrategy<DIM>. Periodic boundary conditions are not supported yet.
The user must perform the following steps to use CellPoissonHypreSolver<DIM>:
Sample parameters for initialization from database (and their default values):
* print_solver_info = FALSE // Whether to print some data for debugging * max_iterations = 10 // Max iterations used by Hypre * relative_residual_tol = 1.0e-8 // Residual tolerance used by Hypre * num_pre_relax_steps = 1 // # of presmoothing steps used by Hypre * num_post_relax_steps = 1 // # of postsmoothing steps used by Hypre * use_smg = FALSE // Whether to use hypre's smg solver * // (alternative is the pfmg solver) *
SAMRAI::solv::CellPoissonHypreSolver< DIM >::CellPoissonHypreSolver | ( | const std::string & | object_name, | |
tbox::Pointer< tbox::Database > | database = NULL | |||
) |
SAMRAI::solv::CellPoissonHypreSolver< DIM >::~CellPoissonHypreSolver | ( | ) |
The Poisson destructor releases all internally managed data.
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::initializeSolverState | ( | tbox::Pointer< hier::PatchHierarchy< DIM > > | hierarchy, | |
int | ln = 0 | |||
) |
Initialize to a given hierarchy.
Initializer Poisson solver for a patch level in a hierarchy.
hierarchy | Hierarchy | |
ln | Level number |
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::deallocateSolverState | ( | ) |
Reset to an uninitialized state.
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::setMatrixCoefficients | ( | const PoissonSpecifications & | spec | ) |
Set the matrix coefficients.
For information describing the Poisson equation parameters, see the light-weight PoissonSpecifications class where you set the values of C and D.
This method must be called before solveSystem().
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::setSolnIdDepth | ( | const int | depth | ) | [inline] |
Set default depth of the solution data involved in the solve.
If the solution data has multiple depths, the solver uses just one depth at a time. The default depth is the first depth. Use this function to change it. This is not used to set the depth of the data (which is not controled by this class) but the depth used in the solve.
Changing the depth after setting up the matrix is permissible, as the solution data does not affect the matrix.
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::setRhsIdDepth | ( | const int | depth | ) | [inline] |
Set default depth of the rhs data involved in the solve.
If the rhs data has multiple depths, the solver uses just one depth at a time. The default depth is the first depth. Use this function to change it. This is not used to set the depth of the data (which is not controled by this class) but the depth used in the solve.
Changing the depth after setting up the matrix is permissible, as the rhs data does not affect the matrix.
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::setStoppingCriteria | ( | const int | max_iterations = 10 , |
|
const double | relative_residual_tol = 1.0e-6 | |||
) | [inline] |
Set the stopping criteria (max iterations and residual tolerance) for the linear solver.
max_iterations | gives the maximum number of iterations | |
relative_residual_tol | the maximum error tolerance |
int SAMRAI::solv::CellPoissonHypreSolver< DIM >::solveSystem | ( | const int | u, | |
const int | f, | |||
bool | homogeneous_bc = false | |||
) |
Solve the linear system Au=f.
The solution u and the right hand side f are specified via patch indices on the patch hierarchy.
Member functions getNumberOfIterations() return the iterations from the solver. Note that the matrix coefficients and boundary condition object must have been set up before this function is called. As long as the matrix coefficients do not change, this routine may be called repeatedly to solve any number of linear systems (with the right-hand side varying). If the boundary conditions or matrix coefficients are changed then function setMatrixCoefficients() must be called again.
When computing the matrix coefficients in setMatrixCoefficients(), the inhomogeneous portion of the boundary condition (constant terms, independent of u and thus having no effect on the matrix) are saved and added to the source term, f, before performing the matrix solve. In some situations, it may be useful to not add the inhomogeneous portion to f. The flag argument homoegneous_bc
is used for this. (This is a sort of optimization, to avoid having to re-call setMatrixCoefficients() to change the inhomogeneous portion.)
u | Descriptor of cell-centered unknown variable. | |
f | Descriptor of cell-centered source variable. | |
homogeneous_bc | Whether homogeneous boundary conditions are assumed. |
int SAMRAI::solv::CellPoissonHypreSolver< DIM >::getNumberOfIterations | ( | ) | const [inline] |
Return the number of iterations taken by the solver to converge.
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::setNumPreRelaxSteps | ( | const int | steps | ) | [inline] |
Set the number of pre-relax steps used by the Hypre solve.
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::setNumPostRelaxSteps | ( | const int | steps | ) | [inline] |
Set the number of post-relax steps used by the Hypre solve.
double SAMRAI::solv::CellPoissonHypreSolver< DIM >::getRelativeResidualNorm | ( | ) | const [inline] |
Return the final residual norm returned by the Hypre solve.
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::setUseSMG | ( | bool | use_smg | ) | [inline] |
Set whether to use Hypre's PFMG algorithm instead of the SMG algorithm.
The flag is used to select which of HYPRE's linear solver algorithms to use if true, the semicoarsening multigrid algorithm is used, and if false, the "PF" multigrid algorithm is used. By default, the SMG algorithm is used.
Changing the algorithm must be done before setting up the matrix coefficients.
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::setBoundaries | ( | const std::string & | boundary_type, | |
const int | fluxes = -1 , |
|||
const int | flags = -1 , |
|||
int * | bdry_types = NULL | |||
) | [inline] |
Specify boundary condition directly, without using a RobinBcCoefStrategy<DIM> object.
Use either setBoundaries() or setPhysicalBcCoefObject(), but not both.
A SimpleCelBcCoef object is used to interpret and implement the specified boundary conditions. See SimpleCellRobinBcCoefs<DIM>::setBoundaries() for an explanation of the arguments.
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::setPhysicalBcCoefObject | ( | const RobinBcCoefStrategy< DIM > * | physical_bc_coef_strategy, | |
const tbox::Pointer< hier::Variable< DIM > > | variable = NULL | |||
) | [inline] |
Specify boundary condition through the use of a Robin boundary condition object.
Use either setBoundaries() or setPhysicalBcCoefObject(), but not both.
The Robin boundary condition object is used when setting the matrix coefficient and when solving the system. If your boundary conditions are fixed values at ghost cell centers, use the GhostCellRobinBcCoefs<DIM> implementation of the RobinBcCoefStrategy<DIM> strategy.
physical_bc_coef_strategy | tbox::Pointer a concrete implementation of the Robin bc strategy. | |
variable | hier::Variable pointer to be passed to RobinBcCoefStrategy<DIM>::setBcCoefs(), but otherwise unused by this class. |
void SAMRAI::solv::CellPoissonHypreSolver< DIM >::setPrintSolverInfo | ( | const bool | ) | [inline] |
Set the flag for printing solver information.
This optional function is used primarily for debugging.
If set true, it will print the HYPRE matrix information to the following files:
If this method is not called, or the flag is set false, no printing will occur.