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-2016 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(), 0.0),
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  << "\n patch type '" << p.type()
148  << "' either q or h and Ta were not found '"
149  << "\n for patch " << p.name()
150  << " of field " << internalField().name()
151  << " in file " << internalField().objectPath()
152  << exit(FatalError);
153  }
154 
155  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
156 
157  if (dict.found("QrPrevious"))
158  {
159  QrPrevious_ = scalarField("QrPrevious", dict, p.size());
160  }
161 
162  if (dict.found("refValue"))
163  {
164  // Full restart
165  refValue() = scalarField("refValue", dict, p.size());
166  refGrad() = scalarField("refGradient", dict, p.size());
167  valueFraction() = scalarField("valueFraction", dict, p.size());
168  }
169  else
170  {
171  // Start from user entered data. Assume fixedValue.
172  refValue() = *this;
173  refGrad() = 0.0;
174  valueFraction() = 1.0;
175  }
176 }
177 
178 
181 (
183 )
184 :
185  mixedFvPatchScalarField(tppsf),
186  temperatureCoupledBase(tppsf),
187  mode_(tppsf.mode_),
188  q_(tppsf.q_),
189  h_(tppsf.h_),
190  Ta_(tppsf.Ta_),
191  QrPrevious_(tppsf.QrPrevious_),
192  QrRelaxation_(tppsf.QrRelaxation_),
193  QrName_(tppsf.QrName_),
194  thicknessLayers_(tppsf.thicknessLayers_),
195  kappaLayers_(tppsf.kappaLayers_)
196 {}
197 
198 
201 (
204 )
205 :
206  mixedFvPatchScalarField(tppsf, iF),
207  temperatureCoupledBase(patch(), tppsf),
208  mode_(tppsf.mode_),
209  q_(tppsf.q_),
210  h_(tppsf.h_),
211  Ta_(tppsf.Ta_),
212  QrPrevious_(tppsf.QrPrevious_),
213  QrRelaxation_(tppsf.QrRelaxation_),
214  QrName_(tppsf.QrName_),
215  thicknessLayers_(tppsf.thicknessLayers_),
216  kappaLayers_(tppsf.kappaLayers_)
217 {}
218 
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
223 (
224  const fvPatchFieldMapper& m
225 )
226 {
227  mixedFvPatchScalarField::autoMap(m);
228  q_.autoMap(m);
229  h_.autoMap(m);
230  Ta_.autoMap(m);
231  QrPrevious_.autoMap(m);
232 }
233 
234 
236 (
237  const fvPatchScalarField& ptf,
238  const labelList& addr
239 )
240 {
241  mixedFvPatchScalarField::rmap(ptf, addr);
242 
244  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
245 
246  q_.rmap(tiptf.q_, addr);
247  h_.rmap(tiptf.h_, addr);
248  Ta_.rmap(tiptf.Ta_, addr);
249  QrPrevious_.rmap(tiptf.QrPrevious_, addr);
250 }
251 
252 
254 {
255  if (updated())
256  {
257  return;
258  }
259 
260  const scalarField Tp(*this);
261  scalarField hp(patch().size(), 0.0);
262 
263  scalarField Qr(Tp.size(), 0.0);
264  if (QrName_ != "none")
265  {
266  Qr = patch().lookupPatchField<volScalarField, scalar>(QrName_);
267 
268  Qr = QrRelaxation_*Qr + (1.0 - QrRelaxation_)*QrPrevious_;
269  QrPrevious_ = Qr;
270  }
271 
272  switch (mode_)
273  {
274  case fixedHeatFlux:
275  {
276  refGrad() = (q_ + Qr)/kappa(Tp);
277  refValue() = 0.0;
278  valueFraction() = 0.0;
279 
280  break;
281  }
282  case fixedHeatTransferCoeff:
283  {
284  scalar totalSolidRes = 0.0;
285  if (thicknessLayers_.size() > 0)
286  {
287  forAll(thicknessLayers_, iLayer)
288  {
289  const scalar l = thicknessLayers_[iLayer];
290  if (kappaLayers_[iLayer] > 0.0)
291  {
292  totalSolidRes += l/kappaLayers_[iLayer];
293  }
294  }
295  }
296  hp = 1.0/(1.0/h_ + totalSolidRes);
297 
298  Qr /= Tp;
299  refGrad() = 0.0;
300  refValue() = hp*Ta_/(hp - Qr);
301  valueFraction() =
302  (hp - Qr)/((hp - Qr) + kappa(Tp)*patch().deltaCoeffs());
303 
304  break;
305  }
306  default:
307  {
309  << "Illegal heat flux mode " << operationModeNames[mode_]
310  << exit(FatalError);
311  }
312  }
313 
314  mixedFvPatchScalarField::updateCoeffs();
315 
316  if (debug)
317  {
318  scalar Q = gSum(kappa(Tp)*patch().magSf()*snGrad());
319 
320  Info<< patch().boundaryMesh().mesh().name() << ':'
321  << patch().name() << ':'
322  << this->internalField().name() << " :"
323  << " heat transfer rate:" << Q
324  << " walltemperature "
325  << " min:" << gMin(*this)
326  << " max:" << gMax(*this)
327  << " avg:" << gAverage(*this)
328  << endl;
329  }
330 }
331 
332 
334 (
335  Ostream& os
336 ) const
337 {
340 
341  QrPrevious_.writeEntry("QrPrevious", os);
342  os.writeKeyword("Qr")<< QrName_ << token::END_STATEMENT << nl;
343  os.writeKeyword("relaxation")<< QrRelaxation_
344  << token::END_STATEMENT << nl;
345 
346  switch (mode_)
347  {
348 
349  case fixedHeatFlux:
350  {
351  q_.writeEntry("q", os);
352  break;
353  }
354  case fixedHeatTransferCoeff:
355  {
356  h_.writeEntry("h", os);
357  Ta_.writeEntry("Ta", os);
358  thicknessLayers_.writeEntry("thicknessLayers", os);
359  kappaLayers_.writeEntry("kappaLayers", os);
360  break;
361  }
362  default:
363  {
365  << "Illegal heat flux mode " << operationModeNames[mode_]
366  << abort(FatalError);
367  }
368  }
369 }
370 
371 
372 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
373 
374 namespace Foam
375 {
377  (
380  );
381 }
382 
383 // ************************************************************************* //
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void write(Ostream &) const
Write.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Type gMin(const FieldField< Field, Type > &f)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const word & name() const
Return name.
Definition: fvPatch.H:149
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
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
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.
virtual label size() const
Return size.
Definition: fvPatch.H:161
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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
Common functions used in temperature coupled boundaries.
Type gAverage(const FieldField< Field, Type > &f)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:396
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.
volScalarField Qr(IOobject("Qr", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("Qr", dimMass/pow3(dimTime), 0.0))
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
This boundary condition supplies a heat flux condition for temperature on an external wall...
runTime 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:451