34template<
class Type,
class Limiter>
35void Foam::fv::cellLimitedGrad<Type, Limiter>::limitGradient
37 const Field<scalar>& limiter,
45template<
class Type,
class Limiter>
46void Foam::fv::cellLimitedGrad<Type, Limiter>::limitGradient
48 const Field<vector>& limiter,
64template<
class Type,
class Limiter>
78 basicGradScheme_().calcGrad(
g, vsf);
85 const labelUList& owner =
mesh.owner();
86 const labelUList& neighbour =
mesh.neighbour();
88 const volVectorField&
C =
mesh.C();
89 const surfaceVectorField& Cf =
mesh.Cf();
96 const label own = owner[facei];
97 const label nei = neighbour[facei];
99 const Type& vsfOwn = vsf[own];
100 const Type& vsfNei = vsf[nei];
102 maxVsf[own] = max(maxVsf[own], vsfNei);
103 minVsf[own] = min(minVsf[own], vsfNei);
105 maxVsf[nei] = max(maxVsf[nei], vsfOwn);
106 minVsf[nei] = min(minVsf[nei], vsfOwn);
114 const fvPatchField<Type>& psf = bsf[patchi];
119 const Field<Type> psfNei(psf.patchNeighbourField());
123 const label own = pOwner[pFacei];
124 const Type& vsfNei = psfNei[pFacei];
126 maxVsf[own] =
max(maxVsf[own], vsfNei);
127 minVsf[own] =
min(minVsf[own], vsfNei);
134 const label own = pOwner[pFacei];
135 const Type& vsfNei = psf[pFacei];
137 maxVsf[own] =
max(maxVsf[own], vsfNei);
138 minVsf[own] =
min(minVsf[own], vsfNei);
148 const Field<Type> maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
160 const label own = owner[facei];
161 const label nei = neighbour[facei];
169 (Cf[facei] -
C[own]) &
g[own]
178 (Cf[facei] -
C[nei]) &
g[nei]
185 const vectorField& pCf = Cf.boundaryField()[patchi];
189 const label own = pOwner[pFacei];
196 ((pCf[pFacei] -
C[own]) &
g[own])
206 Info<<
"gradient limiter for: " << vsf.
name()
207 <<
" min = " <<
limits.min()
208 <<
" max = " <<
limits.max()
209 <<
" average: " << avg <<
endl;
212 limitGradient(limiter,
g);
213 g.correctBoundaryConditions();
214 gaussGrad<Type>::correctBoundaryConditions(vsf,
g);
218template<
class Type,
class Limiter>
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
const uniformDimensionedVectorField & g
Graphite solid properties.
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Generic GeometricField class.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
const fileName & instance() const noexcept
Read access to instance path component.
void size(const label n)
Older name for setAddressableSize.
Mesh data needed to do the Finite Volume discretisation.
virtual bool coupled() const
True if the patch field is coupled.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
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.
void limitFace(Type &limiter, const Type &maxDelta, const Type &minDelta, const Type &extrapolate) const
static void correctBoundaryConditions(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > &)
Correct the boundary values of the gradient using the patchField snGrad functions.
const fvMesh & mesh() const
Return const reference to mesh.
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
A traits class, which is primarily used for primitives and vector-space.
Tensor of scalars, i.e. Tensor<scalar>.
A class for managing temporary objects.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Mesh data needed to do the Finite Volume discretisation.
A class for handling words, derived from Foam::string.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Field< vector > vectorField
Specialisation of Field<T> for vector.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
UList< label > labelUList
A UList of labels.
tmp< areaScalarField > limiter(const areaScalarField &phi)
#define forAll(list, i)
Loop across all elements in list.