Loading...
Searching...
No Matches
multiphaseMixture.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) 2021 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 "multiphaseMixture.H"
30#include "Time.H"
31#include "subCycle.H"
32#include "MULES.H"
33#include "surfaceInterpolate.H"
34#include "fvcGrad.H"
35#include "fvcSnGrad.H"
36#include "fvcDiv.H"
37#include "fvcFlux.H"
38#include "unitConversion.H"
39#include "alphaContactAngleFvPatchScalarField.H"
40
41// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42
43void Foam::multiphaseMixture::calcAlphas()
44{
45 scalar level = 0.0;
46 alphas_ == 0.0;
47
48 for (const phase& ph : phases_)
49 {
50 alphas_ += level * ph;
51 level += 1.0;
52 }
53}
54
55
56// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57
59(
60 const volVectorField& U,
62)
63:
65 (
67 (
68 "transportProperties",
69 U.time().constant(),
70 U.db(),
71 IOobject::READ_MODIFIED,
72 IOobject::NO_WRITE,
73 IOobject::REGISTER
74 )
75 ),
76
77 phases_(lookup("phases"), phase::iNew(U, phi)),
78
79 mesh_(U.mesh()),
80 U_(U),
81 phi_(phi),
82
83 rhoPhi_
84 (
86 (
87 "rhoPhi",
88 mesh_.time().timeName(),
89 mesh_
90 ),
91 mesh_,
93 ),
94
95 alphas_
96 (
98 (
99 "alphas",
100 mesh_.time().timeName(),
101 mesh_,
102 IOobject::NO_READ,
103 IOobject::AUTO_WRITE,
104 IOobject::REGISTER
105 ),
106 mesh_,
108 ),
109
110 nu_
111 (
113 (
114 "nu",
115 mesh_.time().timeName(),
116 mesh_
117 ),
118 mu()/rho()
119 ),
120
121 sigmas_(lookup("sigmas")),
122 dimSigma_(1, 0, -2, 0, 0),
123 deltaN_
124 (
125 "deltaN",
126 1e-8/cbrt(average(mesh_.V()))
127 )
128{
129 rhoPhi_.setOriented();
130
131 calcAlphas();
132 alphas_.write();
133}
134
135
136// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
137
140{
141 auto iter = phases_.cbegin();
142
143 tmp<volScalarField> trho = iter()*iter().rho();
144 volScalarField& rho = trho.ref();
145
146 for (++iter; iter != phases_.cend(); ++iter)
147 {
148 rho += iter()*iter().rho();
149 }
150
151 return trho;
152}
153
154
155Foam::tmp<Foam::scalarField>
156Foam::multiphaseMixture::rho(const label patchi) const
157{
158 auto iter = phases_.cbegin();
159
160 tmp<scalarField> trho = iter().boundaryField()[patchi]*iter().rho().value();
161 scalarField& rho = trho.ref();
162
163 for (++iter; iter != phases_.cend(); ++iter)
164 {
165 rho += iter().boundaryField()[patchi]*iter().rho().value();
166 }
167
168 return trho;
169}
170
171
172Foam::tmp<Foam::volScalarField>
174{
175 auto iter = phases_.cbegin();
176
177 tmp<volScalarField> tmu = iter()*iter().rho()*iter().nu();
178 volScalarField& mu = tmu.ref();
179
180 for (++iter; iter != phases_.cend(); ++iter)
181 {
182 mu += iter()*iter().rho()*iter().nu();
183 }
184
185 return tmu;
186}
187
188
189Foam::tmp<Foam::scalarField>
190Foam::multiphaseMixture::mu(const label patchi) const
191{
192 auto iter = phases_.cbegin();
193
194 tmp<scalarField> tmu =
195 (
196 iter().boundaryField()[patchi]
197 *iter().rho().value()
198 *iter().nu(patchi)
199 );
200
201 scalarField& mu = tmu.ref();
202
203 for (++iter; iter != phases_.cend(); ++iter)
204 {
205 mu +=
206 (
207 iter().boundaryField()[patchi]
208 *iter().rho().value()
209 *iter().nu(patchi)
210 );
211 }
212
213 return tmu;
214}
215
216
217Foam::tmp<Foam::surfaceScalarField>
219{
220 auto iter = phases_.cbegin();
221
223 fvc::interpolate(iter())*iter().rho()*fvc::interpolate(iter().nu());
224 surfaceScalarField& muf = tmuf.ref();
225
226 for (++iter; iter != phases_.cend(); ++iter)
227 {
228 muf +=
229 fvc::interpolate(iter())*iter().rho()*fvc::interpolate(iter().nu());
230 }
231
232 return tmuf;
233}
234
235
236Foam::tmp<Foam::volScalarField>
238{
239 return nu_;
240}
241
242
243Foam::tmp<Foam::scalarField>
244Foam::multiphaseMixture::nu(const label patchi) const
245{
246 return nu_.boundaryField()[patchi];
247}
248
249
250Foam::tmp<Foam::surfaceScalarField>
252{
253 return muf()/fvc::interpolate(rho());
254}
255
256
257Foam::tmp<Foam::surfaceScalarField>
259{
260 auto tstf = surfaceScalarField::New
261 (
262 "surfaceTensionForce",
264 mesh_,
265 dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), Zero)
266 );
267
268 surfaceScalarField& stf = tstf.ref();
269 stf.setOriented();
270
271 forAllConstIters(phases_, iter1)
272 {
273 const phase& alpha1 = iter1();
274
275 auto iter2 = iter1;
276
277 for (++iter2; iter2 != phases_.cend(); ++iter2)
278 {
279 const phase& alpha2 = iter2();
280
281 auto sigma = sigmas_.cfind(interfacePair(alpha1, alpha2));
282
283 if (!sigma.good())
284 {
286 << "Cannot find interface " << interfacePair(alpha1, alpha2)
287 << " in list of sigma values"
288 << exit(FatalError);
289 }
290
291 stf += dimensionedScalar("sigma", dimSigma_, *sigma)
293 (
296 );
297 }
298 }
299
300 return tstf;
301}
302
303
305{
306 correct();
307
308 const Time& runTime = mesh_.time();
309
310 volScalarField& alpha = phases_.first();
311
312 const dictionary& alphaControls = mesh_.solverDict("alpha");
313 label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
314 scalar cAlpha(alphaControls.get<scalar>("cAlpha"));
315
316 if (nAlphaSubCycles > 1)
317 {
318 surfaceScalarField rhoPhiSum
319 (
320 mesh_.newIOobject("rhoPhiSum"),
321 mesh_,
322 dimensionedScalar(rhoPhi_.dimensions(), Zero)
323 );
324
325 dimensionedScalar totalDeltaT = runTime.deltaT();
326
327 for
328 (
330 !(++alphaSubCycle).end();
331 )
332 {
333 solveAlphas(cAlpha);
334 rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi_;
335 }
336
337 rhoPhi_ = rhoPhiSum;
338 }
339 else
340 {
341 solveAlphas(cAlpha);
342 }
343
344 // Update the mixture kinematic viscosity
345 nu_ = mu()/rho();
346}
347
348
350{
351 for (phase& ph : phases_)
352 {
353 ph.correct();
354 }
355}
356
357
358Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseMixture::nHatfv
359(
360 const volScalarField& alpha1,
362) const
363{
364 /*
365 // Cell gradient of alpha
366 volVectorField gradAlpha =
367 alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
368
369 // Interpolated face-gradient of alpha
370 surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
371 */
372
373 surfaceVectorField gradAlphaf
374 (
377 );
378
379 // Face unit interface normal
380 return gradAlphaf/(mag(gradAlphaf) + deltaN_);
381}
382
383
384Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseMixture::nHatf
385(
386 const volScalarField& alpha1,
388) const
389{
390 // Face unit interface normal flux
391 return nHatfv(alpha1, alpha2) & mesh_.Sf();
392}
393
394
395// Correction for the boundary condition on the unit normal nHat on
396// walls to produce the correct contact angle.
397
398// The dynamic contact angle is calculated from the component of the
399// velocity on the direction of the interface, parallel to the wall.
400
401void Foam::multiphaseMixture::correctContactAngle
402(
403 const phase& alpha1,
404 const phase& alpha2,
406) const
407{
410
411 const fvBoundaryMesh& boundary = mesh_.boundary();
412
413 forAll(boundary, patchi)
414 {
416 {
419
420 correctBoundaryContactAngle(acap, patchi, alpha1, alpha2, nHatb);
421 }
422 else if (isA<alphaContactAngleFvPatchScalarField>(gb2f[patchi]))
423 {
426
427 correctBoundaryContactAngle(acap, patchi, alpha2, alpha1, nHatb);
428 }
429 }
430}
431
432
433void Foam::multiphaseMixture::correctBoundaryContactAngle
434(
436 label patchi,
437 const phase& alpha1,
438 const phase& alpha2,
440) const
441{
442 const fvBoundaryMesh& boundary = mesh_.boundary();
443
444 vectorField& nHatPatch = nHatb[patchi];
445
446 vectorField AfHatPatch
447 (
448 mesh_.Sf().boundaryField()[patchi]
449 /mesh_.magSf().boundaryField()[patchi]
450 );
451
452 const auto tp = acap.thetaProps().cfind(interfacePair(alpha1, alpha2));
453
454 if (!tp.good())
455 {
457 << "Cannot find interface " << interfacePair(alpha1, alpha2)
458 << "\n in table of theta properties for patch "
459 << acap.patch().name()
460 << exit(FatalError);
461 }
462
463 const bool matched = (tp.key().first() == alpha1.name());
464
465 const scalar theta0 = degToRad(tp().theta0(matched));
466 scalarField theta(boundary[patchi].size(), theta0);
467
468 const scalar uTheta = tp().uTheta();
469
470 // Calculate the dynamic contact angle if required
471 if (uTheta > SMALL)
472 {
473 const scalar thetaA = degToRad(tp().thetaA(matched));
474 const scalar thetaR = degToRad(tp().thetaR(matched));
475
476 // Calculated the component of the velocity parallel to the wall
477 vectorField Uwall
478 (
479 U_.boundaryField()[patchi].patchInternalField()
480 - U_.boundaryField()[patchi]
481 );
482 Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
483
484 // Find the direction of the interface parallel to the wall
485 vectorField nWall
486 (
487 nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
488 );
489
490 // Normalise nWall
491 nWall /= (mag(nWall) + SMALL);
492
493 // Calculate Uwall resolved normal to the interface parallel to
494 // the interface
495 scalarField uwall(nWall & Uwall);
496
497 theta += (thetaA - thetaR)*tanh(uwall/uTheta);
498 }
499
500
501 // Reset nHatPatch to correspond to the contact angle
502
503 scalarField a12(nHatPatch & AfHatPatch);
504
505 scalarField b1(cos(theta));
506
507 scalarField b2(nHatPatch.size());
508
509 forAll(b2, facei)
510 {
511 b2[facei] = cos(acos(a12[facei]) - theta[facei]);
512 }
513
514 scalarField det(1.0 - a12*a12);
515
516 scalarField a((b1 - a12*b2)/det);
517 scalarField b((b2 - a12*b1)/det);
518
519 nHatPatch = a*AfHatPatch + b*nHatPatch;
520
521 nHatPatch /= (mag(nHatPatch) + deltaN_.value());
522}
523
524
525Foam::tmp<Foam::volScalarField> Foam::multiphaseMixture::K
526(
527 const phase& alpha1,
528 const phase& alpha2
529) const
530{
531 tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
532
533 correctContactAngle(alpha1, alpha2, tnHatfv.ref().boundaryFieldRef());
534
535 // Simple expression for curvature
536 return -fvc::div(tnHatfv & mesh_.Sf());
537}
538
539
540Foam::tmp<Foam::volScalarField>
542{
543 auto tnearInt = volScalarField::New
544 (
545 "nearInterface",
547 mesh_,
549 );
550
551 for (const phase& ph : phases_)
552 {
553 tnearInt.ref() = max(tnearInt(), pos0(ph - 0.01)*pos0(0.99 - ph));
554 }
555
556 return tnearInt;
557}
558
559
560void Foam::multiphaseMixture::solveAlphas
561(
562 const scalar cAlpha
563)
564{
565 static label nSolves(-1);
566 ++nSolves;
567
568 const word alphaScheme("div(phi,alpha)");
569 const word alpharScheme("div(phirb,alpha)");
570
571 surfaceScalarField phic(mag(phi_/mesh_.magSf()));
572 phic = min(cAlpha*phic, max(phic));
573
574 PtrList<surfaceScalarField> alphaPhiCorrs(phases_.size());
575 int phasei = 0;
576
577 for (phase& alpha : phases_)
578 {
579 alphaPhiCorrs.set
580 (
581 phasei,
583 (
584 "phi" + alpha.name() + "Corr",
586 (
587 phi_,
588 alpha,
589 alphaScheme
590 )
591 )
592 );
593
594 surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phasei];
595
596 for (phase& alpha2 : phases_)
597 {
598 if (&alpha2 == &alpha) continue;
599
601
602 alphaPhiCorr += fvc::flux
603 (
605 alpha,
607 );
608 }
609
611 (
612 1.0/mesh_.time().deltaTValue(),
614 alpha,
615 phi_,
616 alphaPhiCorr,
617 zeroField(),
618 zeroField(),
619 oneField(),
620 zeroField(),
621 true
622 );
623
624 ++phasei;
625 }
626
627 MULES::limitSum(alphaPhiCorrs);
628
629 rhoPhi_ = Zero;
630
631 volScalarField sumAlpha
632 (
633 mesh_.newIOobject("sumAlpha"),
634 mesh_,
636 );
637
638 phasei = 0;
639
640 for (phase& alpha : phases_)
641 {
642 surfaceScalarField& alphaPhi = alphaPhiCorrs[phasei];
643 alphaPhi += upwind<scalar>(mesh_, phi_).flux(alpha);
644
646 (
648 alpha,
650 );
651
652 rhoPhi_ += alphaPhi*alpha.rho();
653
654 Info<< alpha.name() << " volume fraction, min, max = "
655 << alpha.weightedAverage(mesh_.V()).value()
656 << ' ' << min(alpha).value()
657 << ' ' << max(alpha).value()
658 << endl;
659
660 sumAlpha += alpha;
661
662 ++phasei;
663 }
664
665 Info<< "Phase-sum volume fraction, min, max = "
666 << sumAlpha.weightedAverage(mesh_.V()).value()
667 << ' ' << min(sumAlpha).value()
668 << ' ' << max(sumAlpha).value()
669 << endl;
670
671 // Correct the sum of the phase-fractions to avoid 'drift'
672 volScalarField sumCorr(1.0 - sumAlpha);
673 for (phase& alpha : phases_)
674 {
675 alpha += alpha*sumCorr;
676 }
677
678 calcAlphas();
679}
680
681
683{
685 {
686 bool readOK = true;
687
688 PtrList<entry> phaseData(lookup("phases"));
689 label phasei = 0;
690
691 for (phase& ph : phases_)
692 {
693 readOK &= ph.read(phaseData[phasei++].dict());
694 }
695
696 readEntry("sigmas", sigmas_);
697
698 return readOK;
699 }
700
701 return false;
702}
703
704
705// ************************************************************************* //
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
const volScalarField & alpha2
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weights, const label comm=UPstream::worldComm) const
Return the global weighted average.
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
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
T & first()
Access first element of the list, position [0].
Definition UList.H:957
Contact-angle boundary condition for multi-phase interface-capturing simulations. Used in conjunction...
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.
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.
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
tmp< surfaceScalarField > surfaceTensionForce() const
void correct()
Correct the mixture properties.
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
tmp< volScalarField > nu() const
Return the kinematic laminar viscosity.
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
tmp< volScalarField > rho() const
Return the mixture density.
multiphaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
void solve()
Solve for the mixture phase-fractions.
bool read()
Read base transportProperties dictionary.
tmp< surfaceScalarField > nuf() const
Return the face-interpolated dynamic laminar viscosity.
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 phase.H:53
Lookup type of boundary radiation properties.
Definition lookup.H:60
bool end() const
Return true if the number of sub-cycles has been reached.
Perform a subCycleTime on a field.
Definition subCycle.H:137
A class for managing temporary objects.
Definition tmp.H:75
virtual bool read()=0
Read transportProperties dictionary.
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
thermo correct()
faceListList boundary
dynamicFvMesh & mesh
engineTime & runTime
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
word alpharScheme("div(phirb,alpha)")
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
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
label phasei
Definition pEqn.H:27
const expr V(m.psi().mesh().V())
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
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)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
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.
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)
dimensionedScalar acos(const dimensionedScalar &ds)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
volScalarField & nu
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.