CoEulerDdtScheme< Type > Class Template Reference

Courant number limited first-order Euler implicit/explicit ddt. More...

Inheritance diagram for CoEulerDdtScheme< Type >:
Collaboration diagram for CoEulerDdtScheme< Type >:

Public Types

typedef ddtScheme< Type >::fluxFieldType fluxFieldType
 
- Public Types inherited from ddtScheme< Type >
typedef SurfaceField< typename flux< Type >::typefluxFieldType
 

Public Member Functions

 TypeName ("CoEuler")
 Runtime type information. More...
 
 CoEulerDdtScheme (const fvMesh &mesh, Istream &is)
 Construct from mesh and Istream. More...
 
 CoEulerDdtScheme (const CoEulerDdtScheme &)=delete
 Disallow default bitwise copy construction. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
virtual tmp< VolField< Type > > fvcDdt (const dimensioned< Type > &)
 
virtual tmp< VolField< Type > > fvcDdt (const VolField< Type > &)
 
virtual tmp< VolField< Type > > fvcDdt (const dimensionedScalar &, const VolField< Type > &)
 
virtual tmp< VolField< Type > > fvcDdt (const volScalarField &, const VolField< Type > &)
 
virtual tmp< VolField< Type > > fvcDdt (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &psi)
 
virtual tmp< fvMatrix< Type > > fvmDdt (const VolField< Type > &)
 
virtual tmp< fvMatrix< Type > > fvmDdt (const dimensionedScalar &, const VolField< Type > &)
 
virtual tmp< fvMatrix< Type > > fvmDdt (const volScalarField &, const VolField< Type > &)
 
virtual tmp< fvMatrix< Type > > fvmDdt (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &psi)
 
virtual tmp< fluxFieldTypefvcDdtUfCorr (const VolField< Type > &U, const SurfaceField< Type > &Uf)
 
virtual tmp< fluxFieldTypefvcDdtPhiCorr (const VolField< Type > &U, const fluxFieldType &phi)
 
virtual tmp< fluxFieldTypefvcDdtUfCorr (const volScalarField &rho, const VolField< Type > &U, const SurfaceField< Type > &rhoUf)
 
virtual tmp< fluxFieldTypefvcDdtPhiCorr (const volScalarField &rho, const VolField< Type > &U, const fluxFieldType &phi)
 
virtual tmp< fluxFieldTypefvcDdtUfCorr (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &U, const SurfaceField< Type > &Uf)
 
virtual tmp< fluxFieldTypefvcDdtPhiCorr (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &U, const fluxFieldType &phi)
 
virtual tmp< surfaceScalarFieldmeshPhi (const VolField< Type > &)
 
virtual tmp< scalarFieldmeshPhi (const VolField< Type > &, const label patchi)
 
void operator= (const CoEulerDdtScheme &)=delete
 Disallow default bitwise assignment. More...
 
tmp< surfaceScalarFieldfvcDdtUfCorr (const VolField< scalar > &U, const SurfaceField< scalar > &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)
 
tmp< surfaceScalarFieldfvcDdtUfCorr (const volScalarField &alpha, const volScalarField &rho, const volScalarField &U, const surfaceScalarField &Uf)
 
tmp< surfaceScalarFieldfvcDdtPhiCorr (const volScalarField &alpha, 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...
 
 ddtScheme (const ddtScheme &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~ddtScheme ()
 Destructor. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
virtual tmp< SurfaceField< Type > > fvcDdt (const SurfaceField< Type > &)
 
virtual tmp< surfaceScalarFieldfvcDdtPhiCoeff (const VolField< Type > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
 
virtual tmp< surfaceScalarFieldfvcDdtPhiCoeff (const VolField< Type > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr, const volScalarField &rho)
 
virtual tmp< surfaceScalarFieldfvcDdtPhiCoeff (const VolField< Type > &U, const fluxFieldType &phi)
 
virtual tmp< surfaceScalarFieldfvcDdtPhiCoeff (const VolField< Type > &U, const fluxFieldType &phi, const volScalarField &rho)
 
void operator= (const ddtScheme &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from refCount
int count () const
 Return the current reference count. More...
 
bool unique () const
 Return true if the reference count is 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...
 
void operator= (const refCount &)=delete
 Disallow bitwise assignment. 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 refCount
 refCount ()
 Construct null initialising count to 0. More...
 
 refCount (const refCount &)=delete
 Disallow copy. More...
 
- Protected Attributes inherited from ddtScheme< Type >
const fvMeshmesh_
 

Detailed Description

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

Courant number limited first-order Euler implicit/explicit ddt.

The time-step is adjusted locally so that the local Courant number does not exceed the specified value.

This scheme should only be used for steady-state computations using transient codes where local time-stepping is preferable to under-relaxation for transport consistency reasons.

Source files

Definition at line 61 of file CoEulerDdtScheme.H.

Member Typedef Documentation

◆ fluxFieldType

Definition at line 169 of file CoEulerDdtScheme.H.

Constructor & Destructor Documentation

◆ CoEulerDdtScheme() [1/2]

CoEulerDdtScheme ( const fvMesh mesh,
Istream is 
)
inline

Construct from mesh and Istream.

Definition at line 96 of file CoEulerDdtScheme.H.

◆ CoEulerDdtScheme() [2/2]

CoEulerDdtScheme ( const CoEulerDdtScheme< Type > &  )
delete

Disallow default bitwise copy construction.

Member Function Documentation

◆ TypeName()

TypeName ( "CoEuler"  )

Runtime type information.

◆ mesh()

const fvMesh& mesh ( ) const
inline

Return mesh reference.

Definition at line 111 of file CoEulerDdtScheme.H.

References ddtScheme< Type >::mesh().

Here is the call graph for this function:

◆ fvcDdt() [1/5]

◆ fvcDdt() [2/5]

tmp< VolField< Type > > fvcDdt ( const VolField< Type > &  vf)
virtual

◆ fvcDdt() [3/5]

tmp< VolField< Type > > fvcDdt ( const dimensionedScalar rho,
const VolField< Type > &  vf 
)
virtual

◆ fvcDdt() [4/5]

tmp< VolField< Type > > fvcDdt ( const volScalarField rho,
const VolField< Type > &  vf 
)
virtual

◆ fvcDdt() [5/5]

tmp< VolField< Type > > fvcDdt ( const volScalarField alpha,
const volScalarField rho,
const VolField< Type > &  psi 
)
virtual

◆ fvmDdt() [1/4]

tmp< fvMatrix< Type > > fvmDdt ( const VolField< Type > &  vf)
virtual

◆ fvmDdt() [2/4]

tmp< fvMatrix< Type > > fvmDdt ( const dimensionedScalar rho,
const VolField< Type > &  vf 
)
virtual

◆ fvmDdt() [3/4]

tmp< fvMatrix< Type > > fvmDdt ( const volScalarField rho,
const VolField< Type > &  vf 
)
virtual

◆ fvmDdt() [4/4]

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

◆ fvcDdtUfCorr() [1/6]

tmp< typename CoEulerDdtScheme< Type >::fluxFieldType > fvcDdtUfCorr ( const VolField< Type > &  U,
const SurfaceField< Type > &  Uf 
)
virtual

Implements ddtScheme< Type >.

Definition at line 557 of file CoEulerDdtScheme.C.

References Foam::fvc::dotInterpolate(), Foam::fvc::interpolate(), IOobject::name(), Foam::compressible::New(), OldTimeField< FieldType >::oldTime(), and U.

Here is the call graph for this function:

◆ fvcDdtPhiCorr() [1/6]

tmp< typename CoEulerDdtScheme< Type >::fluxFieldType > fvcDdtPhiCorr ( const VolField< Type > &  U,
const fluxFieldType phi 
)
virtual

Implements ddtScheme< Type >.

Definition at line 581 of file CoEulerDdtScheme.C.

References Foam::fvc::dotInterpolate(), Foam::fvc::interpolate(), IOobject::name(), Foam::compressible::New(), OldTimeField< FieldType >::oldTime(), and U.

Here is the call graph for this function:

◆ fvcDdtUfCorr() [2/6]

◆ fvcDdtPhiCorr() [2/6]

◆ fvcDdtUfCorr() [3/6]

virtual tmp<fluxFieldType> fvcDdtUfCorr ( const volScalarField alpha,
const volScalarField rho,
const VolField< Type > &  U,
const SurfaceField< Type > &  Uf 
)
inlinevirtual

Implements ddtScheme< Type >.

Definition at line 197 of file CoEulerDdtScheme.H.

References NotImplemented, and GeometricField< Type, PatchField, GeoMesh >::null().

Here is the call graph for this function:

◆ fvcDdtPhiCorr() [3/6]

virtual tmp<fluxFieldType> fvcDdtPhiCorr ( const volScalarField alpha,
const volScalarField rho,
const VolField< Type > &  U,
const fluxFieldType phi 
)
inlinevirtual

Implements ddtScheme< Type >.

Definition at line 209 of file CoEulerDdtScheme.H.

References NotImplemented, and GeometricField< Type, PatchField, GeoMesh >::null().

Here is the call graph for this function:

◆ meshPhi() [1/2]

tmp< surfaceScalarField > meshPhi ( const VolField< Type > &  )
virtual

Implements ddtScheme< Type >.

Definition at line 744 of file CoEulerDdtScheme.C.

References Foam::dimTime, Foam::dimVolume, and GeometricField< Type, PatchField, GeoMesh >::New().

Here is the call graph for this function:

◆ meshPhi() [2/2]

tmp< scalarField > meshPhi ( const VolField< Type > &  ,
const label  patchi 
)
virtual

Implements ddtScheme< Type >.

Definition at line 759 of file CoEulerDdtScheme.C.

References boundary(), and patchi.

Here is the call graph for this function:

◆ operator=()

void operator= ( const CoEulerDdtScheme< Type > &  )
delete

Disallow default bitwise assignment.

◆ fvcDdtUfCorr() [4/6]

tmp< surfaceScalarField > fvcDdtUfCorr ( const VolField< scalar > &  U,
const SurfaceField< scalar > &  Uf 
)

◆ fvcDdtPhiCorr() [4/6]

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

◆ fvcDdtUfCorr() [5/6]

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

◆ fvcDdtPhiCorr() [5/6]

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

◆ fvcDdtUfCorr() [6/6]

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

◆ fvcDdtPhiCorr() [6/6]

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

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