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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "MULES.H"
30#include "upwind.H"
31#include "fvcSurfaceIntegrate.H"
32#include "localEulerDdtScheme.H"
33#include "slicedSurfaceFields.H"
34#include "wedgeFvPatch.H"
35#include "syncTools.H"
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39template<class RdeltaTType, class RhoType, class SpType, class SuType>
41(
42 const RdeltaTType& rDeltaT,
43 const RhoType& rho,
45 const surfaceScalarField& phiPsi,
46 const SpType& Sp,
47 const SuType& Su
48)
49{
50 Info<< "MULES: Solving for " << psi.name() << endl;
51
52 const fvMesh& mesh = psi.mesh();
53
54 scalarField& psiIf = psi;
55 const scalarField& psi0 = psi.oldTime();
56
57 psiIf = 0.0;
58 fvc::surfaceIntegrate(psiIf, phiPsi);
59
60 if (mesh.moving())
61 {
62 psiIf =
63 (
64 mesh.Vsc0()().field()*rho.oldTime().field()
65 *psi0*rDeltaT/mesh.Vsc()().field()
66 + Su.field()
67 - psiIf
68 )/(rho.field()*rDeltaT - Sp.field());
69 }
70 else
71 {
72 psiIf =
73 (
74 rho.oldTime().field()*psi0*rDeltaT
75 + Su.field()
76 - psiIf
77 )/(rho.field()*rDeltaT - Sp.field());
78 }
79
80 psi.correctBoundaryConditions();
81}
82
83
84template<class RhoType>
86(
87 const RhoType& rho,
89 const surfaceScalarField& phiPsi
90)
91{
92 explicitSolve(rho, psi, phiPsi, zeroField(), zeroField());
93}
94
95
96template<class RhoType, class SpType, class SuType>
98(
99 const RhoType& rho,
101 const surfaceScalarField& phiPsi,
102 const SpType& Sp,
103 const SuType& Su
104)
105{
106 const fvMesh& mesh = psi.mesh();
107
108 if (fv::localEulerDdt::enabled(mesh))
109 {
110 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
111 explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
112 }
113 else
114 {
115 const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
116 explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
117 }
118}
119
120
121template<class RhoType, class PsiMaxType, class PsiMinType>
123(
124 const RhoType& rho,
126 const surfaceScalarField& phiBD,
127 surfaceScalarField& phiPsi,
128 const PsiMaxType& psiMax,
129 const PsiMinType& psiMin
130)
131{
133 (
134 rho,
135 psi,
136 phiBD,
137 phiPsi,
138 zeroField(),
139 zeroField(),
140 psiMax,
141 psiMin
142 );
143}
144
145
146template
147<
148 class RhoType,
149 class SpType,
150 class SuType,
151 class PsiMaxType,
152 class PsiMinType
153>
155(
156 const RhoType& rho,
158 const surfaceScalarField& phi,
159 surfaceScalarField& phiPsi,
160 const SpType& Sp,
161 const SuType& Su,
162 const PsiMaxType& psiMax,
163 const PsiMinType& psiMin
164)
165{
166 const fvMesh& mesh = psi.mesh();
167
168 psi.correctBoundaryConditions();
169
170 if (fv::localEulerDdt::enabled(mesh))
171 {
172 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
173 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, false);
174 explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
175 }
176 else
177 {
178 const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
179 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, false);
180 explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
181 }
182}
183
184
185template
186<
187 class RdeltaTType,
188 class RhoType,
189 class SpType,
190 class SuType,
191 class PsiMaxType,
192 class PsiMinType
193>
195(
196 scalarField& allLambda,
197 const RdeltaTType& rDeltaT,
198 const RhoType& rho,
199 const volScalarField& psi,
200 const surfaceScalarField& phiBD,
201 const surfaceScalarField& phiCorr,
202 const SpType& Sp,
203 const SuType& Su,
204 const PsiMaxType& psiMax,
205 const PsiMinType& psiMin
206)
207{
208 const scalarField& psiIf = psi;
209 const volScalarField::Boundary& psiBf = psi.boundaryField();
210
211 const fvMesh& mesh = psi.mesh();
212
213 const dictionary& MULEScontrols = mesh.solverDict(psi.name());
214
215 const label nLimiterIter
216 (
217 MULEScontrols.getOrDefault<label>("nLimiterIter", 3)
218 );
219
220 const scalar smoothLimiter
221 (
222 MULEScontrols.getOrDefault<scalar>("smoothLimiter", 0)
223 );
224
225 const scalar extremaCoeff
226 (
227 MULEScontrols.getOrDefault<scalar>("extremaCoeff", 0)
228 );
229
230 const scalar boundaryExtremaCoeff
231 (
232 MULEScontrols.getOrDefault<scalar>
233 (
234 "boundaryExtremaCoeff",
235 extremaCoeff
236 )
237 );
238
239 const scalar boundaryDeltaExtremaCoeff
240 (
241 max(boundaryExtremaCoeff - extremaCoeff, 0)
242 );
243
244 const scalarField& psi0 = psi.oldTime();
245
246 const labelUList& owner = mesh.owner();
247 const labelUList& neighb = mesh.neighbour();
249 const scalarField& V = tVsc();
250
251 const scalarField& phiBDIf = phiBD;
252 const surfaceScalarField::Boundary& phiBDBf =
253 phiBD.boundaryField();
254
255 const scalarField& phiCorrIf = phiCorr;
256 const surfaceScalarField::Boundary& phiCorrBf =
257 phiCorr.boundaryField();
258
260 (
262 (
263 "lambda",
264 mesh.time().timeName(),
265 mesh,
266 IOobject::NO_READ,
267 IOobject::NO_WRITE,
268 false
269 ),
270 mesh,
271 dimless,
272 allLambda,
273 false // Use slices for the couples
274 );
275
276 scalarField& lambdaIf = lambda;
277 surfaceScalarField::Boundary& lambdaBf = lambda.boundaryFieldRef();
278
279 scalarField psiMaxn(psiIf.size());
280 scalarField psiMinn(psiIf.size());
281
282 psiMaxn = psiMin;
283 psiMinn = psiMax;
284
285 scalarField sumPhiBD(psiIf.size(), Zero);
286
287 scalarField sumPhip(psiIf.size(), Zero);
288 scalarField mSumPhim(psiIf.size(), Zero);
289
290 forAll(phiCorrIf, facei)
291 {
292 const label own = owner[facei];
293 const label nei = neighb[facei];
294
295 psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
296 psiMinn[own] = min(psiMinn[own], psiIf[nei]);
297
298 psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
299 psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
300
301 sumPhiBD[own] += phiBDIf[facei];
302 sumPhiBD[nei] -= phiBDIf[facei];
303
304 const scalar phiCorrf = phiCorrIf[facei];
305
306 if (phiCorrf > 0)
307 {
308 sumPhip[own] += phiCorrf;
309 mSumPhim[nei] += phiCorrf;
310 }
311 else
312 {
313 mSumPhim[own] -= phiCorrf;
314 sumPhip[nei] -= phiCorrf;
315 }
316 }
317
318 forAll(phiCorrBf, patchi)
319 {
320 const fvPatchScalarField& psiPf = psiBf[patchi];
321 const scalarField& phiBDPf = phiBDBf[patchi];
322 const scalarField& phiCorrPf = phiCorrBf[patchi];
323
324 const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
325
326 if (psiPf.coupled())
327 {
328 const scalarField psiPNf(psiPf.patchNeighbourField());
329
330 forAll(phiCorrPf, pFacei)
331 {
332 const label pfCelli = pFaceCells[pFacei];
333
334 psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
335 psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
336 }
337 }
338 else if (psiPf.fixesValue())
339 {
340 forAll(phiCorrPf, pFacei)
341 {
342 const label pfCelli = pFaceCells[pFacei];
343
344 psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
345 psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
346 }
347 }
348 else
349 {
350 // Add the optional additional allowed boundary extrema
351 if (boundaryDeltaExtremaCoeff > 0)
352 {
353 forAll(phiCorrPf, pFacei)
354 {
355 const label pfCelli = pFaceCells[pFacei];
356
357 const scalar extrema =
358 boundaryDeltaExtremaCoeff
359 *(psiMax[pfCelli] - psiMin[pfCelli]);
360
361 psiMaxn[pfCelli] += extrema;
362 psiMinn[pfCelli] -= extrema;
363 }
364 }
365 }
366
367 forAll(phiCorrPf, pFacei)
368 {
369 const label pfCelli = pFaceCells[pFacei];
370
371 sumPhiBD[pfCelli] += phiBDPf[pFacei];
372
373 const scalar phiCorrf = phiCorrPf[pFacei];
374
375 if (phiCorrf > 0)
376 {
377 sumPhip[pfCelli] += phiCorrf;
378 }
379 else
380 {
381 mSumPhim[pfCelli] -= phiCorrf;
382 }
383 }
384 }
385
386 psiMaxn = min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
387 psiMinn = max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
388
389 if (smoothLimiter > SMALL)
390 {
391 psiMaxn =
392 min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
393 psiMinn =
394 max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
395 }
396
397 if (mesh.moving())
398 {
400
401 psiMaxn =
402 V
403 *(
404 (rho.field()*rDeltaT - Sp.field())*psiMaxn
405 - Su.field()
406 )
407 - (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
408 + sumPhiBD;
409
410 psiMinn =
411 V
412 *(
413 Su.field()
414 - (rho.field()*rDeltaT - Sp.field())*psiMinn
415 )
416 + (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
417 - sumPhiBD;
418 }
419 else
420 {
421 psiMaxn =
422 V
423 *(
424 (rho.field()*rDeltaT - Sp.field())*psiMaxn
425 - Su.field()
426 - (rho.oldTime().field()*rDeltaT)*psi0
427 )
428 + sumPhiBD;
429
430 psiMinn =
431 V
432 *(
433 Su.field()
434 - (rho.field()*rDeltaT - Sp.field())*psiMinn
435 + (rho.oldTime().field()*rDeltaT)*psi0
436 )
437 - sumPhiBD;
438 }
439
440 scalarField sumlPhip(psiIf.size());
441 scalarField mSumlPhim(psiIf.size());
442
443 for (int j=0; j<nLimiterIter; j++)
444 {
445 sumlPhip = 0;
446 mSumlPhim = 0;
447
448 forAll(lambdaIf, facei)
449 {
450 const label own = owner[facei];
451 const label nei = neighb[facei];
452
453 scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
454
455 if (lambdaPhiCorrf > 0)
456 {
457 sumlPhip[own] += lambdaPhiCorrf;
458 mSumlPhim[nei] += lambdaPhiCorrf;
459 }
460 else
461 {
462 mSumlPhim[own] -= lambdaPhiCorrf;
463 sumlPhip[nei] -= lambdaPhiCorrf;
464 }
465 }
466
467 forAll(lambdaBf, patchi)
468 {
469 scalarField& lambdaPf = lambdaBf[patchi];
470 const scalarField& phiCorrfPf = phiCorrBf[patchi];
471
472 const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
473
474 forAll(lambdaPf, pFacei)
475 {
476 const label pfCelli = pFaceCells[pFacei];
477 const scalar lambdaPhiCorrf =
478 lambdaPf[pFacei]*phiCorrfPf[pFacei];
479
480 if (lambdaPhiCorrf > 0)
481 {
482 sumlPhip[pfCelli] += lambdaPhiCorrf;
483 }
484 else
485 {
486 mSumlPhim[pfCelli] -= lambdaPhiCorrf;
487 }
488 }
489 }
490
491 forAll(sumlPhip, celli)
492 {
493 sumlPhip[celli] =
494 max(min
495 (
496 (sumlPhip[celli] + psiMaxn[celli])
497 /(mSumPhim[celli] + ROOTVSMALL),
498 1.0), 0.0
499 );
500
501 mSumlPhim[celli] =
502 max(min
503 (
504 (mSumlPhim[celli] + psiMinn[celli])
505 /(sumPhip[celli] + ROOTVSMALL),
506 1.0), 0.0
507 );
508 }
509
510 const scalarField& lambdam = sumlPhip;
511 const scalarField& lambdap = mSumlPhim;
512
513 forAll(lambdaIf, facei)
514 {
515 if (phiCorrIf[facei] > 0)
516 {
517 lambdaIf[facei] = min
518 (
519 lambdaIf[facei],
520 min(lambdap[owner[facei]], lambdam[neighb[facei]])
521 );
522 }
523 else
524 {
525 lambdaIf[facei] = min
526 (
527 lambdaIf[facei],
528 min(lambdam[owner[facei]], lambdap[neighb[facei]])
529 );
530 }
531 }
532
533 forAll(lambdaBf, patchi)
534 {
535 fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
536 const scalarField& phiCorrfPf = phiCorrBf[patchi];
537 const fvPatchScalarField& psiPf = psiBf[patchi];
538
539 if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
540 {
541 lambdaPf = 0;
542 }
543 else if (psiPf.coupled())
544 {
545 const labelList& pFaceCells =
546 mesh.boundary()[patchi].faceCells();
547
548 forAll(lambdaPf, pFacei)
549 {
550 const label pfCelli = pFaceCells[pFacei];
551
552 if (phiCorrfPf[pFacei] > 0)
553 {
554 lambdaPf[pFacei] =
555 min(lambdaPf[pFacei], lambdap[pfCelli]);
556 }
557 else
558 {
559 lambdaPf[pFacei] =
560 min(lambdaPf[pFacei], lambdam[pfCelli]);
561 }
562 }
563 }
564 }
565
566 syncTools::syncFaceList(mesh, allLambda, minEqOp<scalar>());
567 }
568}
569
570
571template
572<
573 class RdeltaTType,
574 class RhoType,
575 class SpType,
576 class SuType,
577 class PsiMaxType,
578 class PsiMinType
579>
581(
582 const RdeltaTType& rDeltaT,
583 const RhoType& rho,
584 const volScalarField& psi,
585 const surfaceScalarField& phi,
586 surfaceScalarField& phiPsi,
587 const SpType& Sp,
588 const SuType& Su,
589 const PsiMaxType& psiMax,
590 const PsiMinType& psiMin,
591 const bool returnCorr
592)
593{
594 const fvMesh& mesh = psi.mesh();
595
597
599 const surfaceScalarField::Boundary& phiPsiBf = phiPsi.boundaryField();
600
601 forAll(phiBDBf, patchi)
602 {
603 fvsPatchScalarField& phiBDPf = phiBDBf[patchi];
604
605 if (!phiBDPf.coupled())
606 {
607 phiBDPf = phiPsiBf[patchi];
608 }
609 }
610
611 surfaceScalarField& phiCorr = phiPsi;
612 phiCorr -= phiBD;
613
614 scalarField allLambda(mesh.nFaces(), 1.0);
615
617 (
619 (
620 "lambda",
621 mesh.time().timeName(),
622 mesh,
623 IOobject::NO_READ,
624 IOobject::NO_WRITE,
625 false
626 ),
627 mesh,
628 dimless,
629 allLambda,
630 false // Use slices for the couples
631 );
632
633 limiter
634 (
635 allLambda,
636 rDeltaT,
637 rho,
638 psi,
639 phiBD,
640 phiCorr,
641 Sp,
642 Su,
643 psiMax,
644 psiMin
645 );
646
647 if (returnCorr)
648 {
649 phiCorr *= lambda;
650 }
651 else
652 {
653 phiPsi = phiBD + lambda*phiCorr;
654 }
655}
656
657
658template
659<
660 class RhoType,
661 class SpType,
662 class SuType,
663 class PsiMaxType,
664 class PsiMinType
665>
667(
668 const RhoType& rho,
669 const volScalarField& psi,
670 const surfaceScalarField& phi,
671 surfaceScalarField& phiPsi,
672 const SpType& Sp,
673 const SuType& Su,
674 const PsiMaxType& psiMax,
675 const PsiMinType& psiMin,
676 const bool rtnCorr
677)
678{
679 const fvMesh& mesh = psi.mesh();
680
681 if (fv::localEulerDdt::enabled(mesh))
682 {
683 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
684 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, rtnCorr);
685 }
686 else
687 {
688 const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
689 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, rtnCorr);
690 }
691}
692
693
694template<class SurfaceScalarFieldList>
695void Foam::MULES::limitSum(SurfaceScalarFieldList& phiPsiCorrs)
696{
697 {
698 UPtrList<scalarField> phiPsiCorrsInternal(phiPsiCorrs.size());
699 forAll(phiPsiCorrs, phasei)
700 {
701 phiPsiCorrsInternal.set(phasei, &phiPsiCorrs[phasei]);
702 }
703
704 limitSum(phiPsiCorrsInternal);
705 }
706
707 const surfaceScalarField::Boundary& bfld =
708 phiPsiCorrs[0].boundaryField();
709
710 forAll(bfld, patchi)
711 {
712 if (bfld[patchi].coupled())
713 {
714 UPtrList<scalarField> phiPsiCorrsPatch(phiPsiCorrs.size());
715 forAll(phiPsiCorrs, phasei)
716 {
717 phiPsiCorrsPatch.set
718 (
719 phasei,
720 &phiPsiCorrs[phasei].boundaryFieldRef()[patchi]
721 );
722 }
723
724 limitSum(phiPsiCorrsPatch);
725 }
726 }
727}
728
729
730template<class SurfaceScalarFieldList>
732(
733 const SurfaceScalarFieldList& alphas,
734 SurfaceScalarFieldList& phiPsiCorrs,
735 const labelHashSet& fixed
736)
737{
738 {
739 UPtrList<const scalarField> alphasInternal(alphas.size());
740 forAll(alphas, phasei)
741 {
742 alphasInternal.set(phasei, &alphas[phasei]);
743 }
744 UPtrList<scalarField> phiPsiCorrsInternal(phiPsiCorrs.size());
745 forAll(phiPsiCorrs, phasei)
746 {
747 phiPsiCorrsInternal.set(phasei, &phiPsiCorrs[phasei]);
748 }
749
750 limitSum(alphasInternal, phiPsiCorrsInternal, fixed);
751 }
752
753 const surfaceScalarField::Boundary& bfld =
754 phiPsiCorrs[0].boundaryField();
755
756 forAll(bfld, patchi)
757 {
758 if (bfld[patchi].coupled())
759 {
760 UPtrList<const scalarField> alphasPatch(alphas.size());
761 forAll(alphas, phasei)
762 {
763 alphasPatch.set
764 (
765 phasei,
766 &alphas[phasei].boundaryField()[patchi]
767 );
768 }
769 UPtrList<scalarField> phiPsiCorrsPatch(phiPsiCorrs.size());
770 forAll(phiPsiCorrs, phasei)
771 {
772 phiPsiCorrsPatch.set
773 (
774 phasei,
775 &phiPsiCorrs[phasei].boundaryFieldRef()[patchi]
776 );
777 }
778
779 limitSum(alphasPatch, phiPsiCorrsPatch, fixed);
780 }
781 }
782}
783// ************************************************************************* //
MULES: Multidimensional universal limiter for explicit solution.
Y[inertIndex] max(0.0)
surfaceScalarField & phi
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
const T * set(const label i) const
Definition: UPtrList.H:248
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:337
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:453
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:350
virtual bool coupled() const
Return true if this patch field is coupled.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
A class for managing temporary objects.
Definition: tmp.H:65
Upwind differencing scheme class.
Definition: upwind.H:59
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:56
const volScalarField & psi
rDeltaTY field()
pressureControl limit(p)
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
dynamicFvMesh & mesh
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
zeroField Su
Definition: alphaSuSp.H:1
zeroField Sp
Definition: alphaSuSp.H:2
label phasei
Definition: pEqn.H:27
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
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)
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)
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:34
const dimensionSet dimless
Dimensionless.
IOstream & fixed(IOstream &io)
Definition: IOstream.H:458
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:38
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333