Second-order backward-differencing ddt using the current and two previous time-step values. More...
Public Types | |
typedef ddtScheme< Type >::fluxFieldType | fluxFieldType |
Public Types inherited from ddtScheme< Type > | |
typedef SurfaceField< typename flux< Type >::type > | fluxFieldType |
Public Member Functions | |
TypeName ("backward") | |
Runtime type information. More... | |
backwardDdtScheme (const fvMesh &mesh) | |
Construct from mesh. More... | |
backwardDdtScheme (const fvMesh &mesh, Istream &is) | |
Construct from mesh and Istream. More... | |
backwardDdtScheme (const backwardDdtScheme &)=delete | |
Disallow default bitwise copy construction. More... | |
const fvMesh & | mesh () const |
Return mesh reference. More... | |
virtual tmp< VolField< Type > > | fvcDdt (const dimensioned< Type > &) |
virtual tmp< VolField< Type > > | fvcDdt (const VolField< Type > &) |
virtual tmp< VolField< Type > > | fvcDdt (const dimensionedScalar &, const VolField< Type > &) |
virtual tmp< VolField< Type > > | fvcDdt (const volScalarField &, const VolField< Type > &) |
virtual tmp< VolField< Type > > | fvcDdt (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &psi) |
virtual tmp< fvMatrix< Type > > | fvmDdt (const VolField< Type > &) |
virtual tmp< fvMatrix< Type > > | fvmDdt (const dimensionedScalar &, const VolField< Type > &) |
virtual tmp< fvMatrix< Type > > | fvmDdt (const volScalarField &, const VolField< Type > &) |
virtual tmp< fvMatrix< Type > > | fvmDdt (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &psi) |
virtual tmp< fluxFieldType > | fvcDdtUfCorr (const VolField< Type > &U, const SurfaceField< Type > &Uf) |
virtual tmp< fluxFieldType > | fvcDdtPhiCorr (const VolField< Type > &U, const fluxFieldType &phi) |
virtual tmp< fluxFieldType > | fvcDdtUfCorr (const volScalarField &rho, const VolField< Type > &U, const SurfaceField< Type > &rhoUf) |
virtual tmp< fluxFieldType > | fvcDdtPhiCorr (const volScalarField &rho, const VolField< Type > &U, const fluxFieldType &phi) |
virtual tmp< fluxFieldType > | fvcDdtUfCorr (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &U, const SurfaceField< Type > &Uf) |
virtual tmp< fluxFieldType > | fvcDdtPhiCorr (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &U, const fluxFieldType &phi) |
virtual tmp< surfaceScalarField > | meshPhi (const VolField< Type > &) |
virtual tmp< scalarField > | meshPhi (const VolField< Type > &, const label patchi) |
void | operator= (const backwardDdtScheme &)=delete |
Disallow default bitwise assignment. More... | |
tmp< surfaceScalarField > | fvcDdtUfCorr (const VolField< scalar > &U, const SurfaceField< scalar > &Uf) |
tmp< surfaceScalarField > | fvcDdtPhiCorr (const volScalarField &U, const surfaceScalarField &phi) |
tmp< surfaceScalarField > | fvcDdtUfCorr (const volScalarField &rho, const volScalarField &U, const surfaceScalarField &rhoUf) |
tmp< surfaceScalarField > | fvcDdtPhiCorr (const volScalarField &rho, const volScalarField &U, const surfaceScalarField &phi) |
tmp< surfaceScalarField > | fvcDdtUfCorr (const volScalarField &alpha, const volScalarField &rho, const volScalarField &U, const surfaceScalarField &Uf) |
tmp< surfaceScalarField > | fvcDdtPhiCorr (const volScalarField &alpha, const volScalarField &rho, const volScalarField &U, const surfaceScalarField &phi) |
tmp< surfaceScalarField > | fvcDdtUfCorr (const volScalarField &alpha, const volScalarField &rho, const volScalarField &U, const surfaceScalarField &Uf) |
tmp< surfaceScalarField > | fvcDdtPhiCorr (const volScalarField &alpha, const volScalarField &rho, const volScalarField &U, const surfaceScalarField &phi) |
Public Member Functions inherited from ddtScheme< Type > | |
virtual const word & | type () const =0 |
Runtime type information. More... | |
declareRunTimeSelectionTable (tmp, ddtScheme, Istream,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData)) | |
ddtScheme (const fvMesh &mesh) | |
Construct from mesh. More... | |
ddtScheme (const fvMesh &mesh, Istream &) | |
Construct from mesh and Istream. More... | |
ddtScheme (const ddtScheme &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~ddtScheme () |
Destructor. More... | |
const fvMesh & | mesh () const |
Return mesh reference. More... | |
virtual tmp< SurfaceField< Type > > | fvcDdt (const SurfaceField< Type > &) |
virtual tmp< surfaceScalarField > | fvcDdtPhiCoeff (const VolField< Type > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr) |
virtual tmp< surfaceScalarField > | fvcDdtPhiCoeff (const VolField< Type > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr, const volScalarField &rho) |
virtual tmp< surfaceScalarField > | fvcDdtPhiCoeff (const VolField< Type > &U, const fluxFieldType &phi) |
virtual tmp< surfaceScalarField > | fvcDdtPhiCoeff (const VolField< Type > &U, const fluxFieldType &phi, const volScalarField &rho) |
void | operator= (const ddtScheme &)=delete |
Disallow default bitwise assignment. More... | |
Public Member Functions inherited from refCount | |
int | count () const |
Return the current reference count. More... | |
bool | unique () const |
Return true if the reference count is zero. More... | |
void | operator++ () |
Increment the reference count. More... | |
void | operator++ (int) |
Increment the reference count. More... | |
void | operator-- () |
Decrement the reference count. More... | |
void | operator-- (int) |
Decrement the reference count. More... | |
void | operator= (const refCount &)=delete |
Disallow bitwise assignment. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from ddtScheme< Type > | |
static tmp< ddtScheme< Type > > | New (const fvMesh &mesh, Istream &schemeData) |
Return a pointer to a new ddtScheme created on freestore. More... | |
Protected Member Functions inherited from refCount | |
refCount () | |
Construct null initialising count to 0. More... | |
refCount (const refCount &)=delete | |
Disallow copy. More... | |
Protected Attributes inherited from ddtScheme< Type > | |
const fvMesh & | mesh_ |
Second-order backward-differencing ddt using the current and two previous time-step values.
Definition at line 55 of file backwardDdtScheme.H.
typedef ddtScheme<Type>::fluxFieldType fluxFieldType |
Definition at line 172 of file backwardDdtScheme.H.
|
inline |
Construct from mesh.
Definition at line 82 of file backwardDdtScheme.H.
References backwardDdtScheme< Type >::mesh(), polyMesh::moving(), and fvMesh::V00().
|
inline |
Construct from mesh and Istream.
Definition at line 95 of file backwardDdtScheme.H.
References backwardDdtScheme< Type >::mesh(), polyMesh::moving(), and fvMesh::V00().
|
delete |
Disallow default bitwise copy construction.
TypeName | ( | "backward" | ) |
Runtime type information.
|
inline |
Return mesh reference.
Definition at line 114 of file backwardDdtScheme.H.
References ddtScheme< Type >::mesh().
Referenced by backwardDdtScheme< Type >::backwardDdtScheme().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 76 of file backwardDdtScheme.C.
References dimensioned< Type >::dimensions(), Foam::dimTime, dimensioned< Type >::name(), GeometricField< Type, PatchField, GeoMesh >::New(), tmp< T >::ref(), dimensioned< Type >::value(), and Foam::Zero.
Implements ddtScheme< Type >.
Definition at line 139 of file backwardDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::New(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), and dimensioned< Type >::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 203 of file backwardDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::New(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), rho, and dimensioned< Type >::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 268 of file backwardDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::New(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), rho, and dimensioned< Type >::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 336 of file backwardDdtScheme.C.
References alpha(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::New(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), rho, and dimensioned< Type >::value().
Implements ddtScheme< Type >.
Definition at line 416 of file backwardDdtScheme.C.
References lduMatrix::diag(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, Foam::dimVol, GeometricField< Type, PatchField, GeoMesh >::oldTime(), tmp< T >::ref(), and fvMatrix< Type >::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 467 of file backwardDdtScheme.C.
References lduMatrix::diag(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, Foam::dimVol, GeometricField< Type, PatchField, GeoMesh >::oldTime(), tmp< T >::ref(), rho, and fvMatrix< Type >::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 518 of file backwardDdtScheme.C.
References lduMatrix::diag(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, Foam::dimVol, GeometricField< Type, PatchField, GeoMesh >::oldTime(), tmp< T >::ref(), rho, and fvMatrix< Type >::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 572 of file backwardDdtScheme.C.
References alpha(), lduMatrix::diag(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, Foam::dimVol, GeometricField< Type, PatchField, GeoMesh >::oldTime(), tmp< T >::ref(), rho, and fvMatrix< Type >::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 638 of file backwardDdtScheme.C.
References Foam::fvc::interpolate(), IOobject::name(), Foam::compressible::New(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 674 of file backwardDdtScheme.C.
References Foam::fvc::dotInterpolate(), IOobject::name(), Foam::compressible::New(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 708 of file backwardDdtScheme.C.
References Foam::abort(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimVelocity, Foam::FatalError, FatalErrorInFunction, Foam::fvc::interpolate(), IOobject::name(), Foam::compressible::New(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), rho, and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 807 of file backwardDdtScheme.C.
References Foam::abort(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimFlux, Foam::dimVelocity, Foam::fvc::dotInterpolate(), Foam::FatalError, FatalErrorInFunction, IOobject::name(), Foam::compressible::New(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), rho, and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 888 of file backwardDdtScheme.C.
References Foam::abort(), alpha(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimVelocity, Foam::FatalError, FatalErrorInFunction, Foam::fvc::interpolate(), IOobject::name(), Foam::compressible::New(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), rho, and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 950 of file backwardDdtScheme.C.
References Foam::abort(), alpha(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimFlux, Foam::dimVelocity, Foam::fvc::dotInterpolate(), Foam::FatalError, FatalErrorInFunction, Foam::fvc::interpolate(), IOobject::name(), Foam::compressible::New(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), rho, and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1009 of file backwardDdtScheme.C.
References Foam::name(), and GeometricField< Type, PatchField, GeoMesh >::New().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1032 of file backwardDdtScheme.C.
References patchi.
|
delete |
Disallow default bitwise assignment.
tmp< surfaceScalarField > fvcDdtUfCorr | ( | const VolField< scalar > & | U, |
const SurfaceField< scalar > & | Uf | ||
) |
tmp< surfaceScalarField > fvcDdtPhiCorr | ( | const volScalarField & | U, |
const surfaceScalarField & | phi | ||
) |
tmp< surfaceScalarField > fvcDdtUfCorr | ( | const volScalarField & | rho, |
const volScalarField & | U, | ||
const surfaceScalarField & | rhoUf | ||
) |
tmp< surfaceScalarField > fvcDdtPhiCorr | ( | const volScalarField & | rho, |
const volScalarField & | U, | ||
const surfaceScalarField & | phi | ||
) |
tmp< surfaceScalarField > fvcDdtUfCorr | ( | const volScalarField & | alpha, |
const volScalarField & | rho, | ||
const volScalarField & | U, | ||
const surfaceScalarField & | Uf | ||
) |
tmp< surfaceScalarField > fvcDdtPhiCorr | ( | const volScalarField & | alpha, |
const volScalarField & | rho, | ||
const volScalarField & | U, | ||
const surfaceScalarField & | phi | ||
) |
tmp< surfaceScalarField > fvcDdtUfCorr | ( | const volScalarField & | alpha, |
const volScalarField & | rho, | ||
const volScalarField & | U, | ||
const surfaceScalarField & | Uf | ||
) |
tmp< surfaceScalarField > fvcDdtPhiCorr | ( | const volScalarField & | alpha, |
const volScalarField & | rho, | ||
const volScalarField & | U, | ||
const surfaceScalarField & | phi | ||
) |