38#ifndef UpwindFitScheme_H
39#define UpwindFitScheme_H
53template<
class Type,
class Polynomial,
class Stencil>
66 const scalar linearLimitFactor_;
69 const scalar centralWeight_;
96 linearLimitFactor_(readScalar(is)),
111 linearLimitFactor_(readScalar(is)),
166#define makeUpwindFitSurfaceInterpolationTypeScheme\
174typedef UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL> \
175 UpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \
176defineTemplateTypeNameAndDebugWithName \
177 (UpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
179surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
180<UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
181 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
183surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
184<UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
185 add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_;
187#define makeUpwindFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \
189makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
190makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
191makeUpwindFitSurfaceInterpolationTypeScheme \
198makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor) \
199makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Data for the quadratic fit correction interpolation scheme to be used with upwind biased stencil.
const List< scalarList > & owncoeffs() const
Return reference to owner fit coefficients.
const List< scalarList > & neicoeffs() const
Return reference to neighbour fit coefficients.
Upwind biased fit surface interpolation scheme that applies an explicit correction to linear.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
UpwindFitScheme(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
TypeName("UpwindFitScheme")
Runtime type information.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Creates upwind stencil by shifting a centred stencil to upwind and downwind faces and optionally remo...
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.
Mesh data needed to do the Finite Volume discretisation.
Central-differencing interpolation scheme class.
const fvMesh & mesh() const
Return mesh reference.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.