ThermoParcelI.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>
30 :
31  ParcelType::constantProperties(),
32  T0_(this->dict_, 0.0),
33  TMin_(this->dict_, 0.0),
34  TMax_(this->dict_, VGREAT),
35  Cp0_(this->dict_, 0.0),
36  epsilon0_(this->dict_, 0.0),
37  f0_(this->dict_, 0.0)
38 {}
39 
40 
41 template<class ParcelType>
43 (
44  const constantProperties& cp
45 )
46 :
47  ParcelType::constantProperties(cp),
48  T0_(cp.T0_),
49  TMin_(cp.TMin_),
50  TMax_(cp.TMax_),
51  Cp0_(cp.Cp0_),
52  epsilon0_(cp.epsilon0_),
53  f0_(cp.f0_)
54 {}
55 
56 
57 template<class ParcelType>
59 (
60  const dictionary& parentDict
61 )
62 :
63  ParcelType::constantProperties(parentDict),
64  T0_(this->dict_, "T0"),
65  TMin_(this->dict_, "TMin", 200.0),
66  TMax_(this->dict_, "TMax", 5000.0),
67  Cp0_(this->dict_, "Cp0"),
68  epsilon0_(this->dict_, "epsilon0"),
69  f0_(this->dict_, "f0")
70 {}
71 
72 
73 template<class ParcelType>
75 (
76  const polyMesh& mesh,
77  const vector& position,
78  const label celli,
79  const label tetFacei,
80  const label tetPtI
81 )
82 :
83  ParcelType(mesh, position, celli, tetFacei, tetPtI),
84  T_(0.0),
85  Cp_(0.0),
86  Tc_(0.0),
87  Cpc_(0.0)
88 {}
89 
90 
91 template<class ParcelType>
93 (
94  const polyMesh& mesh,
95  const vector& position,
96  const label celli,
97  const label tetFacei,
98  const label tetPtI,
99  const label typeId,
100  const scalar nParticle0,
101  const scalar d0,
102  const scalar dTarget0,
103  const vector& U0,
104  const vector& f0,
105  const vector& angularMomentum0,
106  const vector& torque0,
107  const constantProperties& constProps
108 )
109 :
110  ParcelType
111  (
112  mesh,
113  position,
114  celli,
115  tetFacei,
116  tetPtI,
117  typeId,
118  nParticle0,
119  d0,
120  dTarget0,
121  U0,
122  f0,
123  angularMomentum0,
124  torque0,
125  constProps
126  ),
127  T_(constProps.T0()),
128  Cp_(constProps.Cp0()),
129  Tc_(0.0),
130  Cpc_(0.0)
131 {}
132 
133 
134 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
135 
136 template<class ParcelType>
137 inline Foam::scalar
139 {
140  return T0_.value();
141 }
142 
143 
144 template<class ParcelType>
145 inline Foam::scalar
147 {
148  return TMin_.value();
149 }
150 
151 
152 template<class ParcelType>
153 inline Foam::scalar
155 {
156  return TMax_.value();
157 }
158 
159 
160 template<class ParcelType>
161 inline void
163 {
164  TMax_.setValue(TMax);
165 }
166 
167 
168 template<class ParcelType>
169 inline Foam::scalar
171 {
172  return Cp0_.value();
173 }
174 
175 
176 template<class ParcelType>
177 inline Foam::scalar
179 {
180  return epsilon0_.value();
181 }
182 
183 
184 template<class ParcelType>
185 inline Foam::scalar
187 {
188  return f0_.value();
189 }
190 
191 
192 // * * * * * * * * * * ThermoParcel Member Functions * * * * * * * * * * * * //
193 
194 template<class ParcelType>
195 inline Foam::scalar Foam::ThermoParcel<ParcelType>::T() const
196 {
197  return T_;
198 }
199 
200 
201 template<class ParcelType>
202 inline Foam::scalar Foam::ThermoParcel<ParcelType>::Cp() const
203 {
204  return Cp_;
205 }
206 
207 
208 template<class ParcelType>
209 inline Foam::scalar Foam::ThermoParcel<ParcelType>::hs() const
210 {
211  return Cp_*(T_ - 298.15);
212 }
213 
214 
215 template<class ParcelType>
216 inline Foam::scalar Foam::ThermoParcel<ParcelType>::Tc() const
217 {
218  return Tc_;
219 }
220 
221 
222 template<class ParcelType>
223 inline Foam::scalar Foam::ThermoParcel<ParcelType>::Cpc() const
224 {
225  return Cpc_;
226 }
227 
228 
229 template<class ParcelType>
230 inline Foam::scalar& Foam::ThermoParcel<ParcelType>::T()
231 {
232  return T_;
233 }
234 
235 
236 template<class ParcelType>
238 {
239  return Cp_;
240 }
241 
242 
243 // ************************************************************************* //
scalar Cpc() const
Return const access to carrier specific heat capacity.
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
scalar f0() const
Return const access to the particle scattering factor [].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
scalar Cp_
Specific heat capacity [J/(kg.K)].
Definition: ThermoParcel.H:233
Class to hold thermo particle constant properties.
Definition: ThermoParcel.H:78
scalar TMax() const
Return const access to maximum temperature [K].
void setTMax(const scalar TMax)
Set the maximum temperature [K].
scalar T() const
Return const access to temperature.
scalar TMin() const
Return const access to minimum temperature [K].
const Type & value() const
Return the value.
scalar T0() const
Return const access to the particle initial temperature [K].
scalar T_
Temperature [K].
Definition: ThermoParcel.H:230
ThermoParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
Definition: ThermoParcelI.H:75
scalar Cpc_
Specific heat capacity [J/(kg.K)].
Definition: ThermoParcel.H:242
scalar Tc() const
Return const access to carrier temperature.
void setValue(const Type &value)
Set the value.
scalar Cp0() const
Return const access to the particle specific heat capacity.
scalar hs() const
Return the parcel sensible enthalpy.
scalar Cp() const
Return const access to specific heat capacity.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar epsilon0() const
Return const access to the particle emissivity [].
scalar Tc_
Temperature [K].
Definition: ThermoParcel.H:239