comfort.C
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) 2019 OpenFOAM Foundation
9 Copyright (C) 2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "comfort.H"
30#include "wallFvPatch.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace functionObjects
38{
41}
42}
43
44// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45
46Foam::tmp<Foam::volScalarField> Foam::functionObjects::comfort::magU() const
47{
48 tmp<volScalarField> tmagU = mag(lookupObject<volVectorField>("U"));
49 volScalarField& magU = tmagU.ref();
50
51 // Flag to use the averaged velocity field in the domain.
52 // Consistent with EN ISO 7730 but does not make physical sense
53 if (meanVelocity_)
54 {
55 magU = magU.weightedAverage(mesh_.V());
56 }
57
58 return tmagU;
59}
60
61
62Foam::dimensionedScalar Foam::functionObjects::comfort::Trad() const
63{
64 dimensionedScalar Trad(Trad_);
65
66 // The mean radiation is calculated by the mean wall temperatures
67 // which are summed and divided by the area | only walls are taken into
68 // account. This approach might be correct for a squared room but will
69 // defintely be inconsistent for complex room geometries. The norm does
70 // not provide any information about the calculation of this quantity.
71 if (!TradSet_)
72 {
73 const volScalarField::Boundary& TBf =
74 lookupObject<volScalarField>("T").boundaryField();
75
76 scalar areaIntegral = 0;
77 scalar TareaIntegral = 0;
78
79 forAll(TBf, patchi)
80 {
81 const fvPatchScalarField& pT = TBf[patchi];
82 const fvPatch& pTBf = TBf[patchi].patch();
83 const scalarField& pSf = pTBf.magSf();
84
85 if (isType<wallFvPatch>(pTBf))
86 {
87 areaIntegral += gSum(pSf);
88 TareaIntegral += gSum(pSf*pT);
89 }
90 }
91
92 Trad.value() = TareaIntegral/areaIntegral;
93 }
94
95 // Bounds based on EN ISO 7730
96 if ((Trad.value() < 283.15) || (Trad.value() > 313.15))
97 {
99 << "The calculated mean wall radiation temperature is out of the\n"
100 << "bounds specified in EN ISO 7730:2005\n"
101 << "Valid range is 10 degC < T < 40 degC\n"
102 << "The actual value is: " << (Trad.value() - 273.15) << nl << endl;
103 }
104
105 return Trad;
106}
107
108
109Foam::tmp<Foam::volScalarField> Foam::functionObjects::comfort::pSat() const
110{
111 static const dimensionedScalar kPaToPa(dimPressure, 1000);
112 static const dimensionedScalar A(dimless, 16.6563);
113 static const dimensionedScalar B(dimTemperature, 4030.183);
114 static const dimensionedScalar C(dimTemperature, -38.15);
115
116 tmp<volScalarField> tpSat = volScalarField::New("pSat", mesh_, pSat_);
117
118 // Calculate the saturation pressure if no user input is given
119 if (pSat_.value() == 0)
120 {
121 const auto& T = lookupObject<volScalarField>("T");
122
123 // Equation based on ISO 7730:2006
124 tpSat = kPaToPa*exp(A - B/(T + C));
125 }
126
127 return tpSat;
128}
129
130
131Foam::tmp<Foam::volScalarField> Foam::functionObjects::comfort::Tcloth
132(
133 volScalarField& hc,
134 const dimensionedScalar& metabolicRateSI,
135 const dimensionedScalar& extWorkSI,
136 const volScalarField& T,
137 const dimensionedScalar& Trad
138)
139{
140 const dimensionedScalar factor1(dimTemperature, 308.85);
141
142 const dimensionedScalar factor2
143 (
144 dimTemperature/metabolicRateSI.dimensions(),
145 0.028
146 );
147
148 const dimensionedScalar factor3
149 (
151 3.96e-8
152 );
153
154 // Heat transfer coefficient based on forced convection [W/m^2/K]
155 const volScalarField hcForced
156 (
157 dimensionedScalar(hc.dimensions()/sqrt(dimVelocity), 12.1)
158 *sqrt(magU())
159 );
160
161 // Tcl [K] (surface cloth temperature)
162 tmp<volScalarField> tTcl
163 (
165 (
166 "Tcl",
167 T.mesh(),
169 )
170 );
171 volScalarField& Tcl = tTcl.ref();
172
173 // Initial guess
174 Tcl = T;
175
176 label i = 0;
177
178 Tcl.storePrevIter();
179
180 // Same temperatures as for the radiation
181 const dimensionedScalar Tmin(dimTemperature, 283.15);
182 const dimensionedScalar Tmax(dimTemperature, 313.15);
183
184 // Iterative solving of equation (2)
185 do
186 {
187 Tcl = (Tcl + Tcl.prevIter())/2;
188 Tcl.storePrevIter();
189
190 // Heat transfer coefficient based on natural convection
191 volScalarField hcNatural
192 (
193 dimensionedScalar(hc.dimensions()/pow025(dimTemperature), 2.38)
194 *pow025(mag(Tcl - T))
195 );
196
197 // Set heat transfer coefficient based on equation (3)
198 hc =
199 pos(hcForced - hcNatural)*hcForced
200 + neg0(hcForced - hcNatural)*hcNatural;
201
202 // Calculate surface temperature based on equation (2)
203 Tcl =
204 factor1
205 - factor2*(metabolicRateSI - extWorkSI)
206 - Icl_*factor3*fcl_*(pow4(Tcl) - pow4(Trad))
207 - Icl_*fcl_*hc*(Tcl - T);
208
209 // Make sure that Tcl is in some physical limit (same range as we used
210 // for the radiative estimation - based on ISO EN 7730:2005)
211 Tcl.clip(Tmin, Tmax);
212
213 } while (!converged(Tcl) && i++ < maxClothIter_);
214
215 if (i == maxClothIter_)
216 {
218 << "The surface cloth temperature did not converge within " << i
219 << " iterations" << nl;
220 }
221
222 return tTcl;
223}
224
225
227(
228 const volScalarField& phi
229) const
230{
231 return
232 max(mag(phi.primitiveField() - phi.prevIter().primitiveField()))
233 < tolerance_;
234}
235
236
237// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
238
240(
241 const word& name,
242 const Time& runTime,
243 const dictionary& dict
244)
245:
247 clothing_("clothing", dimless, 0),
248 metabolicRate_("metabolicRate", dimMass/pow3(dimTime), 0.8),
249 extWork_("extWork", dimMass/pow3(dimTime), 0),
250 Trad_("Trad", dimTemperature, 0),
251 relHumidity_("relHumidity", dimless, 0.5),
252 pSat_("pSat", dimPressure, 0),
253 Icl_("Icl", pow3(dimTime)*dimTemperature/dimMass, 0),
254 fcl_("fcl", dimless, 0),
255 tolerance_(1e-4),
256 maxClothIter_(100),
257 TradSet_(false),
258 meanVelocity_(false)
259{
260 read(dict);
261}
262
263
264// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
265
267{
269 {
270 clothing_.readIfPresent(dict);
271 metabolicRate_.readIfPresent(dict);
272 extWork_.readIfPresent(dict);
273 pSat_.readIfPresent(dict);
274 tolerance_ = dict.getOrDefault("tolerance", 1e-4);
275 maxClothIter_ = dict.getOrDefault("maxClothIter", 100);
276 meanVelocity_ = dict.getOrDefault<bool>("meanVelocity", false);
277
278 // Read relative humidity if provided and convert from % to fraction
279 if (dict.found(relHumidity_.name()))
280 {
281 relHumidity_.read(dict);
282 relHumidity_ /= 100;
283 }
284
285 // Read radiation temperature if provided
286 if (dict.found(Trad_.name()))
287 {
288 TradSet_ = true;
289 Trad_.read(dict);
290 }
291
292 Icl_ = dimensionedScalar(Icl_.dimensions(), 0.155)*clothing_;
293
294 fcl_.value() =
295 Icl_.value() <= 0.078
296 ? 1.0 + 1.290*Icl_.value()
297 : 1.05 + 0.645*Icl_.value();
298
299 return true;
300 }
301
302 return false;
303}
304
305
307{
308 // Assign and build fields
309 const dimensionedScalar Trad(this->Trad());
310 const volScalarField pSat(this->pSat());
311
312 const dimensionedScalar metabolicRateSI(58.15*metabolicRate_);
313 const dimensionedScalar extWorkSI(58.15*extWork_);
314
315 const auto& T = lookupObject<volScalarField>("T");
316
317 // Heat transfer coefficient [W/m^2/K]
318 // This field is updated in Tcloth()
320 (
322 (
323 "hc",
324 mesh_.time().timeName(),
325 mesh_
326 ),
327 mesh_,
329 );
330
331 // Calculate the surface temperature of the cloth by an iterative
332 // process using equation (2) from DIN EN ISO 7730 [degC]
333 const volScalarField Tcloth
334 (
335 this->Tcloth
336 (
337 hc,
338 metabolicRateSI,
339 extWorkSI,
340 T,
341 Trad
342 )
343 );
344
345 // Calculate the PMV quantity
346 const dimensionedScalar factor1(pow3(dimTime)/dimMass, 0.303);
347 const dimensionedScalar factor2
348 (
349 dimless/metabolicRateSI.dimensions(),
350 -0.036
351 );
352 const dimensionedScalar factor3(factor1.dimensions(), 0.028);
353 const dimensionedScalar factor4(dimLength/dimTime, 3.05e-3);
354 const dimensionedScalar factor5(dimPressure, 5733);
355 const dimensionedScalar factor6(dimTime/dimLength, 6.99);
356 const dimensionedScalar factor8(metabolicRateSI.dimensions(), 58.15);
357 const dimensionedScalar factor9(dimless/dimPressure, 1.7e-5);
358 const dimensionedScalar factor10(dimPressure, 5867);
359 const dimensionedScalar factor11(dimless/dimTemperature, 0.0014);
360 const dimensionedScalar factor12(dimTemperature, 307.15);
361 const dimensionedScalar factor13
362 (
364 3.96e-8
365 );
366
367 const scalar factor7
368 (
369 // Special treatment of Term4
370 // if metaRate - extWork < factor8, set to zero
371 (metabolicRateSI - extWorkSI).value() < factor8.value() ? 0 : 0.42
372 );
373
374 Info<< "Calculating the predicted mean vote (PMV)" << endl;
375
376 // Equation (1)
378 (
379 // Term1: Thermal sensation transfer coefficient
380 (factor1*exp(factor2*metabolicRateSI) + factor3)
381 *(
382 (metabolicRateSI - extWorkSI)
383
384 // Term2: Heat loss difference through skin
385 - factor4
386 *(
387 factor5
388 - factor6*(metabolicRateSI - extWorkSI)
389 - pSat*relHumidity_
390 )
391
392 // Term3: Heat loss through sweating
393 - factor7*(metabolicRateSI - extWorkSI - factor8)
394
395 // Term4: Heat loss through latent respiration
396 - factor9*metabolicRateSI*(factor10 - pSat*relHumidity_)
397
398 // Term5: Heat loss through dry respiration
399 - factor11*metabolicRateSI*(factor12 - T)
400
401 // Term6: Heat loss through radiation
402 - factor13*fcl_*(pow4(Tcloth) - pow4(Trad))
403
404 // Term7: Heat loss through convection
405 - fcl_*hc*(Tcloth - T)
406 )
407 );
408
409 Info<< "Calculating the predicted percentage of dissatisfaction (PPD)"
410 << endl;
411
412 // Equation (5)
414 100 - 95*exp(-0.03353*pow4(PMV()) - 0.21790*sqr(PMV()));
415
416 Info<< "Calculating the draught rating (DR)\n";
417
418 const dimensionedScalar Umin(dimVelocity, 0.05);
419 const dimensionedScalar Umax(dimVelocity, 0.5);
420 const dimensionedScalar pre(dimless, 0.37);
421 const dimensionedScalar C1(dimVelocity, 3.14);
422
423 // Limit the velocity field to the values given in EN ISO 7733
424 volScalarField Umag(mag(lookupObject<volVectorField>("U")));
425 Umag.clip(Umin, Umax);
426
427 // Calculate the turbulent intensity if turbulent kinetic energy field k
428 // exists
430 (
432 (
433 "TI",
434 mesh_.time().timeName(),
435 mesh_
436 ),
437 mesh_,
439 );
440
441 if (foundObject<volScalarField>("k"))
442 {
443 const auto& k = lookupObject<volScalarField>("k");
444 TI = sqrt(2/3*k)/Umag;
445 }
446
447 // For unit correctness
448 const dimensionedScalar correctUnit
449 (
450 dimensionSet(0, -1.62, 1.62, -1, 0, 0, 0),
451 1
452 );
453
454 // Equation (6)
456 correctUnit*(factor12 - T)*pow(Umag - Umin, 0.62)*(pre*Umag*TI + C1);
457
458 // Calculate the operative temperature
459 tmp<volScalarField> Top = 0.5*(T + Trad);
460
461 // Need modifiable field names:
462 word fieldNamePMV = "PMV";
463 word fieldNamePPD = "PPD";
464 word fieldNameDR = "DR";
465 word fieldNameTop = "Top";
466
467 return
468 (
469 store(fieldNamePMV, PMV)
470 && store(fieldNamePPD, PPD)
471 && store(fieldNameDR, DR)
472 && store(fieldNameTop, Top)
473 );
474}
475
476
478{
479 return
480 (
481 writeObject("PMV")
482 && writeObject("PPD")
483 && writeObject("DR")
484 && writeObject("Top")
485 );
486}
487
488
489// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
label k
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
void clip(const dimensioned< MinMax< Type > > &range)
Clip the field to be bounded within the specified range.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual bool read()
Re-read model coefficients if they have changed.
bool converged() const noexcept
Has the solver converged?
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
const dimensionSet & dimensions() const
Return const reference to dimensions.
const Type & value() const
Return const reference to value.
Abstract base-class for Time/database function objects.
Calculates the thermal comfort quantities predicted mean vote (PMV), predicted percentage of dissatis...
Definition: comfort.H:247
virtual bool write()
Write the PPD and PMV fields.
Definition: comfort.C:477
virtual bool read(const dictionary &)
Read the data needed for the comfort calculation.
Definition: comfort.C:266
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes the magnitude of an input field.
Definition: mag.H:153
Computes the power of an input volScalarField.
Definition: pow.H:217
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
engineTime & runTime
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
const dimensionSet dimPressure
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensionedScalar neg0(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
fvPatchField< scalar > fvPatchScalarField
dimensionedScalar pow025(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & C
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333