Loading...
Searching...
No Matches
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-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::MassTransferPhaseSystem
28
29Description
30 Class for mass transfer between phases
31
32SourceFiles
33 MassTransferPhaseSystem.C
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef Foam_multiphaseInter_MassTransferPhaseSystem_H
38#define Foam_multiphaseInter_MassTransferPhaseSystem_H
39
41#include "HashPtrTable.H"
42#include "interfaceCompositionModel.H"
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46namespace Foam
47{
48
49using namespace multiphaseInter;
50
51/*---------------------------------------------------------------------------*\
52 Class MassTransferPhaseSystem Declaration
53\*---------------------------------------------------------------------------*/
54
55template<class BasePhaseSystem>
57:
58 public BasePhaseSystem
59{
60public:
61
62 // Public typedefs
63
64 typedef
66 <
70 >
72
73
75
76protected:
77
78 // Protected typedefs
79
80 typedef
82 <
86 >
88
89
90 // Protected Data
91
92 //- Overall inter-phase mass transfer rates [Kg/s]
94
95 //- Mass transfer models
97
98
99 // Protected Member Functions
100
101 //- Calculate L between phases
104 const volScalarField& dmdtNetki,
105 const phasePairKey& keyik,
106 const phasePairKey& keyki,
107 const volScalarField& T
108 ) const;
109
110
111public:
112
113 // Constructors
114
115 //- Construct from fvMesh
116 explicit MassTransferPhaseSystem(const fvMesh&);
117
118
119 //- Destructor
120 virtual ~MassTransferPhaseSystem() = default;
121
122
123 // Member Functions
124
125 //- Return total interfacial mass flow rate
126 tmp<volScalarField> dmdt(const phasePairKey& key) const;
127
128
129 // Mass transfer functions
130
131 //- Return the heat transfer matrix
132 // NOTE: Call KSu and KSp with T as variable,if not provided uses dmdt.
134
135 //- Return the volumetric rate transfer matrix
136 // NOTE: Call KSu and KSp with p as variable,if not provided uses dmdt.
138
139 //- Correct/calculates mass sources dmdt for phases
140 // NOTE: Call the kexp() for all the mass transfer models.
141 virtual void correctMassSources(const volScalarField& T);
142
143 //- Calculate mass transfer for alpha's
144 virtual void alphaTransfer(SuSpTable& Su, SuSpTable& Sp);
145
146 //- Calculate mass transfer for species
147 virtual void massSpeciesTransfer
148 (
149 const Foam::phaseModel& phase,
152 const word speciesName
153 );
154
155 //- Add volume change in pEq
156 virtual bool includeVolChange();
158
159
160// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161
162} // End namespace Foam
163
164// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165
166#ifdef NoRepository
168#endif
169
170// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171
172#endif
173
174// ************************************************************************* //
DimensionedField< scalar, volMesh > Internal
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdtTable
tmp< volScalarField > dmdt(const phasePairKey &key) const
Return total interfacial mass flow rate.
virtual void correctMassSources(const volScalarField &T)
Correct/calculates mass sources dmdt for phases.
MassTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual ~MassTransferPhaseSystem()=default
Destructor.
HashTable< autoPtr< interfaceCompositionModel >, phasePairKey, phasePairKey::hash > massTransferModelTable
virtual void massSpeciesTransfer(const Foam::phaseModel &phase, volScalarField::Internal &Su, volScalarField::Internal &Sp, const word speciesName)
Calculate mass transfer for species.
HashTable< volScalarField::Internal > SuSpTable
virtual tmp< fvScalarMatrix > volTransfer(const volScalarField &p)
Return the volumetric rate transfer matrix.
virtual bool includeVolChange()
Add volume change in pEq.
virtual void alphaTransfer(SuSpTable &Su, SuSpTable &Sp)
Calculate mass transfer for alpha's.
virtual tmp< fvScalarMatrix > heatTransfer(const volScalarField &T)
Return the heat transfer matrix.
tmp< volScalarField > calculateL(const volScalarField &dmdtNetki, const phasePairKey &keyik, const phasePairKey &keyki, const volScalarField &T) const
Calculate L between phases.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
virtual volScalarField & p()
Pressure [Pa].
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const phaseModel & phase(const label i) const
Constant access phase model i.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
An ordered or unorder pair of phase names. Typically specified as follows.
phasePairKey::hasher hash
Alternative name for functor.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)