MassTransferPhaseSystem.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::MassTransferPhaseSystem
28 
29 Description
30  Class for mass transfer between phases
31 
32 SourceFiles
33  MassTransferPhaseSystem.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef MassTransferPhaseSystem_H
38 #define MassTransferPhaseSystem_H
39 
40 #include "phaseSystem.H"
41 #include "HashPtrTable.H"
42 #include "interfaceCompositionModel.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class MassTransferPhaseSystem Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 template<class BasePhaseSystem>
55 :
56  public BasePhaseSystem
57 {
58 public:
59 
60  // Public typedefs
61 
62  typedef
63  HashTable
64  <
68  >
70 
71 
73 
74 protected:
75 
76  // Protected typedefs
77 
78  typedef
80  <
84  >
85  dmdtTable;
86 
87 
88  // Protected Data
89 
90  //- Overall inter-phase mass transfer rates [Kg/s]
92 
93  //- Mass transfer models
95 
96 
97  // Protected Member Functions
98 
99  //- Calculate L between phases
101  (
102  const volScalarField& dmdtNetki,
103  const phasePairKey& keyik,
104  const phasePairKey& keyki,
105  const volScalarField& T
106  ) const;
107 
108 
109 public:
110 
111  // Constructors
112 
113  //- Construct from fvMesh
114  explicit MassTransferPhaseSystem(const fvMesh&);
115 
116 
117  //- Destructor
118  virtual ~MassTransferPhaseSystem() = default;
119 
120 
121  // Member Functions
122 
123  //- Return total interfacial mass flow rate
125 
126 
127  // Mass transfer functions
128 
129  //- Return the heat transfer matrix
130  // NOTE: Call KSu and KSp with T as variable,if not provided uses dmdt.
132 
133  //- Return the volumetric rate transfer matrix
134  // NOTE: Call KSu and KSp with p as variable,if not provided uses dmdt.
136 
137  //- Correct/calculates mass sources dmdt for phases
138  // NOTE: Call the kexp() for all the mass transfer models.
139  virtual void correctMassSources(const volScalarField& T);
140 
141  //- Calculate mass transfer for alpha's
142  virtual void alphaTransfer(SuSpTable& Su, SuSpTable& Sp);
143 
144  //- Calculate mass transfer for species
145  virtual void massSpeciesTransfer
146  (
147  const phaseModel& phase,
150  const word speciesName
151  );
152 
153  //- Add volume change in pEq
154  virtual bool includeVolChange();
155 };
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 } // End namespace Foam
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 #ifdef NoRepository
165 # include "MassTransferPhaseSystem.C"
166 #endif
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #endif
171 
172 // ************************************************************************* //
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::MassTransferPhaseSystem::SuSpTable
HashTable< volScalarField::Internal > SuSpTable
Definition: MassTransferPhaseSystem.H:71
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::phasePairKey::hasher
Hashing functor for phasePairKey.
Definition: phasePairKey.H:122
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::MassTransferPhaseSystem::calculateL
tmp< volScalarField > calculateL(const volScalarField &dmdtNetki, const phasePairKey &keyik, const phasePairKey &keyki, const volScalarField &T) const
Calculate L between phases.
Definition: MassTransferPhaseSystem.C:79
Foam::MassTransferPhaseSystem::volTransfer
virtual tmp< fvScalarMatrix > volTransfer(const volScalarField &p)
Return the volumetric rate transfer matrix.
Definition: MassTransferPhaseSystem.C:305
Sp
zeroField Sp
Definition: alphaSuSp.H:2
Foam::glTF::key
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:108
Foam::MassTransferPhaseSystem::MassTransferPhaseSystem
MassTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
Definition: MassTransferPhaseSystem.C:40
Foam::MassTransferPhaseSystem::dmdtTable
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdtTable
Definition: MassTransferPhaseSystem.H:84
MassTransferPhaseSystem.C
Foam::MassTransferPhaseSystem::massSpeciesTransfer
virtual void massSpeciesTransfer(const phaseModel &phase, volScalarField::Internal &Su, volScalarField::Internal &Sp, const word speciesName)
Calculate mass transfer for species.
Definition: MassTransferPhaseSystem.C:709
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::GeometricField< scalar, fvPatchField, volMesh >::Internal
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Definition: GeometricField.H:107
Foam::MassTransferPhaseSystem::dmdt_
dmdtTable dmdt_
Overall inter-phase mass transfer rates [Kg/s].
Definition: MassTransferPhaseSystem.H:90
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::MassTransferPhaseSystem::~MassTransferPhaseSystem
virtual ~MassTransferPhaseSystem()=default
Destructor.
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::phasePairKey
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:65
Foam::MassTransferPhaseSystem::massTransferModels_
massTransferModelTable massTransferModels_
Mass transfer models.
Definition: MassTransferPhaseSystem.H:93
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::MassTransferPhaseSystem
Class for mass transfer between phases.
Definition: MassTransferPhaseSystem.H:53
Foam::MassTransferPhaseSystem::massTransferModelTable
HashTable< autoPtr< interfaceCompositionModel >, phasePairKey, phasePairKey::hash > massTransferModelTable
Definition: MassTransferPhaseSystem.H:68
Foam::MassTransferPhaseSystem::heatTransfer
virtual tmp< fvScalarMatrix > heatTransfer(const volScalarField &T)
Return the heat transfer matrix.
Definition: MassTransferPhaseSystem.C:163
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
Foam::autoPtr< interfaceCompositionModel >
Foam::HashPtrTable
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
Definition: HashPtrTable.H:54
HashPtrTable.H
Foam::MassTransferPhaseSystem::includeVolChange
virtual bool includeVolChange()
Add volume change in pEq.
Definition: MassTransferPhaseSystem.C:731
Foam::MassTransferPhaseSystem::correctMassSources
virtual void correctMassSources(const volScalarField &T)
Correct/calculates mass sources dmdt for phases.
Definition: MassTransferPhaseSystem.C:454
Foam::MassTransferPhaseSystem::dmdt
tmp< volScalarField > dmdt(const phasePairKey &key) const
Return total interfacial mass flow rate.
Definition: MassTransferPhaseSystem.C:134
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::MassTransferPhaseSystem::alphaTransfer
virtual void alphaTransfer(SuSpTable &Su, SuSpTable &Sp)
Calculate mass transfer for alpha's.
Definition: MassTransferPhaseSystem.C:505