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 -------------------------------------------------------------------------------
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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 Foam::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 
62 Foam::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 - 273.15 << nl << endl;
103  }
104 
105  return Trad;
106 }
107 
108 
109 Foam::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 
131 Foam::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 
226 bool Foam::functionObjects::comfort::converged
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()
319  volScalarField hc
320  (
321  IOobject
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)
377  tmp<volScalarField> PMV =
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)
413  tmp<volScalarField> PPD =
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
429  volScalarField TI
430  (
431  IOobject
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  // Workaround
459  word fieldNamePMV = "PMV";
460  word fieldNamePPD = "PPD";
461  word fieldNameDR = "DR";
462 
463  return
464  store(fieldNamePMV, PMV)
465  && store(fieldNamePPD, PPD)
466  && store(fieldNameDR, DR);
467 }
468 
470 {
471  return
472  writeObject("PMV")
473  && writeObject("PPD")
474  && writeObject("DR");
475 }
476 
477 
478 // ************************************************************************* //
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:40
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::dimPressure
const dimensionSet dimPressure
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::neg0
dimensionedScalar neg0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:210
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
wallFvPatch.H
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
Foam::GeometricField::clip
void clip(const dimensioned< MinMax< Type >> &range)
Clip the field to be bounded within the specified range.
Definition: GeometricField.C:1154
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
C
volScalarField & C
Definition: readThermalProperties.H:102
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::GeometricField< scalar, fvPatchField, volMesh >::New
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &ds, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field from name, mesh, dimensions and patch type.
Definition: GeometricFieldNew.C:34
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::functionObjects::comfort::execute
virtual bool execute()
Definition: comfort.C:306
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::functionObjects::comfort::read
virtual bool read(const dictionary &)
Read the data needed for the comfort calculation.
Definition: comfort.C:266
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
comfort.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::functionObjects::comfort::comfort
comfort(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: comfort.C:240
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::functionObjects::comfort::write
virtual bool write()
Write the PPD and PMV fields.
Definition: comfort.C:469
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::DimensionedField::weightedAverage
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
Definition: DimensionedField.C:461
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::functionObjects::fvMeshFunctionObject::mesh_
const fvMesh & mesh_
Reference to the fvMesh.
Definition: fvMeshFunctionObject.H:73
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:420
Foam::GeometricField< scalar, fvPatchField, volMesh >
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:179
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177