#include <boundedBackwardFaDdtScheme.H>


Public Member Functions | |
| TypeName ("boundedBackward") | |
| Runtime type information. | |
| boundedBackwardFaDdtScheme (const boundedBackwardFaDdtScheme &)=delete | |
| No copy construct. | |
| void | operator= (const boundedBackwardFaDdtScheme &)=delete |
| No copy assignment. | |
| boundedBackwardFaDdtScheme (const faMesh &mesh) | |
| Construct from mesh. | |
| boundedBackwardFaDdtScheme (const faMesh &mesh, Istream &is) | |
| Construct from mesh and Istream. | |
| const faMesh & | mesh () const |
| Return mesh reference. | |
| tmp< areaScalarField > | facDdt (const dimensionedScalar) |
| tmp< areaScalarField > | facDdt0 (const dimensionedScalar) |
| tmp< areaScalarField > | facDdt (const areaScalarField &) |
| tmp< areaScalarField > | facDdt0 (const areaScalarField &) |
| tmp< areaScalarField > | facDdt (const dimensionedScalar &, const areaScalarField &) |
| tmp< areaScalarField > | facDdt0 (const dimensionedScalar &, const areaScalarField &) |
| tmp< areaScalarField > | facDdt (const areaScalarField &, const areaScalarField &) |
| tmp< areaScalarField > | facDdt0 (const areaScalarField &, const areaScalarField &) |
| tmp< faScalarMatrix > | famDdt (const areaScalarField &) |
| tmp< faScalarMatrix > | famDdt (const dimensionedScalar &, const areaScalarField &) |
| tmp< faScalarMatrix > | famDdt (const areaScalarField &, const areaScalarField &) |
| Public Member Functions inherited from faDdtScheme< scalar > | |
| virtual const word & | type () const=0 |
| Runtime type information. | |
| declareRunTimeSelectionTable (tmp, faDdtScheme, Istream,(const faMesh &mesh, Istream &schemeData),(mesh, schemeData)) | |
| faDdtScheme (const faDdtScheme &)=delete | |
| No copy construct. | |
| void | operator= (const faDdtScheme &)=delete |
| No copy assignment. | |
| virtual | ~faDdtScheme () |
| Destructor. | |
| const faMesh & | mesh () const noexcept |
| Return mesh reference. | |
| virtual tmp< GeometricField< scalar, faPatchField, areaMesh > > | facDdt (const dimensioned< scalar >)=0 |
| virtual tmp< GeometricField< scalar, faPatchField, areaMesh > > | facDdt0 (const dimensioned< scalar >)=0 |
| virtual tmp< faMatrix< scalar > > | famDdt (const GeometricField< scalar, faPatchField, areaMesh > &)=0 |
| Public Member Functions inherited from refCount | |
| constexpr | refCount () noexcept |
| Default construct, initializing count to 0. | |
| int | use_count () const noexcept |
| Return the current reference count. | |
| bool | unique () const noexcept |
| Return true if the reference count is zero. | |
| void | operator++ () noexcept |
| Increment the reference count. | |
| void | operator++ (int) noexcept |
| Increment the reference count. | |
| void | operator-- () noexcept |
| Decrement the reference count. | |
| void | operator-- (int) noexcept |
| Decrement the reference count. | |
Additional Inherited Members | |
| Static Public Member Functions inherited from faDdtScheme< scalar > | |
| static tmp< faDdtScheme< scalar > > | New (const faMesh &mesh, Istream &schemeData) |
| Return a pointer to a new faDdtScheme created on freestore. | |
| Protected Attributes inherited from faDdtScheme< scalar > | |
| const faMesh & | mesh_ |
| Reference to mesh. | |
Definition at line 53 of file boundedBackwardFaDdtScheme.H.
|
delete |
No copy construct.
References boundedBackwardFaDdtScheme().
Referenced by boundedBackwardFaDdtScheme(), and operator=().


|
inline |
Construct from mesh.
Definition at line 112 of file boundedBackwardFaDdtScheme.H.
References faDdtScheme< scalar >::faDdtScheme(), and mesh().

Construct from mesh and Istream.
Definition at line 120 of file boundedBackwardFaDdtScheme.H.
References faDdtScheme< scalar >::faDdtScheme(), and mesh().

| TypeName | ( | "boundedBackward" | ) |
Runtime type information.
|
delete |
No copy assignment.
References boundedBackwardFaDdtScheme().

|
inline |
Return mesh reference.
Definition at line 131 of file boundedBackwardFaDdtScheme.H.
References faDdtScheme< Type >::mesh().
Referenced by boundedBackwardFaDdtScheme(), boundedBackwardFaDdtScheme(), facDdt(), facDdt(), facDdt(), facDdt(), facDdt0(), facDdt0(), facDdt0(), facDdt0(), famDdt(), famDdt(), and famDdt().


| tmp< areaScalarField > facDdt | ( | const dimensionedScalar | dt | ) |
Definition at line 50 of file boundedBackwardFaDdtScheme.C.
References faPatchFieldBase::calculatedType(), TimeState::deltaT(), dimensioned< Type >::dimensions(), Foam::dimTime, mesh(), dimensioned< Type >::name(), tmp< T >::New(), tmp< T >::ref(), faMesh::S(), faMesh::S0(), faMesh::S00(), faMesh::time(), dimensioned< Type >::value(), and Foam::Zero.

| tmp< areaScalarField > facDdt0 | ( | const dimensionedScalar | dt | ) |
Definition at line 105 of file boundedBackwardFaDdtScheme.C.
References TimeState::deltaT(), mesh(), dimensioned< Type >::name(), tmp< T >::ref(), faMesh::S(), faMesh::S0(), faMesh::S00(), faMesh::time(), and dimensioned< Type >::value().

| tmp< areaScalarField > facDdt | ( | const areaScalarField & | vf | ) |
Definition at line 152 of file boundedBackwardFaDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), TimeState::deltaT(), dimensioned< Type >::dimensions(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::limiter(), Foam::mag(), mesh(), IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), Foam::pos(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), Foam::sqr(), faMesh::time(), and dimensioned< Type >::value().

| tmp< areaScalarField > facDdt0 | ( | const areaScalarField & | vf | ) |
Definition at line 249 of file boundedBackwardFaDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), TimeState::deltaT(), dimensioned< Type >::dimensions(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::limiter(), Foam::mag(), mesh(), IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), Foam::pos(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), Foam::sqr(), faMesh::time(), and dimensioned< Type >::value().

| tmp< areaScalarField > facDdt | ( | const dimensionedScalar & | rho, |
| const areaScalarField & | vf ) |
Definition at line 343 of file boundedBackwardFaDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), TimeState::deltaT(), dimensioned< Type >::dimensions(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::limiter(), Foam::mag(), mesh(), IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), Foam::pos(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), rho, Foam::sqr(), faMesh::time(), and dimensioned< Type >::value().

| tmp< areaScalarField > facDdt0 | ( | const dimensionedScalar & | rho, |
| const areaScalarField & | vf ) |
Definition at line 440 of file boundedBackwardFaDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), TimeState::deltaT(), dimensioned< Type >::dimensions(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::limiter(), Foam::mag(), mesh(), IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), Foam::pos(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), rho, Foam::sqr(), faMesh::time(), and dimensioned< Type >::value().

| tmp< areaScalarField > facDdt | ( | const areaScalarField & | rho, |
| const areaScalarField & | vf ) |
Definition at line 535 of file boundedBackwardFaDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), TimeState::deltaT(), dimensioned< Type >::dimensions(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::limiter(), Foam::mag(), mesh(), IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), Foam::pos(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), rho, Foam::sqr(), faMesh::time(), and dimensioned< Type >::value().

| tmp< areaScalarField > facDdt0 | ( | const areaScalarField & | rho, |
| const areaScalarField & | vf ) |
Definition at line 636 of file boundedBackwardFaDdtScheme.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), TimeState::deltaT(), dimensioned< Type >::dimensions(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::limiter(), Foam::mag(), mesh(), IOobject::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), Foam::pos(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), rho, Foam::sqr(), faMesh::time(), and dimensioned< Type >::value().

| tmp< faScalarMatrix > famDdt | ( | const areaScalarField & | vf | ) |
Definition at line 734 of file boundedBackwardFaDdtScheme.C.
References Foam::dimArea, DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, GeometricField< Type, PatchField, GeoMesh >::internalField(), Foam::limiter(), Foam::mag(), mesh(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), Foam::pos(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), tmp< T >::ref(), and faMesh::S().

| tmp< faScalarMatrix > famDdt | ( | const dimensionedScalar & | rho, |
| const areaScalarField & | vf ) |
Definition at line 805 of file boundedBackwardFaDdtScheme.C.
References Foam::dimArea, DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, GeometricField< Type, PatchField, GeoMesh >::internalField(), Foam::limiter(), Foam::mag(), mesh(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), Foam::pos(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), tmp< T >::ref(), rho, and faMesh::S().

| tmp< faScalarMatrix > famDdt | ( | const areaScalarField & | rho, |
| const areaScalarField & | vf ) |
Definition at line 876 of file boundedBackwardFaDdtScheme.C.
References Foam::dimArea, DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimTime, GeometricField< Type, PatchField, GeoMesh >::internalField(), Foam::limiter(), Foam::mag(), mesh(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), Foam::pos(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), tmp< T >::ref(), rho, and faMesh::S().
