steadyStateDdtScheme< Type > Class Template Reference

SteadyState implicit/explicit ddt which returns 0. More...

Inheritance diagram for steadyStateDdtScheme< Type >:
[legend]
Collaboration diagram for steadyStateDdtScheme< Type >:
[legend]

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 ("steadyState")
 Runtime type information. More...
 
 steadyStateDdtScheme (const fvMesh &mesh)
 Construct from mesh. More...
 
 steadyStateDdtScheme (const fvMesh &mesh, Istream &is)
 Construct from mesh and Istream. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
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 > &)
 
template<>
tmp< surfaceScalarFieldfvcDdtUfCorr (const GeometricField< scalar, fvPatchField, volMesh > &U, const GeometricField< scalar, fvsPatchField, surfaceMesh > &Uf)
 
template<>
tmp< surfaceScalarFieldfvcDdtPhiCorr (const volScalarField &U, const surfaceScalarField &phi)
 
template<>
tmp< surfaceScalarFieldfvcDdtUfCorr (const volScalarField &rho, const volScalarField &U, const surfaceScalarField &Uf)
 
template<>
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. 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 &is)
 Construct from mesh and Istream. More...
 
virtual ~ddtScheme ()=default
 Destructor. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
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. More...
 
int count () const noexcept
 Return the current reference count. More...
 
bool unique () const noexcept
 Return true if the reference count is zero. More...
 
void operator++ () noexcept
 Increment the reference count. More...
 
void operator++ (int) noexcept
 Increment the reference count. More...
 
void operator-- () noexcept
 Decrement the reference count. More...
 
void operator-- (int) noexcept
 Decrement the reference count. More...
 
- 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. More...
 
- Static Public Attributes inherited from ddtSchemeBase
static bool experimentalDdtCorr
 
- Protected Member Functions inherited from ddtScheme< Type >
 ddtScheme (const ddtScheme &)=delete
 No copy construct. More...
 
void operator= (const ddtScheme &)=delete
 No copy assignment. More...
 
- Protected Attributes inherited from ddtScheme< Type >
const fvMeshmesh_
 
scalar ddtPhiCoeff_
 Input for fvcDdtPhiCoeff. More...
 

Detailed Description

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

SteadyState implicit/explicit ddt which returns 0.

Source files

Definition at line 59 of file steadyStateDdtScheme.H.

Member Typedef Documentation

◆ fluxFieldType

Definition at line 154 of file steadyStateDdtScheme.H.

Constructor & Destructor Documentation

◆ steadyStateDdtScheme() [1/2]

steadyStateDdtScheme ( const fvMesh mesh)
inline

Construct from mesh.

Definition at line 81 of file steadyStateDdtScheme.H.

◆ steadyStateDdtScheme() [2/2]

steadyStateDdtScheme ( const fvMesh mesh,
Istream is 
)
inline

Construct from mesh and Istream.

Definition at line 87 of file steadyStateDdtScheme.H.

Member Function Documentation

◆ TypeName()

TypeName ( "steadyState"  )

Runtime type information.

◆ mesh()

const fvMesh& mesh ( ) const
inline

Return mesh reference.

Definition at line 96 of file steadyStateDdtScheme.H.

References ddtScheme< Type >::mesh().

Here is the call graph for this function:

◆ fvcDdt() [1/5]

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

Implements ddtScheme< Type >.

Definition at line 47 of file steadyStateDdtScheme.C.

References dimensioned< Type >::dimensions(), Foam::dimTime, mesh, dimensioned< Type >::name(), Foam::New(), and Foam::Zero.

Here is the call graph for this function:

◆ fvcDdt() [2/5]

tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt ( const GeometricField< Type, fvPatchField, volMesh > &  vf)
virtual

Implements ddtScheme< Type >.

Definition at line 68 of file steadyStateDdtScheme.C.

References Foam::dimTime, mesh, Foam::New(), and Foam::Zero.

Here is the call graph for this function:

◆ fvcDdt() [3/5]

tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt ( const dimensionedScalar rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)
virtual

Implements ddtScheme< Type >.

Definition at line 89 of file steadyStateDdtScheme.C.

References Foam::dimTime, mesh, Foam::New(), rho, and Foam::Zero.

Here is the call graph for this function:

◆ fvcDdt() [4/5]

tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt ( const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)
virtual

Implements ddtScheme< Type >.

Definition at line 111 of file steadyStateDdtScheme.C.

References Foam::dimTime, mesh, Foam::New(), rho, and Foam::Zero.

Here is the call graph for this function:

◆ fvcDdt() [5/5]

tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt ( const volScalarField alpha,
const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  psi 
)
virtual

Implements ddtScheme< Type >.

Definition at line 133 of file steadyStateDdtScheme.C.

References Foam::constant::atomic::alpha, Foam::dimTime, mesh, dimensioned< Type >::name(), Foam::New(), rho, and Foam::Zero.

Here is the call graph for this function:

◆ fvmDdt() [1/4]

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

Implements ddtScheme< Type >.

Definition at line 156 of file steadyStateDdtScheme.C.

References Foam::dimTime, and Foam::dimVol.

◆ fvmDdt() [2/4]

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

Implements ddtScheme< Type >.

Definition at line 176 of file steadyStateDdtScheme.C.

References Foam::dimTime, Foam::dimVol, and rho.

◆ fvmDdt() [3/4]

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

Implements ddtScheme< Type >.

Definition at line 197 of file steadyStateDdtScheme.C.

References Foam::dimTime, Foam::dimVol, and rho.

◆ fvmDdt() [4/4]

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

Implements ddtScheme< Type >.

Definition at line 218 of file steadyStateDdtScheme.C.

References Foam::constant::atomic::alpha, dimensioned< Type >::dimensions(), Foam::dimTime, Foam::dimVol, and rho.

Here is the call graph for this function:

◆ fvcDdtUfCorr() [1/4]

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

Implements ddtScheme< Type >.

Definition at line 240 of file steadyStateDdtScheme.C.

References Foam::dimArea, Foam::dimTime, mesh, tmp< T >::ref(), U, Uf, and Foam::Zero.

Here is the call graph for this function:

◆ fvcDdtPhiCorr() [1/4]

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

Implements ddtScheme< Type >.

Definition at line 272 of file steadyStateDdtScheme.C.

References Foam::dimTime, mesh, phi, tmp< T >::ref(), U, and Foam::Zero.

Here is the call graph for this function:

◆ fvcDdtUfCorr() [2/4]

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

Implements ddtScheme< Type >.

Definition at line 304 of file steadyStateDdtScheme.C.

References Foam::dimArea, Foam::dimTime, mesh, tmp< T >::ref(), rho, U, Uf, and Foam::Zero.

Here is the call graph for this function:

◆ fvcDdtPhiCorr() [2/4]

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

Implements ddtScheme< Type >.

Definition at line 339 of file steadyStateDdtScheme.C.

References Foam::dimTime, mesh, phi, tmp< T >::ref(), rho, U, and Foam::Zero.

Here is the call graph for this function:

◆ meshPhi()

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

Implements ddtScheme< Type >.

Definition at line 373 of file steadyStateDdtScheme.C.

References Foam::dimTime, Foam::dimVolume, mesh, tmp< T >::New(), IOobject::NO_READ, IOobject::NO_WRITE, timeName, and Foam::Zero.

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 
)

◆ fvcDdtPhiCorr() [3/4]

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

◆ fvcDdtUfCorr() [4/4]

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

◆ fvcDdtPhiCorr() [4/4]

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

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