Public Types | Public Member Functions | Private Member Functions | Private Attributes
SLTSDdtScheme< Type > Class Template Reference

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...

Inheritance diagram for SLTSDdtScheme< Type >:
Inheritance graph
[legend]
Collaboration diagram for SLTSDdtScheme< Type >:
Collaboration graph
[legend]

Public Types

typedef ddtScheme< Type >::fluxFieldType fluxFieldType
 
- Public Types inherited from ddtScheme< Type >
typedef GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMeshfluxFieldType
 

Public Member Functions

 TypeName ("SLTS")
 Runtime type information. More...
 
 SLTSDdtScheme (const fvMesh &mesh, Istream &is)
 Construct from mesh and Istream. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const dimensioned< Type > &)
 
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const GeometricField< Type, fvPatchField, volMesh > &)
 
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &)
 
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const volScalarField &, const GeometricField< Type, fvPatchField, volMesh > &)
 
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt (const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &psi)
 
tmp< fvMatrix< Type > > fvmDdt (const GeometricField< Type, fvPatchField, volMesh > &)
 
tmp< fvMatrix< Type > > fvmDdt (const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &)
 
tmp< fvMatrix< Type > > fvmDdt (const volScalarField &, const GeometricField< Type, fvPatchField, volMesh > &)
 
tmp< fvMatrix< Type > > fvmDdt (const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &psi)
 
tmp< fluxFieldTypefvcDdtUfCorr (const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
 
tmp< fluxFieldTypefvcDdtPhiCorr (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
 
tmp< fluxFieldTypefvcDdtUfCorr (const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
 
tmp< fluxFieldTypefvcDdtPhiCorr (const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
 
tmp< surfaceScalarFieldmeshPhi (const GeometricField< Type, fvPatchField, volMesh > &)
 
tmp< surfaceScalarFieldfvcDdtUfCorr (const GeometricField< scalar, fvPatchField, volMesh > &U, const GeometricField< scalar, fvsPatchField, surfaceMesh > &Uf)
 
tmp< surfaceScalarFieldfvcDdtPhiCorr (const volScalarField &U, const surfaceScalarField &phi)
 
tmp< surfaceScalarFieldfvcDdtUfCorr (const volScalarField &rho, const volScalarField &U, const surfaceScalarField &Uf)
 
tmp< surfaceScalarFieldfvcDdtPhiCorr (const volScalarField &rho, const volScalarField &U, const surfaceScalarField &phi)
 
- Public Member Functions inherited from ddtScheme< Type >
virtual const wordtype () const =0
 Runtime type information. More...
 
 declareRunTimeSelectionTable (tmp, ddtScheme, Istream,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
 
 ddtScheme (const fvMesh &mesh)
 Construct from mesh. More...
 
 ddtScheme (const fvMesh &mesh, Istream &)
 Construct from mesh and Istream. More...
 
virtual ~ddtScheme ()
 Destructor. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
 
tmp< surfaceScalarFieldfvcDdtPhiCoeff (const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
 
- Public Member Functions inherited from refCount
 refCount ()
 Construct null with zero count. More...
 
int count () const
 Return the reference count. More...
 
bool okToDelete () const
 Return true if the reference count is zero. More...
 
void resetRefCount ()
 Reset the reference count to zero. More...
 
void operator++ ()
 Increment the reference count. More...
 
void operator++ (int)
 Increment the reference count. More...
 
void operator-- ()
 Decrement the reference count. More...
 
void operator-- (int)
 Decrement the reference count. More...
 

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< volScalarFieldSLrDeltaT () 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 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. More...
 
- Protected Member Functions inherited from ddtScheme< Type >
 ddtScheme (const ddtScheme &)
 Disallow copy construct. More...
 
void operator= (const ddtScheme &)
 Disallow default bitwise assignment. More...
 
- Protected Attributes inherited from ddtScheme< Type >
const fvMeshmesh_
 

Detailed Description

template<class Type>
class Foam::fv::SLTSDdtScheme< Type >

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.

See also
Foam::fv::CoEulerDdtScheme
Source files

Definition at line 63 of file SLTSDdtScheme.H.

Member Typedef Documentation

◆ fluxFieldType

Definition at line 174 of file SLTSDdtScheme.H.

Constructor & Destructor Documentation

◆ SLTSDdtScheme() [1/2]

SLTSDdtScheme ( const SLTSDdtScheme< Type > &  )
private

Disallow default bitwise copy construct.

◆ SLTSDdtScheme() [2/2]

SLTSDdtScheme ( const fvMesh mesh,
Istream is 
)
inline

Construct from mesh and Istream.

Definition at line 104 of file SLTSDdtScheme.H.

Member Function Documentation

◆ operator=()

void operator= ( const SLTSDdtScheme< Type > &  )
private

Disallow default bitwise assignment.

◆ relaxedDiag()

void relaxedDiag ( scalarField rD,
const surfaceScalarField phi 
) const
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.

Here is the call graph for this function:

◆ SLrDeltaT()

tmp< volScalarField > SLrDeltaT
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().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "SLTS"  )

Runtime type information.

◆ mesh()

const fvMesh& mesh ( ) const
inline

Return mesh reference.

Definition at line 116 of file SLTSDdtScheme.H.

References ddtScheme< Type >::mesh().

Here is the call graph for this function:

◆ fvcDdt() [1/5]

tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt ( const dimensioned< Type > &  dt)
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().

Here is the call graph for this function:

◆ fvcDdt() [2/5]

tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt ( const GeometricField< Type, fvPatchField, volMesh > &  vf)
virtual

Implements ddtScheme< Type >.

Definition at line 215 of file SLTSDdtScheme.C.

References GeometricField::boundaryField(), GeometricField::internalField(), mesh, GeometricField::oldTime(), and timeName.

Here is the call graph for this function:

◆ fvcDdt() [3/5]

tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt ( const dimensionedScalar rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)
virtual

Implements ddtScheme< Type >.

Definition at line 266 of file SLTSDdtScheme.C.

References GeometricField::boundaryField(), GeometricField::internalField(), mesh, GeometricField::oldTime(), rho, and timeName.

Here is the call graph for this function:

◆ fvcDdt() [4/5]

tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt ( const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)
virtual

Implements ddtScheme< Type >.

Definition at line 318 of file SLTSDdtScheme.C.

References GeometricField::boundaryField(), GeometricField::internalField(), mesh, GeometricField::oldTime(), rho, and timeName.

Here is the call graph for this function:

◆ fvcDdt() [5/5]

tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt ( const volScalarField alpha,
const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  psi 
)
virtual

◆ fvmDdt() [1/4]

tmp< fvMatrix< Type > > fvmDdt ( const GeometricField< Type, fvPatchField, volMesh > &  vf)
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().

Here is the call graph for this function:

◆ fvmDdt() [2/4]

tmp< fvMatrix< Type > > fvmDdt ( const dimensionedScalar rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)
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().

Here is the call graph for this function:

◆ fvmDdt() [3/4]

tmp< fvMatrix< Type > > fvmDdt ( const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)
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().

Here is the call graph for this function:

◆ fvmDdt() [4/4]

tmp< fvMatrix< Type > > fvmDdt ( const volScalarField alpha,
const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  psi 
)
virtual

◆ fvcDdtUfCorr() [1/4]

tmp< typename SLTSDdtScheme< Type >::fluxFieldType > fvcDdtUfCorr ( const GeometricField< Type, fvPatchField, volMesh > &  U,
const GeometricField< Type, fvsPatchField, surfaceMesh > &  Uf 
)
virtual

Implements ddtScheme< Type >.

Definition at line 598 of file SLTSDdtScheme.C.

References Foam::fvc::interpolate(), mesh, timeName, U, and Uf.

Here is the call graph for this function:

◆ fvcDdtPhiCorr() [1/4]

tmp< typename SLTSDdtScheme< Type >::fluxFieldType > fvcDdtPhiCorr ( const GeometricField< Type, fvPatchField, volMesh > &  U,
const fluxFieldType phi 
)
virtual

Implements ddtScheme< Type >.

Definition at line 635 of file SLTSDdtScheme.C.

References Foam::fvc::interpolate(), mesh, phi, timeName, and U.

Here is the call graph for this function:

◆ fvcDdtUfCorr() [2/4]

tmp< typename SLTSDdtScheme< Type >::fluxFieldType > fvcDdtUfCorr ( const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  U,
const GeometricField< Type, fvsPatchField, surfaceMesh > &  Uf 
)
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.

Here is the call graph for this function:

◆ fvcDdtPhiCorr() [2/4]

tmp< typename SLTSDdtScheme< Type >::fluxFieldType > fvcDdtPhiCorr ( const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  U,
const fluxFieldType phi 
)
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.

Here is the call graph for this function:

◆ meshPhi()

tmp< surfaceScalarField > meshPhi ( const GeometricField< Type, fvPatchField, volMesh > &  )
virtual

◆ fvcDdtUfCorr() [3/4]

tmp< surfaceScalarField > fvcDdtUfCorr ( const GeometricField< scalar, fvPatchField, volMesh > &  U,
const GeometricField< scalar, fvsPatchField, surfaceMesh > &  Uf 
)

◆ fvcDdtPhiCorr() [3/4]

tmp< surfaceScalarField > fvcDdtPhiCorr ( const volScalarField U,
const surfaceScalarField phi 
)

◆ fvcDdtUfCorr() [4/4]

tmp< surfaceScalarField > fvcDdtUfCorr ( const volScalarField rho,
const volScalarField U,
const surfaceScalarField Uf 
)

◆ fvcDdtPhiCorr() [4/4]

tmp< surfaceScalarField > fvcDdtPhiCorr ( const volScalarField rho,
const volScalarField U,
const surfaceScalarField phi 
)

Field Documentation

◆ phiName_

word phiName_
private

Name of the flux field used to calculate the local time-step.

Definition at line 70 of file SLTSDdtScheme.H.

◆ rhoName_

word rhoName_
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.

◆ alpha_

scalar alpha_
private

Under-relaxation factor.

Definition at line 77 of file SLTSDdtScheme.H.


The documentation for this class was generated from the following files: