Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as the previous time-step ddt. More...
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 |
![]() | |
typedef GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > | fluxFieldType |
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 tmp< ddtScheme< Type > > | New (const fvMesh &mesh, Istream &schemeData) |
Return a pointer to a new ddtScheme created on freestore. More... | |
![]() | |
ddtScheme (const ddtScheme &) | |
Disallow copy construct. More... | |
void | operator= (const ddtScheme &) |
Disallow default bitwise assignment. More... | |
![]() | |
const fvMesh & | mesh_ |
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.
cnCoeff = 1.0/(1.0 + ocCoeff);
Definition at line 94 of file CrankNicolsonDdtScheme.H.
typedef ddtScheme<Type>::fluxFieldType fluxFieldType |
Definition at line 287 of file CrankNicolsonDdtScheme.H.
|
private |
Disallow default bitwise copy construct.
|
inline |
Construct from mesh.
Definition at line 196 of file CrankNicolsonDdtScheme.H.
|
inline |
Construct from mesh and Istream.
Definition at line 203 of file CrankNicolsonDdtScheme.H.
References Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, and CrankNicolsonDdtScheme< Type >::ocCoeff_.
|
private |
Disallow default bitwise assignment.
|
private |
Definition at line 106 of file CrankNicolsonDdtScheme.C.
References IOobject::AUTO_WRITE, Foam::dimTime, mesh, IOobject::MUST_READ, Foam::name(), IOobject::NO_READ, Time::startTime(), Time::startTimeIndex(), regIOobject::store(), TimeState::timeIndex(), timeName, Time::timeName(), and dimensioned::value().
|
private |
Check if the ddt0 needs to be evaluated for this time-step.
Definition at line 187 of file CrankNicolsonDdtScheme.C.
References mesh.
|
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.
|
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.
|
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.
|
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.
|
private |
Return ddt0 multiplied by the off-centreing coefficient.
Definition at line 255 of file CrankNicolsonDdtScheme.C.
TypeName | ( | "CrankNicolson" | ) |
Runtime type information.
|
inline |
Return mesh reference.
Definition at line 223 of file CrankNicolsonDdtScheme.H.
References ddtScheme< Type >::mesh().
|
inline |
Return the off-centreing coefficient.
Definition at line 229 of file CrankNicolsonDdtScheme.H.
References CrankNicolsonDdtScheme< Type >::ocCoeff_.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 285 of file CrankNicolsonDdtScheme.C.
References dimensioned::dimensions(), Foam::dimTime, mesh, dimensioned::name(), and timeName.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 347 of file CrankNicolsonDdtScheme.C.
References GeometricField::boundaryField(), dimensioned::dimensions(), Foam::fv::ff(), GeometricField::internalField(), mesh, GeometricField::oldTime(), timeName, and dimensioned::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 436 of file CrankNicolsonDdtScheme.C.
References GeometricField::boundaryField(), dimensioned::dimensions(), Foam::fv::ff(), GeometricField::internalField(), mesh, GeometricField::oldTime(), rho, timeName, and dimensioned::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 526 of file CrankNicolsonDdtScheme.C.
References GeometricField::boundaryField(), dimensioned::dimensions(), Foam::fv::ff(), GeometricField::internalField(), mesh, GeometricField::oldTime(), rho, timeName, and dimensioned::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 626 of file CrankNicolsonDdtScheme.C.
References Foam::constant::atomic::alpha, GeometricField::boundaryField(), dimensioned::dimensions(), Foam::fv::ff(), GeometricField::internalField(), mesh, dimensioned::name(), GeometricField::oldTime(), rho, timeName, and dimensioned::value().
|
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().
|
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().
|
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().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1010 of file CrankNicolsonDdtScheme.C.
References Foam::constant::atomic::alpha, dimensioned::dimensions(), Foam::dimTime, Foam::dimVol, Foam::fv::ff(), mesh, dimensioned::name(), GeometricField::oldTime(), rho, and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1121 of file CrankNicolsonDdtScheme.C.
References Foam::fvc::interpolate(), mesh, timeName, U, and Uf.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 1182 of file CrankNicolsonDdtScheme.C.
References Foam::fvc::interpolate(), mesh, phi, timeName, and U.
|
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.
|
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.
|
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.
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 | ||
) |
|
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().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.