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-2016 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 vector& position,
67  const label celli,
68  const label tetFacei,
69  const label tetPtI
70 )
71 :
72  ParcelType(mesh, position, 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  const label tetFacei,
86  const label tetPtI,
87  const label typeId,
88  const scalar nParticle0,
89  const scalar d0,
90  const scalar dTarget0,
91  const vector& U0,
92  const vector& f0,
93  const vector& angularMomentum0,
94  const vector& torque0,
95  const scalarField& Y0,
96  const constantProperties& constProps
97 )
98 :
99  ParcelType
100  (
101  mesh,
102  position,
103  celli,
104  tetFacei,
105  tetPtI,
106  typeId,
107  nParticle0,
108  d0,
109  dTarget0,
110  U0,
111  f0,
112  angularMomentum0,
113  torque0,
114  constProps
115  ),
116  mass0_(0.0),
117  Y_(Y0),
118  pc_(0.0)
119 {
120  // Set initial parcel mass
121  mass0_ = this->mass();
122 }
123 
124 
125 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
126 
127 template<class ParcelType>
128 inline Foam::scalar
130 {
131  return pMin_.value();
132 }
133 
134 
135 template<class ParcelType>
136 inline bool
138 {
139  return constantVolume_.value();
140 }
141 
142 
143 // * * * * * * * * * * ThermoParcel Member Functions * * * * * * * * * * * * //
144 
145 template<class ParcelType>
146 inline Foam::scalar Foam::ReactingParcel<ParcelType>::mass0() const
147 {
148  return mass0_;
149 }
150 
151 
152 template<class ParcelType>
154 {
155  return Y_;
156 }
157 
158 
159 template<class ParcelType>
160 inline Foam::scalar Foam::ReactingParcel<ParcelType>::pc() const
161 {
162  return pc_;
163 }
164 
165 
166 template<class ParcelType>
168 {
169  return pc_;
170 }
171 
172 
173 template<class ParcelType>
175 {
176  return mass0_;
177 }
178 
179 
180 template<class ParcelType>
182 {
183  return Y_;
184 }
185 
186 
187 // ************************************************************************* //
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:137
scalarField Y_
Mass fractions of mixture [].
const Type & value() const
Return the value.
const scalarField & Y() const
Return const access to mass fractions of mixture [].
bool constantVolume() const
Return const access to the constant volume flag.
scalar mass0() const
Return const access to initial mass [kg].
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].
ReactingParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
scalar pc() const
Return the owner cell pressure [Pa].
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar mass0_
Initial mass [kg].