33 #ifndef PureUpwindFitScheme_H
34 #define PureUpwindFitScheme_H
49 template<
class Type,
class Polynomial,
class Stencil>
59 const scalar linearLimitFactor_;
62 const scalar centralWeight_;
157 #define makePureUpwindFitSurfaceInterpolationTypeScheme\
165 typedef PureUpwindFitScheme<TYPE, POLYNOMIAL, STENCIL> \
166 PureUpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \
167 defineTemplateTypeNameAndDebugWithName \
168 (PureUpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
170 surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
171 <PureUpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
172 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
174 surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
175 <PureUpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
176 add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_;
178 #define makePureUpwindFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \
180 makePureUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
181 makePureUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
182 makePureUpwindFitSurfaceInterpolationTypeScheme \
189 makePureUpwindFitSurfaceInterpolationTypeScheme \
196 makePureUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
static Type & New(const word &name, const Mesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Upwind biased fit surface interpolation scheme that applies an explicit correction to upwind.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
TypeName("PureUpwindFitScheme")
Runtime type information.
void operator=(const PureUpwindFitScheme &)=delete
Disallow default bitwise assignment.
PureUpwindFitScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
virtual tmp< SurfaceField< Type > > correction(const VolField< Type > &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.
const List< scalarList > & neicoeffs() const
Return reference to neighbour fit coefficients.
const List< scalarList > & owncoeffs() const
Return reference to owner fit coefficients.
Creates upwind stencil by shifting a centred stencil to upwind and downwind faces and optionally remo...
tmp< SurfaceField< Type > > weightedSum(const surfaceScalarField &phi, const VolField< Type > &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.
const surfaceScalarField & faceFlux_
const fvMesh & mesh() const
Return mesh reference.
A class for managing temporary objects.
Upwind interpolation scheme class.
A class for handling words, derived from string.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.