Go to the documentation of this file.
38 #ifndef PureUpwindFitScheme_H
39 #define PureUpwindFitScheme_H
54 template<
class Type,
class Polynomial,
class Stencil>
64 const scalar linearLimitFactor_;
67 const scalar centralWeight_;
97 linearLimitFactor_(readScalar(is)),
111 linearLimitFactor_(readScalar(is)),
162 #define makePureUpwindFitSurfaceInterpolationTypeScheme\
170 typedef PureUpwindFitScheme<TYPE, POLYNOMIAL, STENCIL> \
171 PureUpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \
172 defineTemplateTypeNameAndDebugWithName \
173 (PureUpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
175 surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
176 <PureUpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
177 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
179 surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
180 <PureUpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
181 add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_;
183 #define makePureUpwindFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \
185 makePureUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
186 makePureUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
187 makePureUpwindFitSurfaceInterpolationTypeScheme \
194 makePureUpwindFitSurfaceInterpolationTypeScheme \
201 makePureUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Upwind differencing scheme class.
static const Type & New(const Mesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
const List< scalarList > & neicoeffs() const
Return reference to neighbour fit coefficients.
Creates upwind stencil by shifting a centred stencil to upwind and downwind faces and optionally remo...
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
const Type & lookupObject(const word &name, const bool recursive=false) const
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Data for the quadratic fit correction interpolation scheme to be used with upwind biased stencil.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > weightedSum(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &fld, const List< List< scalar >> &ownWeights, const List< List< scalar >> &neiWeights) const
Sum vol field contributions to create face values.
const List< scalarList > & owncoeffs() const
Return reference to owner fit coefficients.
Mesh data needed to do the Finite Volume discretisation.
Upwind biased fit surface interpolation scheme that applies an explicit correction to upwind.
const surfaceScalarField & faceFlux_
TypeName("PureUpwindFitScheme")
Runtime type information.
const fvMesh & mesh() const
Return mesh reference.