45#ifndef outletStabilised_H
46#define outletStabilised_H
68 public surfaceInterpolationScheme<Type>
79 outletStabilised(
const outletStabilised&) =
delete;
82 void operator=(
const outletStabilised&) =
delete;
100 surfaceInterpolationScheme<Type>(
mesh),
110 surfaceInterpolationScheme<Type>::
New(
mesh, faceFlux_, is)
123 surfaceInterpolationScheme<Type>(
mesh),
127 surfaceInterpolationScheme<Type>::
New(
mesh, faceFlux, is)
160 : mesh_.boundary()[patchi].faceCells()
163 for (
const label facei :
cells[celli])
165 if (mesh_.isInternalFace(facei))
168 w[facei] =
pos0(faceFlux_[facei]);
181 return tScheme_().corrected();
186 virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
195 tScheme_().correction(vf);
216 : mesh_.boundary()[patchi].faceCells()
219 for (
const label facei :
cells[celli])
221 if (mesh_.isInternalFace(facei))
Generic GeometricField class.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Base class for direction-mixed boundary conditions.
Mesh data needed to do the Finite Volume discretisation.
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Outlet-stabilised interpolation scheme which applies upwind differencing to the faces of the cells ad...
outletStabilised(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
TypeName("outletStabilised")
Runtime type information.
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
outletStabilised(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
const fvMesh & mesh() const
Return mesh reference.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
A class for managing temporary objects.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
A class for handling words, derived from Foam::string.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
dimensionedScalar pos0(const dimensionedScalar &ds)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
List< cell > cellList
List of cell.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
static constexpr const zero Zero
Global zero (0).
#define forAll(list, i)
Loop across all elements in list.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.