Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29
30template<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
43template<class ParcelType>
45(
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
59template<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
75template<class ParcelType>
77(
78 const polyMesh& mesh,
80 const label celli,
81 const label tetFacei,
82 const label tetPti
83)
84:
85 ParcelType(mesh, coordinates, celli, tetFacei, tetPti),
86 T_(0.0),
87 Cp_(0.0)
88{}
89
90
91template<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{}
103
104
105template<class ParcelType>
107(
108 const polyMesh& mesh,
110 const label celli,
111 const label tetFacei,
112 const label tetPti,
113 const label typeId,
114 const scalar nParticle0,
115 const scalar d0,
116 const scalar dTarget0,
117 const vector& U0,
118 const vector& f0,
119 const vector& angularMomentum0,
120 const vector& torque0,
121 const constantProperties& constProps
122)
123:
124 ParcelType
125 (
126 mesh,
128 celli,
129 tetFacei,
130 tetPti,
131 typeId,
132 nParticle0,
133 d0,
134 dTarget0,
135 U0,
136 f0,
137 angularMomentum0,
138 torque0,
139 constProps
140 ),
141 T_(constProps.T0()),
142 Cp_(constProps.Cp0())
144
145
146// * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
148template<class ParcelType>
149inline Foam::scalar
152 return T0_.value();
153}
154
155
156template<class ParcelType>
157inline Foam::scalar
160 return TMin_.value();
161}
163
164template<class ParcelType>
165inline Foam::scalar
168 return TMax_.value();
170
171
172template<class ParcelType>
173inline void
176 TMax_.setValue(TMax);
177}
178
179
180template<class ParcelType>
181inline Foam::scalar
184 return Cp0_.value();
185}
186
187
188template<class ParcelType>
189inline Foam::scalar
192 return epsilon0_.value();
193}
194
195
196template<class ParcelType>
197inline Foam::scalar
199{
200 return f0_.value();
201}
202
203
204// * * * * * * * * * * ThermoParcel Member Functions * * * * * * * * * * * * //
205
206template<class ParcelType>
207inline Foam::scalar Foam::ThermoParcel<ParcelType>::T() const
208{
209 return T_;
210}
211
212
213template<class ParcelType>
214inline Foam::scalar Foam::ThermoParcel<ParcelType>::Cp() const
215{
216 return Cp_;
217}
218
219
220template<class ParcelType>
221inline Foam::scalar Foam::ThermoParcel<ParcelType>::hs() const
222{
223 return Cp_*(T_ - 298.15);
224}
225
226
227template<class ParcelType>
229{
230 return T_;
231}
232
233
234template<class ParcelType>
235inline Foam::scalar& Foam::ThermoParcel<ParcelType>::Cp()
236{
237 return Cp_;
238}
239
240
241// ************************************************************************* //
Class to hold thermo particle constant properties.
void setTMax(const scalar TMax)
Set the maximum temperature [K].
scalar TMin() const
Return const access to minimum temperature [K].
scalar epsilon0() const
Return const access to the particle emissivity [].
scalar TMax() const
Return const access to maximum 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)].
scalar hs() const
Return the parcel sensible enthalpy.
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
scalar T_
Temperature [K].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition particleI.H:116
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition POSIX.C:1065
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Vector< scalar > vector
Definition vector.H:57