MULESlimiter.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-2025 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 "MULES.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "wedgeFvPatch.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 template
34 <
35  class RdeltaTType,
36  class RhoType,
37  class SpType,
38  class PsiMaxType,
39  class PsiMinType
40 >
42 (
43  const control& controls,
45  const RdeltaTType& rDeltaT,
46  const RhoType& rho,
47  const volScalarField& psi,
48  const scalarField& SuCorr,
49  const surfaceScalarField& phiBD,
50  const surfaceScalarField& phiCorr,
51  const SpType& Sp,
52  const PsiMaxType& psiMax,
53  const PsiMinType& psiMin
54 )
55 {
56  const scalarField& psiIf = psi;
57  const volScalarField::Boundary& psiBf = psi.boundaryField();
58 
59  const fvMesh& mesh = psi.mesh();
60 
61  const scalar boundaryDeltaExtremaCoeff
62  (
63  max(controls.boundaryExtremaCoeff - controls.extremaCoeff, 0)
64  );
65 
66  const labelUList& owner = mesh.owner();
67  const labelUList& neighb = mesh.neighbour();
69  const scalarField& V = tVsc();
70 
71  const surfaceScalarField::Boundary& phiBDBf = phiBD.boundaryField();
72 
73  const scalarField& phiCorrIf = phiCorr;
74  const surfaceScalarField::Boundary& phiCorrBf = phiCorr.boundaryField();
75 
76  scalarField& lambdaIf = lambda;
77  surfaceScalarField::Boundary& lambdaBf = lambda.boundaryFieldRef();
78 
79  scalarField psiMaxn(psiIf.size());
80  scalarField psiMinn(psiIf.size());
81 
82  scalarField phiCorrNorm;
83  if (controls.tol != 0)
84  {
85  phiCorrNorm = (V*(rho.primitiveField()*rDeltaT - Sp.primitiveField()));
86  }
87 
88  scalarField sumPhip(psiIf.size(), 0.0);
89  scalarField mSumPhim(psiIf.size(), 0.0);
90 
91  if (controls.globalBounds)
92  {
93  psiMaxn = psiMax;
94  psiMinn = psiMin;
95 
96  forAll(phiCorrIf, facei)
97  {
98  const label own = owner[facei];
99  const label nei = neighb[facei];
100 
101  const scalar phiCorrf = phiCorrIf[facei];
102 
103  if (phiCorrf > 0)
104  {
105  sumPhip[own] += phiCorrf;
106  mSumPhim[nei] += phiCorrf;
107  }
108  else
109  {
110  mSumPhim[own] -= phiCorrf;
111  sumPhip[nei] -= phiCorrf;
112  }
113  }
114 
115  forAll(phiCorrBf, patchi)
116  {
117  const scalarField& phiCorrPf = phiCorrBf[patchi];
118  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
119 
120  forAll(phiCorrPf, pFacei)
121  {
122  const label pfCelli = pFaceCells[pFacei];
123  const scalar phiCorrf = phiCorrPf[pFacei];
124 
125  if (phiCorrf > 0)
126  {
127  sumPhip[pfCelli] += phiCorrf;
128  }
129  else
130  {
131  mSumPhim[pfCelli] -= phiCorrf;
132  }
133  }
134  }
135  }
136  else
137  {
138  psiMaxn = psiMin;
139  psiMinn = psiMax;
140 
141  forAll(phiCorrIf, facei)
142  {
143  const label own = owner[facei];
144  const label nei = neighb[facei];
145 
146  psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
147  psiMinn[own] = min(psiMinn[own], psiIf[nei]);
148 
149  psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
150  psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
151 
152  const scalar phiCorrf = phiCorrIf[facei];
153 
154  if (phiCorrf > 0)
155  {
156  sumPhip[own] += phiCorrf;
157  mSumPhim[nei] += phiCorrf;
158  }
159  else
160  {
161  mSumPhim[own] -= phiCorrf;
162  sumPhip[nei] -= phiCorrf;
163  }
164  }
165 
166  forAll(phiCorrBf, patchi)
167  {
168  const fvPatchScalarField& psiPf = psiBf[patchi];
169  const scalarField& phiCorrPf = phiCorrBf[patchi];
170 
171  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
172 
173  if (psiPf.coupled())
174  {
175  const scalarField psiPNf(psiPf.patchNeighbourField());
176 
177  forAll(phiCorrPf, pFacei)
178  {
179  const label pfCelli = pFaceCells[pFacei];
180 
181  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
182  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
183  }
184  }
185  else if (psiPf.fixesValue())
186  {
187  forAll(phiCorrPf, pFacei)
188  {
189  const label pfCelli = pFaceCells[pFacei];
190 
191  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
192  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
193  }
194  }
195  else
196  {
197  // Add the optional additional allowed boundary extrema
198  if (boundaryDeltaExtremaCoeff > 0)
199  {
200  forAll(phiCorrPf, pFacei)
201  {
202  const label pfCelli = pFaceCells[pFacei];
203 
204  const scalar extrema =
205  boundaryDeltaExtremaCoeff
206  *(psiMax[pfCelli] - psiMin[pfCelli]);
207 
208  psiMaxn[pfCelli] += extrema;
209  psiMinn[pfCelli] -= extrema;
210  }
211  }
212  }
213 
214  forAll(phiCorrPf, pFacei)
215  {
216  const label pfCelli = pFaceCells[pFacei];
217  const scalar phiCorrf = phiCorrPf[pFacei];
218 
219  if (phiCorrf > 0)
220  {
221  sumPhip[pfCelli] += phiCorrf;
222  }
223  else
224  {
225  mSumPhim[pfCelli] -= phiCorrf;
226  }
227  }
228  }
229 
230  if (controls.extremaCoeff > 0)
231  {
232  psiMaxn = min
233  (
234  psiMaxn + controls.extremaCoeff*(psiMax - psiMin),
235  psiMax
236  );
237 
238  psiMinn = max
239  (
240  psiMinn - controls.extremaCoeff*(psiMax - psiMin),
241  psiMin
242  );
243  }
244  else
245  {
246  psiMaxn = min(psiMaxn, psiMax);
247  psiMinn = max(psiMinn, psiMin);
248  }
249 
250  if (controls.smoothingCoeff > small)
251  {
252  psiMaxn = min
253  (
254  controls.smoothingCoeff*psiIf
255  + (1.0 - controls.smoothingCoeff)*psiMaxn,
256  psiMax
257  );
258 
259  psiMinn = max
260  (
261  controls.smoothingCoeff*psiIf
262  + (1.0 - controls.smoothingCoeff)*psiMinn,
263  psiMin
264  );
265  }
266  }
267 
268 
269  psiMaxn =
270  V*((rho.primitiveField()*rDeltaT - Sp.primitiveField())*psiMaxn)
271  - SuCorr;
272 
273  psiMinn =
274  SuCorr
275  - V*((rho.primitiveField()*rDeltaT - Sp.primitiveField())*psiMinn);
276 
277  scalarField sumlPhip(psiIf.size());
278  scalarField mSumlPhim(psiIf.size());
279 
280  // Allocate storage for lambda0 on coupled patches
281  // for optional convergence test
283  if (controls.tol != 0)
284  {
285  forAll(lambdaBf, patchi)
286  {
287  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
288 
289  if (lambdaPf.coupled())
290  {
291  lambdaBf0.set
292  (
293  patchi,
295  (
296  mesh.boundary()[patchi],
297  surfaceScalarField::Internal::null()
298  )
299  );
300  }
301  }
302  }
303 
304  for (int j=0; j<controls.nIter; j++)
305  {
306  // Convergence test parameter
307  scalar maxDeltaLambdaPhiCorrRes = 0;
308 
309  // Sum limited positive and negative fluxes
310  // Not needed for first iteration
311  if (j > 0)
312  {
313  sumlPhip = 0;
314  mSumlPhim = 0;
315 
316  forAll(lambdaIf, facei)
317  {
318  const label own = owner[facei];
319  const label nei = neighb[facei];
320 
321  const scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
322 
323  if (lambdaPhiCorrf > 0)
324  {
325  sumlPhip[own] += lambdaPhiCorrf;
326  mSumlPhim[nei] += lambdaPhiCorrf;
327  }
328  else
329  {
330  mSumlPhim[own] -= lambdaPhiCorrf;
331  sumlPhip[nei] -= lambdaPhiCorrf;
332  }
333  }
334 
335  forAll(lambdaBf, patchi)
336  {
337  scalarField& lambdaPf = lambdaBf[patchi];
338  const scalarField& phiCorrfPf = phiCorrBf[patchi];
339 
340  const labelList& pFaceCells =
341  mesh.boundary()[patchi].faceCells();
342 
343  forAll(lambdaPf, pFacei)
344  {
345  const label pfCelli = pFaceCells[pFacei];
346  const scalar lambdaPhiCorrf =
347  lambdaPf[pFacei]*phiCorrfPf[pFacei];
348 
349  if (lambdaPhiCorrf > 0)
350  {
351  sumlPhip[pfCelli] += lambdaPhiCorrf;
352  }
353  else
354  {
355  mSumlPhim[pfCelli] -= lambdaPhiCorrf;
356  }
357  }
358  }
359  }
360 
361  // Reuse storage of sumlPhip and mSumlPhim for lambdam and lambdap
362  scalarField& lambdam = sumlPhip;
363  scalarField& lambdap = mSumlPhim;
364 
365  if (j == 0)
366  {
367  forAll(lambdam, celli)
368  {
369  lambdam[celli] =
370  max(min
371  (
372  psiMaxn[celli]/(mSumPhim[celli] + rootVSmall),
373  1.0), 0.0
374  );
375 
376  lambdap[celli] =
377  max(min
378  (
379  psiMinn[celli]/(sumPhip[celli] + rootVSmall),
380  1.0), 0.0
381  );
382  }
383  }
384  else
385  {
386  forAll(lambdam, celli)
387  {
388  lambdam[celli] =
389  max(min
390  (
391  (sumlPhip[celli] + psiMaxn[celli])
392  /(mSumPhim[celli] + rootVSmall),
393  1.0), 0.0
394  );
395 
396  lambdap[celli] =
397  max(min
398  (
399  (mSumlPhim[celli] + psiMinn[celli])
400  /(sumPhip[celli] + rootVSmall),
401  1.0), 0.0
402  );
403  }
404  }
405 
406  forAll(lambdaIf, facei)
407  {
408  const scalar lambdaIf0 = lambdaIf[facei];
409 
410  if (phiCorrIf[facei] > 0)
411  {
412  lambdaIf[facei] =
413  min(lambdap[owner[facei]], lambdam[neighb[facei]]);
414  }
415  else
416  {
417  lambdaIf[facei] =
418  min(lambdam[owner[facei]], lambdap[neighb[facei]]);
419  }
420 
421  if (controls.tol > 0)
422  {
423  const scalar phiCorrRes =
424  mag(phiCorrIf[facei])
425  /min(phiCorrNorm[owner[facei]], phiCorrNorm[neighb[facei]]);
426 
427  if (phiCorrRes > controls.tol)
428  {
429  maxDeltaLambdaPhiCorrRes = max
430  (
431  maxDeltaLambdaPhiCorrRes,
432  mag(lambdaIf[facei] - lambdaIf0)*phiCorrRes
433  );
434  }
435  }
436  }
437 
438  // Take minimum of value across coupled patches
439  forAll(lambdaBf, patchi)
440  {
441  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
442  const scalarField& phiCorrfPf = phiCorrBf[patchi];
443  const fvPatchScalarField& psiPf = psiBf[patchi];
444 
445  if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
446  {
447  lambdaPf = 0;
448  }
449  else if (psiPf.coupled())
450  {
451  const labelList& pFaceCells =
452  mesh.boundary()[patchi].faceCells();
453 
454  if (controls.tol > 0)
455  {
456  lambdaBf0[patchi] = lambdaPf;
457  }
458 
459  forAll(lambdaPf, pFacei)
460  {
461  const label pfCelli = pFaceCells[pFacei];
462 
463  if (phiCorrfPf[pFacei] > 0)
464  {
465  lambdaPf[pFacei] = lambdap[pfCelli];
466  }
467  else
468  {
469  lambdaPf[pFacei] = lambdam[pfCelli];
470  }
471  }
472  }
473  else
474  {
475  const labelList& pFaceCells =
476  mesh.boundary()[patchi].faceCells();
477 
478  const scalarField& phiBDPf = phiBDBf[patchi];
479 
480  forAll(lambdaPf, pFacei)
481  {
482  // Limit outlet faces only
483  if ((phiBDPf[pFacei] + phiCorrfPf[pFacei]) > small*small)
484  {
485  const label pfCelli = pFaceCells[pFacei];
486  const scalar lambdaPf0 = lambdaPf[pFacei];
487 
488  if (phiCorrfPf[pFacei] > 0)
489  {
490  lambdaPf[pFacei] = lambdap[pfCelli];
491  }
492  else
493  {
494  lambdaPf[pFacei] = lambdam[pfCelli];
495  }
496 
497  if (controls.tol > 0)
498  {
499  const scalar phiCorrRes =
500  mag(phiCorrfPf[pFacei])/phiCorrNorm[pfCelli];
501 
502  if (phiCorrRes > controls.tol)
503  {
504  maxDeltaLambdaPhiCorrRes = max
505  (
506  maxDeltaLambdaPhiCorrRes,
507  mag(lambdaPf[pFacei] - lambdaPf0)*phiCorrRes
508  );
509  }
510  }
511  }
512  }
513  }
514  }
515 
516  // Take minimum value of limiter across coupled patches
517  surfaceScalarField::Boundary lambdaNbrBf
518  (
519  surfaceScalarField::Internal::null(),
520  lambdaBf.boundaryNeighbourField()
521  );
522 
523  forAll(lambdaBf, patchi)
524  {
525  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
526 
527  if (lambdaPf.coupled())
528  {
529  const fvsPatchScalarField& lambdaNbrPf = lambdaNbrBf[patchi];
530  lambdaPf = min(lambdaPf, lambdaNbrPf);
531 
532  if (controls.tol > 0)
533  {
534  const fvsPatchScalarField& lambdaPf0 = lambdaBf0[patchi];
535  const scalarField& phiCorrfPf = phiCorrBf[patchi];
536 
537  const labelList& pFaceCells =
538  mesh.boundary()[patchi].faceCells();
539 
540  forAll(lambdaPf, pFacei)
541  {
542  const scalar phiCorrRes =
543  mag(phiCorrfPf[pFacei])
544  /phiCorrNorm[pFaceCells[pFacei]];
545 
546  if (phiCorrRes > controls.tol)
547  {
548  maxDeltaLambdaPhiCorrRes = max
549  (
550  maxDeltaLambdaPhiCorrRes,
551  mag(lambdaPf[pFacei] - lambdaPf0[pFacei])
552  *phiCorrRes
553  );
554  }
555  }
556  }
557  }
558  }
559 
560  // Optional convergence test
561  if (controls.tol != 0)
562  {
563  reduce(maxDeltaLambdaPhiCorrRes, maxOp<scalar>());
564 
565  if (debug)
566  {
567  Info<< "MULES: maxDeltaLambdaPhiCorrRes "
568  << maxDeltaLambdaPhiCorrRes << endl;
569  }
570 
571  if (maxDeltaLambdaPhiCorrRes < controls.tol) break;
572  }
573  }
574 }
575 
576 
577 // ************************************************************************* //
MULES: Multidimensional universal limiter for explicit solution.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricBoundaryField class.
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
Foam::calculatedFvsPatchField.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:471
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:962
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:477
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:327
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:340
virtual tmp< Field< Type > > patchNeighbourField(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:445
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:86
virtual bool coupled() const
Return true if this patch field is coupled.
A class for managing temporary objects.
Definition: tmp.H:55
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
const volScalarField & psi
dimensionedScalar lambda(viscosity->lookup("lambda"))
void limiter(const control &controls, surfaceScalarField &lambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const scalarField &SuCorr, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const PsiMaxType &psiMax, const PsiMinType &psiMin)
Definition: MULESlimiter.C:42
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Switch globalBounds
Optional switch to select global bounds only.
Definition: MULES.H:80
scalar smoothingCoeff
Optional coefficient to reduce the allowed range of the solution to.
Definition: MULES.H:103
scalar boundaryExtremaCoeff
Optional coefficient to relax the local boundedness constraint.
Definition: MULES.H:95
label nIter
Optional maximum number of limiter iterations.
Definition: MULES.H:70
scalar extremaCoeff
Optional coefficient to relax the local boundedness constraint.
Definition: MULES.H:90
scalar tol
Optional limiter convergence tolerance.
Definition: MULES.H:76
Foam::surfaceFields.