Loading...
Searching...
No Matches
inhomogeneousMixture.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::inhomogeneousMixture
29
30Group
31 grpReactionThermophysicalMixtures
32
33Description
34 The inhomogeneous mixture contains species ("ft", "b").
35
36SourceFiles
37 inhomogeneousMixture.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef inhomogeneousMixture_H
42#define inhomogeneousMixture_H
43
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51/*---------------------------------------------------------------------------*\
52 Class inhomogeneousMixture Declaration
53\*---------------------------------------------------------------------------*/
54
55template<class ThermoType>
56class inhomogeneousMixture
57:
59{
60 // Private Data
61
62 dimensionedScalar stoicRatio_;
63
64 ThermoType fuel_;
65 ThermoType oxidant_;
66 ThermoType products_;
67
68 mutable ThermoType mixture_;
69
70 //- Mixture fraction
71 volScalarField& ft_;
72
73 //- Regress variable
75
76 //- No copy construct
77 inhomogeneousMixture
78 (
79 const inhomogeneousMixture<ThermoType>&
80 ) = delete;
81
82
83public:
84
85 //- The type of thermodynamics this mixture is instantiated for
86 typedef ThermoType thermoType;
87
88
89 // Constructors
90
91 //- Construct from dictionary, mesh and phase name
92 inhomogeneousMixture
93 (
95 const fvMesh& mesh,
96 const word& phaseName
97 );
98
99
100 //- Destructor
101 virtual ~inhomogeneousMixture() = default;
102
103
104 // Member functions
105
106 //- Return the instantiated type name
107 static word typeName()
109 return "inhomogeneousMixture<" + ThermoType::typeName() + '>';
110 }
111
112 const dimensionedScalar& stoicRatio() const
113 {
114 return stoicRatio_;
115 }
117 const ThermoType& mixture(const scalar, const scalar) const;
118
119 const ThermoType& cellMixture(const label celli) const
120 {
121 return mixture(ft_[celli], b_[celli]);
122 }
123
124 const ThermoType& cellVolMixture
125 (
126 const scalar p,
127 const scalar T,
128 const label celli
129 ) const
130 {
131 return mixture(ft_[celli], b_[celli]);
132 }
134 const ThermoType& patchFaceMixture
135 (
136 const label patchi,
137 const label facei
138 ) const
139 {
140 return mixture
141 (
142 ft_.boundaryField()[patchi][facei],
143 b_.boundaryField()[patchi][facei]
144 );
145 }
146
147 const ThermoType& patchFaceVolMixture
148 (
149 const scalar p,
150 const scalar T,
151 const label patchi,
152 const label facei
153 ) const
154 {
155 return mixture
157 ft_.boundaryField()[patchi][facei],
158 b_.boundaryField()[patchi][facei]
159 );
160 }
161
162 const ThermoType& cellReactants(const label celli) const
163 {
164 return mixture(ft_[celli], 1);
165 }
166
167 const ThermoType& patchFaceReactants
168 (
169 const label patchi,
170 const label facei
171 ) const
172 {
173 return mixture
174 (
175 ft_.boundaryField()[patchi][facei],
177 );
178 }
179
180 const ThermoType& cellProducts(const label celli) const
181 {
182 return mixture(ft_[celli], 0);
183 }
184
185 const ThermoType& patchFaceProducts
186 (
187 const label patchi,
188 const label facei
189 ) const
190 {
191 return mixture
192 (
193 ft_.boundaryField()[patchi][facei],
195 );
196 }
197
198 //- Read dictionary
199 void read(const dictionary&);
200
201 //- Return thermo based on index
202 const ThermoType& getLocalThermo(const label speciei) const;
203};
204
205
206// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207
208} // End namespace Foam
209
210// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211
212#ifdef NoRepository
213 #include "inhomogeneousMixture.C"
214#endif
215
216// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217
218#endif
219
220// ************************************************************************* //
basicCombustionMixture(const dictionary &thermoDict, const wordList &specieNames, const fvMesh &mesh, const word &phaseName)
Construct from dictionary, specie names, mesh and phase name.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const ThermoType & cellMixture(const label celli) const
const ThermoType & cellReactants(const label celli) const
const ThermoType & cellProducts(const label celli) const
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
const ThermoType & cellVolMixture(const scalar p, const scalar T, const label celli) const
const ThermoType & patchFaceProducts(const label patchi, const label facei) const
static word typeName()
Return the instantiated type name.
const ThermoType & patchFaceMixture(const label patchi, const label facei) const
const dimensionedScalar & stoicRatio() const
virtual ~inhomogeneousMixture()=default
Destructor.
const ThermoType & patchFaceVolMixture(const scalar p, const scalar T, const label patchi, const label facei) const
const ThermoType & patchFaceReactants(const label patchi, const label facei) const
void read(const dictionary &)
Read dictionary.
const ThermoType & getLocalThermo(const label speciei) const
Return thermo based on index.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
const dictionary & thermoDict
Definition EEqn.H:16
dynamicFvMesh & mesh
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture