00001 #ifdef CH_LANG_CC
00002
00003
00004
00005
00006
00007
00008
00009 #endif
00010
00011 #ifndef _PETSCSOLVER_H_
00012 #define _PETSCSOLVER_H_
00013
00014 #include "LinearSolver.H"
00015 #include "parstream.H"
00016 #include "CH_Timer.H"
00017 #include "LevelData.H"
00018
00019 #ifdef CH_USE_PETSC
00020 #include "petsc.h"
00021 #include "petscmat.h"
00022 #include "petscksp.h"
00023 #include "petscviewer.h"
00024 #endif
00025
00026 #include "NamespaceHeader.H"
00027
00028
00029
00030
00031
00032
00033
00034
00035 template <class T>
00036 class PetscSolver : public LinearSolver<T>
00037 {
00038 public:
00039 PetscSolver();
00040 virtual ~PetscSolver();
00041 void destroy();
00042
00043
00044
00045
00046 virtual void setHomogeneous(bool a_homogeneous)
00047 {
00048 m_homogeneous = a_homogeneous;
00049 }
00050
00051
00052
00053 #ifdef CH_USE_PETSC
00054 virtual void setFunctionAndJacobian( PetscErrorCode (*f)(SNES,Vec,Vec,void*),
00055 #if PETSC_VERSION_GE(3,5,0)
00056 PetscErrorCode (*j)(SNES,Vec,Mat,Mat,void*)
00057 #else
00058 PetscErrorCode (*j)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)
00059 #endif
00060 )
00061 {
00062 m_function = f;
00063 m_jacobian = j;
00064 }
00065 static PetscErrorCode ksp_monitor_pout(KSP ksp, PetscInt it, PetscReal rnorm ,void *ctx)
00066 {
00067 pout() << " KSP:: iteration = " << it << " residual norm = " << rnorm << std::endl;
00068 return 0;
00069 }
00070 #endif
00071
00072 virtual void define( Real a_dx, bool a_homogeneous = true );
00073 virtual void define( LinearOp<T> *, bool a_homogeneous = true );
00074
00075
00076
00077
00078 virtual void solve( T& a_phi, const T& a_rhs );
00079 int solve_private(T& a_phi, const T& a_rhs);
00080
00081
00082
00083
00084
00085
00086
00087 virtual void solve_mfree( T& a_phi,
00088 const T& a_rhs,
00089 LinearOp<T> *a_op );
00090
00091 int solve_mfree_private(T& a_phi,
00092 const T& a_rhs,
00093 LinearOp<T> *a_op );
00094
00095 LinearOp<T> *m_op_mfree;
00096 T m_phi_mfree;
00097 T m_Lphi_mfree;
00098 bool m_mfree_homogeneous;
00099
00100 #ifdef CH_USE_PETSC
00101 static PetscErrorCode apply_mfree(Mat A, Vec x, Vec f);
00102 #endif
00103
00104
00105
00106
00107
00108 Real computeResidual( );
00109
00110
00111
00112 virtual int applyOp(T & a_phi, const T & a_rhs );
00113
00114
00115
00116 void setInitialGuessNonzero( bool b = true )
00117 {
00118 m_nz_init_guess = b;
00119 }
00120
00121
00122
00123
00124 void setOptionsPrefix( const char prefix[] )
00125 {
00126 #ifdef CH_USE_PETSC
00127 if (m_ksp)
00128 {
00129 KSPSetOptionsPrefix(m_ksp,prefix);
00130 }
00131 else
00132 {
00133 strcpy(m_prestring, prefix);
00134 }
00135 #endif
00136 }
00137 #ifdef CH_USE_PETSC
00138 KSP getKSP() { return m_ksp; }
00139 #endif
00140
00141
00142
00143
00144 virtual int getNNZPerRow() const
00145 {
00146 return 9;
00147 }
00148
00149
00150
00151 virtual bool supportNNZExact() const
00152 {
00153 return false;
00154 }
00155
00156
00157
00158
00159 virtual void rhsOperation(const T& a_rhs)
00160 {
00161 }
00162
00163
00164
00165
00166 bool m_homogeneous;
00167
00168
00169
00170
00171 public:
00172 Real m_dx;
00173
00174
00175
00176
00177
00178
00179
00180 virtual Real addBCdiagValue(const IntVect& a_iv, const IntVect& a_jv, const T& a_rhs, const DataIndex& a_datInd, const Real coeff = 1)
00181 {
00182 return 0;
00183 }
00184 public:
00185 int resetOperator()
00186 {
00187 #ifdef CH_USE_PETSC
00188 if (m_defined == 2)
00189 {
00190 m_defined = 1;
00191 }
00192 #endif
00193 return 0;
00194 }
00195
00196
00197
00198
00199
00200 Real normInfinity( const T& a_phi ) const;
00201
00202 protected:
00203
00204
00205
00206
00207
00208 virtual Real addBCrhsValue(const IntVect& a_iv, const T& a_phi, const DataIndex& a_datInd, const Real& coeff = 1)
00209 {
00210 return 0.0;
00211 }
00212
00213
00214 #ifdef CH_USE_PETSC
00215 public:
00216 Mat m_mat;
00217 void *m_ctx;
00218
00219 Vec m_xx, m_rr, m_bb;
00220 SNES m_snes;
00221 KSP m_ksp;
00222 protected:
00223 PetscInt m_defined;
00224 PetscErrorCode (*m_function)(SNES,Vec,Vec,void*);
00225 #if PETSC_VERSION_GE(3,5,0)
00226 PetscErrorCode (*m_jacobian)(SNES,Vec,Mat,Mat,void*);
00227 #else
00228 PetscErrorCode (*m_jacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
00229 #endif
00230 #endif
00231 protected:
00232 char m_prestring[32];
00233 bool m_null;
00234 bool m_nz_init_guess;
00235 LevelData<BaseFab<bool> > m_bccode;
00236 public:
00237 int setup_solver( const T& a_phi );
00238 int create_mat_vec( const T& a_phi );
00239 #ifdef CH_USE_PETSC
00240 LevelData<BaseFab<PetscInt> > m_gids;
00241 PetscInt m_gid0;
00242 PetscErrorCode putPetscInChombo( T& a_phi, const Vec xx );
00243 PetscErrorCode putChomboInPetsc( Vec out, const T& a_phi );
00244 virtual int formMatrix( Mat, const T* =0, PetscInt my0=0, PetscInt nloc=0, PetscInt *d=0, PetscInt *o=0) = 0;
00245 #endif
00246 virtual BaseFab<Real>& getRegFab(T& a_fab, const DataIndex& a_datInd) = 0;
00247 const virtual BaseFab<Real>& getRegFab(const T& a_fab, const DataIndex& a_datInd) const = 0;
00248 const virtual BaseFab<Real>& getRegFab(const T& a_fab, const DataIndex& a_datInd, Box& a_box) const = 0;
00249 void setNull( bool n = true );
00250 };
00251
00252
00253
00254
00255 template <class T>
00256 class PetscSolverFAB : public PetscSolver<T>
00257 {
00258 public:
00259 PetscSolverFAB() : PetscSolver<LevelData<FArrayBox> >()
00260 {
00261 }
00262
00263 BaseFab<Real>& getRegFab(LevelData<FArrayBox>& a_fab, const DataIndex& a_datInd)
00264 {
00265 BaseFab<Real>& fabref = static_cast<BaseFab<Real>& >(a_fab[a_datInd]);
00266 return fabref;
00267 }
00268
00269 const BaseFab<Real>& getRegFab(const LevelData<FArrayBox>& a_fab, const DataIndex& a_datInd) const
00270 {
00271 const BaseFab<Real>& fabref = static_cast<const BaseFab<Real>& >(a_fab[a_datInd]);
00272 return fabref;
00273 }
00274
00275 const BaseFab<Real>& getRegFab(const LevelData<FArrayBox>& a_fab, const DataIndex& a_datInd, Box& a_box) const
00276 {
00277 const BaseFab<Real>& fabref = static_cast<const BaseFab<Real>& >(a_fab[a_datInd]);
00278 a_box = fabref.box();
00279 return fabref;
00280 }
00281 };
00282
00283
00284
00285
00286
00287
00288 template <class T>
00289 class PetscSolverPoisson : public PetscSolverFAB<T>
00290 {
00291 public:
00292 PetscSolverPoisson();
00293 virtual void define( Real a_dx, bool a_homogeneous = true );
00294 virtual void define( LinearOp<T> *a_op, bool a_homogeneous = true )
00295 {
00296 PetscSolver<T>::define(a_op,a_homogeneous);
00297 }
00298 #ifdef CH_USE_PETSC
00299 virtual int formMatrix(Mat, const LevelData<FArrayBox>* =0,PetscInt my0=0, PetscInt nloc=0, PetscInt *d=0, PetscInt *o=0 );
00300 #endif
00301 Real m_alpha;
00302 Real m_beta;
00303 };
00304
00305
00306
00307
00308
00309
00310 template <class T>
00311 class PetscSolverViscousTensor : public PetscSolverFAB<LevelData<FArrayBox> >
00312 {
00313 public:
00314 PetscSolverViscousTensor();
00315 virtual void define( LinearOp<T>* a_operator, bool a_homogeneous = true );
00316
00317
00318
00319 virtual int getNNZPerRow() const
00320 {
00321 return 13;
00322 }
00323 #ifdef CH_USE_PETSC
00324 virtual int formMatrix(Mat, const LevelData<FArrayBox>* =0, PetscInt my0=0, PetscInt nloc=0, PetscInt *d=0, PetscInt *o=0 );
00325 #endif
00326 protected:
00327 Real m_alpha, m_beta;
00328 Real m_dxCrse;
00329 LevelData<FArrayBox> *m_a;
00330 LevelData<FluxBox> *m_eta;
00331 LevelData<FluxBox> *m_lamb;
00332 public:
00333 void define( Real alpha, Real beta, LevelData<FArrayBox> *a, LevelData<FluxBox> *eta, LevelData<FluxBox> *lam)
00334 {
00335 m_alpha = alpha;
00336 m_beta = beta;
00337 m_a = a;
00338 m_eta = eta;
00339 m_lamb = lam;
00340 }
00341 };
00342
00343 #include "NamespaceFooter.H"
00344
00345 #ifdef CH_USE_PETSC
00346 #ifndef CH_EXPLICIT_TEMPLATES
00347 #include "PetscSolverI.H"
00348 #endif // CH_EXPLICIT_TEMPLATES
00349 #endif
00350
00351 #endif