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-2023 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 #include "ThermoParcel.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class ParcelType>
32 :
33  ParcelType::constantProperties(),
34  T0_(this->dict_, 0.0),
35  TMin_(this->dict_, 0.0),
36  TMax_(this->dict_, vGreat),
37  Cp0_(this->dict_, 0.0),
38  epsilon0_(this->dict_, 0.0),
39  f0_(this->dict_, 0.0)
40 {}
41 
42 
43 template<class ParcelType>
45 (
46  const constantProperties& cp
47 )
48 :
49  ParcelType::constantProperties(cp),
50  T0_(cp.T0_),
51  TMin_(cp.TMin_),
52  TMax_(cp.TMax_),
53  Cp0_(cp.Cp0_),
54  epsilon0_(cp.epsilon0_),
55  f0_(cp.f0_)
56 {}
57 
58 
59 template<class ParcelType>
61 (
62  const dictionary& parentDict
63 )
64 :
65  ParcelType::constantProperties(parentDict),
66  T0_(this->dict_, "T0"),
67  TMin_(this->dict_, "TMin", 200.0),
68  TMax_(this->dict_, "TMax", 5000.0),
69  Cp0_(this->dict_, "Cp0"),
70  epsilon0_(this->dict_, "epsilon0"),
71  f0_(this->dict_, "f0")
72 {}
73 
74 
75 template<class ParcelType>
77 (
78  const polyMesh& mesh,
79  const barycentric& coordinates,
80  const label celli,
81  const label tetFacei,
82  const label tetPti,
83  const label facei
84 )
85 :
86  ParcelType(mesh, coordinates, celli, tetFacei, tetPti, facei),
87  T_(0.0),
88  Cp_(0.0)
89 {}
90 
91 
92 template<class ParcelType>
94 (
95  const polyMesh& mesh,
96  const vector& position,
97  const label celli,
98  label& nLocateBoundaryHits
99 )
100 :
101  ParcelType(mesh, position, celli, nLocateBoundaryHits),
102  T_(0.0),
103  Cp_(0.0)
104 {}
105 
106 
107 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
108 
109 template<class ParcelType>
110 inline Foam::scalar
112 {
113  return T0_.value();
114 }
115 
116 
117 template<class ParcelType>
118 inline Foam::scalar
120 {
121  return TMin_.value();
122 }
123 
124 
125 template<class ParcelType>
126 inline Foam::scalar
128 {
129  return TMax_.value();
130 }
131 
132 
133 template<class ParcelType>
134 inline void
136 {
137  TMax_.setValue(TMax);
138 }
139 
140 
141 template<class ParcelType>
142 inline Foam::scalar
144 {
145  return Cp0_.value();
146 }
147 
148 
149 template<class ParcelType>
150 inline Foam::scalar
152 {
153  return epsilon0_.value();
154 }
155 
156 
157 template<class ParcelType>
158 inline Foam::scalar
160 {
161  return f0_.value();
162 }
163 
164 
165 // * * * * * * * * * * ThermoParcel Member Functions * * * * * * * * * * * * //
166 
167 template<class ParcelType>
168 inline Foam::scalar Foam::ThermoParcel<ParcelType>::T() const
169 {
170  return T_;
171 }
172 
173 
174 template<class ParcelType>
175 inline Foam::scalar Foam::ThermoParcel<ParcelType>::Cp() const
176 {
177  return Cp_;
178 }
179 
180 
181 template<class ParcelType>
182 inline Foam::scalar& Foam::ThermoParcel<ParcelType>::T()
183 {
184  return T_;
185 }
186 
187 
188 template<class ParcelType>
190 {
191  return Cp_;
192 }
193 
194 
195 // ************************************************************************* //
Class to hold thermo particle constant properties.
Definition: ThermoParcel.H:90
void setTMax(const scalar TMax)
Set the maximum temperature [K].
scalar f0() const
Return const access to the particle scattering factor [].
scalar TMin() const
Return const access to minimum temperature [K].
scalar epsilon0() const
Return const access to the particle emissivity [].
scalar Cp0() const
Return const access to the particle specific heat capacity.
scalar TMax() const
Return const access to maximum temperature [K].
scalar T0() const
Return const access to the particle initial temperature [K].
scalar Cp() const
Return const access to specific heat capacity.
scalar T() const
Return const access to temperature.
scalar Cp_
Specific heat capacity [J/kg/K].
Definition: ThermoParcel.H:271
scalar T_
Temperature [K].
Definition: ThermoParcel.H:268
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
Definition: ThermoParcelI.H:77
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
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
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:753