Chombo + EB
3.2
|
#include <MultilevelLinearOp.H>
Public Member Functions | |
MultilevelLinearOp () | |
default constructor (sets defaults for preconditioner) More... | |
virtual void | define (const Vector< DisjointBoxLayout > &a_vectGrids, const Vector< int > &a_refRatios, const Vector< ProblemDomain > &a_domains, const Vector< RealVect > &a_vectDx, RefCountedPtr< AMRLevelOpFactory< LevelData< T > > > &a_opFactory, int a_lBase) |
define function – opFactory should be able to define all required levels More... | |
virtual | ~MultilevelLinearOp () |
virtual void | residual (Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_phi, const Vector< LevelData< T > * > &a_rhs, bool a_homogeneous=false) |
compute residual using AMRLevelOps AMRResidual functionality More... | |
virtual void | preCond (Vector< LevelData< T > * > &a_cor, const Vector< LevelData< T > * > &a_residual) |
Apply preconditioner. More... | |
virtual void | applyOp (Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_phi, bool a_homogeneous=false) |
virtual void | create (Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_rhs) |
virtual void | clear (Vector< LevelData< T > * > &a_lhs) |
virtual void | assign (Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_rhs) |
virtual Real | dotProduct (const Vector< LevelData< T > * > &a_1, const Vector< LevelData< T > * > &a_2) |
virtual void | incr (Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_x, Real a_scale) |
virtual void | axby (Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_x, const Vector< LevelData< T > * > &a_y, Real a_a, Real a_b) |
virtual void | scale (Vector< LevelData< T > * > &a_lhs, const Real &a_scale) |
virtual Real | norm (const Vector< LevelData< T > * > &a_rhs, int a_ord) |
virtual void | setToZero (Vector< LevelData< T > * > &a_lhs) |
virtual void | write (const Vector< LevelData< T > * > *a_data, const char *a_filename) |
Public Member Functions inherited from LinearOp< Vector< LevelData< T > * > > | |
virtual | ~LinearOp () |
virtual void | assignLocal (Vector< LevelData< T > *> &a_lhs, const Vector< LevelData< T > *> &a_rhs) |
virtual void | mDotProduct (const Vector< LevelData< T > *> &a_1, const int a_sz, const Vector< LevelData< T > *> a_2[], Real a_mdots[]) |
virtual Real | dx () const |
Public Attributes | |
bool | m_use_multigrid_preconditioner |
if true, preCond(...) function uses multigrid v-cycle(s) More... | |
int | m_num_mg_iterations |
number of mg v-cycles for each preCond call (if using mg preconditioner) More... | |
int | m_num_mg_smooth |
number of smoothings in the multigrid preconditioner More... | |
int | m_preCondSolverDepth |
max depth for multigrid preconditioner (if using mg preconditioner) More... | |
Vector< DisjointBoxLayout > | m_vectGrids |
Vector< int > | m_refRatios |
Vector< ProblemDomain > | m_domains |
Vector< RealVect > | m_vectDx |
Vector< RefCountedPtr< AMRLevelOp< LevelData< T > > > > | m_vectOperators |
vector of level operators More... | |
int | m_lBase |
AMRMultiGrid< LevelData< T > > | m_preCondSolver |
Solver which is used by the preCond function. More... | |
LinearSolver< LevelData< T > > * | m_precondBottomSolverPtr |
bottom solver for m_preCondSolver More... | |
bool | m_isPrecondSolverInitialized |
has the preconditioner solver been initialized? More... | |
Operator class for Linear solver for solving L(phi) = rhs for the case where <T> is a multilevel (AMR) data structure; as such, T is actually a Vector<T>. The MultilevelLinearOp is also defined with a Vector<AMRLevelOp> which applies the operator in an AMR way. For now, this is hardwired for T = LevelData<FArrayBox>, but we can hopefully get back redo this once I get the kinks worked out...
MultilevelLinearOp< T >::MultilevelLinearOp | ( | ) |
default constructor (sets defaults for preconditioner)
|
virtual |
References CH_TIME.
|
virtual |
define function – opFactory should be able to define all required levels
define function
References CH_assert, CH_TIME, BiCGStabSolver< T >::m_verbosity, Max(), Vector< T >::resize(), and Vector< T >::size().
|
virtual |
compute residual using AMRLevelOps AMRResidual functionality
a_lhs = L(a_phi) - a_rhs. If a_homogeneous is true, evaluate the operator using homogeneous boundary conditions.
Say you are solving L(phi) = rhs. Make a_lhs = L(a_phi) - a_rhs. If a_homogeneous is true, evaluate the operator using homogeneous boundary conditions.
Implements LinearOp< Vector< LevelData< T > * > >.
|
virtual |
Apply preconditioner.
Given the current state of the residual the correction, apply your preconditioner to a_cor.
Implements LinearOp< Vector< LevelData< T > * > >.
|
virtual |
In the context of solving L(phi) = rhs, set a_lhs = L(a_phi). If a_homogeneous is true, evaluate the operator using homogeneous boundary conditions.
Implements LinearOp< Vector< LevelData< T > * > >.
References CH_TIME.
|
virtual |
Create data holder a_lhs that mirrors a_rhs. You do not need to copy the data of a_rhs, just make a holder the same size.
Implements LinearOp< Vector< LevelData< T > * > >.
|
virtual |
Clean up data holder before it goes out of scope. This is necessary because create calls new.
Reimplemented from LinearOp< Vector< LevelData< T > * > >.
References CH_TIME.
|
virtual |
Set a_lhs equal to a_rhs.
Implements LinearOp< Vector< LevelData< T > * > >.
|
virtual |
Compute and return the dot product of a_1 and a_2. In most contexts, this means return the sum over all data points of a_1*a_2.
Implements LinearOp< Vector< LevelData< T > * > >.
References LayoutIterator::begin(), CH_assert, CH_TIME, Box::coarsen(), BoxLayout::dataIterator(), getBoxes(), Box::isEmpty(), BoxLayout::layoutIterator(), BoxLayoutData< T >::nComp(), LayoutIterator::ok(), scale(), and Vector< T >::size().
|
virtual |
Increment by scaled amount (a_lhs += a_scale*a_x).
Implements LinearOp< Vector< LevelData< T > * > >.
References CH_TIME.
|
virtual |
Set input to a scaled sum (a_lhs = a_a*a_x + a_b*a_y).
Implements LinearOp< Vector< LevelData< T > * > >.
References CH_TIME.
|
virtual |
Multiply the input by a given scale (a_lhs *= a_scale).
Implements LinearOp< Vector< LevelData< T > * > >.
References CH_TIME.
|
virtual |
Return the norm of a_rhs. a_ord == 0 max norm, a_ord == 1 sum(abs(a_rhs)), else, L(a_ord) norm.
Implements LinearOp< Vector< LevelData< T > * > >.
References CH_TIME, and MayDay::Error().
|
virtual |
Set a_lhs to zero.
Implements LinearOp< Vector< LevelData< T > * > >.
|
virtual |
Debugging aid for solvers. Print out a "T" to a file named "filename" default implementation is to print out a message saying "LinearOp::write not implemented"
Reimplemented from LinearOp< Vector< LevelData< T > * > >.
References MayDay::Warning(), and writeVectorLevelName().
bool MultilevelLinearOp< T >::m_use_multigrid_preconditioner |
if true, preCond(...) function uses multigrid v-cycle(s)
int MultilevelLinearOp< T >::m_num_mg_iterations |
number of mg v-cycles for each preCond call (if using mg preconditioner)
int MultilevelLinearOp< T >::m_num_mg_smooth |
number of smoothings in the multigrid preconditioner
int MultilevelLinearOp< T >::m_preCondSolverDepth |
max depth for multigrid preconditioner (if using mg preconditioner)
Vector<DisjointBoxLayout> MultilevelLinearOp< T >::m_vectGrids |
Vector<int> MultilevelLinearOp< T >::m_refRatios |
Vector<ProblemDomain> MultilevelLinearOp< T >::m_domains |
Vector<RealVect> MultilevelLinearOp< T >::m_vectDx |
Vector<RefCountedPtr<AMRLevelOp<LevelData<T> > > > MultilevelLinearOp< T >::m_vectOperators |
vector of level operators
int MultilevelLinearOp< T >::m_lBase |
AMRMultiGrid<LevelData<T> > MultilevelLinearOp< T >::m_preCondSolver |
Solver which is used by the preCond function.
if m_use_multigrid_preconditioner is true, this solver is used to apply a AMR multigrid v-cycle as the preconditioner
LinearSolver<LevelData<T> >* MultilevelLinearOp< T >::m_precondBottomSolverPtr |
bottom solver for m_preCondSolver
bool MultilevelLinearOp< T >::m_isPrecondSolverInitialized |
has the preconditioner solver been initialized?