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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "multiphaseMixture.H"
27 #include "alphaContactAngleFvPatchScalarField.H"
28 #include "Time.H"
29 #include "subCycle.H"
30 #include "MULES.H"
31 #include "surfaceInterpolate.H"
32 #include "fvcGrad.H"
33 #include "fvcSnGrad.H"
34 #include "fvcDiv.H"
35 #include "fvcFlux.H"
36 
37 // * * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * //
38 
39 const Foam::scalar Foam::multiphaseMixture::convertToRad =
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 {
47  scalar level = 0.0;
48  alphas_ == 0.0;
49 
50  forAllIter(PtrDictionary<phase>, phases_, iter)
51  {
52  alphas_ += level*iter();
53  level += 1.0;
54  }
55 
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const volVectorField& U,
65  const surfaceScalarField& phi
66 )
67 :
68  IOdictionary
69  (
70  IOobject
71  (
72  "transportProperties",
73  U.time().constant(),
74  U.db(),
75  IOobject::MUST_READ_IF_MODIFIED,
76  IOobject::NO_WRITE
77  )
78  ),
79 
80  phases_(lookup("phases"), phase::iNew(U, phi)),
81 
82  mesh_(U.mesh()),
83  U_(U),
84  phi_(phi),
85 
86  rhoPhi_
87  (
88  IOobject
89  (
90  "rhoPhi",
91  mesh_.time().timeName(),
92  mesh_,
93  IOobject::NO_READ,
94  IOobject::NO_WRITE
95  ),
96  mesh_,
97  dimensionedScalar("rhoPhi", dimMass/dimTime, 0.0)
98  ),
99 
100  alphas_
101  (
102  IOobject
103  (
104  "alphas",
105  mesh_.time().timeName(),
106  mesh_,
107  IOobject::NO_READ,
108  IOobject::AUTO_WRITE
109  ),
110  mesh_,
111  dimensionedScalar("alphas", dimless, 0.0),
112  zeroGradientFvPatchScalarField::typeName
113  ),
114 
115  nu_
116  (
117  IOobject
118  (
119  "nu",
120  mesh_.time().timeName(),
121  mesh_
122  ),
123  mu()/rho()
124  ),
125 
126  sigmas_(lookup("sigmas")),
127  dimSigma_(1, 0, -2, 0, 0),
128  deltaN_
129  (
130  "deltaN",
131  1e-8/pow(average(mesh_.V()), 1.0/3.0)
132  )
133 {
134  calcAlphas();
135  alphas_.write();
136 }
137 
138 
139 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
140 
143 {
144  PtrDictionary<phase>::const_iterator iter = phases_.begin();
145 
146  tmp<volScalarField> trho = iter()*iter().rho();
147 
148  for (++iter; iter != phases_.end(); ++iter)
149  {
150  trho() += iter()*iter().rho();
151  }
152 
153  return trho;
154 }
155 
156 
159 {
160  PtrDictionary<phase>::const_iterator iter = phases_.begin();
161 
162  tmp<scalarField> trho = iter().boundaryField()[patchi]*iter().rho().value();
163 
164  for (++iter; iter != phases_.end(); ++iter)
165  {
166  trho() += iter().boundaryField()[patchi]*iter().rho().value();
167  }
168 
169  return trho;
170 }
171 
172 
175 {
176  PtrDictionary<phase>::const_iterator iter = phases_.begin();
177 
178  tmp<volScalarField> tmu = iter()*iter().rho()*iter().nu();
179 
180  for (++iter; iter != phases_.end(); ++iter)
181  {
182  tmu() += iter()*iter().rho()*iter().nu();
183  }
184 
185  return tmu;
186 }
187 
188 
191 {
192  PtrDictionary<phase>::const_iterator iter = phases_.begin();
193 
194  tmp<scalarField> tmu =
195  iter().boundaryField()[patchi]
196  *iter().rho().value()
197  *iter().nu(patchi);
198 
199  for (++iter; iter != phases_.end(); ++iter)
200  {
201  tmu() +=
202  iter().boundaryField()[patchi]
203  *iter().rho().value()
204  *iter().nu(patchi);
205  }
206 
207  return tmu;
208 }
209 
210 
213 {
214  PtrDictionary<phase>::const_iterator iter = phases_.begin();
215 
216  tmp<surfaceScalarField> tmuf =
217  fvc::interpolate(iter())*iter().rho()*fvc::interpolate(iter().nu());
218 
219  for (++iter; iter != phases_.end(); ++iter)
220  {
221  tmuf() +=
222  fvc::interpolate(iter())*iter().rho()*fvc::interpolate(iter().nu());
223  }
224 
225  return tmuf;
226 }
227 
228 
231 {
232  return nu_;
233 }
234 
235 
238 {
239  return nu_.boundaryField()[patchi];
240 }
241 
242 
245 {
246  return muf()/fvc::interpolate(rho());
247 }
248 
249 
252 {
253  tmp<surfaceScalarField> tstf
254  (
256  (
257  IOobject
258  (
259  "surfaceTensionForce",
260  mesh_.time().timeName(),
261  mesh_
262  ),
263  mesh_,
265  (
266  "surfaceTensionForce",
267  dimensionSet(1, -2, -2, 0, 0),
268  0.0
269  )
270  )
271  );
272 
273  surfaceScalarField& stf = tstf();
274 
275  forAllConstIter(PtrDictionary<phase>, phases_, iter1)
276  {
277  const phase& alpha1 = iter1();
278 
279  PtrDictionary<phase>::const_iterator iter2 = iter1;
280  ++iter2;
281 
282  for (; iter2 != phases_.end(); ++iter2)
283  {
284  const phase& alpha2 = iter2();
285 
286  sigmaTable::const_iterator sigma =
287  sigmas_.find(interfacePair(alpha1, alpha2));
288 
289  if (sigma == sigmas_.end())
290  {
292  << "Cannot find interface " << interfacePair(alpha1, alpha2)
293  << " in list of sigma values"
294  << exit(FatalError);
295  }
296 
297  stf += dimensionedScalar("sigma", dimSigma_, sigma())
299  (
302  );
303  }
304  }
305 
306  return tstf;
307 }
308 
309 
311 {
312  correct();
313 
314  const Time& runTime = mesh_.time();
315 
316  volScalarField& alpha = phases_.first();
317 
318  const dictionary& alphaControls = mesh_.solverDict("alpha");
319  label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")));
320  scalar cAlpha(readScalar(alphaControls.lookup("cAlpha")));
321 
322  if (nAlphaSubCycles > 1)
323  {
324  surfaceScalarField rhoPhiSum
325  (
326  IOobject
327  (
328  "rhoPhiSum",
329  runTime.timeName(),
330  mesh_
331  ),
332  mesh_,
333  dimensionedScalar("0", rhoPhi_.dimensions(), 0)
334  );
335 
336  dimensionedScalar totalDeltaT = runTime.deltaT();
337 
338  for
339  (
340  subCycle<volScalarField> alphaSubCycle(alpha, nAlphaSubCycles);
341  !(++alphaSubCycle).end();
342  )
343  {
344  solveAlphas(cAlpha);
345  rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi_;
346  }
347 
348  rhoPhi_ = rhoPhiSum;
349  }
350  else
351  {
352  solveAlphas(cAlpha);
353  }
354 
355  // Update the mixture kinematic viscosity
356  nu_ = mu()/rho();
357 }
358 
359 
361 {
362  forAllIter(PtrDictionary<phase>, phases_, iter)
363  {
364  iter().correct();
365  }
366 }
367 
368 
370 (
371  const volScalarField& alpha1,
372  const volScalarField& alpha2
373 ) const
374 {
375  /*
376  // Cell gradient of alpha
377  volVectorField gradAlpha =
378  alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
379 
380  // Interpolated face-gradient of alpha
381  surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
382  */
383 
384  surfaceVectorField gradAlphaf
385  (
388  );
389 
390  // Face unit interface normal
391  return gradAlphaf/(mag(gradAlphaf) + deltaN_);
392 }
393 
394 
396 (
397  const volScalarField& alpha1,
398  const volScalarField& alpha2
399 ) const
400 {
401  // Face unit interface normal flux
402  return nHatfv(alpha1, alpha2) & mesh_.Sf();
403 }
404 
405 
406 // Correction for the boundary condition on the unit normal nHat on
407 // walls to produce the correct contact angle.
408 
409 // The dynamic contact angle is calculated from the component of the
410 // velocity on the direction of the interface, parallel to the wall.
411 
413 (
414  const phase& alpha1,
415  const phase& alpha2,
416  surfaceVectorField::GeometricBoundaryField& nHatb
417 ) const
418 {
419  const volScalarField::GeometricBoundaryField& gbf
420  = alpha1.boundaryField();
421 
422  const fvBoundaryMesh& boundary = mesh_.boundary();
423 
425  {
426  if (isA<alphaContactAngleFvPatchScalarField>(gbf[patchi]))
427  {
428  const alphaContactAngleFvPatchScalarField& acap =
429  refCast<const alphaContactAngleFvPatchScalarField>(gbf[patchi]);
430 
431  vectorField& nHatPatch = nHatb[patchi];
432 
433  vectorField AfHatPatch
434  (
435  mesh_.Sf().boundaryField()[patchi]
436  /mesh_.magSf().boundaryField()[patchi]
437  );
438 
439  alphaContactAngleFvPatchScalarField::thetaPropsTable::
440  const_iterator tp =
441  acap.thetaProps().find(interfacePair(alpha1, alpha2));
442 
443  if (tp == acap.thetaProps().end())
444  {
446  << "Cannot find interface " << interfacePair(alpha1, alpha2)
447  << "\n in table of theta properties for patch "
448  << acap.patch().name()
449  << exit(FatalError);
450  }
451 
452  bool matched = (tp.key().first() == alpha1.name());
453 
454  scalar theta0 = convertToRad*tp().theta0(matched);
455  scalarField theta(boundary[patchi].size(), theta0);
456 
457  scalar uTheta = tp().uTheta();
458 
459  // Calculate the dynamic contact angle if required
460  if (uTheta > SMALL)
461  {
462  scalar thetaA = convertToRad*tp().thetaA(matched);
463  scalar thetaR = convertToRad*tp().thetaR(matched);
464 
465  // Calculated the component of the velocity parallel to the wall
466  vectorField Uwall
467  (
468  U_.boundaryField()[patchi].patchInternalField()
469  - U_.boundaryField()[patchi]
470  );
471  Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
472 
473  // Find the direction of the interface parallel to the wall
474  vectorField nWall
475  (
476  nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
477  );
478 
479  // Normalise nWall
480  nWall /= (mag(nWall) + SMALL);
481 
482  // Calculate Uwall resolved normal to the interface parallel to
483  // the interface
484  scalarField uwall(nWall & Uwall);
485 
486  theta += (thetaA - thetaR)*tanh(uwall/uTheta);
487  }
488 
489 
490  // Reset nHatPatch to correspond to the contact angle
491 
492  scalarField a12(nHatPatch & AfHatPatch);
493 
494  scalarField b1(cos(theta));
495 
496  scalarField b2(nHatPatch.size());
497 
498  forAll(b2, facei)
499  {
500  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
501  }
502 
503  scalarField det(1.0 - a12*a12);
504 
505  scalarField a((b1 - a12*b2)/det);
506  scalarField b((b2 - a12*b1)/det);
507 
508  nHatPatch = a*AfHatPatch + b*nHatPatch;
509 
510  nHatPatch /= (mag(nHatPatch) + deltaN_.value());
511  }
512  }
513 }
514 
515 
517 (
518  const phase& alpha1,
519  const phase& alpha2
520 ) const
521 {
522  tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
523 
524  correctContactAngle(alpha1, alpha2, tnHatfv().boundaryField());
525 
526  // Simple expression for curvature
527  return -fvc::div(tnHatfv & mesh_.Sf());
528 }
529 
530 
533 {
534  tmp<volScalarField> tnearInt
535  (
536  new volScalarField
537  (
538  IOobject
539  (
540  "nearInterface",
541  mesh_.time().timeName(),
542  mesh_
543  ),
544  mesh_,
545  dimensionedScalar("nearInterface", dimless, 0.0)
546  )
547  );
548 
549  forAllConstIter(PtrDictionary<phase>, phases_, iter)
550  {
551  tnearInt() = max(tnearInt(), pos(iter() - 0.01)*pos(0.99 - iter()));
552  }
553 
554  return tnearInt;
555 }
556 
557 
559 (
560  const scalar cAlpha
561 )
562 {
563  static label nSolves=-1;
564  nSolves++;
565 
566  word alphaScheme("div(phi,alpha)");
567  word alpharScheme("div(phirb,alpha)");
568 
569  surfaceScalarField phic(mag(phi_/mesh_.magSf()));
570  phic = min(cAlpha*phic, max(phic));
571 
572  PtrList<surfaceScalarField> alphaPhiCorrs(phases_.size());
573  int phasei = 0;
574 
575  forAllIter(PtrDictionary<phase>, phases_, iter)
576  {
577  phase& alpha = iter();
578 
579  alphaPhiCorrs.set
580  (
581  phasei,
583  (
584  "phi" + alpha.name() + "Corr",
585  fvc::flux
586  (
587  phi_,
588  alpha,
589  alphaScheme
590  )
591  )
592  );
593 
594  surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phasei];
595 
596  forAllIter(PtrDictionary<phase>, phases_, iter2)
597  {
598  phase& alpha2 = iter2();
599 
600  if (&alpha2 == &alpha) continue;
601 
603 
604  alphaPhiCorr += fvc::flux
605  (
607  alpha,
609  );
610  }
611 
613  (
614  1.0/mesh_.time().deltaT().value(),
615  geometricOneField(),
616  alpha,
617  phi_,
618  alphaPhiCorr,
619  zeroField(),
620  zeroField(),
621  1,
622  0,
623  true
624  );
625 
626  phasei++;
627  }
628 
629  MULES::limitSum(alphaPhiCorrs);
630 
631  rhoPhi_ = dimensionedScalar("0", dimensionSet(1, 0, -1, 0, 0), 0);
632 
633  volScalarField sumAlpha
634  (
635  IOobject
636  (
637  "sumAlpha",
638  mesh_.time().timeName(),
639  mesh_
640  ),
641  mesh_,
642  dimensionedScalar("sumAlpha", dimless, 0)
643  );
644 
645  phasei = 0;
646 
647  forAllIter(PtrDictionary<phase>, phases_, iter)
648  {
649  phase& alpha = iter();
650 
651  surfaceScalarField& alphaPhi = alphaPhiCorrs[phasei];
652  alphaPhi += upwind<scalar>(mesh_, phi_).flux(alpha);
653 
655  (
656  geometricOneField(),
657  alpha,
658  alphaPhi,
659  zeroField(),
660  zeroField()
661  );
662 
663  rhoPhi_ += alphaPhi*alpha.rho();
664 
665  Info<< alpha.name() << " volume fraction, min, max = "
666  << alpha.weightedAverage(mesh_.V()).value()
667  << ' ' << min(alpha).value()
668  << ' ' << max(alpha).value()
669  << endl;
670 
671  sumAlpha += alpha;
672 
673  phasei++;
674  }
675 
676  Info<< "Phase-sum volume fraction, min, max = "
677  << sumAlpha.weightedAverage(mesh_.V()).value()
678  << ' ' << min(sumAlpha).value()
679  << ' ' << max(sumAlpha).value()
680  << endl;
681 
682  calcAlphas();
683 }
684 
685 
687 {
688  if (transportModel::read())
689  {
690  bool readOK = true;
691 
692  PtrList<entry> phaseData(lookup("phases"));
693  label phasei = 0;
694 
695  forAllIter(PtrDictionary<phase>, phases_, iter)
696  {
697  readOK &= iter().read(phaseData[phasei++].dict());
698  }
699 
700  lookup("sigmas") >> sigmas_;
701 
702  return readOK;
703  }
704  else
705  {
706  return false;
707  }
708 }
709 
710 
711 // ************************************************************************* //
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::multiphaseMixture::convertToRad
static const scalar convertToRad
Conversion factor for degrees into radians.
Definition: multiphaseMixture.H:157
Foam::multiphaseMixture::calcAlphas
void calcAlphas()
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::MULES::explicitSolve
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Definition: MULESTemplates.C:38
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::multiphaseMixture::K
tmp< volScalarField > K(const phase &alpha1, const phase &alpha2) const
subCycle.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::multiphaseMixture::multiphaseMixture
multiphaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
Foam::multiphaseMixture::correct
void correct()
Correct the mixture properties.
fvcSnGrad.H
Calculate the snGrad of the given volField.
boundary
faceListList boundary(nPatches)
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
fvcDiv.H
Calculate the divergence of the given field.
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:235
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
tp
volScalarField tp(IOobject("tp", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar("tp", dimensionSet(0, 0, 1, 0, 0, 0, 0), 0))
Foam::multiphaseMixture::nu
tmp< volScalarField > nu() const
Return the kinematic laminar viscosity.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
phir
surfaceScalarField phir(phic *interface.nHatf())
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
U
U
Definition: pEqn.H:46
constant
Constant dispersed-phase particle diameter model.
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:49
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:56
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::transportModel::read
virtual bool read()=0
Read transportProperties dictionary.
Definition: transportModel.C:50
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
alphaPhi
surfaceScalarField alphaPhi(IOobject("alphaPhi", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), phi *fvc::interpolate(alpha1))
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::multiphaseMixture::nuf
tmp< surfaceScalarField > nuf() const
Return the face-interpolated dynamic laminar viscosity.
nAlphaSubCycles
label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")))
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
Foam::Info
messageStream Info
Foam::multiphaseMixture::nHatfv
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
correct
fvOptions correct(rho)
Foam::MULES::limitSum
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:51
Foam::multiphaseMixture::surfaceTensionForce
tmp< surfaceScalarField > surfaceTensionForce() const
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::multiphaseMixture::read
bool read()
Read base transportProperties dictionary.
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
alpha2
alpha2
Definition: alphaEqn.H:112
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
boundaryField
cellIbMask *cellIbMaskExt *faceIbMask *cellIbMask boundaryField().evaluateCoupled()
Foam::multiphaseMixture::solve
void solve()
Solve for the mixture phase-fractions.
surfaceInterpolate.H
Surface Interpolation.
Foam::MULES::limit
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin, const bool returnCorr)
Definition: MULESTemplates.C:526
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::multiphaseMixture::rho
tmp< volScalarField > rho() const
Return the mixture density.
phic
phic
Definition: alphaEqnsSubCycle.H:5
alpha1
volScalarField & alpha1
Definition: createFields.H:15
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
trho
tmp< volScalarField > trho
Definition: setRegionSolidFields.H:4
Foam::multiphaseMixture::mu
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
alpharScheme
word alpharScheme("div(phirb,alpha)")
phasei
label phasei
Definition: pEqn.H:37
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::multiphaseMixture::solveAlphas
void solveAlphas(const scalar cAlpha)
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:259
fvcGrad.H
Calculate the gradient of the given field.
fvcFlux.H
Calculate the face-flux of the given field.
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::multiphaseMixture::nearInterface
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
multiphaseMixture.H
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::det
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:60
Foam::fvc::flux
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcFlux.C:45
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:55
MULES.H
MULES: Multidimensional universal limiter for explicit solution.
Foam::multiphaseMixture::nHatf
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Foam::multiphaseMixture::correctContactAngle
void correctContactAngle(const phase &alpha1, const phase &alpha2, surfaceVectorField::GeometricBoundaryField &nHatb) const
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::multiphaseMixture::alphas_
volScalarField alphas_
Definition: multiphaseMixture.H:143
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:335
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::multiphaseMixture::phases_
PtrDictionary< phase > phases_
Dictionary of phases.
Definition: multiphaseMixture.H:136
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256
Foam::multiphaseMixture::muf
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190
alphaControls
const dictionary & alphaControls
Definition: alphaControls.H:1