Go to the documentation of this file.
48 const fvMesh&
mesh = vsf.mesh();
50 tmp<volVectorField> tGrad = basicGradScheme_().calcGrad(vsf,
name);
70 label own = owner[facei];
71 label nei = neighbour[facei];
73 scalar vsfOwn = vsf[own];
74 scalar vsfNei = vsf[nei];
76 maxVsf[own] =
max(maxVsf[own], vsfNei);
77 minVsf[own] =
min(minVsf[own], vsfNei);
79 maxVsf[nei] =
max(maxVsf[nei], vsfOwn);
80 minVsf[nei] =
min(minVsf[nei], vsfOwn);
84 const volScalarField::GeometricBoundaryField& bsf = vsf.boundaryField();
94 const scalarField psfNei(psf.patchNeighbourField());
98 label own = pOwner[pFacei];
99 scalar vsfNei = psfNei[pFacei];
101 maxVsf[own] =
max(maxVsf[own], vsfNei);
102 minVsf[own] =
min(minVsf[own], vsfNei);
109 label own = pOwner[pFacei];
110 scalar vsfNei = psf[pFacei];
112 maxVsf[own] =
max(maxVsf[own], vsfNei);
113 minVsf[own] =
min(minVsf[own], vsfNei);
123 const scalarField maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
137 label own = owner[facei];
138 label nei = neighbour[facei];
146 (Cf[facei] -
C[own]) &
g[own]
155 (Cf[facei] -
C[nei]) &
g[nei]
166 label own = pOwner[pFacei];
173 (pCf[pFacei] -
C[own]) &
g[own]
180 Info<<
"gradient limiter for: " << vsf.name()
187 g.correctBoundaryConditions();
224 label own = owner[facei];
225 label nei = neighbour[facei];
227 const vector& vsfOwn = vsf[own];
228 const vector& vsfNei = vsf[nei];
230 maxVsf[own] =
max(maxVsf[own], vsfNei);
231 minVsf[own] =
min(minVsf[own], vsfNei);
233 maxVsf[nei] =
max(maxVsf[nei], vsfOwn);
234 minVsf[nei] =
min(minVsf[nei], vsfOwn);
251 label own = pOwner[pFacei];
252 const vector& vsfNei = psfNei[pFacei];
254 maxVsf[own] =
max(maxVsf[own], vsfNei);
255 minVsf[own] =
min(minVsf[own], vsfNei);
262 label own = pOwner[pFacei];
263 const vector& vsfNei = psf[pFacei];
265 maxVsf[own] =
max(maxVsf[own], vsfNei);
266 minVsf[own] =
min(minVsf[own], vsfNei);
276 const vectorField maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
290 label own = owner[facei];
291 label nei = neighbour[facei];
299 (Cf[facei] -
C[own]) &
g[own]
308 (Cf[facei] -
C[nei]) &
g[nei]
319 label own = pOwner[pFacei];
326 ((pCf[pFacei] -
C[own]) &
g[own])
333 Info<<
"gradient limiter for: " << vsf.name()
351 g.correctBoundaryConditions();
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from string.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
Type gAverage(const FieldField< Field, Type > &f)
Tensor< scalar > tensor
Tensor of scalars.
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 scalar psiMax, const scalar psiMin)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
makeFvGradScheme(cellLimitedGrad)
const dimensionedVector & g
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual bool coupled() const
Return true if this patch field is coupled.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
InternalField & internalField()
Return internal field.
volVectorField vectorField(fieldObject, mesh)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Mesh data needed to do the Finite Volume discretisation.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
U correctBoundaryConditions()
volScalarField scalarField(fieldObject, mesh)
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
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.
Type gMin(const FieldField< Field, Type > &f)
UList< label > labelUList
Graphite solid properties.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Generic GeometricField class.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)
word name(const complex &)
Return a string representation of a complex.