Stabilised local time-step first-order Euler implicit/explicit ddt. The time-step is adjusted locally so that an advective equations remains diagonally dominant. 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_ |
Stabilised local time-step first-order Euler implicit/explicit ddt. The time-step is adjusted locally so that an advective equations remains diagonally dominant.
This scheme should only be used for steady-state computations using transient codes where local time-stepping is preferably to under-relaxation for transport consistency reasons.
Definition at line 64 of file SLTSDdtScheme.H.
typedef ddtScheme<Type>::fluxFieldType fluxFieldType |
Definition at line 175 of file SLTSDdtScheme.H.
|
inline |
Definition at line 105 of file SLTSDdtScheme.H.
TypeName | ( | "SLTS" | ) |
|
inline |
Definition at line 117 of file SLTSDdtScheme.H.
References ddtScheme< Type >::mesh().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 151 of file SLTSDdtScheme.C.
References dimensioned::dimensions(), Foam::dimTime, mesh, dimensioned::name(), GeometricField::primitiveField(), tmp::ref(), timeName, dimensioned::value(), and Foam::Zero.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 200 of file SLTSDdtScheme.C.
References GeometricField::boundaryField(), mesh, GeometricField::oldTime(), GeometricField::primitiveField(), and timeName.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 251 of file SLTSDdtScheme.C.
References GeometricField::boundaryField(), mesh, GeometricField::oldTime(), GeometricField::primitiveField(), rho, and timeName.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 303 of file SLTSDdtScheme.C.
References GeometricField::boundaryField(), mesh, GeometricField::oldTime(), GeometricField::primitiveField(), rho, and timeName.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 358 of file SLTSDdtScheme.C.
References Foam::constant::atomic::alpha, GeometricField::boundaryField(), dimensioned::dimensions(), mesh, dimensioned::name(), GeometricField::oldTime(), GeometricField::primitiveField(), rho, and timeName.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 427 of file SLTSDdtScheme.C.
References Foam::dimTime, Foam::dimVol, mesh, GeometricField::oldTime(), tmp::ref(), and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 462 of file SLTSDdtScheme.C.
References Foam::dimTime, Foam::dimVol, mesh, GeometricField::oldTime(), tmp::ref(), rho, and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 499 of file SLTSDdtScheme.C.
References Foam::dimTime, Foam::dimVol, mesh, GeometricField::oldTime(), tmp::ref(), rho, and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 538 of file SLTSDdtScheme.C.
References Foam::constant::atomic::alpha, dimensioned::dimensions(), Foam::dimTime, Foam::dimVol, mesh, GeometricField::oldTime(), tmp::ref(), rho, and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 581 of file SLTSDdtScheme.C.
References Foam::fvc::dotInterpolate(), Foam::fvc::interpolate(), mesh, timeName, U, and Uf.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 614 of file SLTSDdtScheme.C.
References Foam::fvc::dotInterpolate(), Foam::fvc::interpolate(), mesh, phi, timeName, and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 646 of file SLTSDdtScheme.C.
References Foam::abort(), Foam::dimDensity, Foam::dimVelocity, Foam::fvc::dotInterpolate(), Foam::FatalError, FatalErrorInFunction, Foam::fvc::interpolate(), mesh, rho, timeName, U, and Uf.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 731 of file SLTSDdtScheme.C.
References Foam::abort(), Foam::dimArea, Foam::dimVelocity, Foam::fvc::dotInterpolate(), Foam::FatalError, FatalErrorInFunction, mesh, phi, rho, timeName, and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 821 of file SLTSDdtScheme.C.
References Foam::dimTime, Foam::dimVolume, mesh, IOobject::NO_READ, IOobject::NO_WRITE, tmp::ref(), 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.