Go to the documentation of this file.
81 public surfaceInterpolationScheme<Type>,
82 public blendedSchemeBase<Type>
144 if (
Co1_ < 0 || Co2_ < 0 || Co1_ >=
Co2_)
147 <<
"coefficients = " <<
Co1_ <<
" and " <<
Co2_
148 <<
" should be > 0 and Co2 > Co1"
175 if (
Co1_ < 0 || Co2_ < 0 || Co1_ >=
Co2_)
178 <<
"coefficients = " <<
Co1_ <<
" and " <<
Co2_
179 <<
" should be > 0 and Co2 > Co1"
201 mesh.objectRegistry::template lookupObject<volScalarField>
209 <<
"dimensions of faceFlux are not correct"
213 return tmp<surfaceScalarField>
217 vf.name() +
"BlendingFactor",
238 tmp<surfaceScalarField>
241 const GeometricField<Type, fvPatchField, volMesh>& vf
248 + (scalar(1.0) - bf)*
tScheme2_().weights(vf);
264 + (scalar(1.0) - bf)*
tScheme2_().interpolate(vf);
316 return tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
A class for handling words, derived from string.
tmp< surfaceInterpolationScheme< Type > > tScheme1_
Scheme 1.
dimensionedScalar deltaT() const
Return time step.
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
A class for managing temporary objects.
const dimensionSet dimVelocity
const dimensionSet dimDensity
const scalar Co2_
Courant number above which scheme2 is used.
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.
const surfaceScalarField & faceFlux_
The face-flux used to compute the face Courant number.
CoBlended(const CoBlended &)
Disallow default bitwise copy construct.
void operator=(const CoBlended &)
Disallow default bitwise assignment.
dimensioned< scalar > mag(const dimensioned< Type > &)
const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
TypeName("CoBlended")
Runtime type information.
const dimensionSet dimArea(sqr(dimLength))
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
Mesh data needed to do the Finite Volume discretisation.
const scalar Co1_
Courant number below which scheme1 is used.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< surfaceInterpolationScheme< Type > > tScheme2_
Scheme 2.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Two-scheme Courant number based blending differencing scheme.
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.