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 -------------------------------------------------------------------------------
11 License
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 
46 void 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 
59 void 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",
74  fvc::flux
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.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(),
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);
158 
160  (
161  geometricOneField(),
162  phase,
163  alphaPhi
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 
197 Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseSystem::nHatfv
198 (
199  const volScalarField& alpha1,
200  const volScalarField& alpha2
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 
223 Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::nHatf
224 (
225  const volScalarField& alpha1,
226  const volScalarField& alpha2
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 
240 void Foam::multiphaseSystem::correctContactAngle
241 (
242  const phaseModel& phase1,
243  const phaseModel& phase2,
245 ) const
246 {
247  const volScalarField::Boundary& gbf
248  = phase1.boundaryField();
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.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 
355 Foam::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,
375  const surfaceScalarField& phi
376 )
377 :
379  (
380  IOobject
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  (
397  IOobject
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 
481 Foam::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 
514 Foam::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  tmp<volScalarField> tCvm
540  (
541  new volScalarField
542  (
543  IOobject
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.good())
564  {
565  tCvm.ref() += iterCvm()*phase2.rho()*phase2;
566  }
567  else
568  {
569  iterCvm = Cvms_.cfind(interfacePair(phase2, phase));
570 
571  if (iterCvm.good())
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 {
587  tmp<volVectorField> tSvm
588  (
589  new volVectorField
590  (
591  IOobject
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.good())
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.good())
625  {
626  tSvm.ref() += Cvm()*phase.rho()*phase2*phase2.DDtU();
627  }
628  }
629  }
630 
631  volVectorField::Boundary& SvmBf =
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(),
669  dm.residualPhaseFraction()
670  )
671  *dm.K
672  (
673  max
674  (
675  mag(dm.phase1().U() - dm.phase2().U()),
676  dm.residualSlip()
677  )
678  )
679  ).ptr();
680 
681  volScalarField::Boundary& Kbf = Kptr->boundaryFieldRef();
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,
708  const dragCoeffFields& dragCoeffs
709 ) const
710 {
711  tmp<volScalarField> tdragCoeff
712  (
713  new volScalarField
714  (
715  IOobject
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();
732  dragCoeffFields::const_iterator dcIter = dragCoeffs.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  (
763  IOobject
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.good())
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  (
810  new volScalarField
811  (
812  IOobject
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  {
855  volScalarField& alpha = phase;
856 
857  alpha0s.set
858  (
859  phasei,
860  new volScalarField(alpha.oldTime())
861  );
862 
863  alphaPhiSums.set
864  (
865  phasei,
867  (
868  IOobject
869  (
870  "phiSum" + alpha.name(),
871  runTime.timeName(),
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  {
906  volScalarField& alpha = phase;
907 
908  phase.alphaPhi() = alphaPhiSums[phasei];
909 
910  // Correct the time index of the field
911  // to correspond to the global time
912  alpha.timeIndex() = runTime.timeIndex();
913 
914  // Reset the old-time field value
915  alpha.oldTime() = alpha0s[phasei];
916  alpha.oldTime().timeIndex() = runTime.timeIndex();
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 // ************************************************************************* //
tmp< volScalarField > nu() const
Return the mixture laminar viscosity.
faceListList boundary
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
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)
dimensionedScalar tanh(const dimensionedScalar &ds)
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:233
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensionedScalar acos(const dimensionedScalar &ds)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
label phasei
Definition: pEqn.H:27
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual bool read()
Read object.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
bool found(const Key &key) const
Same as contains()
Definition: HashTable.H:1354
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
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:489
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
Unit conversion functions.
Type & refCast(U &obj)
A dynamic_cast (for references). Generates a FatalError on failed casts and uses the virtual type() m...
Definition: typeInfo.H:159
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< volScalarField > trho
autoPtr< multiphaseSystem::dragCoeffFields > dragCoeffs(fluid.dragCoeffs())
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Calculate the snGrad of the given volField.
const volScalarField & alpha2
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
dimensionedScalar det(const dimensionedSphericalTensor &dt)
const dimensionSet dimless
Dimensionless.
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.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
phaseModel & phase2
CGAL::Exact_predicates_exact_constructions_kernel K
const surfaceScalarField & phi() const
Definition: phaseModel.H:218
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
const volVectorField & DDtU() const
Definition: phaseModel.H:208
Area-weighted average a surfaceField creating a volField.
tmp< volScalarField > rho() const
Return the mixture density.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const surfaceScalarField & alphaPhi() const
Definition: phaseModel.H:228
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const volVectorField & U() const
Definition: phaseModel.H:198
word timeName
Definition: getTimeIndex.H:3
tmp< volScalarField > Cvm(const phaseModel &phase) const
Return the virtual-mass coefficient for the given phase.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
const dictionary & alphaControls
Definition: alphaControls.H:1
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:27
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Calculate the gradient of the given field.
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Time & time() const noexcept
Return time registry.
Calculate the face-flux of the given field.
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
dimensionedScalar cbrt(const dimensionedScalar &ds)
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
label timeIndex() const noexcept
Return the current time index.
Definition: TimeStateI.H:43
tmp< volScalarField > dragCoeff(const phaseModel &phase, const dragCoeffFields &dragCoeffs) const
Return the sum of the drag coefficients for the given phase.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
void solve()
Solve for the mixture phase-fractions.
Calculate the divergence of the given field.
dimensionedScalar pos0(const dimensionedScalar &ds)
tmp< volVectorField > Svm(const phaseModel &phase) const
Return the virtual-mass source for the given phase.
phaseModel & phase1
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
bool read()
Read base transportProperties dictionary.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
static autoPtr< dragModel > New(const dictionary &dict, const phaseModel &phase1, const phaseModel &phase2)
Definition: dragModel.C:57
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.
const dimensionSet dimDensity
const word & name() const noexcept
Return const reference to name.
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
const dimensionedScalar mu
Atomic mass unit.
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
U
Definition: pEqn.H:72
const dimensionedScalar & rho() const
Definition: phaseModel.H:193
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:50
const word & name() const
Definition: phaseModel.H:166
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
Field< vector > vectorField
Specialisation of Field<T> for vector.
MULES: Multidimensional universal limiter for explicit solution.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:61
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
autoPtr< dragCoeffFields > dragCoeffs() const
Return the drag coefficients for all of the interfaces.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const volScalarField & alpha1