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-2020 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 
27 #include "volFields.H"
30 
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  template<>
38  const char*
39  NamedEnum
40  <
42  3
43  >::names[] =
44  {
45  "power",
46  "flux",
47  "coefficient"
48  };
49 }
50 
51 const Foam::NamedEnum
52 <
54  3
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
62 (
63  const fvPatch& p,
65 )
66 :
67  mixedFvPatchScalarField(p, iF),
68  temperatureCoupledBase(patch()),
69  mode_(fixedHeatFlux),
70  Q_(0),
71  Ta_(),
72  relaxation_(1),
73  emissivity_(0),
74  qrRelaxation_(1),
75  qrName_("undefined-qr"),
76  thicknessLayers_(),
77  kappaLayers_()
78 {
79  refValue() = 0;
80  refGrad() = 0;
81  valueFraction() = 1;
82 }
83 
84 
87 (
88  const fvPatch& p,
90  const dictionary& dict
91 )
92 :
93  mixedFvPatchScalarField(p, iF),
94  temperatureCoupledBase(patch(), dict),
95  mode_(operationModeNames.read(dict.lookup("mode"))),
96  Q_(0),
97  Ta_(),
98  relaxation_(dict.lookupOrDefault<scalar>("relaxation", 1)),
99  emissivity_(dict.lookupOrDefault<scalar>("emissivity", 0)),
100  qrRelaxation_(dict.lookupOrDefault<scalar>("qrRelaxation", 1)),
101  qrName_(dict.lookupOrDefault<word>("qr", "none")),
102  thicknessLayers_(),
103  kappaLayers_()
104 {
105  switch (mode_)
106  {
107  case fixedPower:
108  {
109  dict.lookup("Q") >> Q_;
110 
111  break;
112  }
113  case fixedHeatFlux:
114  {
115  q_ = scalarField("q", dict, p.size());
116 
117  break;
118  }
119  case fixedHeatTransferCoeff:
120  {
121  h_ = scalarField("h", dict, p.size());
122  Ta_ = Function1<scalar>::New("Ta", dict);
123 
124  if (dict.found("thicknessLayers"))
125  {
126  dict.lookup("thicknessLayers") >> thicknessLayers_;
127  dict.lookup("kappaLayers") >> kappaLayers_;
128  }
129 
130  break;
131  }
132  }
133 
134  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
135 
136  if (qrName_ != "none")
137  {
138  if (dict.found("qrPrevious"))
139  {
140  qrPrevious_ = scalarField("qrPrevious", dict, p.size());
141  }
142  else
143  {
144  qrPrevious_.setSize(p.size(), 0);
145  }
146  }
147 
148  if (dict.found("refValue"))
149  {
150  // Full restart
151  refValue() = scalarField("refValue", dict, p.size());
152  refGrad() = scalarField("refGradient", dict, p.size());
153  valueFraction() = scalarField("valueFraction", dict, p.size());
154  }
155  else
156  {
157  // Start from user entered data. Assume fixedValue.
158  refValue() = *this;
159  refGrad() = 0;
160  valueFraction() = 1;
161  }
162 }
163 
164 
167 (
169  const fvPatch& p,
171  const fvPatchFieldMapper& mapper
172 )
173 :
174  mixedFvPatchScalarField(ptf, p, iF, mapper),
175  temperatureCoupledBase(patch(), ptf),
176  mode_(ptf.mode_),
177  Q_(ptf.Q_),
178  Ta_(ptf.Ta_, false),
179  relaxation_(ptf.relaxation_),
180  emissivity_(ptf.emissivity_),
181  qrRelaxation_(ptf.qrRelaxation_),
182  qrName_(ptf.qrName_),
183  thicknessLayers_(ptf.thicknessLayers_),
184  kappaLayers_(ptf.kappaLayers_)
185 {
186  switch (mode_)
187  {
188  case fixedPower:
189  {
190  break;
191  }
192  case fixedHeatFlux:
193  {
194  mapper(q_, ptf.q_);
195  break;
196  }
197  case fixedHeatTransferCoeff:
198  {
199  mapper(h_, ptf.h_);
200  break;
201  }
202  }
203 
204  if (qrName_ != "none")
205  {
206  mapper(qrPrevious_, ptf.qrPrevious_);
207  }
208 }
209 
210 
213 (
216 )
217 :
218  mixedFvPatchScalarField(tppsf, iF),
219  temperatureCoupledBase(patch(), 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 
235 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
236 
238 (
239  const fvPatchFieldMapper& m
240 )
241 {
242  mixedFvPatchScalarField::autoMap(m);
243 
244  switch (mode_)
245  {
246  case fixedPower:
247  {
248  break;
249  }
250  case fixedHeatFlux:
251  {
252  m(q_, q_);
253 
254  break;
255  }
256  case fixedHeatTransferCoeff:
257  {
258  m(h_, h_);
259 
260  break;
261  }
262  }
263 
264  if (qrName_ != "none")
265  {
266  m(qrPrevious_, qrPrevious_);
267  }
268 }
269 
270 
272 (
273  const fvPatchScalarField& ptf,
274  const labelList& addr
275 )
276 {
277  mixedFvPatchScalarField::rmap(ptf, addr);
278 
280  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
281 
282  switch (mode_)
283  {
284  case fixedPower:
285  {
286  break;
287  }
288  case fixedHeatFlux:
289  {
290  q_.rmap(tiptf.q_, addr);
291 
292  break;
293  }
294  case fixedHeatTransferCoeff:
295  {
296  h_.rmap(tiptf.h_, addr);
297 
298  break;
299  }
300  }
301 
302  if (qrName_ != "none")
303  {
304  qrPrevious_.rmap(tiptf.qrPrevious_, addr);
305  }
306 }
307 
308 
310 {
311  if (updated())
312  {
313  return;
314  }
315 
316  const scalarField& Tp(*this);
317 
318  // Store current valueFraction and refValue for relaxation
319  const scalarField valueFraction0(valueFraction());
320  const scalarField refValue0(refValue());
321 
322  scalarField qr(Tp.size(), 0);
323  if (qrName_ != "none")
324  {
325  qr =
326  qrRelaxation_
327  *patch().lookupPatchField<volScalarField, scalar>(qrName_)
328  + (1 - qrRelaxation_)*qrPrevious_;
329 
330  qrPrevious_ = qr;
331  }
332 
333  switch (mode_)
334  {
335  case fixedPower:
336  {
337  refGrad() = (Q_/gSum(patch().magSf()) + qr)/kappa(*this);
338  refValue() = Tp;
339  valueFraction() = 0;
340 
341  break;
342  }
343  case fixedHeatFlux:
344  {
345  refGrad() = (q_ + qr)/kappa(*this);
346  refValue() = Tp;
347  valueFraction() = 0;
348 
349  break;
350  }
351  case fixedHeatTransferCoeff:
352  {
353  scalar totalSolidRes = 0;
354  if (thicknessLayers_.size())
355  {
356  forAll(thicknessLayers_, iLayer)
357  {
358  const scalar l = thicknessLayers_[iLayer];
359  if (kappaLayers_[iLayer] > 0)
360  {
361  totalSolidRes += l/kappaLayers_[iLayer];
362  }
363  }
364  }
365 
366  const scalar Ta = Ta_->value(this->db().time().timeOutputValue());
367 
368  const scalarField hp
369  (
370  1
371  /(
372  1
373  /(
374  (emissivity_ > 0)
375  ? (
376  h_
377  + emissivity_*sigma.value()
378  *((pow3(Ta) + pow3(Tp)) + Ta*Tp*(Ta + Tp))
379  )()
380  : h_
381  ) + totalSolidRes
382  )
383  );
384 
385  const scalarField hpTa(hp*Ta);
386 
387  const scalarField kappaDeltaCoeffs
388  (
389  this->kappa(*this)*patch().deltaCoeffs()
390  );
391 
392  refGrad() = 0;
393 
394  forAll(Tp, i)
395  {
396  if (qr[i] < 0)
397  {
398  const scalar hpmqr = hp[i] - qr[i]/Tp[i];
399 
400  refValue()[i] = hpTa[i]/hpmqr;
401  valueFraction()[i] = hpmqr/(hpmqr + kappaDeltaCoeffs[i]);
402  }
403  else
404  {
405  refValue()[i] = (hpTa[i] + qr[i])/hp[i];
406  valueFraction()[i] = hp[i]/(hp[i] + kappaDeltaCoeffs[i]);
407  }
408  }
409 
410  break;
411  }
412  }
413 
414  valueFraction() =
415  relaxation_*valueFraction()
416  + (1 - relaxation_)*valueFraction0;
417 
418  refValue() = relaxation_*refValue() + (1 - relaxation_)*refValue0;
419 
420  mixedFvPatchScalarField::updateCoeffs();
421 
422  if (debug)
423  {
424  const scalar Q = gSum(kappa(*this)*patch().magSf()*snGrad());
425 
426  Info<< patch().boundaryMesh().mesh().name() << ':'
427  << patch().name() << ':'
428  << this->internalField().name() << " :"
429  << " heat transfer rate:" << Q
430  << " walltemperature "
431  << " min:" << gMin(*this)
432  << " max:" << gMax(*this)
433  << " avg:" << gAverage(*this)
434  << endl;
435  }
436 }
437 
438 
440 (
441  Ostream& os
442 ) const
443 {
445 
446  writeEntry(os, "mode", operationModeNames[mode_]);
448 
449  switch (mode_)
450  {
451  case fixedPower:
452  {
453  writeEntry(os, "Q", Q_);
454 
455  break;
456  }
457  case fixedHeatFlux:
458  {
459  writeEntry(os, "q", q_);
460 
461  break;
462  }
463  case fixedHeatTransferCoeff:
464  {
465  writeEntry(os, "h", h_);
466  writeEntry(os, Ta_());
467 
468  if (relaxation_ < 1)
469  {
470  writeEntry(os, "relaxation", relaxation_);
471  }
472 
473  if (emissivity_ > 0)
474  {
475  writeEntry(os, "emissivity", emissivity_);
476  }
477 
478  if (thicknessLayers_.size())
479  {
480  writeEntry(os, "thicknessLayers", thicknessLayers_);
481  writeEntry(os, "kappaLayers", kappaLayers_);
482  }
483 
484  break;
485  }
486  }
487 
488  writeEntry(os, "qr", qrName_);
489 
490  if (qrName_ != "none")
491  {
492  writeEntry(os, "qrRelaxation", qrRelaxation_);
493 
494  writeEntry(os, "qrPrevious", qrPrevious_);
495  }
496 
497  writeEntry(os, "refValue", refValue());
498  writeEntry(os, "refGradient", refGrad());
499  writeEntry(os, "valueFraction", valueFraction());
500  writeEntry(os, "value", *this);
501 }
502 
503 
504 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
505 
506 namespace Foam
507 {
509  (
512  );
513 }
514 
515 // ************************************************************************* //
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:643
#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:156
Type gMin(const FieldField< Field, Type > &f)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
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:260
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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:156
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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:275
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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
Namespace for OpenFOAM.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844