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 Copyright (C) 2020 OpenCFD Ltd.
10 Copyright (C) 2020 Henning Scheufler
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "fvcSmooth.H"
31#include "volFields.H"
32#include "FaceCellWave.H"
33#include "smoothData.H"
34#include "sweepData.H"
35#include "fvMatrices.H"
36#include "fvcVolumeIntegrate.H"
37#include "fvmLaplacian.H"
38#include "fvmSup.H"
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
44(
46 const scalar coeff
47)
48{
49 const fvMesh& mesh = field.mesh();
50 scalar maxRatio = 1 + coeff;
51
52 DynamicList<label> changedFaces(mesh.nFaces()/100 + 100);
53 DynamicList<smoothData> changedFacesInfo(changedFaces.size());
54
55 const labelUList& owner = mesh.owner();
56 const labelUList& neighbour = mesh.neighbour();
57
58 forAll(owner, facei)
59 {
60 const label own = owner[facei];
61 const label nbr = neighbour[facei];
62
63 // Check if owner value much larger than neighbour value or vice versa
64 if (field[own] > maxRatio*field[nbr])
65 {
66 changedFaces.append(facei);
67 changedFacesInfo.append(smoothData(field[own]));
68 }
69 else if (field[nbr] > maxRatio*field[own])
70 {
71 changedFaces.append(facei);
72 changedFacesInfo.append(smoothData(field[nbr]));
73 }
74 }
75
76 // Insert all faces of coupled patches - FaceCellWave will correct them
77 forAll(mesh.boundaryMesh(), patchi)
78 {
79 const polyPatch& patch = mesh.boundaryMesh()[patchi];
80
81 if (patch.coupled())
82 {
83 forAll(patch, patchFacei)
84 {
85 label facei = patch.start() + patchFacei;
86 label own = mesh.faceOwner()[facei];
87
88 changedFaces.append(facei);
89 changedFacesInfo.append(smoothData(field[own]));
90 }
91 }
92 }
93
94 changedFaces.shrink();
95 changedFacesInfo.shrink();
96
97 // Set initial field on cells
98 List<smoothData> cellData(mesh.nCells());
99
100 forAll(field, celli)
101 {
102 cellData[celli] = field[celli];
103 }
104
105 // Set initial field on faces
106 List<smoothData> faceData(mesh.nFaces());
107
109 td.maxRatio = maxRatio;
110
111 // Propagate information over whole domain
113 (
114 mesh,
115 changedFaces,
116 changedFacesInfo,
117 faceData,
118 cellData,
119 mesh.globalData().nTotalCells(), // max iterations
120 td
121 );
122
123 forAll(field, celli)
124 {
125 field[celli] = cellData[celli].value();
126 }
127
128 field.correctBoundaryConditions();
129}
130
131
133(
135 const volScalarField& alpha,
136 const label nLayers,
137 const scalar alphaDiff,
138 const scalar alphaMax,
139 const scalar alphaMin
140)
141{
142 const fvMesh& mesh = field.mesh();
143
144 DynamicList<label> changedFaces(mesh.nFaces()/100 + 100);
145 DynamicList<smoothData> changedFacesInfo(changedFaces.size());
146
147 // Set initial field on cells
148 List<smoothData> cellData(mesh.nCells());
149
150 forAll(field, celli)
151 {
152 cellData[celli] = field[celli];
153 }
154
155 // Set initial field on faces
156 List<smoothData> faceData(mesh.nFaces());
157
158 const labelUList& owner = mesh.owner();
159 const labelUList& neighbour = mesh.neighbour();
160
161 forAll(owner, facei)
162 {
163 const label own = owner[facei];
164 const label nbr = neighbour[facei];
165
166 if (mag(alpha[own] - alpha[nbr]) > alphaDiff)
167 {
168 changedFaces.append(facei);
169 changedFacesInfo.append
170 (
171 smoothData(max(field[own], field[nbr]))
172 );
173 }
174 }
175
176 // Insert all faces of coupled patches - FaceCellWave will correct them
177 forAll(mesh.boundaryMesh(), patchi)
178 {
179 const polyPatch& patch = mesh.boundaryMesh()[patchi];
180
181 if (patch.coupled())
182 {
183 forAll(patch, patchFacei)
184 {
185 label facei = patch.start() + patchFacei;
186 label own = mesh.faceOwner()[facei];
187
188 const scalarField alphapn
189 (
190 alpha.boundaryField()[patchi].patchNeighbourField()
191 );
192
193 if (mag(alpha[own] - alphapn[patchFacei]) > alphaDiff)
194 {
195 changedFaces.append(facei);
196 changedFacesInfo.append(smoothData(field[own]));
197 }
198 }
199 }
200 }
201
202 changedFaces.shrink();
203 changedFacesInfo.shrink();
204
206 td.maxRatio = 1.0;
207
208 // Propagate information over whole domain
210 (
211 mesh,
212 faceData,
213 cellData,
214 td
215 );
216
217 smoothData.setFaceInfo(changedFaces, changedFacesInfo);
218
219 smoothData.iterate(nLayers);
220
221 forAll(field, celli)
222 {
223 field[celli] = cellData[celli].value();
224 }
225
226 field.correctBoundaryConditions();
227}
228
229
231(
233 const volScalarField& alpha,
234 const label nLayers,
235 const scalar alphaDiff
236)
237{
238 const fvMesh& mesh = field.mesh();
239
240 DynamicList<label> changedFaces(mesh.nFaces()/100 + 100);
241 DynamicList<sweepData> changedFacesInfo(changedFaces.size());
242
243 // Set initial field on cells
244 List<sweepData> cellData(mesh.nCells());
245
246 // Set initial field on faces
247 List<sweepData> faceData(mesh.nFaces());
248
249 const labelUList& owner = mesh.owner();
250 const labelUList& neighbour = mesh.neighbour();
251 const vectorField& Cf = mesh.faceCentres();
252
253 forAll(owner, facei)
254 {
255 const label own = owner[facei];
256 const label nbr = neighbour[facei];
257
258 if (mag(alpha[own] - alpha[nbr]) > alphaDiff)
259 {
260 changedFaces.append(facei);
261 changedFacesInfo.append
262 (
263 sweepData(max(field[own], field[nbr]), Cf[facei])
264 );
265 }
266 }
267
268 // Insert all faces of coupled patches - FaceCellWave will correct them
269 forAll(mesh.boundaryMesh(), patchi)
270 {
271 const polyPatch& patch = mesh.boundaryMesh()[patchi];
272
273 if (patch.coupled())
274 {
275 forAll(patch, patchFacei)
276 {
277 label facei = patch.start() + patchFacei;
278 label own = mesh.faceOwner()[facei];
279
280 const scalarField alphapn
281 (
282 alpha.boundaryField()[patchi].patchNeighbourField()
283 );
284
285 if (mag(alpha[own] - alphapn[patchFacei]) > alphaDiff)
286 {
287 changedFaces.append(facei);
288 changedFacesInfo.append
289 (
290 sweepData(field[own], Cf[facei])
291 );
292 }
293 }
294 }
295 }
296
297 changedFaces.shrink();
298 changedFacesInfo.shrink();
299
300 // Propagate information over whole domain
302 (
303 mesh,
304 faceData,
305 cellData
306 );
307
308 sweepData.setFaceInfo(changedFaces, changedFacesInfo);
309
310 sweepData.iterate(nLayers);
311
312 forAll(field, celli)
313 {
314 if (cellData[celli].valid(sweepData.data()))
315 {
316 field[celli] = max(field[celli], cellData[celli].value());
317 }
318 }
319
320 field.correctBoundaryConditions();
321}
322
323
325(
326 volScalarField& mDotOut,
327 const volScalarField& mDotIn,
328 const volScalarField& alpha1,
329 const volScalarField& alpha2,
330 const dimensionedScalar& D,
331 const scalar cutoff
332)
333{
334 const fvMesh& mesh = alpha1.mesh();
335
336 volScalarField mDotSmear
337 (
339 (
340 "mDotSmear",
341 mesh.time().timeName(),
342 mesh,
343 IOobject::NO_READ,
344 IOobject::NO_WRITE
345 ),
346 mesh,
347 dimensionedScalar(mDotOut.dimensions(), Zero),
349 );
350
351 //- Smearing of source term field
352 fvScalarMatrix mSourceEqn
353 (
354 fvm::Sp(scalar(1), mDotSmear)
355 - fvm::laplacian(D, mDotSmear)
356 ==
357 mDotIn
358 );
359
360 mSourceEqn.solve();
361
362 // Cut cells with cutoff < alpha1 < 1-cutoff and rescale remaining
363 // source term field
364
365 dimensionedScalar intvDotLiquid("intvDotLiquid", dimMass/dimTime, 0.0);
366 dimensionedScalar intvDotVapor ("intvDotVapor", dimMass/dimTime, 0.0);
367
368 const scalarField& Vol = mesh.V();
369
370 forAll(mesh.C(), celli)
371 {
372 if (alpha1[celli] < cutoff)
373 {
374 intvDotVapor.value() +=
375 alpha2[celli]*mDotSmear[celli]*Vol[celli];
376 }
377 else if (alpha1[celli] > 1.0 - cutoff)
378 {
379 intvDotLiquid.value() +=
380 alpha1[celli]*mDotSmear[celli]*Vol[celli];
381 }
382 }
383
384 reduce(intvDotVapor.value(), sumOp<scalar>());
385 reduce(intvDotLiquid.value(), sumOp<scalar>());
386
387 //- Calculate Nl and Nv
388 dimensionedScalar Nl ("Nl", dimless, Zero);
389 dimensionedScalar Nv ("Nv", dimless, Zero);
390
391 const dimensionedScalar intmSource0(fvc::domainIntegrate(mDotIn));
392
393 if (intvDotVapor.value() > VSMALL)
394 {
395 Nv = intmSource0/intvDotVapor;
396 }
397 if (intvDotLiquid.value() > VSMALL)
398 {
399 Nl = intmSource0/intvDotLiquid;
400 }
401
402 //- Set source terms in cells with alpha1 < cutoff or alpha1 > 1-cutoff
403 forAll(mesh.C(), celli)
404 {
405 if (alpha1[celli] < cutoff)
406 {
407 mDotOut[celli] = Nv.value()*(1 - alpha1[celli])*mDotSmear[celli];
408 }
409 else if (alpha1[celli] > 1.0 - cutoff)
410 {
411 //mDotOut[celli] = 0;
412 mDotOut[celli] = -Nl.value()*alpha1[celli]*mDotSmear[celli];
413 }
414 else
415 {
416 mDotOut[celli] = 0;
417 }
418 }
419}
420// ************************************************************************* //
Y[inertIndex] max(0.0)
reduce(hasMovingMesh, orOp< bool >())
const volScalarField & alpha1
const volScalarField & alpha2
const dimensionSet & dimensions() const
Return dimensions.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:434
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:81
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const Type & value() const
Return const reference to value.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Class used to pass additional data in.
Definition: smoothData.H:81
scalar maxRatio
Cut off distance.
Definition: smoothData.H:85
Helper class used by the fvc::smooth and fvc::spread functions.
Definition: smoothData.H:57
Helper class used by fvc::sweep function.
Definition: sweepData.H:57
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
rDeltaTY field()
dynamicFvMesh & mesh
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Volume integrate volField creating a volField.
Calculate the matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
void sweep(volScalarField &field, const volScalarField &alpha, const label nLayers, const scalar alphaDiff=0.2)
Definition: fvcSmooth.C:231
void spreadSource(volScalarField &mDotOut, const volScalarField &mDotIn, const volScalarField &alpha1, const volScalarField &alpha2, const dimensionedScalar &D, const scalar cutoff)
Definition: fvcSmooth.C:325
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:44
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:133
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
volScalarField & alpha
const dimensionedScalar & D
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333