75#ifndef cellCoBlended_H
76#define cellCoBlended_H
80#include "surfaceInterpolate.H"
96 public surfaceInterpolationScheme<Type>,
120 cellCoBlended(
const cellCoBlended&) =
delete;
123 void operator=(
const cellCoBlended&) =
delete;
143 surfaceInterpolationScheme<Type>(
mesh),
144 Co1_(readScalar(is)),
147 surfaceInterpolationScheme<Type>::
New(
mesh, is)
149 Co2_(readScalar(is)),
152 surfaceInterpolationScheme<Type>::
New(
mesh, is)
159 if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
162 <<
"coefficients = " << Co1_ <<
" and " << Co2_
163 <<
" should be > 0 and Co2 > Co1"
178 Co1_(readScalar(is)),
183 Co2_(readScalar(is)),
186 surfaceInterpolationScheme<Type>::
New(
mesh, faceFlux, is)
190 if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
193 <<
"coefficients = " << Co1_ <<
" and " << Co2_
194 <<
" should be > 0 and Co2 > Co1"
216 mesh.objectRegistry::template lookupObject<volScalarField>
224 <<
"dimensions of faceFlux are not correct"
233 mesh.time().timeName(),
246 Co.primitiveFieldRef() =
248 Co.correctBoundaryConditions();
254 vf.
name() +
"BlendingFactor",
267 tmp<surfaceScalarField>
270 const GeometricField<Type, fvPatchField, volMesh>& vf
276 bf*tScheme1_().weights(vf)
277 + (scalar(1) - bf)*tScheme2_().weights(vf);
283 tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
286 const GeometricField<Type, fvPatchField, volMesh>& vf
292 bf*tScheme1_().interpolate(vf)
293 + (scalar(1) - bf)*tScheme2_().interpolate(vf);
300 return tScheme1_().corrected() || tScheme2_().corrected();
306 virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
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.
Two-scheme cell-based 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.
cellCoBlended(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.
TypeName("cellCoBlended")
Runtime type information.
cellCoBlended(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
Mesh data needed to do the Finite Volume discretisation.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
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.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
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.
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimless
Dimensionless.
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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 ...
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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.