42#ifndef backwardDdtScheme_H
43#define backwardDdtScheme_H
69 scalar deltaT_()
const;
72 scalar deltaT0_()
const;
76 template<
class GeoField>
77 scalar deltaT0_(
const GeoField&)
const;
Generic GeometricField class.
bool good() const noexcept
True if next operation might succeed.
bool eof() const noexcept
True if end of input seen.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
Second-order backward-differencing ddt using the current and two previous time-step values.
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
TypeName("backward")
Runtime type information.
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
backwardDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
const fvMesh & mesh() const
Return mesh reference.
backwardDdtScheme(const fvMesh &mesh)
Construct from mesh.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
Abstract base class for ddt schemes.
scalar ddtPhiCoeff_
Input for fvcDdtPhiCoeff.
const fvMesh & mesh() const
Return mesh reference.
bool moving() const noexcept
Is mesh moving.
A class for managing temporary objects.
const volScalarField & psi
autoPtr< surfaceVectorField > Uf
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.