IBAMR  IBAMR version 0.19.
Public Member Functions | List of all members
SAMRAI::solv::CellPoissonFACOps< DIM > Class Template Reference

FAC operator class to solve Poisson's equation on a SAMR grid, using cell-centered, second-order finite-volume method, with Robin boundary conditions. More...

#include <CellPoissonFACOps.h>

Inheritance diagram for SAMRAI::solv::CellPoissonFACOps< DIM >:
Inheritance graph
[legend]

Public Member Functions

 CellPoissonFACOps (const std::string &object_name=std::string(), tbox::Pointer< tbox::Database > database=NULL)
 Constructor. More...
 
 ~CellPoissonFACOps ()
 Destructor. More...
 
void setPoissonSpecifications (const PoissonSpecifications &spec)
 Set the scalar Poisson equation specifications. More...
 
void enableLogging (bool enable_logging)
 Enable logging. More...
 
Functions for setting solver mathematic algorithm controls
void setSmoothingChoice (const std::string &smoothing_choice)
 Set the choice of smoothing algorithms. More...
 
void setCoarsestLevelSolverChoice (const std::string &choice)
 Set coarse level solver. More...
 
void setCoarsestLevelSolverTolerance (double tol)
 Set tolerance for coarse level solve. More...
 
void setCoarsestLevelSolverMaxIterations (int max_iterations)
 Set max iterations for coarse level solve. More...
 
void setCoarseFineDiscretization (const std::string &coarsefine_method)
 Set the coarse-fine boundary discretization method. More...
 
void setProlongationMethod (const std::string &prolongation_method)
 Set the name of the prolongation method. More...
 
Functions for setting patch data indices and coefficients
void setFluxId (int flux_id)
 Set the scratch patch data index for the flux. More...
 
void setPhysicalBcCoefObject (const RobinBcCoefStrategy< DIM > *physical_bc_coef)
 Provide an implementation for getting the physical bc coefficients. More...
 
Functions for checking validity and correctness of state.
void checkInputPatchDataIndices () const
 Check validity and correctness of input patch data indices. More...
 
void computeVectorWeights (tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy, int weight_id, int coarsest_ln=-1, int finest_ln=-1) const
 Set weight appropriate for computing vector norms. More...
 
void setPreconditioner (const FACPreconditioner< DIM > *preconditioner)
 Set the FAC preconditioner that will be using this object. More...
 
void computeFluxOnPatch (const hier::Patch< DIM > &patch, const hier::IntVector< DIM > &ratio_to_coarser_level, const pdat::CellData< DIM, double > &w_data, pdat::SideData< DIM, double > &Dgradw_data) const
 function to compute flux, using general diffusion coefficient data. More...
 
virtual void restrictSolution (const SAMRAIVectorReal< DIM, double > &source, SAMRAIVectorReal< DIM, double > &dest, int dest_ln)
 Restrict the solution quantity to the specified level from the next finer level. More...
 
virtual void restrictResidual (const SAMRAIVectorReal< DIM, double > &source, SAMRAIVectorReal< DIM, double > &dest, int dest_ln)
 Restrict the residual quantity to the specified level from the next finer level. More...
 
virtual void prolongErrorAndCorrect (const SAMRAIVectorReal< DIM, double > &source, SAMRAIVectorReal< DIM, double > &dest, int dest_ln)
 Prolong the error quantity to the specified level from the next coarser level and apply the correction to the fine-level error. More...
 
virtual void smoothError (SAMRAIVectorReal< DIM, double > &error, const SAMRAIVectorReal< DIM, double > &residual, int ln, int num_sweeps)
 Perform a given number of relaxations on the error. More...
 
virtual int solveCoarsestLevel (SAMRAIVectorReal< DIM, double > &error, const SAMRAIVectorReal< DIM, double > &residual, int coarsest_ln)
 Solve the residual equation Ae=r on the coarsest level in the FAC iteration. More...
 
virtual void computeCompositeResidualOnLevel (SAMRAIVectorReal< DIM, double > &residual, const SAMRAIVectorReal< DIM, double > &solution, const SAMRAIVectorReal< DIM, double > &rhs, int ln, bool error_equation_indicator)
 Compute composite grid residual on a single level. More...
 
virtual double computeResidualNorm (const SAMRAIVectorReal< DIM, double > &residual, int fine_ln, int coarse_ln)
 Compute the norm of the residual quantity. More...
 
virtual void initializeOperatorState (const SAMRAIVectorReal< DIM, double > &solution, const SAMRAIVectorReal< DIM, double > &rhs)
 Compute hierarchy-dependent data if any is required. More...
 
virtual void deallocateOperatorState ()
 Remove all hierarchy-dependent data. More...
 
virtual void postprocessOneCycle (int fac_cycle_num, const SAMRAIVectorReal< DIM, double > &current_soln, const SAMRAIVectorReal< DIM, double > &residual)
 Regular call back routine to be called after each FAC cycle. More...
 

Private Attributes

Various refine and coarsen objects used internally.
tbox::Pointer< xfer::RefineOperator< DIM > > d_prolongation_refine_operator
 Error prolongation (refinement) operator. More...
 
tbox::Pointer< xfer::RefineAlgorithm< DIM > > d_prolongation_refine_algorithm
 
tbox::Array< tbox::Pointer< xfer::RefineSchedule< DIM > > > d_prolongation_refine_schedules
 
tbox::Pointer< xfer::CoarsenOperator< DIM > > d_urestriction_coarsen_operator
 Solution restriction (coarsening) operator. More...
 
tbox::Pointer< xfer::CoarsenAlgorithm< DIM > > d_urestriction_coarsen_algorithm
 
tbox::Array< tbox::Pointer< xfer::CoarsenSchedule< DIM > > > d_urestriction_coarsen_schedules
 
tbox::Pointer< xfer::CoarsenOperator< DIM > > d_rrestriction_coarsen_operator
 Residual restriction (coarsening) operator. More...
 
tbox::Pointer< xfer::CoarsenAlgorithm< DIM > > d_rrestriction_coarsen_algorithm
 
tbox::Array< tbox::Pointer< xfer::CoarsenSchedule< DIM > > > d_rrestriction_coarsen_schedules
 
tbox::Pointer< xfer::CoarsenOperator< DIM > > d_flux_coarsen_operator
 Coarsen operator for outerflux-to-flux. More...
 
tbox::Pointer< xfer::CoarsenAlgorithm< DIM > > d_flux_coarsen_algorithm
 
tbox::Array< tbox::Pointer< xfer::CoarsenSchedule< DIM > > > d_flux_coarsen_schedules
 
tbox::Pointer< xfer::RefineOperator< DIM > > d_ghostfill_refine_operator
 Refine operator for cell-like data from coarser level. More...
 
tbox::Pointer< xfer::RefineAlgorithm< DIM > > d_ghostfill_refine_algorithm
 
tbox::Array< tbox::Pointer< xfer::RefineSchedule< DIM > > > d_ghostfill_refine_schedules
 
tbox::Pointer< xfer::RefineOperator< DIM > > d_ghostfill_nocoarse_refine_operator
 Refine operator for cell-like data from same level. More...
 
tbox::Pointer< xfer::RefineAlgorithm< DIM > > d_ghostfill_nocoarse_refine_algorithm
 
tbox::Array< tbox::Pointer< xfer::RefineSchedule< DIM > > > d_ghostfill_nocoarse_refine_schedules
 
CartesianRobinBcHelper< DIM > d_bc_helper
 Utility object employed in setting ghost cells and providing xfer::RefinePatchStrategy<DIM> implementation. More...
 
Non-essential objects used in outputs and debugging.
bool d_enable_logging
 Logging flag. More...
 
const FACPreconditioner< DIM > * d_preconditioner
 Preconditioner using this object. More...
 
tbox::Pointer< math::HierarchyCellDataOpsReal< DIM, double > > d_hopscell
 Hierarchy cell operator used in debugging. More...
 
tbox::Pointer< math::HierarchySideDataOpsReal< DIM, double > > d_hopsside
 Hierarchy side operator used in debugging. More...
 
tbox::Pointer< tbox::Timert_restrict_solution
 Timers for performance measurement. More...
 
tbox::Pointer< tbox::Timert_restrict_residual
 
tbox::Pointer< tbox::Timert_prolong
 
tbox::Pointer< tbox::Timert_smooth_error
 
tbox::Pointer< tbox::Timert_solve_coarsest
 
tbox::Pointer< tbox::Timert_compute_composite_residual
 
tbox::Pointer< tbox::Timert_compute_residual_norm
 

Private workhorse functions.

std::string d_object_name
 Object name. More...
 
tbox::Pointer< hier::PatchHierarchy< DIM > > d_hierarchy
 Reference hierarchy. More...
 
int d_ln_min
 Coarsest level for solve. More...
 
int d_ln_max
 Finest level for solve. More...
 
tbox::Array< hier::CoarseFineBoundary< DIM > > d_cf_boundary
 Description of coarse-fine boundaries. More...
 
void smoothErrorByRedBlack (SAMRAIVectorReal< DIM, double > &error, const SAMRAIVectorReal< DIM, double > &residual, int ln, int num_sweeps, double residual_tolerance=-1.0)
 Red-black Gauss-Seidel error smoothing on a level. More...
 
int solveCoarsestLevel_HYPRE (SAMRAIVectorReal< DIM, double > &error, const SAMRAIVectorReal< DIM, double > &residual, int ln)
 Solve the coarsest level using HYPRE. More...
 
void ewingFixFlux (const hier::Patch< DIM > &patch, const pdat::CellData< DIM, double > &soln_data, pdat::SideData< DIM, double > &flux_data, const hier::IntVector< DIM > &ratio_to_coarser) const
 Fix flux per Ewing's coarse-fine boundary treatment. More...
 
void computeResidualOnPatch (const hier::Patch< DIM > &patch, const pdat::SideData< DIM, double > &flux_data, const pdat::CellData< DIM, double > &soln_data, const pdat::CellData< DIM, double > &rhs_data, pdat::CellData< DIM, double > &residual_data) const
 AMR-unaware function to compute residual on a single patch, with variable scalar field. More...
 
void redOrBlackSmoothingOnPatch (const hier::Patch< DIM > &patch, const pdat::SideData< DIM, double > &flux_data, const pdat::CellData< DIM, double > &rhs_data, pdat::CellData< DIM, double > &soln_data, char red_or_black, double *p_maxres=NULL) const
 AMR-unaware function to red or black smoothing on a single patch, for variable diffusion coefficient and variable scalar field. More...
 
void xeqScheduleProlongation (int dst_id, int src_id, int scr_id, int dest_ln)
 Execute a refinement schedule for prolonging cell data. More...
 
void xeqScheduleURestriction (int dst_id, int src_id, int dest_ln)
 Execute schedule for restricting solution to the specified level or reregister an existing one. More...
 
void xeqScheduleRRestriction (int dst_id, int src_id, int dest_ln)
 Execute schedule for restricting residual to the specified level or reregister an existing one. More...
 
void xeqScheduleFluxCoarsen (int dst_id, int src_id, int dest_ln)
 Execute schedule for coarsening flux to the specified level or reregister an existing one. More...
 
void xeqScheduleGhostFill (int dst_id, int dest_ln)
 Execute schedule for filling ghosts on the specified level or reregister an existing one. More...
 
void xeqScheduleGhostFillNoCoarse (int dst_id, int dest_ln)
 Execute schedule for filling ghosts on the specified level or reregister an existing one. This version does not get data from coarser levels. More...
 
int registerCellScratch () const
 Return the patch data index for cell scratch data. More...
 
int registerFluxScratch () const
 Return the patch data index for flux scratch data. More...
 
int registerOfluxScratch () const
 Return the patch data index for outerflux scratch data. More...
 
static void freeVariables ()
 Free static variables at shutdown time. More...
 

Private state variables for solution process.

PoissonSpecifications d_poisson_spec
 Scalar Poisson equations specifications. More...
 
std::string d_smoothing_choice
 Smoothing choice. More...
 
std::string d_coarse_solver_choice
 Coarse level solver. More...
 
std::string d_cf_discretization
 Coarse-fine discretization method. More...
 
std::string d_prolongation_method
 Coarse-fine discretization method. More...
 
double d_coarse_solver_tolerance
 Tolerance specified to coarse solver. More...
 
int d_coarse_solver_max_iterations
 Coarse level solver iteration limit. More...
 
double d_residual_tolerance_during_smoothing
 Residual tolerance to govern smoothing. More...
 
int d_flux_id
 Id of the flux. More...
 
const RobinBcCoefStrategy< DIM > * d_physical_bc_coef
 Externally provided physical boundary condition object. More...
 
tbox::Pointer< hier::VariableContextd_context
 Default context of internally maintained hierarchy data. More...
 
int d_cell_scratch_id
 ID of the solution-like scratch data. More...
 
int d_flux_scratch_id
 ID of the side-centered scratch data. More...
 
int d_oflux_scratch_id
 ID of the outerside-centered scratch data. More...
 
static tbox::Pointer< pdat::CellVariable< DIM, double > > s_cell_scratch_var
 
static tbox::Pointer< pdat::SideVariable< DIM, double > > s_flux_scratch_var
 
static tbox::Pointer< pdat::OutersideVariable< DIM, double > > s_oflux_scratch_var
 

Detailed Description

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

This class provides operators that are used by the FAC preconditioner FACPreconditioner<DIM>. It is used to solve the scalar Poisson's equation using a cell-centered second-order finite-volume discretization. It is designed to provide all operations specific to the scalar Poisson's equation,

\[ \nabla \cdot D \nabla u + C u = f \]

(see PoissonSpecifications) where

You are left to provide the source function, initial guess, etc., by specifying them in specific forms.

This class provides:

  1. 5-point (second order), cell-centered stencil operations for the discrete Laplacian.
  2. Red-black Gauss-Seidel smoothing.
  3. Provisions for working Robin boundary conditions (see RobinBcCoefStrategy).

This class is meant to provide the Poisson-specific operator used by the FAC preconditioner, FACPreconditioner<DIM>. To use the preconditioner with this class, you will have to provide:

  1. The solution vector SAMRAIVectorReal<DIM>, with appropriate norm weighting for the cell-centered AMR mesh. This class provides the function computeVectorWeights() to help with computing the appropriate weights. Since this is for a scalar equation, only the first depth of the first component of the vectors are used. All other parts are ignored.
  2. The source vector SAMRAIVectorReal<DIM> for f.
  3. A PoissonSpecifications objects to specify the cell-centered scalar field C and the side-centered diffusion coefficients D
  4. The boundary condition specifications in terms of the coefficients \( \alpha \), \( \beta \) and \( \gamma \) in the Robin formula \( \alpha u + \beta u_n = \gamma \) applied on the boundary faces. See RobinBcCoefStrategy<DIM>.

This class allocates and deallocates only its own scratch data. Other data that it manipuates are passed in as function arguments. Hence, it owns none of the solution vectors, error vectors, diffusion coefficient data, or any such things.

Input Examples

* coarse_solver_choice = "hypre"    // see setCoarsestLevelSolverChoice()
* coarse_solver_tolerance = 1e-14   // see setCoarsestLevelSolverTolerance()
* coarse_solver_max_iterations = 10 // see setCoarsestLevelSolverMaxIterations()
* smoothing_choice = "redblack"     // see setSmoothingChoice()
* cf_discretization = "Ewing"       // see setCoarseFineDiscretization()
* prolongation_method = "LINEAR_REFINE" // see setProlongationMethod()
* hypre_solver = { ... }            // tbox::Database for initializing Hypre solver
* 

Constructor & Destructor Documentation

◆ CellPoissonFACOps()

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

If you want standard output and logging, pass in valid pointers for those streams.

Parameters
object_nameOjbect name
databaseInput database

◆ ~CellPoissonFACOps()

template<int DIM>
SAMRAI::solv::CellPoissonFACOps< DIM >::~CellPoissonFACOps ( )

Deallocate internal data.

Member Function Documentation

◆ setPoissonSpecifications()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setPoissonSpecifications ( const PoissonSpecifications spec)

◆ enableLogging()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::enableLogging ( bool  enable_logging)

By default, logging is disabled. The logging flag is propagated to the major components used by this class.

◆ setSmoothingChoice()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setSmoothingChoice ( const std::string &  smoothing_choice)

Current smoothing choices are:

  • "redblack": Red-black Gauss-Seidel smoothing.

◆ setCoarsestLevelSolverChoice()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setCoarsestLevelSolverChoice ( const std::string &  choice)

Select from these:

  • "redblack" (red-black smoothing until convergence–very slow!)
  • "hypre" (only if the HYPRE library is available).

◆ setCoarsestLevelSolverTolerance()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setCoarsestLevelSolverTolerance ( double  tol)

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

◆ setCoarsestLevelSolverMaxIterations()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setCoarsestLevelSolverMaxIterations ( int  max_iterations)

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

◆ setCoarseFineDiscretization()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setCoarseFineDiscretization ( const std::string &  coarsefine_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", Ewing's coarse-fine discretization is used (a constant refinement is performed, and the flux is later corrected to result in Ewing's scheme). For a reference to Ewing's discretization 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
coarsefine_methodString selecting the coarse-fine discretization method.

◆ setProlongationMethod()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setProlongationMethod ( const std::string &  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
prolongation_methodString selecting the coarse-fine discretization method.

◆ setFluxId()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setFluxId ( int  flux_id)

The use of this function is optional. The patch data index should be a pdat::SideData<DIM> type of variable. If the flux id is -1 (the default initial value), scratch space for the flux is allocated as needed and immediately deallocated afterward, level by level. If you have space preallocated for flux and you would like that to be used, set flux id to the patch data index of that space.

◆ setPhysicalBcCoefObject()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setPhysicalBcCoefObject ( const RobinBcCoefStrategy< DIM > *  physical_bc_coef)

If your solution is fixed at the physical boundary ghost cell centers AND those cells have the correct values before entering solveSystem(), you may use a GhostCellRobinBcCoefs<DIM> object.

If your solution is not fixed at the ghost cell centers, the ghost cell values will change as the interior cell values change. In those cases, the flexible Robin boundary conditions are applied. You must call this function to provide the implementation for determining the boundary condition coefficients.

Parameters
physical_bc_coeftbox::Pointer to an object that can set the Robin bc coefficients.

◆ checkInputPatchDataIndices()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::checkInputPatchDataIndices ( ) const

Descriptors checked:

  1. Diffusion coefficient (see setDiffcoefId())
  2. Flux (see setFluxId())
  3. Source (see setScalarFieldId())

◆ computeVectorWeights()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::computeVectorWeights ( tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy,
int  weight_id,
int  coarsest_ln = -1,
int  finest_ln = -1 
) const

If you this function to set the weights used when you SAMRAIVectorReal<DIM>::addComponent, you can use the vector norm functions of SAMRAIVectorReal<DIM>, and the weights will be used to blank out coarse grid regions under fine grids.

The weights computed are specific to the cell-centered discretization used by this class. The weight is equal to the cell volume if the cell has not been refined, and zero if it has.

This function is state-independent. All inputs are in the argument list.

Parameters
hierarchyHierarchy configuration to compute weights for
weight_idhier::Patch data index of the weight
coarsest_lnCoarsest level number. Must be included in hierarchy. Must not be greater than finest_ln. Default to 0.
finest_lnFinest level number. Must be included in hierarchy. Must not be less than coarsest_ln. Default to finest level in hierarchy.

◆ setPreconditioner()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setPreconditioner ( const FACPreconditioner< DIM > *  preconditioner)

The FAC preconditioner is accessed to get convergence data during the cycle postprocessing step. It is optional.

◆ computeFluxOnPatch()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::computeFluxOnPatch ( const hier::Patch< DIM > &  patch,
const hier::IntVector< DIM > &  ratio_to_coarser_level,
const pdat::CellData< DIM, double > &  w_data,
pdat::SideData< DIM, double > &  Dgradw_data 
) const

Recall that this solver class discretizes the PDE

\[ \nabla \cdot D \nabla u + C u = f \]

on an AMR grid. This member function allows users of this solver class to compute gradient terms,

\[ D \nabla w \]

, in their code in a manner consistent with the solver discretization. In particular, when solving PDE systems, it may be necessary to discretize the gradient operator appearing in equations not treated by the solver class in the same way as those treated by this class. These funtions allow users to do this easily. The divergence operator used in this solver is the standard sum of centered differences involving flux terms on the cell sides computed by these routines.

Note that the patch must exist on a level in an AMR hierarchy so that the discretization can be computed properly at the coarse-fine interface. Poisson coefficients C and D must exist on the patch, if they are variable. Also, calling this function does not affect the internal solver state in any way. However, the solver must be fully initialized before it is called and care should be exercised to pass arguments so that the solver solution quantity and other internal solver quantities are not adversely affected.

Parameters
patchpatch on which computation will take place
ratio_to_coarser_levelrefinement ratio from coarser level to level on which patch lives; if current patch level is level zero, this is ignored
w_datacell-centered data
Dgradw_dataside-centered flux data (i.e., D (grad w))

◆ restrictSolution()

template<int DIM>
virtual void SAMRAI::solv::CellPoissonFACOps< DIM >::restrictSolution ( const SAMRAIVectorReal< DIM, double > &  source,
SAMRAIVectorReal< DIM, double > &  dest,
int  dest_ln 
)
virtual

Restrict the residual data to level dest_ln in the destination vector d, from level dest_ln+1 in the source vector s.

Can assume:

  1. dest_ln is not the finest level in the range being solved.
  2. corresponding solution has been computed on level dest_ln+1.
  3. the source and destination residual vectors (s and d) may or may not be the same. (This function must work in either case.)

Upon return from this function, the solution on the refined region of the coarse level will represent the coarsened version of the fine solution in a manner that is consistent with the linear system approximation on the composite grid. This function must not change the solution values anywhere except on level dest_ln of the destination vector.

The source and destination vectors may be the same.

Parameters
sourcesource solution
destdestination solution
dest_lndestination level number

Implements SAMRAI::solv::FACOperatorStrategy< DIM >.

◆ restrictResidual()

template<int DIM>
virtual void SAMRAI::solv::CellPoissonFACOps< DIM >::restrictResidual ( const SAMRAIVectorReal< DIM, double > &  source,
SAMRAIVectorReal< DIM, double > &  dest,
int  dest_ln 
)
virtual

Restrict the residual data to level dest_ln in the destination vector d, from level dest_ln+1 in the source vector s.

Can assume:

  1. dest_ln is not the finest level in the range being solved.
  2. correspnding residual has been computed on level dest_ln+1.
  3. the source and destination residual vectors (s and d) may or may not be the same. (This function must work in either case.)

Upon return from this function, the residual on the refined region of the coarse level will represent the coarsened version of the fine residual in a manner that is consistent with the linear system approximation on the composite grid. This function must not change the residual values anywhere except on level dest_ln of the destination vector.

The source and destination vectors may be the same.

Parameters
sourcesource residual
destdestination residual
dest_lndestination level number

Implements SAMRAI::solv::FACOperatorStrategy< DIM >.

◆ prolongErrorAndCorrect()

template<int DIM>
virtual void SAMRAI::solv::CellPoissonFACOps< DIM >::prolongErrorAndCorrect ( const SAMRAIVectorReal< DIM, double > &  source,
SAMRAIVectorReal< DIM, double > &  dest,
int  dest_ln 
)
virtual

On the part of the coarse level that does not overlap the fine level, the error is the corection to Au=f.

On the part of the coarse level that does overlap the fine level, the error is the corection to Ae=r of the fine level.

This function should apply the coarse-level correction to the fine level, that is

\[ e^{fine} \leftarrow e^{fine} + I^{fine}_{coarse} e^{coarse} \]

Note: You probably have to store the refined error in a temporary location before adding it to the current error.

The array of boundary information contains a description of the coarse-fine level boundary for each patch on the level; the boundary information for patch N is obtained as the N-th element in the array, coarse_fine_boundary[N].

Upon return from this function, the error on the fine level must represent the correction to the solution on that level. Also, this function must not change the error values on the coarse level.

The source and destination vectors may be the same.

Parameters
sourcesource error vector
destdestination error vector
dest_lndestination level number of data transfer

Implements SAMRAI::solv::FACOperatorStrategy< DIM >.

◆ smoothError()

template<int DIM>
virtual void SAMRAI::solv::CellPoissonFACOps< DIM >::smoothError ( SAMRAIVectorReal< DIM, double > &  error,
const SAMRAIVectorReal< DIM, double > &  residual,
int  ln,
int  num_sweeps 
)
virtual

Relax the residual equation Ae=r by applying the given number of smoothing sweeps on the specified level. The relaxation may ignore the possible existence of finer levels on a given level.

The array of boundary information contains a description of the coarse-fine level boundary for each patch on the level; the boundary information for patch N is obtained as the N-th element in the array, coarse_fine_boundary[N].

May assume:

  • If intermediate data from level l+1 is needed (for example, to match flux at coarse-fine boundaries), that data is already computed and stored on level l+1.
  • The error in the next finer level has been computed and stored there.

Steps for each iteration.

  1. Fill ghost boundaries
  2. Compute intermediate data (if needed) and coarsen intermediate data stored in level l+1 (if needed).
  3. Perform relaxation step (update e toward a better approximation).

Final step before leaving function.

  • If needed, compute and store intermediate data for next coarser level l-1.
Parameters
errorerror vector
residualresidual vector
lnlevel number
num_sweepsnumber of sweeps

Implements SAMRAI::solv::FACOperatorStrategy< DIM >.

◆ solveCoarsestLevel()

template<int DIM>
virtual int SAMRAI::solv::CellPoissonFACOps< DIM >::solveCoarsestLevel ( SAMRAIVectorReal< DIM, double > &  error,
const SAMRAIVectorReal< DIM, double > &  residual,
int  coarsest_ln 
)
virtual

Here e is the given error quantity and r is the given residual quantity. The array of boundary information contains a description of the coarse-fine level boundary for each patch on the level; the boundary information for patch N is obtained as the N-th element in the array, coarse_fine_boundary[N].

This routine must fill boundary values for given solution quantity on all patches on the specified level before the solve is performed.

Parameters
errorerror vector
residualresidual vector
coarsest_lncoarsest level number
Returns
0 if solver converged to specified level, nonzero otherwise.

Implements SAMRAI::solv::FACOperatorStrategy< DIM >.

◆ computeCompositeResidualOnLevel()

template<int DIM>
virtual void SAMRAI::solv::CellPoissonFACOps< DIM >::computeCompositeResidualOnLevel ( SAMRAIVectorReal< DIM, double > &  residual,
const SAMRAIVectorReal< DIM, double > &  solution,
const SAMRAIVectorReal< DIM, double > &  rhs,
int  ln,
bool  error_equation_indicator 
)
virtual

For the specified level number ln, compute the composite residual r=f-Au, where f is the right hand side and u is the solution. Note that the composite residual is not a one-level residual. It must take into account the composite grid stencil around the coarse-fine grid interface.

May assume:

  • Composite residual on next finer level l+1, has been computed already.
  • If any intermediately computed data is needed from level l+1, it has been done and stored on that level.
  • Residual computations for the original equation and the error equations will not be intermingled within one FAC cycle.

Steps:

  1. Fill boundary ghosts.
  2. If needed, coarsen intermediate data from level l+1.
  3. Compute residual \( r^l \leftarrow f - A u^l \).

Final step before leaving function:

  • If any intermediately computed data is needed in at level l-1, it must be computed and stored before leaving this function.

Important: Do not restrict residual from finer levels. (However, you must write the function restrictResidual() to do this.)

Important: This function must also work when the right-hand-side and the residual are identical. In that case, it should effectively do \( r \leftarrow r - A u \).

Parameters
residualresidual vector
solutionsolution vector
rhssource (right hand side) vector
lnlevel number
error_equation_indicatorflag stating whether u is an error vector or a solution vector

Implements SAMRAI::solv::FACOperatorStrategy< DIM >.

◆ computeResidualNorm()

template<int DIM>
virtual double SAMRAI::solv::CellPoissonFACOps< DIM >::computeResidualNorm ( const SAMRAIVectorReal< DIM, double > &  residual,
int  fine_ln,
int  coarse_ln 
)
virtual

Compute norm of the given residual on the given range of hierarchy levels. The residual vector is computed already and you should not change it. The only purpose of this function to allow you to choose how to define the norm.

The norm value is used during the FAC iteration to determine convergence of the composite grid linear system.

Residual values that lie under a finer level should not be counted.

Parameters
residualresidual vector
fine_lnfinest level number
coarse_lncoarsest level number
Returns
norm value of residual vector, which should be non-negative

Implements SAMRAI::solv::FACOperatorStrategy< DIM >.

◆ initializeOperatorState()

template<int DIM>
virtual void SAMRAI::solv::CellPoissonFACOps< DIM >::initializeOperatorState ( const SAMRAIVectorReal< DIM, double > &  solution,
const SAMRAIVectorReal< DIM, double > &  rhs 
)
virtual

This function is called when the hierarchy configuration changes. If you maintain any hierarchy-dependent data in your implementation (for example, caching communication schedules or computing coarse-fine boundaries), use this function to update that data.

If you do not maintain such data, this function may be empty.

Note that although the vector arguments given to other methods in this class may not necessarily be the same as those given to this method, there will be similarities, including:

  • hierarchy configuration (hierarchy pointer and level range)
  • number, type and alignment of vector component data
  • ghost cell width of data in the solution (or solution-like) vector
Parameters
solutionsolution vector u
rhsright hand side vector f

The default implementation does nothing.

Reimplemented from SAMRAI::solv::FACOperatorStrategy< DIM >.

◆ deallocateOperatorState()

template<int DIM>
virtual void SAMRAI::solv::CellPoissonFACOps< DIM >::deallocateOperatorState ( )
virtual

Remove all hierarchy-dependent data set by initializeOperatorState().

See also
initializeOperatorState

Reimplemented from SAMRAI::solv::FACOperatorStrategy< DIM >.

◆ postprocessOneCycle()

template<int DIM>
virtual void SAMRAI::solv::CellPoissonFACOps< DIM >::postprocessOneCycle ( int  fac_cycle_num,
const SAMRAIVectorReal< DIM, double > &  current_soln,
const SAMRAIVectorReal< DIM, double > &  residual 
)
virtual

This function is called after each FAC cycle. It allows you to monitor the progress and do other things. You should not modify the solution vector in the argument.

The default implementation does nothing.

Parameters
fac_cycle_numFAC cycle number completed
current_solncurrent solution
residualresidual based on the current solution

Reimplemented from SAMRAI::solv::FACOperatorStrategy< DIM >.

◆ smoothErrorByRedBlack()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::smoothErrorByRedBlack ( SAMRAIVectorReal< DIM, double > &  error,
const SAMRAIVectorReal< DIM, double > &  residual,
int  ln,
int  num_sweeps,
double  residual_tolerance = -1.0 
)
private

Smoothes on the residual equation \( Ae=r \) on a level.

Parameters
errorerror vector
residualresidual vector
lnlevel number
num_sweepsnumber of sweeps
residual_tolerancethe maximum residual considered to be converged

◆ solveCoarsestLevel_HYPRE()

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::solveCoarsestLevel_HYPRE ( SAMRAIVectorReal< DIM, double > &  error,
const SAMRAIVectorReal< DIM, double > &  residual,
int  ln 
)
private

◆ ewingFixFlux()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::ewingFixFlux ( const hier::Patch< DIM > &  patch,
const pdat::CellData< DIM, double > &  soln_data,
pdat::SideData< DIM, double > &  flux_data,
const hier::IntVector< DIM > &  ratio_to_coarser 
) const
private

Ewing's coarse-fine boundary treatment can be implemented using a constant refinement into the fine-grid ghost boundary, naively computing the flux using the constant-refined data then fixing up the flux to correct the error.

To use this function

  1. you must use constant refinement to fill the fine level ghost cells
  2. the flux must first be computed and stored
Parameters
patchpatch
soln_datacell-centered solution data
flux_dataside-centered flux data
diffcoef_dataside-centered diffusion coefficient data
cfbcoarse-fine boundary object for the level in which patch resides
ratio_to_coarserRefinement ratio to the next coarser level.

◆ computeResidualOnPatch()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::computeResidualOnPatch ( const hier::Patch< DIM > &  patch,
const pdat::SideData< DIM, double > &  flux_data,
const pdat::CellData< DIM, double > &  soln_data,
const pdat::CellData< DIM, double > &  rhs_data,
pdat::CellData< DIM, double > &  residual_data 
) const
private
Parameters
patchpatch
flux_dataside-centered flux data
soln_datacell-centered solution data
rhs_datacell-centered rhs data
residual_datacell-centered residual data

◆ redOrBlackSmoothingOnPatch()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::redOrBlackSmoothingOnPatch ( const hier::Patch< DIM > &  patch,
const pdat::SideData< DIM, double > &  flux_data,
const pdat::CellData< DIM, double > &  rhs_data,
pdat::CellData< DIM, double > &  soln_data,
char  red_or_black,
double p_maxres = NULL 
) const
private
Parameters
patchpatch
flux_dataside-centered flux data
rhs_datacell-centered rhs data
scalar_field_datacell-centered scalar field data
soln_datacell-centered solution data
red_or_blackred-black switch. Set to 'r' or 'b'.
p_maxresmax residual output. Set to NULL to avoid computing.

◆ xeqScheduleProlongation()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::xeqScheduleProlongation ( int  dst_id,
int  src_id,
int  scr_id,
int  dest_ln 
)
private

General notes regarding internal objects for communication: We maintain objects to support caching schedules to improve efficiency. Communication is needed in 5 distinct tasks.

  1. Prolongation
  2. Restriction
  3. Flux coarsening. Changing the coarse grid flux to the composite grid flux by coarsening the fine grid flux at the coarse-fine boundaries.
  4. Fill boundary data from other patches in the same level and physical boundary condition.
  5. Fill boundary data from same level, coarser levels and physical boundary condition.

For each task, we maintain a refine or coarsen operator, and a array of communication schedules (one for each destination level).

The 5 member functions named xeqSchedule... execute communication schedules appropriate for five specific tasks. They use a cached schedule if possible or create and cache a new schedule if needed. These functions and the data they manipulate are as follows:

  1. xeqScheduleProlongation(): d_prolongation_refine_operator d_prolongation_refine_schedules
  2. xeqScheduleURestriction(): d_restriction_coarsen_operator, d_urestriction_coarsen_schedules.
  3. xeqScheduleRRestriction(): d_restriction_coarsen_operator, d_rrestriction_coarsen_schedules.
  4. xeqScheduleFluxCoarsen(): d_flux_coarsen_operator, d_flux_coarsen_schedules.
  5. xeqScheduleGhostFill(): d_ghostfill_refine_operator, d_ghostfill_refine_schedules.
  6. xeqScheduleGhostFillNoCoarse(): d_ghostfill_nocoarse_refine_operator, d_ghostfill_nocoarse_refine_schedules.
Returns
refinement schedule for prolongation

◆ xeqScheduleURestriction()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::xeqScheduleURestriction ( int  dst_id,
int  src_id,
int  dest_ln 
)
private

See general notes for xeqScheduleProlongation().

Returns
coarsening schedule for restriction

◆ xeqScheduleRRestriction()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::xeqScheduleRRestriction ( int  dst_id,
int  src_id,
int  dest_ln 
)
private

See general notes for xeqScheduleProlongation().

Returns
coarsening schedule for restriction

◆ xeqScheduleFluxCoarsen()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::xeqScheduleFluxCoarsen ( int  dst_id,
int  src_id,
int  dest_ln 
)
private

See general notes for xeqScheduleProlongation().

Returns
coarsening schedule for setting composite grid flux at coarse-fine boundaries.

◆ xeqScheduleGhostFill()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::xeqScheduleGhostFill ( int  dst_id,
int  dest_ln 
)
private

See general notes for xeqScheduleProlongation().

Returns
refine schedule for filling ghost data from coarser level and physical bc.

◆ xeqScheduleGhostFillNoCoarse()

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::xeqScheduleGhostFillNoCoarse ( int  dst_id,
int  dest_ln 
)
private

See general notes for xeqScheduleProlongation().

This function is used for the bottom solve level, since it does not access data from any coarser level. (Ghost data obtained from coarser level must have been placed there before solve begins!)

Returns
refine schedule for filling ghost data from same level and physical bc.

◆ registerCellScratch()

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::registerCellScratch ( ) const
private

◆ registerFluxScratch()

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::registerFluxScratch ( ) const
private

◆ registerOfluxScratch()

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::registerOfluxScratch ( ) const
private

◆ freeVariables()

template<int DIM>
static void SAMRAI::solv::CellPoissonFACOps< DIM >::freeVariables ( )
staticprivate

Member Data Documentation

◆ d_object_name

template<int DIM>
std::string SAMRAI::solv::CellPoissonFACOps< DIM >::d_object_name
private

◆ d_hierarchy

template<int DIM>
tbox::Pointer< hier::PatchHierarchy<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_hierarchy
private

This variable is non-null between the initializeOperatorState() and deallocateOperatorState() calls. It is not truly needed, because the hierarchy is obtainable through variables in most function argument lists. We use it to enforce working on one hierarchy at a time.

◆ d_ln_min

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::d_ln_min
private

◆ d_ln_max

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::d_ln_max
private

◆ d_cf_boundary

template<int DIM>
tbox::Array< hier::CoarseFineBoundary<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_cf_boundary
private

There is one coarse-fine boundary object for each level. d_coarse_fine_boundary[i] is the description of the coarse-fine boundary between level i and level i-1. The coarse-fine boundary does not exist at the coarsest level, although the hier::CoarseFineBoundary<DIM> object still exists (it should not contain any boxes).

This array is initialized in initializeOperatorState() and deallocated in deallocateOperatorState(). When allocated, it is allocated for the index range [0,d_ln_max], though the range [0,d_ln_min-1] is not used. This is okay because hier::CoarseFineBoundary<DIM> is a light object before it is set for a level.

◆ d_poisson_spec

template<int DIM>
PoissonSpecifications SAMRAI::solv::CellPoissonFACOps< DIM >::d_poisson_spec
private

◆ d_smoothing_choice

template<int DIM>
std::string SAMRAI::solv::CellPoissonFACOps< DIM >::d_smoothing_choice
private

◆ d_coarse_solver_choice

template<int DIM>
std::string SAMRAI::solv::CellPoissonFACOps< DIM >::d_coarse_solver_choice
private

◆ d_cf_discretization

template<int DIM>
std::string SAMRAI::solv::CellPoissonFACOps< DIM >::d_cf_discretization
private

◆ d_prolongation_method

template<int DIM>
std::string SAMRAI::solv::CellPoissonFACOps< DIM >::d_prolongation_method
private

The name of the refinement operator used to prolong the coarse grid correction.

See also
setProlongationMethod()

◆ d_coarse_solver_tolerance

template<int DIM>
double SAMRAI::solv::CellPoissonFACOps< DIM >::d_coarse_solver_tolerance
private

◆ d_coarse_solver_max_iterations

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::d_coarse_solver_max_iterations
private

◆ d_residual_tolerance_during_smoothing

template<int DIM>
double SAMRAI::solv::CellPoissonFACOps< DIM >::d_residual_tolerance_during_smoothing
private

When we use one of the internal error smoothing functions and want to terminate the smoothing sweeps at a certain level of residual, this will be set to > 0. If it is < 0, the smoothing function effectively ignores it.

This variable is needed because some coarse-level solver simply runs the smoothing function until convergence. It sets this variable to > 0, calls the smoothing function, then resets it to < 0.

◆ d_flux_id

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::d_flux_id
private

If set to -1, create and delete storage space on the fly. Else, user has provided space for flux.

See also
setFluxId

◆ d_physical_bc_coef

template<int DIM>
const RobinBcCoefStrategy<DIM>* SAMRAI::solv::CellPoissonFACOps< DIM >::d_physical_bc_coef
private

◆ s_cell_scratch_var

template<int DIM>
tbox::Pointer<pdat::CellVariable<DIM,double> > SAMRAI::solv::CellPoissonFACOps< DIM >::s_cell_scratch_var
staticprivate

◆ s_flux_scratch_var

template<int DIM>
tbox::Pointer<pdat::SideVariable<DIM,double> > SAMRAI::solv::CellPoissonFACOps< DIM >::s_flux_scratch_var
staticprivate

◆ s_oflux_scratch_var

template<int DIM>
tbox::Pointer<pdat::OutersideVariable<DIM,double> > SAMRAI::solv::CellPoissonFACOps< DIM >::s_oflux_scratch_var
staticprivate

◆ d_context

template<int DIM>
tbox::Pointer<hier::VariableContext> SAMRAI::solv::CellPoissonFACOps< DIM >::d_context
private

◆ d_cell_scratch_id

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::d_cell_scratch_id
private

Set in constructor and never changed. Corresponds to a pdat::CellVariable<DIM,double> named d_object_name+"::cell_scratch". Scratch data is allocated and removed as needed to reduce memory usage.

◆ d_flux_scratch_id

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::d_flux_scratch_id
private

Set in constructor and never changed. Corresponds to a pdat::SideVariable<DIM,double> named d_object_name+"::flux_scratch".

This data is allocated only as needed and deallocated immediately after use.

◆ d_oflux_scratch_id

template<int DIM>
int SAMRAI::solv::CellPoissonFACOps< DIM >::d_oflux_scratch_id
private

Set in constructor and never changed. Corresponds to a pdat::OutersideVariable<DIM,double> named d_object_name+"::oflux_scratch".

◆ d_prolongation_refine_operator

template<int DIM>
tbox::Pointer<xfer::RefineOperator<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_prolongation_refine_operator
private

◆ d_prolongation_refine_algorithm

template<int DIM>
tbox::Pointer<xfer::RefineAlgorithm<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_prolongation_refine_algorithm
private

◆ d_prolongation_refine_schedules

template<int DIM>
tbox::Array<tbox::Pointer<xfer::RefineSchedule<DIM> > > SAMRAI::solv::CellPoissonFACOps< DIM >::d_prolongation_refine_schedules
private

◆ d_urestriction_coarsen_operator

template<int DIM>
tbox::Pointer<xfer::CoarsenOperator<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_urestriction_coarsen_operator
private

◆ d_urestriction_coarsen_algorithm

template<int DIM>
tbox::Pointer<xfer::CoarsenAlgorithm<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_urestriction_coarsen_algorithm
private

◆ d_urestriction_coarsen_schedules

template<int DIM>
tbox::Array<tbox::Pointer<xfer::CoarsenSchedule<DIM> > > SAMRAI::solv::CellPoissonFACOps< DIM >::d_urestriction_coarsen_schedules
private

◆ d_rrestriction_coarsen_operator

template<int DIM>
tbox::Pointer<xfer::CoarsenOperator<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_rrestriction_coarsen_operator
private

◆ d_rrestriction_coarsen_algorithm

template<int DIM>
tbox::Pointer<xfer::CoarsenAlgorithm<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_rrestriction_coarsen_algorithm
private

◆ d_rrestriction_coarsen_schedules

template<int DIM>
tbox::Array<tbox::Pointer<xfer::CoarsenSchedule<DIM> > > SAMRAI::solv::CellPoissonFACOps< DIM >::d_rrestriction_coarsen_schedules
private

◆ d_flux_coarsen_operator

template<int DIM>
tbox::Pointer<xfer::CoarsenOperator<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_flux_coarsen_operator
private

◆ d_flux_coarsen_algorithm

template<int DIM>
tbox::Pointer<xfer::CoarsenAlgorithm<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_flux_coarsen_algorithm
private

◆ d_flux_coarsen_schedules

template<int DIM>
tbox::Array<tbox::Pointer<xfer::CoarsenSchedule<DIM> > > SAMRAI::solv::CellPoissonFACOps< DIM >::d_flux_coarsen_schedules
private

◆ d_ghostfill_refine_operator

template<int DIM>
tbox::Pointer<xfer::RefineOperator<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_ghostfill_refine_operator
private

◆ d_ghostfill_refine_algorithm

template<int DIM>
tbox::Pointer<xfer::RefineAlgorithm<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_ghostfill_refine_algorithm
private

◆ d_ghostfill_refine_schedules

template<int DIM>
tbox::Array<tbox::Pointer<xfer::RefineSchedule<DIM> > > SAMRAI::solv::CellPoissonFACOps< DIM >::d_ghostfill_refine_schedules
private

◆ d_ghostfill_nocoarse_refine_operator

template<int DIM>
tbox::Pointer<xfer::RefineOperator<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_ghostfill_nocoarse_refine_operator
private

◆ d_ghostfill_nocoarse_refine_algorithm

template<int DIM>
tbox::Pointer<xfer::RefineAlgorithm<DIM> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_ghostfill_nocoarse_refine_algorithm
private

◆ d_ghostfill_nocoarse_refine_schedules

template<int DIM>
tbox::Array<tbox::Pointer<xfer::RefineSchedule<DIM> > > SAMRAI::solv::CellPoissonFACOps< DIM >::d_ghostfill_nocoarse_refine_schedules
private

◆ d_bc_helper

template<int DIM>
CartesianRobinBcHelper<DIM> SAMRAI::solv::CellPoissonFACOps< DIM >::d_bc_helper
private

Since this class deals only in scalar variables having Robin boundary conditions, we take advantage of the corresponding implementation in CartesianRobinBcHelper<DIM>. Whenever we need an implementation of xfer::RefinePatchStrategy<DIM>, this object is used. Note that in the code, before we use this object to set ghost cell values, directly or indirectly by calling xfer::RefineSchedule::fillData(), we must tell d_bc_helper the patch data index we want to set and whether we are setting data with homogeneous boundary condition.

◆ d_enable_logging

template<int DIM>
bool SAMRAI::solv::CellPoissonFACOps< DIM >::d_enable_logging
private

◆ d_preconditioner

template<int DIM>
const FACPreconditioner<DIM>* SAMRAI::solv::CellPoissonFACOps< DIM >::d_preconditioner
private

◆ d_hopscell

template<int DIM>
tbox::Pointer<math::HierarchyCellDataOpsReal<DIM,double> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_hopscell
private

◆ d_hopsside

template<int DIM>
tbox::Pointer<math::HierarchySideDataOpsReal<DIM,double> > SAMRAI::solv::CellPoissonFACOps< DIM >::d_hopsside
private

◆ t_restrict_solution

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::solv::CellPoissonFACOps< DIM >::t_restrict_solution
private

◆ t_restrict_residual

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::solv::CellPoissonFACOps< DIM >::t_restrict_residual
private

◆ t_prolong

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::solv::CellPoissonFACOps< DIM >::t_prolong
private

◆ t_smooth_error

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::solv::CellPoissonFACOps< DIM >::t_smooth_error
private

◆ t_solve_coarsest

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::solv::CellPoissonFACOps< DIM >::t_solve_coarsest
private

◆ t_compute_composite_residual

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::solv::CellPoissonFACOps< DIM >::t_compute_composite_residual
private

◆ t_compute_residual_norm

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::solv::CellPoissonFACOps< DIM >::t_compute_residual_norm
private

The documentation for this class was generated from the following file: