62 <<
"Discretisation scheme not specified"
64 <<
"Valid schemes are :" <<
nl
65 << MeshConstructorTablePtr_->sortedToc()
69 const word schemeName(schemeData);
71 auto* ctorPtr = MeshConstructorTable(schemeName);
80 *MeshConstructorTablePtr_
84 return ctorPtr(
mesh, schemeData);
96 const word& snGradName
108 snGradName +
"("+vf.
name()+
')',
131 deltaCoeffs[facei]*(vf[neighbour[facei]] - vf[owner[facei]]);
143 ssfbf[patchi] = pvf.
snGrad(tdeltaCoeffs().boundaryField()[patchi]);
147 ssfbf[patchi] = pvf.
snGrad();
160 const word& sndGradName
163 return snGrad(vf, vf.
mesh().nonOrthDeltaCoeffs(), sndGradName);
176 snGrad(vf, deltaCoeffs(vf))
void setOriented(bool on=true) noexcept
Set the oriented flag: on/off.
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic GeometricField class.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
const fileName & instance() const noexcept
Read access to instance path component.
bool eof() const noexcept
True if end of input seen.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Mesh data needed to do the Finite Volume discretisation.
virtual bool coupled() const
True if the patch field is coupled.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > sndGrad(const GeometricField< Type, fvPatchField, volMesh > &, const word &snGradName="sndGrad")
Return the sndGrad of the given cell field.
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.
const fvMesh & mesh() const
Return const reference to mesh.
virtual bool corrected() const noexcept
Return true if this scheme uses an explicit correction.
static tmp< snGradScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &) const =0
Return the interpolation weighting factors for the given field.
A class for managing temporary objects.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
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.
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for finite-volume.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
errorManipArg< error, int > exit(error &err, const int errNo=1)
UList< label > labelUList
A UList of labels.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.