Go to the documentation of this file.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 backwardDdtScheme_H
00037 #define backwardDdtScheme_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 backwardDdtScheme
00057 :
00058 public fv::ddtScheme<Type>
00059 {
00060
00061
00062
00063 scalar deltaT_() const;
00064
00065
00066 scalar deltaT0_() const;
00067
00068
00069
00070 template<class GeoField>
00071 scalar deltaT0_(const GeoField&) const;
00072
00073
00074 backwardDdtScheme(const backwardDdtScheme&);
00075
00076
00077 void operator=(const backwardDdtScheme&);
00078
00079
00080 public:
00081
00082
00083 TypeName("backward");
00084
00085
00086
00087
00088
00089 backwardDdtScheme(const fvMesh& mesh)
00090 :
00091 ddtScheme<Type>(mesh)
00092 {}
00093
00094
00095 backwardDdtScheme(const fvMesh& mesh, Istream& is)
00096 :
00097 ddtScheme<Type>(mesh, is)
00098 {}
00099
00100
00101
00102
00103
00104 const fvMesh& mesh() const
00105 {
00106 return fv::ddtScheme<Type>::mesh();
00107 }
00108
00109 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00110 (
00111 const dimensioned<Type>&
00112 );
00113
00114 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00115 (
00116 const GeometricField<Type, fvPatchField, volMesh>&
00117 );
00118
00119 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00120 (
00121 const dimensionedScalar&,
00122 const GeometricField<Type, fvPatchField, volMesh>&
00123 );
00124
00125 tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00126 (
00127 const volScalarField&,
00128 const GeometricField<Type, fvPatchField, volMesh>&
00129 );
00130
00131 tmp<fvMatrix<Type> > fvmDdt
00132 (
00133 GeometricField<Type, fvPatchField, volMesh>&
00134 );
00135
00136 tmp<fvMatrix<Type> > fvmDdt
00137 (
00138 const dimensionedScalar&,
00139 GeometricField<Type, fvPatchField, volMesh>&
00140 );
00141
00142 tmp<fvMatrix<Type> > fvmDdt
00143 (
00144 const volScalarField&,
00145 GeometricField<Type, fvPatchField, volMesh>&
00146 );
00147
00148 typedef typename ddtScheme<Type>::fluxFieldType fluxFieldType;
00149
00150 tmp<fluxFieldType> fvcDdtPhiCorr
00151 (
00152 const volScalarField& rA,
00153 const GeometricField<Type, fvPatchField, volMesh>& U,
00154 const fluxFieldType& phi
00155 );
00156
00157 tmp<fluxFieldType> fvcDdtPhiCorr
00158 (
00159 const volScalarField& rA,
00160 const volScalarField& rho,
00161 const GeometricField<Type, fvPatchField, volMesh>& U,
00162 const fluxFieldType& phi
00163 );
00164
00165
00166 tmp<surfaceScalarField> meshPhi
00167 (
00168 const GeometricField<Type, fvPatchField, volMesh>&
00169 );
00170 };
00171
00172
00173 template<>
00174 tmp<surfaceScalarField> backwardDdtScheme<scalar>::fvcDdtPhiCorr
00175 (
00176 const volScalarField& rA,
00177 const volScalarField& U,
00178 const surfaceScalarField& phi
00179 );
00180
00181
00182 template<>
00183 tmp<surfaceScalarField> backwardDdtScheme<scalar>::fvcDdtPhiCorr
00184 (
00185 const volScalarField& rA,
00186 const volScalarField& rho,
00187 const volScalarField& U,
00188 const surfaceScalarField& phi
00189 );
00190
00191
00192
00193
00194 }
00195
00196
00197
00198 }
00199
00200
00201
00202 #ifdef NoRepository
00203 # include <finiteVolume/backwardDdtScheme.C>
00204 #endif
00205
00206
00207
00208 #endif
00209
00210