CMULESTemplates.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) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "CMULES.H"
30 #include "fvcSurfaceIntegrate.H"
31 #include "localEulerDdtScheme.H"
32 #include "slicedSurfaceFields.H"
33 #include "wedgeFvPatch.H"
34 #include "syncTools.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 template<class RdeltaTType, class RhoType, class SpType, class SuType>
40 (
41  const RdeltaTType& rDeltaT,
42  const RhoType& rho,
44  const surfaceScalarField& phiCorr,
45  const SpType& Sp,
46  const SuType& Su
47 )
48 {
49  Info<< "MULES: Correcting " << psi.name() << endl;
50 
51  const fvMesh& mesh = psi.mesh();
52 
53  scalarField psiIf(psi.size(), Zero);
54  fvc::surfaceIntegrate(psiIf, phiCorr);
55 
56  if (mesh.moving())
57  {
58  psi.primitiveFieldRef() =
59  (
60  rho.field()*psi.primitiveField()*rDeltaT
61  + Su.field()
62  - psiIf
63  )/(rho.field()*rDeltaT - Sp.field());
64  }
65  else
66  {
67  psi.primitiveFieldRef() =
68  (
69  rho.field()*psi.primitiveField()*rDeltaT
70  + Su.field()
71  - psiIf
72  )/(rho.field()*rDeltaT - Sp.field());
73  }
74 
75  psi.correctBoundaryConditions();
76 }
77 
78 
79 template<class RhoType>
81 (
82  const RhoType& rho,
84  const surfaceScalarField& phiCorr
85 )
86 {
87  correct(rho, psi, phiCorr, zeroField(), zeroField());
88 }
89 
90 
91 template<class RhoType, class SpType, class SuType>
93 (
94  const RhoType& rho,
96  const surfaceScalarField& phiCorr,
97  const SpType& Sp,
98  const SuType& Su
99 )
100 {
101  const fvMesh& mesh = psi.mesh();
102 
103  if (fv::localEulerDdt::enabled(mesh))
104  {
105  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
106  correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
107  }
108  else
109  {
110  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
111  correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
112  }
113 }
114 
115 
116 template<class RhoType, class PsiMaxType, class PsiMinType>
118 (
119  const RhoType& rho,
121  const surfaceScalarField& phi,
122  surfaceScalarField& phiCorr,
123  const PsiMaxType& psiMax,
124  const PsiMinType& psiMin
125 )
126 {
127  correct(rho, psi, phi, phiCorr, zeroField(), zeroField(), psiMax, psiMin);
128 }
129 
130 
131 template
132 <
133  class RhoType,
134  class SpType,
135  class SuType,
136  class PsiMaxType,
137  class PsiMinType
138 >
140 (
141  const RhoType& rho,
143  const surfaceScalarField& phi,
144  surfaceScalarField& phiCorr,
145  const SpType& Sp,
146  const SuType& Su,
147  const PsiMaxType& psiMax,
148  const PsiMinType& psiMin
149 )
150 {
151  const fvMesh& mesh = psi.mesh();
152 
153  if (fv::localEulerDdt::enabled(mesh))
154  {
155  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
156 
157  limitCorr
158  (
159  rDeltaT,
160  rho,
161  psi,
162  phi,
163  phiCorr,
164  Sp,
165  Su,
166  psiMax,
167  psiMin
168  );
169 
170  correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
171  }
172  else
173  {
174  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
175 
176  limitCorr
177  (
178  rDeltaT,
179  rho,
180  psi,
181  phi,
182  phiCorr,
183  Sp,
184  Su,
185  psiMax,
186  psiMin
187  );
188 
189  correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
190  }
191 }
192 
193 
194 template
195 <
196  class RdeltaTType,
197  class RhoType,
198  class SpType,
199  class SuType,
200  class PsiMaxType,
201  class PsiMinType
202 >
204 (
205  scalarField& allLambda,
206  const RdeltaTType& rDeltaT,
207  const RhoType& rho,
208  const volScalarField& psi,
209  const surfaceScalarField& phi,
210  const surfaceScalarField& phiCorr,
211  const SpType& Sp,
212  const SuType& Su,
213  const PsiMaxType& psiMax,
214  const PsiMinType& psiMin
215 )
216 {
217  const scalarField& psiIf = psi;
218  const volScalarField::Boundary& psiBf = psi.boundaryField();
219 
220  const fvMesh& mesh = psi.mesh();
221 
222  const dictionary& MULEScontrols = mesh.solverDict(psi.name());
223 
224  const label nLimiterIter
225  (
226  MULEScontrols.get<label>("nLimiterIter")
227  );
228 
229  const scalar smoothLimiter
230  (
231  MULEScontrols.getOrDefault<scalar>("smoothLimiter", 0)
232  );
233 
234  const scalar extremaCoeff
235  (
236  MULEScontrols.getOrDefault<scalar>("extremaCoeff", 0)
237  );
238 
239  const scalar boundaryExtremaCoeff
240  (
241  MULEScontrols.getOrDefault<scalar>
242  (
243  "boundaryExtremaCoeff",
244  extremaCoeff
245  )
246  );
247 
248  const scalar boundaryDeltaExtremaCoeff
249  (
250  max(boundaryExtremaCoeff - extremaCoeff, 0)
251  );
252 
253  const labelUList& owner = mesh.owner();
254  const labelUList& neighb = mesh.neighbour();
255  tmp<volScalarField::Internal> tVsc = mesh.Vsc();
256  const scalarField& V = tVsc();
257 
258  const surfaceScalarField::Boundary& phiBf =
259  phi.boundaryField();
260 
261  const scalarField& phiCorrIf = phiCorr;
262  const surfaceScalarField::Boundary& phiCorrBf =
263  phiCorr.boundaryField();
264 
266  (
267  IOobject
268  (
269  "lambda",
270  mesh.time().timeName(),
271  mesh,
272  IOobject::NO_READ,
273  IOobject::NO_WRITE,
274  false
275  ),
276  mesh,
277  dimless,
278  allLambda,
279  false // Use slices for the couples
280  );
281 
282  scalarField& lambdaIf = lambda;
283  surfaceScalarField::Boundary& lambdaBf =
284  lambda.boundaryFieldRef();
285 
286  scalarField psiMaxn(psiIf.size());
287  scalarField psiMinn(psiIf.size());
288 
289  psiMaxn = psiMin;
290  psiMinn = psiMax;
291 
292  scalarField sumPhip(psiIf.size(), Zero);
293  scalarField mSumPhim(psiIf.size(), Zero);
294 
295  forAll(phiCorrIf, facei)
296  {
297  const label own = owner[facei];
298  const label nei = neighb[facei];
299 
300  psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
301  psiMinn[own] = min(psiMinn[own], psiIf[nei]);
302 
303  psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
304  psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
305 
306  const scalar phiCorrf = phiCorrIf[facei];
307 
308  if (phiCorrf > 0)
309  {
310  sumPhip[own] += phiCorrf;
311  mSumPhim[nei] += phiCorrf;
312  }
313  else
314  {
315  mSumPhim[own] -= phiCorrf;
316  sumPhip[nei] -= phiCorrf;
317  }
318  }
319 
320  forAll(phiCorrBf, patchi)
321  {
322  const fvPatchScalarField& psiPf = psiBf[patchi];
323  const scalarField& phiCorrPf = phiCorrBf[patchi];
324 
325  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
326 
327  if (psiPf.coupled())
328  {
329  const scalarField psiPNf(psiPf.patchNeighbourField());
330 
331  forAll(phiCorrPf, pFacei)
332  {
333  label pfCelli = pFaceCells[pFacei];
334 
335  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
336  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
337  }
338  }
339  else if (psiPf.fixesValue())
340  {
341  forAll(phiCorrPf, pFacei)
342  {
343  const label pfCelli = pFaceCells[pFacei];
344 
345  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
346  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
347  }
348  }
349  else
350  {
351  // Add the optional additional allowed boundary extrema
352  if (boundaryDeltaExtremaCoeff > 0)
353  {
354  forAll(phiCorrPf, pFacei)
355  {
356  const label pfCelli = pFaceCells[pFacei];
357 
358  const scalar extrema =
359  boundaryDeltaExtremaCoeff
360  *(psiMax[pfCelli] - psiMin[pfCelli]);
361 
362  psiMaxn[pfCelli] += extrema;
363  psiMinn[pfCelli] -= extrema;
364  }
365  }
366  }
367 
368  forAll(phiCorrPf, pFacei)
369  {
370  const label pfCelli = pFaceCells[pFacei];
371 
372  const scalar phiCorrf = phiCorrPf[pFacei];
373 
374  if (phiCorrf > 0)
375  {
376  sumPhip[pfCelli] += phiCorrf;
377  }
378  else
379  {
380  mSumPhim[pfCelli] -= phiCorrf;
381  }
382  }
383  }
384 
385  psiMaxn = min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
386  psiMinn = max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
387 
388  if (smoothLimiter > SMALL)
389  {
390  psiMaxn =
391  min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
392  psiMinn =
393  max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
394  }
395 
396  psiMaxn =
397  V
398  *(
399  (rho.field()*rDeltaT - Sp.field())*psiMaxn
400  - Su.field()
401  - rho.field()*psi.primitiveField()*rDeltaT
402  );
403 
404  psiMinn =
405  V
406  *(
407  Su.field()
408  - (rho.field()*rDeltaT - Sp.field())*psiMinn
409  + rho.field()*psi.primitiveField()*rDeltaT
410  );
411 
412  scalarField sumlPhip(psiIf.size());
413  scalarField mSumlPhim(psiIf.size());
414 
415  for (int j=0; j<nLimiterIter; j++)
416  {
417  sumlPhip = 0;
418  mSumlPhim = 0;
419 
420  forAll(lambdaIf, facei)
421  {
422  const label own = owner[facei];
423  const label nei = neighb[facei];
424 
425  const scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
426 
427  if (lambdaPhiCorrf > 0)
428  {
429  sumlPhip[own] += lambdaPhiCorrf;
430  mSumlPhim[nei] += lambdaPhiCorrf;
431  }
432  else
433  {
434  mSumlPhim[own] -= lambdaPhiCorrf;
435  sumlPhip[nei] -= lambdaPhiCorrf;
436  }
437  }
438 
439  forAll(lambdaBf, patchi)
440  {
441  scalarField& lambdaPf = lambdaBf[patchi];
442  const scalarField& phiCorrfPf = phiCorrBf[patchi];
443 
444  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
445 
446  forAll(lambdaPf, pFacei)
447  {
448  label pfCelli = pFaceCells[pFacei];
449 
450  scalar lambdaPhiCorrf = lambdaPf[pFacei]*phiCorrfPf[pFacei];
451 
452  if (lambdaPhiCorrf > 0)
453  {
454  sumlPhip[pfCelli] += lambdaPhiCorrf;
455  }
456  else
457  {
458  mSumlPhim[pfCelli] -= lambdaPhiCorrf;
459  }
460  }
461  }
462 
463  forAll(sumlPhip, celli)
464  {
465  sumlPhip[celli] =
466  max(min
467  (
468  (sumlPhip[celli] + psiMaxn[celli])
469  /(mSumPhim[celli] + ROOTVSMALL),
470  1.0), 0.0
471  );
472 
473  mSumlPhim[celli] =
474  max(min
475  (
476  (mSumlPhim[celli] + psiMinn[celli])
477  /(sumPhip[celli] + ROOTVSMALL),
478  1.0), 0.0
479  );
480  }
481 
482  const scalarField& lambdam = sumlPhip;
483  const scalarField& lambdap = mSumlPhim;
484 
485  forAll(lambdaIf, facei)
486  {
487  if (phiCorrIf[facei] > 0)
488  {
489  lambdaIf[facei] = min
490  (
491  lambdaIf[facei],
492  min(lambdap[owner[facei]], lambdam[neighb[facei]])
493  );
494  }
495  else
496  {
497  lambdaIf[facei] = min
498  (
499  lambdaIf[facei],
500  min(lambdam[owner[facei]], lambdap[neighb[facei]])
501  );
502  }
503  }
504 
505 
506  forAll(lambdaBf, patchi)
507  {
508  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
509  const scalarField& phiCorrfPf = phiCorrBf[patchi];
510  const fvPatchScalarField& psiPf = psiBf[patchi];
511 
512  if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
513  {
514  lambdaPf = 0;
515  }
516  else if (psiPf.coupled())
517  {
518  const labelList& pFaceCells =
519  mesh.boundary()[patchi].faceCells();
520 
521  forAll(lambdaPf, pFacei)
522  {
523  const label pfCelli = pFaceCells[pFacei];
524 
525  if (phiCorrfPf[pFacei] > 0)
526  {
527  lambdaPf[pFacei] =
528  min(lambdaPf[pFacei], lambdap[pfCelli]);
529  }
530  else
531  {
532  lambdaPf[pFacei] =
533  min(lambdaPf[pFacei], lambdam[pfCelli]);
534  }
535  }
536  }
537  else
538  {
539  const labelList& pFaceCells =
540  mesh.boundary()[patchi].faceCells();
541  const scalarField& phiPf = phiBf[patchi];
542 
543  forAll(lambdaPf, pFacei)
544  {
545  // Limit outlet faces only
546  if ((phiPf[pFacei] + phiCorrfPf[pFacei]) > SMALL*SMALL)
547  {
548  const label pfCelli = pFaceCells[pFacei];
549 
550  if (phiCorrfPf[pFacei] > 0)
551  {
552  lambdaPf[pFacei] =
553  min(lambdaPf[pFacei], lambdap[pfCelli]);
554  }
555  else
556  {
557  lambdaPf[pFacei] =
558  min(lambdaPf[pFacei], lambdam[pfCelli]);
559  }
560  }
561  }
562  }
563  }
564 
565  syncTools::syncFaceList(mesh, allLambda, minEqOp<scalar>());
566  }
567 }
568 
569 
570 template
571 <
572  class RdeltaTType,
573  class RhoType,
574  class SpType,
575  class SuType,
576  class PsiMaxType,
577  class PsiMinType
578 >
580 (
581  const RdeltaTType& rDeltaT,
582  const RhoType& rho,
583  const volScalarField& psi,
584  const surfaceScalarField& phi,
585  surfaceScalarField& phiCorr,
586  const SpType& Sp,
587  const SuType& Su,
588  const PsiMaxType& psiMax,
589  const PsiMinType& psiMin
590 )
591 {
592  const fvMesh& mesh = psi.mesh();
593 
594  scalarField allLambda(mesh.nFaces(), 1.0);
595 
597  (
598  IOobject
599  (
600  "lambda",
601  mesh.time().timeName(),
602  mesh,
603  IOobject::NO_READ,
604  IOobject::NO_WRITE,
605  false
606  ),
607  mesh,
608  dimless,
609  allLambda,
610  false // Use slices for the couples
611  );
612 
614  (
615  allLambda,
616  rDeltaT,
617  rho,
618  psi,
619  phi,
620  phiCorr,
621  Sp,
622  Su,
623  psiMax,
624  psiMin
625  );
626 
627  phiCorr *= lambda;
628 }
629 
630 
631 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
CMULES.H
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::vtk::Tools::zeroField
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Definition: foamVtkToolsTemplates.C:327
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Sp
zeroField Sp
Definition: alphaSuSp.H:2
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::fvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:332
rho
rho
Definition: readInitialConditions.H:88
syncTools.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::minEqOp
Definition: ops.H:81
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::fvPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:345
localEulerDdtScheme.H
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
correct
fvOptions correct(rho)
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
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
fvcSurfaceIntegrate.H
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
Foam::MULES::correct
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
Definition: CMULESTemplates.C:40
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::fvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:448
Foam::MULES::limiterCorr
void limiterCorr(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
Definition: CMULESTemplates.C:204
Foam::MULES::limitCorr
void limitCorr(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
Definition: CMULESTemplates.C:580
Foam::List< label >
Foam::fvc::surfaceIntegrate
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcSurfaceIntegrate.C:46
Foam::UList< label >
slicedSurfaceFields.H
Foam::SlicedGeometricField
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
Definition: slicedSurfaceFieldsFwd.H:58
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::GeometricField< scalar, fvPatchField, volMesh >
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
wedgeFvPatch.H