Go to the documentation of this file.
34 template<
class Type,
class Limiter>
45 template<
class Type,
class Limiter>
64 template<
class Type,
class Limiter>
76 const GeometricField<Type, fvPatchField, volMesh>& vsf,
80 const fvMesh&
mesh = vsf.mesh();
86 > tGrad = basicGradScheme_().calcGrad(vsf,
name);
106 Field<Type> maxVsf(vsf.primitiveField());
107 Field<Type> minVsf(vsf.primitiveField());
111 const label own = owner[facei];
112 const label nei = neighbour[facei];
114 const Type& vsfOwn = vsf[own];
115 const Type& vsfNei = vsf[nei];
117 maxVsf[own] =
max(maxVsf[own], vsfNei);
118 minVsf[own] =
min(minVsf[own], vsfNei);
120 maxVsf[nei] =
max(maxVsf[nei], vsfOwn);
121 minVsf[nei] =
min(minVsf[nei], vsfOwn);
125 const typename GeometricField<Type, fvPatchField, volMesh>::Boundary& bsf =
130 const fvPatchField<Type>& psf = bsf[patchi];
135 const Field<Type> psfNei(psf.patchNeighbourField());
139 const label own = pOwner[pFacei];
140 const Type& vsfNei = psfNei[pFacei];
142 maxVsf[own] =
max(maxVsf[own], vsfNei);
143 minVsf[own] =
min(minVsf[own], vsfNei);
150 const label own = pOwner[pFacei];
151 const Type& vsfNei = psf[pFacei];
153 maxVsf[own] =
max(maxVsf[own], vsfNei);
154 minVsf[own] =
min(minVsf[own], vsfNei);
164 const Field<Type> maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
172 Field<Type>
limiter(vsf.primitiveField().size(), pTraits<Type>::one);
176 const label own = owner[facei];
177 const label nei = neighbour[facei];
185 (Cf[facei] -
C[own]) &
g[own]
194 (Cf[facei] -
C[nei]) &
g[nei]
201 const vectorField& pCf = Cf.boundaryField()[patchi];
205 const label own = pOwner[pFacei];
212 ((pCf[pFacei] -
C[own]) &
g[own])
219 Info<<
"gradient limiter for: " << vsf.name()
226 g.correctBoundaryConditions();
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
A class for managing temporary objects.
Type gAverage(const FieldField< Field, Type > &f)
Mesh data needed to do the Finite Volume discretisation.
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Ostream & endl(Ostream &os)
label min(const labelHashSet &set, label minValue=labelMax)
Field< vector > vectorField
Specialisation of Field<T> for vector.
cellMask correctBoundaryConditions()
label max(const labelHashSet &set, label maxValue=labelMin)
const uniformDimensionedVectorField & g
GeometricField< vector, fvPatchField, volMesh > volVectorField
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
word name(const expressions::valueTypeCode typeCode)
Type gMin(const FieldField< Field, Type > &f)
UList< label > labelUList
A UList of labels.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Generic GeometricField class.
tmp< areaScalarField > limiter(const areaScalarField &phi)
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Type gMax(const FieldField< Field, Type > &f)
cellLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.