Data Structures | Public Types | Public Member Functions | Private Member Functions | Private Attributes
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 >:
Inheritance graph
[legend]
Collaboration diagram for CrankNicolsonDdtScheme< Type >:
Collaboration graph
[legend]

Data Structures

class  DDt0Field
 Class to store the ddt0 fields on the objectRegistry for use in the. More...
 

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...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
scalar ocCoeff () const
 Return the off-centreing coefficient. 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 > &)
 
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. 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...
 
virtual ~ddtScheme ()
 Destructor. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
 
tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
 
- Public Member Functions inherited from refCount
 refCount ()
 Construct null with zero count. More...
 
int count () const
 Return the reference count. More...
 
bool okToDelete () const
 Return true if the reference count is zero. More...
 
void resetRefCount ()
 Reset the reference count to 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...
 

Private Member Functions

 CrankNicolsonDdtScheme (const CrankNicolsonDdtScheme &)
 Disallow default bitwise copy construct. More...
 
void operator= (const CrankNicolsonDdtScheme &)
 Disallow default bitwise assignment. More...
 
template<class GeoField >
DDt0Field< GeoField > & ddt0_ (const word &name, const dimensionSet &dims)
 
template<class GeoField >
bool evaluate (const DDt0Field< GeoField > &ddt0) const
 Check if the ddt0 needs to be evaluated for this time-step. More...
 
template<class GeoField >
scalar coef_ (const DDt0Field< GeoField > &) const
 Return the coefficient for Euler scheme for the first time-step. More...
 
template<class GeoField >
scalar coef0_ (const DDt0Field< GeoField > &) const
 Return the old time-step coefficient for Euler scheme for the. More...
 
template<class GeoField >
dimensionedScalar rDtCoef_ (const DDt0Field< GeoField > &) const
 Return the reciprocal time-step coefficient for Euler for the. More...
 
template<class GeoField >
dimensionedScalar rDtCoef0_ (const DDt0Field< GeoField > &) const
 Return the reciprocal old time-step coefficient for Euler for the. More...
 
template<class GeoField >
tmp< GeoField > offCentre_ (const GeoField &ddt0) const
 Return ddt0 multiplied by the off-centreing coefficient. More...
 

Private Attributes

scalar ocCoeff_
 Off-centering coefficient, 1 -> CN, less than one blends with EI. 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 ddtScheme< Type >
 ddtScheme (const ddtScheme &)
 Disallow copy construct. More...
 
void operator= (const ddtScheme &)
 Disallow default bitwise assignment. 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 in the range [0,1] following the scheme name e.g.

ddtSchemes
{
    default         CrankNicolson 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 94 of file CrankNicolsonDdtScheme.H.

Member Typedef Documentation

◆ fluxFieldType

Definition at line 287 of file CrankNicolsonDdtScheme.H.

Constructor & Destructor Documentation

◆ CrankNicolsonDdtScheme() [1/3]

CrankNicolsonDdtScheme ( const CrankNicolsonDdtScheme< Type > &  )
private

Disallow default bitwise copy construct.

◆ CrankNicolsonDdtScheme() [2/3]

CrankNicolsonDdtScheme ( const fvMesh mesh)
inline

Construct from mesh.

Definition at line 196 of file CrankNicolsonDdtScheme.H.

◆ CrankNicolsonDdtScheme() [3/3]

CrankNicolsonDdtScheme ( const fvMesh mesh,
Istream is 
)
inline

Construct from mesh and Istream.

Definition at line 203 of file CrankNicolsonDdtScheme.H.

References Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, and CrankNicolsonDdtScheme< Type >::ocCoeff_.

Here is the call graph for this function:

Member Function Documentation

◆ operator=()

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

Disallow default bitwise assignment.

◆ ddt0_()

CrankNicolsonDdtScheme< Type >::DDt0Field< GeoField > & ddt0_ ( const word name,
const dimensionSet dims 
)
private

◆ evaluate()

bool evaluate ( const DDt0Field< GeoField > &  ddt0) const
private

Check if the ddt0 needs to be evaluated for this time-step.

Definition at line 187 of file CrankNicolsonDdtScheme.C.

References mesh.

◆ coef_()

scalar coef_ ( const DDt0Field< GeoField > &  ddt0) const
private

Return the coefficient for Euler scheme for the first time-step.

for and CN thereafter

Definition at line 197 of file CrankNicolsonDdtScheme.C.

References mesh, CrankNicolsonDdtScheme< Type >::DDt0Field< GeoField >::startTimeIndex(), and timeIndex.

Here is the call graph for this function:

◆ coef0_()

scalar coef0_ ( const DDt0Field< GeoField > &  ddt0) const
private

Return the old time-step coefficient for Euler scheme for the.

second time-step and for CN thereafter

Definition at line 215 of file CrankNicolsonDdtScheme.C.

References mesh, CrankNicolsonDdtScheme< Type >::DDt0Field< GeoField >::startTimeIndex(), and timeIndex.

Here is the call graph for this function:

◆ rDtCoef_()

dimensionedScalar rDtCoef_ ( const DDt0Field< GeoField > &  ddt0) const
private

Return the reciprocal time-step coefficient for Euler for the.

first time-step and CN thereafter

Definition at line 233 of file CrankNicolsonDdtScheme.C.

References mesh.

◆ rDtCoef0_()

dimensionedScalar rDtCoef0_ ( const DDt0Field< GeoField > &  ddt0) const
private

Return the reciprocal old time-step coefficient for Euler for the.

second time-step and CN thereafter

Definition at line 244 of file CrankNicolsonDdtScheme.C.

References mesh.

◆ offCentre_()

tmp< GeoField > offCentre_ ( const GeoField &  ddt0) const
private

Return ddt0 multiplied by the off-centreing coefficient.

Definition at line 255 of file CrankNicolsonDdtScheme.C.

◆ TypeName()

TypeName ( "CrankNicolson"  )

Runtime type information.

◆ mesh()

const fvMesh& mesh ( ) const
inline

Return mesh reference.

Definition at line 223 of file CrankNicolsonDdtScheme.H.

References ddtScheme< Type >::mesh().

Here is the call graph for this function:

◆ ocCoeff()

scalar ocCoeff ( ) const
inline

Return the off-centreing coefficient.

Definition at line 229 of file CrankNicolsonDdtScheme.H.

References CrankNicolsonDdtScheme< Type >::ocCoeff_.

◆ fvcDdt() [1/5]

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

Implements ddtScheme< Type >.

Definition at line 285 of file CrankNicolsonDdtScheme.C.

References dimensioned::dimensions(), Foam::dimTime, mesh, dimensioned::name(), and timeName.

Here is the call graph for this function:

◆ fvcDdt() [2/5]

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

◆ fvcDdt() [3/5]

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

◆ fvcDdt() [4/5]

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

◆ fvcDdt() [5/5]

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

◆ fvmDdt() [1/4]

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

Implements ddtScheme< Type >.

Definition at line 756 of file CrankNicolsonDdtScheme.C.

References Foam::dimTime, Foam::dimVol, Foam::fv::ff(), mesh, GeometricField::oldTime(), and fvMatrix::source().

Here is the call 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 838 of file CrankNicolsonDdtScheme.C.

References Foam::dimTime, Foam::dimVol, Foam::fv::ff(), mesh, GeometricField::oldTime(), rho, and fvMatrix::source().

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 920 of file CrankNicolsonDdtScheme.C.

References Foam::dimTime, Foam::dimVol, Foam::fv::ff(), mesh, GeometricField::oldTime(), rho, and fvMatrix::source().

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

◆ fvcDdtUfCorr() [1/4]

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

Implements ddtScheme< Type >.

Definition at line 1121 of file CrankNicolsonDdtScheme.C.

References Foam::fvc::interpolate(), mesh, timeName, U, and Uf.

Here is the call graph for this function:

◆ fvcDdtPhiCorr() [1/4]

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

Implements ddtScheme< Type >.

Definition at line 1182 of file CrankNicolsonDdtScheme.C.

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

Here is the call graph for this function:

◆ fvcDdtUfCorr() [2/4]

tmp< typename CrankNicolsonDdtScheme< 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 1243 of file CrankNicolsonDdtScheme.C.

References Foam::abort(), Foam::fvc::ddtCorr(), Foam::dimVelocity, Foam::FatalError, FatalErrorInFunction, Foam::fvc::interpolate(), mesh, rho, timeName, U, and Uf.

Here is the call graph for this function:

◆ fvcDdtPhiCorr() [2/4]

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

Implements ddtScheme< Type >.

Definition at line 1338 of file CrankNicolsonDdtScheme.C.

References Foam::abort(), Foam::fvc::ddtCorr(), Foam::dimArea, Foam::dimVelocity, Foam::FatalError, FatalErrorInFunction, Foam::fvc::interpolate(), mesh, phi, rho, timeName, and U.

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 1432 of file CrankNicolsonDdtScheme.C.

References Foam::dimVolume, mesh, Foam::name(), IOobject::NO_READ, IOobject::NO_WRITE, 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 
)

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

Field Documentation

◆ ocCoeff_

scalar ocCoeff_
private

Off-centering coefficient, 1 -> CN, less than one blends with EI.

Definition at line 140 of file CrankNicolsonDdtScheme.H.

Referenced by CrankNicolsonDdtScheme< Type >::CrankNicolsonDdtScheme(), and CrankNicolsonDdtScheme< Type >::ocCoeff().


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