00001 #ifdef CH_LANG_CC
00002
00003
00004
00005
00006
00007
00008
00009 #endif
00010
00011
00012
00013
00014
00015 #ifndef _MULTILEVELLINEAROP_H_
00016 #define _MULTILEVELLINEAROP_H_
00017
00018 #include <cmath>
00019
00020 #include "RefCountedPtr.H"
00021 #include "Vector.H"
00022 #include "RealVect.H"
00023 #include "LevelData.H"
00024 #include "LinearSolver.H"
00025 #include "AMRMultiGrid.H"
00026 #include "BiCGStabSolver.H"
00027
00028 #include "NamespaceHeader.H"
00029
00031
00039 template <class T>
00040 class MultilevelLinearOp : public LinearOp< Vector<LevelData<T>* > >
00041 {
00042
00043 public:
00045 MultilevelLinearOp();
00046
00048
00050 void define(const Vector<DisjointBoxLayout>& a_vectGrids,
00051 const Vector<int>& a_refRatios,
00052 const Vector<ProblemDomain>& a_domains,
00053 const Vector<RealVect>& a_vectDx,
00054 RefCountedPtr<AMRLevelOpFactory<LevelData<T> > >& a_opFactory,
00055 int a_lBase);
00056
00057
00059 virtual ~MultilevelLinearOp();
00060
00062
00067 virtual void residual( Vector<LevelData<T>* >& a_lhs,
00068 const Vector<LevelData<T>* >& a_phi,
00069 const Vector<LevelData<T>* >& a_rhs,
00070 bool a_homogeneous = false);
00071
00073
00077 virtual void preCond(Vector<LevelData<T>* >& a_cor,
00078 const Vector<LevelData<T>* >& a_residual);
00079
00081
00086 virtual void applyOp(Vector<LevelData<T>* >& a_lhs,
00087 const Vector<LevelData<T>* >& a_phi,
00088 bool a_homogeneous = false);
00089
00091
00095 virtual void create(Vector<LevelData<T>* >& a_lhs,
00096 const Vector<LevelData<T>* >& a_rhs);
00097
00099
00103 virtual void clear(Vector<LevelData<T>* >& a_lhs);
00104
00105
00107
00110 virtual void assign(Vector<LevelData<T>* >& a_lhs,
00111 const Vector<LevelData<T>* >& a_rhs);
00112
00114
00118 virtual Real dotProduct(const Vector<LevelData<T>* >& a_1,
00119 const Vector<LevelData<T>* >& a_2);
00120
00122
00125 virtual void incr (Vector<LevelData<T>* >& a_lhs,
00126 const Vector<LevelData<T>* >& a_x,
00127 Real a_scale);
00128
00130
00133 virtual void axby(Vector<LevelData<T>* >& a_lhs,
00134 const Vector<LevelData<T>* >& a_x,
00135 const Vector<LevelData<T>* >& a_y,
00136 Real a_a,
00137 Real a_b);
00138
00140
00143 virtual void scale(Vector<LevelData<T>* >& a_lhs,
00144 const Real& a_scale);
00145
00147
00151 virtual Real norm(const Vector<LevelData<T>* >& a_rhs,
00152 int a_ord);
00153
00155
00158 virtual void setToZero(Vector<LevelData<T>* >& a_lhs);
00159
00160
00161 virtual void write(const Vector<LevelData<T>* >* a_data,
00162 const char* a_filename);
00163
00165 bool m_use_multigrid_preconditioner;
00166
00168 int m_num_mg_iterations;
00169
00171 int m_num_mg_smooth;
00172
00174 int m_preCondSolverDepth;
00175
00177 Vector<DisjointBoxLayout> m_vectGrids;
00178
00180 Vector<int> m_refRatios;
00181
00183 Vector<ProblemDomain> m_domains;
00184
00185
00186 Vector<RealVect> m_vectDx;
00187
00189 Vector<RefCountedPtr<AMRLevelOp<LevelData<T> > > > m_vectOperators;
00190
00192 int m_lBase;
00193
00195
00199 AMRMultiGrid<LevelData<T> > m_preCondSolver;
00200
00202 LinearSolver<LevelData<T> >* m_precondBottomSolverPtr;
00203
00205 bool m_isPrecondSolverInitialized;
00206 };
00207
00208 #include "NamespaceFooter.H"
00209
00210 #include "MultilevelLinearOpI.H"
00211 #endif