boundedDdtScheme< Type > Class Template Reference

Bounded form of the selected ddt scheme. More...

Inheritance diagram for boundedDdtScheme< Type >:
Collaboration diagram for boundedDdtScheme< 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 ("bounded")
 Runtime type information. More...
 
 boundedDdtScheme (const fvMesh &mesh, Istream &is)
 Construct from mesh and Istream. More...
 
 boundedDdtScheme (const boundedDdtScheme &)=delete
 Disallow default bitwise copy construction. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const dimensioned< Type > &)
 
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const GeometricField< Type, fvPatchField, volMesh > &)
 
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &)
 
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const volScalarField &, const GeometricField< Type, fvPatchField, volMesh > &)
 
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &psi)
 
virtual tmp< fvMatrix< Type > > fvmDdt (const GeometricField< Type, fvPatchField, volMesh > &)
 
virtual tmp< fvMatrix< Type > > fvmDdt (const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &)
 
virtual tmp< fvMatrix< Type > > fvmDdt (const volScalarField &, const GeometricField< Type, fvPatchField, volMesh > &)
 
virtual tmp< fvMatrix< Type > > fvmDdt (const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &psi)
 
virtual tmp< fluxFieldTypefvcDdtUfCorr (const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
 
virtual tmp< fluxFieldTypefvcDdtPhiCorr (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
 
virtual tmp< fluxFieldTypefvcDdtUfCorr (const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
 
virtual tmp< fluxFieldTypefvcDdtPhiCorr (const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
 
virtual tmp< surfaceScalarFieldmeshPhi (const GeometricField< Type, fvPatchField, volMesh > &)
 
void operator= (const boundedDdtScheme &)=delete
 Disallow default bitwise assignment. More...
 
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 &)
 Construct from mesh and Istream. More...
 
 ddtScheme (const ddtScheme &)
 Disallow default bitwise copy construction. More...
 
virtual ~ddtScheme ()
 Destructor. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > fvcDdt (const GeometricField< Type, fvsPatchField, surfaceMesh > &)
 
virtual tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
 
virtual tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr, const volScalarField &rho)
 
virtual tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
 
virtual tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &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...
 

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 initializing count to 0. More...
 
- Protected Attributes inherited from ddtScheme< Type >
const fvMeshmesh_
 

Detailed Description

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

Bounded form of the selected ddt scheme.

Boundedness is achieved by subtracting ddt(phi)*vf or Sp(ddt(rho), vf) which is non-conservative if ddt(rho) != 0 but conservative otherwise.

Can be used for the ddt of bounded scalar properties to improve stability if insufficient convergence of the pressure equation causes temporary divergence of the flux field.

Source files

Definition at line 61 of file boundedDdtScheme.H.

Member Typedef Documentation

◆ fluxFieldType

Definition at line 153 of file boundedDdtScheme.H.

Constructor & Destructor Documentation

◆ boundedDdtScheme() [1/2]

boundedDdtScheme ( const fvMesh mesh,
Istream is 
)
inline

Construct from mesh and Istream.

Definition at line 79 of file boundedDdtScheme.H.

◆ boundedDdtScheme() [2/2]

boundedDdtScheme ( const boundedDdtScheme< Type > &  )
delete

Disallow default bitwise copy construction.

Member Function Documentation

◆ TypeName()

TypeName ( "bounded"  )

Runtime type information.

◆ mesh()

const fvMesh& mesh ( ) const
inline

Return mesh reference.

Definition at line 95 of file boundedDdtScheme.H.

References alpha(), boundedDdtScheme< Type >::fvcDdt(), boundedDdtScheme< Type >::fvmDdt(), ddtScheme< Type >::mesh(), psi, and rho.

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 boundedDdtScheme.C.

Referenced by boundedDdtScheme< Type >::fvcDdt(), and boundedDdtScheme< Type >::mesh().

Here is the caller 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 58 of file boundedDdtScheme.C.

References boundedDdtScheme< Type >::fvcDdt().

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 69 of file boundedDdtScheme.C.

References boundedDdtScheme< Type >::fvcDdt().

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 81 of file boundedDdtScheme.C.

References Foam::fvc::ddt(), and boundedDdtScheme< Type >::fvcDdt().

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 93 of file boundedDdtScheme.C.

References Foam::fvc::ddt(), and boundedDdtScheme< Type >::fvmDdt().

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 106 of file boundedDdtScheme.C.

Referenced by boundedDdtScheme< Type >::fvcDdt(), boundedDdtScheme< Type >::fvmDdt(), and boundedDdtScheme< Type >::mesh().

Here is the caller graph for this function:

◆ fvmDdt() [2/4]

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

Implements ddtScheme< Type >.

Definition at line 117 of file boundedDdtScheme.C.

References boundedDdtScheme< Type >::fvmDdt().

Here is the call graph for this function:

◆ fvmDdt() [3/4]

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

Implements ddtScheme< Type >.

Definition at line 129 of file boundedDdtScheme.C.

References Foam::fvc::ddt(), boundedDdtScheme< Type >::fvmDdt(), and Foam::fvm::Sp().

Here is the call graph for this function:

◆ 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 141 of file boundedDdtScheme.C.

References Foam::fvc::ddt(), boundedDdtScheme< Type >::fvcDdtUfCorr(), and Foam::fvm::Sp().

Here is the call graph for this function:

◆ fvcDdtUfCorr() [1/4]

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

Implements ddtScheme< Type >.

Definition at line 156 of file boundedDdtScheme.C.

References boundedDdtScheme< Type >::fvcDdtPhiCorr().

Referenced by boundedDdtScheme< Type >::fvcDdtPhiCorr(), and boundedDdtScheme< Type >::fvmDdt().

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

◆ fvcDdtPhiCorr() [1/4]

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

Implements ddtScheme< Type >.

Definition at line 168 of file boundedDdtScheme.C.

References boundedDdtScheme< Type >::fvcDdtUfCorr().

Referenced by boundedDdtScheme< Type >::fvcDdtUfCorr().

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

◆ fvcDdtUfCorr() [2/4]

tmp< typename boundedDdtScheme< 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 180 of file boundedDdtScheme.C.

References boundedDdtScheme< Type >::fvcDdtPhiCorr().

Here is the call graph for this function:

◆ fvcDdtPhiCorr() [2/4]

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

Implements ddtScheme< Type >.

Definition at line 193 of file boundedDdtScheme.C.

References boundedDdtScheme< Type >::meshPhi().

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 205 of file boundedDdtScheme.C.

Referenced by boundedDdtScheme< Type >::fvcDdtPhiCorr().

Here is the caller graph for this function:

◆ operator=()

void operator= ( const boundedDdtScheme< Type > &  )
delete

Disallow default bitwise assignment.

◆ 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: