ReactingParcelI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
28 template<class ParcelType>
29 inline
31 :
32  ParcelType::constantProperties(),
33  pMin_(this->dict_, 0.0),
34  constantVolume_(this->dict_, false)
35 {}
36 
37 
38 template<class ParcelType>
40 (
41  const constantProperties& cp
42 )
43 :
44  ParcelType::constantProperties(cp),
45  pMin_(cp.pMin_),
46  constantVolume_(cp.constantVolume_)
47 {}
48 
49 
50 template<class ParcelType>
52 (
53  const dictionary& parentDict
54 )
55 :
56  ParcelType::constantProperties(parentDict),
57  pMin_(this->dict_, "pMin", 1000.0),
58  constantVolume_(this->dict_, word("constantVolume"))
59 {}
60 
61 
62 template<class ParcelType>
64 (
65  const polyMesh& mesh,
66  const barycentric& coordinates,
67  const label celli,
68  const label tetFacei,
69  const label tetPti
70 )
71 :
72  ParcelType(mesh, coordinates, celli, tetFacei, tetPti),
73  mass0_(0.0),
74  Y_(0)
75 {}
76 
77 
78 template<class ParcelType>
80 (
81  const polyMesh& mesh,
82  const vector& position,
83  const label celli
84 )
85 :
86  ParcelType(mesh, position, celli),
87  mass0_(0.0),
88  Y_(0)
89 {}
90 
91 
92 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
93 
94 template<class ParcelType>
95 inline Foam::scalar
97 {
98  return pMin_.value();
99 }
100 
101 
102 template<class ParcelType>
103 inline bool
105 {
106  return constantVolume_.value();
107 }
108 
109 
110 // * * * * * * * * * * ThermoParcel Member Functions * * * * * * * * * * * * //
111 
112 template<class ParcelType>
113 inline Foam::scalar Foam::ReactingParcel<ParcelType>::mass0() const
114 {
115  return mass0_;
116 }
117 
118 
119 template<class ParcelType>
121 {
122  return Y_;
123 }
124 
125 
126 template<class ParcelType>
128 {
129  return mass0_;
130 }
131 
132 
133 template<class ParcelType>
135 {
136  return Y_;
137 }
138 
139 
140 // ************************************************************************* //
Class to hold reacting parcel constant properties.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
scalarField Y_
Mass fractions of mixture [].
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
bool constantVolume() const
Return const access to the constant volume flag.
const Type & value() const
Return the value.
A class for handling words, derived from string.
Definition: word.H:59
scalar pMin() const
Return const access to the minimum pressure.
scalar mass0() const
Return const access to initial mass [kg].
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar mass0_
Initial mass [kg].
const scalarField & Y() const
Return const access to mass fractions of mixture [].