15 #ifndef _AMRMULTIGRID_H_ 16 #define _AMRMULTIGRID_H_ 32 #include "NamespaceHeader.H" 72 return this->
norm(a_coarResid, 0);
95 virtual void AMRResidual(T& a_residual,
const T& a_phiFine,
const T& a_phi,
96 const T& a_phiCoarse,
const T& a_rhs,
97 bool a_homogeneousDomBC,
105 virtual void AMRResidualNF(T& a_residual,
const T& a_phi,
const T& a_phiCoarse,
106 const T& a_rhs,
bool a_homogeneousBC) = 0;
113 virtual void AMRResidualNC(T& a_residual,
const T& a_phiFine,
const T& a_phi,
114 const T& a_rhs,
bool a_homogeneousBC,
122 const T& a_phiFine,
const T& a_phi,
123 const T& a_phiCoarse,
124 bool a_homogeneousDomBC,
134 const T& a_phiCoarse,
135 bool a_homogeneousBC) = 0;
145 bool a_homogeneousBC,
150 virtual void AMRRestrict(T& a_resCoarse,
const T& a_residual,
const T& a_correction,
151 const T& a_coarseCorrection,
bool a_skip_res) = 0;
155 virtual void AMRProlong(T& a_correction,
const T& a_coarseCorrection) = 0;
160 const T& a_coarseCorrection) = 0;
167 const int& a_refRat) = 0;
182 this->
assign(a_lhs, a_rhs);
192 return this->
norm(a_phi, 0);
196 virtual void AMRProlongS(T& a_correction,
const T& a_coarseCorrection,
197 T& a_temp,
const Copier& a_copier)
204 T& a_temp,
const Copier& a_copier,
205 const Copier& a_cornerCopier,
208 AMRProlong( a_correction, a_coarseCorrection );
211 virtual void AMRRestrictS(T& a_resCoarse,
const T& a_residual,
const T& a_correction,
212 const T& a_coarseCorrection, T& scratch,
bool a_skip_res =
false )
214 AMRRestrict(a_resCoarse, a_residual, a_correction, a_coarseCorrection, a_skip_res);
251 virtual int refToFiner(
const ProblemDomain& a_indexSpace)
const =0;
280 virtual void recordResiduals(
const Vector<T*>& a_residuals,
292 virtual void recordCorrections(
const Vector<T*>& a_corrections,
337 m_inspectors.push_back(a_inspector);
346 int l_max,
int l_base,
bool a_zeroPhi=
true,
347 bool forceHomogeneous =
false );
354 int l_max,
int l_base,
bool a_zeroPhi=
true,
355 bool forceHomogeneous =
false);
360 int l_max,
int l_base,
bool a_zeroPhi=
true,
361 bool forceHomogeneous =
false);
364 int l_max,
int l_base);
366 virtual void AMRVCycle(
Vector<T*>& a_correction,
368 int l,
int l_max,
int l_base);
370 void setMGCycle(
int a_numMG);
373 int l_max,
int l_base);
376 int l_max,
int l_base);
381 int m_pre, m_post, m_bottom, m_numMG;
408 bool a_homogeneousBC=
false,
409 bool a_computeNorm=
true);
425 bool a_homogeneousBC=
false);
451 void setSolverParameters(
const int& a_pre,
455 const int& a_iterMax,
458 const Real& a_normThresh);
466 void setBottomSolver(
int l_max,
int l_base);
468 void setBottomSolverEpsCushion(
Real a_bottomSolverEpsCushion);
471 void getInfo()
const;
474 void relax(T& phi, T& R,
int depth,
int nRelax = 2);
476 void computeAMRResidualLevel(
Vector<T*>& a_resid,
479 int l_max,
int l_base,
int ilev,
480 bool a_homogeneousBC);
519 for (
int ilev = a_lbase; ilev <= a_lmax; ilev++)
523 m_op[a_lbase]->outputAMR(outputData, a_name);
532 const int& a_iterMax,
535 const Real& a_normThresh)
537 m_solverParamsSet =
true;
543 m_normThresh = a_normThresh;
544 m_iterMax = a_iterMax;
545 for (
int img = 0; img < m_mg.size(); img++)
547 m_mg[img]->m_pre = a_pre;
548 m_mg[img]->m_post = a_post;
549 m_mg[img]->m_bottom = a_bottom;
552 m_bottomSolverEpsCushion = 1.0;
559 for (
int iop = 0; iop <
m_op.size(); iop++)
565 for (
int img = 0; img < m_mg.size(); img++)
578 for (
int iop = 0; iop <
m_op.size(); iop++)
592 for (
int img = 0; img < m_mg.size(); img++)
594 retval[img] = m_mg[img]->getAllOperators();
613 m_convergenceMetric(0.),
614 m_bottomSolverEpsCushion(1.0),
615 m_bottomSolver(NULL),
626 for (
int ilev = 0; ilev <
m_op.
size(); ilev++)
628 m_mg[ilev]->m_numMG = a_numMG;
629 m_mg[ilev]->m_cycle = a_numMG;
636 CH_TIME(
"AMRMultiGrid::relax");
638 if (
m_op[depth]->refToCoarser() > 2)
643 int intermediateDepth = 0;
644 int r =
m_op[depth]->refToCoarser();
650 int tmp =
m_mg[depth]->m_depth;
651 m_mg[depth]->m_depth = intermediateDepth;
656 m_mg[depth]->cycle(0, a_correction, a_residual);
657 m_mg[depth]->m_depth = tmp;
663 m_op[depth]->relax(a_correction, a_residual, a_numSmooth);
678 return *(
m_op[level]);
688 bool a_homogeneousBC,
691 CH_TIME(
"AMRMultiGrid::computeAMRResidual");
695 for (
int ilev = l_base; ilev <= l_max; ilev++)
701 l_max, l_base, ilev, a_homogeneousBC);
706 localNorm =
m_op[ilev]->localMaxNorm(*a_resid[ilev]);
711 localNorm =
m_op[ilev]->localMaxNorm(*a_resid[ilev]);
713 const int normType = 0;
715 localNorm =
m_op[ilev]->norm(*a_resid[ilev],normType);
716 rnorm =
Max(localNorm, rnorm);
724 int result = MPI_Allreduce(&rnorm, &recv, 1,
MPI_CH_REAL,
725 MPI_MAX, Chombo_MPI::comm);
726 if (result != MPI_SUCCESS)
729 MayDay::Error(
"sorry, but I had a communcation error on norm");
783 bool a_homogeneousBC)
785 CH_TIME(
"AMRMultiGrid<T>::computeAMROperator");
787 for(
int ilev=l_base; ilev<=l_max; ilev++)
793 m_op[l_max]->AMROperatorNF(*(a_lphi[l_max]),*(a_phi[l_max]),
794 *(a_phi[l_max-1]), a_homogeneousBC);
796 else if (ilev == l_base)
800 m_op[l_base]->AMROperatorNC(*(a_lphi[l_base]), *(a_phi[l_base+1]),
802 a_homogeneousBC,
m_op[l_base+1]);
806 m_op[l_base]->AMROperator(*a_lphi[l_base], *a_phi[l_base+1], *a_phi[l_base],
808 a_homogeneousBC,
m_op[l_base+1]);
813 m_op[ilev]->AMROperator(*a_lphi[ilev], *a_phi[ilev+1], *a_phi[ilev],
815 a_homogeneousBC,
m_op[ilev+1]);
823 m_op[l_max]->applyOp(*a_lphi[l_max], *a_phi[l_max], a_homogeneousBC);
827 m_op[l_max]->AMROperatorNF(*(a_lphi[l_max]), *(a_phi[l_max]),
840 int l_max,
int l_base,
int ilev,
841 bool a_homogeneousBC)
843 CH_TIME(
"AMRMultiGrid<T>::computeAMRResidualLevel");
851 m_op[l_max]->AMRResidualNF(*(a_resid[l_max]), *(a_phi[l_max]),
852 *(a_phi[l_max-1]), *(a_rhs[l_max]),
855 else if (ilev == l_base)
859 m_op[l_base]->AMRResidualNC(*(a_resid[l_base]), *(a_phi[l_base+1]),
860 *(a_phi[l_base]), *(a_rhs[l_base]),
861 a_homogeneousBC,
m_op[l_base+1]);
865 m_op[l_base]->AMRResidual(*a_resid[l_base], *a_phi[l_base+1], *a_phi[l_base],
866 *a_phi[l_base-1], *a_rhs[l_base],
867 a_homogeneousBC,
m_op[l_base+1]);
872 m_op[ilev]->AMRResidual(*a_resid[ilev], *a_phi[ilev+1], *a_phi[ilev],
873 *a_phi[ilev-1], *a_rhs[ilev],
874 a_homogeneousBC,
m_op[ilev+1]);
882 m_op[l_max]->residual(*a_resid[l_max], *a_phi[l_max], *a_rhs[l_max],a_homogeneousBC);
886 m_op[l_max]->AMRResidualNF(*(a_resid[l_max]), *(a_phi[l_max]),
887 *(a_phi[l_max-1]), *(a_rhs[l_max]),
896 int l_max,
int l_base,
bool a_zeroPhi,
897 bool a_forceHomogeneous)
899 CH_TIME(
"AMRMultiGrid::solve");
900 init(a_phi, a_rhs, l_max, l_base);
901 solveNoInit(a_phi, a_rhs, l_max, l_base, a_zeroPhi, a_forceHomogeneous);
903 revert(a_phi, a_rhs, l_max, l_base);
909 int l_max,
int l_base,
bool a_zeroPhi,
910 bool a_forceHomogeneous)
912 CH_TIME(
"AMRMultiGrid::solveNoInit");
918 for (
int ilev = lowlim; ilev <= l_max; ilev++)
920 uberResidual[ilev] =
new T();
923 m_op[ilev]->create(*uberResidual[ilev], *a_rhs[ilev]);
926 solveNoInitResid(a_phi, uberResidual, a_rhs, l_max, l_base, a_zeroPhi, a_forceHomogeneous);
927 for (
int i = lowlim; i <= l_max; i++)
930 delete uberResidual[i];
937 int l_max,
int l_base,
bool a_zeroPhi,
938 bool a_forceHomogeneous)
940 CH_TIMERS(
"AMRMultiGrid::solveNo-InitResid");
941 CH_TIMER(
"AMRMultiGrid::AMRVcycle", vtimer);
955 bool outputIntermediates =
false;
957 for (
int ilev = lowlim; ilev <= l_max; ilev++)
959 uberCorrection[ilev] =
new T();
960 m_op[ilev]->create(*uberCorrection[ilev], *a_phi[ilev]);
963 m_op[ilev]->create(*uberResidual[ilev], *a_rhs[ilev]);
964 m_op[ilev]->setToZero(*(uberResidual[ilev]));
966 m_op[ilev]->setToZero(*(uberCorrection[ilev]));
970 for (
int ilev = l_base; ilev <=l_max; ++ilev)
972 m_op[ilev]->setToZero(*(a_phi[ilev]));
978 CH_TIME(
"Initial AMR Residual");
998 << setiosflags(ios::showpoint)
999 << setiosflags(ios::scientific);
1000 pout() <<
" AMRMultiGrid:: iteration = " << iter <<
", residual norm = " << rnorm << std::endl;
1006 bool goHang = iter <
m_imin || rnorm <(1-
m_hang)*norm_last;
1008 while (goMin || (goIter && goRedu && goHang && goNorm))
1013 m_inspectors[i]->recordResiduals(uberResidual, l_base, l_max, iter);
1015 if (outputIntermediates)
1017 char strresname[100];
1018 sprintf(strresname,
"amrmg.res.iter.%03d", iter);
1019 string nameres(strresname);
1020 outputAMR(uberResidual, nameres, l_max, l_base);
1027 AMRVCycle(uberCorrection, uberResidual, l_max, l_max, l_base);
1031 sprintf(charname,
"resid_iter%03d.%dd.hdf5", iter,
SpaceDim);
1032 string sname(charname);
1036 m_inspectors[i]->recordCorrections(uberCorrection, l_base, l_max, iter);
1039 for (
int ilev = l_base; ilev <= l_max; ilev++)
1041 if (outputIntermediates)
1043 char strcorname[100];
1044 sprintf(strcorname,
"amrmg.phi.iter.%03d", iter);
1045 string namecor(strcorname);
1046 outputAMR(a_phi, namecor, l_max, l_base);
1049 m_op[ilev]->incr(*(a_phi[ilev]), *(uberCorrection[ilev]), 1.0);
1051 if (outputIntermediates)
1053 char strcorname[100];
1054 sprintf(strcorname,
"amrmg.cor.iter.%03d", iter);
1055 string namecor(strcorname);
1056 outputAMR(uberCorrection, namecor, l_max, l_base);
1059 m_op[ilev]->setToZero(*(uberCorrection[ilev]));
1066 if (
m_op[0]->orderOfAccuracy()>2)
1068 for (
int ilev=l_max; ilev>l_base; ilev--)
1070 m_op[ilev]->enforceCFConsistency(*a_phi[ilev-1], *a_phi[ilev]);
1076 rnorm =
computeAMRResidual(uberResidual, a_phi, a_rhs, l_max, l_base, a_forceHomogeneous,
true);
1080 pout() <<
" AMRMultiGrid:: iteration = " << iter <<
", residual norm = " << rnorm;
1083 pout() <<
", rate = " << norm_last/rnorm;
1085 pout() << std::endl;
1099 m_exitStatus = int(!goRedu) + int(!goIter)*2 + int(!goHang)*4 + int(!goNorm)*8;
1102 pout() <<
" AMRMultiGrid:: iteration = " << iter <<
", residual norm = " << rnorm << std::endl;
1106 if (!goIter && goRedu && goNorm)
1108 pout() <<
" AMRMultiGrid:: WARNING: Exit because max iteration count exceeded" << std::endl;
1110 if (!goHang && goRedu && goNorm)
1112 pout() <<
" AMRMultiGrid:: WARNING: Exit because of solver hang" << std::endl;
1119 for (
int i = lowlim; i <= l_max; i++)
1122 delete uberCorrection[i];
1128 int l_max,
int l_base)
1130 CH_TIME(
"AMRMultiGrid::relaxOnly");
1133 init(a_phi, a_rhs, l_max, l_base);
1139 for (
int ilev = 0; ilev <
m_op.
size(); ilev++)
1141 uberResidual[ilev] =
new T();
1142 m_op[ilev]->create(*uberResidual[ilev], *a_rhs[ilev]);
1146 m_op[0]->residual(*uberResidual[0], *a_phi[0], *a_rhs[0],
true);
1154 pout() <<
" AMRMultiGrid::relaxOnly iteration = " << iter <<
", residual norm = " << rnorm << std::endl;
1159 bool goHang = iter <
m_imin || rnorm <(1-
m_hang)*norm_last;
1160 while (goIter && goRedu && goHang && goNorm)
1163 m_op[0]->relax(*a_phi[0], *a_rhs[0], 1);
1167 m_op[0]->residual(*uberResidual[0], *a_phi[0], *a_rhs[0],
true);
1168 rnorm =
m_op[0]->norm(*uberResidual[0], 2);
1171 pout() <<
" AMRMultiGrid::relaxOnly iteration = " << iter <<
", residual norm = " << rnorm
1172 <<
", rate = " << norm_last/rnorm << std::endl;
1180 for (
int i = 0; i <
m_op.
size(); i++)
1183 delete uberResidual[i];
1190 for (
int i = 0; i <
m_op.
size(); i++)
1213 int l_max,
int l_base)
1216 CH_TIME(
"AMRMultiGrid::init");
1219 for (
int i = l_base; i <= l_max; i++)
1229 int intermediateDepth = 0;
1234 intermediateDepth++;
1236 m_mg[i]->m_depth = intermediateDepth;
1238 m_mg[i]->init(*a_phi[i], *a_rhs[i]);
1253 int l_max,
int l_base)
1255 CH_TIME(
"AMRMultiGrid::revert");
1256 for (
int i = l_base; i <= l_max; i++)
1260 m_mg[i]->m_depth =
m_mg[i]->m_defaultDepth;
1269 for (
int ilev = l_base; ilev <= l_max; ilev++)
1287 CH_assert((a_bottomSolverEpsCushion > 0.) &&
1288 (a_bottomSolverEpsCushion <= 1.));
1298 CH_TIME(
"AMRMultiGrid::define");
1306 m_resC. resize( a_maxAMRLevels, NULL);
1312 for (
int i = 0; i < a_maxAMRLevels; i++)
1323 if (i < a_maxAMRLevels-1)
1334 for(
unsigned int i=0; i<
m_mg.
size(); ++i)
1336 pout()<<
"MG Level "<<i<<
" depth "<<
m_mg[i]->m_depth
1337 <<
" cycles "<<
m_mg[i]->m_cycle<<
"domain "<<
m_mg[i]->m_topLevelDomain<<
"\n";
1344 int ilev,
int l_max,
int l_base)
1348 for (
int level = l_base; level <= l_max; level++)
1350 m_op[level]->assignLocal(*
m_residual[level], *a_uberResidual[level]);
1355 if (l_max == l_base)
1358 m_mg[l_base]->oneCycle(*(a_uberCorrection[ilev]), *(a_uberResidual[ilev]));
1360 else if (ilev == l_base)
1380 l_max, l_base, ilev-1,
1388 *(a_uberCorrection[ilev]));
1396 for (
int img = 0; img <
m_numMG; img++)
1398 AMRVCycle(a_uberCorrection, a_uberResidual, ilev-1, l_max, l_base);
1410 T& dCorr = *(a_uberCorrection[ilev]);
1411 m_op[ilev]->setToZero(dCorr);
1421 #include "NamespaceFooter.H" std::ostream & pout()
Use this in place of std::cout for program output.
int m_bottom
Definition: AMRMultiGrid.H:381
virtual AMRLevelOp< T > * AMRnewOp(const ProblemDomain &a_indexSpace)=0
virtual int refToFiner(const ProblemDomain &a_indexSpace) const =0
#define CH_TIMERS(name)
Definition: CH_Timer.H:133
virtual Real AMRNorm(const T &a_coarResid, const T &a_fineResid, const int &a_refRat, const int &a_ord)
Definition: AMRMultiGrid.H:67
bool m_solverParamsSet
Definition: AMRMultiGrid.H:379
A reference-counting handle class.
Definition: RefCountedPtr.H:173
virtual void AMRUpdateResidual(T &a_residual, const T &a_correction, const T &a_coarseCorrection)=0
virtual void AMROperator(T &a_LofPhi, const T &a_phiFine, const T &a_phi, const T &a_phiCoarse, bool a_homogeneousDomBC, AMRLevelOp< T > *a_finerOp)=0
virtual void AMRResidualNC(T &a_residual, const T &a_phiFine, const T &a_phi, const T &a_rhs, bool a_homogeneousBC, AMRLevelOp< T > *a_finerOp)=0
#define CH_assert(cond)
Definition: CHArray.H:37
virtual void buildCopier(Copier &a_copier, const T &a_lhs, const T &a_rhs)
Definition: AMRMultiGrid.H:176
A class to facilitate interaction with physical boundary conditions.
Definition: ProblemDomain.H:141
void revert(const Vector< T *> &a_phi, const Vector< T *> &a_rhs, int l_max, int l_base)
Definition: AMRMultiGrid.H:1252
virtual void define(const ProblemDomain &a_coarseDomain, AMRLevelOpFactory< T > &a_factory, LinearSolver< T > *a_bottomSolver, int a_numLevels)
Definition: AMRMultiGrid.H:1293
virtual ~AMRMultiGridInspector()
Destructor.
Definition: AMRMultiGrid.H:269
virtual void solve(Vector< T *> &a_phi, const Vector< T *> &a_rhs, int l_max, int l_base, bool a_zeroPhi=true, bool forceHomogeneous=false)
Definition: AMRMultiGrid.H:895
int m_iterMax
Definition: AMRMultiGrid.H:380
virtual void AMROperatorNF(T &a_LofPhi, const T &a_phi, const T &a_phiCoarse, bool a_homogeneousBC)=0
LinearSolver< T > * m_bottomSolver
Definition: AMRMultiGrid.H:492
AMRLevelOp< T > & levelOp(int level)
Definition: AMRMultiGrid.H:676
Real m_initial_rnorm
initial residual from this solve
Definition: AMRMultiGrid.H:394
A strange but true thing to make copying from one boxlayoutdata to another fast.
Definition: Copier.H:145
int m_verbosity
Definition: AMRMultiGrid.H:380
#define CH_START(tpointer)
Definition: CH_Timer.H:145
Vector< T * > m_residual
Definition: AMRMultiGrid.H:485
AMRMultiGridInspector()
Base class constructor. This must be called by all subclasses.
Definition: AMRMultiGrid.H:264
virtual void AMRRestrictS(T &a_resCoarse, const T &a_residual, const T &a_correction, const T &a_coarseCorrection, T &scratch, bool a_skip_res=false)
Definition: AMRMultiGrid.H:211
virtual Real norm(const T &a_rhs, int a_ord)=0
int m_pre
Definition: AMRMultiGrid.H:381
void addInspector(RefCountedPtr< AMRMultiGridInspector< T > > &a_inspector)
Definition: AMRMultiGrid.H:334
Definition: NoOpSolver.H:20
Vector< Copier > m_reverseCopier
Definition: AMRMultiGrid.H:488
virtual void AMRProlongS(T &a_correction, const T &a_coarseCorrection, T &a_temp, const Copier &a_copier)
Definition: AMRMultiGrid.H:196
virtual void createCoarsened(T &a_lhs, const T &a_rhs, const int &a_refRat)=0
Real m_eps
Definition: AMRMultiGrid.H:378
Definition: AMRMultiGrid.H:308
void clear()
Definition: AMRMultiGrid.H:1188
Vector< MultiGrid< T > * > m_mg
Definition: AMRMultiGrid.H:483
void computeAMRResidualLevel(Vector< T *> &a_resid, Vector< T *> &a_phi, const Vector< T *> &a_rhs, int l_max, int l_base, int ilev, bool a_homogeneousBC)
Definition: AMRMultiGrid.H:837
Real m_convergenceMetric
Definition: AMRMultiGrid.H:391
virtual void AMRResidualNF(T &a_residual, const T &a_phi, const T &a_phiCoarse, const T &a_rhs, bool a_homogeneousBC)=0
const int SpaceDim
Definition: SPACE.H:38
void append(const Vector< T > &invec)
Definition: Vector.H:330
Definition: AMRMultiGrid.H:39
#define MPI_CH_REAL
Definition: REAL.H:34
virtual void assignCopier(T &a_lhs, const T &a_rhs, const Copier &a_copier)
Definition: AMRMultiGrid.H:180
virtual void setToZero(T &a_lhs)=0
virtual void AMRProlongS_2(T &a_correction, const T &a_coarseCorrection, T &a_temp, const Copier &a_copier, const Copier &a_cornerCopier, const AMRLevelOp< LevelData< FArrayBox > > *a_crsOp)
Definition: AMRMultiGrid.H:203
#define CH_TIMER(name, tpointer)
Definition: CH_Timer.H:63
void resize(unsigned int isize)
Definition: Vector.H:346
virtual void dumpStuff(Vector< T *> data, string filename)
Definition: AMRMultiGrid.H:58
Vector< MGLevelOp< T > *> getOperatorsOp()
Definition: AMRMultiGrid.H:575
MGLevelOp & operator=(const MGLevelOp &)
virtual void AMRVCycle(Vector< T *> &a_correction, Vector< T *> &a_residual, int l, int l_max, int l_base)
Definition: AMRMultiGrid.H:1342
virtual void AMRResidual(T &a_residual, const T &a_phiFine, const T &a_phi, const T &a_phiCoarse, const T &a_rhs, bool a_homogeneousDomBC, AMRLevelOp< T > *a_finerOp)=0
int m_numMG
Definition: AMRMultiGrid.H:381
void push_back(const T &in)
Definition: Vector.H:295
void clear()
Definition: Vector.H:180
Vector< MGLevelOp< T > *> getAllOperators()
Definition: AMRMultiGrid.H:556
virtual int refToCoarser()=0
#define CH_TIME(name)
Definition: CH_Timer.H:82
Definition: MultiGrid.H:294
const char * name(const FArrayBox &a_dummySpecializationArg)
Definition: CH_HDF5.H:907
Real m_bottomSolverEpsCushion
Definition: AMRMultiGrid.H:397
int m_maxDepth
max no. of coarsenings – -1 (default) means coarsen as far as possible
Definition: AMRMultiGrid.H:385
double Real
Definition: REAL.H:33
Definition: MultiGrid.H:30
virtual void outputAMR(Vector< T *> &a_rhs, string &a_name)
Definition: AMRMultiGrid.H:81
AMRMultiGrid()
Definition: AMRMultiGrid.H:600
Real m_hang
Definition: AMRMultiGrid.H:378
int m_iterMin
Definition: AMRMultiGrid.H:380
Vector< RefCountedPtr< AMRMultiGridInspector< T > > > m_inspectors
Definition: AMRMultiGrid.H:501
size_t size() const
Definition: Vector.H:192
Vector< T * > m_correction
Definition: AMRMultiGrid.H:484
Vector< T * > m_resC
Definition: AMRMultiGrid.H:486
AMRLevelOp()
Constructor.
Definition: AMRMultiGrid.H:53
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 enforceCFConsistency(T &a_coarseCorrection, const T &a_correction)
This routine is for operators with orderOfAccuracy()>2.
Definition: AMRMultiGrid.H:223
Definition: MultiGrid.H:331
virtual void clear(T &a_lhs)
Definition: LinearSolver.H:70
Definition: AMRMultiGrid.H:259
virtual void outputLevel(T &a_rhs, string &a_name)
Definition: AMRMultiGrid.H:76
void outputAMR(Vector< T *> &a_data, string &a_name, int a_lmax, int a_lbase)
Definition: AMRMultiGrid.H:516
int m_exitStatus
Definition: AMRMultiGrid.H:380
void getInfo() const
Dump out the state of the solver. AMR levels, MG levels, bottom solver, box counts, etc.
Definition: AMRMultiGrid.H:1331
ProblemDomain & refine(int a_refinement_ratio)
Refine this problem domain.
virtual void AMROperatorNC(T &a_LofPhi, const T &a_phiFine, const T &a_phi, bool a_homogeneousBC, AMRLevelOp< T > *a_finerOp)=0
virtual void solveNoInitResid(Vector< T *> &a_phi, Vector< T *> &a_finalResid, const Vector< T *> &a_rhs, int l_max, int l_base, bool a_zeroPhi=true, bool forceHomogeneous=false)
use if you want final residual
Definition: AMRMultiGrid.H:935
virtual void create(T &a_lhs, const T &a_rhs)=0
virtual void AMRProlong(T &a_correction, const T &a_coarseCorrection)=0
Definition: LinearSolver.H:156
virtual void relax(T &a_correction, const T &a_residual, int a_iterations)=0
void relax(T &phi, T &R, int depth, int nRelax=2)
Definition: AMRMultiGrid.H:634
int m_post
Definition: AMRMultiGrid.H:381
int m_imin
Definition: AMRMultiGrid.H:380
virtual void solveNoInit(Vector< T *> &a_phi, const Vector< T *> &a_rhs, int l_max, int l_base, bool a_zeroPhi=true, bool forceHomogeneous=false)
Definition: AMRMultiGrid.H:907
Real computeAMRResidual(Vector< T *> &a_resid, Vector< T *> &a_phi, const Vector< T *> &a_rhs, int l_max, int l_base, bool a_homogeneousBC=false, bool a_computeNorm=true)
Definition: AMRMultiGrid.H:683
void computeAMROperator(Vector< T *> &a_lph, Vector< T *> &a_phi, int l_max, int l_base, bool a_homogeneousBC=false)
Definition: AMRMultiGrid.H:779
#define CH_STOP(tpointer)
Definition: CH_Timer.H:150
Vector< char > m_hasInitBeenCalled
Definition: AMRMultiGrid.H:494
void setSolverParameters(const int &a_pre, const int &a_post, const int &a_bottom, const int &a_numMG, const int &a_iterMax, const Real &a_eps, const Real &a_hang, const Real &a_normThresh)
Definition: AMRMultiGrid.H:528
Vector< Vector< MGLevelOp< T > *> > getOperatorsMG()
Definition: AMRMultiGrid.H:588
T Max(const T &a_a, const T &a_b)
Definition: Misc.H:39
virtual void AMRRestrict(T &a_resCoarse, const T &a_residual, const T &a_correction, const T &a_coarseCorrection, bool a_skip_res)=0
virtual void init(const Vector< T *> &a_phi, const Vector< T *> &a_rhs, int l_max, int l_base)
Definition: AMRMultiGrid.H:1212
virtual void assign(T &a_lhs, const T &a_rhs)=0
Vector< Copier > m_resCopier
Definition: AMRMultiGrid.H:487
virtual Real localMaxNorm(const T &a_phi)
Definition: AMRMultiGrid.H:190
MGLevelOp< T > * m_op
Definition: MultiGrid.H:76
virtual void dumpLevel(T &a_data, string name)
Definition: AMRMultiGrid.H:48
Definition: AMRMultiGrid.H:233
NoOpSolver< T > m_nosolve
Definition: AMRMultiGrid.H:490
virtual void dumpAMR(Vector< T *> &a_data, string name)
Definition: AMRMultiGrid.H:44
Real m_normThresh
Definition: AMRMultiGrid.H:378
virtual ~AMRLevelOpFactory()
Definition: AMRMultiGrid.H:236
Vector< AMRLevelOp< T > *> & getAMROperators()
Definition: AMRMultiGrid.H:434
virtual ~AMRLevelOp()
Destructor.
Definition: AMRMultiGrid.H:63
virtual void zeroCovered(T &a_lhs, T &a_rhs, const Copier &a_copier)
Definition: AMRMultiGrid.H:185
void setMGCycle(int a_numMG)
Definition: AMRMultiGrid.H:624
void relaxOnlyHomogeneous(Vector< T *> &a_phi, const Vector< T *> &a_rhs, int l_max, int l_base)
Definition: AMRMultiGrid.H:1127
void setBottomSolverEpsCushion(Real a_bottomSolverEpsCushion)
Definition: AMRMultiGrid.H:1285
virtual unsigned int orderOfAccuracy(void) const
Definition: AMRMultiGrid.H:217
Vector< AMRLevelOp< T > * > m_op
Definition: AMRMultiGrid.H:482
virtual ~AMRMultiGrid()
Definition: AMRMultiGrid.H:669
void setBottomSolver(int l_max, int l_base)
set up bottom solver for internal MG solver
Definition: AMRMultiGrid.H:1266