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