phaseSystem.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2025 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 "phaseSystem.H"
28 #include "surfaceInterpolate.H"
29 #include "fvcDdt.H"
30 #include "localEulerDdtScheme.H"
31 #include "fvcDiv.H"
32 #include "fvcGrad.H"
33 #include "fvcSnGrad.H"
34 #include "fvCorrectPhi.H"
35 #include "fvcMeshPhi.H"
38 #include "correctContactAngle.H"
42 #include "pressureReference.H"
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
49 }
50 
51 
52 const Foam::word Foam::phaseSystem::propertiesName("phaseProperties");
53 
54 
55 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
56 
58 {
60  (
62  (
63  "sumAlphaMoving",
65  calculatedFvPatchScalarField::typeName
66  )
67  );
68 
69  for
70  (
71  label movingPhasei=1;
72  movingPhasei<movingPhaseModels_.size();
73  movingPhasei++
74  )
75  {
76  sumAlphaMoving.ref() += movingPhaseModels_[movingPhasei];
77  }
78 
79  return sumAlphaMoving;
80 }
81 
82 
84 {
85  // Calculate the mean velocity difference with respect to Um0
86  // from the current velocity of the moving phases
87  volVectorField dUm(Um0);
88 
89  forAll(movingPhaseModels_, movingPhasei)
90  {
91  dUm -=
92  movingPhaseModels_[movingPhasei]
93  *movingPhaseModels_[movingPhasei].U();
94  }
95 
96  forAll(movingPhaseModels_, movingPhasei)
97  {
98  movingPhaseModels_[movingPhasei].URef() += dUm;
99  }
100 }
101 
102 
104 (
105  const PtrList<surfaceScalarField>& alphafs,
106  const surfaceScalarField& phim0
107 )
108 {
109  // Calculate the mean flux difference with respect to phim0
110  // from the current flux of the moving phases
111  surfaceScalarField dphim(phim0);
112 
113  forAll(movingPhaseModels_, movingPhasei)
114  {
115  dphim -=
116  alphafs[movingPhaseModels_[movingPhasei].index()]
117  *movingPhaseModels_[movingPhasei].phi();
118  }
119 
120  forAll(movingPhaseModels_, movingPhasei)
121  {
122  movingPhaseModels_[movingPhasei].phiRef() += dphim;
123  }
124 }
125 
126 
128 (
129  const volScalarField& alpha1,
130  const volScalarField& alpha2
131 ) const
132 {
133  /*
134  // Cell gradient of alpha
135  volVectorField gradAlpha =
136  alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
137 
138  // Interpolated face-gradient of alpha
139  surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
140  */
141 
142  surfaceVectorField gradAlphaf
143  (
145  - fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2))
146  );
147 
148  // Face unit interface normal
149  return gradAlphaf/(mag(gradAlphaf) + deltaN_);
150 }
151 
152 
154 (
155  const volScalarField& alpha1,
156  const volScalarField& alpha2
157 ) const
158 {
159  // Face unit interface normal flux
160  return nHatfv(alpha1, alpha2) & mesh_.Sf();
161 }
162 
163 
165 (
166  const phaseModel& phase1,
167  const phaseModel& phase2
168 ) const
169 {
170  tmp<surfaceVectorField> tnHatfv = nHatfv(phase1, phase2);
171 
173  (
174  phase1,
175  phase2,
176  phase1.U()().boundaryField(),
177  deltaN_,
178  tnHatfv.ref().boundaryFieldRef()
179  );
180 
181  // Simple expression for curvature
182  return -fvc::div(tnHatfv & mesh_.Sf());
183 }
184 
185 
186 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
187 
189 (
190  const fvMesh& mesh
191 )
192 :
194  (
195  IOobject
196  (
197  propertiesName,
198  mesh.time().constant(),
199  mesh,
200  IOobject::MUST_READ_IF_MODIFIED,
201  IOobject::NO_WRITE
202  )
203  ),
204 
205  mesh_(mesh),
206 
207  pimple_(mesh_.lookupObject<pimpleNoLoopControl>("solutionControl")),
208 
209  MRF_(mesh_),
210 
211  referencePhaseName_(lookupOrDefault("referencePhase", word::null)),
212 
213  phaseModels_
214  (
215  lookup("phases"),
216  phaseModel::iNew(*this, referencePhaseName_)
217  ),
218 
219  phi_
220  (
221  IOobject
222  (
223  "phi",
224  mesh.time().name(),
225  mesh
226  ),
227  mesh,
229  ),
230 
231  dpdt_
232  (
233  IOobject
234  (
235  "dpdt",
236  mesh.time().name(),
237  mesh
238  ),
239  mesh,
241  ),
242 
243  cAlphas_
244  (
245  found("interfaceCompression")
246  ? generateInterfacialValues<scalar>
247  (
248  *this,
249  subDict("interfaceCompression")
250  )
251  : cAlphaTable()
252  ),
253 
254  deltaN_
255  (
256  "deltaN",
257  1e-8/pow(average(mesh_.V()), 1.0/3.0)
258  ),
259 
260  surfaceTensionCoefficientModels_
261  (
263  (
264  *this,
266  )
267  )
268 {
269  // Groupings
270  label movingPhasei = 0;
271  label stationaryPhasei = 0;
272  label thermalPhasei = 0;
273  label multicomponentPhasei = 0;
274  forAll(phaseModels_, phasei)
275  {
276  phaseModel& phase = phaseModels_[phasei];
277  movingPhasei += !phase.stationary();
278  stationaryPhasei += phase.stationary();
279  thermalPhasei += !phase.isothermal();
280  multicomponentPhasei += !phase.pure();
281  }
282  movingPhaseModels_.resize(movingPhasei);
283  stationaryPhaseModels_.resize(stationaryPhasei);
284  thermalPhaseModels_.resize(thermalPhasei);
285  multicomponentPhaseModels_.resize(multicomponentPhasei);
286 
287  movingPhasei = 0;
288  stationaryPhasei = 0;
289  thermalPhasei = 0;
290  multicomponentPhasei = 0;
291  forAll(phaseModels_, phasei)
292  {
293  phaseModel& phase = phaseModels_[phasei];
294  if (!phase.stationary())
295  {
296  movingPhaseModels_.set(movingPhasei++, &phase);
297  }
298  if (phase.stationary())
299  {
300  stationaryPhaseModels_.set(stationaryPhasei++, &phase);
301  }
302  if (!phase.isothermal())
303  {
304  thermalPhaseModels_.set(thermalPhasei++, &phase);
305  }
306  if (!phase.pure())
307  {
308  multicomponentPhaseModels_.set(multicomponentPhasei++, &phase);
309  }
310  }
311 
312  forAll(movingPhaseModels_, movingPhasei)
313  {
314  phi_ +=
315  fvc::interpolate(movingPhaseModels_[movingPhasei])
316  *movingPhaseModels_[movingPhasei].phi();
317  }
318 
319  // Write phi
321 
322  // Update motion fields
324 
325  // Set the optional reference phase fraction from the other phases
327  {
328  phaseModel* referencePhasePtr = &phases()[referencePhaseName_];
329  volScalarField& referenceAlpha = *referencePhasePtr;
330 
331  referenceAlpha = 1;
332 
333  forAll(phaseModels_, phasei)
334  {
335  if (&phaseModels_[phasei] != referencePhasePtr)
336  {
337  referenceAlpha -= phaseModels_[phasei];
338  }
339  }
340  }
341 
342  forAll(phaseModels_, phasei)
343  {
344  const volScalarField& alphai = phases()[phasei];
345  mesh_.schemes().setFluxRequired(alphai.name());
346  }
347 
348  // Check for and warn about the type entry being present
349  if (found("type"))
350  {
352  << "The phase system type entry - type - in "
353  << relativeObjectPath() << " is no longer used"
354  << endl;
355  }
356 
357  // Check for and warn/error about phase change models being in the old
358  // location in constant/phaseProperties
359  const wordList modelEntries
360  ({
361  "phaseTransfer",
362  "saturationTemperature",
363  "interfaceComposition",
364  "diffusiveMassTransfer"
365  });
366 
367  OStringStream modelEntriesString;
368  forAll(modelEntries, i)
369  {
370  modelEntriesString<< modelEntries[i];
371  if (i < modelEntries.size() - 2) modelEntriesString<< ", ";
372  if (i == modelEntries.size() - 2) modelEntriesString<< " and ";
373  }
374 
375  label warnOrError = 0; // 1 == warn, 2 = error
376  forAll(modelEntries, i)
377  {
378  if (!found(modelEntries[i])) continue;
379 
380  warnOrError =
381  isDict(modelEntries[i]) && !subDict(modelEntries[i]).empty()
382  ? 2
383  : 1;
384  }
385 
386  OStringStream msg;
387  if (warnOrError != 0)
388  {
389  msg << "Phase change model entries - "
390  << modelEntriesString.str().c_str() << " - in "
391  << relativeObjectPath() << " are no longer used. These models "
392  << "are now specified as fvModels.";
393  }
394  if (warnOrError == 1)
395  {
397  << msg.str().c_str() << endl;
398  }
399  if (warnOrError == 2)
400  {
402  << msg.str().c_str() << exit(FatalIOError);
403  }
404 }
405 
406 
407 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
408 
410 {}
411 
412 
413 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
414 
416 {
417  nAlphaSubCyclesPtr = Function1<scalar>::New
418  (
419  dict.found("nAlphaSubCycles")
420  ? "nAlphaSubCycles"
421  : "nSubCycles",
422  dimless,
423  dimless,
424  dict
425  );
426 
428  (
429  {"nCorrectors", "nAlphaCorr"},
430  1
431  );
432 
433  MULESCorr = dict.lookupOrDefault<Switch>("MULESCorr", false);
434 
435  MULESCorrRelax = dict.lookupOrDefault<scalar>("MULESCorrRelax", 0.5);
436 
437  vDotResidualAlpha =
438  dict.lookupOrDefault("vDotResidualAlpha", 1e-4);
439 
440  MULES.read(dict);
441 
442  clip = dict.lookupOrDefault<Switch>("clip", true);
443 }
444 
445 
447 {
448  nAlphaSubCycles = ceil(nAlphaSubCyclesPtr->value(CoNum));
449 }
450 
451 
453 {
454  tmp<volScalarField> rho(movingPhaseModels_[0]*movingPhaseModels_[0].rho());
455 
456  for
457  (
458  label movingPhasei=1;
459  movingPhasei<movingPhaseModels_.size();
460  movingPhasei++
461  )
462  {
463  rho.ref() +=
464  movingPhaseModels_[movingPhasei]
465  *movingPhaseModels_[movingPhasei].rho();
466  }
467 
468  if (stationaryPhaseModels_.empty())
469  {
470  return rho;
471  }
472  else
473  {
474  return rho/sumAlphaMoving();
475  }
476 }
477 
478 
480 {
481  tmp<volVectorField> U(movingPhaseModels_[0]*movingPhaseModels_[0].U());
482 
483  for
484  (
485  label movingPhasei=1;
486  movingPhasei<movingPhaseModels_.size();
487  movingPhasei++
488  )
489  {
490  U.ref() +=
491  movingPhaseModels_[movingPhasei]
492  *movingPhaseModels_[movingPhasei].U();
493  }
494 
495  if (stationaryPhaseModels_.empty())
496  {
497  return U;
498  }
499  else
500  {
501  return U/sumAlphaMoving();
502  }
503 }
504 
505 
508 {
509  if (surfaceTensionCoefficientModels_.found(key))
510  {
511  return surfaceTensionCoefficientModels_[key]->sigma();
512  }
513  else
514  {
515  return volScalarField::New
516  (
517  surfaceTensionCoefficientModel::typeName + ":sigma",
518  mesh_,
520  );
521  }
522 }
523 
524 
527 {
528  if (surfaceTensionCoefficientModels_.found(key))
529  {
530  return surfaceTensionCoefficientModels_[key]->sigma(patchi);
531  }
532  else
533  {
534  return tmp<scalarField>
535  (
536  new scalarField(mesh_.boundary()[patchi].size(), 0)
537  );
538  }
539 }
540 
541 
544 {
545  tmp<volScalarField> tnearInt
546  (
548  (
549  "nearInterface",
550  mesh_,
552  )
553  );
554 
555  forAll(phases(), phasei)
556  {
557  tnearInt.ref() = max
558  (
559  tnearInt(),
560  pos0(phases()[phasei] - 0.01)*pos0(0.99 - phases()[phasei])
561  );
562  }
563 
564  return tnearInt;
565 }
566 
567 
569 (
570  const phaseModel& phase1
571 ) const
572 {
573  tmp<surfaceScalarField> tSurfaceTension
574  (
576  (
577  "surfaceTension",
578  mesh_,
579  dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), 0)
580  )
581  );
582 
583  forAll(phases(), phasej)
584  {
585  const phaseModel& phase2 = phases()[phasej];
586 
587  if (&phase2 != &phase1)
588  {
589  const phaseInterface interface(phase1, phase2);
590 
591  if (cAlphas_.found(interface))
592  {
593  tSurfaceTension.ref() +=
594  fvc::interpolate(sigma(interface)*K(phase1, phase2))
595  *(
596  fvc::interpolate(phase2)*fvc::snGrad(phase1)
597  - fvc::interpolate(phase1)*fvc::snGrad(phase2)
598  );
599  }
600  }
601  }
602 
603  return tSurfaceTension;
604 }
605 
606 
608 {
609  forAll(phaseModels_, phasei)
610  {
611  if (!phaseModels_[phasei].incompressible())
612  {
613  return false;
614  }
615  }
616 
617  return true;
618 }
619 
620 
622 {
623  forAll(phaseModels_, phasei)
624  {
625  phaseModels_[phasei].correct();
626  }
627 }
628 
629 
631 (
633 )
634 {
635  forAll(movingPhaseModels_, movingPhasei)
636  {
637  phaseModel& phase = movingPhaseModels_[movingPhasei];
638  const volScalarField& alpha = phase;
639  volScalarField& rho = phase.rho();
640 
641  volScalarField source
642  (
644  (
645  IOobject::groupName("source", phase.name()),
646  mesh_,
648  )
649  );
650 
651  if (fvModels().addsSupToField(rho.name()))
652  {
653  source += fvModels().source(alpha, rho)&rho;
654  }
655 
656  if (dmdts.set(phase.index()))
657  {
658  source.internalFieldRef() += dmdts[phase.index()];
659  }
660 
661  phase.correctContinuityError(source);
662  }
663 }
664 
665 
667 {
668  bool updateDpdt = false;
669 
670  forAll(phaseModels_, phasei)
671  {
672  phaseModels_[phasei].correctKinematics();
673 
674  // Update the pressure time-derivative if required
675  if (!updateDpdt && phaseModels_[phasei].thermo().dpdt())
676  {
677  dpdt_ = fvc::ddt(phaseModels_[phasei].fluidThermo().p());
678  updateDpdt = true;
679  }
680  }
681 }
682 
683 
685 {
686  forAll(phaseModels_, phasei)
687  {
688  phaseModels_[phasei].correctThermo();
689  }
690 }
691 
692 
694 {
695  forAll(phaseModels_, phasei)
696  {
697  phaseModels_[phasei].correctReactions();
698  }
699 }
700 
701 
703 {
704  forAll(phaseModels_, phasei)
705  {
706  phaseModels_[phasei].correctSpecies();
707  }
708 }
709 
710 
712 {
713  forAll(phaseModels_, phasei)
714  {
715  phaseModels_[phasei].predictMomentumTransport();
716  }
717 }
718 
719 
721 {
722  forAll(phaseModels_, phasei)
723  {
724  phaseModels_[phasei].predictThermophysicalTransport();
725  }
726 }
727 
728 
730 {
731  forAll(phaseModels_, phasei)
732  {
733  phaseModels_[phasei].correctMomentumTransport();
734  }
735 }
736 
737 
739 {
740  forAll(phaseModels_, phasei)
741  {
742  phaseModels_[phasei].correctThermophysicalTransport();
743  }
744 }
745 
746 
748 {
749  if (mesh_.changing())
750  {
751  MRF_.update();
752 
753  // forAll(phaseModels_, phasei)
754  // {
755  // phaseModels_[phasei].meshUpdate();
756  // }
757  }
758 }
759 
760 
762 {
763  forAll(movingPhaseModels_, movingPhasei)
764  {
765  phaseModel& phase = movingPhaseModels_[movingPhasei];
766 
767  tmp<volVectorField> tU(phase.U());
768  const volVectorField::Boundary& UBf = tU().boundaryField();
769 
771  (
772  MRF_.relative(mesh_.Sf().boundaryField() & UBf)
773  );
774 
776 
777  forAll(mesh_.boundary(), patchi)
778  {
779  if
780  (
781  isA<fixedValueFvsPatchScalarField>(phiBf[patchi])
782  && !isA<movingWallVelocityFvPatchVectorField>(UBf[patchi])
783  && !isA<movingWallSlipVelocityFvPatchVectorField>(UBf[patchi])
784  )
785  {
786  phiBf[patchi] == phiRelBf[patchi];
787  }
788  }
789  }
790 }
791 
792 
794 (
795  const volScalarField& p_rgh,
796  const autoPtr<volScalarField>& divU,
799 )
800 {
801  // Update the phase fluxes from the phase face-velocity and make relative
802  forAll(movingPhaseModels_, movingPhasei)
803  {
804  phaseModel& phase = movingPhaseModels_[movingPhasei];
805  phase.phiRef() = mesh_.Sf() & phase.UfRef();
806  MRF_.makeRelative(phase.phiRef());
807  fvc::makeRelative(phase.phiRef(), phase.U());
808  }
809 
810  forAll(movingPhaseModels_, movingPhasei)
811  {
812  phaseModel& phase = movingPhaseModels_[movingPhasei];
813 
816 
817  forAll(Ubf, patchi)
818  {
819  if (Ubf[patchi].fixesValue())
820  {
821  Ubf[patchi].initEvaluate();
822  }
823  }
824 
825  forAll(Ubf, patchi)
826  {
827  if (Ubf[patchi].fixesValue())
828  {
829  Ubf[patchi].evaluate();
830  UfBf[patchi] = Ubf[patchi];
831  }
832  }
833  }
834 
835  // Correct fixed-flux BCs to be consistent with the velocity BCs
836  correctBoundaryFlux();
837 
838  phi_ = Zero;
839  PtrList<surfaceScalarField> alphafs(phaseModels_.size());
840  forAll(movingPhaseModels_, movingPhasei)
841  {
842  phaseModel& phase = movingPhaseModels_[movingPhasei];
843  const label phasei = phase.index();
844  const volScalarField& alpha = phase;
845 
846  alphafs.set(phasei, fvc::interpolate(alpha).ptr());
847 
848  // Calculate absolute flux
849  // from the mapped surface velocity
850  phi_ += alphafs[phasei]*(mesh_.Sf() & phase.UfRef());
851  }
852 
853  if (incompressible())
854  {
856  (
857  phi_,
858  movingPhaseModels_[0].U(),
859  p_rgh,
861  divU,
863  pimple
864  );
865  }
866  else
867  {
869  (
871  (
872  "psi",
873  mesh_,
875  )
876  );
877 
878  forAll(phases(), phasei)
879  {
880  phaseModel& phase = phases()[phasei];
881  const volScalarField& alpha = phase;
882 
883  psi += alpha*phase.fluidThermo().psi()/phase.rho();
884  }
885 
887  (
888  phi_,
889  p_rgh,
890  psi,
892  divU(),
893  pimple
894  );
895  }
896 
897  // Make the flux relative to the mesh motion
898  MRF_.makeRelative(phi_);
899  fvc::makeRelative(phi_, movingPhaseModels_[0].U());
900 
901  setMixturePhi(alphafs, phi_);
902 }
903 
904 
906 {
907  if (regIOobject::read())
908  {
909  bool readOK = true;
910 
911  forAll(phaseModels_, phasei)
912  {
913  readOK &= phaseModels_[phasei].read();
914  }
915 
916  // models ...
917 
918  return readOK;
919  }
920  else
921  {
922  return false;
923  }
924 }
925 
926 
928 {
930  {
931  return fv::localEulerDdt::localRDeltaT(vf.mesh())*vf;
932  }
933  else
934  {
935  return vf/vf.mesh().time().deltaT();
936  }
937 }
938 
939 
941 {
943  {
945  }
946  else
947  {
948  return sf/sf.mesh().time().deltaT();
949  }
950 }
951 
952 
953 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
const Mesh & mesh() const
Return mesh.
Generic field type.
Definition: FieldField.H:77
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Generic GeometricBoundaryField class.
void evaluate()
Evaluate boundary conditions.
Generic GeometricField class.
Internal & internalFieldRef()
Return a reference to the dimensioned internal field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
An STL-conforming hash table.
Definition: HashTable.H:127
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
writeOption writeOpt() const
Definition: IOobject.H:367
fileName relativeObjectPath() const
Return complete relativePath + object name.
Definition: IOobject.H:420
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
Output to memory buffer stream.
Definition: OStringStream.H:52
string str() const
Return the string.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
bool set(const label) const
Is element set.
Definition: UPtrListI.H:87
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void resize(const label)
Reset size of UPtrList. This can only be used to set the size.
Definition: UPtrListI.H:71
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
T lookupOrDefaultBackwardsCompatible(const wordList &, const T &) const
Find and return a T, trying a list of keywords in sequence,.
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:803
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
Dimension set for the base types.
Definition: dimensionSet.H:125
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:56
virtual const volScalarField & psi() const =0
Compressibility [s^2/m^2].
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1799
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:447
static const surfaceScalarField & localRDeltaTf(const fvMesh &mesh)
Return the reciprocal of the local face time-step.
Definition: localEulerDdt.C:58
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:46
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:37
Non-orthogonal solution control class. Provides non-orthogonal-loop control methods.
Word-pair based class used for keying interface models in hash tables.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
virtual bool pure() const =0
Return whether the phase is pure (i.e., not multi-component)
virtual surfaceScalarField & phiRef()=0
Access the volumetric flux.
virtual bool stationary() const =0
Return whether the phase is stationary.
virtual const rhoFluidThermo & fluidThermo() const =0
Return the thermophysical model.
label index() const
Return the index of the phase.
Definition: phaseModel.C:157
virtual tmp< volVectorField > U() const =0
Return the velocity.
virtual bool isothermal() const =0
Return whether the phase is isothermal.
virtual volVectorField & URef()=0
Access the velocity.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
virtual const volScalarField & rho() const =0
Return the density field.
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
Definition: phaseModel.C:205
virtual surfaceVectorField & UfRef()=0
Access the face velocity.
Class to represent a system of phases.
Definition: phaseSystem.H:74
phaseModelPartialList stationaryPhaseModels_
Stationary phase models.
Definition: phaseSystem.H:175
void correctKinematics()
Correct the kinematics.
Definition: phaseSystem.C:666
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
Definition: phaseSystem.C:543
void correctPhi(const volScalarField &p_rgh, const autoPtr< volScalarField > &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pimple)
...
Definition: phaseSystem.C:794
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:252
void correctReactions()
Correct the reactions.
Definition: phaseSystem.C:693
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:156
phaseModelPartialList multicomponentPhaseModels_
Multi-component phase models.
Definition: phaseSystem.H:181
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:169
phaseModelPartialList movingPhaseModels_
Moving phase models.
Definition: phaseSystem.H:172
void setMixtureU(const volVectorField &Um)
Re-normalise the velocity of the phases.
Definition: phaseSystem.C:83
void correctSpecies()
Correct the species mass fractions.
Definition: phaseSystem.C:702
tmp< surfaceScalarField > surfaceTension(const phaseModel &) const
Return the surface tension force.
Definition: phaseSystem.C:569
void correctThermo()
Correct the thermodynamics.
Definition: phaseSystem.C:684
void correct()
Correct the fluid properties other than those listed below.
Definition: phaseSystem.C:621
tmp< volScalarField > sigma(const phaseInterfaceKey &key) const
Return the surface tension coefficient for an interface.
Definition: phaseSystem.C:507
void correctContinuityError(const PtrList< volScalarField::Internal > &dmdts)
Correct the continuity errors.
Definition: phaseSystem.C:631
void predictThermophysicalTransport()
Predict the energy transport e.g. alphat.
Definition: phaseSystem.C:720
tmp< volScalarField > K(const phaseModel &alpha1, const phaseModel &alpha2) const
Curvature of interface between two phases.
Definition: phaseSystem.C:165
tmp< volVectorField > U() const
Return the mixture velocity.
Definition: phaseSystem.C:479
void correctMomentumTransport()
Correct the momentumTransport.
Definition: phaseSystem.C:729
void correctBoundaryFlux()
Correct fixed-flux BCs to be consistent with the velocity BCs.
Definition: phaseSystem.C:761
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:104
void correctThermophysicalTransport()
Correct the energy transport e.g. alphat.
Definition: phaseSystem.C:738
phaseModelPartialList thermalPhaseModels_
Thermal phase models.
Definition: phaseSystem.H:178
void meshUpdate()
Update the fluid properties for mesh changes.
Definition: phaseSystem.C:747
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:184
tmp< volScalarField > sumAlphaMoving() const
Return the sum of the phase fractions of the moving phases.
Definition: phaseSystem.C:57
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Definition: phaseSystem.C:189
void setMixturePhi(const PtrList< surfaceScalarField > &alphafs, const surfaceScalarField &phim)
Re-normalise the flux of the phases.
Definition: phaseSystem.C:104
virtual ~phaseSystem()
Destructor.
Definition: phaseSystem.C:409
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases dotted with face areas.
Definition: phaseSystem.C:154
tmp< volScalarField > rho() const
Return the mixture density.
Definition: phaseSystem.C:452
word referencePhaseName_
Name of optional reference phase which is not solved for.
Definition: phaseSystem.H:166
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases.
Definition: phaseSystem.C:128
void predictMomentumTransport()
Predict the momentumTransport.
Definition: phaseSystem.C:711
virtual bool read()
Read base phaseProperties dictionary.
Definition: phaseSystem.C:905
bool incompressible() const
Return incompressibility.
Definition: phaseSystem.C:607
Pimple no-loop control class. Implements various option flags, but leaves loop controls to the deriva...
Provides controls for the pressure reference in closed-volume simulations.
virtual bool read()
Read object.
Abstract base-class for interface surface-tension models which can be used when interface compression...
static const dimensionSet dimSigma
Coefficient dimensions.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
scalar CoNum
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
Flux correction functions to ensure continuity.
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Calculate the snGrad of the given volField.
label patchi
volScalarField sf(fieldObject, mesh)
const volScalarField & psi
K
Definition: pEqn.H:75
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
void correctPhi(surfaceScalarField &phi, const volVectorField &U, const volScalarField &p, const autoPtr< volScalarField > &rAU, const autoPtr< volScalarField > &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pcorrControl)
Definition: fvCorrectPhi.C:32
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:89
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.
void generateInterfacialModels(PtrList< phaseInterface > &interfaces, PtrList< ModelType > &models, const phaseSystem &fluid, const dictionary &dict, const wordHashSet &ignoreKeys, const phaseInterface &interface, const Args &... args)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:106
dimensionedScalar pos0(const dimensionedScalar &ds)
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
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimVolumetricFlux
const dimensionSet dimless
dimensioned< Type > average(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
SurfaceField< scalar > surfaceScalarField
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:927
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTime
void correctContactAngle(const volScalarField &alpha1, const volScalarField &alpha2, const volVectorField::Boundary &Ubf, const dimensionedScalar &deltaN, surfaceVectorField::Boundary &nHatbf)
Correct the contact angle for the two volume fraction fields.
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
Foam::HashTable< ValueType, Foam::phaseInterfaceKey, Foam::phaseInterfaceKey::hash > generateInterfacialValues(const phaseSystem &fluid, const dictionary &dict, const wordHashSet &ignoreKeys=wordHashSet())
const dimensionSet dimDensity
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31
void correct(const scalar CoNum)
Correct nAlphaSubCycles for the current Courant number.
Definition: phaseSystem.C:446
void read(const dictionary &dict)
Read the alpha and MULES controls from dict.
Definition: phaseSystem.C:415