Go to the documentation of this file.
41 scalar maxRatio = 1 +
coeff;
51 const label own = owner[facei];
52 const label nbr = neighbour[facei];
55 if (field[own] > maxRatio*field[nbr])
57 changedFaces.
append(facei);
60 else if (field[nbr] > maxRatio*field[own])
62 changedFaces.
append(facei);
79 changedFaces.
append(facei);
93 cellData[celli] = field[celli];
110 mesh.globalData().nTotalCells(),
116 field[celli] = cellData[celli].value();
128 const scalar alphaDiff,
130 const scalar alphaMin
143 cellData[celli] = field[celli];
154 const label own = owner[facei];
155 const label nbr = neighbour[facei];
159 changedFaces.
append(facei);
181 alpha.boundaryField()[
patchi].patchNeighbourField()
184 if (
mag(
alpha[own] - alphapn[patchFacei]) > alphaDiff)
186 changedFaces.
append(facei);
194 changedFacesInfo.
shrink();
208 smoothData.setFaceInfo(changedFaces, changedFacesInfo);
214 field[celli] = cellData[celli].value();
226 const scalar alphaDiff
246 const label own = owner[facei];
247 const label nbr = neighbour[facei];
251 changedFaces.
append(facei);
273 alpha.boundaryField()[
patchi].patchNeighbourField()
276 if (
mag(
alpha[own] - alphapn[patchFacei]) > alphaDiff)
278 changedFaces.
append(facei);
289 changedFacesInfo.
shrink();
299 sweepData.setFaceInfo(changedFaces, changedFacesInfo);
305 if (cellData[celli].valid(
sweepData.data()))
307 field[celli] =
max(field[celli], cellData[celli].value());
Class used to pass additional data in.
Helper class used by fvc::sweep function.
#define forAll(list, i)
Loop across all elements in list.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
dimensioned< scalar > mag(const dimensioned< Type > &)
Helper class used by the fvc::smooth and fvc::spread functions.
void spread(volScalarField &field, const volScalarField &alpha, const label nLayers, const scalar alphaDiff=0.2, const scalar alphaMax=0.99, const scalar alphaMin=0.01)
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.
A patch is a list of labels that address the faces in the global face list.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Mesh data needed to do the Finite Volume discretisation.
scalar maxRatio
Cut off distance.
label start() const
Return start label of this patch in the polyMesh face list.
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
void correctBoundaryConditions()
Correct boundary field.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Wave propagation of information through grid. Every iteration information goes through one layer of c...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Generic GeometricField class.
void sweep(volScalarField &field, const volScalarField &alpha, const label nLayers, const scalar alphaDiff=0.2)
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
void smooth(volScalarField &field, const scalar coeff)