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

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 <source/solvers/poisson/CellPoissonFACOps.h>

Inheritance diagram for SAMRAI::solv::CellPoissonFACOps< DIM >:

Inheritance graph
[legend]
List of all members.

Public Member Functions

 CellPoissonFACOps (const string &object_name=string(), tbox::Pointer< tbox::Database > database=NULL)
 Constructor.
 ~CellPoissonFACOps (void)
 Destructor.
void setPoissonSpecifications (const PoissonSpecifications &spec)
 Set the scalar Poisson equation specifications.
void enableLogging (bool enable_logging)
 Enable logging.
void setSmoothingChoice (const string &smoothing_choice)
 Set the choice of smoothing algorithms.
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 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 setUseSMG (bool use_smg)
 Set whether to use Hypre's PFMG algorithm instead of the SMG algorithm.
void setFluxId (int flux_id)
 Set the scratch patch data index for the flux.
void setPhysicalBcCoefObject (const RobinBcCoefStrategy< DIM > *physical_bc_coef)
 Provide an implementation for getting the physical bc coefficients.
void checkInputPatchDataIndices () const
 Check validity and correctness of input patch data indices.
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.
void setPreconditioner (const FACPreconditioner< DIM > *preconditioner)
 Set the FAC preconditioner that will be using this object.
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.
name FACOperatorStrategy<
DIM > virtual virtuals void 
restrictSolution (const SAMRAIVectorReal< DIM, double > &s, SAMRAIVectorReal< DIM, double > &d, int dest_ln)
 Restrict the solution quantity to the specified level from the next finer level.
virtual void restrictResidual (const SAMRAIVectorReal< DIM, double > &s, SAMRAIVectorReal< DIM, double > &d, int dest_ln)
 Restrict the residual quantity to the specified level from the next finer level.
virtual void prolongErrorAndCorrect (const SAMRAIVectorReal< DIM, double > &s, SAMRAIVectorReal< DIM, double > &d, 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.
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.
virtual void solveCoarsestLevel (SAMRAIVectorReal< DIM, double > &error, const SAMRAIVectorReal< DIM, double > &residual, int ln)
 Solve the residual equation Ae=r on the coarsest level in the FAC iteration.
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.
virtual double computeResidualNorm (const SAMRAIVectorReal< DIM, double > &residual, int fine_ln, int coarse_ln)
 Compute the norm of the residual quantity.
virtual void initializeOperatorState (const SAMRAIVectorReal< DIM, double > &solution, const SAMRAIVectorReal< DIM, double > &rhs)
 Compute hierarchy-dependent data if any is required.
virtual void deallocateOperatorState ()
 Remove all hierarchy-dependent data.
virtual void postprocessOneCycle (int iteration_num, const SAMRAIVectorReal< DIM, double > &current_soln, const SAMRAIVectorReal< DIM, double > &residual)
 Regular call back routine to be called after each FAC cycle.

Detailed Description

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

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

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

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

Constructor.

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

Parameters:
object_name Ojbect name
database Input database

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

Destructor.

Deallocate internal data.


Member Function Documentation

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setPoissonSpecifications const PoissonSpecifications spec  )  [inline]
 

Set the scalar Poisson equation specifications.

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::enableLogging bool  enable_logging  )  [inline]
 

Enable logging.

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

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setSmoothingChoice const string &  smoothing_choice  )  [inline]
 

Set the choice of smoothing algorithms.

Current smoothing choices are:

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

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

Set coarse level solver.

Select from these:

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

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< 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::CellPoissonFACOps< 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::CellPoissonFACOps< 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", 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:
method String selecting the coarse-fine discretization method.

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< 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:
prolongation_method String selecting the coarse-fine discretization method.

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setUseSMG bool  use_smg  )  [inline]
 

Set whether to use Hypre's PFMG algorithm instead of the SMG algorithm.

This flag affects the Hypre solver (used to solve the coarsest level). 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 initializing the solver state and must NOT be done while the state is initialized (the program will exit), as that would corrupt the state.

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::setFluxId int  flux_id  )  [inline]
 

Set the scratch patch data index for the flux.

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.

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

Provide an implementation for getting the physical bc coefficients.

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_coef tbox::Pointer to an object that can set the Robin bc coefficients.

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

Check validity and correctness of input patch data indices.

Descriptors checked:

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

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
 

Set weight appropriate for computing vector norms.

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:
hierarchy Hierarchy configuration to compute weights for
weight_id hier::Patch data index of the weight
coarsest_ln Coarsest level number. Must be included in hierarchy. Must not be greater than finest_ln. Default to 0.
finest_ln Finest level number. Must be included in hierarchy. Must not be less than coarsest_ln. Default to finest level in hierarchy.

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

Set the FAC preconditioner that will be using this object.

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

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
 

function to compute flux, using general diffusion coefficient data.

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:
patch patch on which computation will take place
ratio_to_coarser_level refinement ratio from coarser level to level on which patch lives; if current patch level is level zero, this is ignored
w_data cell-centered data
Dgradw_data side-centered flux data (i.e., D (grad w))

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::restrictSolution const SAMRAIVectorReal< DIM, double > &  s,
SAMRAIVectorReal< DIM, double > &  d,
int  dest_ln
[virtual]
 

Restrict the solution quantity to the specified level from the next finer level.

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:
source source solution
dest destination solution
dest_ln destination level number

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

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::restrictResidual const SAMRAIVectorReal< DIM, double > &  s,
SAMRAIVectorReal< DIM, double > &  d,
int  dest_ln
[virtual]
 

Restrict the residual quantity to the specified level from the next finer level.

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:
source source residual
dest destination residual
dest_ln destination level number

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

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::prolongErrorAndCorrect const SAMRAIVectorReal< DIM, double > &  s,
SAMRAIVectorReal< DIM, double > &  d,
int  dest_ln
[virtual]
 

Prolong the error quantity to the specified level from the next coarser level and apply the correction to the fine-level error.

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:
source source error vector
dest destination error vector
dest_ln destination level number of data transfer

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

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

Perform a given number of relaxations on the error.

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:
error error vector
residual residual vector
ln level number
num_sweeps number of sweeps

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

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::solveCoarsestLevel SAMRAIVectorReal< DIM, double > &  error,
const SAMRAIVectorReal< DIM, double > &  residual,
int  ln
[virtual]
 

Solve the residual equation Ae=r on the coarsest level in the FAC iteration.

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:
error error vector
residual residual vector
coarsest_ln coarsest level number

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

template<int DIM>
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]
 

Compute composite grid residual on a single level.

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:
residual residual vector
solution solution vector
rhs source (right hand side) vector
ln level number
error_equation_indicator flag stating whether u is an error vector or a solution vector

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

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

Compute the norm of the residual quantity.

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:
residual residual vector
fine_ln finest level number
coarse_ln coarsest level number
Returns:
norm value of residual vector, which should be non-negative

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

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

Compute hierarchy-dependent data if any is required.

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:
solution solution vector u
rhs right hand side vector f
The default implementation does nothing.

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

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

Remove all hierarchy-dependent data.

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

See also:
initializeOperatorState

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

template<int DIM>
void SAMRAI::solv::CellPoissonFACOps< DIM >::postprocessOneCycle int  iteration_num,
const SAMRAIVectorReal< DIM, double > &  current_soln,
const SAMRAIVectorReal< DIM, double > &  residual
[virtual]
 

Regular call back routine to be called after each FAC cycle.

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_num FAC cycle number completed
current_soln current solution
residual residual based on the current solution

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


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