Main Page | Namespace List | Class Hierarchy | Class List | File List | Namespace Members | Class Members | File Members | Related Pages

SAMRAI::solv::CellPoissonFACSolver< DIM > Class Template Reference

Class for solving scalar Poisson's equation on SAMR grid, wrapping up lower-level components (FAC cycling, Poisson equation operations and boundary conditions) in a single high-level interface. More...

#include <source/solvers/poisson/CellPoissonFACSolver.h>

List of all members.

Public Member Functions

 CellPoissonFACSolver (const string &object_name, tbox::Pointer< tbox::Database > database=tbox::Pointer< tbox::Database >())
 Construct a solver.
 ~CellPoissonFACSolver (void)
 Destructor.
void enableLogging (bool logging)
 Enable logging.
bool solveSystem (const int solution, const int rhs, tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy, int coarse_ln=-1, int fine_ln=-1)
 Solve Poisson's equation, assuming an uninitialized solver state.
bool solveSystem (const int solution, const int rhs)
 Solve Poisson's equation using the current solver state set by initializeSolverState().
void setBoundaries (const string &boundary_type, const int fluxes=-1, const int flags=-1, int *bdry_types=NULL)
 Specify the boundary conditions that are to be used at the physical domain boundary.
void setBcObject (const RobinBcCoefStrategy< DIM > *bc_object)
 Override internal implementation to set boundary condition coefficients with user-provided implementation.
name Specifying PDE parameters
void 
setDPatchDataId (int id)
 Set the patch data index for variable D.
void setDConstant (double scalar)
 Set the scalar value variable D.
void setCPatchDataId (int id)
 Set the scalar value variable C.
void setCConstant (double scalar)
 Set the patch data index for variable C.
name Functions for setting
solver mathematic algorithm
controls void 
setCoarsestLevelSolverChoice (const string &choice)
 Set coarse level solver.
void setCoarsestLevelSolverTolerance (double tol)
 Set tolerance for coarse level solve.
void setCoarsestLevelSolverMaxIterations (int max_iterations)
 Set max iterations for coarse level solve.
void setUseSMG (bool use_smg)
 Set whether to use HYPRe's PFMG algorithm instead of the SMG algorithm.
void setCoarseFineDiscretization (const string &coarsefine_method)
 Set the coarse-fine boundary discretization method.
void setProlongationMethod (const string &prolongation_method)
 Set the name of the prolongation method.
void setPresmoothingSweeps (int num_pre_sweeps)
 Set the number of pre-smoothing sweeps during FAC iteration process.
void setPostsmoothingSweeps (int num_post_sweeps)
 Set the number of post-smoothing sweeps during FAC iteration process.
void setMaxCycles (int max_cycles)
 Set the max number of iterations (cycles) to use per solve.
void setResidualTolerance (double residual_tol)
 Set the residual tolerance for stopping.
void initializeSolverState (const int solution, const int rhs, tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy, const int coarse_level=-1, const int fine_level=-1)
 Prepare the solver's internal state for solving.
void deallocateSolverState ()
 Remove the solver's internal state data.
int getNumberIterations () const
 Return FAC iteration count from last (or current if there is one) FAC iteration process.
void getConvergenceFactors (double &avg_factor, double &final_factor) const
 Get average convergance rate and convergence rate of the last (or current if there is one) FAC solve.
double getResidualNorm () const
 Return residual norm from the just-completed FAC iteration.


Detailed Description

template<int DIM>
class SAMRAI::solv::CellPoissonFACSolver< DIM >

Class for solving scalar Poisson's equation on SAMR grid, wrapping up lower-level components (FAC cycling, Poisson equation operations and boundary conditions) in a single high-level interface.

Note: this class provides a backward-compatible interface to the soon-to-be obsolete PoissonHierarchySolver<DIM> class. Although this class hides the lower-level components (FAC cycling, Poisson equation operations and boundary conditions), it is perfectly acceptable to use those lower-level components directly.

We solve the equation div(D grad(u)) + Cu = f where D is a side-centered array and C is a cell-centered array. u and f are also cell-centered. Boundary conditions supported are Dirichlet, Neumann and mixed (Dirichlet on some faces and Neumann on others).

This class is a wrapper, providing a single class that coordinates three major components: the FAC solver, the cell-centered Poisson FAC operator and a default Robin bc coefficient implelemtation. It is perfectly acceptable to use those classes outside of this class.

The underlying solver is an FAC solver using cell-centered discretization. The difference scheme is second-order central-difference. On coarse-fine boundaries within the solution levels, the composite grid operator uses, by default, the discretization method of Ewing, Lazarov and Vassilevski ("Local Refinement Techniques for Elliptic Problems on Cell-Centered Grids, I. Error Analysis", Mathematics of Computation, Vol. 56, No. 194, April 1991, pp. 437-461).

Typical use of this class is:

  1. Construct a CellPoissonFACSolver<DIM> object, providing it the hierarchy and range of levels participating in the solve.
  2. Set the parameters C and D using the functions named setC... and setD... By default, D=1 and C=0 everywhere.
  3. Call setBoundaries() to state the types boundary conditions, along with supplemental data for setting those boundary conditions.
  4. Call initializeSolverState() to set up information internal to the solver. This is step is not required but will save setup costs if you are making multiple solves. This commits the object to the current hierarchy state and the specific types of boundary conditions you selected, It does NOT commit to the specific values of the boundary condition. A hierarchy change (through adaption or other means) invalidates the state, thus you must reinitialize or deallocateSolverState() the state before another solve.
  5. Solve the equation with solveSystem(). You provide the patch data indices for the solution u and the right hand side f. u must have at least one ghost cell and where a Dirichlet boundary condition applies, those cells must be set to the value on the boundary. If only Neumann boundary conditions are used, the ghost cell values do not matter.
  6. Call deallocateSolverState() to free up internal resources, if initializeSolverState() was called before the solve.

After the solve, information on the solve can be obtained by calling one of these functions:

Finer solver controls can be set using the functions in this class.

Object of this class can be set using input databases. The following parameters can be set. Each is shown with its default value in the case where hypre is used.

 * enable_logging = TRUE // Bool flag to switch logging on/off
 * max_cycles = 10       // Integer number of max FAC cycles to use
 * residual_tol = 1.e-6  // Residual tolerance to solve for
 * num_pre_sweeps = 1    // Number of presmoothing sweeps to use
 * num_post_sweeps = 1   // Number of postsmoothing sweeps to use
 * coarse_fine_discretization = "Ewing" // Name of coarse-fine discretization
 * prolongation_method = "CONSTANT_REFINE" // Name of prolongation method
 * coarse_solver_choice = "hypre"  // Name of coarse level solver
 * coarse_solver_tolerance = 1e-10 // Coarse level tolerance
 * coarse_solver_max_iterations = 20 // Coarse level max iterations
 * use_smg = "FALSE"     // Whether to use hypre's smg solver
 *                       // (alternative is the pfmg solver)
 * 


Constructor & Destructor Documentation

template<int DIM>
SAMRAI::solv::CellPoissonFACSolver< DIM >::CellPoissonFACSolver const string &  object_name,
tbox::Pointer< tbox::Database database = tbox::Pointertbox::Database >()
 

Construct a solver.

If the database is not NULL, initial settings will be set using the database. The solver is uninitialized until initializeSolverState() is called.

Parameters:
object_name Name of object used in outputs
database tbox::Database for initialization (may be NULL)

template<int DIM>
SAMRAI::solv::CellPoissonFACSolver< DIM >::~CellPoissonFACSolver void   ) 
 

Destructor.


Member Function Documentation

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::enableLogging bool  logging  ) 
 

Enable logging.

To disable, pass in false.

template<int DIM>
bool SAMRAI::solv::CellPoissonFACSolver< DIM >::solveSystem const int  solution,
const int  rhs,
tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy,
int  coarse_ln = -1,
int  fine_ln = -1
 

Solve Poisson's equation, assuming an uninitialized solver state.

Here, u is the "solution" patch data index and f is the right hand side patch data index. The return value is true if the solver converged and false otherwise.

This function is a wrapper. It simply initializes the solver state, call the solveSystem(const int,const int) for the initialized solver then deallocates the solver state.

Upon return from this function, solution will contain the result of the solve.

See initializeSolverState() for opportunities to save overhead when using multiple consecutive solves.

See also:
solveSystem(const int,const int)
Parameters:
solution hier::Patch data index for solution u
rhs hier::Patch data index for right hand side f
hierarchy The patch hierarchy to solve on
coarse_ln The coarsest level in the solve.
fine_ln The finest level in the solve.
Returns:
whether solver converged to specified level
See also:
initializeSolverState

template<int DIM>
bool SAMRAI::solv::CellPoissonFACSolver< DIM >::solveSystem const int  solution,
const int  rhs
 

Solve Poisson's equation using the current solver state set by initializeSolverState().

When the solver state has been initialized, this function may be called repeadedly with different values on the rhs. There is some cost savings for multiple solves when this is done.

Before calling this function, the solution and right-hand-side quantities should be set properly by the user on all patch interiors on the range of levels covered by the FAC iteration. All data for these patch data index should be allocated. Thus, the user is responsible for managing the storage for the solution and right-hand-side.

Returns:
whether solver converged to specified level
See also:
solveSystem( const int, const int, tbox::Pointer< hier::PatchHierarchy<DIM> >, int, int);

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setBoundaries const string &  boundary_type,
const int  fluxes = -1,
const int  flags = -1,
int *  bdry_types = NULL
 

Specify the boundary conditions that are to be used at the physical domain boundary.

This method is used to set up the default SimpleCellRobinBcCoefs<DIM> object for specifying boundary conditions. Note that you may alternatively provide your own implementation of the Robin boundary condition coefficients using the setBcObject() method.

The boundary conditions specified as the string argument "boundary_type." The boundary type argument can be "Dirichlet", "Neumann", or "Mixed".

If using Dirichlet boundary conditions, then before the solver is called, the storage for the unknown u must have a layer of ghost cells at least one cell wide that includes the Dirichlet boundary values.

If using Neumann boundary conditions, then before the solver is called, the outerface boundary flux data must be set for the Neumann conditions. The fluxes argument gives the patch data index of this flux data.

The mixed boundary type is for a mixture of Dirichlet and Neumann boundary conditions are used at the physical domain boundary. The fluxes argument gives the patch data index of the outerface data that specifies the flux data for the Neumann conditions. The flags array is an outerface data array of integer flags that specifies whether Dirichlet (flag == zero) or Neumann (flag == one) conditions are to be used at a particular cell boundary face. Note that the flag data must be set before the matrix entries can be computed and the flux data must be set before the solver is called. The bdry_types argument can be used if the boundary conditions are mixed but one or more of the faces of the physical boundary are entirely either Dirichlet or Neumann boundaries. The bdry_types argument should be an array of 2*DIM integers, specifying the boundary conditions on each side of the physical domain. It should be ordered {x_lo, x_hi, y_lo, y_hi, z_lo, z_hi}, with the values for each face being 0 for Dirichlet conditions, 1 for Neumann conditions, and 2 for mixed boundary conditions. The bdry_type argument is never required, but if used it can sometimes make the PoissonHYPRESolver class more efficient.

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setBcObject const RobinBcCoefStrategy< DIM > *  bc_object  ) 
 

Override internal implementation to set boundary condition coefficients with user-provided implementation.

This function is used to override the default internal object for setting Robin boundary condition coefficients. You should override when you need to avoid the limitations of the SimpleCellRobinBcCoefs<DIM> class or you prefer to use your own implementation.

Note that an important limitation of the SimpleCellRobinBcCoefs<DIM> class is the inability to support linear interpolation in the prolongation step.

Once the boundary condition object is overwritten by this method, you must no longer call the setBoundaries() method.

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setDPatchDataId int  id  )  [inline]
 

Set the patch data index for variable D.

In addition, disregard any previous D specified by setDConstant().

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setDConstant double  scalar  )  [inline]
 

Set the scalar value variable D.

In addition, disregard any previous D specified by setDPatchDataId().

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setCPatchDataId int  id  )  [inline]
 

Set the scalar value variable C.

In addition, disregard any previous C specified by setCConstant().

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setCConstant double  scalar  )  [inline]
 

Set the patch data index for variable C.

In addition, disregard any previous C specified by setCConstant().

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setCoarsestLevelSolverChoice const string &  choice  )  [inline]
 

Set coarse level solver.

Select from these:

  • "redblack"
  • "hypre" (only if the HYPRE library is available).

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setCoarsestLevelSolverTolerance double  tol  )  [inline]
 

Set tolerance for coarse level solve.

If the coarse level solver requires a tolerance (currently, they all do), the specified value is used.

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setCoarsestLevelSolverMaxIterations int  max_iterations  )  [inline]
 

Set max iterations for coarse level solve.

If the coarse level solver requires a max iteration limit (currently, they all do), the specified value is used.

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< 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.

This setting has effect only when HYPRe is chosen for the coarsest level solver. See setCoarsestLevelSolverChoice().

Changing the algorithm must be done before setting up the matrix coefficients.

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setCoarseFineDiscretization const string &  coarsefine_method  )  [inline]
 

Set the coarse-fine boundary discretization method.

Specify the op_name string which will be passed to xfer::Geometry<DIM>::lookupRefineOperator() to get the operator for setting fine grid ghost cells from the coarse grid. Note that chosing this operator implicitly choses the discretization method at the coarse-fine boundary.

There is one important instance where this string is not passed to xfer::Geometry<DIM>::lookupRefineOperator(). If this variable is set to "Ewing", a constant refinement method is used along with Ewing's correction. For a reference to the correction method, see "Local Refinement Techniques for Elliptic Problems on Cell-Centered Grids, I. Error Analysis", Mathematics of Computation, Vol. 56, No. 194, April 1991, pp. 437-461.

Parameters:
method String selecting the coarse-fine discretization method.

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setProlongationMethod const string &  prolongation_method  )  [inline]
 

Set the name of the prolongation method.

Specify the op_name string which will be passed to xfer::Geometry<DIM>::lookupRefineOperator() to get the operator for prolonging the coarse-grid correction.

By default, "CONSTANT_REFINE" is used. "LINEAR_REFINE" seems to to lead to faster convergence, but it does NOT satisfy the Galerkin condition.

Prolonging using linear refinement requires a Robin bc coefficient implementation that is capable of delivering coefficients for non-hierarchy data, because linear refinement requires boundary conditions to be set on temporary levels.

Parameters:
method String selecting the coarse-fine discretization method.

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setPresmoothingSweeps int  num_pre_sweeps  )  [inline]
 

Set the number of pre-smoothing sweeps during FAC iteration process.

Presmoothing is applied during the fine-to-coarse phase of the iteration. The default is to use one sweep.

Parameters:
num_pre_sweeps Number of presmoothing sweeps

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setPostsmoothingSweeps int  num_post_sweeps  )  [inline]
 

Set the number of post-smoothing sweeps during FAC iteration process.

Postsmoothing is applied during the coarse-to-fine phase of the iteration. The default is to use one sweep.

Parameters:
num_post_sweeps Number of postsmoothing sweeps

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setMaxCycles int  max_cycles  )  [inline]
 

Set the max number of iterations (cycles) to use per solve.

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::setResidualTolerance double  residual_tol  )  [inline]
 

Set the residual tolerance for stopping.

If you want the prescribed maximum number of cycles to always be taken, set the residual tolerance to a negative number.

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::initializeSolverState const int  solution,
const int  rhs,
tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy,
const int  coarse_level = -1,
const int  fine_level = -1
 

Prepare the solver's internal state for solving.

In the interest of efficiency, this class may prepare and cache some hierarchy-dependent objects. Though it is not required, initializing the solver state makes for greater efficiency when you are doing multiple solves on the same system of equation. If you do not initialize the state, it is initialized and deallocated each time you call solveSystem(const int, const int). The state must be reinitialized if the hierarchy or a boundary condition type changes.

To unset the data set in this function, see deallocateSolverState().

The solution and rhs patch data indices in the argument list are used to determine the form of the data you plan to use in the solve. They need not be the same data you solve on, but they should be similar. Both must represent cell-centered double data. The solution must have at least one ghost cell width, though this is not checked in the initialize phase, because data is not required yet.

Parameters:
solution solution patch data index for u
rhs right hand side patch data index for f
hierarchy The patch hierarchy to solve on
coarse_ln The coarsest level in the solve
fine_ln The finest level in the solve

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::deallocateSolverState  ) 
 

Remove the solver's internal state data.

Remove all hierarchy-dependent data set by initializeSolverState. It is safe to call deallocateSolverState() even state is already deallocated, but nothing is done in that case.

See also:
initializeSolverState()

template<int DIM>
int SAMRAI::solv::CellPoissonFACSolver< DIM >::getNumberIterations  )  const [inline]
 

Return FAC iteration count from last (or current if there is one) FAC iteration process.

template<int DIM>
void SAMRAI::solv::CellPoissonFACSolver< DIM >::getConvergenceFactors double &  avg_factor,
double &  final_factor
const [inline]
 

Get average convergance rate and convergence rate of the last (or current if there is one) FAC solve.

Parameters:
avg_factor average convergence factor over current FAC cycles
final_factor convergence factor of the last FAC cycle

template<int DIM>
double SAMRAI::solv::CellPoissonFACSolver< DIM >::getResidualNorm  )  const [inline]
 

Return residual norm from the just-completed FAC iteration.

The norm return value is computed as the maximum norm over all patch levels involved in the solve. The value corresponds to the norm applied in the user-defined residual computation.

The latest computed norm is the one returned.


The documentation for this class was generated from the following files:
Generated on Fri Dec 2 11:32:42 2005 for SAMRAI by  doxygen 1.4.2