Surface gradient scheme with full explicit non-orthogonal correction. More...
#include <correctedSnGrad.H>


Public Member Functions | |
| TypeName ("corrected") | |
| Runtime type information. | |
| correctedSnGrad (const fvMesh &mesh) | |
| Construct from mesh. | |
| correctedSnGrad (const fvMesh &mesh, Istream &) | |
| Construct from mesh and data stream. | |
| virtual | ~correctedSnGrad ()=default |
| Destructor. | |
| virtual tmp< surfaceScalarField > | deltaCoeffs (const GeometricField< Type, fvPatchField, volMesh > &) const |
| Return the interpolation weighting factors for the given field. | |
| virtual bool | corrected () const noexcept |
| Return true if this scheme uses an explicit correction. | |
| tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > | fullGradCorrection (const GeometricField< Type, fvPatchField, volMesh > &) const |
| Return the explicit correction to the correctedSnGrad for the given field using the gradient of the field. | |
| virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > | correction (const GeometricField< Type, fvPatchField, volMesh > &) const |
| Return the explicit correction to the correctedSnGrad for the given field using the gradients of the field components. | |
| tmp< surfaceScalarField > | correction (const volScalarField &vsf) const |
| tmp< surfaceVectorField > | correction (const volVectorField &vvf) const |
| Foam::tmp< Foam::surfaceScalarField > | correction (const volScalarField &vsf) const |
| Foam::tmp< Foam::surfaceVectorField > | correction (const volVectorField &vvf) const |
| Public Member Functions inherited from snGradScheme< Type > | |
| virtual const word & | type () const =0 |
| Runtime type information. | |
| declareRunTimeSelectionTable (tmp, snGradScheme, Mesh,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData)) | |
| snGradScheme (const fvMesh &mesh) | |
| Construct from mesh. | |
| virtual | ~snGradScheme ()=default |
| Destructor. | |
| const fvMesh & | mesh () const |
| Return const reference to mesh. | |
| virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > | snGrad (const GeometricField< Type, fvPatchField, volMesh > &) const |
| Return the snGrad of the given cell field with explicit correction. | |
| tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > | snGrad (const tmp< GeometricField< Type, fvPatchField, volMesh > > &) const |
| Return the snGrad of the given tmp cell field with explicit correction. | |
| 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 snGradScheme< Type > | |
| static tmp< snGradScheme< Type > > | New (const fvMesh &mesh, Istream &schemeData) |
| Return new tmp interpolation scheme. | |
| static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > | snGrad (const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const word &snGradName="snGrad") |
| Return the snGrad of the given cell field by using the given deltaCoeffs. | |
| static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > | sndGrad (const GeometricField< Type, fvPatchField, volMesh > &, const word &snGradName="sndGrad") |
| Return the sndGrad of the given cell field. | |
Surface gradient scheme with full explicit non-orthogonal correction.
system/fvSchemes: snGradSchemes
{
snGrad(<term>) corrected;
}
Definition at line 65 of file correctedSnGrad.H.
|
inline |
Construct from mesh.
Definition at line 90 of file correctedSnGrad.H.
References snGradScheme< Type >::mesh().
Referenced by correction().


Construct from mesh and data stream.
Definition at line 98 of file correctedSnGrad.H.
References snGradScheme< Type >::mesh().

|
virtualdefault |
Destructor.
| TypeName | ( | "corrected" | ) |
Runtime type information.
|
inlinevirtual |
Return the interpolation weighting factors for the given field.
Implements snGradScheme< Type >.
Definition at line 115 of file correctedSnGrad.H.
References snGradScheme< Type >::mesh(), and surfaceInterpolation::nonOrthDeltaCoeffs().

|
inlinevirtualnoexcept |
Return true if this scheme uses an explicit correction.
Reimplemented from snGradScheme< Type >.
Definition at line 126 of file correctedSnGrad.H.
References Foam::noexcept.
| Foam::tmp< Foam::GeometricField< Type, Foam::fvsPatchField, Foam::surfaceMesh > > fullGradCorrection | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) | const |
Return the explicit correction to the correctedSnGrad for the given field using the gradient of the field.
Definition at line 33 of file correctedSnGrad.C.
References snGradScheme< Type >::mesh(), mesh, IOobject::name(), gradScheme< Type >::New(), and tmp< T >::ref().
Referenced by correction(), correction(), and correction().


|
virtual |
Return the explicit correction to the correctedSnGrad for the given field using the gradients of the field components.
Reimplemented from snGradScheme< Type >.
Definition at line 59 of file correctedSnGrad.C.
References GeometricField< Type, PatchField, GeoMesh >::component(), correctedSnGrad(), DimensionedField< Type, GeoMesh >::dimensions(), fullGradCorrection(), IOobject::instance(), snGradScheme< Type >::mesh(), mesh, IOobject::name(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, tmp< T >::ref(), GeometricField< Type, PatchField, GeoMesh >::replace(), and DimensionedField< Type, GeoMesh >::setOriented().

| tmp< surfaceScalarField > correction | ( | const volScalarField & | vsf | ) | const |
| tmp< surfaceVectorField > correction | ( | const volVectorField & | vvf | ) | const |
| Foam::tmp< Foam::surfaceScalarField > correction | ( | const volScalarField & | vsf | ) | const |
Definition at line 30 of file correctedSnGrads.C.
References fullGradCorrection().

| Foam::tmp< Foam::surfaceVectorField > correction | ( | const volVectorField & | vvf | ) | const |
Definition at line 41 of file correctedSnGrads.C.
References fullGradCorrection().
