phaseSystem.H
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-2021 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 Class
25  Foam::phaseSystem
26 
27 Description
28  Class to represent a system of phases and model interfacial transfers
29  between them.
30 
31 SourceFiles
32  phaseSystem.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef phaseSystem_H
37 #define phaseSystem_H
38 
39 #include "IOdictionary.H"
40 
41 #include "phaseModel.H"
42 #include "phasePair.H"
43 #include "orderedPhasePair.H"
44 #include "HashPtrTable.H"
45 #include "PtrListDictionary.H"
46 
47 #include "IOMRFZoneList.H"
48 #include "fvModels.H"
49 #include "fvConstraints.H"
50 
51 #include "volFields.H"
52 #include "surfaceFields.H"
53 #include "fvMatricesFwd.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 
60 class blendingMethod;
61 template<class modelType> class BlendedInterfacialModel;
62 class surfaceTensionModel;
63 class aspectRatioModel;
64 class pressureReference;
65 class nonOrthogonalSolutionControl;
66 
67 /*---------------------------------------------------------------------------*\
68  Class phaseSystem Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 class phaseSystem
72 :
73  public IOdictionary
74 {
75 public:
76 
77  // Public Typedefs
78 
80 
82 
84 
86 
88 
89  typedef
90  HashTable<autoPtr<phasePair>, phasePairKey, phasePairKey::hash>
92 
93  typedef
95  dmdtfTable;
96 
97  typedef
99  <
100  HashPtrTable<volScalarField>,
101  phasePairKey,
103  >
104  dmidtfTable;
105 
106 
107 
108 protected:
109 
110  // Protected typedefs
111 
112  typedef
114  dictTable;
115 
116  typedef
117  HashTable<autoPtr<blendingMethod>, word, word::hash>
119 
120  typedef
121  HashTable
122  <
123  autoPtr<surfaceTensionModel>,
124  phasePairKey,
126  >
128 
129  typedef
130  HashTable
131  <
132  autoPtr<aspectRatioModel>,
133  phasePairKey,
135  >
137 
139  cAlphaTable;
140 
141 
142  // Protected data
143 
144  //- Reference to the mesh
145  const fvMesh& mesh_;
146 
147  //- Name of optional reference phase which is not solved for
148  // but obtained from the sum of the other phases
149  word referencePhaseName_;
150 
151  //- Phase models
152  phaseModelList phaseModels_;
153 
154  //- Moving phase models
155  phaseModelPartialList movingPhaseModels_;
156 
157  //- Stationary phase models
158  phaseModelPartialList stationaryPhaseModels_;
159 
160  //- Anisothermal phase models
161  phaseModelPartialList anisothermalPhaseModels_;
162 
163  //- Multi-component phase models
164  phaseModelPartialList multiComponentPhaseModels_;
165 
166  //- Phase pairs
167  phasePairTable phasePairs_;
168 
169  //- Total volumetric flux
171 
172  //- Rate of change of pressure
174 
175  //- Optional MRF zones
177 
178  //- Blending methods
179  blendingMethodTable blendingMethods_;
180 
181  //- Interface compression coefficients
182  cAlphaTable cAlphas_;
183 
184  //- Stabilisation for normalisation of the interface normal
186 
187 
188  // Sub Models
189 
190  //- Surface tension models
192 
193  //- Aspect ratio models
195 
196 
197  //- Flag to indicate that returned lists of fields are "complete"; i.e.,
198  // that an absence of force is returned as a constructed list of zeros,
199  // rather than a null pointer
200  static const bool fillFields_ = false;
201 
202 
203  // Protected member functions
204 
205  //- Calculate and return the mixture flux
207  (
208  const phaseModelList& phaseModels
209  ) const;
210 
211  //- Generate pairs
212  void generatePairs(const dictTable& modelDicts);
213 
214  //- Generate pairs and sub-model tables
215  template<class modelType>
216  void createSubModels
217  (
218  const dictTable& modelDicts,
219  HashTable
220  <
222  phasePairKey,
224  >& models
225  );
226 
227  //- Generate pairs and sub-model tables
228  template<class modelType>
230  (
231  const word& modelName,
232  HashTable
233  <
235  phasePairKey,
237  >& models,
238  const bool correctFixedFluxBCs = true
239  );
240 
241  //- Generate pairs and blended sub-model tables
242  template<class modelType>
244  (
245  const word& modelName,
246  HashTable
247  <
249  phasePairKey,
251  >& models,
252  const bool correctFixedFluxBCs = true
253  );
254 
255  //- Generate pairs and two-sided sub-model tables
256  template<class modelType>
258  (
259  const word& modelName,
260  HashTable
261  <
263  phasePairKey,
265  >& models,
266  const bool correctFixedFluxBCs = true
267  );
268 
269  //- Check that mass transfer is supported across the given interface
270  template<class modelType>
271  void validateMassTransfer(const phasePair& key) const;
272 
273  //- Return the sum of the phase fractions of the moving phases
275 
276  //- Re-normalise the velocity of the phases
277  // around the specified mixture mean
278  void setMixtureU(const volVectorField& Um);
279 
280 
281  // Functions required for interface compression
282 
283  //- Normal to interface between two phases
284  // Used for interface compression
286  (
287  const volScalarField& alpha1,
288  const volScalarField& alpha2
289  ) const;
290 
291  //- Normal to interface between two phases dotted with face areas
292  // Used for interface compression
294  (
295  const volScalarField& alpha1,
296  const volScalarField& alpha2
297  ) const;
298 
299  //- Correction for the boundary condition on the unit normal nHat on
300  // walls to produce the correct contact angle.
301  //
302  // The dynamic contact angle is calculated from the component of
303  // the velocity on the direction of the interface, parallel to the
304  // wall. Used for interface compression.
306  (
307  const phaseModel& alpha1,
308  const phaseModel& alpha2,
309  surfaceVectorField::Boundary& nHatb
310  ) const;
311 
312  //- Curvature of interface between two phases
313  // Used for interface compression
315  (
316  const phaseModel& alpha1,
317  const phaseModel& alpha2
318  ) const;
319 
320 
321  // Functions required by twoPhaseSystem
322 
323  //- Return the drag coefficient for phase pair
324  virtual tmp<volScalarField> Kd(const phasePairKey& key) const = 0;
325 
326  //- Return the virtual mass coefficient for phase pair
327  virtual tmp<volScalarField> Vm(const phasePairKey& key) const = 0;
328 
329  //- Return the face drag coefficient for phase pair
331  (
332  const phasePairKey& key
333  ) const = 0;
334 
335 
336 public:
337 
338  //- Runtime type information
339  TypeName("phaseSystem");
340 
341  //- Default name of the phase properties dictionary
342  static const word propertiesName;
343 
344 
345  // Declare runtime construction
346 
348  (
349  autoPtr,
350  phaseSystem,
351  dictionary,
352  (
353  const fvMesh& mesh
354  ),
355  (mesh)
356  );
357 
358 
359  // Constructors
360 
361  //- Construct from fvMesh
362  phaseSystem(const fvMesh& mesh);
363 
364 
365  // Selectors
366 
368  (
369  const fvMesh& mesh
370  );
371 
372 
373  //- Destructor
374  virtual ~phaseSystem();
375 
376 
377  // Member Functions
378 
379  // Access
380 
381  //- Return the mesh
382  inline const fvMesh& mesh() const;
383 
384  //- Return the phase models
385  inline const phaseModelList& phases() const;
386 
387  //- Access the phase models
388  inline phaseModelList& phases();
389 
390  //- Return the models for phases that are moving
391  inline const phaseModelPartialList& movingPhases() const;
392 
393  //- Access the models for phases that are moving
394  inline phaseModelPartialList& movingPhases();
395 
396  //- Return the models for phases that are stationary
397  inline const phaseModelPartialList& stationaryPhases() const;
398 
399  //- Access the models for phases that are stationary
400  inline phaseModelPartialList& stationaryPhases();
401 
402  //- Return the models for phases that have variable temperature
403  inline const phaseModelPartialList& anisothermalPhases() const;
404 
405  //- Access the models for phases that have variable temperature
406  inline phaseModelPartialList& anisothermalPhases();
407 
408  //- Return the models for phases that have multiple species
409  inline const phaseModelPartialList& multiComponentPhases() const;
410 
411  //- Access the models for phases that have multiple species
412  inline phaseModelPartialList& multiComponentPhases();
413 
414  //- Return the phase pairs
415  inline const phasePairTable& phasePairs() const;
416 
417  //- Return the phase not given as an argument in a two-phase system
418  // An error is generated if the system is not two-phase
419  inline const phaseModel& otherPhase(const phaseModel& phase) const;
420 
421  //- Return the mixture flux
422  inline const surfaceScalarField& phi() const;
423 
424  //- Access the mixture flux
425  inline surfaceScalarField& phi();
426 
427  //- Return the rate of change of the pressure
428  inline const volScalarField& dpdt() const;
429 
430  //- Access the rate of change of the pressure
431  inline volScalarField& dpdt();
432 
433  //- Return MRF zones
434  inline const IOMRFZoneList& MRF() const;
435 
436  //- Access the fvModels
437  inline Foam::fvModels& fvModels(fvMesh& mesh);
438 
439  //- Access the fvModels
440  inline const Foam::fvModels& fvModels() const;
441 
442  //- Access the fvConstraints
444 
445  //- Access the fvConstraints
446  inline const Foam::fvConstraints& fvConstraints() const;
447 
448 
449  // Sub-model lookup
450 
451  //- Check availability of a sub model for a given phase pair
452  template<class modelType>
453  bool foundSubModel(const phasePair& key) const;
454 
455  //- Return a sub model between a phase pair
456  template<class modelType>
457  const modelType& lookupSubModel(const phasePair& key) const;
458 
459  //- Check availability of a sub model between two phases
460  template<class modelType>
461  bool foundSubModel
462  (
463  const phaseModel& dispersed,
464  const phaseModel& continuous
465  ) const;
466 
467  //- Return a sub model between two phases
468  template<class modelType>
469  const modelType& lookupSubModel
470  (
471  const phaseModel& dispersed,
472  const phaseModel& continuous
473  ) const;
474 
475  //- Check availability of a blended sub model for a given phase pair
476  template<class modelType>
477  bool foundBlendedSubModel(const phasePair& key) const;
478 
479  //- Return a blended sub model between a phase pair
480  template<class modelType>
482  lookupBlendedSubModel(const phasePair& key) const;
483 
484 
485  // Field construction
486 
487  //- Fill up gaps in a phase-indexed list of fields with zeros
488  template
489  <
490  class Type,
491  template<class> class PatchField,
492  class GeoMesh
493  >
494  void fillFields
495  (
496  const word& name,
497  const dimensionSet& dims,
499  ) const;
500 
501  //- Fill up gaps in a phase-indexed table of fields with zeros
502  template
503  <
504  class Type,
505  template<class> class PatchField,
506  class GeoMesh
507  >
508  void fillFields
509  (
510  const word& name,
511  const dimensionSet& dims,
513  fieldTable
514  ) const;
515 
516 
517  // Properties
518 
519  //- Return the mixture density
520  tmp<volScalarField> rho() const;
521 
522  //- Return the mixture velocity
523  tmp<volVectorField> U() const;
524 
525  //- Return the aspect-ratio for a pair
526  tmp<volScalarField> E(const phasePairKey& key) const;
527 
528  //- Return the surface tension coefficient for a pair
529  tmp<volScalarField> sigma(const phasePairKey& key) const;
530 
531  //- Return the surface tension coefficient for a pair on a patch
533  (
534  const phasePairKey& key,
535  label patchi
536  ) const;
537 
538  //- Indicator of the proximity of the interface
539  // Field values are 1 near and 0 away for the interface.
541 
542  //- Return the mass transfer rate for an interface
543  virtual tmp<volScalarField> dmdtf(const phasePairKey& key) const;
544 
545  //- Return the mass transfer rates for each phase
546  virtual PtrList<volScalarField> dmdts() const;
547 
548  //- Return the mass transfer pressure implicit coefficients
549  // for each phase
550  virtual PtrList<volScalarField> d2mdtdps() const;
551 
552  //- Return incompressibility
553  bool incompressible() const;
554 
555 
556  // Transfers
557 
558  //- Return the momentum transfer matrices for the cell-based
559  // algorithm
561 
562  //- Return the momentum transfer matrices for the face-based
563  // algorithm
565 
566  //- Return the implicit force coefficients for the face-based
567  // algorithm
568  virtual PtrList<surfaceScalarField> AFfs() const = 0;
569 
570  //- Return the force fluxes for the cell-based algorithm
572  (
574  ) = 0;
575 
576  //- Return the force fluxes for the face-based algorithm
578  (
580  ) = 0;
581 
582  //- Return the force fluxes for the cell-based algorithm
584  (
586  ) const = 0;
587 
588  //- Return the force fluxes for the face-based algorithm
590  (
592  ) const = 0;
593 
594  //- Return the explicit part of the drag force
596  (
598  ) const = 0;
599 
600  //- Returns true if the phase pressure is treated implicitly
601  // in the phase fraction equation
602  virtual bool implicitPhasePressure(const phaseModel& phase) const;
603 
604  //- Returns true if the phase pressure is treated implicitly
605  // in the phase fraction equation for any phase
606  virtual bool implicitPhasePressure() const;
607 
608  //- Return the phase diffusivity
609  // divided by the momentum central coefficient
611  (
614  ) const = 0;
615 
616  //- Solve the drag system for the new velocities and fluxes
617  virtual void partialElimination
618  (
623  ) = 0;
624 
625  //- Solve the drag system for the new fluxes
626  virtual void partialEliminationf
627  (
631  ) = 0;
632 
633  //- Re-normalise the flux of the phases
634  // around the specified mixture mean
635  void setMixturePhi
636  (
638  const surfaceScalarField& phim
639  );
640 
641  //- Return the flux corrections for the cell-based algorithm
643  (
645  const bool includeVirtualMass = false
646  ) const = 0;
647 
648  //- Return the heat transfer matrices
649  virtual autoPtr<heatTransferTable> heatTransfer() const = 0;
650 
651  //- Return the specie transfer matrices
652  virtual autoPtr<specieTransferTable> specieTransfer() const = 0;
653 
654  //- Return the surface tension force
656  (
657  const phaseModel& phase
658  ) const;
659 
660 
661  // Evolution
662 
663  //- Solve for the phase fractions
664  virtual void solve
665  (
668  );
669 
670  //- Correct the fluid properties other than those listed below
671  virtual void correct();
672 
673  //- Correct the continuity errors
674  virtual void correctContinuityError();
675 
676  //- Correct the kinematics
677  virtual void correctKinematics();
678 
679  //- Correct the thermodynamics
680  virtual void correctThermo();
681 
682  //- Correct the reactions
683  virtual void correctReactions();
684 
685  //- Correct the species mass fractions
686  virtual void correctSpecies();
687 
688  //- Correct the turbulence
689  virtual void correctTurbulence();
690 
691  //- Correct the energy transport e.g. alphat
692  virtual void correctEnergyTransport();
693 
694  //- Update the fluid properties for mesh changes
695  virtual void meshUpdate();
696 
697  //- Correct fixed-flux BCs to be consistent with the velocity BCs
698  void correctBoundaryFlux();
699 
700  void correctPhi
701  (
702  const volScalarField& p_rgh,
703  const tmp<volScalarField>& divU,
706  );
707 
708 
709  // IO
710 
711  //- Read base phaseProperties dictionary
712  virtual bool read();
713 };
714 
715 
716 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
717 
720 
721 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
722 
723 } // End namespace Foam
724 
725 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
726 
727 #include "phaseSystemI.H"
728 
729 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
730 
731 #ifdef NoRepository
732  #include "phaseSystemTemplates.C"
733 #endif
734 
735 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
736 
737 #endif
738 
739 // ************************************************************************* //
Foam::surfaceFields.
const BlendedInterfacialModel< modelType > & lookupBlendedSubModel(const phasePair &key) const
Return a blended sub model between a phase pair.
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:341
Hashing function class, shared by all the derived classes.
Definition: string.H:92
pressureReference & pressureReference
volScalarField & p_rgh
virtual PtrList< surfaceScalarField > AFfs() const =0
Return the implicit force coefficients for the face-based.
void generatePairs(const dictTable &modelDicts)
Generate pairs.
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &phiKdPhifs)=0
Solve the drag system for the new fluxes.
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:49
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
Definition: phaseSystem.H:190
tmp< volScalarField > E(const phasePairKey &key) const
Return the aspect-ratio for a pair.
virtual void correctTurbulence()
Correct the turbulence.
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
virtual tmp< volScalarField > Kd(const phasePairKey &key) const =0
Return the drag coefficient for phase pair.
pimpleNoLoopControl & pimple
const Foam::fvModels & fvModels() const
Access the fvModels.
Definition: phaseSystemI.H:170
virtual PtrList< surfaceScalarField > phiFfs(const PtrList< surfaceScalarField > &rAUfs)=0
Return the force fluxes for the face-based algorithm.
virtual autoPtr< specieTransferTable > specieTransfer() const =0
Return the specie transfer matrices.
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases dotted with face areas.
virtual PtrList< volScalarField > d2mdtdps() const
Return the mass transfer pressure implicit coefficients.
const Foam::fvConstraints & fvConstraints() const
Access the fvConstraints.
Definition: phaseSystemI.H:182
declareRunTimeSelectionTable(autoPtr, phaseSystem, dictionary,(const fvMesh &mesh),(mesh))
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition: phaseSystem.H:90
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const dimensionedScalar deltaN_
Stabilisation for normalisation of the interface normal.
Definition: phaseSystem.H:184
PtrList< surfaceScalarField > alphafs(phases.size())
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);if(fluid.found("pMin")){ IOWarningInFunction(fluid)<< "Pressure limits, pMin and pMax, are now read from "<< pimple.dict().name()<< endl;}pressureReference pressureReference(p, p_rgh, pimple.dict(), fluid.incompressible());if(fluid.incompressible()){ p=p_rgh+fluid.rho() *gh;}if(p_rgh.needReference() &&fluid.incompressible()){ p+=dimensionedScalar("p", p.dimensions(), pressureReference.refValue() - getRefCellValue(p, pressureReference.refCell()));}p_rgh=p - fluid.rho() *gh;mesh.setFluxRequired(p_rgh.name());PtrList< volScalarField > rAUs
Definition: createFields.H:67
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:78
void setMixturePhi(const PtrList< surfaceScalarField > &alphafs, const surfaceScalarField &phim)
Re-normalise the flux of the phases.
phasePairTable phasePairs_
Phase pairs.
Definition: phaseSystem.H:166
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
bool foundBlendedSubModel(const phasePair &key) const
Check availability of a blended sub model for a given phase pair.
virtual autoPtr< heatTransferTable > heatTransfer() const =0
Return the heat transfer matrices.
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
volScalarField dpdt_
Rate of change of pressure.
Definition: phaseSystem.H:172
volScalarField & alpha1(mixture.alpha1())
TypeName("phaseSystem")
Runtime type information.
virtual autoPtr< momentumTransferTable > momentumTransferf()=0
Return the momentum transfer matrices for the face-based.
tmp< volScalarField > K(const phaseModel &alpha1, const phaseModel &alpha2) const
Curvature of interface between two phases.
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:80
phaseModelPartialList stationaryPhaseModels_
Stationary phase models.
Definition: phaseSystem.H:157
void correctPhi(const volScalarField &p_rgh, const tmp< volScalarField > &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pimple)
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:82
virtual PtrList< surfaceScalarField > ddtCorrByAs(const PtrList< volScalarField > &rAUs, const bool includeVirtualMass=false) const =0
Return the flux corrections for the cell-based algorithm.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:50
Provides controls for the pressure reference in closed-volume simulations.
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
HashTable< autoPtr< surfaceTensionModel >, phasePairKey, phasePairKey::hash > surfaceTensionModelTable
Definition: phaseSystem.H:126
virtual autoPtr< momentumTransferTable > momentumTransfer()=0
Return the momentum transfer matrices for the cell-based.
virtual tmp< volScalarField > Vm(const phasePairKey &key) const =0
Return the virtual mass coefficient for phase pair.
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
void correctBoundaryFlux()
Correct fixed-flux BCs to be consistent with the velocity BCs.
HashPtrTable< HashPtrTable< volScalarField >, phasePairKey, phasePairKey::hash > dmidtfTable
Definition: phaseSystem.H:103
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:151
alpha2
Definition: alphaEqn.H:115
phaseModelPartialList multiComponentPhaseModels_
Multi-component phase models.
Definition: phaseSystem.H:163
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const phaseModelPartialList & stationaryPhases() const
Return the models for phases that are stationary.
Definition: phaseSystemI.H:63
virtual void solve(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs)
Solve for the phase fractions.
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
Return the surface tension force.
const phasePairTable & phasePairs() const
Return the phase pairs.
Definition: phaseSystemI.H:105
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:35
virtual void correctContinuityError()
Correct the continuity errors.
virtual PtrList< volVectorField > KdUByAs(const PtrList< volScalarField > &rAUs) const =0
Return the explicit part of the drag force.
virtual void correctKinematics()
Correct the kinematics.
Dimension set for the base types.
Definition: dimensionSet.H:120
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:70
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdtfTable
Definition: phaseSystem.H:94
tmp< volVectorField > U() const
Return the mixture velocity.
Finite volume models.
Definition: fvModels.H:60
const IOMRFZoneList & MRF() const
Return MRF zones.
Definition: phaseSystemI.H:158
const phaseModelPartialList & anisothermalPhases() const
Return the models for phases that have variable temperature.
Definition: phaseSystemI.H:77
virtual void correctThermo()
Correct the thermodynamics.
PtrListDictionary< phaseModel > phaseModelList
Definition: phaseSystem.H:84
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:112
word referencePhaseName_
Name of optional reference phase which is not solved for.
Definition: phaseSystem.H:148
UPtrList< phaseModel > phaseModelPartialList
Definition: phaseSystem.H:86
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:144
static autoPtr< phaseSystem > New(const fvMesh &mesh)
virtual PtrList< surfaceScalarField > phiFs(const PtrList< volScalarField > &rAUs)=0
Return the force fluxes for the cell-based algorithm.
bool foundSubModel(const phasePair &key) const
Check availability of a sub model for a given phase pair.
An STL-conforming hash table.
Definition: HashTable.H:61
void validateMassTransfer(const phasePair &key) const
Check that mass transfer is supported across the given interface.
phaseModelPartialList anisothermalPhaseModels_
Anisothermal phase models.
Definition: phaseSystem.H:160
cAlphaTable cAlphas_
Interface compression coefficients.
Definition: phaseSystem.H:181
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases.
phaseModelPartialList movingPhaseModels_
Moving phase models.
Definition: phaseSystem.H:154
virtual void correct()
Correct the fluid properties other than those listed below.
tmp< volScalarField > rho() const
Return the mixture density.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
tmp< volScalarField > sumAlphaMoving() const
Return the sum of the phase fractions of the moving phases.
aspectRatioModelTable aspectRatioModels_
Aspect ratio models.
Definition: phaseSystem.H:193
IOMRFZoneList MRF_
Optional MRF zones.
Definition: phaseSystem.H:175
PtrList< surfaceScalarField > rAUfs
Definition: createFields.H:68
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
volScalarField sf(fieldObject, mesh)
const volScalarField & dpdt() const
Return the rate of change of the pressure.
Definition: phaseSystemI.H:146
bool incompressible() const
Return incompressibility.
virtual tmp< surfaceScalarField > Kdf(const phasePairKey &key) const =0
Return the face drag coefficient for phase pair.
tmp< volScalarField > byDt(const volScalarField &vf)
virtual void correctReactions()
Correct the reactions.
virtual void partialElimination(const PtrList< volScalarField > &rAUs, const PtrList< volVectorField > &KdUByAs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &phiKdPhis)=0
Solve the drag system for the new velocities and fluxes.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
void fillFields(const word &name, const dimensionSet &dims, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fieldList) const
Fill up gaps in a phase-indexed list of fields with zeros.
Forward declarations of fvMatrix specialisations.
label patchi
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: phaseSystemI.H:134
virtual ~phaseSystem()
Destructor.
Finite volume constraints.
Definition: fvConstraints.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
virtual bool read()
Read base phaseProperties dictionary.
blendingMethodTable blendingMethods_
Blending methods.
Definition: phaseSystem.H:178
HashTable< autoPtr< blendingMethod >, word, word::hash > blendingMethodTable
Definition: phaseSystem.H:117
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models, const bool correctFixedFluxBCs=true)
Generate pairs and sub-model tables.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual PtrList< surfaceScalarField > DByAfs(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs) const =0
Return the phase diffusivity.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
HashTable< autoPtr< aspectRatioModel >, phasePairKey, phasePairKey::hash > aspectRatioModelTable
Definition: phaseSystem.H:135
HashTable< scalar, phasePairKey, phasePairKey::hash > cAlphaTable
Definition: phaseSystem.H:138
void setMixtureU(const volVectorField &Um)
Re-normalise the velocity of the phases.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
void correctContactAngle(const phaseModel &alpha1, const phaseModel &alpha2, surfaceVectorField::Boundary &nHatb) const
Correction for the boundary condition on the unit normal nHat on.
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:169
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
virtual PtrList< surfaceScalarField > phiKdPhifs(const PtrList< surfaceScalarField > &rAUfs) const =0
Return the force fluxes for the face-based algorithm.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual PtrList< surfaceScalarField > phiKdPhis(const PtrList< volScalarField > &rAUs) const =0
Return the force fluxes for the cell-based algorithm.
const phaseModelPartialList & multiComponentPhases() const
Return the models for phases that have multiple species.
Definition: phaseSystemI.H:91
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:28
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Non-orthogonal solution control class. Provides non-orthogonal-loop control methods.
virtual void meshUpdate()
Update the fluid properties for mesh changes.
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries...
Definition: IOMRFZoneList.H:64
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Definition: phaseSystem.H:113
virtual bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
virtual void correctSpecies()
Correct the species mass fractions.
Namespace for OpenFOAM.
virtual tmp< volScalarField > dmdtf(const phasePairKey &key) const
Return the mass transfer rate for an interface.