Go to the documentation of this file.
33 #ifndef UpwindFitScheme_H
34 #define UpwindFitScheme_H
48 template<
class Type,
class Polynomial,
class Stencil>
161 #define makeUpwindFitSurfaceInterpolationTypeScheme\
169 typedef UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL> \
170 UpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \
171 defineTemplateTypeNameAndDebugWithName \
172 (UpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
174 surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
175 <UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL> > \
176 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
178 surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
179 <UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL> > \
180 add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_;
182 #define makeUpwindFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \
184 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
185 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
186 makeUpwindFitSurfaceInterpolationTypeScheme \
193 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor) \
194 makeUpwindFitSurfaceInterpolationTypeScheme(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.
const scalar linearLimitFactor_
Factor the fit is allowed to deviate from linear.
A class for managing temporary objects.
const scalar centralWeight_
Weights for central stencil.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
static const UpwindFitData< Polynomial > & New(const fvMesh &mesh)
UpwindFitScheme(const UpwindFitScheme &)
Disallow default bitwise copy construct.
void operator=(const UpwindFitScheme &)
Disallow default bitwise assignment.
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...
TypeName("UpwindFitScheme")
Runtime type information.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
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.
const surfaceScalarField & faceFlux_
Reference to the surface flux used to choose upwind direction.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
const fvMesh & mesh() const
Return mesh reference.
Central-differencing interpolation scheme class.
Generic GeometricField class.
Upwind biased fit surface interpolation scheme that applies an explicit correction to linear.