Go to the documentation of this file.
40 #ifndef cellLimitedGrad_H
41 #define cellLimitedGrad_H
101 ) <<
"coefficient = " <<
k_
102 <<
" should be >= 0 and <= 1"
113 const Type& maxDelta,
114 const Type& minDelta,
115 const Type& extrapolate
138 const scalar& maxDelta,
139 const scalar& minDelta,
140 const scalar& extrapolate
143 if (extrapolate > maxDelta + VSMALL)
147 else if (extrapolate < minDelta - VSMALL)
158 const Type& maxDelta,
159 const Type& minDelta,
160 const Type& extrapolate
163 for (
direction cmpt=0; cmpt<Type::nComponents; cmpt++)
168 maxDelta.component(cmpt),
169 minDelta.component(cmpt),
170 extrapolate.component(cmpt)
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
cellLimitedGrad(const cellLimitedGrad &)
Disallow default bitwise copy construct.
const scalar k_
Limiter coefficient.
A class for handling words, derived from string.
A class for managing temporary objects.
Mesh data needed to do the Finite Volume discretisation.
typeOfRank< typename pTraits< arg1 >::cmptType, int(pTraits< arg1 >::rank)+int(pTraits< arg2 >::rank) >::type type
void limiter(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
static tmp< gradScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new gradScheme created on freestore.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
void operator=(const cellLimitedGrad &)
Disallow default bitwise assignment.
TypeName("cellLimited")
RunTime type information.
const fvMesh & mesh() const
Return mesh reference.
tmp< fv::gradScheme< Type > > basicGradScheme_
Abstract base class for gradient schemes.
Mesh data needed to do the Finite Volume discretisation.
errorManipArg< error, int > exit(error &err, const int errNo=1)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
cellLimitedGrad(const fvMesh &mesh, Istream &schemeData)
Construct from mesh and schemeData.
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Generic GeometricField class.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
static void limitFace(Type &limiter, const Type &maxDelta, const Type &minDelta, const Type &extrapolate)
word name(const complex &)
Return a string representation of a complex.
cellLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.