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-2020 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 "fvOptions.H"
49 
50 #include "volFields.H"
51 #include "surfaceFields.H"
52 #include "fvMatricesFwd.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 
59 class blendingMethod;
60 template<class modelType> class BlendedInterfacialModel;
61 class surfaceTensionModel;
62 class aspectRatioModel;
63 
64 /*---------------------------------------------------------------------------*\
65  Class phaseSystem Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class phaseSystem
69 :
70  public IOdictionary
71 {
72 public:
73 
74  // Public Typedefs
75 
77 
79 
81 
83 
85 
86  typedef
87  HashTable<autoPtr<phasePair>, phasePairKey, phasePairKey::hash>
89 
90  typedef
92  dmdtfTable;
93 
94  typedef
96  <
97  HashPtrTable<volScalarField>,
98  phasePairKey,
100  >
101  dmidtfTable;
102 
103 
104 
105 protected:
106 
107  // Protected typedefs
108 
109  typedef
111  dictTable;
112 
113  typedef
114  HashTable<autoPtr<blendingMethod>, word, word::hash>
116 
117  typedef
118  HashTable
119  <
120  autoPtr<surfaceTensionModel>,
121  phasePairKey,
123  >
125 
126  typedef
127  HashTable
128  <
129  autoPtr<aspectRatioModel>,
130  phasePairKey,
132  >
134 
136  cAlphaTable;
137 
138 
139  // Protected data
140 
141  //- Reference to the mesh
142  const fvMesh& mesh_;
143 
144  //- Name of optional reference phase which is not solved for
145  // but obtained from the sum of the other phases
146  word referencePhaseName_;
147 
148  //- Phase models
149  phaseModelList phaseModels_;
150 
151  //- Moving phase models
152  phaseModelPartialList movingPhaseModels_;
153 
154  //- Stationary phase models
155  phaseModelPartialList stationaryPhaseModels_;
156 
157  //- Anisothermal phase models
158  phaseModelPartialList anisothermalPhaseModels_;
159 
160  //- Multi-component phase models
161  phaseModelPartialList multiComponentPhaseModels_;
162 
163  //- Phase pairs
164  phasePairTable phasePairs_;
165 
166  //- Total volumetric flux
168 
169  //- Rate of change of pressure
171 
172  //- Optional MRF zones
174 
175  //- Blending methods
176  blendingMethodTable blendingMethods_;
177 
178  //- Interface compression coefficients
179  cAlphaTable cAlphas_;
180 
181  //- Stabilisation for normalisation of the interface normal
183 
184 
185  // Sub Models
186 
187  //- Surface tension models
189 
190  //- Aspect ratio models
192 
193 
194  //- Flag to indicate that returned lists of fields are "complete"; i.e.,
195  // that an absence of force is returned as a constructed list of zeros,
196  // rather than a null pointer
197  static const bool fillFields_ = false;
198 
199 
200  // Protected member functions
201 
202  //- Calculate and return the mixture flux
204  (
205  const phaseModelList& phaseModels
206  ) const;
207 
208  //- Generate pairs
209  void generatePairs(const dictTable& modelDicts);
210 
211  //- Generate pairs and sub-model tables
212  template<class modelType>
213  void createSubModels
214  (
215  const dictTable& modelDicts,
216  HashTable
217  <
219  phasePairKey,
221  >& models
222  );
223 
224  //- Generate pairs and sub-model tables
225  template<class modelType>
227  (
228  const word& modelName,
229  HashTable
230  <
232  phasePairKey,
234  >& models,
235  const bool correctFixedFluxBCs = true
236  );
237 
238  //- Generate pairs and blended sub-model tables
239  template<class modelType>
241  (
242  const word& modelName,
243  HashTable
244  <
246  phasePairKey,
248  >& models,
249  const bool correctFixedFluxBCs = true
250  );
251 
252  //- Generate pairs and two-sided sub-model tables
253  template<class modelType>
255  (
256  const word& modelName,
257  HashTable
258  <
260  phasePairKey,
262  >& models,
263  const bool correctFixedFluxBCs = true
264  );
265 
266  //- Return the sum of the phase fractions of the moving phases
268 
269  //- Re-normalise the velocity of the phases
270  // around the specified mixture mean
271  void setMixtureU(const volVectorField& Um);
272 
273  //- Re-normalise the flux of the phases
274  // around the specified mixture mean
275  void setMixturePhi
276  (
278  const surfaceScalarField& phim
279  );
280 
281 
282  // Functions required for interface compression
283 
284  //- Normal to interface between two phases
285  // Used for interface compression
287  (
288  const volScalarField& alpha1,
289  const volScalarField& alpha2
290  ) const;
291 
292  //- Normal to interface between two phases dotted with face areas
293  // Used for interface compression
295  (
296  const volScalarField& alpha1,
297  const volScalarField& alpha2
298  ) const;
299 
300  //- Correction for the boundary condition on the unit normal nHat on
301  // walls to produce the correct contact angle.
302  //
303  // The dynamic contact angle is calculated from the component of the
304  // velocity on the direction of the interface, parallel to the wall.
305  // Used for interface compression
307  (
308  const phaseModel& alpha1,
309  const phaseModel& alpha2,
310  surfaceVectorField::Boundary& nHatb
311  ) const;
312 
313 
314  //- Curvature of interface between two phases
315  // Used for interface compression
317  (
318  const phaseModel& alpha1,
319  const phaseModel& alpha2
320  ) const;
321 
322 
323  // Functions required by twoPhaseSystem
324 
325  //- Return the drag coefficient for phase pair
326  virtual tmp<volScalarField> Kd(const phasePairKey& key) const = 0;
327 
328  //- Return the virtual mass coefficient for phase pair
329  virtual tmp<volScalarField> Vm(const phasePairKey& key) const = 0;
330 
331  //- Return the face drag coefficient for phase pair
333  (
334  const phasePairKey& key
335  ) const = 0;
336 
337 
338 public:
339 
340  //- Runtime type information
341  TypeName("phaseSystem");
342 
343  //- Default name of the phase properties dictionary
344  static const word propertiesName;
345 
346 
347  // Declare runtime construction
348 
350  (
351  autoPtr,
352  phaseSystem,
353  dictionary,
354  (
355  const fvMesh& mesh
356  ),
357  (mesh)
358  );
359 
360 
361  // Constructors
362 
363  //- Construct from fvMesh
364  phaseSystem(const fvMesh& mesh);
365 
366 
367  // Selectors
368 
370  (
371  const fvMesh& mesh
372  );
373 
374 
375  //- Destructor
376  virtual ~phaseSystem();
377 
378 
379  // Member Functions
380 
381  // Access
382 
383  //- Return the mesh
384  inline const fvMesh& mesh() const;
385 
386  //- Return the phase models
387  inline const phaseModelList& phases() const;
388 
389  //- Access the phase models
390  inline phaseModelList& phases();
391 
392  //- Return the models for phases that are moving
393  inline const phaseModelPartialList& movingPhases() const;
394 
395  //- Access the models for phases that are moving
396  inline phaseModelPartialList& movingPhases();
397 
398  //- Return the models for phases that are stationary
399  inline const phaseModelPartialList& stationaryPhases() const;
400 
401  //- Access the models for phases that are stationary
402  inline phaseModelPartialList& stationaryPhases();
403 
404  //- Return the models for phases that have variable temperature
405  inline const phaseModelPartialList& anisothermalPhases() const;
406 
407  //- Access the models for phases that have variable temperature
408  inline phaseModelPartialList& anisothermalPhases();
409 
410  //- Return the models for phases that have multiple species
411  inline const phaseModelPartialList& multiComponentPhases() const;
412 
413  //- Access the models for phases that have multiple species
414  inline phaseModelPartialList& multiComponentPhases();
415 
416  //- Return the phase pairs
417  inline const phasePairTable& phasePairs() const;
418 
419  //- Return the phase not given as an argument in a two-phase system
420  // An error is generated if the system is not two-phase
421  inline const phaseModel& otherPhase(const phaseModel& phase) const;
422 
423  //- Return the mixture flux
424  inline const surfaceScalarField& phi() const;
425 
426  //- Access the mixture flux
427  inline surfaceScalarField& phi();
428 
429  //- Return the rate of change of the pressure
430  inline const volScalarField& dpdt() const;
431 
432  //- Access the rate of change of the pressure
433  inline volScalarField& dpdt();
434 
435  //- Return MRF zones
436  inline const IOMRFZoneList& MRF() const;
437 
438  //- Access the fvOptions
439  inline fv::options& fvOptions() const;
440 
441 
442  // Sub-model lookup
443 
444  //- Check availability of a sub model for a given phase pair
445  template<class modelType>
446  bool foundSubModel(const phasePair& key) const;
447 
448  //- Return a sub model between a phase pair
449  template<class modelType>
450  const modelType& lookupSubModel(const phasePair& key) const;
451 
452  //- Check availability of a sub model between two phases
453  template<class modelType>
454  bool foundSubModel
455  (
456  const phaseModel& dispersed,
457  const phaseModel& continuous
458  ) const;
459 
460  //- Return a sub model between two phases
461  template<class modelType>
462  const modelType& lookupSubModel
463  (
464  const phaseModel& dispersed,
465  const phaseModel& continuous
466  ) const;
467 
468  //- Check availability of a blended sub model for a given phase pair
469  template<class modelType>
470  bool foundBlendedSubModel(const phasePair& key) const;
471 
472  //- Return a blended sub model between a phase pair
473  template<class modelType>
475  lookupBlendedSubModel(const phasePair& key) const;
476 
477 
478  // Field construction
479 
480  //- Fill up gaps in a phase-indexed list of fields with zeros
481  template
482  <
483  class Type,
484  template<class> class PatchField,
485  class GeoMesh
486  >
487  void fillFields
488  (
489  const word& name,
490  const dimensionSet& dims,
492  ) const;
493 
494  //- Fill up gaps in a phase-indexed table of fields with zeros
495  template
496  <
497  class Type,
498  template<class> class PatchField,
499  class GeoMesh
500  >
501  void fillFields
502  (
503  const word& name,
504  const dimensionSet& dims,
506  fieldTable
507  ) const;
508 
509 
510  // Properties
511 
512  //- Return the mixture density
513  tmp<volScalarField> rho() const;
514 
515  //- Return the mixture velocity
516  tmp<volVectorField> U() const;
517 
518  //- Return the aspect-ratio for a pair
519  tmp<volScalarField> E(const phasePairKey& key) const;
520 
521  //- Return the surface tension coefficient for a pair
522  tmp<volScalarField> sigma(const phasePairKey& key) const;
523 
524  //- Return the surface tension coefficient for a pair on a patch
526  (
527  const phasePairKey& key,
528  label patchi
529  ) const;
530 
531  //- Indicator of the proximity of the interface
532  // Field values are 1 near and 0 away for the interface.
534 
535  //- Return the mass transfer rate for an interface
536  virtual tmp<volScalarField> dmdtf(const phasePairKey& key) const;
537 
538  //- Return the mass transfer rates for each phase
539  virtual PtrList<volScalarField> dmdts() const;
540 
541  //- Return incompressibility
542  bool incompressible() const;
543 
544 
545  // Transfers
546 
547  //- Return the momentum transfer matrices for the cell-based
548  // algorithm
550 
551  //- Return the momentum transfer matrices for the face-based
552  // algorithm
554 
555  //- Return the implicit force coefficients for the face-based
556  // algorithm
557  virtual PtrList<surfaceScalarField> AFfs() const = 0;
558 
559  //- Return the force fluxes for the cell-based algorithm
561  (
563  ) = 0;
564 
565  //- Return the force fluxes for the face-based algorithm
567  (
569  ) = 0;
570 
571  //- Return the force fluxes for the cell-based algorithm
573  (
575  ) const = 0;
576 
577  //- Return the force fluxes for the face-based algorithm
579  (
581  ) const = 0;
582 
583  //- Return the explicit part of the drag force
585  (
587  ) const = 0;
588 
589  //- Returns true if the phase pressure is treated implicitly
590  // in the phase fraction equation
591  virtual bool implicitPhasePressure(const phaseModel& phase) const;
592 
593  //- Returns true if the phase pressure is treated implicitly
594  // in the phase fraction equation for any phase
595  virtual bool implicitPhasePressure() const;
596 
597  //- Return the phase diffusivity
598  // divided by the momentum central coefficient
600  (
603  ) const = 0;
604 
605  //- Solve the drag system for the new velocities and fluxes
606  virtual void partialElimination
607  (
612  ) = 0;
613 
614  //- Solve the drag system for the new fluxes
615  virtual void partialEliminationf
616  (
620  ) = 0;
621 
622  //- Return the flux corrections for the cell-based algorithm
624  (
626  const bool includeVirtualMass = false
627  ) const = 0;
628 
629  //- Return the heat transfer matrices
630  virtual autoPtr<heatTransferTable> heatTransfer() const = 0;
631 
632  //- Return the specie transfer matrices
633  virtual autoPtr<specieTransferTable> specieTransfer() const = 0;
634 
635  //- Return the surface tension force
637  (
638  const phaseModel& phase
639  ) const;
640 
641 
642  // Evolution
643 
644  //- Solve for the phase fractions
645  virtual void solve
646  (
649  );
650 
651  //- Correct the fluid properties other than those listed below
652  virtual void correct();
653 
654  //- Correct the continuity errors
655  virtual void correctContinuityError();
656 
657  //- Correct the kinematics
658  virtual void correctKinematics();
659 
660  //- Correct the thermodynamics
661  virtual void correctThermo();
662 
663  //- Correct the reactions
664  virtual void correctReactions();
665 
666  //- Correct the species mass fractions
667  virtual void correctSpecies();
668 
669  //- Correct the turbulence
670  virtual void correctTurbulence();
671 
672  //- Correct the energy transport e.g. alphat
673  virtual void correctEnergyTransport();
674 
675 
676  // IO
677 
678  //- Read base phaseProperties dictionary
679  virtual bool read();
680 };
681 
682 
683 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
684 
687 
688 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
689 
690 } // End namespace Foam
691 
692 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
693 
694 #include "phaseSystemI.H"
695 
696 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
697 
698 #ifdef NoRepository
699  #include "phaseSystemTemplates.C"
700 #endif
701 
702 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
703 
704 #endif
705 
706 // ************************************************************************* //
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:343
Hashing function class, shared by all the derived classes.
Definition: string.H:92
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.
fv::options & fvOptions() const
Access the fvOptions.
Definition: phaseSystemI.H:164
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:49
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
Definition: phaseSystem.H:187
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.
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.
declareRunTimeSelectionTable(autoPtr, phaseSystem, dictionary,(const fvMesh &mesh),(mesh))
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition: phaseSystem.H:87
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const dimensionedScalar deltaN_
Stabilisation for normalisation of the interface normal.
Definition: phaseSystem.H:181
PtrList< surfaceScalarField > alphafs(phases.size())
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:75
void setMixturePhi(const PtrList< surfaceScalarField > &alphafs, const surfaceScalarField &phim)
Re-normalise the flux of the phases.
phasePairTable phasePairs_
Phase pairs.
Definition: phaseSystem.H:163
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:169
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:77
phaseModelPartialList stationaryPhaseModels_
Stationary phase models.
Definition: phaseSystem.H:154
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:79
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 specialization for hashing pointers.
Definition: HashPtrTable.H:50
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:123
Finite-volume options.
Definition: fvOptions.H:52
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.
HashPtrTable< HashPtrTable< volScalarField >, phasePairKey, phasePairKey::hash > dmidtfTable
Definition: phaseSystem.H:100
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:148
alpha2
Definition: alphaEqn.H:115
phaseModelPartialList multiComponentPhaseModels_
Multi-component phase models.
Definition: phaseSystem.H:160
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
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
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);label pRefCell=0;scalar pRefValue=0.0;if(fluid.incompressible()){ p=max(p_rgh+fluid.rho() *gh, pMin);if(p_rgh.needReference()) { setRefCell(p, p_rgh, pimple.dict(), pRefCell, pRefValue);p+=dimensionedScalar("p", p.dimensions(), pRefValue - getRefCellValue(p, pRefCell));p_rgh=p - fluid.rho() *gh;}}mesh.setFluxRequired(p_rgh.name());PtrList< volScalarField > rAUs
Definition: createFields.H:67
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdtfTable
Definition: phaseSystem.H:91
tmp< volVectorField > U() const
Return the mixture velocity.
volScalarField & alpha1(mixture.alpha1())
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:81
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:145
UPtrList< phaseModel > phaseModelPartialList
Definition: phaseSystem.H:83
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:141
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.
const word & name() const
Name function is needed to disambiguate those inherited.
An STL-conforming hash table.
Definition: HashTable.H:61
phaseModelPartialList anisothermalPhaseModels_
Anisothermal phase models.
Definition: phaseSystem.H:157
cAlphaTable cAlphas_
Interface compression coefficients.
Definition: phaseSystem.H:178
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases.
phaseModelPartialList movingPhaseModels_
Moving phase models.
Definition: phaseSystem.H:151
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:190
IOMRFZoneList MRF_
Optional MRF zones.
Definition: phaseSystem.H:172
PtrList< surfaceScalarField > rAUfs
Definition: createFields.H:68
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 specializations.
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.
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:175
HashTable< autoPtr< blendingMethod >, word, word::hash > blendingMethodTable
Definition: phaseSystem.H:114
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:132
HashTable< scalar, phasePairKey, phasePairKey::hash > cAlphaTable
Definition: phaseSystem.H:135
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:166
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
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries...
Definition: IOMRFZoneList.H:66
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Definition: phaseSystem.H:110
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.