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
00036
00037 #ifndef EulerDdtScheme_H
00038 #define EulerDdtScheme_H
00039
00040 #include "ddtScheme.H"
00041
00042
00043
00044 namespace Foam
00045 {
00046
00047
00048
00049 namespace fv
00050 {
00051
00052
00053
00054
00055
00056 template<class Type>
00057 class EulerDdtScheme
00058 :
00059 public ddtScheme<Type>
00060 {
00061
00062
00063
00064 EulerDdtScheme(const EulerDdtScheme&);
00065
00066
00067 void operator=(const EulerDdtScheme&);
00068
00069
00070 public:
00071
00072
00073 TypeName("Euler");
00074
00075
00076
00077
00078
00079 EulerDdtScheme(const fvMesh& mesh)
00080 :
00081 ddtScheme<Type>(mesh)
00082 {}
00083
00084
00085 EulerDdtScheme(const fvMesh& mesh, Istream& is)
00086 :
00087 ddtScheme<Type>(mesh, is)
00088 {}
00089
00090
00091
00092
00093
00094 const fvMesh& mesh() const
00095 {
00096 return fv::ddtScheme<Type>::mesh();
00097 }
00098
00099 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00100 (
00101 const dimensioned<Type>&
00102 );
00103
00104 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00105 (
00106 const GeometricField<Type, fvPatchField, volMesh>&
00107 );
00108
00109 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00110 (
00111 const dimensionedScalar&,
00112 const GeometricField<Type, fvPatchField, volMesh>&
00113 );
00114
00115 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00116 (
00117 const volScalarField&,
00118 const GeometricField<Type, fvPatchField, volMesh>&
00119 );
00120
00121 tmp<fvMatrix<Type> > fvmDdt
00122 (
00123 GeometricField<Type, fvPatchField, volMesh>&
00124 );
00125
00126 tmp<fvMatrix<Type> > fvmDdt
00127 (
00128 const dimensionedScalar&,
00129 GeometricField<Type, fvPatchField, volMesh>&
00130 );
00131
00132 tmp<fvMatrix<Type> > fvmDdt
00133 (
00134 const volScalarField&,
00135 GeometricField<Type, fvPatchField, volMesh>&
00136 );
00137
00138 typedef typename ddtScheme<Type>::fluxFieldType fluxFieldType;
00139
00140 tmp<fluxFieldType> fvcDdtPhiCorr
00141 (
00142 const volScalarField& rA,
00143 const GeometricField<Type, fvPatchField, volMesh>& U,
00144 const fluxFieldType& phi
00145 );
00146
00147 tmp<fluxFieldType> fvcDdtPhiCorr
00148 (
00149 const volScalarField& rA,
00150 const volScalarField& rho,
00151 const GeometricField<Type, fvPatchField, volMesh>& U,
00152 const fluxFieldType& phi
00153 );
00154
00155 tmp<surfaceScalarField> meshPhi
00156 (
00157 const GeometricField<Type, fvPatchField, volMesh>&
00158 );
00159 };
00160
00161
00162 template<>
00163 tmp<surfaceScalarField> EulerDdtScheme<scalar>::fvcDdtPhiCorr
00164 (
00165 const volScalarField& rA,
00166 const volScalarField& U,
00167 const surfaceScalarField& phi
00168 );
00169
00170
00171 template<>
00172 tmp<surfaceScalarField> EulerDdtScheme<scalar>::fvcDdtPhiCorr
00173 (
00174 const volScalarField& rA,
00175 const volScalarField& rho,
00176 const volScalarField& U,
00177 const surfaceScalarField& phi
00178 );
00179
00180
00181
00182
00183 }
00184
00185
00186
00187 }
00188
00189
00190
00191 #ifdef NoRepository
00192 # include "EulerDdtScheme.C"
00193 #endif
00194
00195
00196
00197 #endif
00198
00199