MultiComponentPhaseModel.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::MultiComponentPhaseModel
28
29Description
30 Class which represents a phase with multiple species. Returns the species'
31 mass fractions, and their governing equations.
32
33SourceFiles
34 MultiComponentPhaseModel.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef MultiComponentPhaseModel_H
39#define MultiComponentPhaseModel_H
40
41#include "phaseModel.H"
42#include "hashedWordList.H"
43#include "rhoReactionThermo.H"
44#include "basicThermo.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51/*---------------------------------------------------------------------------*\
52 Class phaseModel Declaration
53\*---------------------------------------------------------------------------*/
54
55template<class BasePhaseModel, class phaseThermo>
57:
58 public BasePhaseModel
59{
60protected:
61
62 // Protected data
63
64 //- Species table
66
67 //- Inert species index
68 label inertIndex_;
69
70 //- Thermophysical model
72
73 //- Ptr list of volumetric fractions for species
75
76 //- Add diffusion transport on Yi's Eq
77 bool addDiffusion_;
78
79 //- Schmidt number
80 scalar Sct_;
81
82
83 // Protected functions
84
85 //- Transfor volume fraction into mass fractions
87
88 //- Transfor mass fraction into volume fractions
90
91
92public:
93
94 // Constructors
95
97 (
99 const word& phaseName
100 );
101
102
103 //- Destructor
104 virtual ~MultiComponentPhaseModel() = default;
105
106
107 // Member Functions
108
109 // Access
110
111 //- Species table
112 const hashedWordList& species() const
113 {
114 return species_;
115 }
116
117
118 // Thermo
119
120 //- Access to thermo
121 virtual const phaseThermo& thermo() const;
122
123 //- Access non-const thermo
124 virtual phaseThermo& thermo();
125
126 //- Correct phase thermo
127 virtual void correct();
128
129 //- Solve species fraction equation
130 virtual void solveYi
131 (
134 );
135
136 //- Constant access the species mass fractions
137 virtual const PtrList<volScalarField>& Y() const;
138
139 //- Access the species mass fractions
140 virtual PtrList<volScalarField>& Y();
141
142 //- Return inert species index
143 label inertIndex() const;
144};
145
146
147// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148
149} // End namespace Foam
150
151// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152
153#ifdef NoRepository
154# include "MultiComponentPhaseModel.C"
155#endif
156
157// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158
159#endif
160
161// ************************************************************************* //
twoPhaseSystem & fluid
Class which represents a phase with multiple species. Returns the species' mass fractions,...
void calculateMassFractions()
Transfor volume fraction into mass fractions.
virtual void solveYi(PtrList< volScalarField::Internal > &, PtrList< volScalarField::Internal > &)
Solve species fraction equation.
virtual void correct()
Correct phase thermo.
label inertIndex() const
Return inert species index.
autoPtr< phaseThermo > thermoPtr_
Thermophysical model.
PtrList< volScalarField > X_
Ptr list of volumetric fractions for species.
virtual const phaseThermo & thermo() const
Access to thermo.
void calculateVolumeFractions()
Transfor mass fraction into volume fractions.
label inertIndex_
Inert species index.
virtual ~MultiComponentPhaseModel()=default
Destructor.
bool addDiffusion_
Add diffusion transport on Yi's Eq.
hashedWordList species_
Species table.
const hashedWordList & species() const
Species table.
virtual const PtrList< volScalarField > & Y() const
Constant access the species mass fractions.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
A class for handling words, derived from Foam::string.
Definition: word.H:68
Namespace for OpenFOAM.