#include <source/solvers/FAC/FACOperatorStrategy.h>
Inheritance diagram for SAMRAI::solv::FACOperatorStrategy< DIM >:
Public Member Functions | |
FACOperatorStrategy () | |
Empty constructor. | |
virtual | ~FACOperatorStrategy () |
Virtual destructor. | |
Operator-dependent virtual methods | |
virtual void | restrictSolution (const SAMRAIVectorReal< DIM, double > &source, SAMRAIVectorReal< DIM, double > &dest, int dest_ln)=0 |
Restrict the solution quantity to the specified level from the next finer level. | |
virtual void | restrictResidual (const SAMRAIVectorReal< DIM, double > &source, SAMRAIVectorReal< DIM, double > &dest, int dest_ln)=0 |
Restrict the residual quantity to the specified level from the next finer level. | |
virtual void | prolongErrorAndCorrect (const SAMRAIVectorReal< DIM, double > &source, SAMRAIVectorReal< DIM, double > &dest, int dest_ln)=0 |
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)=0 |
Perform a given number of relaxations on the error. | |
virtual int | solveCoarsestLevel (SAMRAIVectorReal< DIM, double > &error, const SAMRAIVectorReal< DIM, double > &residual, int coarsest_ln)=0 |
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)=0 |
Compute composite grid residual on a single level. | |
virtual double | computeResidualNorm (const SAMRAIVectorReal< DIM, double > &residual, int fine_ln, int coarse_ln)=0 |
Compute the norm of the residual quantity. | |
virtual void | postprocessOneCycle (int fac_cycle_num, const SAMRAIVectorReal< DIM, double > ¤t_soln, const SAMRAIVectorReal< DIM, double > &residual) |
Regular call back routine to be called after each FAC cycle. | |
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. |
The FACPreconditioner<DIM> constructor accepts a concrete implementation of this interface and calls the concrete implementations of the virtual functions declared herein during the solution process.
All vector arguments in these interfaces are guaranteed to be either the vectors given to in FACPreconditioner<DIM>::solveSystem() or FACPreconditioner<DIM>::initializeSolverState() or vectors cloned from them.
SAMRAI::solv::FACOperatorStrategy< DIM >::FACOperatorStrategy | ( | ) |
Empty constructor.
SAMRAI::solv::FACOperatorStrategy< DIM >::~FACOperatorStrategy | ( | ) | [virtual] |
Virtual destructor.
virtual void SAMRAI::solv::FACOperatorStrategy< DIM >::restrictSolution | ( | const SAMRAIVectorReal< DIM, double > & | source, | |
SAMRAIVectorReal< DIM, double > & | dest, | |||
int | dest_ln | |||
) | [pure 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:
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.
source | source solution | |
dest | destination solution | |
dest_ln | destination level number |
Implemented in SAMRAI::solv::CellPoissonFACOps< DIM >.
virtual void SAMRAI::solv::FACOperatorStrategy< DIM >::restrictResidual | ( | const SAMRAIVectorReal< DIM, double > & | source, | |
SAMRAIVectorReal< DIM, double > & | dest, | |||
int | dest_ln | |||
) | [pure 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:
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.
source | source residual | |
dest | destination residual | |
dest_ln | destination level number |
Implemented in SAMRAI::solv::CellPoissonFACOps< DIM >.
virtual void SAMRAI::solv::FACOperatorStrategy< DIM >::prolongErrorAndCorrect | ( | const SAMRAIVectorReal< DIM, double > & | source, | |
SAMRAIVectorReal< DIM, double > & | dest, | |||
int | dest_ln | |||
) | [pure 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
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.
source | source error vector | |
dest | destination error vector | |
dest_ln | destination level number of data transfer |
Implemented in SAMRAI::solv::CellPoissonFACOps< DIM >.
virtual void SAMRAI::solv::FACOperatorStrategy< DIM >::smoothError | ( | SAMRAIVectorReal< DIM, double > & | error, | |
const SAMRAIVectorReal< DIM, double > & | residual, | |||
int | ln, | |||
int | num_sweeps | |||
) | [pure 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:
Steps for each iteration.
Final step before leaving function.
error | error vector | |
residual | residual vector | |
ln | level number | |
num_sweeps | number of sweeps |
Implemented in SAMRAI::solv::CellPoissonFACOps< DIM >.
virtual int SAMRAI::solv::FACOperatorStrategy< DIM >::solveCoarsestLevel | ( | SAMRAIVectorReal< DIM, double > & | error, | |
const SAMRAIVectorReal< DIM, double > & | residual, | |||
int | coarsest_ln | |||
) | [pure 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.
error | error vector | |
residual | residual vector | |
coarsest_ln | coarsest level number |
Implemented in SAMRAI::solv::CellPoissonFACOps< DIM >.
virtual void SAMRAI::solv::FACOperatorStrategy< DIM >::computeCompositeResidualOnLevel | ( | SAMRAIVectorReal< DIM, double > & | residual, | |
const SAMRAIVectorReal< DIM, double > & | solution, | |||
const SAMRAIVectorReal< DIM, double > & | rhs, | |||
int | ln, | |||
bool | error_equation_indicator | |||
) | [pure 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:
Steps:
Final step before leaving 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 .
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 |
Implemented in SAMRAI::solv::CellPoissonFACOps< DIM >.
virtual double SAMRAI::solv::FACOperatorStrategy< DIM >::computeResidualNorm | ( | const SAMRAIVectorReal< DIM, double > & | residual, | |
int | fine_ln, | |||
int | coarse_ln | |||
) | [pure 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.
residual | residual vector | |
fine_ln | finest level number | |
coarse_ln | coarsest level number |
Implemented in SAMRAI::solv::CellPoissonFACOps< DIM >.
void SAMRAI::solv::FACOperatorStrategy< DIM >::postprocessOneCycle | ( | int | fac_cycle_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.
fac_cycle_num | FAC cycle number completed | |
current_soln | current solution | |
residual | residual based on the current solution |
Reimplemented in SAMRAI::solv::CellPoissonFACOps< DIM >.
void SAMRAI::solv::FACOperatorStrategy< 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:
solution | solution vector u | |
rhs | right hand side vector f |
Reimplemented in SAMRAI::solv::CellPoissonFACOps< DIM >.
void SAMRAI::solv::FACOperatorStrategy< DIM >::deallocateOperatorState | ( | ) | [virtual] |
Remove all hierarchy-dependent data.
Remove all hierarchy-dependent data set by initializeOperatorState().
Reimplemented in SAMRAI::solv::CellPoissonFACOps< DIM >.