Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as the previous time-step ddt. More...
Public Types | |
typedef ddtScheme< Type >::fluxFieldType | fluxFieldType |
![]() | |
typedef GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > | fluxFieldType |
Additional Inherited Members | |
![]() | |
static tmp< ddtScheme< Type > > | New (const fvMesh &mesh, Istream &schemeData) |
![]() | |
static bool | experimentalDdtCorr |
![]() | |
ddtScheme (const ddtScheme &)=delete | |
void | operator= (const ddtScheme &)=delete |
![]() | |
const fvMesh & | mesh_ |
scalar | ddtPhiCoeff_ |
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 288 of file CrankNicolsonDdtScheme.H.
CrankNicolsonDdtScheme | ( | const fvMesh & | mesh | ) |
Definition at line 279 of file CrankNicolsonDdtScheme.C.
References CrankNicolsonDdtScheme< Type >::mesh(), polyMesh::moving(), and fvMesh::V00().
CrankNicolsonDdtScheme | ( | const fvMesh & | mesh, |
Istream & | is | ||
) |
Definition at line 295 of file CrankNicolsonDdtScheme.C.
References dict, Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, token::isNumber(), mesh, token::number(), ocCoeff, and Istream::putBack().
TypeName | ( | "CrankNicolson" | ) |
|
inline |
Definition at line 224 of file CrankNicolsonDdtScheme.H.
References ddtScheme< Type >::mesh().
Referenced by CrankNicolsonDdtScheme< Type >::CrankNicolsonDdtScheme(), and CrankNicolsonDdtScheme< Type >::ocCoeff().
|
inline |
Definition at line 230 of file CrankNicolsonDdtScheme.H.
References CrankNicolsonDdtScheme< Type >::mesh().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 341 of file CrankNicolsonDdtScheme.C.
References dimensioned::dimensions(), Foam::dimTime, Foam::stringOps::evaluate(), mesh, dimensioned::name(), tmp::ref(), timeName, and Foam::Zero.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 398 of file CrankNicolsonDdtScheme.C.
References GeometricField::boundaryField(), Foam::stringOps::evaluate(), Foam::fv::ff(), mesh, GeometricField::oldTime(), timeName, and dimensioned::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 485 of file CrankNicolsonDdtScheme.C.
References GeometricField::boundaryField(), dimensioned::dimensions(), Foam::stringOps::evaluate(), Foam::fv::ff(), mesh, GeometricField::oldTime(), GeometricField::primitiveField(), rho, timeName, and dimensioned::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 575 of file CrankNicolsonDdtScheme.C.
References GeometricField::boundaryField(), dimensioned::dimensions(), Foam::stringOps::evaluate(), Foam::fv::ff(), mesh, GeometricField::oldTime(), GeometricField::primitiveField(), rho, timeName, and dimensioned::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 675 of file CrankNicolsonDdtScheme.C.
References Foam::constant::atomic::alpha, GeometricField::boundaryField(), dimensioned::dimensions(), Foam::stringOps::evaluate(), Foam::fv::ff(), mesh, dimensioned::name(), GeometricField::oldTime(), GeometricField::primitiveField(), rho, timeName, and dimensioned::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 805 of file CrankNicolsonDdtScheme.C.
References Foam::dimTime, Foam::dimVol, Foam::stringOps::evaluate(), Foam::fv::ff(), mesh, GeometricField::oldTime(), and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 888 of file CrankNicolsonDdtScheme.C.
References Foam::dimTime, Foam::dimVol, Foam::stringOps::evaluate(), Foam::fv::ff(), mesh, GeometricField::oldTime(), rho, and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 970 of file CrankNicolsonDdtScheme.C.
References Foam::dimTime, Foam::dimVol, Foam::stringOps::evaluate(), Foam::fv::ff(), mesh, GeometricField::oldTime(), rho, and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1060 of file CrankNicolsonDdtScheme.C.
References Foam::constant::atomic::alpha, dimensioned::dimensions(), Foam::dimTime, Foam::dimVol, Foam::stringOps::evaluate(), Foam::fv::ff(), mesh, dimensioned::name(), GeometricField::oldTime(), rho, and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1171 of file CrankNicolsonDdtScheme.C.
References Foam::stringOps::evaluate(), Foam::fvc::interpolate(), mesh, timeName, U, and Uf.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1232 of file CrankNicolsonDdtScheme.C.
References Foam::fvc::dotInterpolate(), Foam::stringOps::evaluate(), mesh, phi, timeName, and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1295 of file CrankNicolsonDdtScheme.C.
References Foam::abort(), ddtCorr, Foam::dimVelocity, Foam::stringOps::evaluate(), Foam::FatalError, FatalErrorInFunction, Foam::fvc::interpolate(), mesh, rho, timeName, U, and Uf.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1452 of file CrankNicolsonDdtScheme.C.
References Foam::abort(), ddtCorr, Foam::dimArea, Foam::dimVelocity, Foam::fvc::dotInterpolate(), Foam::stringOps::evaluate(), Foam::FatalError, FatalErrorInFunction, mesh, phi, rho, timeName, and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1597 of file CrankNicolsonDdtScheme.C.
References Foam::dimVolume, Foam::stringOps::evaluate(), mesh, Foam::name(), IOobject::NO_READ, IOobject::NO_WRITE, phi, and timeName.
CrankNicolsonDdtScheme<Type>::template DDt0Field<GeoField>& ddt0_ | ( | const word & | name, |
const dimensionSet & | dims | ||
) |
Definition at line 103 of file CrankNicolsonDdtScheme.C.
References IOobject::AUTO_WRITE, Foam::dimTime, mesh, IOobject::MUST_READ, Foam::name(), IOobject::NO_READ, runTime, regIOobject::store(), timeName, and Foam::Zero.
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 | ||
) |
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.