73#include "surfaceInterpolate.H"
111 CoBlended(
const CoBlended&) =
delete;
114 void operator=(
const CoBlended&) =
delete;
134 surfaceInterpolationScheme<Type>(
mesh),
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"
168 surfaceInterpolationScheme<Type>(
mesh),
169 Co1_(readScalar(is)),
172 surfaceInterpolationScheme<Type>::
New(
mesh, faceFlux, is)
174 Co2_(readScalar(is)),
177 surfaceInterpolationScheme<Type>::
New(
mesh, faceFlux, is)
181 if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
184 <<
"coefficients = " << Co1_ <<
" and " << Co2_
185 <<
" should be > 0 and Co2 > Co1"
196 const GeometricField<Type, fvPatchField, volMesh>& vf
200 tmp<surfaceScalarField> tUflux = faceFlux_;
207 mesh.objectRegistry::template lookupObject<volScalarField>
215 <<
"dimensions of faceFlux are not correct"
219 return tmp<surfaceScalarField>
223 vf.name() +
"BlendingFactor",
228 mesh.time().deltaT()*
mesh.deltaCoeffs()
249 bf*tScheme1_().weights(vf)
250 + (scalar(1) - bf)*tScheme2_().weights(vf);
256 tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
259 const GeometricField<Type, fvPatchField, volMesh>& vf
265 bf*tScheme1_().interpolate(vf)
266 + (scalar(1) - bf)*tScheme2_().interpolate(vf);
273 return tScheme1_().corrected() || tScheme2_().corrected();
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.
CoBlended(const fvMesh &mesh, Istream &is)
Construct from mesh 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 noexcept
Return dimensions.
Generic GeometricField class.
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)....
blendedSchemeBase()=default
Constructor.
Mesh data needed to do the Finite Volume discretisation.
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.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
#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.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
errorManipArg< error, int > exit(error &err, const int errNo=1)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.