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 |
Private Member Functions | |
SLTSDdtScheme (const SLTSDdtScheme &) | |
Disallow default bitwise copy construct. More... | |
void | operator= (const SLTSDdtScheme &) |
Disallow default bitwise assignment. More... | |
void | relaxedDiag (scalarField &rD, const surfaceScalarField &phi) const |
Calculate a relaxed diagonal from the given flux field. More... | |
tmp< volScalarField > | SLrDeltaT () const |
Return the reciprocal of the stabilised local time-step. More... | |
Private Attributes | |
word | phiName_ |
Name of the flux field used to calculate the local time-step. More... | |
word | rhoName_ |
Name of the density field used to obtain the volumetric flux. More... | |
scalar | alpha_ |
Under-relaxation factor. 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_ |
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 63 of file SLTSDdtScheme.H.
typedef ddtScheme<Type>::fluxFieldType fluxFieldType |
Definition at line 174 of file SLTSDdtScheme.H.
|
private |
Disallow default bitwise copy construct.
|
inline |
Construct from mesh and Istream.
Definition at line 104 of file SLTSDdtScheme.H.
|
private |
Disallow default bitwise assignment.
|
private |
Calculate a relaxed diagonal from the given flux field.
Definition at line 44 of file SLTSDdtScheme.C.
References Foam::diag(), polyPatch::faceCells(), forAll, mesh, fvPatch::patch(), fvsPatchField::patch(), patchi, and phi.
|
private |
Return the reciprocal of the stabilised local time-step.
Definition at line 90 of file SLTSDdtScheme.C.
References Foam::abort(), GeometricField::correctBoundaryConditions(), Foam::dimless, Foam::dimTime, Foam::FatalError, FatalErrorInFunction, GeometricField::internalField(), Foam::max(), mesh, phi, rho, trDeltaT, and dimensioned::value().
TypeName | ( | "SLTS" | ) |
Runtime type information.
|
inline |
Return mesh reference.
Definition at line 116 of file SLTSDdtScheme.H.
References ddtScheme< Type >::mesh().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 156 of file SLTSDdtScheme.C.
References dimensioned::dimensions(), Foam::dimTime, GeometricField::internalField(), mesh, dimensioned::name(), timeName, and dimensioned::value().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 215 of file SLTSDdtScheme.C.
References GeometricField::boundaryField(), GeometricField::internalField(), mesh, GeometricField::oldTime(), and timeName.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 266 of file SLTSDdtScheme.C.
References GeometricField::boundaryField(), GeometricField::internalField(), mesh, GeometricField::oldTime(), rho, and timeName.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 318 of file SLTSDdtScheme.C.
References GeometricField::boundaryField(), GeometricField::internalField(), mesh, GeometricField::oldTime(), rho, and timeName.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 373 of file SLTSDdtScheme.C.
References Foam::constant::atomic::alpha, GeometricField::boundaryField(), dimensioned::dimensions(), GeometricField::internalField(), mesh, dimensioned::name(), GeometricField::oldTime(), rho, and timeName.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 442 of file SLTSDdtScheme.C.
References Foam::dimTime, Foam::dimVol, Foam::endl(), Foam::gMax(), Foam::gMin(), Foam::Info, internalField(), mesh, GeometricField::oldTime(), and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 480 of file SLTSDdtScheme.C.
References Foam::dimTime, Foam::dimVol, internalField(), mesh, GeometricField::oldTime(), rho, and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 517 of file SLTSDdtScheme.C.
References Foam::dimTime, Foam::dimVol, internalField(), mesh, GeometricField::oldTime(), rho, and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 556 of file SLTSDdtScheme.C.
References Foam::constant::atomic::alpha, dimensioned::dimensions(), Foam::dimTime, Foam::dimVol, internalField(), mesh, GeometricField::oldTime(), rho, and fvMatrix::source().
|
virtual |
Implements ddtScheme< Type >.
Definition at line 598 of file SLTSDdtScheme.C.
References Foam::fvc::interpolate(), mesh, timeName, U, and Uf.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 635 of file SLTSDdtScheme.C.
References Foam::fvc::interpolate(), mesh, phi, timeName, and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 667 of file SLTSDdtScheme.C.
References Foam::abort(), Foam::dimDensity, Foam::dimVelocity, Foam::FatalError, FatalErrorInFunction, Foam::fvc::interpolate(), mesh, rho, timeName, U, and Uf.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 734 of file SLTSDdtScheme.C.
References Foam::abort(), Foam::dimArea, Foam::dimVelocity, Foam::FatalError, FatalErrorInFunction, Foam::fvc::interpolate(), mesh, phi, rho, timeName, and U.
|
virtual |
Implements ddtScheme< Type >.
Definition at line 795 of file SLTSDdtScheme.C.
References Foam::dimTime, Foam::dimVolume, mesh, IOobject::NO_READ, IOobject::NO_WRITE, 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 |
Name of the flux field used to calculate the local time-step.
Definition at line 70 of file SLTSDdtScheme.H.
|
private |
Name of the density field used to obtain the volumetric flux.
from the mass flux if required
Definition at line 74 of file SLTSDdtScheme.H.
|
private |
Under-relaxation factor.
Definition at line 77 of file SLTSDdtScheme.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.