Loading...
Searching...
No Matches
backwardDdtScheme< Type > Class Template Reference

Second-order backward-differencing ddt using the current and two previous time-step values. More...

#include <backwardDdtScheme.H>

Inheritance diagram for backwardDdtScheme< Type >:
Collaboration diagram for backwardDdtScheme< Type >:

Public Types

typedef ddtScheme< Type >::fluxFieldType fluxFieldType
Public Types inherited from ddtScheme< Type >
typedef GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMeshfluxFieldType

Public Member Functions

 TypeName ("backward")
 Runtime type information.
 backwardDdtScheme (const fvMesh &mesh)
 Construct from mesh.
 backwardDdtScheme (const fvMesh &mesh, Istream &is)
 Construct from mesh and Istream.
const fvMeshmesh () const
 Return mesh reference.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const GeometricField< Type, fvPatchField, volMesh > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const volScalarField &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &psi)
tmp< fvMatrix< Type > > fvmDdt (const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fvMatrix< Type > > fvmDdt (const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fvMatrix< Type > > fvmDdt (const volScalarField &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fvMatrix< Type > > fvmDdt (const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &psi)
tmp< fluxFieldTypefvcDdtUfCorr (const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
tmp< fluxFieldTypefvcDdtPhiCorr (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
tmp< fluxFieldTypefvcDdtUfCorr (const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
tmp< fluxFieldTypefvcDdtPhiCorr (const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
tmp< surfaceScalarFieldmeshPhi (const GeometricField< Type, fvPatchField, volMesh > &)
tmp< surfaceScalarFieldfvcDdtUfCorr (const GeometricField< scalar, fvPatchField, volMesh > &U, const GeometricField< scalar, fvsPatchField, surfaceMesh > &Uf)
tmp< surfaceScalarFieldfvcDdtPhiCorr (const volScalarField &U, const surfaceScalarField &phi)
tmp< surfaceScalarFieldfvcDdtUfCorr (const volScalarField &rho, const volScalarField &U, const surfaceScalarField &Uf)
tmp< surfaceScalarFieldfvcDdtPhiCorr (const volScalarField &rho, const volScalarField &U, const surfaceScalarField &phi)
Public Member Functions inherited from ddtScheme< Type >
virtual const wordtype () const =0
 Runtime type information.
 declareRunTimeSelectionTable (tmp, ddtScheme, Istream,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
 ddtScheme (const fvMesh &mesh)
 Construct from mesh.
 ddtScheme (const fvMesh &mesh, Istream &is)
 Construct from mesh and Istream.
virtual ~ddtScheme ()=default
 Destructor.
const fvMeshmesh () const
 Return mesh reference.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > fvcDdt (const GeometricField< Type, fvsPatchField, surfaceMesh > &)
tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
tmp< surfaceScalarFieldfvcDdtPhiCoeffExperimental (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr, const volScalarField &rho)
tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &rhoU, const fluxFieldType &phi, const volScalarField &rho)
Public Member Functions inherited from refCount
constexpr refCount () noexcept
 Default construct, initializing count to 0.
int use_count () const noexcept
 Return the current reference count.
bool unique () const noexcept
 Return true if the reference count is zero.
void operator++ () noexcept
 Increment the reference count.
void operator++ (int) noexcept
 Increment the reference count.
void operator-- () noexcept
 Decrement the reference count.
void operator-- (int) noexcept
 Decrement the reference count.
Public Member Functions inherited from ddtSchemeBase
 ddtSchemeBase ()
 ddtSchemeBase ()

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.
Static Public Attributes inherited from ddtSchemeBase
static bool experimentalDdtCorr
 Flag to use experimental ddtCorr from org version Default is off for backwards compatibility.
Protected Member Functions inherited from ddtScheme< Type >
 ddtScheme (const ddtScheme &)=delete
 No copy construct.
void operator= (const ddtScheme &)=delete
 No copy assignment.
Protected Attributes inherited from ddtScheme< Type >
const fvMeshmesh_
scalar ddtPhiCoeff_
 Input for fvcDdtPhiCoeff.

Detailed Description

template<class Type>
class Foam::fv::backwardDdtScheme< Type >

Second-order backward-differencing ddt using the current and two previous time-step values.

Source files

Definition at line 57 of file backwardDdtScheme.H.

Member Typedef Documentation

◆ fluxFieldType

template<class Type>
typedef ddtScheme<Type>::fluxFieldType fluxFieldType

Definition at line 194 of file backwardDdtScheme.H.

Constructor & Destructor Documentation

◆ backwardDdtScheme() [1/2]

template<class Type>
backwardDdtScheme ( const fvMesh & mesh)
inline

Construct from mesh.

Definition at line 105 of file backwardDdtScheme.H.

References ddtScheme< Type >::ddtScheme(), and mesh().

Here is the call graph for this function:

◆ backwardDdtScheme() [2/2]

template<class Type>
backwardDdtScheme ( const fvMesh & mesh,
Istream & is )
inline

Construct from mesh and Istream.

Definition at line 113 of file backwardDdtScheme.H.

References ddtScheme< Type >::ddtScheme(), IOstream::eof(), IOstream::good(), if(), and mesh().

Here is the call graph for this function:

Member Function Documentation

◆ TypeName()

template<class Type>
TypeName ( "backward" )

Runtime type information.

◆ mesh()

template<class Type>
const fvMesh & mesh ( ) const
inline

Return mesh reference.

Definition at line 136 of file backwardDdtScheme.H.

References ddtScheme< Type >::mesh().

Referenced by backwardDdtScheme(), backwardDdtScheme(), fvcDdt(), fvcDdt(), fvcDdt(), fvcDdt(), fvcDdt(), fvcDdtPhiCorr(), fvcDdtPhiCorr(), fvcDdtUfCorr(), fvcDdtUfCorr(), fvmDdt(), fvmDdt(), fvmDdt(), and fvmDdt().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ fvcDdt() [1/5]

template<class Type>
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt ( const dimensioned< Type > & dt)
virtual

◆ fvcDdt() [2/5]

◆ fvcDdt() [3/5]

◆ fvcDdt() [4/5]

◆ fvcDdt() [5/5]

◆ fvmDdt() [1/4]

template<class Type>
tmp< fvMatrix< Type > > fvmDdt ( const GeometricField< Type, fvPatchField, volMesh > & vf)
virtual

◆ fvmDdt() [2/4]

template<class Type>
tmp< fvMatrix< Type > > fvmDdt ( const dimensionedScalar & rho,
const GeometricField< Type, fvPatchField, volMesh > & vf )
virtual

◆ fvmDdt() [3/4]

template<class Type>
tmp< fvMatrix< Type > > fvmDdt ( const volScalarField & rho,
const GeometricField< Type, fvPatchField, volMesh > & vf )
virtual

◆ fvmDdt() [4/4]

template<class Type>
tmp< fvMatrix< Type > > fvmDdt ( const volScalarField & alpha,
const volScalarField & rho,
const GeometricField< Type, fvPatchField, volMesh > & psi )
virtual

◆ fvcDdtUfCorr() [1/4]

template<class Type>
tmp< typename backwardDdtScheme< Type >::fluxFieldType > fvcDdtUfCorr ( const GeometricField< Type, fvPatchField, volMesh > & U,
const GeometricField< Type, fvsPatchField, surfaceMesh > & Uf )
virtual

Implements ddtScheme< Type >.

Definition at line 687 of file backwardDdtScheme.C.

References ddtScheme< Type >::fvcDdtPhiCoeff(), Foam::fvc::interpolate(), mesh(), mesh, timeName, U, and Uf.

Here is the call graph for this function:

◆ fvcDdtPhiCorr() [1/4]

template<class Type>
tmp< typename backwardDdtScheme< Type >::fluxFieldType > fvcDdtPhiCorr ( const GeometricField< Type, fvPatchField, volMesh > & U,
const fluxFieldType & phi )
virtual

Implements ddtScheme< Type >.

Definition at line 731 of file backwardDdtScheme.C.

References Foam::fvc::dotInterpolate(), ddtScheme< Type >::fvcDdtPhiCoeff(), mesh(), mesh, phi, timeName, and U.

Here is the call graph for this function:

◆ fvcDdtUfCorr() [2/4]

template<class Type>
tmp< typename backwardDdtScheme< Type >::fluxFieldType > fvcDdtUfCorr ( const volScalarField & rho,
const GeometricField< Type, fvPatchField, volMesh > & U,
const GeometricField< Type, fvsPatchField, surfaceMesh > & Uf )
virtual

◆ fvcDdtPhiCorr() [2/4]

template<class Type>
tmp< typename backwardDdtScheme< Type >::fluxFieldType > fvcDdtPhiCorr ( const volScalarField & rho,
const GeometricField< Type, fvPatchField, volMesh > & U,
const fluxFieldType & phi )
virtual

◆ meshPhi()

template<class Type>
tmp< surfaceScalarField > meshPhi ( const GeometricField< Type, fvPatchField, volMesh > & vf)
virtual

Implements ddtScheme< Type >.

Definition at line 982 of file backwardDdtScheme.C.

References mesh, name, IOobjectOption::NO_READ, IOobjectOption::NO_REGISTER, IOobjectOption::NO_WRITE, oldTime(), phi, and timeName.

Here is the call graph for this function:

◆ fvcDdtUfCorr() [3/4]

tmp< surfaceScalarField > fvcDdtUfCorr ( const GeometricField< scalar, fvPatchField, volMesh > & U,
const GeometricField< scalar, fvsPatchField, surfaceMesh > & Uf )

References U, and Uf.

◆ fvcDdtPhiCorr() [3/4]

tmp< surfaceScalarField > fvcDdtPhiCorr ( const volScalarField & U,
const surfaceScalarField & phi )

References phi, and U.

◆ fvcDdtUfCorr() [4/4]

tmp< surfaceScalarField > fvcDdtUfCorr ( const volScalarField & rho,
const volScalarField & U,
const surfaceScalarField & Uf )

References rho, U, and Uf.

◆ fvcDdtPhiCorr() [4/4]

tmp< surfaceScalarField > fvcDdtPhiCorr ( const volScalarField & rho,
const volScalarField & U,
const surfaceScalarField & phi )

References phi, rho, and U.


The documentation for this class was generated from the following files: