CMULESTemplates.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) 2013-2022 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 "CMULES.H"
27 #include "fvcSurfaceIntegrate.H"
28 #include "localEulerDdtScheme.H"
29 #include "slicedSurfaceFields.H"
30 #include "wedgeFvPatch.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 template<class RdeltaTType, class RhoType, class SpType, class SuType>
36 (
37  const RdeltaTType& rDeltaT,
38  const RhoType& rho,
39  volScalarField& psi,
40  const surfaceScalarField& phiCorr,
41  const SpType& Sp,
42  const SuType& Su
43 )
44 {
45  Info<< "MULES: Correcting " << psi.name() << endl;
46 
47  const fvMesh& mesh = psi.mesh();
48 
49  scalarField psiIf(psi.size(), 0);
50  fvc::surfaceIntegrate(psiIf, phiCorr);
51 
52  if (mesh.moving())
53  {
54  psi.primitiveFieldRef() =
55  (
56  rho.field()*psi.primitiveField()*rDeltaT
57  + Su.field()
58  - psiIf
59  )/(rho.field()*rDeltaT - Sp.field());
60  }
61  else
62  {
63  psi.primitiveFieldRef() =
64  (
65  rho.field()*psi.primitiveField()*rDeltaT
66  + Su.field()
67  - psiIf
68  )/(rho.field()*rDeltaT - Sp.field());
69  }
70 
72 }
73 
74 
75 template<class RhoType>
77 (
78  const RhoType& rho,
79  volScalarField& psi,
80  const surfaceScalarField& phiCorr
81 )
82 {
83  correct(rho, psi, phiCorr, zeroField(), zeroField());
84 }
85 
86 
87 template<class RhoType, class SpType, class SuType>
89 (
90  const RhoType& rho,
91  volScalarField& psi,
92  const surfaceScalarField& phiCorr,
93  const SpType& Sp,
94  const SuType& Su
95 )
96 {
97  const fvMesh& mesh = psi.mesh();
98 
99  if (fv::localEulerDdt::enabled(mesh))
100  {
101  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
102  correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
103  }
104  else
105  {
106  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
107  correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
108  }
109 }
110 
111 
112 template<class RhoType, class PsiMaxType, class PsiMinType>
114 (
115  const RhoType& rho,
116  volScalarField& psi,
117  const surfaceScalarField& phi,
118  surfaceScalarField& phiCorr,
119  const PsiMaxType& psiMax,
120  const PsiMinType& psiMin
121 )
122 {
123  correct(rho, psi, phi, phiCorr, zeroField(), zeroField(), psiMax, psiMin);
124 }
125 
126 
127 template
128 <
129  class RhoType,
130  class SpType,
131  class SuType,
132  class PsiMaxType,
133  class PsiMinType
134 >
136 (
137  const RhoType& rho,
138  volScalarField& psi,
139  const surfaceScalarField& phi,
140  surfaceScalarField& phiCorr,
141  const SpType& Sp,
142  const SuType& Su,
143  const PsiMaxType& psiMax,
144  const PsiMinType& psiMin
145 )
146 {
147  const fvMesh& mesh = psi.mesh();
148 
149  if (fv::localEulerDdt::enabled(mesh))
150  {
151  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
152 
153  limitCorr
154  (
155  rDeltaT,
156  rho,
157  psi,
158  phi,
159  phiCorr,
160  Sp,
161  Su,
162  psiMax,
163  psiMin
164  );
165 
166  correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
167  }
168  else
169  {
170  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
171 
172  limitCorr
173  (
174  rDeltaT,
175  rho,
176  psi,
177  phi,
178  phiCorr,
179  Sp,
180  Su,
181  psiMax,
182  psiMin
183  );
184 
185  correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
186  }
187 }
188 
189 
190 template
191 <
192  class RdeltaTType,
193  class RhoType,
194  class SpType,
195  class SuType,
196  class PsiMaxType,
197  class PsiMinType
198 >
200 (
201  surfaceScalarField& lambda,
202  const RdeltaTType& rDeltaT,
203  const RhoType& rho,
204  const volScalarField& psi,
205  const surfaceScalarField& phi,
206  const surfaceScalarField& phiCorr,
207  const SpType& Sp,
208  const SuType& Su,
209  const PsiMaxType& psiMax,
210  const PsiMinType& psiMin
211 )
212 {
213  const scalarField& psiIf = psi;
214  const volScalarField::Boundary& psiBf = psi.boundaryField();
215 
216  const fvMesh& mesh = psi.mesh();
217 
218  const dictionary& MULEScontrols = mesh.solution().solverDict(psi.name());
219 
220  const label nLimiterIter
221  (
222  MULEScontrols.lookup<label>("nLimiterIter")
223  );
224 
225  const scalar smoothLimiter
226  (
227  MULEScontrols.lookupOrDefault<scalar>("smoothLimiter", 0)
228  );
229 
230  const scalar extremaCoeff
231  (
232  MULEScontrols.lookupOrDefault<scalar>("extremaCoeff", 0)
233  );
234 
235  const scalar boundaryExtremaCoeff
236  (
237  MULEScontrols.lookupOrDefault<scalar>
238  (
239  "boundaryExtremaCoeff",
240  extremaCoeff
241  )
242  );
243 
244  const scalar boundaryDeltaExtremaCoeff
245  (
246  max(boundaryExtremaCoeff - extremaCoeff, 0)
247  );
248 
249  const labelUList& owner = mesh.owner();
250  const labelUList& neighb = mesh.neighbour();
251  tmp<volScalarField::Internal> tVsc = mesh.Vsc();
252  const scalarField& V = tVsc();
253 
254  const surfaceScalarField::Boundary& phiBf =
255  phi.boundaryField();
256 
257  const scalarField& phiCorrIf = phiCorr;
258  const surfaceScalarField::Boundary& phiCorrBf =
259  phiCorr.boundaryField();
260 
261  scalarField& lambdaIf = lambda;
262  surfaceScalarField::Boundary& lambdaBf =
263  lambda.boundaryFieldRef();
264 
265  scalarField psiMaxn(psiIf.size());
266  scalarField psiMinn(psiIf.size());
267 
268  psiMaxn = psiMin;
269  psiMinn = psiMax;
270 
271  scalarField sumPhip(psiIf.size(), 0.0);
272  scalarField mSumPhim(psiIf.size(), 0.0);
273 
274  forAll(phiCorrIf, facei)
275  {
276  const label own = owner[facei];
277  const label nei = neighb[facei];
278 
279  psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
280  psiMinn[own] = min(psiMinn[own], psiIf[nei]);
281 
282  psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
283  psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
284 
285  const scalar phiCorrf = phiCorrIf[facei];
286 
287  if (phiCorrf > 0)
288  {
289  sumPhip[own] += phiCorrf;
290  mSumPhim[nei] += phiCorrf;
291  }
292  else
293  {
294  mSumPhim[own] -= phiCorrf;
295  sumPhip[nei] -= phiCorrf;
296  }
297  }
298 
299  forAll(phiCorrBf, patchi)
300  {
301  const fvPatchScalarField& psiPf = psiBf[patchi];
302  const scalarField& phiCorrPf = phiCorrBf[patchi];
303 
304  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
305 
306  if (psiPf.coupled())
307  {
308  const scalarField psiPNf(psiPf.patchNeighbourField());
309 
310  forAll(phiCorrPf, pFacei)
311  {
312  label pfCelli = pFaceCells[pFacei];
313 
314  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
315  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
316  }
317  }
318  else if (psiPf.fixesValue())
319  {
320  forAll(phiCorrPf, pFacei)
321  {
322  const label pfCelli = pFaceCells[pFacei];
323 
324  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
325  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
326  }
327  }
328  else
329  {
330  // Add the optional additional allowed boundary extrema
331  if (boundaryDeltaExtremaCoeff > 0)
332  {
333  forAll(phiCorrPf, pFacei)
334  {
335  const label pfCelli = pFaceCells[pFacei];
336 
337  const scalar extrema =
338  boundaryDeltaExtremaCoeff
339  *(psiMax[pfCelli] - psiMin[pfCelli]);
340 
341  psiMaxn[pfCelli] += extrema;
342  psiMinn[pfCelli] -= extrema;
343  }
344  }
345  }
346 
347  forAll(phiCorrPf, pFacei)
348  {
349  const label pfCelli = pFaceCells[pFacei];
350 
351  const scalar phiCorrf = phiCorrPf[pFacei];
352 
353  if (phiCorrf > 0)
354  {
355  sumPhip[pfCelli] += phiCorrf;
356  }
357  else
358  {
359  mSumPhim[pfCelli] -= phiCorrf;
360  }
361  }
362  }
363 
364  psiMaxn = min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
365  psiMinn = max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
366 
367  if (smoothLimiter > small)
368  {
369  psiMaxn =
370  min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
371  psiMinn =
372  max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
373  }
374 
375  psiMaxn =
376  V
377  *(
378  (rho.field()*rDeltaT - Sp.field())*psiMaxn
379  - Su.field()
380  - rho.field()*psi.primitiveField()*rDeltaT
381  );
382 
383  psiMinn =
384  V
385  *(
386  Su.field()
387  - (rho.field()*rDeltaT - Sp.field())*psiMinn
388  + rho.field()*psi.primitiveField()*rDeltaT
389  );
390 
391  scalarField sumlPhip(psiIf.size());
392  scalarField mSumlPhim(psiIf.size());
393 
394  for (int j=0; j<nLimiterIter; j++)
395  {
396  sumlPhip = 0;
397  mSumlPhim = 0;
398 
399  forAll(lambdaIf, facei)
400  {
401  const label own = owner[facei];
402  const label nei = neighb[facei];
403 
404  const scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
405 
406  if (lambdaPhiCorrf > 0)
407  {
408  sumlPhip[own] += lambdaPhiCorrf;
409  mSumlPhim[nei] += lambdaPhiCorrf;
410  }
411  else
412  {
413  mSumlPhim[own] -= lambdaPhiCorrf;
414  sumlPhip[nei] -= lambdaPhiCorrf;
415  }
416  }
417 
418  forAll(lambdaBf, patchi)
419  {
420  scalarField& lambdaPf = lambdaBf[patchi];
421  const scalarField& phiCorrfPf = phiCorrBf[patchi];
422 
423  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
424 
425  forAll(lambdaPf, pFacei)
426  {
427  label pfCelli = pFaceCells[pFacei];
428 
429  scalar lambdaPhiCorrf = lambdaPf[pFacei]*phiCorrfPf[pFacei];
430 
431  if (lambdaPhiCorrf > 0)
432  {
433  sumlPhip[pfCelli] += lambdaPhiCorrf;
434  }
435  else
436  {
437  mSumlPhim[pfCelli] -= lambdaPhiCorrf;
438  }
439  }
440  }
441 
442  forAll(sumlPhip, celli)
443  {
444  sumlPhip[celli] =
445  max(min
446  (
447  (sumlPhip[celli] + psiMaxn[celli])
448  /(mSumPhim[celli] + rootVSmall),
449  1.0), 0.0
450  );
451 
452  mSumlPhim[celli] =
453  max(min
454  (
455  (mSumlPhim[celli] + psiMinn[celli])
456  /(sumPhip[celli] + rootVSmall),
457  1.0), 0.0
458  );
459  }
460 
461  const scalarField& lambdam = sumlPhip;
462  const scalarField& lambdap = mSumlPhim;
463 
464  forAll(lambdaIf, facei)
465  {
466  if (phiCorrIf[facei] > 0)
467  {
468  lambdaIf[facei] = min
469  (
470  lambdaIf[facei],
471  min(lambdap[owner[facei]], lambdam[neighb[facei]])
472  );
473  }
474  else
475  {
476  lambdaIf[facei] = min
477  (
478  lambdaIf[facei],
479  min(lambdam[owner[facei]], lambdap[neighb[facei]])
480  );
481  }
482  }
483 
484 
485  forAll(lambdaBf, patchi)
486  {
487  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
488  const scalarField& phiCorrfPf = phiCorrBf[patchi];
489  const fvPatchScalarField& psiPf = psiBf[patchi];
490 
491  if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
492  {
493  lambdaPf = 0;
494  }
495  else if (psiPf.coupled())
496  {
497  const labelList& pFaceCells =
498  mesh.boundary()[patchi].faceCells();
499 
500  forAll(lambdaPf, pFacei)
501  {
502  const label pfCelli = pFaceCells[pFacei];
503 
504  if (phiCorrfPf[pFacei] > 0)
505  {
506  lambdaPf[pFacei] =
507  min(lambdaPf[pFacei], lambdap[pfCelli]);
508  }
509  else
510  {
511  lambdaPf[pFacei] =
512  min(lambdaPf[pFacei], lambdam[pfCelli]);
513  }
514  }
515  }
516  else
517  {
518  const labelList& pFaceCells =
519  mesh.boundary()[patchi].faceCells();
520  const scalarField& phiPf = phiBf[patchi];
521 
522  forAll(lambdaPf, pFacei)
523  {
524  // Limit outlet faces only
525  if ((phiPf[pFacei] + phiCorrfPf[pFacei]) > small*small)
526  {
527  const label pfCelli = pFaceCells[pFacei];
528 
529  if (phiCorrfPf[pFacei] > 0)
530  {
531  lambdaPf[pFacei] =
532  min(lambdaPf[pFacei], lambdap[pfCelli]);
533  }
534  else
535  {
536  lambdaPf[pFacei] =
537  min(lambdaPf[pFacei], lambdam[pfCelli]);
538  }
539  }
540  }
541  }
542  }
543 
544  // Take minimum of value across coupled patches
545  surfaceScalarField::Boundary lambdaNbrBf
546  (
547  surfaceScalarField::Internal::null(),
548  lambdaBf.boundaryNeighbourField()
549  );
550  forAll(lambdaBf, patchi)
551  {
552  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
553  const fvsPatchScalarField& lambdaNbrPf = lambdaNbrBf[patchi];
554  if (lambdaPf.coupled())
555  {
556  lambdaPf = min(lambdaPf, lambdaNbrPf);
557  }
558  }
559  }
560 }
561 
562 
563 template
564 <
565  class RdeltaTType,
566  class RhoType,
567  class SpType,
568  class SuType,
569  class PsiMaxType,
570  class PsiMinType
571 >
573 (
574  const RdeltaTType& rDeltaT,
575  const RhoType& rho,
576  const volScalarField& psi,
577  const surfaceScalarField& phi,
578  surfaceScalarField& phiCorr,
579  const SpType& Sp,
580  const SuType& Su,
581  const PsiMaxType& psiMax,
582  const PsiMinType& psiMin
583 )
584 {
585  const fvMesh& mesh = psi.mesh();
586 
588  (
589  IOobject
590  (
591  "lambda",
592  mesh.time().timeName(),
593  mesh,
594  IOobject::NO_READ,
595  IOobject::NO_WRITE,
596  false
597  ),
598  mesh,
600  );
601 
603  (
604  lambda,
605  rDeltaT,
606  rho,
607  psi,
608  phi,
609  phiCorr,
610  Sp,
611  Su,
612  psiMax,
613  psiMin
614  );
615 
616  phiCorr *= lambda;
617 }
618 
619 
620 // ************************************************************************* //
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void limiterCorr(surfaceScalarField &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)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const word & name() const
Return name.
Definition: IOobject.H:315
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
bool moving() const
Is mesh moving.
Definition: polyMesh.H:525
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:310
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensionedScalar lambda(viscosity->lookup("lambda"))
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const fvSolution & solution() const
Return the fvSchemes.
Definition: fvMesh.C:1683
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
const dimensionSet dimless
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:335
fvMesh & mesh
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:423
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
virtual bool coupled() const
Return true if this patch field is coupled.
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:50
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
const volScalarField & psi
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:323
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
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
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
const Mesh & mesh() const
Return mesh.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:417
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
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)
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:53
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:433
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800