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 #ifndef CrankNicholsonDdtScheme_H
00037 #define CrankNicholsonDdtScheme_H
00038
00039 #include <finiteVolume/ddtScheme.H>
00040
00041
00042
00043 namespace Foam
00044 {
00045
00046
00047
00048 namespace fv
00049 {
00050
00051
00052
00053
00054
00055 template<class Type>
00056 class CrankNicholsonDdtScheme
00057 :
00058 public fv::ddtScheme<Type>
00059 {
00060
00061
00062
00063
00064
00065 template<class GeoField>
00066 class DDt0Field
00067 :
00068 public GeoField
00069 {
00070 label startTimeIndex_;
00071
00072 public:
00073
00074
00075 DDt0Field
00076 (
00077 const IOobject& io,
00078 const fvMesh& mesh
00079 );
00080
00081
00082
00083 DDt0Field
00084 (
00085 const IOobject& io,
00086 const fvMesh& mesh,
00087 const dimensioned<typename GeoField::value_type>& dimType
00088 );
00089
00090
00091 label startTimeIndex() const;
00092
00093
00094 GeoField& operator()();
00095
00096
00097 void operator=(const GeoField& gf);
00098 };
00099
00100
00101
00102 scalar ocCoeff_;
00103
00104
00105
00106
00107
00108 CrankNicholsonDdtScheme(const CrankNicholsonDdtScheme&);
00109
00110
00111 void operator=(const CrankNicholsonDdtScheme&);
00112
00113 template<class GeoField>
00114 DDt0Field<GeoField>& ddt0_
00115 (
00116 const word& name,
00117 const dimensionSet& dims
00118 );
00119
00120
00121 template<class GeoField>
00122 bool evaluate(const DDt0Field<GeoField>& ddt0) const;
00123
00124
00125
00126 template<class GeoField>
00127 scalar coef_(const DDt0Field<GeoField>&) const;
00128
00129
00130
00131 template<class GeoField>
00132 scalar coef0_(const DDt0Field<GeoField>&) const;
00133
00134
00135
00136 template<class GeoField>
00137 dimensionedScalar rDtCoef_(const DDt0Field<GeoField>&) const;
00138
00139
00140
00141 template<class GeoField>
00142 dimensionedScalar rDtCoef0_(const DDt0Field<GeoField>&) const;
00143
00144
00145 template<class GeoField>
00146 tmp<GeoField> offCentre_(const GeoField& ddt0) const;
00147
00148
00149 public:
00150
00151
00152 TypeName("CrankNicholson");
00153
00154
00155
00156
00157
00158 CrankNicholsonDdtScheme(const fvMesh& mesh)
00159 :
00160 ddtScheme<Type>(mesh),
00161 ocCoeff_(1.0)
00162 {}
00163
00164
00165 CrankNicholsonDdtScheme(const fvMesh& mesh, Istream& is)
00166 :
00167 ddtScheme<Type>(mesh, is),
00168 ocCoeff_(readScalar(is))
00169 {
00170 if (ocCoeff_ < 0 || ocCoeff_ > 1)
00171 {
00172 FatalIOErrorIn
00173 (
00174 "CrankNicholsonDdtScheme(const fvMesh& mesh, Istream& is)",
00175 is
00176 ) << "coefficient = " << ocCoeff_
00177 << " should be >= 0 and <= 1"
00178 << exit(FatalIOError);
00179 }
00180 }
00181
00182
00183
00184
00185
00186 const fvMesh& mesh() const
00187 {
00188 return fv::ddtScheme<Type>::mesh();
00189 }
00190
00191 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00192 (
00193 const dimensioned<Type>&
00194 );
00195
00196 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00197 (
00198 const GeometricField<Type, fvPatchField, volMesh>&
00199 );
00200
00201 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00202 (
00203 const dimensionedScalar&,
00204 const GeometricField<Type, fvPatchField, volMesh>&
00205 );
00206
00207 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00208 (
00209 const volScalarField&,
00210 const GeometricField<Type, fvPatchField, volMesh>&
00211 );
00212
00213 tmp<fvMatrix<Type> > fvmDdt
00214 (
00215 GeometricField<Type, fvPatchField, volMesh>&
00216 );
00217
00218 tmp<fvMatrix<Type> > fvmDdt
00219 (
00220 const dimensionedScalar&,
00221 GeometricField<Type, fvPatchField, volMesh>&
00222 );
00223
00224 tmp<fvMatrix<Type> > fvmDdt
00225 (
00226 const volScalarField&,
00227 GeometricField<Type, fvPatchField, volMesh>&
00228 );
00229
00230 typedef typename ddtScheme<Type>::fluxFieldType fluxFieldType;
00231
00232 tmp<fluxFieldType> fvcDdtPhiCorr
00233 (
00234 const volScalarField& rA,
00235 const GeometricField<Type, fvPatchField, volMesh>& U,
00236 const fluxFieldType& phi
00237 );
00238
00239 tmp<fluxFieldType> fvcDdtPhiCorr
00240 (
00241 const volScalarField& rA,
00242 const volScalarField& rho,
00243 const GeometricField<Type, fvPatchField, volMesh>& U,
00244 const fluxFieldType& phi
00245 );
00246
00247
00248 tmp<surfaceScalarField> meshPhi
00249 (
00250 const GeometricField<Type, fvPatchField, volMesh>&
00251 );
00252 };
00253
00254
00255 template<>
00256 tmp<surfaceScalarField> CrankNicholsonDdtScheme<scalar>::fvcDdtPhiCorr
00257 (
00258 const volScalarField& rA,
00259 const volScalarField& U,
00260 const surfaceScalarField& phi
00261 );
00262
00263
00264 template<>
00265 tmp<surfaceScalarField> CrankNicholsonDdtScheme<scalar>::fvcDdtPhiCorr
00266 (
00267 const volScalarField& rA,
00268 const volScalarField& rho,
00269 const volScalarField& U,
00270 const surfaceScalarField& phi
00271 );
00272
00273
00274
00275
00276 }
00277
00278
00279
00280 }
00281
00282
00283
00284 #ifdef NoRepository
00285 # include <finiteVolume/CrankNicholsonDdtScheme.C>
00286 #endif
00287
00288
00289
00290 #endif
00291
00292