ThermoParcelI.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-2021 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 {}
87 
88 
89 template<class ParcelType>
91 (
92  const polyMesh& mesh,
93  const vector& position,
94  const label celli
95 )
96 :
97  ParcelType(mesh, position, celli),
98  T_(0.0),
99  Cp_(0.0)
100 {}
101 
102 
103 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
104 
105 template<class ParcelType>
106 inline Foam::scalar
108 {
109  return T0_.value();
110 }
111 
112 
113 template<class ParcelType>
114 inline Foam::scalar
116 {
117  return TMin_.value();
118 }
119 
120 
121 template<class ParcelType>
122 inline Foam::scalar
124 {
125  return TMax_.value();
126 }
127 
128 
129 template<class ParcelType>
130 inline void
132 {
133  TMax_.setValue(TMax);
134 }
135 
136 
137 template<class ParcelType>
138 inline Foam::scalar
140 {
141  return Cp0_.value();
142 }
143 
144 
145 template<class ParcelType>
146 inline Foam::scalar
148 {
149  return epsilon0_.value();
150 }
151 
152 
153 template<class ParcelType>
154 inline Foam::scalar
156 {
157  return f0_.value();
158 }
159 
160 
161 // * * * * * * * * * * ThermoParcel Member Functions * * * * * * * * * * * * //
162 
163 template<class ParcelType>
164 inline Foam::scalar Foam::ThermoParcel<ParcelType>::T() const
165 {
166  return T_;
167 }
168 
169 
170 template<class ParcelType>
171 inline Foam::scalar Foam::ThermoParcel<ParcelType>::Cp() const
172 {
173  return Cp_;
174 }
175 
176 
177 template<class ParcelType>
178 inline Foam::scalar& Foam::ThermoParcel<ParcelType>::T()
179 {
180  return T_;
181 }
182 
183 
184 template<class ParcelType>
186 {
187  return Cp_;
188 }
189 
190 
191 // ************************************************************************* //
scalar T0() const
Return const access to the particle initial temperature [K].
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:156
scalar Cp_
Specific heat capacity [J/kg/K].
Definition: ThermoParcel.H:271
scalar T() const
Return const access to temperature.
Class to hold thermo particle constant properties.
Definition: ThermoParcel.H:87
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:268
scalar Cp0() const
Return const access to the particle specific heat capacity.
void setValue(const Type &value)
Set the value.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
scalar TMax() const
Return const access to maximum temperature [K].
scalar epsilon0() const
Return const access to the particle emissivity [].