multiphaseSystem.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-2018 OpenFOAM Foundation
9 Copyright (C) 2020-2022 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 "multiphaseSystem.H"
30#include "alphaContactAngleFvPatchScalarField.H"
32#include "Time.H"
33#include "subCycle.H"
34#include "MULES.H"
35#include "surfaceInterpolate.H"
36#include "fvcGrad.H"
37#include "fvcSnGrad.H"
38#include "fvcDiv.H"
39#include "fvcFlux.H"
40#include "fvcAverage.H"
41
42#include "unitConversion.H"
43
44// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45
46void Foam::multiphaseSystem::calcAlphas()
47{
48 scalar level = 0.0;
49 alphas_ == 0.0;
50
51 for (const phaseModel& phase : phases_)
52 {
53 alphas_ += level * phase;
54 level += 1.0;
55 }
56}
57
58
59void Foam::multiphaseSystem::solveAlphas()
60{
61 PtrList<surfaceScalarField> alphaPhiCorrs(phases_.size());
62 int phasei = 0;
63
64 for (phaseModel& phase : phases_)
65 {
66 volScalarField& alpha1 = phase;
67
68 alphaPhiCorrs.set
69 (
70 phasei,
72 (
73 "phi" + alpha1.name() + "Corr",
75 (
76 phi_,
77 phase,
78 "div(phi," + alpha1.name() + ')'
79 )
80 )
81 );
82
83 surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phasei];
84
85 for (phaseModel& phase2 : phases_)
86 {
88
89 if (&phase2 == &phase) continue;
90
91 surfaceScalarField phir(phase.phi() - phase2.phi());
92
93 const auto cAlpha = cAlphas_.cfind(interfacePair(phase, phase2));
94
95 if (cAlpha.found())
96 {
98 (
99 (mag(phi_) + mag(phir))/mesh_.magSf()
100 );
101
102 phir += min(cAlpha()*phic, max(phic))*nHatf(phase, phase2);
103 }
104
105 const word phirScheme
106 (
107 "div(phir," + alpha2.name() + ',' + alpha1.name() + ')'
108 );
109
110 alphaPhiCorr += fvc::flux
111 (
112 -fvc::flux(-phir, phase2, phirScheme),
113 phase,
114 phirScheme
115 );
116 }
117
118 phase.correctInflowOutflow(alphaPhiCorr);
119
121 (
122 1.0/mesh_.time().deltaT().value(),
123 geometricOneField(),
124 phase,
125 phi_,
126 alphaPhiCorr,
127 zeroField(),
128 zeroField(),
129 oneField(),
130 zeroField(),
131 true
132 );
133
134 ++phasei;
135 }
136
137 MULES::limitSum(alphaPhiCorrs);
138
139 volScalarField sumAlpha
140 (
141 IOobject
142 (
143 "sumAlpha",
144 mesh_.time().timeName(),
145 mesh_
146 ),
147 mesh_,
148 dimensionedScalar("sumAlpha", dimless, 0)
149 );
150
151 phasei = 0;
152
153 for (phaseModel& phase : phases_)
154 {
155 surfaceScalarField& alphaPhi = alphaPhiCorrs[phasei];
156 alphaPhi += upwind<scalar>(mesh_, phi_).flux(phase);
157 phase.correctInflowOutflow(alphaPhi);
158
160 (
161 geometricOneField(),
162 phase,
164 );
165
166 phase.alphaPhi() = alphaPhi;
167
168 Info<< phase.name() << " volume fraction, min, max = "
169 << phase.weightedAverage(mesh_.V()).value()
170 << ' ' << min(phase).value()
171 << ' ' << max(phase).value()
172 << endl;
173
174 sumAlpha += phase;
175
176 ++phasei;
177 }
178
179 Info<< "Phase-sum volume fraction, min, max = "
180 << sumAlpha.weightedAverage(mesh_.V()).value()
181 << ' ' << min(sumAlpha).value()
182 << ' ' << max(sumAlpha).value()
183 << endl;
184
185 // Correct the sum of the phase-fractions to avoid 'drift'
186 volScalarField sumCorr(1.0 - sumAlpha);
187 for (phaseModel& phase : phases_)
188 {
189 volScalarField& alpha = phase;
190 alpha += alpha*sumCorr;
191 }
192
193 calcAlphas();
194}
195
196
197Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseSystem::nHatfv
198(
199 const volScalarField& alpha1,
201) const
202{
203 /*
204 // Cell gradient of alpha
205 volVectorField gradAlpha =
206 alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
207
208 // Interpolated face-gradient of alpha
209 surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
210 */
211
212 surfaceVectorField gradAlphaf
213 (
216 );
217
218 // Face unit interface normal
219 return gradAlphaf/(mag(gradAlphaf) + deltaN_);
220}
221
222
224(
225 const volScalarField& alpha1,
227) const
228{
229 // Face unit interface normal flux
230 return nHatfv(alpha1, alpha2) & mesh_.Sf();
231}
232
233
234// Correction for the boundary condition on the unit normal nHat on
235// walls to produce the correct contact angle.
236
237// The dynamic contact angle is calculated from the component of the
238// velocity on the direction of the interface, parallel to the wall.
239
240void Foam::multiphaseSystem::correctContactAngle
241(
242 const phaseModel& phase1,
243 const phaseModel& phase2,
245) const
246{
247 const volScalarField::Boundary& gbf
249
250 const fvBoundaryMesh& boundary = mesh_.boundary();
251
252 forAll(boundary, patchi)
253 {
254 if
255 (
256 isA<multiphaseEuler::alphaContactAngleFvPatchScalarField>
257 (
258 gbf[patchi]
259 )
260 )
261 {
262 const auto& acap =
263 refCast
264 <
265 const multiphaseEuler::alphaContactAngleFvPatchScalarField
266 >
267 (
268 gbf[patchi]
269 );
270
271 vectorField& nHatPatch = nHatb[patchi];
272
273 vectorField AfHatPatch
274 (
275 mesh_.Sf().boundaryField()[patchi]
276 /mesh_.magSf().boundaryField()[patchi]
277 );
278
279 const auto tp =
280 acap.thetaProps().cfind(interfacePair(phase1, phase2));
281
282 if (!tp.found())
283 {
285 << "Cannot find interface " << interfacePair(phase1, phase2)
286 << "\n in table of theta properties for patch "
287 << acap.patch().name()
288 << exit(FatalError);
289 }
290
291 bool matched = (tp.key().first() == phase1.name());
292
293 const scalar theta0 = degToRad(tp().theta0(matched));
294 scalarField theta(boundary[patchi].size(), theta0);
295
296 scalar uTheta = tp().uTheta();
297
298 // Calculate the dynamic contact angle if required
299 if (uTheta > SMALL)
300 {
301 const scalar thetaA = degToRad(tp().thetaA(matched));
302 const scalar thetaR = degToRad(tp().thetaR(matched));
303
304 // Calculated the component of the velocity parallel to the wall
305 vectorField Uwall
306 (
307 phase1.U().boundaryField()[patchi].patchInternalField()
308 - phase1.U().boundaryField()[patchi]
309 );
310 Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
311
312 // Find the direction of the interface parallel to the wall
313 vectorField nWall
314 (
315 nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
316 );
317
318 // Normalise nWall
319 nWall /= (mag(nWall) + SMALL);
320
321 // Calculate Uwall resolved normal to the interface parallel to
322 // the interface
323 scalarField uwall(nWall & Uwall);
324
325 theta += (thetaA - thetaR)*tanh(uwall/uTheta);
326 }
327
328
329 // Reset nHatPatch to correspond to the contact angle
330
331 scalarField a12(nHatPatch & AfHatPatch);
332
333 scalarField b1(cos(theta));
334
335 scalarField b2(nHatPatch.size());
336
337 forAll(b2, facei)
338 {
339 b2[facei] = cos(acos(a12[facei]) - theta[facei]);
340 }
341
342 scalarField det(1.0 - a12*a12);
343
344 scalarField a((b1 - a12*b2)/det);
345 scalarField b((b2 - a12*b1)/det);
346
347 nHatPatch = a*AfHatPatch + b*nHatPatch;
348
349 nHatPatch /= (mag(nHatPatch) + deltaN_.value());
350 }
351 }
352}
353
354
356(
357 const phaseModel& phase1,
358 const phaseModel& phase2
359) const
360{
361 tmp<surfaceVectorField> tnHatfv = nHatfv(phase1, phase2);
362
363 correctContactAngle(phase1, phase2, tnHatfv.ref().boundaryFieldRef());
364
365 // Simple expression for curvature
366 return -fvc::div(tnHatfv & mesh_.Sf());
367}
368
369
370// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
371
373(
374 const volVectorField& U,
376)
377:
379 (
381 (
382 "transportProperties",
383 U.time().constant(),
384 U.db(),
385 IOobject::MUST_READ_IF_MODIFIED,
386 IOobject::NO_WRITE
387 )
388 ),
389
390 phases_(lookup("phases"), phaseModel::iNew(U.mesh())),
391
392 mesh_(U.mesh()),
393 phi_(phi),
394
395 alphas_
396 (
398 (
399 "alphas",
400 mesh_.time().timeName(),
401 mesh_,
402 IOobject::NO_READ,
403 IOobject::AUTO_WRITE
404 ),
405 mesh_,
406 dimensionedScalar("alphas", dimless, 0.0)
407 ),
408
409 sigmas_(lookup("sigmas")),
410 dimSigma_(1, 0, -2, 0, 0),
411 cAlphas_(lookup("interfaceCompression")),
412 Cvms_(lookup("virtualMass")),
413 deltaN_
414 (
415 "deltaN",
416 1e-8/cbrt(average(mesh_.V()))
417 )
418{
419 calcAlphas();
420 alphas_.write();
421
422 interfaceDictTable dragModelsDict(lookup("drag"));
423
424 forAllConstIters(dragModelsDict, iter)
425 {
426 dragModels_.set
427 (
428 iter.key(),
430 (
431 iter(),
432 *phases_.lookup(iter.key().first()),
433 *phases_.lookup(iter.key().second())
434 ).ptr()
435 );
436 }
437
438 for (const phaseModel& phase1 : phases_)
439 {
440 for (const phaseModel& phase2 : phases_)
441 {
442 if (&phase2 == &phase1)
443 {
444 continue;
445 }
446
447 const interfacePair key(phase1, phase2);
448
449 if (sigmas_.found(key) && !cAlphas_.found(key))
450 {
452 << "Compression coefficient not specified for phase pair ("
453 << phase1.name() << ' ' << phase2.name()
454 << ") for which a surface tension coefficient is specified"
455 << endl;
456 }
457 }
458 }
459}
460
461
462// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
463
465{
466 auto iter = phases_.cbegin();
467
468 tmp<volScalarField> trho = iter()*iter().rho();
469 volScalarField& rho = trho.ref();
470
471 for (++iter; iter != phases_.cend(); ++iter)
472 {
473 rho += iter()*iter().rho();
474 }
475
476 return trho;
477}
478
479
481Foam::multiphaseSystem::rho(const label patchi) const
482{
483 auto iter = phases_.cbegin();
484
485 tmp<scalarField> trho = iter().boundaryField()[patchi]*iter().rho().value();
486 scalarField& rho = trho.ref();
487
488 for (++iter; iter != phases_.cend(); ++iter)
489 {
490 rho += iter().boundaryField()[patchi]*iter().rho().value();
491 }
492
493 return trho;
494}
495
496
498{
499 auto iter = phases_.cbegin();
500
501 tmp<volScalarField> tmu = iter()*(iter().rho()*iter().nu());
502 volScalarField& mu = tmu.ref();
503
504 for (++iter; iter != phases_.cend(); ++iter)
505 {
506 mu += iter()*(iter().rho()*iter().nu());
507 }
508
509 return tmu/rho();
510}
511
512
514Foam::multiphaseSystem::nu(const label patchi) const
515{
516 auto iter = phases_.cbegin();
517
518 tmp<scalarField> tmu =
519 iter().boundaryField()[patchi]
520 *(iter().rho().value()*iter().nu().value());
521 scalarField& mu = tmu.ref();
522
523 for (++iter; iter != phases_.cend(); ++iter)
524 {
525 mu +=
526 iter().boundaryField()[patchi]
527 *(iter().rho().value()*iter().nu().value());
528 }
529
530 return tmu/rho(patchi);
531}
532
533
535(
536 const phaseModel& phase
537) const
538{
540 (
542 (
544 (
545 "Cvm",
546 mesh_.time().timeName(),
547 mesh_
548 ),
549 mesh_,
551 )
552 );
553
554 for (const phaseModel& phase2 : phases_)
555 {
556 if (&phase2 == &phase)
557 {
558 continue;
559 }
560
561 auto iterCvm = Cvms_.cfind(interfacePair(phase, phase2));
562
563 if (iterCvm.found())
564 {
565 tCvm.ref() += iterCvm()*phase2.rho()*phase2;
566 }
567 else
568 {
569 iterCvm = Cvms_.cfind(interfacePair(phase2, phase));
570
571 if (iterCvm.found())
572 {
573 tCvm.ref() += iterCvm()*phase.rho()*phase2;
574 }
575 }
576 }
577
578 return tCvm;
579}
580
581
583(
584 const phaseModel& phase
585) const
586{
588 (
590 (
592 (
593 "Svm",
594 mesh_.time().timeName(),
595 mesh_
596 ),
597 mesh_,
599 (
600 "Svm",
601 dimensionSet(1, -2, -2, 0, 0),
602 Zero
603 )
604 )
605 );
606
607 for (const phaseModel& phase2 : phases_)
608 {
609 if (&phase2 == &phase)
610 {
611 continue;
612 }
613
614 auto Cvm = Cvms_.cfind(interfacePair(phase, phase2));
615
616 if (Cvm.found())
617 {
618 tSvm.ref() += Cvm()*phase2.rho()*phase2*phase2.DDtU();
619 }
620 else
621 {
622 Cvm = Cvms_.cfind(interfacePair(phase2, phase));
623
624 if (Cvm.found())
625 {
626 tSvm.ref() += Cvm()*phase.rho()*phase2*phase2.DDtU();
627 }
628 }
629 }
630
632 tSvm.ref().boundaryFieldRef();
633
634 // Remove virtual mass at fixed-flux boundaries
635 forAll(phase.phi().boundaryField(), patchi)
636 {
637 if
638 (
639 isA<fixedValueFvsPatchScalarField>
640 (
641 phase.phi().boundaryField()[patchi]
642 )
643 )
644 {
645 SvmBf[patchi] = Zero;
646 }
647 }
648
649 return tSvm;
650}
651
652
655{
656 autoPtr<dragCoeffFields> dragCoeffsPtr(new dragCoeffFields);
657
658 forAllConstIters(dragModels_, iter)
659 {
660 const multiphaseEuler::dragModel& dm = *iter();
661
662 volScalarField* Kptr =
663 (
664 max
665 (
666 // fvc::average(dm.phase1()*dm.phase2()),
667 // fvc::average(dm.phase1())*fvc::average(dm.phase2()),
668 dm.phase1()*dm.phase2(),
670 )
671 *dm.K
672 (
673 max
674 (
675 mag(dm.phase1().U() - dm.phase2().U()),
676 dm.residualSlip()
677 )
678 )
679 ).ptr();
680
682
683 // Remove drag at fixed-flux boundaries
684 forAll(dm.phase1().phi().boundaryField(), patchi)
685 {
686 if
687 (
688 isA<fixedValueFvsPatchScalarField>
689 (
690 dm.phase1().phi().boundaryField()[patchi]
691 )
692 )
693 {
694 Kbf[patchi] = 0.0;
695 }
696 }
697
698 dragCoeffsPtr().set(iter.key(), Kptr);
699 }
700
701 return dragCoeffsPtr;
702}
703
704
706(
707 const phaseModel& phase,
709) const
710{
711 tmp<volScalarField> tdragCoeff
712 (
714 (
716 (
717 "dragCoeff",
718 mesh_.time().timeName(),
719 mesh_
720 ),
721 mesh_,
723 (
724 "dragCoeff",
725 dimensionSet(1, -3, -1, 0, 0),
726 0
727 )
728 )
729 );
730
731 dragModelTable::const_iterator dmIter = dragModels_.begin();
733 for
734 (
735 ;
736 dmIter.good() && dcIter.good();
737 ++dmIter, ++dcIter
738 )
739 {
740 if
741 (
742 &phase == &dmIter()->phase1()
743 || &phase == &dmIter()->phase2()
744 )
745 {
746 tdragCoeff.ref() += *dcIter();
747 }
748 }
749
750 return tdragCoeff;
751}
752
753
755(
756 const phaseModel& phase1
757) const
758{
759 tmp<surfaceScalarField> tSurfaceTension
760 (
762 (
764 (
765 "surfaceTension",
766 mesh_.time().timeName(),
767 mesh_
768 ),
769 mesh_,
771 (
772 "surfaceTension",
773 dimensionSet(1, -2, -2, 0, 0),
774 0
775 )
776 )
777 );
778 tSurfaceTension.ref().setOriented();
779
780 for (const phaseModel& phase2 : phases_)
781 {
782 if (&phase2 == &phase1)
783 {
784 continue;
785 }
786
787 const auto sigma = sigmas_.cfind(interfacePair(phase1, phase2));
788
789 if (sigma.found())
790 {
791 tSurfaceTension.ref() +=
792 dimensionedScalar("sigma", dimSigma_, *sigma)
794 (
797 );
798 }
799 }
800
801 return tSurfaceTension;
802}
803
804
807{
808 tmp<volScalarField> tnearInt
809 (
811 (
813 (
814 "nearInterface",
815 mesh_.time().timeName(),
816 mesh_
817 ),
818 mesh_,
819 dimensionedScalar("nearInterface", dimless, 0.0)
820 )
821 );
822
823 for (const phaseModel& phase : phases_)
824 {
825 tnearInt.ref() =
826 max(tnearInt(), pos0(phase - 0.01)*pos0(0.99 - phase));
827 }
828
829 return tnearInt;
830}
831
832
834{
835 for (phaseModel& phase : phases_)
836 {
837 phase.correct();
838 }
839
840 const Time& runTime = mesh_.time();
841
842 const dictionary& alphaControls = mesh_.solverDict("alpha");
843 label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
844
845 if (nAlphaSubCycles > 1)
846 {
847 dimensionedScalar totalDeltaT = runTime.deltaT();
848
849 PtrList<volScalarField> alpha0s(phases_.size());
850 PtrList<surfaceScalarField> alphaPhiSums(phases_.size());
851
852 label phasei = 0;
853 for (phaseModel& phase : phases_)
854 {
856
857 alpha0s.set
858 (
859 phasei,
861 );
862
863 alphaPhiSums.set
864 (
865 phasei,
867 (
869 (
870 "phiSum" + alpha.name(),
872 mesh_
873 ),
874 mesh_,
875 dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0)
876 )
877 );
878
879 ++phasei;
880 }
881
882 for
883 (
884 subCycleTime alphaSubCycle
885 (
886 const_cast<Time&>(runTime),
888 );
889 !(++alphaSubCycle).end();
890 )
891 {
892 solveAlphas();
893
894 label phasei = 0;
895 for (const phaseModel& phase : phases_)
896 {
897 alphaPhiSums[phasei] += phase.alphaPhi()/nAlphaSubCycles;
898
899 ++phasei;
900 }
901 }
902
903 phasei = 0;
904 for (phaseModel& phase : phases_)
905 {
907
908 phase.alphaPhi() = alphaPhiSums[phasei];
909
910 // Correct the time index of the field
911 // to correspond to the global time
913
914 // Reset the old-time field value
915 alpha.oldTime() = alpha0s[phasei];
917
918 ++phasei;
919 }
920 }
921 else
922 {
923 solveAlphas();
924 }
925}
926
927
929{
930 if (regIOobject::read())
931 {
932 bool readOK = true;
933
934 PtrList<entry> phaseData(lookup("phases"));
935 label phasei = 0;
936
937 for (phaseModel& phase : phases_)
938 {
939 readOK &= phase.read(phaseData[phasei++].dict());
940 }
941
942 lookup("sigmas") >> sigmas_;
943 lookup("interfaceCompression") >> cAlphas_;
944 lookup("virtualMass") >> Cvms_;
945
946 return readOK;
947 }
948 else
949 {
950 return false;
951 }
952}
953
954
955// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
MULES: Multidimensional universal limiter for explicit solution.
surfaceScalarField & phi
const volScalarField & alpha1
phaseModel & phase1
const volScalarField & alpha2
phaseModel & phase2
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
label timeIndex() const
Return the time index of the field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
Definition: HashPtrTable.H:68
bool set(const Key &key, T *ptr)
Assign a new entry, overwriting existing entries.
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:330
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
label timeIndex() const noexcept
Return current time index.
Definition: TimeStateI.H:37
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:55
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
void set(T *p) noexcept
Deprecated(2018-02) Identical to reset().
Definition: autoPtr.H:288
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
const surfaceScalarField & nHatf() const
areaScalarField & surfaceTension()
Return surface tension field.
const phaseModel & phase2() const
Definition: dragModel.H:119
const dimensionedScalar & residualSlip() const
Definition: dragModel.H:129
virtual tmp< volScalarField > K(const volScalarField &Ur) const =0
The drag function K used in the momentum eq.
const phaseModel & phase1() const
Definition: dragModel.H:114
const dimensionedScalar & residualPhaseFraction() const
Definition: dragModel.H:124
Name pair for the interface.
Incompressible multi-phase mixture with built in solution for the phase fractions with interface comp...
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
autoPtr< dragCoeffFields > dragCoeffs() const
Return the drag coefficients for all of the interfaces.
tmp< volScalarField > nu() const
Return the mixture laminar viscosity.
tmp< volScalarField > rho() const
Return the mixture density.
tmp< volScalarField > dragCoeff(const phaseModel &phase, const dragCoeffFields &dragCoeffs) const
Return the sum of the drag coefficients for the given phase.
tmp< volVectorField > Svm(const phaseModel &phase) const
Return the virtual-mass source for the given phase.
void solve()
Solve for the mixture phase-fractions.
bool read()
Read base transportProperties dictionary.
const Time & time() const noexcept
Return time registry.
constant condensation/saturation model.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
const volVectorField & U() const
Definition: phaseModel.H:176
const surfaceScalarField & phi() const
Definition: phaseModel.H:196
const dimensionedScalar & rho() const
Definition: phaseModel.H:171
const volVectorField & DDtU() const
Definition: phaseModel.H:186
const word & name() const
Definition: phaseModel.H:144
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
void correct()
Correct the phase properties.
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition: phase.H:140
bool read(const dictionary &phaseDict)
Read base transportProperties dictionary.
Lookup type of boundary radiation properties.
Definition: lookup.H:66
virtual bool write(const bool valid=true) const
Write using setting from DB.
virtual bool read()
Read object.
A class for managing sub-cycling times.
Definition: subCycleTime.H:51
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
virtual tmp< volScalarField > Cvm() const
Virtual mass coefficient.
U
Definition: pEqn.H:72
const volScalarField & mu
faceListList boundary
dynamicFvMesh & mesh
engineTime & runTime
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Area-weighted average a surfaceField creating a volField.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the gradient of the given field.
Calculate the snGrad of the given volField.
word timeName
Definition: getTimeIndex.H:3
#define WarningInFunction
Report a warning using Foam::Warning.
autoPtr< multiphaseSystem::dragCoeffFields > dragCoeffs(fluid.dragCoeffs())
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 limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:34
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:34
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensionedScalar det(const dimensionedSphericalTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
volScalarField & alpha
dictionary dict
tmp< volScalarField > trho
volScalarField & b
Definition: createFields.H:27
volScalarField & e
Definition: createFields.H:11
const dictionary & alphaControls
Definition: alphaControls.H:1
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
Unit conversion functions.