Go to the documentation of this file.
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"
238 extrapolatedCalculatedFvPatchScalarField::typeName
246 Co.primitiveFieldRef() =
248 Co.correctBoundaryConditions();
254 vf.name() +
"BlendingFactor",
280 bf*tScheme1_().weights(vf)
281 + (scalar(1) - bf)*tScheme2_().weights(vf);
296 bf*tScheme1_().interpolate(vf)
297 + (scalar(1) - bf)*tScheme2_().interpolate(vf);
304 return tScheme1_().corrected() || tScheme2_().corrected();
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimVelocity
const dimensionSet dimDensity
TypeName("cellCoBlended")
Runtime type information.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
static word timeName(const scalar t, const int precision=precision_)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
scalar deltaTValue() const noexcept
Return time step value.
const dimensionSet dimArea(sqr(dimLength))
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Type & lookupObject(const word &name, const bool recursive=false) const
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
Mesh data needed to do the Finite Volume discretisation.
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< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Two-scheme cell-based Courant number based blending differencing scheme.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
Abstract base class for surface interpolation schemes.
const Time & time() const
Return the top-level database.
const fvMesh & mesh() const
Return mesh reference.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Base class for blended schemes to provide access to the blending factor surface field.
const dimensionSet dimless
Dimensionless.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.