externalWallHeatFluxTemperatureFvPatchScalarField.C
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-2019 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 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
31 
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  template<>
39  const char*
40  NamedEnum
41  <
43  3
44  >::names[] =
45  {
46  "power",
47  "flux",
48  "coefficient"
49  };
50 }
51 
52 const Foam::NamedEnum
53 <
55  3
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
63 (
64  const fvPatch& p,
66 )
67 :
68  mixedFvPatchScalarField(p, iF),
69  temperatureCoupledBase(patch(), "undefined", "undefined", "undefined-K"),
70  mode_(fixedHeatFlux),
71  Q_(0),
72  Ta_(),
73  relaxation_(1),
74  emissivity_(0),
75  qrRelaxation_(1),
76  qrName_("undefined-qr"),
77  thicknessLayers_(),
78  kappaLayers_()
79 {
80  refValue() = 0;
81  refGrad() = 0;
82  valueFraction() = 1;
83 }
84 
85 
88 (
89  const fvPatch& p,
91  const dictionary& dict
92 )
93 :
94  mixedFvPatchScalarField(p, iF),
95  temperatureCoupledBase(patch(), dict),
96  mode_(operationModeNames.read(dict.lookup("mode"))),
97  Q_(0),
98  Ta_(),
99  relaxation_(dict.lookupOrDefault<scalar>("relaxation", 1)),
100  emissivity_(dict.lookupOrDefault<scalar>("emissivity", 0)),
101  qrRelaxation_(dict.lookupOrDefault<scalar>("qrRelaxation", 1)),
102  qrName_(dict.lookupOrDefault<word>("qr", "none")),
103  thicknessLayers_(),
104  kappaLayers_()
105 {
106  switch (mode_)
107  {
108  case fixedPower:
109  {
110  dict.lookup("Q") >> Q_;
111 
112  break;
113  }
114  case fixedHeatFlux:
115  {
116  q_ = scalarField("q", dict, p.size());
117 
118  break;
119  }
120  case fixedHeatTransferCoeff:
121  {
122  h_ = scalarField("h", dict, p.size());
123  Ta_ = Function1<scalar>::New("Ta", dict);
124 
125  if (dict.found("thicknessLayers"))
126  {
127  dict.lookup("thicknessLayers") >> thicknessLayers_;
128  dict.lookup("kappaLayers") >> kappaLayers_;
129  }
130 
131  break;
132  }
133  }
134 
135  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
136 
137  if (qrName_ != "none")
138  {
139  if (dict.found("qrPrevious"))
140  {
141  qrPrevious_ = scalarField("qrPrevious", dict, p.size());
142  }
143  else
144  {
145  qrPrevious_.setSize(p.size(), 0);
146  }
147  }
148 
149  if (dict.found("refValue"))
150  {
151  // Full restart
152  refValue() = scalarField("refValue", dict, p.size());
153  refGrad() = scalarField("refGradient", dict, p.size());
154  valueFraction() = scalarField("valueFraction", dict, p.size());
155  }
156  else
157  {
158  // Start from user entered data. Assume fixedValue.
159  refValue() = *this;
160  refGrad() = 0;
161  valueFraction() = 1;
162  }
163 }
164 
165 
168 (
170  const fvPatch& p,
172  const fvPatchFieldMapper& mapper
173 )
174 :
175  mixedFvPatchScalarField(ptf, p, iF, mapper),
176  temperatureCoupledBase(patch(), ptf),
177  mode_(ptf.mode_),
178  Q_(ptf.Q_),
179  Ta_(ptf.Ta_, false),
180  relaxation_(ptf.relaxation_),
181  emissivity_(ptf.emissivity_),
182  qrRelaxation_(ptf.qrRelaxation_),
183  qrName_(ptf.qrName_),
184  thicknessLayers_(ptf.thicknessLayers_),
185  kappaLayers_(ptf.kappaLayers_)
186 {
187  switch (mode_)
188  {
189  case fixedPower:
190  {
191  break;
192  }
193  case fixedHeatFlux:
194  {
195  mapper(q_, ptf.q_);
196  break;
197  }
198  case fixedHeatTransferCoeff:
199  {
200  mapper(h_, ptf.h_);
201  break;
202  }
203  }
204 
205  if (qrName_ != "none")
206  {
207  mapper(qrPrevious_, ptf.qrPrevious_);
208  }
209 }
210 
211 
214 (
216 )
217 :
218  mixedFvPatchScalarField(tppsf),
219  temperatureCoupledBase(tppsf),
220  mode_(tppsf.mode_),
221  Q_(tppsf.Q_),
222  q_(tppsf.q_),
223  h_(tppsf.h_),
224  Ta_(tppsf.Ta_, false),
225  relaxation_(tppsf.relaxation_),
226  emissivity_(tppsf.emissivity_),
227  qrPrevious_(tppsf.qrPrevious_),
228  qrRelaxation_(tppsf.qrRelaxation_),
229  qrName_(tppsf.qrName_),
230  thicknessLayers_(tppsf.thicknessLayers_),
231  kappaLayers_(tppsf.kappaLayers_)
232 {}
233 
234 
237 (
240 )
241 :
242  mixedFvPatchScalarField(tppsf, iF),
243  temperatureCoupledBase(patch(), tppsf),
244  mode_(tppsf.mode_),
245  Q_(tppsf.Q_),
246  q_(tppsf.q_),
247  h_(tppsf.h_),
248  Ta_(tppsf.Ta_, false),
249  relaxation_(tppsf.relaxation_),
250  emissivity_(tppsf.emissivity_),
251  qrPrevious_(tppsf.qrPrevious_),
252  qrRelaxation_(tppsf.qrRelaxation_),
253  qrName_(tppsf.qrName_),
254  thicknessLayers_(tppsf.thicknessLayers_),
255  kappaLayers_(tppsf.kappaLayers_)
256 {}
257 
258 
259 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
260 
262 (
263  const fvPatchFieldMapper& m
264 )
265 {
266  mixedFvPatchScalarField::autoMap(m);
267 
268  switch (mode_)
269  {
270  case fixedPower:
271  {
272  break;
273  }
274  case fixedHeatFlux:
275  {
276  m(q_, q_);
277 
278  break;
279  }
280  case fixedHeatTransferCoeff:
281  {
282  m(h_, h_);
283 
284  break;
285  }
286  }
287 
288  if (qrName_ != "none")
289  {
290  m(qrPrevious_, qrPrevious_);
291  }
292 }
293 
294 
296 (
297  const fvPatchScalarField& ptf,
298  const labelList& addr
299 )
300 {
301  mixedFvPatchScalarField::rmap(ptf, addr);
302 
304  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
305 
306  switch (mode_)
307  {
308  case fixedPower:
309  {
310  break;
311  }
312  case fixedHeatFlux:
313  {
314  q_.rmap(tiptf.q_, addr);
315 
316  break;
317  }
318  case fixedHeatTransferCoeff:
319  {
320  h_.rmap(tiptf.h_, addr);
321 
322  break;
323  }
324  }
325 
326  if (qrName_ != "none")
327  {
328  qrPrevious_.rmap(tiptf.qrPrevious_, addr);
329  }
330 }
331 
332 
334 {
335  if (updated())
336  {
337  return;
338  }
339 
340  const scalarField& Tp(*this);
341 
342  // Store current valueFraction and refValue for relaxation
343  const scalarField valueFraction0(valueFraction());
344  const scalarField refValue0(refValue());
345 
346  scalarField qr(Tp.size(), 0);
347  if (qrName_ != "none")
348  {
349  qr =
350  qrRelaxation_
351  *patch().lookupPatchField<volScalarField, scalar>(qrName_)
352  + (1 - qrRelaxation_)*qrPrevious_;
353 
354  qrPrevious_ = qr;
355  }
356 
357  switch (mode_)
358  {
359  case fixedPower:
360  {
361  refGrad() = (Q_/gSum(patch().magSf()) + qr)/kappa(Tp);
362  refValue() = Tp;
363  valueFraction() = 0;
364 
365  break;
366  }
367  case fixedHeatFlux:
368  {
369  refGrad() = (q_ + qr)/kappa(Tp);
370  refValue() = Tp;
371  valueFraction() = 0;
372 
373  break;
374  }
375  case fixedHeatTransferCoeff:
376  {
377  scalar totalSolidRes = 0;
378  if (thicknessLayers_.size())
379  {
380  forAll(thicknessLayers_, iLayer)
381  {
382  const scalar l = thicknessLayers_[iLayer];
383  if (kappaLayers_[iLayer] > 0)
384  {
385  totalSolidRes += l/kappaLayers_[iLayer];
386  }
387  }
388  }
389  scalarField hp(1/(1/h_ + totalSolidRes));
390 
391  const scalar Ta = Ta_->value(this->db().time().timeOutputValue());
392  scalarField hpTa(hp*Ta);
393 
394  if (emissivity_ > 0)
395  {
396  // Evaluate the radiative flux to the environment
397  // from the surface temperature ...
398  if (totalSolidRes > 0)
399  {
400  // ... including the effect of the solid wall thermal
401  // resistance
402  scalarField TpLambda(h_/(h_ + 1/totalSolidRes));
403  scalarField Ts(TpLambda*Tp + (1 - TpLambda)*Ta);
404  scalarField lambdaTa4(pow4((1 - TpLambda)*Ta));
405 
406  hp += emissivity_*sigma.value()*(pow4(Ts) - lambdaTa4)/Tp;
407  hpTa += emissivity_*sigma.value()*(lambdaTa4 + pow4(Ta));
408  }
409  else
410  {
411  // ... if there is no solid wall thermal resistance use
412  // the current wall temperature
413  hp += emissivity_*sigma.value()*pow3(Tp);
414  hpTa += emissivity_*sigma.value()*pow4(Ta);
415  }
416  }
417 
418  const scalarField kappaDeltaCoeffs
419  (
420  this->kappa(Tp)*patch().deltaCoeffs()
421  );
422 
423  refGrad() = 0;
424 
425  forAll(Tp, i)
426  {
427  if (qr[i] < 0)
428  {
429  const scalar hpmqr = hp[i] - qr[i]/Tp[i];
430 
431  refValue()[i] = hpTa[i]/hpmqr;
432  valueFraction()[i] = hpmqr/(hpmqr + kappaDeltaCoeffs[i]);
433  }
434  else
435  {
436  refValue()[i] = (hpTa[i] + qr[i])/hp[i];
437  valueFraction()[i] = hp[i]/(hp[i] + kappaDeltaCoeffs[i]);
438  }
439  }
440 
441  break;
442  }
443  }
444 
445  valueFraction() =
446  relaxation_*valueFraction()
447  + (1 - relaxation_)*valueFraction0;
448 
449  refValue() = relaxation_*refValue() + (1 - relaxation_)*refValue0;
450 
451  mixedFvPatchScalarField::updateCoeffs();
452 
453  if (debug)
454  {
455  const scalar Q = gSum(kappa(Tp)*patch().magSf()*snGrad());
456 
457  Info<< patch().boundaryMesh().mesh().name() << ':'
458  << patch().name() << ':'
459  << this->internalField().name() << " :"
460  << " heat transfer rate:" << Q
461  << " walltemperature "
462  << " min:" << gMin(*this)
463  << " max:" << gMax(*this)
464  << " avg:" << gAverage(*this)
465  << endl;
466  }
467 }
468 
469 
471 (
472  Ostream& os
473 ) const
474 {
476 
477  writeEntry(os, "mode", operationModeNames[mode_]);
479 
480  switch (mode_)
481  {
482  case fixedPower:
483  {
484  writeEntry(os, "Q", Q_);
485 
486  break;
487  }
488  case fixedHeatFlux:
489  {
490  writeEntry(os, "q", q_);
491 
492  break;
493  }
494  case fixedHeatTransferCoeff:
495  {
496  writeEntry(os, "h", h_);
497  writeEntry(os, Ta_());
498 
499  if (relaxation_ < 1)
500  {
501  writeEntry(os, "relaxation", relaxation_);
502  }
503 
504  if (emissivity_ > 0)
505  {
506  writeEntry(os, "emissivity", emissivity_);
507  }
508 
509  if (thicknessLayers_.size())
510  {
511  writeEntry(os, "thicknessLayers", thicknessLayers_);
512  writeEntry(os, "kappaLayers", kappaLayers_);
513  }
514 
515  break;
516  }
517  }
518 
519  writeEntry(os, "qr", qrName_);
520 
521  if (qrName_ != "none")
522  {
523  writeEntry(os, "qrRelaxation", qrRelaxation_);
524 
525  writeEntry(os, "qrPrevious", qrPrevious_);
526  }
527 
528  writeEntry(os, "refValue", refValue());
529  writeEntry(os, "refGradient", refGrad());
530  writeEntry(os, "valueFraction", valueFraction());
531  writeEntry(os, "value", *this);
532 }
533 
534 
535 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
536 
537 namespace Foam
538 {
540  (
543  );
544 }
545 
546 // ************************************************************************* //
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
Type gMin(const FieldField< Field, Type > &f)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const Type & value() const
Return const reference to value.
virtual label size() const
Return size.
Definition: fvPatch.H:155
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Type gMax(const FieldField< Field, Type > &f)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar pow3(const dimensionedScalar &ds)
void setSize(const label)
Reset size of List.
Definition: List.C:281
Common functions used in temperature coupled boundaries.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Type gAverage(const FieldField< Field, Type > &f)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:295
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow4(const dimensionedScalar &ds)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
messageStream Info
externalWallHeatFluxTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
This boundary condition applies a heat flux condition to temperature on an external wall in one of th...
void write(Ostream &) const
Write.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
static autoPtr< Function1< Type > > New(const word &entryName, const dictionary &dict)
Selector.
Definition: Function1New.C:32
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583