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 | Copyright (C) 2015 OpenCFD Ltd.
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  if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
141  {
142  dict.lookup("kappaLayers") >> kappaLayers_;
143 
144  if (thicknessLayers_.size() != kappaLayers_.size())
145  {
147  << "\n number of layers for thicknessLayers and "
148  << "kappaLayers must be the same"
149  << "\n for patch " << p.name()
150  << " of field " << dimensionedInternalField().name()
151  << " in file " << dimensionedInternalField().objectPath()
152  << exit(FatalIOError);
153  }
154  }
155  }
156  else
157  {
159  << "\n patch type '" << p.type()
160  << "' either q or h and Ta were not found '"
161  << "\n for patch " << p.name()
162  << " of field " << dimensionedInternalField().name()
163  << " in file " << dimensionedInternalField().objectPath()
164  << exit(FatalError);
165  }
166 
167  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
168 
169  if (dict.found("QrPrevious"))
170  {
171  QrPrevious_ = scalarField("QrPrevious", dict, p.size());
172  }
173 
174  if (dict.found("refValue"))
175  {
176  // Full restart
177  refValue() = scalarField("refValue", dict, p.size());
178  refGrad() = scalarField("refGradient", dict, p.size());
179  valueFraction() = scalarField("valueFraction", dict, p.size());
180  }
181  else
182  {
183  // Start from user entered data. Assume fixedValue.
184  refValue() = *this;
185  refGrad() = 0.0;
186  valueFraction() = 1.0;
187  }
188 }
189 
190 
193 (
195 )
196 :
197  mixedFvPatchScalarField(tppsf),
198  temperatureCoupledBase(tppsf),
199  mode_(tppsf.mode_),
200  q_(tppsf.q_),
201  h_(tppsf.h_),
202  Ta_(tppsf.Ta_),
203  QrPrevious_(tppsf.QrPrevious_),
204  QrRelaxation_(tppsf.QrRelaxation_),
205  QrName_(tppsf.QrName_),
206  thicknessLayers_(tppsf.thicknessLayers_),
207  kappaLayers_(tppsf.kappaLayers_)
208 {}
209 
210 
213 (
216 )
217 :
218  mixedFvPatchScalarField(tppsf, iF),
219  temperatureCoupledBase(patch(), tppsf),
220  mode_(tppsf.mode_),
221  q_(tppsf.q_),
222  h_(tppsf.h_),
223  Ta_(tppsf.Ta_),
224  QrPrevious_(tppsf.QrPrevious_),
225  QrRelaxation_(tppsf.QrRelaxation_),
226  QrName_(tppsf.QrName_),
227  thicknessLayers_(tppsf.thicknessLayers_),
228  kappaLayers_(tppsf.kappaLayers_)
229 {}
230 
231 
232 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
233 
235 (
236  const fvPatchFieldMapper& m
237 )
238 {
239  mixedFvPatchScalarField::autoMap(m);
240  q_.autoMap(m);
241  h_.autoMap(m);
242  Ta_.autoMap(m);
243 }
244 
245 
247 (
248  const fvPatchScalarField& ptf,
249  const labelList& addr
250 )
251 {
252  mixedFvPatchScalarField::rmap(ptf, addr);
253 
255  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
256 
257  q_.rmap(tiptf.q_, addr);
258  h_.rmap(tiptf.h_, addr);
259  Ta_.rmap(tiptf.Ta_, addr);
260 }
261 
262 
264 {
265  if (updated())
266  {
267  return;
268  }
269 
270  const scalarField Tp(*this);
271  scalarField hp(patch().size(), 0.0);
272 
273  scalarField Qr(Tp.size(), 0.0);
274  if (QrName_ != "none")
275  {
276  Qr = patch().lookupPatchField<volScalarField, scalar>(QrName_);
277 
279  QrPrevious_ = Qr;
280  }
281 
282  switch (mode_)
283  {
284  case fixedHeatFlux:
285  {
286  refGrad() = (q_ + Qr)/kappa(Tp);
287  refValue() = 0.0;
288  valueFraction() = 0.0;
289 
290  break;
291  }
293  {
294  scalar totalSolidRes = 0.0;
295  if (thicknessLayers_.size() > 0)
296  {
297  forAll (thicknessLayers_, iLayer)
298  {
299  const scalar l = thicknessLayers_[iLayer];
300  if (kappaLayers_[iLayer] > 0.0)
301  {
302  totalSolidRes += l/kappaLayers_[iLayer];
303  }
304  }
305  }
306  hp = 1.0/(1.0/h_ + totalSolidRes);
307 
308  Qr /= Tp;
309  refGrad() = 0.0;
310  refValue() = hp*Ta_/(hp - Qr);
311  valueFraction() =
312  (hp - Qr)/((hp - Qr) + kappa(Tp)*patch().deltaCoeffs());
313 
314  break;
315  }
316  default:
317  {
319  << "Illegal heat flux mode " << operationModeNames[mode_]
320  << exit(FatalError);
321  }
322  }
323 
324  mixedFvPatchScalarField::updateCoeffs();
325 
326  if (debug)
327  {
328  scalar Q = gSum(kappa(Tp)*patch().magSf()*snGrad());
329 
330  Info<< patch().boundaryMesh().mesh().name() << ':'
331  << patch().name() << ':'
332  << this->dimensionedInternalField().name() << " :"
333  << " heat transfer rate:" << Q
334  << " wall temperature "
335  << " min:" << gMin(*this)
336  << " max:" << gMax(*this)
337  << " avg:" << gAverage(*this)
338  << endl;
339  }
340 }
341 
342 
344 (
345  Ostream& os
346 ) const
347 {
350 
351  QrPrevious_.writeEntry("QrPrevious", os);
352  os.writeKeyword("Qr")<< QrName_ << token::END_STATEMENT << nl;
353  os.writeKeyword("relaxation")<< QrRelaxation_
354  << token::END_STATEMENT << nl;
355 
356  switch (mode_)
357  {
358 
359  case fixedHeatFlux:
360  {
361  q_.writeEntry("q", os);
362  break;
363  }
364  case fixedHeatTransferCoeff:
365  {
366  h_.writeEntry("h", os);
367  Ta_.writeEntry("Ta", os);
368  thicknessLayers_.writeEntry("thicknessLayers", os);
369  kappaLayers_.writeEntry("kappaLayers", os);
370  break;
371  }
372  default:
373  {
375  << "Illegal heat flux mode " << operationModeNames[mode_]
376  << abort(FatalError);
377  }
378  }
379 }
380 
381 
382 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
383 
384 namespace Foam
385 {
387  (
390  );
391 }
392 
393 // ************************************************************************* //
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::kappaLayers_
scalarList kappaLayers_
Conductivity of layers.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:219
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::q_
scalarField q_
Heat flux / [W/m2].
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:198
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::mode_
operationMode mode_
Operation mode.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:195
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::fixedHeatFlux
@ fixedHeatFlux
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:182
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:235
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::h_
scalarField h_
Heat transfer coefficient / [W/m2K].
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:201
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:247
Foam::temperatureCoupledBase
Common functions for use in temperature coupled boundaries.
Definition: temperatureCoupledBase.H:104
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::Ta_
scalarField Ta_
Ambient temperature / [K].
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:204
Foam::Q
This function object calculates and outputs the second invariant of the velocity gradient tensor [1/s...
Definition: Q.H:119
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
externalWallHeatFluxTemperatureFvPatchScalarField.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Qr
volScalarField Qr(IOobject("Qr", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("Qr", dimMass/pow3(dimTime), 0.0))
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::fixedHeatTransferCoeff
@ fixedHeatTransferCoeff
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:183
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::operationMode
operationMode
Operation mode enumeration.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:180
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::thicknessLayers_
scalarList thicknessLayers_
Thickness of layers.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:216
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::externalWallHeatFluxTemperatureFvPatchScalarField
externalWallHeatFluxTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:62
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::write
void write(Ostream &) const
Write.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:344
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::QrPrevious_
scalarField QrPrevious_
Chache Qr for relaxation.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:207
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::temperatureCoupledBase::kappa
tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Definition: temperatureCoupledBase.C:101
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::externalWallHeatFluxTemperatureFvPatchScalarField
This boundary condition supplies a heat flux condition for temperature on an external wall....
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:170
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:263
Foam::temperatureCoupledBase::write
void write(Ostream &) const
Write.
Definition: temperatureCoupledBase.C:225
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::operationModeNames
static const NamedEnum< operationMode, 3 > operationModeNames
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:187
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::QrRelaxation_
scalar QrRelaxation_
Relaxation for Qr.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:210
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::QrName_
const word QrName_
Name of the radiative heat flux.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:213
write
Tcoeff write()
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
mappedPatchBase.H