40 #ifndef cellLimitedGrad_H
41 #define cellLimitedGrad_H
58 template<
class Type,
class Limiter>
103 if (k_ < 0 || k_ > 1)
108 ) <<
"coefficient = " << k_
109 <<
" should be >= 0 and <= 1"
123 const scalar maxDelta,
124 const scalar minDelta,
125 const scalar extrapolate
131 const Type& maxDelta,
132 const Type& minDelta,
133 const Type& extrapolate
155 template<
class Type,
class Limiter>
159 const scalar maxDelta,
160 const scalar minDelta,
161 const scalar extrapolate
166 if (extrapolate > small)
168 r = maxDelta/extrapolate;
170 else if (extrapolate < -small)
172 r = minDelta/extrapolate;
183 template<
class Type,
class Limiter>
187 const Type& maxDelta,
188 const Type& minDelta,
189 const Type& extrapolate
192 for (
direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Mesh data needed to do the Finite Volume discretisation.
cellLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
virtual tmp< VolField< typename outerProduct< vector, Type >::type > > calcGrad(const VolField< Type > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad.
void limitFace(Type &limiter, const Type &maxDelta, const Type &minDelta, const Type &extrapolate) const
void limitFaceCmpt(scalar &limiter, const scalar maxDelta, const scalar minDelta, const scalar extrapolate) const
TypeName("cellLimited")
RunTime type information.
void operator=(const cellLimitedGrad &)=delete
Disallow default bitwise assignment.
cellLimitedGrad(const fvMesh &mesh, Istream &schemeData)
Construct from mesh and schemeData.
Abstract base class for gradient schemes.
const fvMesh & mesh() const
Return mesh reference.
static tmp< gradScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new gradScheme created on freestore.
A class for managing temporary objects.
A class for handling words, derived from string.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
void limiter(surfaceScalarField &lambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
errorManipArg< error, int > exit(error &err, const int errNo=1)
word name(const bool)
Return a word representation of a bool.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
label & setComponent(label &l, const direction)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.