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) 2015-2018 OpenFOAM Foundation
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
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47
48/*---------------------------------------------------------------------------*\
49 Class MultiComponentPhaseModel Declaration
50\*---------------------------------------------------------------------------*/
51
52template<class BasePhaseModel>
53class MultiComponentPhaseModel
54:
55 public BasePhaseModel
56{
57protected:
58
59 // Protected data
60
61 //- Turbulent Schmidt number
63
64 //- Residual phase fraction
66
67 //- Inert species index
68 label inertIndex_;
69
70 //- Pointer list to active species
72
73
74public:
75
76 // Constructors
77
79 (
80 const phaseSystem& fluid,
81 const word& phaseName,
82 const label index
83 );
84
85
86 //- Destructor
88
89
90 // Member Functions
91
92 //- Correct the thermodynamics
93 virtual void correctThermo();
94
95 // Species
96
97 //- Return whether the phase is pure (i.e., not multi-component)
98 virtual bool pure() const;
99
100 //- Return the species fraction equation
102
103 //- Return the species mass fractions
104 virtual const PtrList<volScalarField>& Y() const;
105
106 //- Return a species mass fraction by name
107 virtual const volScalarField& Y(const word& name) const;
108
109 //- Access the species mass fractions
110 virtual PtrList<volScalarField>& YRef();
111
112 //- Return the active species mass fractions
113 virtual const UPtrList<volScalarField>& YActive() const;
114
115 //- Access the active species mass fractions
117};
118
119
120// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121
122} // End namespace Foam
123
124// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125
126#ifdef NoRepository
127 #include "MultiComponentPhaseModel.C"
128#endif
129
130// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131
132#endif
133
134// ************************************************************************* //
twoPhaseSystem & fluid
Class which represents a phase with multiple species. Returns the species' mass fractions,...
virtual ~MultiComponentPhaseModel()
Destructor.
virtual void correctThermo()
Correct the thermodynamics.
virtual const PtrList< volScalarField > & Y() const
Return the species mass fractions.
virtual const UPtrList< volScalarField > & YActive() const
Return the active species mass fractions.
virtual PtrList< volScalarField > & YRef()
Access the species mass fractions.
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
virtual bool pure() const
Return whether the phase is pure (i.e., not multi-component)
dimensionedScalar residualAlpha_
Residual phase fraction.
label inertIndex_
Inert species index.
dimensionedScalar Sct_
Turbulent Schmidt number.
virtual ~MultiComponentPhaseModel()=default
Destructor.
UPtrList< volScalarField > YActive_
Pointer list to active species.
virtual UPtrList< volScalarField > & YActiveRef()
Access the active species mass fractions.
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
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:76
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
Namespace for OpenFOAM.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59