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));
134 label own = owner[facei];
135 label nei = neighbour[facei];
164 label own = pOwner[pFacei];
176 g.correctBoundaryConditions();
213 label own = owner[facei];
214 label nei = neighbour[facei];
216 const vector& vsfOwn = vsf[own];
217 const vector& vsfNei = vsf[nei];
219 maxVsf[own] =
max(maxVsf[own], vsfNei);
220 minVsf[own] =
min(minVsf[own], vsfNei);
222 maxVsf[nei] =
max(maxVsf[nei], vsfOwn);
223 minVsf[nei] =
min(minVsf[nei], vsfOwn);
240 label own = pOwner[pFacei];
241 const vector& vsfNei = psfNei[pFacei];
243 maxVsf[own] =
max(maxVsf[own], vsfNei);
244 minVsf[own] =
min(minVsf[own], vsfNei);
251 label own = pOwner[pFacei];
252 const vector& vsfNei = psf[pFacei];
254 maxVsf[own] =
max(maxVsf[own], vsfNei);
255 minVsf[own] =
min(minVsf[own], vsfNei);
265 const vectorField maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
276 label own = owner[facei];
277 label nei = neighbour[facei];
306 label own = pOwner[pFacei];
318 g.correctBoundaryConditions();
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
makeFvGradScheme(cellMDLimitedGrad)
A class for handling words, derived from string.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
const dimensionedVector & g
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()
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.
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...
UList< label > labelUList
Graphite solid properties.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Generic GeometricField class.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.