fvcSmooth.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvcSmooth.H"
27 #include "volFields.H"
28 #include "FaceCellWave.H"
29 #include "smoothData.H"
30 #include "sweepData.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
35 (
36  volScalarField& field,
37  const scalar coeff
38 )
39 {
40  const fvMesh& mesh = field.mesh();
41  scalar maxRatio = 1 + coeff;
42 
43  DynamicList<label> changedFaces(mesh.nFaces()/100 + 100);
44  DynamicList<smoothData> changedFacesInfo(changedFaces.size());
45 
46  const labelUList& owner = mesh.owner();
47  const labelUList& neighbour = mesh.neighbour();
48 
49  forAll(owner, facei)
50  {
51  const label own = owner[facei];
52  const label nbr = neighbour[facei];
53 
54  // Check if owner value much larger than neighbour value or vice versa
55  if (field[own] > maxRatio*field[nbr])
56  {
57  changedFaces.append(facei);
58  changedFacesInfo.append(smoothData(field[own]));
59  }
60  else if (field[nbr] > maxRatio*field[own])
61  {
62  changedFaces.append(facei);
63  changedFacesInfo.append(smoothData(field[nbr]));
64  }
65  }
66 
67  // Insert all faces of coupled patches - FaceCellWave will correct them
68  forAll(mesh.boundaryMesh(), patchi)
69  {
70  const polyPatch& patch = mesh.boundaryMesh()[patchi];
71 
72  if (patch.coupled())
73  {
74  forAll(patch, patchFacei)
75  {
76  label facei = patch.start() + patchFacei;
77  label own = mesh.faceOwner()[facei];
78 
79  changedFaces.append(facei);
80  changedFacesInfo.append(smoothData(field[own]));
81  }
82  }
83  }
84 
85  changedFaces.shrink();
86  changedFacesInfo.shrink();
87 
88  // Set initial field on cells
89  List<smoothData> cellData(mesh.nCells());
90 
91  forAll(field, celli)
92  {
93  cellData[celli] = field[celli];
94  }
95 
96  // Set initial field on faces
97  List<smoothData> faceData(mesh.nFaces());
98 
100  td.maxRatio = maxRatio;
101 
102  // Propagate information over whole domain
104  (
105  mesh,
106  changedFaces,
107  changedFacesInfo,
108  faceData,
109  cellData,
110  mesh.globalData().nTotalCells(), // max iterations
111  td
112  );
113 
114  forAll(field, celli)
115  {
116  field[celli] = cellData[celli].value();
117  }
118 
120 }
121 
122 
124 (
125  volScalarField& field,
126  const volScalarField& alpha,
127  const label nLayers,
128  const scalar alphaDiff,
129  const scalar alphaMax,
130  const scalar alphaMin
131 )
132 {
133  const fvMesh& mesh = field.mesh();
134 
135  DynamicList<label> changedFaces(mesh.nFaces()/100 + 100);
136  DynamicList<smoothData> changedFacesInfo(changedFaces.size());
137 
138  // Set initial field on cells
139  List<smoothData> cellData(mesh.nCells());
140 
141  forAll(field, celli)
142  {
143  cellData[celli] = field[celli];
144  }
145 
146  // Set initial field on faces
147  List<smoothData> faceData(mesh.nFaces());
148 
149  const labelUList& owner = mesh.owner();
150  const labelUList& neighbour = mesh.neighbour();
151 
152  forAll(owner, facei)
153  {
154  const label own = owner[facei];
155  const label nbr = neighbour[facei];
156 
157  if (mag(alpha[own] - alpha[nbr]) > alphaDiff)
158  {
159  changedFaces.append(facei);
160  changedFacesInfo.append
161  (
162  smoothData(max(field[own], field[nbr]))
163  );
164  }
165  }
166 
167  // Insert all faces of coupled patches - FaceCellWave will correct them
168  forAll(mesh.boundaryMesh(), patchi)
169  {
170  const polyPatch& patch = mesh.boundaryMesh()[patchi];
171 
172  if (patch.coupled())
173  {
174  forAll(patch, patchFacei)
175  {
176  label facei = patch.start() + patchFacei;
177  label own = mesh.faceOwner()[facei];
178 
179  const scalarField alphapn
180  (
181  alpha.boundaryField()[patchi].patchNeighbourField()
182  );
183 
184  if (mag(alpha[own] - alphapn[patchFacei]) > alphaDiff)
185  {
186  changedFaces.append(facei);
187  changedFacesInfo.append(smoothData(field[own]));
188  }
189  }
190  }
191  }
192 
193  changedFaces.shrink();
194  changedFacesInfo.shrink();
195 
197  td.maxRatio = 1.0;
198 
199  // Propagate information over whole domain
201  (
202  mesh,
203  faceData,
204  cellData,
205  td
206  );
207 
208  smoothData.setFaceInfo(changedFaces, changedFacesInfo);
209 
210  smoothData.iterate(nLayers);
211 
212  forAll(field, celli)
213  {
214  field[celli] = cellData[celli].value();
215  }
216 
218 }
219 
220 
221 void Foam::fvc::sweep
222 (
223  volScalarField& field,
224  const volScalarField& alpha,
225  const label nLayers,
226  const scalar alphaDiff
227 )
228 {
229  const fvMesh& mesh = field.mesh();
230 
231  DynamicList<label> changedFaces(mesh.nFaces()/100 + 100);
232  DynamicList<sweepData> changedFacesInfo(changedFaces.size());
233 
234  // Set initial field on cells
235  List<sweepData> cellData(mesh.nCells());
236 
237  // Set initial field on faces
238  List<sweepData> faceData(mesh.nFaces());
239 
240  const labelUList& owner = mesh.owner();
241  const labelUList& neighbour = mesh.neighbour();
242  const vectorField& Cf = mesh.faceCentres();
243 
244  forAll(owner, facei)
245  {
246  const label own = owner[facei];
247  const label nbr = neighbour[facei];
248 
249  if (mag(alpha[own] - alpha[nbr]) > alphaDiff)
250  {
251  changedFaces.append(facei);
252  changedFacesInfo.append
253  (
254  sweepData(max(field[own], field[nbr]), Cf[facei])
255  );
256  }
257  }
258 
259  // Insert all faces of coupled patches - FaceCellWave will correct them
260  forAll(mesh.boundaryMesh(), patchi)
261  {
262  const polyPatch& patch = mesh.boundaryMesh()[patchi];
263 
264  if (patch.coupled())
265  {
266  forAll(patch, patchFacei)
267  {
268  label facei = patch.start() + patchFacei;
269  label own = mesh.faceOwner()[facei];
270 
271  const scalarField alphapn
272  (
273  alpha.boundaryField()[patchi].patchNeighbourField()
274  );
275 
276  if (mag(alpha[own] - alphapn[patchFacei]) > alphaDiff)
277  {
278  changedFaces.append(facei);
279  changedFacesInfo.append
280  (
281  sweepData(field[own], Cf[facei])
282  );
283  }
284  }
285  }
286  }
287 
288  changedFaces.shrink();
289  changedFacesInfo.shrink();
290 
291  // Propagate information over whole domain
293  (
294  mesh,
295  faceData,
296  cellData
297  );
298 
299  sweepData.setFaceInfo(changedFaces, changedFacesInfo);
300 
301  sweepData.iterate(nLayers);
302 
303  forAll(field, celli)
304  {
305  if (cellData[celli].valid(sweepData.data()))
306  {
307  field[celli] = max(field[celli], cellData[celli].value());
308  }
309  }
310 
312 }
313 
314 
315 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual label iterate(const label maxIter)
Iterate until no changes or maxIter reached. Returns actual.
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)
Definition: fvcSmooth.C:124
Class used to pass additional data in.
Definition: smoothData.H:53
void setFaceInfo(const labelList &changedFaces, const List< Type > &changedFacesInfo)
Set initial changed faces.
Definition: FaceCellWave.C:317
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
label nFaces() const
label nCells() const
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:76
label nTotalCells() const
Return total number of cells in decomposed mesh.
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:35
Helper class used by fvc::sweep function.
Definition: sweepData.H:47
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:284
Helper class used by the fvc::smooth and fvc::spread functions.
Definition: smoothData.H:47
void sweep(volScalarField &field, const volScalarField &alpha, const label nLayers, const scalar alphaDiff=0.2)
Definition: fvcSmooth.C:222
dynamicFvMesh & mesh
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:319
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
const Mesh & mesh() const
Return mesh.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:278
const vectorField & faceCentres() const
label patchi
const TrackingData & data() const
Additional data to be passed into container.
Definition: FaceCellWave.H:340
scalar maxRatio
Cut off distance.
Definition: smoothData.H:58
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:309
void correctBoundaryConditions()
Correct boundary field.
dimensioned< scalar > mag(const dimensioned< Type > &)
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...