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-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 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 "phaseInterface.H"
43 #include "phaseInterfaceKey.H"
44 #include "HashPtrTable.H"
45 #include "PtrListDictionary.H"
46 #include "hashedWordList.H"
47 
48 #include "pimpleNoLoopControl.H"
49 
50 #include "IOMRFZoneList.H"
51 #include "fvModels.H"
52 #include "fvConstraints.H"
53 
54 #include "volFields.H"
55 #include "surfaceFields.H"
56 #include "fvMatricesFwd.H"
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace Foam
61 {
62 
63 class interfaceSurfaceTensionModel;
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
91  <
95  >
96  dmdtfTable;
97 
98  typedef
100  <
104  >
105  dmidtfTable;
106 
107 
108 
109 protected:
110 
111  // Protected typedefs
112 
113  typedef
114  HashTable
115  <
119  >
121 
123  cAlphaTable;
124 
125 
126  // Protected data
127 
128  //- Reference to the mesh
129  const fvMesh& mesh_;
130 
131  //- Reference to pimpleNoLoopControl
133 
134  //- Optional MRF zones
136 
137  //- Name of optional reference phase which is not solved for
138  // but obtained from the sum of the other phases
140 
141  //- Phase models
143 
144  //- Moving phase models
146 
147  //- Stationary phase models
149 
150  //- Anisothermal phase models
152 
153  //- Multi-component phase models
155 
156  //- Total volumetric flux
158 
159  //- Rate of change of pressure
161 
162  //- Interface compression coefficients
164 
165  //- Stabilisation for normalisation of the interface normal
167 
168 
169  // Sub Models
170 
171  //- Surface tension models
173 
174 
175  //- Flag to indicate that returned lists of fields are "complete"; i.e.,
176  // that an absence of force is returned as a constructed list of zeros,
177  // rather than a null pointer
178  static const bool fillFields_ = false;
179 
180 
181  // Protected member functions
182 
183  //- Calculate and return the mixture flux
185  (
186  const phaseModelList& phaseModels
187  ) const;
188 
189  //- Return the sum of the phase fractions of the moving phases
191 
192  //- Re-normalise the velocity of the phases
193  // around the specified mixture mean
194  void setMixtureU(const volVectorField& Um);
195 
196 
197  // Functions required for interface compression
198 
199  //- Normal to interface between two phases
200  // Used for interface compression
202  (
203  const volScalarField& alpha1,
204  const volScalarField& alpha2
205  ) const;
206 
207  //- Normal to interface between two phases dotted with face areas
208  // Used for interface compression
210  (
211  const volScalarField& alpha1,
212  const volScalarField& alpha2
213  ) const;
214 
215  //- Curvature of interface between two phases
216  // Used for interface compression
218  (
219  const phaseModel& alpha1,
220  const phaseModel& alpha2
221  ) const;
222 
223 
224  // Sub-model construction
225 
226  //- Return the dictionary containing interfacial model or value
227  // settings for the given name. Performs backwards compatibility
228  // conversions.
229  template<class Type>
230  dictionary interfacialDict(const word& name) const;
231 
232 
233 public:
234 
235  //- Runtime type information
236  TypeName("phaseSystem");
237 
238  //- Default name of the phase properties dictionary
239  static const word propertiesName;
240 
241 
242  // Declare runtime construction
243 
245  (
246  autoPtr,
247  phaseSystem,
248  dictionary,
249  (
250  const fvMesh& mesh
251  ),
252  (mesh)
253  );
254 
255 
256  // Constructors
257 
258  //- Construct from fvMesh
259  phaseSystem(const fvMesh& mesh);
260 
261 
262  // Selectors
263 
265  (
266  const fvMesh& mesh
267  );
268 
269 
270  //- Destructor
271  virtual ~phaseSystem();
272 
273 
274  // Member Functions
275 
276  // Access
277 
278  //- Return the mesh
279  inline const fvMesh& mesh() const;
280 
281  //- Return pimpleNoLoopControl
282  inline const pimpleNoLoopControl& pimple() const;
283 
284  //- Return the phase models
285  inline const phaseModelList& phases() const;
286 
287  //- Access the phase models
288  inline phaseModelList& phases();
289 
290  //- Return the models for phases that are moving
291  inline const phaseModelPartialList& movingPhases() const;
292 
293  //- Access the models for phases that are moving
295 
296  //- Return the models for phases that are stationary
297  inline const phaseModelPartialList& stationaryPhases() const;
298 
299  //- Access the models for phases that are stationary
301 
302  //- Return the models for phases that have variable temperature
303  inline const phaseModelPartialList& anisothermalPhases() const;
304 
305  //- Access the models for phases that have variable temperature
307 
308  //- Return the models for phases that have multiple species
309  inline const phaseModelPartialList& multicomponentPhases() const;
310 
311  //- Access the models for phases that have multiple species
313 
314  //- Return the phase not given as an argument in a two-phase system
315  // An error is generated if the system is not two-phase
316  inline const phaseModel& otherPhase(const phaseModel& phase) const;
317 
318  //- Return the mixture flux
319  inline const surfaceScalarField& phi() const;
320 
321  //- Access the mixture flux
322  inline surfaceScalarField& phi();
323 
324  //- Return the rate of change of the pressure
325  inline const volScalarField& dpdt() const;
326 
327  //- Access the rate of change of the pressure
328  inline volScalarField& dpdt();
329 
330  //- Return MRF zones
331  inline const IOMRFZoneList& MRF() const;
332 
333  //- Access the fvModels
335 
336  //- Access the fvModels
337  inline const Foam::fvModels& fvModels() const;
338 
339  //- Access the fvConstraints
341 
342  //- Access the fvConstraints
343  inline const Foam::fvConstraints& fvConstraints() const;
344 
345 
346  // Sub-model construction
347 
348  //- Return the model name. This is the same as the model's typename
349  // but without "Model" on the end.
350  template<class ModelType>
351  word modelName() const;
352 
353  //- Generate interfacial-model lists
354  template<class ModelType, class ... InterfaceTypes>
356  (
357  const dictionary& dict,
358  const phaseInterface& interface,
359  PtrList<phaseInterface>& interfaces,
360  PtrList<ModelType>& models
361  ) const;
362 
363  //- Generate interfacial-model tables
364  template<class ModelType>
366  (
367  const dictionary& dict,
368  HashTable
369  <
373  >& models
374  ) const;
375 
376  //- Generate interfacial-model tables
377  template<class ModelType>
379  (
380  HashTable
381  <
385  >& models
386  ) const;
387 
388  //- Generate interfacial-model tables
389  template<class ValueType>
391  (
392  const dictionary& dict,
393  HashTable
394  <
395  ValueType,
398  >& values
399  ) const;
400 
401  //- Generate interfacial-model tables
402  template<class ValueType>
404  (
405  const word& valueName,
406  HashTable
407  <
408  ValueType,
411  >& values
412  ) const;
413 
414  //- Return the dictionary from which to construct a low-level
415  // sub-model. Checks that there is just one sub-dictionary then
416  // returns it.
417  template<class ModelType>
418  static const dictionary& modelSubDict(const dictionary& dict);
419 
420  //- Check that mass transfer is supported across the given interface
421  template<class ModelType>
422  void validateMassTransfer(const phaseInterface& interface) const;
423 
424 
425  // Sub-model lookup
426 
427  //- Check availability of a sub model for a given interface
428  template<class ModelType>
429  bool foundInterfacialModel(const phaseInterface& interface) const;
430 
431  //- Return a sub model for an interface
432  template<class ModelType>
433  const ModelType& lookupInterfacialModel
434  (
435  const phaseInterface& interface
436  ) const;
437 
438 
439  // Field construction
440 
441  //- Fill up gaps in a phase-indexed list of fields with zeros
442  template
443  <
444  class Type,
445  template<class> class PatchField,
446  class GeoMesh
447  >
448  void fillFields
449  (
450  const word& name,
451  const dimensionSet& dims,
453  ) const;
454 
455  //- Fill up gaps in a phase-indexed table of fields with zeros
456  template
457  <
458  class Type,
459  template<class> class PatchField,
460  class GeoMesh
461  >
462  void fillFields
463  (
464  const word& name,
465  const dimensionSet& dims,
467  fieldTable
468  ) const;
469 
470 
471  // Properties
472 
473  //- Return the mixture density
474  tmp<volScalarField> rho() const;
475 
476  //- Return the mixture velocity
477  tmp<volVectorField> U() const;
478 
479  //- Return the surface tension coefficient for an interface
480  tmp<volScalarField> sigma(const phaseInterfaceKey& key) const;
481 
482  //- Return the surface tension coefficient for an interface on a
483  // patch
485  (
486  const phaseInterfaceKey& key,
487  const label patchi
488  ) const;
489 
490  //- Indicator of the proximity of the interface
491  // Field values are 1 near and 0 away for the interface.
493 
494  //- Stabilisation for normalisation of the interface normal
495  inline const dimensionedScalar& deltaN() const;
496 
497  //- Return the mass transfer rate for an interface
498  virtual tmp<volScalarField> dmdtf
499  (
500  const phaseInterfaceKey& key
501  ) const;
502 
503  //- Return the mass transfer rates for each phase
504  virtual PtrList<volScalarField> dmdts() const;
505 
506  //- Return the mass transfer pressure implicit coefficients
507  // for each phase
508  virtual PtrList<volScalarField> d2mdtdps() const;
509 
510  //- Return incompressibility
511  bool incompressible() const;
512 
513 
514  // Transfers
515 
516  //- Return the momentum transfer matrices for the cell-based
517  // algorithm
519 
520  //- Return the momentum transfer matrices for the face-based
521  // algorithm
523 
524  //- Return the implicit force coefficients for the face-based
525  // algorithm
526  virtual PtrList<surfaceScalarField> KdVmfs() const = 0;
527 
528  //- Return the force fluxes for the cell-based algorithm
529  virtual PtrList<surfaceScalarField> Fs() const = 0;
530 
531  //- Return the force fluxes for the face-based algorithm
532  virtual PtrList<surfaceScalarField> Ffs() const = 0;
533 
534  //- Return the force fluxes for the cell-based algorithm
535  virtual PtrList<surfaceScalarField> KdPhis() const = 0;
536 
537  //- Return the force fluxes for the face-based algorithm
538  virtual PtrList<surfaceScalarField> KdPhifs() const = 0;
539 
540  //- Return the implicit part of the drag force
541  virtual PtrList<volScalarField> Kds() const = 0;
542 
543  //- Return the explicit part of the drag force
544  virtual PtrList<volVectorField> KdUs() const = 0;
545 
546  //- Returns true if the phase pressure is treated implicitly
547  // in the phase fraction equation
548  virtual bool implicitPhasePressure(const phaseModel& phase) const;
549 
550  //- Returns true if the phase pressure is treated implicitly
551  // in the phase fraction equation for any phase
552  virtual bool implicitPhasePressure() const;
553 
554  //- Return the phase diffusivity
555  // divided by the momentum central coefficient
557  (
558  const PtrList<volScalarField>& rAUs,
559  const PtrList<surfaceScalarField>& rAUfs
560  ) const = 0;
561 
562  //- Return the flux corrections for the cell-based algorithm
563  virtual PtrList<surfaceScalarField> ddtCorrs() const = 0;
564 
565  //- Set the cell and faces drag correction fields
566  virtual void dragCorrs
567  (
569  PtrList<surfaceScalarField>& dragCorrf
570  ) const = 0;
571 
572  //- Solve the drag system for the new velocities and fluxes
573  virtual void partialElimination
574  (
575  const PtrList<volScalarField>& rAUs,
577  const PtrList<surfaceScalarField>& alphafs,
578  const PtrList<surfaceScalarField>& rAUfs,
580  ) = 0;
581 
582  //- Solve the drag system for the new fluxes
583  virtual void partialEliminationf
584  (
585  const PtrList<surfaceScalarField>& rAUfs,
586  const PtrList<surfaceScalarField>& alphafs,
588  ) = 0;
589 
590  //- Re-normalise the flux of the phases
591  // around the specified mixture mean
592  void setMixturePhi
593  (
594  const PtrList<surfaceScalarField>& alphafs,
595  const surfaceScalarField& phim
596  );
597 
598  //- Return the heat transfer matrices
599  virtual autoPtr<heatTransferTable> heatTransfer() const = 0;
600 
601  //- Return the specie transfer matrices
602  virtual autoPtr<specieTransferTable> specieTransfer() const = 0;
603 
604  //- Return the surface tension force
606  (
607  const phaseModel& phase
608  ) const;
609 
610 
611  // Evolution
612 
613  //- Solve for the phase fractions
614  virtual void solve
615  (
616  const PtrList<volScalarField>& rAUs,
617  const PtrList<surfaceScalarField>& rAUfs
618  );
619 
620  //- Correct the fluid properties other than those listed below
621  virtual void correct();
622 
623  //- Correct the continuity errors
624  virtual void correctContinuityError();
625 
626  //- Correct the kinematics
627  virtual void correctKinematics();
628 
629  //- Correct the thermodynamics
630  virtual void correctThermo();
631 
632  //- Correct the reactions
633  virtual void correctReactions();
634 
635  //- Correct the species mass fractions
636  virtual void correctSpecies();
637 
638  //- Predict the momentumTransport
639  virtual void predictMomentumTransport();
640 
641  //- Predict the energy transport e.g. alphat
642  virtual void predictThermophysicalTransport();
643 
644  //- Correct the momentumTransport
645  virtual void correctMomentumTransport();
646 
647  //- Correct the energy transport e.g. alphat
648  virtual void correctThermophysicalTransport();
649 
650  //- Update the fluid properties for mesh changes
651  virtual void meshUpdate();
652 
653  //- Correct fixed-flux BCs to be consistent with the velocity BCs
654  void correctBoundaryFlux();
655 
656  void correctPhi
657  (
658  const volScalarField& p_rgh,
659  const autoPtr<volScalarField>& divU,
662  );
663 
664 
665  // IO
666 
667  //- Read base phaseProperties dictionary
668  virtual bool read();
669 };
670 
671 
672 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
673 
676 
677 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
678 
679 } // End namespace Foam
680 
681 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
682 
683 #include "phaseSystemI.H"
684 
685 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
686 
687 #ifdef NoRepository
688  #include "phaseSystemTemplates.C"
689 #endif
690 
691 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
692 
693 #endif
694 
695 // ************************************************************************* //
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Generic GeometricField class.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
An STL-conforming hash table.
Definition: HashTable.H:127
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries,...
Definition: IOMRFZoneList.H:68
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Dimension set for the base types.
Definition: dimensionSet.H:122
Finite volume constraints.
Definition: fvConstraints.H:62
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume models.
Definition: fvModels.H:65
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 ...
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 autoPtr< specieTransferTable > specieTransfer() const =0
Return the specie transfer matrices.
virtual PtrList< surfaceScalarField > KdPhis() const =0
Return the force fluxes for the cell-based algorithm.
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
virtual PtrList< surfaceScalarField > Ffs() const =0
Return the force fluxes for the face-based algorithm.
dictionary interfacialDict(const word &name) const
Return the dictionary containing interfacial model or value.
void correctPhi(const volScalarField &p_rgh, const autoPtr< volScalarField > &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pimple)
Definition: phaseSystem.C:730
bool foundInterfacialModel(const phaseInterface &interface) const
Check availability of a sub model for a given interface.
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:238
virtual PtrList< volScalarField > Kds() const =0
Return the implicit part of the drag force.
void validateMassTransfer(const phaseInterface &interface) const
Check that mass transfer is supported across the given interface.
void generateInterfacialModels(const dictionary &dict, const phaseInterface &interface, PtrList< phaseInterface > &interfaces, PtrList< ModelType > &models) const
Generate interfacial-model lists.
const IOMRFZoneList & MRF() const
Return MRF zones.
Definition: phaseSystemI.H:157
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:111
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
const dimensionedScalar & deltaN() const
Stabilisation for normalisation of the interface normal.
Definition: phaseSystemI.H:187
virtual PtrList< surfaceScalarField > ddtCorrs() const =0
Return the flux corrections for the cell-based algorithm.
virtual void correctSpecies()
Correct the species mass fractions.
Definition: phaseSystem.C:639
const Foam::fvModels & fvModels() const
Access the fvModels.
Definition: phaseSystemI.H:169
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: phaseSystemI.H:133
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseSystem.C:621
const pimpleNoLoopControl & pimple_
Reference to pimpleNoLoopControl.
Definition: phaseSystem.H:131
const volScalarField & dpdt() const
Return the rate of change of the pressure.
Definition: phaseSystemI.H:145
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
const phaseModelPartialList & multicomponentPhases() const
Return the models for phases that have multiple species.
Definition: phaseSystemI.H:97
virtual autoPtr< heatTransferTable > heatTransfer() const =0
Return the heat transfer matrices.
virtual PtrList< surfaceScalarField > KdVmfs() const =0
Return the implicit force coefficients for the face-based.
virtual PtrList< volVectorField > KdUs() const =0
Return the explicit part of the drag force.
phaseModelPartialList anisothermalPhaseModels_
Anisothermal phase models.
Definition: phaseSystem.H:150
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:55
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
HashTable< scalar, phaseInterfaceKey, phaseInterfaceKey::hash > cAlphaTable
Definition: phaseSystem.H:122
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
static autoPtr< phaseSystem > New(const fvMesh &mesh)
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:41
PtrListDictionary< phaseModel > phaseModelList
Definition: phaseSystem.H:84
const phaseModelPartialList & anisothermalPhases() const
Return the models for phases that have variable temperature.
Definition: phaseSystemI.H:83
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &KdPhifs)=0
Solve the drag system for the new fluxes.
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
Definition: phaseSystem.C:468
virtual void dragCorrs(PtrList< volVectorField > &dragCorrs, PtrList< surfaceScalarField > &dragCorrf) const =0
Set the cell and faces drag correction fields.
virtual void solve(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs)
Solve for the phase fractions.
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.
const dimensionedScalar deltaN_
Stabilisation for normalisation of the interface normal.
Definition: phaseSystem.H:165
interfaceSurfaceTensionModelTable interfaceSurfaceTensionModels_
Surface tension models.
Definition: phaseSystem.H:171
const pimpleNoLoopControl & pimple() const
Return pimpleNoLoopControl.
Definition: phaseSystemI.H:34
virtual autoPtr< momentumTransferTable > momentumTransferf()=0
Return the momentum transfer matrices for the face-based.
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:82
virtual void correctThermophysicalTransport()
Correct the energy transport e.g. alphat.
Definition: phaseSystem.C:675
static const bool fillFields_
Flag to indicate that returned lists of fields are "complete"; i.e.,.
Definition: phaseSystem.H:177
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
virtual autoPtr< momentumTransferTable > momentumTransfer()=0
Return the momentum transfer matrices for the cell-based.
const phaseModelPartialList & stationaryPhases() const
Return the models for phases that are stationary.
Definition: phaseSystemI.H:69
word modelName() const
Return the model name. This is the same as the model's typename.
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:156
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:78
virtual PtrList< surfaceScalarField > Fs() const =0
Return the force fluxes for the cell-based algorithm.
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
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:80
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Definition: phaseSystem.C:210
virtual PtrList< surfaceScalarField > KdPhifs() const =0
Return the force fluxes for the face-based algorithm.
void setMixturePhi(const PtrList< surfaceScalarField > &alphafs, const surfaceScalarField &phim)
Re-normalise the flux of the phases.
Definition: phaseSystem.C:125
HashTable< autoPtr< interfaceSurfaceTensionModel >, phaseInterfaceKey, phaseInterfaceKey::hash > interfaceSurfaceTensionModelTable
Definition: phaseSystem.H:119
virtual void partialElimination(const PtrList< volScalarField > &rAUs, const PtrList< volVectorField > &KdUs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &KdPhis)=0
Solve the drag system for the new velocities and fluxes.
virtual ~phaseSystem()
Destructor.
Definition: phaseSystem.C:345
HashPtrTable< HashPtrTable< volScalarField >, phaseInterfaceKey, phaseInterfaceKey::hash > dmidtfTable
Definition: phaseSystem.H:104
void generateInterfacialValues(const dictionary &dict, HashTable< ValueType, phaseInterfaceKey, phaseInterfaceKey::hash > &values) const
Generate interfacial-model tables.
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:28
TypeName("phaseSystem")
Runtime type information.
UPtrList< phaseModel > phaseModelPartialList
Definition: phaseSystem.H:86
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
IOMRFZoneList MRF_
Optional MRF zones.
Definition: phaseSystem.H:134
HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > dmdtfTable
Definition: phaseSystem.H:95
const Foam::fvConstraints & fvConstraints() const
Access the fvConstraints.
Definition: phaseSystemI.H:181
virtual tmp< surfaceScalarField > alphaDByAf(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs) const =0
Return the phase diffusivity.
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
volScalarField dpdt_
Rate of change of pressure.
Definition: phaseSystem.H:159
virtual void predictMomentumTransport()
Predict the momentumTransport.
Definition: phaseSystem.C:648
static const dictionary & modelSubDict(const dictionary &dict)
Return the dictionary from which to construct a low-level.
declareRunTimeSelectionTable(autoPtr, phaseSystem, dictionary,(const fvMesh &mesh),(mesh))
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
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
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.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Forward declarations of fvMatrix specialisations.
label patchi
volScalarField sf(fieldObject, mesh)
Namespace for OpenFOAM.
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
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:853
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dictionary dict
Foam::surfaceFields.