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-2015 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"
30 #include "mappedPatchBase.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  template<>
37  const char*
38  NamedEnum
39  <
41  3
42  >::names[] =
43  {
44  "fixed_heat_flux",
45  "fixed_heat_transfer_coefficient",
46  "unknown"
47  };
48 
49 } // End namespace Foam
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(), "undefined", "undefined", "undefined-K"),
69  mode_(unknown),
70  q_(p.size(), 0.0),
71  h_(p.size(), 0.0),
72  Ta_(p.size(), 0.0),
73  QrPrevious_(p.size()),
74  QrRelaxation_(1),
75  QrName_("undefined-Qr"),
76  thicknessLayers_(),
77  kappaLayers_()
78 {
79  refValue() = 0.0;
80  refGrad() = 0.0;
81  valueFraction() = 1.0;
82 }
83 
84 
87 (
89  const fvPatch& p,
91  const fvPatchFieldMapper& mapper
92 )
93 :
94  mixedFvPatchScalarField(ptf, p, iF, mapper),
95  temperatureCoupledBase(patch(), ptf),
96  mode_(ptf.mode_),
97  q_(ptf.q_, mapper),
98  h_(ptf.h_, mapper),
99  Ta_(ptf.Ta_, mapper),
100  QrPrevious_(ptf.QrPrevious_, mapper),
101  QrRelaxation_(ptf.QrRelaxation_),
102  QrName_(ptf.QrName_),
103  thicknessLayers_(ptf.thicknessLayers_),
104  kappaLayers_(ptf.kappaLayers_)
105 {}
106 
107 
110 (
111  const fvPatch& p,
113  const dictionary& dict
114 )
115 :
116  mixedFvPatchScalarField(p, iF),
117  temperatureCoupledBase(patch(), dict),
118  mode_(unknown),
119  q_(p.size(), 0.0),
120  h_(p.size(), 0.0),
121  Ta_(p.size(), 0.0),
122  QrPrevious_(p.size(), 0.0),
123  QrRelaxation_(dict.lookupOrDefault<scalar>("relaxation", 1)),
124  QrName_(dict.lookupOrDefault<word>("Qr", "none")),
125  thicknessLayers_(),
126  kappaLayers_()
127 {
128  if (dict.found("q") && !dict.found("h") && !dict.found("Ta"))
129  {
130  mode_ = fixedHeatFlux;
131  q_ = scalarField("q", dict, p.size());
132  }
133  else if (dict.found("h") && dict.found("Ta") && !dict.found("q"))
134  {
135  mode_ = fixedHeatTransferCoeff;
136  h_ = scalarField("h", dict, p.size());
137  Ta_ = scalarField("Ta", dict, p.size());
138  if (dict.found("thicknessLayers"))
139  {
140  dict.lookup("thicknessLayers") >> thicknessLayers_;
141  dict.lookup("kappaLayers") >> kappaLayers_;
142  }
143  }
144  else
145  {
147  (
148  "externalWallHeatFluxTemperatureFvPatchScalarField::"
149  "externalWallHeatFluxTemperatureFvPatchScalarField\n"
150  "(\n"
151  " const fvPatch& p,\n"
152  " const DimensionedField<scalar, volMesh>& iF,\n"
153  " const dictionary& dict\n"
154  ")\n"
155  ) << "\n patch type '" << p.type()
156  << "' either q or h and Ta were not found '"
157  << "\n for patch " << p.name()
158  << " of field " << dimensionedInternalField().name()
159  << " in file " << dimensionedInternalField().objectPath()
160  << exit(FatalError);
161  }
162 
163  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
164 
165  if (dict.found("QrPrevious"))
166  {
167  QrPrevious_ = scalarField("QrPrevious", dict, p.size());
168  }
169 
170  if (dict.found("refValue"))
171  {
172  // Full restart
173  refValue() = scalarField("refValue", dict, p.size());
174  refGrad() = scalarField("refGradient", dict, p.size());
175  valueFraction() = scalarField("valueFraction", dict, p.size());
176  }
177  else
178  {
179  // Start from user entered data. Assume fixedValue.
180  refValue() = *this;
181  refGrad() = 0.0;
182  valueFraction() = 1.0;
183  }
184 }
185 
186 
189 (
191 )
192 :
193  mixedFvPatchScalarField(tppsf),
194  temperatureCoupledBase(tppsf),
195  mode_(tppsf.mode_),
196  q_(tppsf.q_),
197  h_(tppsf.h_),
198  Ta_(tppsf.Ta_),
199  QrPrevious_(tppsf.QrPrevious_),
200  QrRelaxation_(tppsf.QrRelaxation_),
201  QrName_(tppsf.QrName_),
202  thicknessLayers_(tppsf.thicknessLayers_),
203  kappaLayers_(tppsf.kappaLayers_)
204 {}
205 
206 
209 (
212 )
213 :
214  mixedFvPatchScalarField(tppsf, iF),
215  temperatureCoupledBase(patch(), tppsf),
216  mode_(tppsf.mode_),
217  q_(tppsf.q_),
218  h_(tppsf.h_),
219  Ta_(tppsf.Ta_),
220  QrPrevious_(tppsf.QrPrevious_),
221  QrRelaxation_(tppsf.QrRelaxation_),
222  QrName_(tppsf.QrName_),
223  thicknessLayers_(tppsf.thicknessLayers_),
224  kappaLayers_(tppsf.kappaLayers_)
225 {}
226 
227 
228 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
229 
231 (
232  const fvPatchFieldMapper& m
233 )
234 {
235  mixedFvPatchScalarField::autoMap(m);
236  q_.autoMap(m);
237  h_.autoMap(m);
238  Ta_.autoMap(m);
239 }
240 
241 
243 (
244  const fvPatchScalarField& ptf,
245  const labelList& addr
246 )
247 {
248  mixedFvPatchScalarField::rmap(ptf, addr);
249 
251  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
252 
253  q_.rmap(tiptf.q_, addr);
254  h_.rmap(tiptf.h_, addr);
255  Ta_.rmap(tiptf.Ta_, addr);
256 }
257 
258 
260 {
261  if (updated())
262  {
263  return;
264  }
265 
266  const scalarField Tp(*this);
267  scalarField hp(patch().size(), 0.0);
268 
269  scalarField Qr(Tp.size(), 0.0);
270  if (QrName_ != "none")
271  {
272  Qr = patch().lookupPatchField<volScalarField, scalar>(QrName_);
273 
274  Qr = QrRelaxation_*Qr + (1.0 - QrRelaxation_)*QrPrevious_;
275  QrPrevious_ = Qr;
276  }
277 
278  switch (mode_)
279  {
280  case fixedHeatFlux:
281  {
282  refGrad() = (q_ + Qr)/kappa(Tp);
283  refValue() = 0.0;
284  valueFraction() = 0.0;
285 
286  break;
287  }
288  case fixedHeatTransferCoeff:
289  {
290  scalar totalSolidRes = 0.0;
291  if (thicknessLayers_.size() > 0)
292  {
293  forAll (thicknessLayers_, iLayer)
294  {
295  const scalar l = thicknessLayers_[iLayer];
296  if (kappaLayers_[iLayer] > 0.0)
297  {
298  totalSolidRes += l/kappaLayers_[iLayer];
299  }
300  }
301  }
302  hp = 1.0/(1.0/h_ + totalSolidRes);
303 
304  Qr /= Tp;
305  refGrad() = 0.0;
306  refValue() = hp*Ta_/(hp - Qr);
307  valueFraction() =
308  (hp - Qr)/((hp - Qr) + kappa(Tp)*patch().deltaCoeffs());
309 
310  break;
311  }
312  default:
313  {
315  (
316  "externalWallHeatFluxTemperatureFvPatchScalarField"
317  "::updateCoeffs()"
318  ) << "Illegal heat flux mode " << operationModeNames[mode_]
319  << exit(FatalError);
320  }
321  }
322 
323  mixedFvPatchScalarField::updateCoeffs();
324 
325  if (debug)
326  {
327  scalar Q = gSum(kappa(Tp)*patch().magSf()*snGrad());
328 
329  Info<< patch().boundaryMesh().mesh().name() << ':'
330  << patch().name() << ':'
331  << this->dimensionedInternalField().name() << " :"
332  << " heat transfer rate:" << Q
333  << " walltemperature "
334  << " min:" << gMin(*this)
335  << " max:" << gMax(*this)
336  << " avg:" << gAverage(*this)
337  << endl;
338  }
339 }
340 
341 
343 (
344  Ostream& os
345 ) const
346 {
349 
350  QrPrevious_.writeEntry("QrPrevious", os);
351  os.writeKeyword("Qr")<< QrName_ << token::END_STATEMENT << nl;
352  os.writeKeyword("relaxation")<< QrRelaxation_
353  << token::END_STATEMENT << nl;
354 
355  switch (mode_)
356  {
357 
358  case fixedHeatFlux:
359  {
360  q_.writeEntry("q", os);
361  break;
362  }
363  case fixedHeatTransferCoeff:
364  {
365  h_.writeEntry("h", os);
366  Ta_.writeEntry("Ta", os);
367  thicknessLayers_.writeEntry("thicknessLayers", os);
368  kappaLayers_.writeEntry("kappaLayers", os);
369  break;
370  }
371  default:
372  {
374  (
375  "void externalWallHeatFluxTemperatureFvPatchScalarField::write"
376  "("
377  "Ostream& os"
378  ") const"
379  ) << "Illegal heat flux mode " << operationModeNames[mode_]
380  << abort(FatalError);
381  }
382  }
383 }
384 
385 
386 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
387 
388 namespace Foam
389 {
391  (
394  );
395 }
396 
397 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:387
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A class for handling words, derived from string.
Definition: word.H:59
void write(Ostream &) const
Write.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::fvPatchFieldMapper.
messageStream Info
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
const word & name() const
Return name.
Definition: fvPatch.H:149
runTime write()
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Type gSum(const FieldField< Field, Type > &f)
externalWallHeatFluxTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
#define forAll(list, i)
Definition: UList.H:421
Macros for easy insertion into run-time selection tables.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Type gMin(const FieldField< Field, Type > &f)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
rDeltaT dimensionedInternalField()
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
Type gAverage(const FieldField< Field, Type > &f)
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Common functions for use in temperature coupled boundaries.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
virtual label size() const
Return size.
Definition: fvPatch.H:161
This function object calculates and outputs the second invariant of the velocity gradient tensor [1/s...
Definition: Q.H:66
This boundary condition supplies a heat flux condition for temperature on an external wall...
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Type gMax(const FieldField< Field, Type > &f)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.