ReactingParcelI.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 | Copyright (C) 2011-2017 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  pc_(0.0)
76 {}
77 
78 
79 template<class ParcelType>
81 (
82  const polyMesh& mesh,
83  const vector& position,
84  const label celli
85 )
86 :
87  ParcelType(mesh, position, celli),
88  mass0_(0.0),
89  Y_(0),
90  pc_(0.0)
91 {}
92 
93 
94 template<class ParcelType>
96 (
97  const polyMesh& mesh,
98  const barycentric& coordinates,
99  const label celli,
100  const label tetFacei,
101  const label tetPti,
102  const label typeId,
103  const scalar nParticle0,
104  const scalar d0,
105  const scalar dTarget0,
106  const vector& U0,
107  const vector& f0,
108  const vector& angularMomentum0,
109  const vector& torque0,
110  const scalarField& Y0,
111  const constantProperties& constProps
112 )
113 :
114  ParcelType
115  (
116  mesh,
117  coordinates,
118  celli,
119  tetFacei,
120  tetPti,
121  typeId,
122  nParticle0,
123  d0,
124  dTarget0,
125  U0,
126  f0,
127  angularMomentum0,
128  torque0,
129  constProps
130  ),
131  mass0_(0.0),
132  Y_(Y0),
133  pc_(0.0)
134 {
135  // Set initial parcel mass
136  mass0_ = this->mass();
137 }
138 
139 
140 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
141 
142 template<class ParcelType>
143 inline Foam::scalar
145 {
146  return pMin_.value();
147 }
148 
149 
150 template<class ParcelType>
151 inline bool
153 {
154  return constantVolume_.value();
155 }
156 
157 
158 // * * * * * * * * * * ThermoParcel Member Functions * * * * * * * * * * * * //
159 
160 template<class ParcelType>
161 inline Foam::scalar Foam::ReactingParcel<ParcelType>::mass0() const
162 {
163  return mass0_;
164 }
165 
166 
167 template<class ParcelType>
169 {
170  return Y_;
171 }
172 
173 
174 template<class ParcelType>
175 inline Foam::scalar Foam::ReactingParcel<ParcelType>::pc() const
176 {
177  return pc_;
178 }
179 
180 
181 template<class ParcelType>
183 {
184  return pc_;
185 }
186 
187 
188 template<class ParcelType>
190 {
191  return mass0_;
192 }
193 
194 
195 template<class ParcelType>
197 {
198  return Y_;
199 }
200 
201 
202 // ************************************************************************* //
Class to hold reacting parcel constant properties.
scalar pc() const
Return the owner cell pressure [Pa].
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:137
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 pc_
Pressure [Pa].
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 [].