Go to the documentation of this file.
70 #ifndef cellCoBlended_H
71 #define cellCoBlended_H
154 if (
Co1_ < 0 || Co2_ < 0 || Co1_ >=
Co2_)
157 <<
"coefficients = " <<
Co1_ <<
" and " <<
Co2_
158 <<
" should be > 0 and Co2 > Co1"
185 if (
Co1_ < 0 || Co2_ < 0 || Co1_ >=
Co2_)
188 <<
"coefficients = " <<
Co1_ <<
" and " <<
Co2_
189 <<
" should be > 0 and Co2 > Co1"
211 mesh.objectRegistry::template lookupObject<volScalarField>
219 <<
"dimensions of faceFlux are not correct"
233 zeroGradientFvPatchScalarField::typeName
249 vf.name() +
"BlendingFactor",
276 + (scalar(1.0) - bf)*
tScheme2_().weights(vf);
292 + (scalar(1.0) - bf)*
tScheme2_().interpolate(vf);
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
A class for handling words, derived from string.
A class for managing temporary objects.
const dimensionSet dimVelocity
void operator=(const cellCoBlended &)
Disallow default bitwise assignment.
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.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< surfaceInterpolationScheme< Type > > tScheme2_
Scheme 2.
tmp< surfaceInterpolationScheme< Type > > tScheme1_
Scheme 1.
const dimensionSet dimArea(sqr(dimLength))
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Pre-declare SubField and related Field type.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
InternalField & internalField()
Return internal field.
scalar deltaTValue() const
Return time step value.
const surfaceScalarField & faceFlux_
The face-flux used to compute the face Courant number.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
Generic dimensioned Type class.
const scalar Co1_
Courant number below which scheme1 is used.
cellCoBlended(const cellCoBlended &)
Disallow default bitwise copy construct.
Mesh data needed to do the Finite Volume discretisation.
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)
void correctBoundaryConditions()
Correct boundary field.
conserve internalField()+
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const scalar Co2_
Courant number above which scheme2 is used.
#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.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Two-scheme cell-based Courant number based blending differencing scheme.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
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.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Generic GeometricField class.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Base class for blended schemes to provide access to the blending factor surface field.