Loading...
Searching...
No Matches
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 {
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
92
93 const auto cAlpha = cAlphas_.cfind(interfacePair(phase, phase2));
94
95 if (cAlpha.good())
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().deltaTValue(),
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 (
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 (
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 {
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
223Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::nHatf
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 (
257 (
258 gbf[patchi]
259 )
260 )
261 {
262 const auto& acap =
263 refCast
264 <
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.good())
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
355Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::K
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
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 }
458 }
459}
460
461
462// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
463
464Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::rho() const
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();
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());
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{
539 auto tCvm = volScalarField::New
540 (
541 "Cvm",
543 mesh_,
545 );
546
547 for (const phaseModel& phase2 : phases_)
548 {
549 if (&phase2 == &phase)
550 {
551 continue;
552 }
553
554 auto iterCvm = Cvms_.cfind(interfacePair(phase, phase2));
555
556 if (iterCvm.good())
557 {
558 tCvm.ref() += iterCvm()*phase2.rho()*phase2;
559 }
560 else
561 {
562 iterCvm = Cvms_.cfind(interfacePair(phase2, phase));
563
564 if (iterCvm.good())
565 {
566 tCvm.ref() += iterCvm()*phase.rho()*phase2;
567 }
569 }
570
571 return tCvm;
572}
573
574
575Foam::tmp<Foam::volVectorField> Foam::multiphaseSystem::Svm
576(
577 const phaseModel& phase
578) const
579{
580 auto tSvm = volVectorField::New
581 (
582 "Svm",
584 mesh_,
586 (
587 dimensionSet(1, -2, -2, 0, 0),
588 Zero
589 )
590 );
591
592 for (const phaseModel& phase2 : phases_)
593 {
594 if (&phase2 == &phase)
595 {
596 continue;
597 }
598
599 auto Cvm = Cvms_.cfind(interfacePair(phase, phase2));
600
601 if (Cvm.good())
602 {
603 tSvm.ref() += Cvm()*phase2.rho()*phase2*phase2.DDtU();
604 }
605 else
606 {
607 Cvm = Cvms_.cfind(interfacePair(phase2, phase));
608
609 if (Cvm.good())
610 {
611 tSvm.ref() += Cvm()*phase.rho()*phase2*phase2.DDtU();
612 }
613 }
614 }
615
617 tSvm.ref().boundaryFieldRef();
618
619 // Remove virtual mass at fixed-flux boundaries
620 forAll(phase.phi().boundaryField(), patchi)
621 {
622 if
623 (
625 (
626 phase.phi().boundaryField()[patchi]
627 )
628 )
629 {
630 SvmBf[patchi] = Zero;
631 }
633
634 return tSvm;
635}
636
637
640{
641 auto dragCoeffsPtr = autoPtr<dragCoeffFields>::New();
642
643 forAllConstIters(dragModels_, iter)
644 {
645 const multiphaseEuler::dragModel& dm = *iter();
646
647 volScalarField* Kptr =
648 (
649 max
650 (
651 // fvc::average(dm.phase1()*dm.phase2()),
652 // fvc::average(dm.phase1())*fvc::average(dm.phase2()),
653 dm.phase1()*dm.phase2(),
654 dm.residualPhaseFraction()
655 )
656 *dm.K
657 (
658 max
659 (
660 mag(dm.phase1().U() - dm.phase2().U()),
661 dm.residualSlip()
662 )
663 )
664 ).ptr();
665
667
668 // Remove drag at fixed-flux boundaries
669 forAll(dm.phase1().phi().boundaryField(), patchi)
670 {
671 if
672 (
674 (
675 dm.phase1().phi().boundaryField()[patchi]
676 )
677 )
678 {
679 Kbf[patchi] = 0.0;
680 }
681 }
682
683 dragCoeffsPtr().set(iter.key(), Kptr);
684 }
685
686 return dragCoeffsPtr;
687}
688
689
691(
692 const phaseModel& phase,
693 const dragCoeffFields& dragCoeffs
694) const
695{
696 auto tdragCoeff = volScalarField::New
697 (
698 "dragCoeff",
700 mesh_,
701 dimensionedScalar(dimensionSet(1, -3, -1, 0, 0), Zero)
702 );
703
704 dragModelTable::const_iterator dmIter = dragModels_.begin();
705 dragCoeffFields::const_iterator dcIter = dragCoeffs.begin();
706 for
707 (
708 ;
709 dmIter.good() && dcIter.good();
710 ++dmIter, ++dcIter
711 )
712 {
713 if
714 (
715 &phase == &dmIter()->phase1()
716 || &phase == &dmIter()->phase2()
717 )
718 {
719 tdragCoeff.ref() += *dcIter();
721 }
722
723 return tdragCoeff;
724}
725
726
727Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::surfaceTension
728(
729 const phaseModel& phase1
730) const
731{
732 auto tSurfaceTension = surfaceScalarField::New
733 (
734 "surfaceTension",
736 mesh_,
737 dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), Zero)
738 );
739 tSurfaceTension.ref().setOriented();
740
741 for (const phaseModel& phase2 : phases_)
742 {
743 if (&phase2 == &phase1)
744 {
745 continue;
746 }
747
748 const auto sigma = sigmas_.cfind(interfacePair(phase1, phase2));
749
750 if (sigma.good())
751 {
752 tSurfaceTension.ref() +=
753 dimensionedScalar("sigma", dimSigma_, *sigma)
755 (
758 );
759 }
761
762 return tSurfaceTension;
763}
764
765
768{
769 auto tnearInt = volScalarField::New
770 (
771 "nearInterface",
773 mesh_,
775 );
776
777 for (const phaseModel& phase : phases_)
778 {
779 tnearInt.ref() =
780 max(tnearInt(), pos0(phase - 0.01)*pos0(0.99 - phase));
781 }
782
783 return tnearInt;
784}
785
786
788{
789 for (phaseModel& phase : phases_)
790 {
791 phase.correct();
792 }
793
794 const Time& runTime = mesh_.time();
795
796 const dictionary& alphaControls = mesh_.solverDict("alpha");
797 label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
798
799 if (nAlphaSubCycles > 1)
800 {
801 dimensionedScalar totalDeltaT = runTime.deltaT();
802
803 PtrList<volScalarField> alpha0s(phases_.size());
804 PtrList<surfaceScalarField> alphaPhiSums(phases_.size());
805
806 label phasei = 0;
807 for (phaseModel& phase : phases_)
808 {
810
811 alpha0s.set
812 (
813 phasei,
815 );
816
817 alphaPhiSums.set
818 (
819 phasei,
821 (
823 (
824 "phiSum" + alpha.name(),
826 mesh_
827 ),
828 mesh_,
829 dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
830 )
831 );
832
833 ++phasei;
834 }
835
836 for
837 (
838 subCycleTime alphaSubCycle
839 (
840 const_cast<Time&>(runTime),
842 );
843 !(++alphaSubCycle).end();
844 )
845 {
846 solveAlphas();
847
848 label phasei = 0;
849 for (const phaseModel& phase : phases_)
850 {
851 alphaPhiSums[phasei] += phase.alphaPhi()/nAlphaSubCycles;
852
853 ++phasei;
854 }
855 }
856
857 phasei = 0;
858 for (phaseModel& phase : phases_)
859 {
861
862 phase.alphaPhi() = alphaPhiSums[phasei];
863
864 // Correct the time index of the field
865 // to correspond to the global time
867
868 // Reset the old-time field value
869 alpha.oldTime() = alpha0s[phasei];
871
872 ++phasei;
873 }
874 }
875 else
876 {
877 solveAlphas();
878 }
879}
880
881
883{
884 if (regIOobject::read())
885 {
886 bool readOK = true;
887
888 PtrList<entry> phaseData(lookup("phases"));
889 label phasei = 0;
890
891 for (phaseModel& phase : phases_)
892 {
893 readOK &= phase.read(phaseData[phasei++].dict());
894 }
895
896 lookup("sigmas") >> sigmas_;
897 lookup("interfaceCompression") >> cAlphas_;
898 lookup("virtualMass") >> Cvms_;
899
900 return readOK;
901 }
902 else
903 {
904 return false;
905 }
906}
907
908
909// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
MULES: Multidimensional universal limiter for explicit solution.
const volScalarField & alpha1
phaseModel & phase1
const volScalarField & alpha2
phaseModel & phase2
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weights, const label comm=UPstream::worldComm) const
Return the global weighted average.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
label timeIndex() const noexcept
The time index of the field.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition List.H:469
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
label timeIndex() const noexcept
Return the current time index.
Definition TimeStateI.H:43
dimensionedScalar deltaT() const
Return time step.
Definition TimeStateI.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream.
Definition dictionary.C:367
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
const phaseModel & phase2() const
Definition dragModel.H:119
const dimensionedScalar & residualSlip() const
Definition dragModel.H:129
static autoPtr< dragModel > New(const dictionary &dict, const phaseModel &phase1, const phaseModel &phase2)
Definition dragModel.C:57
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.
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
HashPtrTable< volScalarField, interfacePair, interfacePair::symmHash > dragCoeffFields
autoPtr< dragCoeffFields > dragCoeffs() const
Return the drag coefficients for all of the interfaces.
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
tmp< volScalarField > nu() const
Return the mixture laminar viscosity.
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
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< volScalarField > Cvm(const phaseModel &phase) const
Return the virtual-mass coefficient 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.
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition oneField.H:50
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
const volVectorField & U() const
Definition phaseModel.H:198
const surfaceScalarField & phi() const
Definition phaseModel.H:218
const dimensionedScalar & rho() const
Definition phaseModel.H:193
const volVectorField & DDtU() const
Definition phaseModel.H:208
const word & name() const
Definition phaseModel.H:166
const surfaceScalarField & phi() const
Return the mixture flux.
tmp< volVectorField > U() const
Return the mixture velocity.
const fvMesh & mesh() const
Return the mesh.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
void correct()
Correct the phase properties.
const word & name() const
Definition phase.H:114
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition phase.H:151
bool read(const dictionary &phaseDict)
Read base transportProperties dictionary.
Lookup type of boundary radiation properties.
Definition lookup.H:60
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
virtual bool read()
Read object.
A class for managing sub-cycling times.
bool end() const
Return true if the number of sub-cycles has been reached.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Upwind differencing scheme class.
Definition upwind.H:55
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition zeroField.H:50
U
Definition pEqn.H:72
faceListList boundary
dynamicFvMesh & mesh
engineTime & runTime
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
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:27
Different types of constants.
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:47
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvcSnGrad.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition fvcFlux.C:27
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
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:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
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:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dimensionedScalar acos(const dimensionedScalar &ds)
volScalarField & alpha
dictionary dict
tmp< volScalarField > trho
volScalarField & b
volScalarField & e
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
const dictionary & alphaControls
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
Unit conversion functions.