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-2023 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"
36 #include "correctContactAngle.H"
37 #include "dragModel.H"
39 #include "pressureReference.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
47 }
48 
49 const Foam::word Foam::phaseSystem::propertiesName("phaseProperties");
50 
51 
52 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53 
55 (
56  const phaseModelList& phaseModels
57 ) const
58 {
60  (
62  (
63  "phi",
64  fvc::interpolate(phaseModels[0])*phaseModels[0].phi()
65  )
66  );
67 
68  for (label phasei=1; phasei<phaseModels.size(); phasei++)
69  {
70  tmpPhi.ref() +=
71  fvc::interpolate(phaseModels[phasei])*phaseModels[phasei].phi();
72  }
73 
74  return tmpPhi;
75 }
76 
77 
79 {
80  tmp<volScalarField> sumAlphaMoving
81  (
83  (
84  "sumAlphaMoving",
85  movingPhaseModels_[0],
86  calculatedFvPatchScalarField::typeName
87  )
88  );
89 
90  for
91  (
92  label movingPhasei=1;
93  movingPhasei<movingPhaseModels_.size();
94  movingPhasei++
95  )
96  {
97  sumAlphaMoving.ref() += movingPhaseModels_[movingPhasei];
98  }
99 
100  return sumAlphaMoving;
101 }
102 
103 
105 {
106  // Calculate the mean velocity difference with respect to Um0
107  // from the current velocity of the moving phases
108  volVectorField dUm(Um0);
109 
110  forAll(movingPhaseModels_, movingPhasei)
111  {
112  dUm -=
113  movingPhaseModels_[movingPhasei]
114  *movingPhaseModels_[movingPhasei].U();
115  }
116 
117  forAll(movingPhaseModels_, movingPhasei)
118  {
119  movingPhaseModels_[movingPhasei].URef() += dUm;
120  }
121 }
122 
123 
125 (
126  const PtrList<surfaceScalarField>& alphafs,
127  const surfaceScalarField& phim0
128 )
129 {
130  // Calculate the mean flux difference with respect to phim0
131  // from the current flux of the moving phases
132  surfaceScalarField dphim(phim0);
133 
134  forAll(movingPhaseModels_, movingPhasei)
135  {
136  dphim -=
137  alphafs[movingPhaseModels_[movingPhasei].index()]
138  *movingPhaseModels_[movingPhasei].phi();
139  }
140 
141  forAll(movingPhaseModels_, movingPhasei)
142  {
143  movingPhaseModels_[movingPhasei].phiRef() += dphim;
144  }
145 }
146 
147 
149 (
150  const volScalarField& alpha1,
151  const volScalarField& alpha2
152 ) const
153 {
154  /*
155  // Cell gradient of alpha
156  volVectorField gradAlpha =
157  alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
158 
159  // Interpolated face-gradient of alpha
160  surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
161  */
162 
163  surfaceVectorField gradAlphaf
164  (
166  - fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2))
167  );
168 
169  // Face unit interface normal
170  return gradAlphaf/(mag(gradAlphaf) + deltaN_);
171 }
172 
173 
175 (
176  const volScalarField& alpha1,
177  const volScalarField& alpha2
178 ) const
179 {
180  // Face unit interface normal flux
181  return nHatfv(alpha1, alpha2) & mesh_.Sf();
182 }
183 
184 
186 (
187  const phaseModel& phase1,
188  const phaseModel& phase2
189 ) const
190 {
191  tmp<surfaceVectorField> tnHatfv = nHatfv(phase1, phase2);
192 
194  (
195  phase1,
196  phase2,
197  phase1.U()().boundaryField(),
198  deltaN_,
199  tnHatfv.ref().boundaryFieldRef()
200  );
201 
202  // Simple expression for curvature
203  return -fvc::div(tnHatfv & mesh_.Sf());
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
208 
210 (
211  const fvMesh& mesh
212 )
213 :
215  (
216  IOobject
217  (
218  propertiesName,
219  mesh.time().constant(),
220  mesh,
221  IOobject::MUST_READ_IF_MODIFIED,
222  IOobject::NO_WRITE
223  )
224  ),
225 
226  mesh_(mesh),
227 
228  pimple_(mesh_.lookupObject<pimpleNoLoopControl>("solutionControl")),
229 
230  MRF_(mesh_),
231 
232  referencePhaseName_(lookupOrDefault("referencePhase", word::null)),
233 
234  phaseModels_
235  (
236  lookup("phases"),
237  phaseModel::iNew(*this, referencePhaseName_)
238  ),
239 
240  phi_("phi", calcPhi(phaseModels_)),
241 
242  dpdt_
243  (
244  IOobject
245  (
246  "dpdt",
247  mesh.time().name(),
248  mesh
249  ),
250  mesh,
252  ),
253 
254  deltaN_
255  (
256  "deltaN",
257  1e-8/pow(average(mesh_.V()), 1.0/3.0)
258  )
259 {
260  // Groupings
261  label movingPhasei = 0;
262  label stationaryPhasei = 0;
263  label anisothermalPhasei = 0;
264  label multicomponentPhasei = 0;
265  forAll(phaseModels_, phasei)
266  {
267  phaseModel& phase = phaseModels_[phasei];
268  movingPhasei += !phase.stationary();
269  stationaryPhasei += phase.stationary();
270  anisothermalPhasei += !phase.isothermal();
271  multicomponentPhasei += !phase.pure();
272  }
273  movingPhaseModels_.resize(movingPhasei);
274  stationaryPhaseModels_.resize(stationaryPhasei);
275  anisothermalPhaseModels_.resize(anisothermalPhasei);
276  multicomponentPhaseModels_.resize(multicomponentPhasei);
277 
278  movingPhasei = 0;
279  stationaryPhasei = 0;
280  anisothermalPhasei = 0;
281  multicomponentPhasei = 0;
282  forAll(phaseModels_, phasei)
283  {
284  phaseModel& phase = phaseModels_[phasei];
285  if (!phase.stationary())
286  {
287  movingPhaseModels_.set(movingPhasei++, &phase);
288  }
289  if (phase.stationary())
290  {
291  stationaryPhaseModels_.set(stationaryPhasei++, &phase);
292  }
293  if (!phase.isothermal())
294  {
295  anisothermalPhaseModels_.set(anisothermalPhasei++, &phase);
296  }
297  if (!phase.pure())
298  {
299  multicomponentPhaseModels_.set(multicomponentPhasei++, &phase);
300  }
301  }
302 
303  // Write phi
305 
306  // Interface compression coefficients
307  if (this->found("interfaceCompression"))
308  {
309  generateInterfacialValues("interfaceCompression", cAlphas_);
310  }
311 
312  // Surface tension models
314 
315  // Update motion fields
317 
318  // Set the optional reference phase fraction from the other phases
320  {
321  phaseModel* referencePhasePtr = &phases()[referencePhaseName_];
322  volScalarField& referenceAlpha = *referencePhasePtr;
323 
324  referenceAlpha = 1;
325 
326  forAll(phaseModels_, phasei)
327  {
328  if (&phaseModels_[phasei] != referencePhasePtr)
329  {
330  referenceAlpha -= phaseModels_[phasei];
331  }
332  }
333  }
334 
335  forAll(phases(), phasei)
336  {
337  const volScalarField& alphai = phases()[phasei];
338  mesh_.schemes().setFluxRequired(alphai.name());
339  }
340 }
341 
342 
343 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
344 
346 {}
347 
348 
349 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
350 
352 {
353  tmp<volScalarField> rho(movingPhaseModels_[0]*movingPhaseModels_[0].rho());
354 
355  for
356  (
357  label movingPhasei=1;
358  movingPhasei<movingPhaseModels_.size();
359  movingPhasei++
360  )
361  {
362  rho.ref() +=
363  movingPhaseModels_[movingPhasei]
364  *movingPhaseModels_[movingPhasei].rho();
365  }
366 
367  if (stationaryPhaseModels_.empty())
368  {
369  return rho;
370  }
371  else
372  {
373  return rho/sumAlphaMoving();
374  }
375 }
376 
377 
379 {
380  tmp<volVectorField> U(movingPhaseModels_[0]*movingPhaseModels_[0].U());
381 
382  for
383  (
384  label movingPhasei=1;
385  movingPhasei<movingPhaseModels_.size();
386  movingPhasei++
387  )
388  {
389  U.ref() +=
390  movingPhaseModels_[movingPhasei]
391  *movingPhaseModels_[movingPhasei].U();
392  }
393 
394  if (stationaryPhaseModels_.empty())
395  {
396  return U;
397  }
398  else
399  {
400  return U/sumAlphaMoving();
401  }
402 }
403 
404 
407 {
408  if (interfaceSurfaceTensionModels_.found(key))
409  {
410  return interfaceSurfaceTensionModels_[key]->sigma();
411  }
412  else
413  {
414  return volScalarField::New
415  (
416  interfaceSurfaceTensionModel::typeName + ":sigma",
417  mesh_,
419  );
420  }
421 }
422 
423 
426 {
427  if (interfaceSurfaceTensionModels_.found(key))
428  {
429  return interfaceSurfaceTensionModels_[key]->sigma(patchi);
430  }
431  else
432  {
433  return tmp<scalarField>
434  (
435  new scalarField(mesh_.boundary()[patchi].size(), 0)
436  );
437  }
438 }
439 
440 
443 {
444  tmp<volScalarField> tnearInt
445  (
447  (
448  "nearInterface",
449  mesh_,
451  )
452  );
453 
454  forAll(phases(), phasei)
455  {
456  tnearInt.ref() = max
457  (
458  tnearInt(),
459  pos0(phases()[phasei] - 0.01)*pos0(0.99 - phases()[phasei])
460  );
461  }
462 
463  return tnearInt;
464 }
465 
466 
468 (
469  const phaseInterfaceKey& key
470 ) const
471 {
472  return volScalarField::New
473  (
474  IOobject::groupName("dmdtf", phaseInterface(*this, key).name()),
475  mesh(),
477  );
478 }
479 
480 
482 {
483  return PtrList<volScalarField>(phaseModels_.size());
484 }
485 
486 
488 {
489  return PtrList<volScalarField>(phaseModels_.size());
490 }
491 
492 
494 {
495  forAll(phaseModels_, phasei)
496  {
497  if (!phaseModels_[phasei].incompressible())
498  {
499  return false;
500  }
501  }
502 
503  return true;
504 }
505 
506 
508 {
509  return false;
510 }
511 
512 
514 {
515  return false;
516 }
517 
518 
520 (
521  const phaseModel& phase1
522 ) const
523 {
524  tmp<surfaceScalarField> tSurfaceTension
525  (
527  (
528  "surfaceTension",
529  mesh_,
530  dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), 0)
531  )
532  );
533 
534  forAll(phases(), phasej)
535  {
536  const phaseModel& phase2 = phases()[phasej];
537 
538  if (&phase2 != &phase1)
539  {
540  const phaseInterface interface(phase1, phase2);
541 
542  if (cAlphas_.found(interface))
543  {
544  tSurfaceTension.ref() +=
545  fvc::interpolate(sigma(interface)*K(phase1, phase2))
546  *(
547  fvc::interpolate(phase2)*fvc::snGrad(phase1)
548  - fvc::interpolate(phase1)*fvc::snGrad(phase2)
549  );
550  }
551  }
552  }
553 
554  return tSurfaceTension;
555 }
556 
557 
559 {
560  forAll(phaseModels_, phasei)
561  {
562  phaseModels_[phasei].correct();
563  }
564 }
565 
566 
568 {
569  const PtrList<volScalarField> dmdts = this->dmdts();
570 
571  forAll(movingPhaseModels_, movingPhasei)
572  {
573  phaseModel& phase = movingPhaseModels_[movingPhasei];
574  const volScalarField& alpha = phase;
575  volScalarField& rho = phase.rho();
576 
577  volScalarField source
578  (
580  (
581  IOobject::groupName("source", phase.name()),
582  mesh_,
584  )
585  );
586 
587  if (fvModels().addsSupToField(rho.name()))
588  {
589  source += fvModels().source(alpha, rho)&rho;
590  }
591 
592  if (dmdts.set(phase.index()))
593  {
594  source += dmdts[phase.index()];
595  }
596 
597  phase.correctContinuityError(source);
598  }
599 }
600 
601 
603 {
604  bool updateDpdt = false;
605 
606  forAll(phaseModels_, phasei)
607  {
608  phaseModels_[phasei].correctKinematics();
609 
610  updateDpdt = updateDpdt || phaseModels_[phasei].thermo().dpdt();
611  }
612 
613  // Update the pressure time-derivative if required
614  if (updateDpdt)
615  {
616  dpdt_ = fvc::ddt(phaseModels_.begin()().thermo().p());
617  }
618 }
619 
620 
622 {
623  forAll(phaseModels_, phasei)
624  {
625  phaseModels_[phasei].correctThermo();
626  }
627 }
628 
629 
631 {
632  forAll(phaseModels_, phasei)
633  {
634  phaseModels_[phasei].correctReactions();
635  }
636 }
637 
638 
640 {
641  forAll(phaseModels_, phasei)
642  {
643  phaseModels_[phasei].correctSpecies();
644  }
645 }
646 
647 
649 {
650  forAll(phaseModels_, phasei)
651  {
652  phaseModels_[phasei].predictMomentumTransport();
653  }
654 }
655 
656 
658 {
659  forAll(phaseModels_, phasei)
660  {
661  phaseModels_[phasei].predictThermophysicalTransport();
662  }
663 }
664 
665 
667 {
668  forAll(phaseModels_, phasei)
669  {
670  phaseModels_[phasei].correctMomentumTransport();
671  }
672 }
673 
674 
676 {
677  forAll(phaseModels_, phasei)
678  {
679  phaseModels_[phasei].correctThermophysicalTransport();
680  }
681 }
682 
683 
685 {
686  if (mesh_.changing())
687  {
688  MRF_.update();
689 
690  // forAll(phaseModels_, phasei)
691  // {
692  // phaseModels_[phasei].meshUpdate();
693  // }
694  }
695 }
696 
697 
699 {
700  forAll(movingPhases(), movingPhasei)
701  {
702  phaseModel& phase = movingPhases()[movingPhasei];
703 
704  tmp<volVectorField> tU(phase.U());
705  const volVectorField::Boundary& UBf = tU().boundaryField();
706 
708  (
709  MRF_.relative(mesh_.Sf().boundaryField() & UBf)
710  );
711 
713 
714  forAll(mesh_.boundary(), patchi)
715  {
716  if
717  (
718  isA<fixedValueFvsPatchScalarField>(phiBf[patchi])
719  && !isA<movingWallVelocityFvPatchVectorField>(UBf[patchi])
720  )
721  {
722  phiBf[patchi] == phiRelBf[patchi];
723  }
724  }
725  }
726 }
727 
728 
730 (
731  const volScalarField& p_rgh,
732  const autoPtr<volScalarField>& divU,
735 )
736 {
737  forAll(movingPhases(), movingPhasei)
738  {
739  phaseModel& phase = movingPhases()[movingPhasei];
740 
743 
744  forAll(Ubf, patchi)
745  {
746  if (Ubf[patchi].fixesValue())
747  {
748  Ubf[patchi].initEvaluate();
749  }
750  }
751 
752  forAll(Ubf, patchi)
753  {
754  if (Ubf[patchi].fixesValue())
755  {
756  Ubf[patchi].evaluate();
757  UfBf[patchi] = Ubf[patchi];
758  }
759  }
760  }
761 
762  // Correct fixed-flux BCs to be consistent with the velocity BCs
763  correctBoundaryFlux();
764 
765  phi_ = Zero;
766  PtrList<surfaceScalarField> alphafs(phaseModels_.size());
767  forAll(movingPhases(), movingPhasei)
768  {
769  phaseModel& phase = movingPhases()[movingPhasei];
770  const label phasei = phase.index();
771  const volScalarField& alpha = phase;
772 
773  alphafs.set(phasei, fvc::interpolate(alpha).ptr());
774 
775  // Calculate absolute flux
776  // from the mapped surface velocity
777  phi_ += alphafs[phasei]*(mesh_.Sf() & phase.UfRef());
778  }
779 
780  if (incompressible())
781  {
783  (
784  phi_,
785  movingPhases()[0].U(),
786  p_rgh,
788  divU,
790  pimple
791  );
792  }
793  else
794  {
796  (
798  (
799  "psi",
800  mesh_,
802  )
803  );
804 
805  forAll(phases(), phasei)
806  {
807  phaseModel& phase = phases()[phasei];
808  const volScalarField& alpha = phase;
809 
810  psi += alpha*phase.thermo().psi()/phase.rho();
811  }
812 
814  (
815  phi_,
816  p_rgh,
817  psi,
819  divU(),
820  pimple
821  );
822  }
823 
824  // Make the flux relative to the mesh motion
825  fvc::makeRelative(phi_, movingPhases()[0].U());
826 
827  setMixturePhi(alphafs, phi_);
828 }
829 
830 
832 {
833  if (regIOobject::read())
834  {
835  bool readOK = true;
836 
837  forAll(phaseModels_, phasei)
838  {
839  readOK &= phaseModels_[phasei].read();
840  }
841 
842  // models ...
843 
844  return readOK;
845  }
846  else
847  {
848  return false;
849  }
850 }
851 
852 
854 {
856  {
857  return fv::localEulerDdt::localRDeltaT(vf.mesh())*vf;
858  }
859  else
860  {
861  return vf/vf.mesh().time().deltaT();
862  }
863 }
864 
865 
867 {
869  {
871  }
872  else
873  {
874  return sf/sf.mesh().time().deltaT();
875  }
876 }
877 
878 
879 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const Mesh & mesh() const
Return mesh.
Generic field type.
Definition: FieldField.H:77
Generic GeometricBoundaryField class.
void evaluate()
Evaluate boundary conditions.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
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:370
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
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:65
bool set(const label) const
Is element set.
Definition: UPtrListI.H:78
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
Dimension set for the base types.
Definition: dimensionSet.H:122
virtual const volScalarField & psi() const =0
Compressibility [s^2/m^2].
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1748
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:503
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
static const dimensionSet dimSigma
Coefficient dimensions.
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 const rhoThermo & thermo() const =0
Return the thermophysical model.
virtual bool stationary() const =0
Return whether the phase is stationary.
label index() const
Return the index of the phase.
Definition: phaseModel.C:121
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:109
virtual const volScalarField & rho() const =0
Return the density field.
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
Definition: phaseModel.C:163
virtual surfaceVectorField & UfRef()=0
Access the face velocity.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
phaseModelPartialList stationaryPhaseModels_
Stationary phase models.
Definition: phaseSystem.H:147
virtual void correctKinematics()
Correct the kinematics.
Definition: phaseSystem.C:602
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
Definition: phaseSystem.C:442
void correctPhi(const volScalarField &p_rgh, const autoPtr< volScalarField > &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pimple)
Definition: phaseSystem.C:730
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:238
void generateInterfacialModels(const dictionary &dict, const phaseInterface &interface, PtrList< phaseInterface > &interfaces, PtrList< ModelType > &models) const
Generate interfacial-model lists.
virtual void correctReactions()
Correct the reactions.
Definition: phaseSystem.C:630
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:128
phaseModelPartialList multicomponentPhaseModels_
Multi-component phase models.
Definition: phaseSystem.H:153
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:141
phaseModelPartialList movingPhaseModels_
Moving phase models.
Definition: phaseSystem.H:144
void setMixtureU(const volVectorField &Um)
Re-normalise the velocity of the phases.
Definition: phaseSystem.C:104
virtual void correctSpecies()
Correct the species mass fractions.
Definition: phaseSystem.C:639
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: phaseSystemI.H:133
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseSystem.C:621
virtual void correct()
Correct the fluid properties other than those listed below.
Definition: phaseSystem.C:558
tmp< volScalarField > sigma(const phaseInterfaceKey &key) const
Return the surface tension coefficient for an interface.
Definition: phaseSystem.C:406
phaseModelPartialList anisothermalPhaseModels_
Anisothermal phase models.
Definition: phaseSystem.H:150
virtual void predictThermophysicalTransport()
Predict the energy transport e.g. alphat.
Definition: phaseSystem.C:657
tmp< volScalarField > K(const phaseModel &alpha1, const phaseModel &alpha2) const
Curvature of interface between two phases.
Definition: phaseSystem.C:186
tmp< volVectorField > U() const
Return the mixture velocity.
Definition: phaseSystem.C:378
virtual void correctMomentumTransport()
Correct the momentumTransport.
Definition: phaseSystem.C:666
void correctBoundaryFlux()
Correct fixed-flux BCs to be consistent with the velocity BCs.
Definition: phaseSystem.C:698
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
Definition: phaseSystem.C:55
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:41
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
Definition: phaseSystem.C:468
interfaceSurfaceTensionModelTable interfaceSurfaceTensionModels_
Surface tension models.
Definition: phaseSystem.H:171
virtual void correctThermophysicalTransport()
Correct the energy transport e.g. alphat.
Definition: phaseSystem.C:675
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
Definition: phaseSystem.C:481
virtual void meshUpdate()
Update the fluid properties for mesh changes.
Definition: phaseSystem.C:684
cAlphaTable cAlphas_
Interface compression coefficients.
Definition: phaseSystem.H:162
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:156
virtual bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
Definition: phaseSystem.C:513
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
Return the surface tension force.
Definition: phaseSystem.C:520
tmp< volScalarField > sumAlphaMoving() const
Return the sum of the phase fractions of the moving phases.
Definition: phaseSystem.C:78
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Definition: phaseSystem.C:210
void setMixturePhi(const PtrList< surfaceScalarField > &alphafs, const surfaceScalarField &phim)
Re-normalise the flux of the phases.
Definition: phaseSystem.C:125
virtual ~phaseSystem()
Destructor.
Definition: phaseSystem.C:345
void generateInterfacialValues(const dictionary &dict, HashTable< ValueType, phaseInterfaceKey, phaseInterfaceKey::hash > &values) const
Generate interfacial-model tables.
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases dotted with face areas.
Definition: phaseSystem.C:175
tmp< volScalarField > rho() const
Return the mixture density.
Definition: phaseSystem.C:351
word referencePhaseName_
Name of optional reference phase which is not solved for.
Definition: phaseSystem.H:138
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases.
Definition: phaseSystem.C:149
virtual void predictMomentumTransport()
Predict the momentumTransport.
Definition: phaseSystem.C:648
virtual PtrList< volScalarField > d2mdtdps() const
Return the mass transfer pressure implicit coefficients.
Definition: phaseSystem.C:487
virtual bool read()
Read base phaseProperties dictionary.
Definition: phaseSystem.C:831
virtual void correctContinuityError()
Correct the continuity errors.
Definition: phaseSystem.C:567
bool incompressible() const
Return incompressibility.
Definition: phaseSystem.C:493
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.
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:181
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
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))
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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:105
dimensionedScalar pos0(const dimensionedScalar &ds)
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:853
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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.
const dimensionSet dimDensity
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47