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 ddtScheme_H
00036 #define ddtScheme_H
00037 
00038 #include <OpenFOAM/tmp.H>
00039 #include <OpenFOAM/dimensionedType.H>
00040 #include <finiteVolume/volFieldsFwd.H>
00041 #include <finiteVolume/surfaceFieldsFwd.H>
00042 #include <OpenFOAM/typeInfo.H>
00043 #include <OpenFOAM/runTimeSelectionTables.H>
00044 
00045 
00046 
00047 namespace Foam
00048 {
00049 
00050 template<class Type>
00051 class fvMatrix;
00052 
00053 class fvMesh;
00054 
00055 
00056 
00057 namespace fv
00058 {
00059 
00060 
00061 
00062 
00063 
00064 template<class Type>
00065 class ddtScheme
00066 :
00067     public refCount
00068 {
00069 
00070 protected:
00071 
00072     
00073 
00074         const fvMesh& mesh_;
00075 
00076 
00077     
00078 
00079         
00080         ddtScheme(const ddtScheme&);
00081 
00082         
00083         void operator=(const ddtScheme&);
00084 
00085 
00086 public:
00087 
00088     
00089     virtual const word& type() const = 0;
00090 
00091 
00092     
00093 
00094         declareRunTimeSelectionTable
00095         (
00096             tmp,
00097             ddtScheme,
00098             Istream,
00099             (const fvMesh& mesh, Istream& schemeData),
00100             (mesh, schemeData)
00101         );
00102 
00103 
00104     
00105 
00106         
00107         ddtScheme(const fvMesh& mesh)
00108         :
00109             mesh_(mesh)
00110         {}
00111 
00112         
00113         ddtScheme(const fvMesh& mesh, Istream&)
00114         :
00115             mesh_(mesh)
00116         {}
00117 
00118 
00119     
00120 
00121         
00122         static tmp<ddtScheme<Type> > New
00123         (
00124             const fvMesh& mesh,
00125             Istream& schemeData
00126         );
00127 
00128 
00129     
00130 
00131         virtual ~ddtScheme();
00132 
00133 
00134     
00135 
00136         
00137         const fvMesh& mesh() const
00138         {
00139             return mesh_;
00140         }
00141 
00142         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00143         (
00144             const dimensioned<Type>&
00145         ) = 0;
00146 
00147         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00148         (
00149             const GeometricField<Type, fvPatchField, volMesh>&
00150         ) = 0;
00151 
00152         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00153         (
00154             const dimensionedScalar&,
00155             const GeometricField<Type, fvPatchField, volMesh>&
00156         ) = 0;
00157 
00158         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00159         (
00160             const volScalarField&,
00161             const GeometricField<Type, fvPatchField, volMesh>&
00162         ) = 0;
00163 
00164         virtual tmp<fvMatrix<Type> > fvmDdt
00165         (
00166             GeometricField<Type, fvPatchField, volMesh>&
00167         ) = 0;
00168 
00169         virtual tmp<fvMatrix<Type> > fvmDdt
00170         (
00171             const dimensionedScalar&,
00172             GeometricField<Type, fvPatchField, volMesh>&
00173         ) = 0;
00174 
00175         virtual tmp<fvMatrix<Type> > fvmDdt
00176         (
00177             const volScalarField&,
00178             GeometricField<Type, fvPatchField, volMesh>&
00179         ) = 0;
00180 
00181 
00182         typedef GeometricField
00183         <
00184             typename flux<Type>::type,
00185             fvsPatchField,
00186             surfaceMesh
00187         > fluxFieldType;
00188 
00189         tmp<surfaceScalarField> fvcDdtPhiCoeff
00190         (
00191             const GeometricField<Type, fvPatchField, volMesh>& U,
00192             const fluxFieldType& phi,
00193             const fluxFieldType& phiCorr
00194         );
00195 
00196         tmp<surfaceScalarField> fvcDdtPhiCoeff
00197         (
00198             const GeometricField<Type, fvPatchField, volMesh>& U,
00199             const fluxFieldType& phi
00200         );
00201 
00202         virtual tmp<fluxFieldType> fvcDdtPhiCorr
00203         (
00204             const volScalarField& rA,
00205             const GeometricField<Type, fvPatchField, volMesh>& U,
00206             const fluxFieldType& phi
00207         ) = 0;
00208 
00209         tmp<surfaceScalarField> fvcDdtPhiCoeff
00210         (
00211             const volScalarField& rho,
00212             const GeometricField<Type, fvPatchField, volMesh>& rhoU,
00213             const fluxFieldType& phi
00214         );
00215 
00216         virtual tmp<fluxFieldType> fvcDdtPhiCorr
00217         (
00218             const volScalarField& rA,
00219             const volScalarField& rho,
00220             const GeometricField<Type, fvPatchField, volMesh>& U,
00221             const fluxFieldType& phi
00222         ) = 0;
00223 
00224 
00225         virtual tmp<surfaceScalarField> meshPhi
00226         (
00227             const GeometricField<Type, fvPatchField, volMesh>&
00228         ) = 0;
00229 };
00230 
00231 
00232 
00233 
00234 } 
00235 
00236 
00237 
00238 } 
00239 
00240 
00241 
00242 
00243 
00244 #define makeFvDdtTypeScheme(SS, Type)                                          \
00245                                                                                \
00246 defineNamedTemplateTypeNameAndDebug(SS<Type>, 0);                              \
00247                                                                                \
00248 ddtScheme<Type>::addIstreamConstructorToTable<SS<Type> >                       \
00249     add##SS##Type##IstreamConstructorToTable_;
00250 
00251 
00252 #define makeFvDdtScheme(SS)                                                    \
00253                                                                                \
00254 makeFvDdtTypeScheme(SS, scalar)                                                \
00255 makeFvDdtTypeScheme(SS, vector)                                                \
00256 makeFvDdtTypeScheme(SS, sphericalTensor)                                       \
00257 makeFvDdtTypeScheme(SS, symmTensor)                                            \
00258 makeFvDdtTypeScheme(SS, tensor)                                                \
00259                                                                                \
00260 template<>                                                                     \
00261 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr                              \
00262 (                                                                              \
00263     const volScalarField& rA,                                                  \
00264     const volScalarField& U,                                                   \
00265     const surfaceScalarField& phi                                              \
00266 )                                                                              \
00267 {                                                                              \
00268     notImplemented(#SS"<scalar>::fvcDdtPhiCorr");                              \
00269     return surfaceScalarField::null();                                         \
00270 }                                                                              \
00271                                                                                \
00272 template<>                                                                     \
00273 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr                              \
00274 (                                                                              \
00275     const volScalarField& rA,                                                  \
00276     const volScalarField& rho,                                                 \
00277     const volScalarField& U,                                                   \
00278     const surfaceScalarField& phi                                              \
00279 )                                                                              \
00280 {                                                                              \
00281     notImplemented(#SS"<scalar>::fvcDdtPhiCorr");                              \
00282     return surfaceScalarField::null();                                         \
00283 }
00284 
00285 
00286 
00287 
00288 #ifdef NoRepository
00289 #   include <finiteVolume/ddtScheme.C>
00290 #endif
00291 
00292 
00293 
00294 #endif
00295 
00296