00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035 #ifndef laplacianScheme_H
00036 #define laplacianScheme_H
00037
00038 #include <OpenFOAM/tmp.H>
00039 #include <finiteVolume/volFieldsFwd.H>
00040 #include <finiteVolume/surfaceFieldsFwd.H>
00041 #include <finiteVolume/linear.H>
00042 #include <finiteVolume/correctedSnGrad.H>
00043 #include <OpenFOAM/typeInfo.H>
00044 #include <OpenFOAM/runTimeSelectionTables.H>
00045
00046
00047
00048 namespace Foam
00049 {
00050
00051 template<class Type>
00052 class fvMatrix;
00053
00054 class fvMesh;
00055
00056
00057
00058 namespace fv
00059 {
00060
00061
00062
00063
00064
00065 template<class Type, class GType>
00066 class laplacianScheme
00067 :
00068 public refCount
00069 {
00070
00071 protected:
00072
00073
00074
00075 const fvMesh& mesh_;
00076 tmp<surfaceInterpolationScheme<GType> > tinterpGammaScheme_;
00077 tmp<snGradScheme<Type> > tsnGradScheme_;
00078
00079
00080
00081
00082
00083 laplacianScheme(const laplacianScheme&);
00084
00085
00086 void operator=(const laplacianScheme&);
00087
00088
00089 public:
00090
00091
00092 virtual const word& type() const = 0;
00093
00094
00095
00096
00097 declareRunTimeSelectionTable
00098 (
00099 tmp,
00100 laplacianScheme,
00101 Istream,
00102 (const fvMesh& mesh, Istream& schemeData),
00103 (mesh, schemeData)
00104 );
00105
00106
00107
00108
00109
00110 laplacianScheme(const fvMesh& mesh)
00111 :
00112 mesh_(mesh),
00113 tinterpGammaScheme_(new linear<GType>(mesh)),
00114 tsnGradScheme_(new correctedSnGrad<Type>(mesh))
00115 {}
00116
00117
00118 laplacianScheme(const fvMesh& mesh, Istream& is)
00119 :
00120 mesh_(mesh),
00121 tinterpGammaScheme_(NULL),
00122 tsnGradScheme_(NULL)
00123 {
00124 tinterpGammaScheme_ = tmp<surfaceInterpolationScheme<GType> >
00125 (
00126 surfaceInterpolationScheme<GType>::New(mesh, is)
00127 );
00128
00129 tsnGradScheme_ = tmp<snGradScheme<Type> >
00130 (
00131 snGradScheme<Type>::New(mesh, is)
00132 );
00133 }
00134
00135
00136 laplacianScheme
00137 (
00138 const fvMesh& mesh,
00139 const tmp<surfaceInterpolationScheme<GType> >& igs,
00140 const tmp<snGradScheme<Type> >& sngs
00141 )
00142 :
00143 mesh_(mesh),
00144 tinterpGammaScheme_(igs),
00145 tsnGradScheme_(sngs)
00146 {}
00147
00148
00149
00150
00151
00152 static tmp<laplacianScheme<Type, GType> > New
00153 (
00154 const fvMesh& mesh,
00155 Istream& schemeData
00156 );
00157
00158
00159
00160
00161 virtual ~laplacianScheme();
00162
00163
00164
00165
00166
00167 const fvMesh& mesh() const
00168 {
00169 return mesh_;
00170 }
00171
00172 virtual tmp<fvMatrix<Type> > fvmLaplacian
00173 (
00174 const GeometricField<GType, fvsPatchField, surfaceMesh>&,
00175 GeometricField<Type, fvPatchField, volMesh>&
00176 ) = 0;
00177
00178 virtual tmp<fvMatrix<Type> > fvmLaplacian
00179 (
00180 const GeometricField<GType, fvPatchField, volMesh>&,
00181 GeometricField<Type, fvPatchField, volMesh>&
00182 );
00183
00184 virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcLaplacian
00185 (
00186 const GeometricField<Type, fvPatchField, volMesh>&
00187 ) = 0;
00188
00189 virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcLaplacian
00190 (
00191 const GeometricField<GType, fvsPatchField, surfaceMesh>&,
00192 const GeometricField<Type, fvPatchField, volMesh>&
00193 ) = 0;
00194
00195 virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcLaplacian
00196 (
00197 const GeometricField<GType, fvPatchField, volMesh>&,
00198 const GeometricField<Type, fvPatchField, volMesh>&
00199 );
00200 };
00201
00202
00203
00204
00205 }
00206
00207
00208
00209 }
00210
00211
00212
00213
00214
00215 #define makeFvLaplacianTypeScheme(SS, Type, GType) \
00216 \
00217 typedef SS<Type, GType> SS##Type##GType; \
00218 defineNamedTemplateTypeNameAndDebug(SS##Type##GType, 0); \
00219 \
00220 laplacianScheme<Type, GType>:: \
00221 addIstreamConstructorToTable<SS<Type, GType> > \
00222 add##SS##Type##GType##IstreamConstructorToTable_;
00223
00224
00225 #define makeFvLaplacianScheme(SS) \
00226 \
00227 makeFvLaplacianTypeScheme(SS, scalar, scalar) \
00228 makeFvLaplacianTypeScheme(SS, scalar, symmTensor) \
00229 makeFvLaplacianTypeScheme(SS, scalar, tensor) \
00230 makeFvLaplacianTypeScheme(SS, vector, scalar) \
00231 makeFvLaplacianTypeScheme(SS, sphericalTensor, scalar) \
00232 makeFvLaplacianTypeScheme(SS, symmTensor, scalar) \
00233 makeFvLaplacianTypeScheme(SS, tensor, scalar) \
00234
00235
00236
00237
00238 #ifdef NoRepository
00239 # include <finiteVolume/laplacianScheme.C>
00240 #endif
00241
00242
00243
00244 #endif
00245
00246