46#ifndef cellMDLimitedGrad_H
47#define cellMDLimitedGrad_H
101 k_(readScalar(schemeData))
103 if (k_ < 0 || k_ > 1)
106 <<
"coefficient = " << k_
107 <<
" should be >= 0 and <= 1"
118 const Type& maxDelta,
119 const Type& minDelta,
143 const scalar& maxDelta,
144 const scalar& minDelta,
148 const scalar extrapolate = dcf &
g;
150 if (extrapolate > maxDelta)
152 g =
g + dcf*(maxDelta - extrapolate)/
magSqr(dcf);
154 else if (extrapolate < minDelta)
156 g =
g + dcf*(minDelta - extrapolate)/
magSqr(dcf);
165 const Type& maxDelta,
166 const Type& minDelta,
170 for (
direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
172 vector gi(
g[cmpt],
g[cmpt+3],
g[cmpt+6]);
177 minDelta.component(cmpt),
const uniformDimensionedVectorField & g
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
const Cmpt & component(const direction) const
const Cmpt & z() const
Access to the vector z component.
const Cmpt & y() const
Access to the vector y component.
const Cmpt & x() const
Access to the vector x component.
Mesh data needed to do the Finite Volume discretisation.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
cellMDLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
static void limitFace(typename outerProduct< vector, Type >::type &g, const Type &maxDelta, const Type &minDelta, const vector &dcf)
cellMDLimitedGrad(const fvMesh &mesh, Istream &schemeData)
Construct from mesh and schemeData.
TypeName("cellMDLimited")
RunTime type information.
Abstract base class for gradient schemes.
const fvMesh & mesh() const
Return const reference to mesh.
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.
Mesh data needed to do the Finite Volume discretisation.
type
Volume classification types.
A class for handling words, derived from Foam::string.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.