CrankNicolsonDdtScheme< Type > Class Template Reference

Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as the previous time-step ddt. More...

Inheritance diagram for CrankNicolsonDdtScheme< Type >:
Collaboration diagram for CrankNicolsonDdtScheme< 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 ("CrankNicolson")
 Runtime type information. More...
 
 CrankNicolsonDdtScheme (const fvMesh &mesh)
 Construct from mesh. More...
 
 CrankNicolsonDdtScheme (const fvMesh &mesh, Istream &is)
 Construct from mesh and Istream. More...
 
 CrankNicolsonDdtScheme (const CrankNicolsonDdtScheme &)=delete
 Disallow default bitwise copy construction. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
scalar ocCoeff () const
 Return the current off-centreing coefficient. 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 CrankNicolsonDdtScheme &)=delete
 Disallow default bitwise assignment. More...
 
template<class GeoField >
CrankNicolsonDdtScheme< Type >::template DDt0Field< GeoField > & ddt0_ (const word &name, const dimensionSet &dims)
 
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::CrankNicolsonDdtScheme< Type >

Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as the previous time-step ddt.

The Crank-Nicolson scheme is often unstable for complex flows in complex geometries and it is necessary to "off-centre" the scheme to stabilize it while retaining greater temporal accuracy than the first-order Euler-implicit scheme. Off-centering is specified via the mandatory coefficient ocCoeff in the range [0,1] following the scheme name e.g.

ddtSchemes
{
    default         CrankNicolson 0.9;
}

or with an optional "ramp" function to transition from the Euler scheme to Crank-Nicolson over a initial period to avoid start-up problems, e.g.

ddtSchemes
{
    default         CrankNicolson
    ocCoeff
    {
        type scale;
        scale linearRamp;
        duration 0.01;
        value 0.9;
    };
}

With a coefficient of 1 the scheme is fully centred and second-order, with a coefficient of 0 the scheme is equivalent to Euler-implicit. A coefficient of 0.9 has been found to be suitable for a range of cases for which higher-order accuracy in time is needed and provides similar accuracy and stability to the backward scheme. However, the backward scheme has been found to be more robust and provides formal second-order accuracy in time.

The advantage of the Crank-Nicolson scheme over backward is that only the new and old-time values are needed, the additional terms relating to the fluxes and sources are evaluated at the mid-point of the time-step which provides the opportunity to limit the fluxes in such a way as to ensure boundedness while maintaining greater accuracy in time compared to the Euler-implicit scheme. This approach is now used with MULES in the interFoam family of solvers. Boundedness cannot be guaranteed with the backward scheme.

Note
The Crank-Nicolson coefficient for the implicit part of the RHS is related to the off-centering coefficient by
    cnCoeff = 1.0/(1.0 + ocCoeff);
See also
Foam::Euler Foam::backward
Source files

Definition at line 110 of file CrankNicolsonDdtScheme.H.

Member Typedef Documentation

◆ fluxFieldType

Definition at line 283 of file CrankNicolsonDdtScheme.H.

Constructor & Destructor Documentation

◆ CrankNicolsonDdtScheme() [1/3]

CrankNicolsonDdtScheme ( const fvMesh mesh)

Construct from mesh.

Definition at line 281 of file CrankNicolsonDdtScheme.C.

References polyMesh::moving(), and fvMesh::V00().

Here is the call graph for this function:

◆ CrankNicolsonDdtScheme() [2/3]

◆ CrankNicolsonDdtScheme() [3/3]

CrankNicolsonDdtScheme ( const CrankNicolsonDdtScheme< Type > &  )
delete

Disallow default bitwise copy construction.

Member Function Documentation

◆ TypeName()

TypeName ( "CrankNicolson"  )

Runtime type information.

◆ mesh()

const fvMesh& mesh ( ) const
inline

◆ ocCoeff()

scalar ocCoeff ( ) const
inline

Return the current off-centreing coefficient.

Definition at line 225 of file CrankNicolsonDdtScheme.H.

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

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

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

◆ fvcDdt() [1/5]

◆ fvcDdt() [2/5]

◆ fvcDdt() [3/5]

◆ fvcDdt() [4/5]

◆ fvcDdt() [5/5]

◆ fvmDdt() [1/4]

◆ fvmDdt() [2/4]

◆ fvmDdt() [3/4]

◆ fvmDdt() [4/4]

◆ fvcDdtUfCorr() [1/4]

◆ fvcDdtPhiCorr() [1/4]

◆ fvcDdtUfCorr() [2/4]

◆ fvcDdtPhiCorr() [2/4]

◆ meshPhi()

◆ operator=()

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

Disallow default bitwise assignment.

◆ ddt0_()

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