11 #ifndef _MULTILEVELLINEAROPI_H_ 12 #define _MULTILEVELLINEAROPI_H_ 16 #include "NamespaceHeader.H" 24 m_use_multigrid_preconditioner =
true;
25 m_num_mg_iterations = 1;
28 m_preCondSolverDepth = -1;
29 m_precondBottomSolverPtr = NULL;
42 CH_TIME(
"MultilevelLinearOp::define");
44 int numLevels = a_vectGrids.
size();
55 m_vectGrids = a_vectGrids;
56 m_refRatios = a_refRatios;
57 m_domains = a_domains;
59 m_vectOperators.
resize(numLevels);
62 int coarsestLevel =
Max(0, m_lBase-1);
63 for (
int level=coarsestLevel; level<numLevels; level++)
66 m_vectOperators[level] = op;
70 if (m_use_multigrid_preconditioner)
75 m_precondBottomSolverPtr = newPtr;
76 int bottomSolverVerbosity = 1;
80 if (m_preCondSolverDepth >= 0)
82 m_preCondSolver.m_maxDepth = m_preCondSolverDepth;
84 m_preCondSolver.define(m_domains[0],
86 m_precondBottomSolverPtr,
91 Real solverEps = 1.0e-7;
92 int maxIterations = 1;
94 Real normThresh = 1.0e-30;
97 m_preCondSolver.setSolverParameters(m_num_mg_smooth,
107 m_preCondSolver.m_verbosity = 1;
111 m_isPrecondSolverInitialized =
false;
120 CH_TIME(
"MultilevelLinearOp::~MultilevelLinearOp");
122 if (m_precondBottomSolverPtr != NULL)
124 delete m_precondBottomSolverPtr;
125 m_precondBottomSolverPtr = NULL;
142 CH_TIME(
"MultilevelLinearOp::residual");
147 int numLevels = a_lhs.size();
150 CH_assert (m_vectOperators.size() >= numLevels);
152 int maxLevel = numLevels -1;
153 for (
int level = m_lBase; level < numLevels; level++)
163 if (level < maxLevel)
166 m_vectOperators[level]->AMRResidual(thisResid,
172 (m_vectOperators[level+1].operator->()));
177 m_vectOperators[level]->AMRResidualNF(thisResid,
188 if (level < maxLevel)
191 m_vectOperators[level]->AMRResidualNC(thisResid,
196 (m_vectOperators[level+1].operator->()));
201 m_vectOperators[level]->residual(thisResid,
221 CH_TIME(
"MultilevelLinearOp::preCond");
224 int numLevels = a_cor.size();
225 CH_assert (a_residual.size() == numLevels);
226 CH_assert (m_vectOperators.size() <= numLevels);
229 if (m_use_multigrid_preconditioner)
231 CH_TIME(
"MultilevelLinearOp::preCond::Multigrid");
235 create(localCorr, a_cor);
236 create(localResid, a_residual);
244 CH_TIME(
"MultilevelLinearOp::preCond::Multigrid::NewStorage");
247 for (
int lev=0; lev<numLevels; lev++)
249 tempCorr[lev] = localCorr[lev];
250 tempResid[lev] =
const_cast<LevelData<T>*
>(localResid[lev]);
254 int finestLevel = numLevels -1;
256 CH_TIME(
"MultilevelLinearOp::preCond::Multigrid::Initialize");
258 if (!m_isPrecondSolverInitialized)
261 m_preCondSolver.init(tempCorr, tempResid, finestLevel,
264 m_preCondSolver.setBottomSolver(finestLevel, m_lBase);
265 m_isPrecondSolverInitialized =
true;
269 for (
int iter = 0; iter<m_num_mg_iterations; iter++)
271 CH_TIME(
"MultilevelLinearOp::preCond::Multigrid::Iterate");
280 bool homogeneous =
true;
281 residual(localResid, a_cor, a_residual, homogeneous);
282 setToZero(localCorr);
284 m_preCondSolver.AMRVCycle(tempCorr, tempResid,
285 finestLevel, finestLevel,
289 incr(a_cor, localCorr, 1.0);
293 CH_TIME(
"MultilevelLinearOp::preCond::Multigrid::Clear");
301 CH_TIME(
"MultilevelLinearOp::preCond::LevelBased");
304 for (
int level = m_lBase; level < numLevels; level++)
306 m_vectOperators[level]->preCond(*a_cor[level], *a_residual[level]);
323 CH_TIME(
"MultilevelLinearOp::applyOp");
325 int numLevels = a_phi.size();
326 int maxLevel = numLevels-1;
328 for (
int level = m_lBase; level<numLevels; level++)
332 if (level < maxLevel)
334 m_vectOperators[level]->AMROperatorNC(*a_lhs[level],
338 m_vectOperators[level+1].operator->());
343 m_vectOperators[level]->applyOp(*a_lhs[level],
351 if (level < maxLevel)
353 m_vectOperators[level]->AMROperator(*a_lhs[level],
358 m_vectOperators[level+1].operator->());
362 m_vectOperators[level]->AMROperatorNF(*a_lhs[level],
382 CH_TIME(
"MultilevelLinearOp::create");
384 a_lhs.resize(a_rhs.size());
386 int coarsestLevel =
Max(0, m_lBase-1);
387 for (
int level =coarsestLevel; level < a_rhs.size(); level++)
390 m_vectOperators[level]->create(*a_lhs[level],
404 CH_TIME(
"MultilevelLinearOp::clear");
406 for (
int level =0; level < a_lhs.size(); level++)
408 if (a_lhs[level] != NULL)
425 CH_TIME(
"MultilevelLinearOp::assign");
427 CH_assert (a_lhs.size() == a_rhs.size());
429 for (
int level = m_lBase; level < a_lhs.size(); level++)
431 if (!(a_lhs[level] ==NULL) && !(a_rhs[level] == NULL))
433 m_vectOperators[level]->assign(*a_lhs[level], *a_rhs[level]);
448 CH_TIME(
"MultilevelLinearOp::dotProduct");
465 for (
int level =m_lBase; level<temp1.
size()-1; level++)
473 int nRefFine = m_refRatios[level];
479 for (levelDit.
begin(); levelDit.
ok(); ++levelDit)
481 const Box& thisBox = levelGrids[levelDit];
482 for (finerLit.
begin(); finerLit.
ok(); ++finerLit)
484 Box testBox = finerGrids[finerLit];
489 temp1Level[levelDit].setVal(0.0, testBox,
490 0, temp1Level.
nComp());
492 temp2Level[levelDit].setVal(0.0, testBox,
493 0, temp2Level.
nComp());
501 for (
int level =m_lBase; level<temp1.
size(); level++)
503 Real levelProd = m_vectOperators[level]->dotProduct(*temp1[level],
508 Real scale = D_TERM(dxLev[0], *dxLev[1], *dxLev[2]);
529 CH_TIME(
"MultilevelLinearOp::incr");
532 for (
int level=m_lBase; level<a_lhs.size(); level++)
534 m_vectOperators[level]->incr(*a_lhs[level], *a_x[level], a_scale);
550 CH_TIME(
"MultilevelLinearOp::axby");
553 for (
int level=m_lBase; level<a_lhs.size(); ++level)
555 m_vectOperators[level]->axby(*a_lhs[level],
571 CH_TIME(
"MultilevelLinearOp::scale");
574 for (
int level =m_lBase; level<a_lhs.size(); level++)
576 m_vectOperators[level]->scale(*a_lhs[level], a_scale);
591 CH_TIME(
"MultilevelLinearOp::norm");
594 int maxLevel = a_rhs.size()-1;
596 for (
int level = m_lBase; level<a_rhs.size(); level++)
600 if (level < maxLevel)
603 int refRatio = m_refRatios[level];
604 normLevel = m_vectOperators[level]->AMRNorm(*a_rhs[level],
614 normLevel = m_vectOperators[level]->AMRNorm(*a_rhs[level],
623 MayDay::Error(
"MultilevelLinearOp::norm -- undefined for order > 2");
628 normLevel = normLevel*normLevel;
633 thisNorm += normLevel;
635 else if (normLevel > thisNorm)
637 thisNorm = normLevel;
643 thisNorm = sqrt(thisNorm);
657 CH_TIME(
"MultilevelLinearOp::setToZero");
661 int coarsestLevel =
Max(0, m_lBase-1);
662 for (
int level = coarsestLevel; level< a_lhs.size(); level++)
664 m_vectOperators[level]->setToZero(*a_lhs[level]);
671 const char* a_filename)
676 MayDay::Warning(
"MultilevelLinearOp<T>::write unimplemented since CH_USE_HDF5 undefined");
680 #include "NamespaceFooter.H" Box & coarsen(int refinement_ratio)
virtual void incr(Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_x, Real a_scale)
Definition: MultilevelLinearOpI.H:525
virtual Real norm(const Vector< LevelData< T > * > &a_rhs, int a_ord)
Definition: MultilevelLinearOpI.H:588
A reference-counting handle class.
Definition: RefCountedPtr.H:173
virtual void write(const Vector< LevelData< T > * > *a_data, const char *a_filename)
Definition: MultilevelLinearOpI.H:670
#define CH_assert(cond)
Definition: CHArray.H:37
int m_verbosity
Definition: BiCGStabSolver.H:78
void begin()
initialize this iterator to the first Box in its Layout
Definition: LayoutIterator.H:122
virtual bool ok() const
return true if this iterator is still in its Layout
Definition: LayoutIterator.H:117
IndexTM< T, N > scale(const IndexTM< T, N > &a_p, T a_s)
Definition: IndexTMI.H:382
Definition: DataIterator.H:190
An Iterator based on a BoxLayout object.
Definition: LayoutIterator.H:35
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)
Definition: MultilevelLinearOpI.H:544
void resize(unsigned int isize)
Definition: Vector.H:346
void getBoxes(Vector< Vector< Box > > &a_boxes, Vector< int > &a_refRat, const Box &a_domain)
virtual ~MultilevelLinearOp()
Definition: MultilevelLinearOpI.H:118
virtual void setToZero(Vector< LevelData< T > * > &a_lhs)
Definition: MultilevelLinearOpI.H:655
#define CH_TIME(name)
Definition: CH_Timer.H:82
MultilevelLinearOp()
default constructor (sets defaults for preconditioner)
Definition: MultilevelLinearOpI.H:21
virtual void scale(Vector< LevelData< T > * > &a_lhs, const Real &a_scale)
Definition: MultilevelLinearOpI.H:568
new code
Definition: BoxLayoutData.H:170
double Real
Definition: REAL.H:33
A BoxLayout that has a concept of disjointedness.
Definition: DisjointBoxLayout.H:30
virtual void clear(Vector< LevelData< T > * > &a_lhs)
Definition: MultilevelLinearOpI.H:402
size_t size() const
Definition: Vector.H:192
static void Error(const char *const a_msg=m_nullString, int m_exitCode=CH_DEFAULT_ERROR_CODE)
Print out message to cerr and exit with the specified exit code.
virtual void assign(Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_rhs)
Definition: MultilevelLinearOpI.H:422
virtual void applyOp(Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_phi, bool a_homogeneous=false)
Definition: MultilevelLinearOpI.H:319
static void Warning(const char *const a_msg=m_nullString)
Print out message to cerr and continue.
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
Definition: MultilevelLinearOpI.H:137
A Rectangular Domain on an Integer Lattice.
Definition: Box.H:465
A Real vector in SpaceDim-dimensional space.
Definition: RealVect.H:41
virtual void preCond(Vector< LevelData< T > * > &a_cor, const Vector< LevelData< T > * > &a_residual)
Apply preconditioner.
Definition: MultilevelLinearOpI.H:217
int nComp() const
Definition: BoxLayoutData.H:306
T Max(const T &a_a, const T &a_b)
Definition: Misc.H:39
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
Definition: MultilevelLinearOpI.H:35
virtual Real dotProduct(const Vector< LevelData< T > * > &a_1, const Vector< LevelData< T > * > &a_2)
Definition: MultilevelLinearOpI.H:445
void writeVectorLevelName(const Vector< LevelData< FArrayBox > *> *a_dataPtr, const Vector< int > *a_refRatios, const char *a_filename)
LayoutIterator layoutIterator() const
Iterator that processes through ALL the boxes in a BoxLayout.
Definition: AMRMultiGrid.H:233
DataIterator dataIterator() const
Parallel iterator.
bool isEmpty() const
{ Comparison Functions}
Definition: Box.H:1846
virtual void create(Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_rhs)
Definition: MultilevelLinearOpI.H:379
Definition: BiCGStabSolver.H:26