Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as the previous time-step ddt. More...
#include <CrankNicolsonDdtScheme.H>


Public Types | |
| typedef ddtScheme< Type >::fluxFieldType | fluxFieldType |
| Public Types inherited from ddtScheme< Type > | |
| typedef GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > | fluxFieldType |
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. | |
| Static Public Attributes inherited from ddtSchemeBase | |
| static bool | experimentalDdtCorr |
| Flag to use experimental ddtCorr from org version Default is off for backwards compatibility. | |
| Protected Member Functions inherited from ddtScheme< Type > | |
| ddtScheme (const ddtScheme &)=delete | |
| No copy construct. | |
| void | operator= (const ddtScheme &)=delete |
| No copy assignment. | |
| Protected Attributes inherited from ddtScheme< Type > | |
| const fvMesh & | mesh_ |
| scalar | ddtPhiCoeff_ |
| Input for fvcDdtPhiCoeff. | |
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.
cnCoeff = 1.0/(1.0 + ocCoeff);
Definition at line 112 of file CrankNicolsonDdtScheme.H.
| typedef ddtScheme<Type>::fluxFieldType fluxFieldType |
Definition at line 334 of file CrankNicolsonDdtScheme.H.
| CrankNicolsonDdtScheme | ( | const fvMesh & | mesh | ) |
Construct from mesh.
Definition at line 269 of file CrankNicolsonDdtScheme.C.
References ddtScheme< Type >::ddtScheme(), and mesh().

Construct from mesh and Istream.
Definition at line 284 of file CrankNicolsonDdtScheme.C.
References ddtScheme< Type >::ddtScheme(), dict, Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, token::isNumber(), mesh(), token::number(), ocCoeff(), and Istream::putBack().

| TypeName | ( | "CrankNicolson" | ) |
|
inline |
Return mesh reference.
Definition at line 268 of file CrankNicolsonDdtScheme.H.
References ddtScheme< Type >::mesh().
Referenced by CrankNicolsonDdtScheme(), CrankNicolsonDdtScheme(), fvcDdt(), fvcDdt(), fvcDdt(), fvcDdt(), fvcDdt(), fvcDdtPhiCorr(), fvcDdtPhiCorr(), fvcDdtUfCorr(), fvcDdtUfCorr(), ocCoeff(), and TypeName().


|
inline |
Return the current off-centreing coefficient.
Definition at line 276 of file CrankNicolsonDdtScheme.H.
References mesh().
Referenced by CrankNicolsonDdtScheme().


|
virtual |
Implements ddtScheme< Type >.
Definition at line 330 of file CrankNicolsonDdtScheme.C.
References dimensioned< Type >::dimensions(), Foam::dimTime, mesh(), mesh, dimensioned< Type >::name(), tmp< T >::ref(), and timeName.

|
virtual |
Implements ddtScheme< Type >.
Definition at line 393 of file CrankNicolsonDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::fv::ff(), mesh(), mesh, IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), tmp< T >::ref(), timeName, and dimensioned< Type >::value().

|
virtual |
Implements ddtScheme< Type >.
Definition at line 493 of file CrankNicolsonDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), dimensioned< Type >::dimensions(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::fv::ff(), mesh(), mesh, IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), tmp< T >::ref(), rho, timeName, and dimensioned< Type >::value().

|
virtual |
Implements ddtScheme< Type >.
Definition at line 596 of file CrankNicolsonDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), dimensioned< Type >::dimensions(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::fv::ff(), mesh(), mesh, IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), tmp< T >::ref(), rho, timeName, and dimensioned< Type >::value().

|
virtual |
Implements ddtScheme< Type >.
Definition at line 709 of file CrankNicolsonDdtScheme.C.
References alpha, GeometricField< Type, PatchField, GeoMesh >::boundaryField(), dimensioned< Type >::dimensions(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::fv::ff(), mesh(), mesh, IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), tmp< T >::ref(), rho, timeName, and dimensioned< Type >::value().

|
virtual |
Implements ddtScheme< Type >.
Definition at line 852 of file CrankNicolsonDdtScheme.C.
References DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, Foam::dimVol, Foam::fv::ff(), mesh, IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), and tmp< T >::ref().

|
virtual |
Implements ddtScheme< Type >.
Definition at line 941 of file CrankNicolsonDdtScheme.C.
References DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, Foam::dimVol, Foam::fv::ff(), mesh, IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), tmp< T >::ref(), and rho.

|
virtual |
Implements ddtScheme< Type >.
Definition at line 1029 of file CrankNicolsonDdtScheme.C.
References DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, Foam::dimVol, Foam::fv::ff(), mesh, IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), tmp< T >::ref(), and rho.

|
virtual |
Implements ddtScheme< Type >.
Definition at line 1125 of file CrankNicolsonDdtScheme.C.
References alpha, DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, Foam::dimVol, Foam::fv::ff(), mesh, IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), tmp< T >::ref(), and rho.

|
virtual |
Implements ddtScheme< Type >.
Definition at line 1242 of file CrankNicolsonDdtScheme.C.
References ddtScheme< Type >::fvcDdtPhiCoeff(), Foam::fvc::interpolate(), mesh(), mesh, timeName, U, and Uf.

|
virtual |
Implements ddtScheme< Type >.
Definition at line 1303 of file CrankNicolsonDdtScheme.C.
References Foam::fvc::dotInterpolate(), ddtScheme< Type >::fvcDdtPhiCoeff(), mesh(), mesh, phi, timeName, and U.

|
virtual |
Implements ddtScheme< Type >.
Definition at line 1366 of file CrankNicolsonDdtScheme.C.
References Foam::abort(), Foam::dimVelocity, Foam::FatalError, FatalErrorInFunction, ddtScheme< Type >::fvcDdtPhiCoeff(), Foam::fvc::interpolate(), mesh(), mesh, GeometricField< Type, PatchField, GeoMesh >::null(), oldTime(), rho, timeName, U, and Uf.

|
virtual |
Implements ddtScheme< Type >.
Definition at line 1523 of file CrankNicolsonDdtScheme.C.
References Foam::abort(), Foam::dimArea, Foam::dimVelocity, Foam::fvc::dotInterpolate(), Foam::FatalError, FatalErrorInFunction, ddtScheme< Type >::fvcDdtPhiCoeff(), mesh(), mesh, GeometricField< Type, PatchField, GeoMesh >::null(), oldTime(), phi, rho, timeName, and U.

|
virtual |
Implements ddtScheme< Type >.
Definition at line 1668 of file CrankNicolsonDdtScheme.C.
References Foam::dimVolume, mesh, name, IOobjectOption::NO_READ, IOobjectOption::NO_REGISTER, IOobjectOption::NO_WRITE, phi, and timeName.
| CrankNicolsonDdtScheme< Type >::template DDt0Field< GeoField > & ddt0_ | ( | const word & | name, |
| const dimensionSet & | dims ) |
Definition at line 95 of file CrankNicolsonDdtScheme.C.
References IOobjectOption::AUTO_WRITE, Foam::dimTime, mesh, IOobjectOption::MUST_READ, Foam::name(), IOobjectOption::NO_READ, IOobjectOption::REGISTER, runTime, regIOobject::store(), and timeName.

| tmp< surfaceScalarField > fvcDdtUfCorr | ( | const GeometricField< scalar, fvPatchField, volMesh > & | U, |
| const GeometricField< scalar, fvsPatchField, surfaceMesh > & | Uf ) |
| tmp< surfaceScalarField > fvcDdtPhiCorr | ( | const volScalarField & | U, |
| const surfaceScalarField & | phi ) |
| tmp< surfaceScalarField > fvcDdtUfCorr | ( | const volScalarField & | rho, |
| const volScalarField & | U, | ||
| const surfaceScalarField & | Uf ) |
| tmp< surfaceScalarField > fvcDdtPhiCorr | ( | const volScalarField & | rho, |
| const volScalarField & | U, | ||
| const surfaceScalarField & | phi ) |