31 template<
class Type,
class Limiter>
34 const Field<scalar>& limiter,
42 template<
class Type,
class Limiter>
45 const Field<vector>& limiter,
61 template<
class Type,
class Limiter>
66 typename Foam::outerProduct<Foam::vector, Type>::type,
73 const GeometricField<Type, fvPatchField, volMesh>& vsf,
77 const fvMesh& mesh = vsf.mesh();
82 <
typename outerProduct<vector, Type>::type, fvPatchField, volMesh>
83 > tGrad = basicGradScheme_().calcGrad(vsf, name);
92 typename outerProduct<vector, Type>::type,
98 const labelUList& neighbour = mesh.neighbour();
103 Field<Type> maxVsf(vsf.primitiveField());
104 Field<Type> minVsf(vsf.primitiveField());
108 label own = owner[facei];
109 label nei = neighbour[facei];
111 const Type& vsfOwn = vsf[own];
112 const Type& vsfNei = vsf[nei];
114 maxVsf[own] =
max(maxVsf[own], vsfNei);
115 minVsf[own] =
min(minVsf[own], vsfNei);
117 maxVsf[nei] =
max(maxVsf[nei], vsfOwn);
118 minVsf[nei] =
min(minVsf[nei], vsfOwn);
122 const typename GeometricField<Type, fvPatchField, volMesh>::Boundary& bsf =
127 const fvPatchField<Type>& psf = bsf[
patchi];
132 const Field<Type> psfNei(psf.patchNeighbourField());
136 label own = pOwner[pFacei];
137 const Type& vsfNei = psfNei[pFacei];
139 maxVsf[own] =
max(maxVsf[own], vsfNei);
140 minVsf[own] =
min(minVsf[own], vsfNei);
147 label own = pOwner[pFacei];
148 const Type& vsfNei = psf[pFacei];
150 maxVsf[own] =
max(maxVsf[own], vsfNei);
151 minVsf[own] =
min(minVsf[own], vsfNei);
161 const Field<Type> maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
169 Field<Type>
limiter(vsf.primitiveField().size(), pTraits<Type>::one);
173 label own = owner[facei];
174 label nei = neighbour[facei];
182 (Cf[facei] - C[own]) & g[own]
191 (Cf[facei] - C[nei]) & g[nei]
202 label own = pOwner[pFacei];
209 ((pCf[pFacei] - C[own]) & g[own])
216 Info<<
"gradient limiter for: " << vsf.name()
217 <<
" max = " <<
gMax(limiter)
218 <<
" min = " <<
gMin(limiter)
222 limitGradient(limiter, g);
223 g.correctBoundaryConditions();
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
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.
rDeltaT correctBoundaryConditions()
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Type gMin(const FieldField< Field, Type > &f)
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
volVectorField vectorField(fieldObject, mesh)
Generic GeometricField class.
GeometricField< vector, fvPatchField, volMesh > volVectorField
UList< label > labelUList
Mesh data needed to do the Finite Volume discretisation.
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 PsiMaxType &psiMax, const PsiMinType &psiMin)
cellLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Type gAverage(const FieldField< Field, Type > &f)
A class for managing temporary objects.
Tensor< scalar > tensor
Tensor of scalars.