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