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