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-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>
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 barycentric& coordinates,
78  const label celli,
79  const label tetFacei,
80  const label tetPti
81 )
82 :
83  ParcelType(mesh, coordinates, 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 )
98 :
99  ParcelType(mesh, position, celli),
100  T_(0.0),
101  Cp_(0.0),
102  Tc_(0.0),
103  Cpc_(0.0)
104 {}
105 
106 
107 template<class ParcelType>
109 (
110  const polyMesh& mesh,
111  const barycentric& coordinates,
112  const label celli,
113  const label tetFacei,
114  const label tetPti,
115  const label typeId,
116  const scalar nParticle0,
117  const scalar d0,
118  const scalar dTarget0,
119  const vector& U0,
120  const vector& f0,
121  const vector& angularMomentum0,
122  const vector& torque0,
123  const constantProperties& constProps
124 )
125 :
126  ParcelType
127  (
128  mesh,
129  coordinates,
130  celli,
131  tetFacei,
132  tetPti,
133  typeId,
134  nParticle0,
135  d0,
136  dTarget0,
137  U0,
138  f0,
139  angularMomentum0,
140  torque0,
141  constProps
142  ),
143  T_(constProps.T0()),
144  Cp_(constProps.Cp0()),
145  Tc_(0.0),
146  Cpc_(0.0)
147 {}
148 
149 
150 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
151 
152 template<class ParcelType>
153 inline Foam::scalar
155 {
156  return T0_.value();
157 }
158 
159 
160 template<class ParcelType>
161 inline Foam::scalar
163 {
164  return TMin_.value();
165 }
166 
167 
168 template<class ParcelType>
169 inline Foam::scalar
171 {
172  return TMax_.value();
173 }
174 
175 
176 template<class ParcelType>
177 inline void
179 {
180  TMax_.setValue(TMax);
181 }
182 
183 
184 template<class ParcelType>
185 inline Foam::scalar
187 {
188  return Cp0_.value();
189 }
190 
191 
192 template<class ParcelType>
193 inline Foam::scalar
195 {
196  return epsilon0_.value();
197 }
198 
199 
200 template<class ParcelType>
201 inline Foam::scalar
203 {
204  return f0_.value();
205 }
206 
207 
208 // * * * * * * * * * * ThermoParcel Member Functions * * * * * * * * * * * * //
209 
210 template<class ParcelType>
211 inline Foam::scalar Foam::ThermoParcel<ParcelType>::T() const
212 {
213  return T_;
214 }
215 
216 
217 template<class ParcelType>
218 inline Foam::scalar Foam::ThermoParcel<ParcelType>::Cp() const
219 {
220  return Cp_;
221 }
222 
223 
224 template<class ParcelType>
225 inline Foam::scalar Foam::ThermoParcel<ParcelType>::hs() const
226 {
227  return Cp_*(T_ - 298.15);
228 }
229 
230 
231 template<class ParcelType>
232 inline Foam::scalar Foam::ThermoParcel<ParcelType>::Tc() const
233 {
234  return Tc_;
235 }
236 
237 
238 template<class ParcelType>
239 inline Foam::scalar Foam::ThermoParcel<ParcelType>::Cpc() const
240 {
241  return Cpc_;
242 }
243 
244 
245 template<class ParcelType>
246 inline Foam::scalar& Foam::ThermoParcel<ParcelType>::T()
247 {
248  return T_;
249 }
250 
251 
252 template<class ParcelType>
254 {
255  return Cp_;
256 }
257 
258 
259 // ************************************************************************* //
scalar T0() const
Return const access to the particle initial temperature [K].
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 TMin() const
Return const access to minimum temperature [K].
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
scalar T() const
Return const access to temperature.
Class to hold thermo particle constant properties.
Definition: ThermoParcel.H:78
void setTMax(const scalar TMax)
Set the maximum temperature [K].
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Definition: ThermoParcelI.H:75
scalar f0() const
Return const access to the particle scattering factor [].
const Type & value() const
Return the value.
scalar Cp() const
Return const access to specific heat capacity.
scalar T_
Temperature [K].
Definition: ThermoParcel.H:230
scalar hs() const
Return the parcel sensible enthalpy.
scalar Cp0() const
Return const access to the particle specific heat capacity.
scalar Cpc_
Specific heat capacity [J/(kg.K)].
Definition: ThermoParcel.H:242
scalar Tc() const
Return const access to carrier temperature.
scalar Cpc() const
Return const access to carrier specific heat capacity.
void setValue(const Type &value)
Set the value.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar TMax() const
Return const access to maximum temperature [K].
scalar Tc_
Temperature [K].
Definition: ThermoParcel.H:239
scalar epsilon0() const
Return const access to the particle emissivity [].