MULESTemplates.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-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "MULES.H"
29 #include "upwind.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& phiPsi,
45  const SpType& Sp,
46  const SuType& Su
47 )
48 {
49  Info<< "MULES: Solving for " << psi.name() << endl;
50 
51  const fvMesh& mesh = psi.mesh();
52 
53  scalarField& psiIf = psi;
54  const scalarField& psi0 = psi.oldTime();
55 
56  psiIf = 0.0;
57  fvc::surfaceIntegrate(psiIf, phiPsi);
58 
59  if (mesh.moving())
60  {
61  psiIf =
62  (
63  mesh.Vsc0()().field()*rho.oldTime().field()
64  *psi0*rDeltaT/mesh.Vsc()().field()
65  + Su.field()
66  - psiIf
67  )/(rho.field()*rDeltaT - Sp.field());
68  }
69  else
70  {
71  psiIf =
72  (
73  rho.oldTime().field()*psi0*rDeltaT
74  + Su.field()
75  - psiIf
76  )/(rho.field()*rDeltaT - Sp.field());
77  }
78 
79  psi.correctBoundaryConditions();
80 }
81 
82 
83 template<class RhoType>
85 (
86  const RhoType& rho,
88  const surfaceScalarField& phiPsi
89 )
90 {
91  explicitSolve(rho, psi, phiPsi, zeroField(), zeroField());
92 }
93 
94 
95 template<class RhoType, class SpType, class SuType>
97 (
98  const RhoType& rho,
100  const surfaceScalarField& phiPsi,
101  const SpType& Sp,
102  const SuType& Su
103 )
104 {
105  const fvMesh& mesh = psi.mesh();
106 
107  if (fv::localEulerDdt::enabled(mesh))
108  {
109  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
110  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
111  }
112  else
113  {
114  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
115  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
116  }
117 }
118 
119 
120 template<class RhoType, class PsiMaxType, class PsiMinType>
122 (
123  const RhoType& rho,
125  const surfaceScalarField& phiBD,
126  surfaceScalarField& phiPsi,
127  const PsiMaxType& psiMax,
128  const PsiMinType& psiMin
129 )
130 {
132  (
133  rho,
134  psi,
135  phiBD,
136  phiPsi,
137  zeroField(),
138  zeroField(),
139  psiMax,
140  psiMin
141  );
142 }
143 
144 
145 template
146 <
147  class RhoType,
148  class SpType,
149  class SuType,
150  class PsiMaxType,
151  class PsiMinType
152 >
154 (
155  const RhoType& rho,
157  const surfaceScalarField& phi,
158  surfaceScalarField& phiPsi,
159  const SpType& Sp,
160  const SuType& Su,
161  const PsiMaxType& psiMax,
162  const PsiMinType& psiMin
163 )
164 {
165  const fvMesh& mesh = psi.mesh();
166 
167  psi.correctBoundaryConditions();
168 
169  if (fv::localEulerDdt::enabled(mesh))
170  {
171  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
172  limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, false);
173  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
174  }
175  else
176  {
177  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
178  limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, false);
179  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
180  }
181 }
182 
183 
184 template
185 <
186  class RdeltaTType,
187  class RhoType,
188  class SpType,
189  class SuType,
190  class PsiMaxType,
191  class PsiMinType
192 >
194 (
195  scalarField& allLambda,
196  const RdeltaTType& rDeltaT,
197  const RhoType& rho,
198  const volScalarField& psi,
199  const surfaceScalarField& phiBD,
200  const surfaceScalarField& phiCorr,
201  const SpType& Sp,
202  const SuType& Su,
203  const PsiMaxType& psiMax,
204  const PsiMinType& psiMin
205 )
206 {
207  const scalarField& psiIf = psi;
208  const volScalarField::Boundary& psiBf = psi.boundaryField();
209 
210  const fvMesh& mesh = psi.mesh();
211 
212  const dictionary& MULEScontrols = mesh.solverDict(psi.name());
213 
214  const label nLimiterIter
215  (
216  MULEScontrols.lookupOrDefault<label>("nLimiterIter", 3)
217  );
218 
219  const scalar smoothLimiter
220  (
221  MULEScontrols.lookupOrDefault<scalar>("smoothLimiter", 0)
222  );
223 
224  const scalar extremaCoeff
225  (
226  MULEScontrols.lookupOrDefault<scalar>("extremaCoeff", 0)
227  );
228 
229  const scalar boundaryExtremaCoeff
230  (
231  MULEScontrols.lookupOrDefault<scalar>
232  (
233  "boundaryExtremaCoeff",
234  extremaCoeff
235  )
236  );
237 
238  const scalar boundaryDeltaExtremaCoeff
239  (
240  max(boundaryExtremaCoeff - extremaCoeff, 0)
241  );
242 
243  const scalarField& psi0 = psi.oldTime();
244 
245  const labelUList& owner = mesh.owner();
246  const labelUList& neighb = mesh.neighbour();
247  tmp<volScalarField::Internal> tVsc = mesh.Vsc();
248  const scalarField& V = tVsc();
249 
250  const scalarField& phiBDIf = phiBD;
251  const surfaceScalarField::Boundary& phiBDBf =
252  phiBD.boundaryField();
253 
254  const scalarField& phiCorrIf = phiCorr;
255  const surfaceScalarField::Boundary& phiCorrBf =
256  phiCorr.boundaryField();
257 
259  (
260  IOobject
261  (
262  "lambda",
263  mesh.time().timeName(),
264  mesh,
265  IOobject::NO_READ,
266  IOobject::NO_WRITE,
267  false
268  ),
269  mesh,
270  dimless,
271  allLambda,
272  false // Use slices for the couples
273  );
274 
275  scalarField& lambdaIf = lambda;
276  surfaceScalarField::Boundary& lambdaBf = lambda.boundaryFieldRef();
277 
278  scalarField psiMaxn(psiIf.size());
279  scalarField psiMinn(psiIf.size());
280 
281  psiMaxn = psiMin;
282  psiMinn = psiMax;
283 
284  scalarField sumPhiBD(psiIf.size(), Zero);
285 
286  scalarField sumPhip(psiIf.size(), Zero);
287  scalarField mSumPhim(psiIf.size(), Zero);
288 
289  forAll(phiCorrIf, facei)
290  {
291  const label own = owner[facei];
292  const label nei = neighb[facei];
293 
294  psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
295  psiMinn[own] = min(psiMinn[own], psiIf[nei]);
296 
297  psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
298  psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
299 
300  sumPhiBD[own] += phiBDIf[facei];
301  sumPhiBD[nei] -= phiBDIf[facei];
302 
303  const scalar phiCorrf = phiCorrIf[facei];
304 
305  if (phiCorrf > 0)
306  {
307  sumPhip[own] += phiCorrf;
308  mSumPhim[nei] += phiCorrf;
309  }
310  else
311  {
312  mSumPhim[own] -= phiCorrf;
313  sumPhip[nei] -= phiCorrf;
314  }
315  }
316 
317  forAll(phiCorrBf, patchi)
318  {
319  const fvPatchScalarField& psiPf = psiBf[patchi];
320  const scalarField& phiBDPf = phiBDBf[patchi];
321  const scalarField& phiCorrPf = phiCorrBf[patchi];
322 
323  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
324 
325  if (psiPf.coupled())
326  {
327  const scalarField psiPNf(psiPf.patchNeighbourField());
328 
329  forAll(phiCorrPf, pFacei)
330  {
331  const label pfCelli = pFaceCells[pFacei];
332 
333  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
334  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
335  }
336  }
337  else if (psiPf.fixesValue())
338  {
339  forAll(phiCorrPf, pFacei)
340  {
341  const label pfCelli = pFaceCells[pFacei];
342 
343  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
344  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
345  }
346  }
347  else
348  {
349  // Add the optional additional allowed boundary extrema
350  if (boundaryDeltaExtremaCoeff > 0)
351  {
352  forAll(phiCorrPf, pFacei)
353  {
354  const label pfCelli = pFaceCells[pFacei];
355 
356  const scalar extrema =
357  boundaryDeltaExtremaCoeff
358  *(psiMax[pfCelli] - psiMin[pfCelli]);
359 
360  psiMaxn[pfCelli] += extrema;
361  psiMinn[pfCelli] -= extrema;
362  }
363  }
364  }
365 
366  forAll(phiCorrPf, pFacei)
367  {
368  const label pfCelli = pFaceCells[pFacei];
369 
370  sumPhiBD[pfCelli] += phiBDPf[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  if (mesh.moving())
397  {
399 
400  psiMaxn =
401  V
402  *(
403  (rho.field()*rDeltaT - Sp.field())*psiMaxn
404  - Su.field()
405  )
406  - (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
407  + sumPhiBD;
408 
409  psiMinn =
410  V
411  *(
412  Su.field()
413  - (rho.field()*rDeltaT - Sp.field())*psiMinn
414  )
415  + (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
416  - sumPhiBD;
417  }
418  else
419  {
420  psiMaxn =
421  V
422  *(
423  (rho.field()*rDeltaT - Sp.field())*psiMaxn
424  - Su.field()
425  - (rho.oldTime().field()*rDeltaT)*psi0
426  )
427  + sumPhiBD;
428 
429  psiMinn =
430  V
431  *(
432  Su.field()
433  - (rho.field()*rDeltaT - Sp.field())*psiMinn
434  + (rho.oldTime().field()*rDeltaT)*psi0
435  )
436  - sumPhiBD;
437  }
438 
439  scalarField sumlPhip(psiIf.size());
440  scalarField mSumlPhim(psiIf.size());
441 
442  for (int j=0; j<nLimiterIter; j++)
443  {
444  sumlPhip = 0;
445  mSumlPhim = 0;
446 
447  forAll(lambdaIf, facei)
448  {
449  const label own = owner[facei];
450  const label nei = neighb[facei];
451 
452  scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
453 
454  if (lambdaPhiCorrf > 0)
455  {
456  sumlPhip[own] += lambdaPhiCorrf;
457  mSumlPhim[nei] += lambdaPhiCorrf;
458  }
459  else
460  {
461  mSumlPhim[own] -= lambdaPhiCorrf;
462  sumlPhip[nei] -= lambdaPhiCorrf;
463  }
464  }
465 
466  forAll(lambdaBf, patchi)
467  {
468  scalarField& lambdaPf = lambdaBf[patchi];
469  const scalarField& phiCorrfPf = phiCorrBf[patchi];
470 
471  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
472 
473  forAll(lambdaPf, pFacei)
474  {
475  const label pfCelli = pFaceCells[pFacei];
476  const scalar lambdaPhiCorrf =
477  lambdaPf[pFacei]*phiCorrfPf[pFacei];
478 
479  if (lambdaPhiCorrf > 0)
480  {
481  sumlPhip[pfCelli] += lambdaPhiCorrf;
482  }
483  else
484  {
485  mSumlPhim[pfCelli] -= lambdaPhiCorrf;
486  }
487  }
488  }
489 
490  forAll(sumlPhip, celli)
491  {
492  sumlPhip[celli] =
493  max(min
494  (
495  (sumlPhip[celli] + psiMaxn[celli])
496  /(mSumPhim[celli] + ROOTVSMALL),
497  1.0), 0.0
498  );
499 
500  mSumlPhim[celli] =
501  max(min
502  (
503  (mSumlPhim[celli] + psiMinn[celli])
504  /(sumPhip[celli] + ROOTVSMALL),
505  1.0), 0.0
506  );
507  }
508 
509  const scalarField& lambdam = sumlPhip;
510  const scalarField& lambdap = mSumlPhim;
511 
512  forAll(lambdaIf, facei)
513  {
514  if (phiCorrIf[facei] > 0)
515  {
516  lambdaIf[facei] = min
517  (
518  lambdaIf[facei],
519  min(lambdap[owner[facei]], lambdam[neighb[facei]])
520  );
521  }
522  else
523  {
524  lambdaIf[facei] = min
525  (
526  lambdaIf[facei],
527  min(lambdam[owner[facei]], lambdap[neighb[facei]])
528  );
529  }
530  }
531 
532  forAll(lambdaBf, patchi)
533  {
534  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
535  const scalarField& phiCorrfPf = phiCorrBf[patchi];
536  const fvPatchScalarField& psiPf = psiBf[patchi];
537 
538  if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
539  {
540  lambdaPf = 0;
541  }
542  else if (psiPf.coupled())
543  {
544  const labelList& pFaceCells =
545  mesh.boundary()[patchi].faceCells();
546 
547  forAll(lambdaPf, pFacei)
548  {
549  const label pfCelli = pFaceCells[pFacei];
550 
551  if (phiCorrfPf[pFacei] > 0)
552  {
553  lambdaPf[pFacei] =
554  min(lambdaPf[pFacei], lambdap[pfCelli]);
555  }
556  else
557  {
558  lambdaPf[pFacei] =
559  min(lambdaPf[pFacei], lambdam[pfCelli]);
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& phiPsi,
586  const SpType& Sp,
587  const SuType& Su,
588  const PsiMaxType& psiMax,
589  const PsiMinType& psiMin,
590  const bool returnCorr
591 )
592 {
593  const fvMesh& mesh = psi.mesh();
594 
595  surfaceScalarField phiBD(upwind<scalar>(psi.mesh(), phi).flux(psi));
596 
597  surfaceScalarField::Boundary& phiBDBf = phiBD.boundaryFieldRef();
598  const surfaceScalarField::Boundary& phiPsiBf = phiPsi.boundaryField();
599 
600  forAll(phiBDBf, patchi)
601  {
602  fvsPatchScalarField& phiBDPf = phiBDBf[patchi];
603 
604  if (!phiBDPf.coupled())
605  {
606  phiBDPf = phiPsiBf[patchi];
607  }
608  }
609 
610  surfaceScalarField& phiCorr = phiPsi;
611  phiCorr -= phiBD;
612 
613  scalarField allLambda(mesh.nFaces(), 1.0);
614 
616  (
617  IOobject
618  (
619  "lambda",
620  mesh.time().timeName(),
621  mesh,
622  IOobject::NO_READ,
623  IOobject::NO_WRITE,
624  false
625  ),
626  mesh,
627  dimless,
628  allLambda,
629  false // Use slices for the couples
630  );
631 
632  limiter
633  (
634  allLambda,
635  rDeltaT,
636  rho,
637  psi,
638  phiBD,
639  phiCorr,
640  Sp,
641  Su,
642  psiMax,
643  psiMin
644  );
645 
646  if (returnCorr)
647  {
648  phiCorr *= lambda;
649  }
650  else
651  {
652  phiPsi = phiBD + lambda*phiCorr;
653  }
654 }
655 
656 
657 template
658 <
659  class RhoType,
660  class SpType,
661  class SuType,
662  class PsiMaxType,
663  class PsiMinType
664 >
666 (
667  const RhoType& rho,
668  const volScalarField& psi,
669  const surfaceScalarField& phi,
670  surfaceScalarField& phiPsi,
671  const SpType& Sp,
672  const SuType& Su,
673  const PsiMaxType& psiMax,
674  const PsiMinType& psiMin,
675  const bool rtnCorr
676 )
677 {
678  const fvMesh& mesh = psi.mesh();
679 
680  if (fv::localEulerDdt::enabled(mesh))
681  {
682  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
683  limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, rtnCorr);
684  }
685  else
686  {
687  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
688  limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, rtnCorr);
689  }
690 }
691 
692 
693 template<class SurfaceScalarFieldList>
694 void Foam::MULES::limitSum(SurfaceScalarFieldList& phiPsiCorrs)
695 {
696  {
697  UPtrList<scalarField> phiPsiCorrsInternal(phiPsiCorrs.size());
698  forAll(phiPsiCorrs, phasei)
699  {
700  phiPsiCorrsInternal.set(phasei, &phiPsiCorrs[phasei]);
701  }
702 
703  limitSum(phiPsiCorrsInternal);
704  }
705 
706  const surfaceScalarField::Boundary& bfld =
707  phiPsiCorrs[0].boundaryField();
708 
709  forAll(bfld, patchi)
710  {
711  if (bfld[patchi].coupled())
712  {
713  UPtrList<scalarField> phiPsiCorrsPatch(phiPsiCorrs.size());
714  forAll(phiPsiCorrs, phasei)
715  {
716  phiPsiCorrsPatch.set
717  (
718  phasei,
719  &phiPsiCorrs[phasei].boundaryFieldRef()[patchi]
720  );
721  }
722 
723  limitSum(phiPsiCorrsPatch);
724  }
725  }
726 }
727 
728 
729 template<class SurfaceScalarFieldList>
731 (
732  const SurfaceScalarFieldList& alphas,
733  SurfaceScalarFieldList& phiPsiCorrs,
734  const labelHashSet& fixed
735 )
736 {
737  {
738  UPtrList<const scalarField> alphasInternal(alphas.size());
739  forAll(alphas, phasei)
740  {
741  alphasInternal.set(phasei, &alphas[phasei]);
742  }
743  UPtrList<scalarField> phiPsiCorrsInternal(phiPsiCorrs.size());
744  forAll(phiPsiCorrs, phasei)
745  {
746  phiPsiCorrsInternal.set(phasei, &phiPsiCorrs[phasei]);
747  }
748 
749  limitSum(alphasInternal, phiPsiCorrsInternal, fixed);
750  }
751 
752  const surfaceScalarField::Boundary& bfld =
753  phiPsiCorrs[0].boundaryField();
754 
755  forAll(bfld, patchi)
756  {
757  if (bfld[patchi].coupled())
758  {
759  UPtrList<const scalarField> alphasPatch(alphas.size());
760  forAll(alphas, phasei)
761  {
762  alphasPatch.set
763  (
764  phasei,
765  &alphas[phasei].boundaryField()[patchi]
766  );
767  }
768  UPtrList<scalarField> phiPsiCorrsPatch(phiPsiCorrs.size());
769  forAll(phiPsiCorrs, phasei)
770  {
771  phiPsiCorrsPatch.set
772  (
773  phasei,
774  &phiPsiCorrs[phasei].boundaryFieldRef()[patchi]
775  );
776  }
777 
778  limitSum(alphasPatch, phiPsiCorrsPatch, fixed);
779  }
780  }
781 }
782 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::MULES::explicitSolve
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Definition: MULESTemplates.C:40
Foam::MULES::limit
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
Definition: MULESTemplates.C:580
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::upwind
Upwind differencing scheme class.
Definition: upwind.H:56
Foam::vtk::Tools::zeroField
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Definition: foamVtkToolsTemplates.C:258
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
phasei
label phasei
Definition: pEqn.H:27
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::fvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:318
Foam::HashSet< label, Hash< label > >
rho
rho
Definition: readInitialConditions.H:96
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:290
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:331
localEulerDdtScheme.H
Foam::label
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:62
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::MULES::limitSum
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:35
upwind.H
Foam::UPtrList
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:63
field
rDeltaTY field()
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::limiter
void limiter(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
Definition: MULESTemplates.C:194
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::fixed
IOstream & fixed(IOstream &io)
Definition: IOstream.H:441
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
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:84
Foam::fvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:434
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
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::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:34
Foam::GeometricField< scalar, fvPatchField, volMesh >
MULES.H
MULES: Multidimensional universal limiter for explicit solution.
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::limiter
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:37
Foam::fvsPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvsPatchField.H:310
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
wedgeFvPatch.H