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