38 #ifndef UpwindFitScheme_H 39 #define UpwindFitScheme_H 53 template<
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\ 174 typedef UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL> \ 175 UpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \ 176 defineTemplateTypeNameAndDebugWithName \ 177 (UpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \ 179 surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \ 180 <UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \ 181 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \ 183 surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \ 184 <UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \ 185 add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_; 187 #define makeUpwindFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \ 189 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \ 190 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \ 191 makeUpwindFitSurfaceInterpolationTypeScheme \ 198 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor) \ 199 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Central-differencing interpolation scheme class.
static const Type & New(const Mesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
const fvMesh & mesh() const
Return mesh reference.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
const List< scalarList > & neicoeffs() const
Return reference to neighbour fit coefficients.
TypeName("UpwindFitScheme")
Runtime type information.
A class for handling words, derived from Foam::string.
Upwind biased fit surface interpolation scheme that applies an explicit correction to linear...
const List< scalarList > & owncoeffs() const
Return reference to owner fit coefficients.
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.
Data for the quadratic fit correction interpolation scheme to be used with upwind biased stencil...
Creates upwind stencil by shifting a centred stencil to upwind and downwind faces and optionally remo...
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Mesh data needed to do the Finite Volume discretisation.
A class for managing temporary objects.