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