75 #ifndef cellCoBlended_H 76 #define cellCoBlended_H 80 #include "surfaceInterpolate.H" 144 Co1_(readScalar(is)),
149 Co2_(readScalar(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)),
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" 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>>
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
scalar deltaTValue() const noexcept
Return time step value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
const word & name() const noexcept
Return the object name.
TypeName("cellCoBlended")
Runtime type information.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
const dimensionSet dimless
Dimensionless.
const fvMesh & mesh() const
Return mesh reference.
const Time & time() const
Return the top-level database.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
const dimensionSet dimVolume(pow3(dimLength))
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
A class for handling words, derived from Foam::string.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Base class for blended schemes to provide access to the blending factor surface field.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
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.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
A class for managing temporary objects.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Abstract base class for surface interpolation schemes.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
Two-scheme cell-based Courant number based blending differencing scheme.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
const dimensionSet & dimensions() const noexcept
Return dimensions.
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)