externalWallHeatFluxTemperatureFvPatchScalarField.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 | Copyright (C) 2011-2017 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_(Function1<scalar>::New("Ta", dict)),
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 
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  q_.map(ptf.q_, mapper);
195 
196  break;
197  }
198  case fixedHeatTransferCoeff:
199  {
200  h_.map(ptf.h_, mapper);
201 
202  break;
203  }
204  }
205 
206  if (qrName_ != "none")
207  {
208  qrPrevious_.map(ptf.qrPrevious_, mapper);
209  }
210 }
211 
212 
215 (
217 )
218 :
219  mixedFvPatchScalarField(tppsf),
220  temperatureCoupledBase(tppsf),
221  mode_(tppsf.mode_),
222  Q_(tppsf.Q_),
223  q_(tppsf.q_),
224  h_(tppsf.h_),
225  Ta_(tppsf.Ta_, false),
226  relaxation_(tppsf.relaxation_),
227  emissivity_(tppsf.emissivity_),
228  qrPrevious_(tppsf.qrPrevious_),
229  qrRelaxation_(tppsf.qrRelaxation_),
230  qrName_(tppsf.qrName_),
231  thicknessLayers_(tppsf.thicknessLayers_),
232  kappaLayers_(tppsf.kappaLayers_)
233 {}
234 
235 
238 (
241 )
242 :
243  mixedFvPatchScalarField(tppsf, iF),
244  temperatureCoupledBase(patch(), tppsf),
245  mode_(tppsf.mode_),
246  Q_(tppsf.Q_),
247  q_(tppsf.q_),
248  h_(tppsf.h_),
249  Ta_(tppsf.Ta_, false),
250  relaxation_(tppsf.relaxation_),
251  emissivity_(tppsf.emissivity_),
252  qrPrevious_(tppsf.qrPrevious_),
253  qrRelaxation_(tppsf.qrRelaxation_),
254  qrName_(tppsf.qrName_),
255  thicknessLayers_(tppsf.thicknessLayers_),
256  kappaLayers_(tppsf.kappaLayers_)
257 {}
258 
259 
260 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
261 
263 (
264  const fvPatchFieldMapper& m
265 )
266 {
267  mixedFvPatchScalarField::autoMap(m);
268 
269  switch (mode_)
270  {
271  case fixedPower:
272  {
273  break;
274  }
275  case fixedHeatFlux:
276  {
277  q_.autoMap(m);
278 
279  break;
280  }
281  case fixedHeatTransferCoeff:
282  {
283  h_.autoMap(m);
284 
285  break;
286  }
287  }
288 
289  if (qrName_ != "none")
290  {
291  qrPrevious_.autoMap(m);
292  }
293 }
294 
295 
297 (
298  const fvPatchScalarField& ptf,
299  const labelList& addr
300 )
301 {
302  mixedFvPatchScalarField::rmap(ptf, addr);
303 
305  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
306 
307  switch (mode_)
308  {
309  case fixedPower:
310  {
311  break;
312  }
313  case fixedHeatFlux:
314  {
315  q_.rmap(tiptf.q_, addr);
316 
317  break;
318  }
319  case fixedHeatTransferCoeff:
320  {
321  h_.rmap(tiptf.h_, addr);
322 
323  break;
324  }
325  }
326 
327  if (qrName_ != "none")
328  {
329  qrPrevious_.rmap(tiptf.qrPrevious_, addr);
330  }
331 }
332 
333 
335 {
336  if (updated())
337  {
338  return;
339  }
340 
341  const scalarField& Tp(*this);
342 
343  scalarField qr(Tp.size(), 0);
344  if (qrName_ != "none")
345  {
346  qr =
347  qrRelaxation_
348  *patch().lookupPatchField<volScalarField, scalar>(qrName_)
349  + (1 - qrRelaxation_)*qrPrevious_;
350 
351  qrPrevious_ = qr;
352  }
353 
354  switch (mode_)
355  {
356  case fixedPower:
357  {
358  refGrad() = (Q_/gSum(patch().magSf()) + qr)/kappa(Tp);
359  refValue() = 0;
360  valueFraction() = 0;
361 
362  break;
363  }
364  case fixedHeatFlux:
365  {
366  refGrad() = (q_ + qr)/kappa(Tp);
367  refValue() = 0;
368  valueFraction() = 0;
369 
370  break;
371  }
372  case fixedHeatTransferCoeff:
373  {
374  scalar totalSolidRes = 0;
375  if (thicknessLayers_.size())
376  {
377  forAll(thicknessLayers_, iLayer)
378  {
379  const scalar l = thicknessLayers_[iLayer];
380  if (kappaLayers_[iLayer] > 0)
381  {
382  totalSolidRes += l/kappaLayers_[iLayer];
383  }
384  }
385  }
386  scalarField hp(1/(1/h_ + totalSolidRes));
387 
388  const scalar Ta = Ta_->value(this->db().time().timeOutputValue());
389  scalarField hpTa(hp*Ta);
390 
391  if (emissivity_ > 0)
392  {
393  // Evaluate the radiative flux to the environment
394  // from the surface temperature ...
395  if (totalSolidRes > 0)
396  {
397  // ... including the effect of the solid wall thermal
398  // resistance
399  scalarField TpLambda(h_/(h_ + 1/totalSolidRes));
400  scalarField Ts(TpLambda*Tp + (1 - TpLambda)*Ta);
401  scalarField lambdaTa4(pow4((1 - TpLambda)*Ta));
402 
403  hp += emissivity_*sigma.value()*(pow4(Ts) - lambdaTa4)/Tp;
404  hpTa += sigma.value()*(emissivity_*lambdaTa4 + pow4(Ta));
405  }
406  else
407  {
408  // ... if there is no solid wall thermal resistance use
409  // the current wall temperature
410  hp += emissivity_*sigma.value()*pow3(Tp);
411  hpTa += sigma.value()*pow4(Ta);
412  }
413  }
414 
415  const scalarField kappaDeltaCoeffs
416  (
417  this->kappa(Tp)*patch().deltaCoeffs()
418  );
419 
420  refGrad() = 0;
421 
422  forAll(Tp, i)
423  {
424  if (qr[i] < 0)
425  {
426  const scalar hpmqr = hp[i] - qr[i]/Tp[i];
427 
428  refValue()[i] = hpTa[i]/hpmqr;
429  valueFraction()[i] = hpmqr/(hpmqr + kappaDeltaCoeffs[i]);
430  }
431  else
432  {
433  refValue()[i] = (hpTa[i] + qr[i])/hp[i];
434  valueFraction()[i] = hp[i]/(hp[i] + kappaDeltaCoeffs[i]);
435  }
436  }
437 
438  break;
439  }
440  }
441 
442  valueFraction() = relaxation_*valueFraction() + (1 - relaxation_);
443  refValue() = relaxation_*refValue() + (1 - relaxation_)*Tp;
444 
445  mixedFvPatchScalarField::updateCoeffs();
446 
447  if (debug)
448  {
449  const scalar Q = gSum(kappa(Tp)*patch().magSf()*snGrad());
450 
451  Info<< patch().boundaryMesh().mesh().name() << ':'
452  << patch().name() << ':'
453  << this->internalField().name() << " :"
454  << " heat transfer rate:" << Q
455  << " walltemperature "
456  << " min:" << gMin(*this)
457  << " max:" << gMax(*this)
458  << " avg:" << gAverage(*this)
459  << endl;
460  }
461 }
462 
463 
465 (
466  Ostream& os
467 ) const
468 {
470 
471  os.writeKeyword("mode")
472  << operationModeNames[mode_] << token::END_STATEMENT << nl;
474 
475  switch (mode_)
476  {
477  case fixedPower:
478  {
479  os.writeKeyword("Q")
480  << Q_ << token::END_STATEMENT << nl;
481 
482  break;
483  }
484  case fixedHeatFlux:
485  {
486  q_.writeEntry("q", os);
487 
488  break;
489  }
490  case fixedHeatTransferCoeff:
491  {
492  h_.writeEntry("h", os);
493  Ta_->writeData(os);
494 
495  if (relaxation_ < 1)
496  {
497  os.writeKeyword("relaxation")
498  << relaxation_ << token::END_STATEMENT << nl;
499  }
500 
501  if (emissivity_ > 0)
502  {
503  os.writeKeyword("emissivity")
504  << emissivity_ << token::END_STATEMENT << nl;
505  }
506 
507  if (thicknessLayers_.size())
508  {
509  thicknessLayers_.writeEntry("thicknessLayers", os);
510  kappaLayers_.writeEntry("kappaLayers", os);
511  }
512 
513  break;
514  }
515  }
516 
517  os.writeKeyword("qr")<< qrName_ << token::END_STATEMENT << nl;
518 
519  if (qrName_ != "none")
520  {
521  os.writeKeyword("qrRelaxation")
522  << qrRelaxation_ << token::END_STATEMENT << nl;
523 
524  qrPrevious_.writeEntry("qrPrevious", os);
525  }
526 
527  refValue().writeEntry("refValue", os);
528  refGrad().writeEntry("refGradient", os);
529  valueFraction().writeEntry("valueFraction", os);
530  writeEntry("value", os);
531 }
532 
533 
534 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
535 
536 namespace Foam
537 {
539  (
542  );
543 }
544 
545 // ************************************************************************* //
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:53
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
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:253
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
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:362
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
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
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
Type gMax(const FieldField< Field, Type > &f)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
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:395
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.
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.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576