Go to the documentation of this file.
33 #ifndef PureUpwindFitScheme_H
34 #define PureUpwindFitScheme_H
49 template<
class Type,
class Polynomial,
class Stencil>
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)
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
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.
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 UpwindFitData< Polynomial > & New(const fvMesh &mesh)
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)....
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.
const List< scalarList > & owncoeffs() const
Return reference to owner fit coefficients.
Mesh data needed to do the Finite Volume discretisation.
PureUpwindFitScheme(const PureUpwindFitScheme &)
Disallow default bitwise copy construct.
Upwind biased fit surface interpolation scheme that applies an explicit correction to upwind.
const surfaceScalarField & faceFlux_
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
TypeName("PureUpwindFitScheme")
Runtime type information.
const scalar linearLimitFactor_
Factor the fit is allowed to deviate from linear.
const scalar centralWeight_
Weights for central stencil.
const fvMesh & mesh() const
Return mesh reference.
Generic GeometricField class.
void operator=(const PureUpwindFitScheme &)
Disallow default bitwise assignment.