Go to the documentation of this file.
40 #ifndef outletStabilised_H
41 #define outletStabilised_H
154 forAll(pFaceCells, pFacei)
156 const cell& pFaceCell =
cells[pFaceCells[pFacei]];
160 label facei = pFaceCell[fi];
213 forAll(pFaceCells, pFacei)
215 const cell& pFaceCell =
cells[pFaceCells[pFacei]];
219 label facei = pFaceCell[fi];
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
A class for handling words, derived from string.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const cellList & cells() const
TypeName("outletStabilised")
Runtime type information.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
void operator=(const outletStabilised &)
Disallow default bitwise assignment.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
tmp< surfaceInterpolationScheme< Type > > tScheme_
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
const fvMesh & mesh_
Hold reference to mesh.
Mesh data needed to do the Finite Volume discretisation.
Outlet-stabilised interpolation scheme which applies upwind differencing to the faces of the cells ad...
This boundary condition provides a base class for 'mixed' type boundary conditions,...
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
const surfaceScalarField & faceFlux_
outletStabilised(const outletStabilised &)
Disallow default bitwise copy construct.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Traits class for primitives.
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Abstract base class for surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Generic GeometricField class.
A cell is defined as a list of faces with extra functionality.
Base class for direction-mixed boundary conditions.
dimensionedScalar pos(const dimensionedScalar &ds)