73#include "surfaceInterpolate.H"
87 public surfaceInterpolationScheme<Type>,
88 public blendedSchemeBase<Type>
114 void operator=(
const CoBlended&) =
delete;
135 Co1_(readScalar(is)),
140 Co2_(readScalar(is)),
150 if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
153 <<
"coefficients = " << Co1_ <<
" and " << Co2_
154 <<
" should be > 0 and Co2 > Co1"
169 Co1_(readScalar(is)),
174 Co2_(readScalar(is)),
181 if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
184 <<
"coefficients = " << Co1_ <<
" and " << Co2_
185 <<
" should be > 0 and Co2 > Co1"
207 mesh.objectRegistry::template lookupObject<volScalarField>
215 <<
"dimensions of faceFlux are not correct"
219 return tmp<surfaceScalarField>
223 vf.
name() +
"BlendingFactor",
244 tmp<surfaceScalarField>
247 const GeometricField<Type, fvPatchField, volMesh>& vf
253 bf*tScheme1_().weights(vf)
254 + (scalar(1) - bf)*tScheme2_().weights(vf);
269 bf*tScheme1_().interpolate(vf)
270 + (scalar(1) - bf)*tScheme2_().interpolate(vf);
277 return tScheme1_().corrected() || tScheme2_().corrected();
322 return tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
Two-scheme Courant number based blending differencing scheme.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
TypeName("CoBlended")
Runtime type information.
CoBlended(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
const dimensionSet & dimensions() const
Return dimensions.
const word & name() const noexcept
Return the object name.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
dimensionedScalar deltaT() const
Return time step.
Base class for blended schemes to provide access to the blending factor surface field.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Abstract base class for surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
virtual const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimArea(sqr(dimLength))
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.