ddtScheme< Type > Class Template Referenceabstract

Abstract base class for ddt schemes. More...

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

Public Types

typedef GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMeshfluxFieldType
 

Public Member Functions

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, fvPatchField, volMesh > > fvcDdt (const dimensioned< Type > &)=0
 
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const GeometricField< Type, fvPatchField, volMesh > &)=0
 
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &)=0
 
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const volScalarField &, const GeometricField< Type, fvPatchField, volMesh > &)=0
 
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &)=0
 
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > fvcDdt (const GeometricField< Type, fvsPatchField, surfaceMesh > &)
 
virtual tmp< fvMatrix< Type > > fvmDdt (const GeometricField< Type, fvPatchField, volMesh > &)=0
 
virtual tmp< fvMatrix< Type > > fvmDdt (const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &)=0
 
virtual tmp< fvMatrix< Type > > fvmDdt (const volScalarField &, const GeometricField< Type, fvPatchField, volMesh > &)=0
 
virtual tmp< fvMatrix< Type > > fvmDdt (const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &vf)=0
 
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)
 
virtual tmp< fluxFieldTypefvcDdtUfCorr (const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)=0
 
virtual tmp< fluxFieldTypefvcDdtPhiCorr (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)=0
 
virtual tmp< fluxFieldTypefvcDdtUfCorr (const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)=0
 
virtual tmp< fluxFieldTypefvcDdtPhiCorr (const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)=0
 
virtual tmp< surfaceScalarFieldmeshPhi (const GeometricField< Type, fvPatchField, volMesh > &)=0
 
- 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 ()
 

Static Public Member Functions

static tmp< ddtScheme< Type > > New (const fvMesh &mesh, Istream &schemeData)
 Return a pointer to a new ddtScheme created on freestore. More...
 

Protected Member Functions

 ddtScheme (const ddtScheme &)=delete
 No copy construct. More...
 
void operator= (const ddtScheme &)=delete
 No copy assignment. More...
 

Protected Attributes

const fvMeshmesh_
 
scalar ddtPhiCoeff_
 Input for fvcDdtPhiCoeff. More...
 

Additional Inherited Members

- Static Public Attributes inherited from ddtSchemeBase
static bool experimentalDdtCorr
 

Detailed Description

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

Abstract base class for ddt schemes.

Source files

Definition at line 83 of file ddtScheme.H.

Member Typedef Documentation

◆ fluxFieldType

Definition at line 230 of file ddtScheme.H.

Constructor & Destructor Documentation

◆ ddtScheme() [1/3]

ddtScheme ( const ddtScheme< Type > &  )
protecteddelete

No copy construct.

◆ ddtScheme() [2/3]

ddtScheme ( const fvMesh mesh)
inline

Construct from mesh.

Definition at line 131 of file ddtScheme.H.

◆ ddtScheme() [3/3]

ddtScheme ( const fvMesh mesh,
Istream is 
)
inline

Construct from mesh and Istream.

Definition at line 138 of file ddtScheme.H.

◆ ~ddtScheme()

virtual ~ddtScheme ( )
virtualdefault

Destructor.

Member Function Documentation

◆ operator=()

void operator= ( const ddtScheme< Type > &  )
protecteddelete

No copy assignment.

◆ type()

virtual const word & type ( ) const
pure virtual

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( tmp  ,
ddtScheme< Type >  ,
Istream  ,
(const fvMesh &mesh, Istream &schemeData)  ,
(mesh, schemeData)   
)

◆ New()

tmp< ddtScheme< Type > > New ( const fvMesh mesh,
Istream schemeData 
)
static

Return a pointer to a new ddtScheme created on freestore.

Definition at line 49 of file ddtScheme.C.

References Foam::endl(), IOstream::eof(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, FatalIOErrorInLookup, InfoInFunction, and mesh.

Here is the call graph for this function:

◆ mesh()

const fvMesh & mesh ( ) const
inline

Return mesh reference.

Definition at line 162 of file ddtScheme.H.

References ddtScheme< Type >::mesh_.

Referenced by backwardDdtScheme< Type >::mesh(), boundedDdtScheme< Type >::mesh(), CoEulerDdtScheme< Type >::mesh(), CrankNicolsonDdtScheme< Type >::mesh(), EulerDdtScheme< Type >::mesh(), localEulerDdtScheme< Type >::mesh(), SLTSDdtScheme< Type >::mesh(), and steadyStateDdtScheme< Type >::mesh().

Here is the caller graph for this function:

◆ fvcDdt() [1/6]

◆ fvcDdt() [2/6]

◆ fvcDdt() [3/6]

◆ fvcDdt() [4/6]

◆ fvcDdt() [5/6]

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

Implemented in backwardDdtScheme< Type >, boundedDdtScheme< Type >, CoEulerDdtScheme< Type >, CrankNicolsonDdtScheme< Type >, EulerDdtScheme< Type >, localEulerDdtScheme< Type >, SLTSDdtScheme< Type >, and steadyStateDdtScheme< Type >.

Definition at line 91 of file ddtScheme.C.

References NotImplemented, and GeometricField< Type, PatchField, GeoMesh >::null().

Here is the call graph for this function:

◆ fvcDdt() [6/6]

tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > fvcDdt ( const GeometricField< Type, fvsPatchField, surfaceMesh > &  sf)
virtual

Reimplemented in EulerDdtScheme< Type >, and localEulerDdtScheme< Type >.

Definition at line 127 of file ddtScheme.C.

References NotImplemented, and GeometricField< Type, PatchField, GeoMesh >::null().

Here is the call graph for this function:

◆ fvmDdt() [1/4]

◆ fvmDdt() [2/4]

◆ fvmDdt() [3/4]

◆ fvmDdt() [4/4]

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

Implemented in backwardDdtScheme< Type >, boundedDdtScheme< Type >, CoEulerDdtScheme< Type >, CrankNicolsonDdtScheme< Type >, EulerDdtScheme< Type >, localEulerDdtScheme< Type >, SLTSDdtScheme< Type >, and steadyStateDdtScheme< Type >.

Definition at line 108 of file ddtScheme.C.

References alpha, DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, Foam::dimVol, Foam::New(), NotImplemented, and rho.

Here is the call graph for this function:

◆ fvcDdtPhiCoeff() [1/4]

tmp< surfaceScalarField > fvcDdtPhiCoeff ( const GeometricField< Type, fvPatchField, volMesh > &  U,
const fluxFieldType phi,
const fluxFieldType phiCorr 
)

Definition at line 143 of file ddtScheme.C.

References boundary, GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), Foam::dimless, Foam::endl(), forAll, Foam::gAverage(), Foam::gMax(), Foam::gMin(), InfoInFunction, Foam::mag(), mesh, Foam::min(), phi, GeometricField< Type, PatchField, GeoMesh >::primitiveField(), tmp< T >::ref(), and U.

Here is the call graph for this function:

◆ fvcDdtPhiCoeffExperimental()

tmp< surfaceScalarField > fvcDdtPhiCoeffExperimental ( const GeometricField< Type, fvPatchField, volMesh > &  U,
const fluxFieldType phi,
const fluxFieldType phiCorr 
)

Definition at line 217 of file ddtScheme.C.

References boundary, GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), Foam::dimless, Foam::endl(), forAll, Foam::gAverage(), Foam::gMax(), Foam::gMin(), InfoInFunction, Foam::mag(), mesh, Foam::min(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), tmp< T >::ref(), DimensionedField< Type, GeoMesh >::setOriented(), and U.

Here is the call graph for this function:

◆ fvcDdtPhiCoeff() [2/4]

tmp< surfaceScalarField > fvcDdtPhiCoeff ( const GeometricField< Type, fvPatchField, volMesh > &  U,
const fluxFieldType phi,
const fluxFieldType phiCorr,
const volScalarField rho 
)

Definition at line 298 of file ddtScheme.C.

References Foam::fvc::interpolate(), phi, rho, and U.

Here is the call graph for this function:

◆ fvcDdtPhiCoeff() [3/4]

tmp< surfaceScalarField > fvcDdtPhiCoeff ( const GeometricField< Type, fvPatchField, volMesh > &  U,
const fluxFieldType phi 
)

Definition at line 324 of file ddtScheme.C.

References Foam::fvc::dotInterpolate(), mesh, phi, and U.

Here is the call graph for this function:

◆ fvcDdtPhiCoeff() [4/4]

tmp< surfaceScalarField > fvcDdtPhiCoeff ( const GeometricField< Type, fvPatchField, volMesh > &  rhoU,
const fluxFieldType phi,
const volScalarField rho 
)

Definition at line 354 of file ddtScheme.C.

References Foam::fvc::dotInterpolate(), Foam::fvc::interpolate(), mesh, phi, and rho.

Here is the call graph for this function:

◆ fvcDdtUfCorr() [1/2]

◆ fvcDdtPhiCorr() [1/2]

◆ fvcDdtUfCorr() [2/2]

◆ fvcDdtPhiCorr() [2/2]

◆ meshPhi()

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

Implemented in backwardDdtScheme< Type >, boundedDdtScheme< Type >, CoEulerDdtScheme< Type >, CrankNicolsonDdtScheme< Type >, EulerDdtScheme< Type >, localEulerDdtScheme< Type >, SLTSDdtScheme< Type >, and steadyStateDdtScheme< Type >.

Referenced by Foam::fvc::meshPhi().

Here is the caller graph for this function:

Member Data Documentation

◆ mesh_

const fvMesh& mesh_
protected

Definition at line 93 of file ddtScheme.H.

Referenced by ddtScheme< Type >::mesh().

◆ ddtPhiCoeff_

scalar ddtPhiCoeff_
protected

Input for fvcDdtPhiCoeff.

If set to -1 (default) the code will calculate the coefficient If > 0 the coupling coeff is set to this value

Definition at line 98 of file ddtScheme.H.

Referenced by backwardDdtScheme< Type >::backwardDdtScheme().


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