filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "mappedPatchBase.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
41 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
42 filmModel() const
43 {
44  HashTable<const filmModelType*> models
45  = db().time().lookupClass<filmModelType>();
46 
47  forAllConstIters(models, iter)
48  {
49  if (iter()->regionMesh().name() == filmRegionName_)
50  {
51  return *iter();
52  }
53  }
54 
55  DynamicList<word> modelNames;
56  forAllConstIters(models, iter)
57  {
58  modelNames.append(iter()->regionMesh().name());
59  }
60 
62  << "Unable to locate film region " << filmRegionName_
63  << ". Available regions include: " << modelNames
64  << abort(FatalError);
65 
66  return **models.begin();
67 }
68 
69 
72 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
73 pyrModel() const
74 {
75  HashTable<const pyrolysisModelType*> models =
77 
78  forAllConstIters(models, iter)
79  {
80  if (iter()->regionMesh().name() == pyrolysisRegionName_)
81  {
82  return *iter();
83  }
84  }
85 
86  DynamicList<word> modelNames;
87  forAllConstIters(models, iter)
88  {
89  modelNames.append(iter()->regionMesh().name());
90  }
91 
93  << "Unable to locate pyrolysis region " << pyrolysisRegionName_
94  << ". Available regions include: " << modelNames
95  << abort(FatalError);
96 
97  return **models.begin();
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
102 
105 (
106  const fvPatch& p,
108 )
109 :
110  mixedFvPatchScalarField(p, iF),
111  temperatureCoupledBase(patch()), // default method (fluidThermo)
112  filmRegionName_("surfaceFilmProperties"),
113  pyrolysisRegionName_("pyrolysisProperties"),
114  TnbrName_("undefined-Tnbr"),
115  qrName_("undefined-qr"),
116  convectiveScaling_(1.0),
117  filmDeltaDry_(0.0),
118  filmDeltaWet_(0.0)
119 {
120  this->refValue() = 0.0;
121  this->refGrad() = 0.0;
122  this->valueFraction() = 1.0;
123 }
124 
125 
128 (
130  const fvPatch& p,
132  const fvPatchFieldMapper& mapper
133 )
134 :
135  mixedFvPatchScalarField(psf, p, iF, mapper),
137  filmRegionName_(psf.filmRegionName_),
138  pyrolysisRegionName_(psf.pyrolysisRegionName_),
139  TnbrName_(psf.TnbrName_),
140  qrName_(psf.qrName_),
141  convectiveScaling_(psf.convectiveScaling_),
142  filmDeltaDry_(psf.filmDeltaDry_),
143  filmDeltaWet_(psf.filmDeltaWet_)
144 {}
145 
146 
149 (
150  const fvPatch& p,
152  const dictionary& dict
153 )
154 :
155  mixedFvPatchScalarField(p, iF),
157  filmRegionName_
158  (
159  dict.getOrDefault<word>("filmRegion", "surfaceFilmProperties")
160  ),
161  pyrolysisRegionName_
162  (
163  dict.getOrDefault<word>("pyrolysisRegion", "pyrolysisProperties")
164  ),
165  TnbrName_(dict.lookup("Tnbr")),
166  qrName_(dict.lookup("qr")),
167  convectiveScaling_(dict.getOrDefault<scalar>("convectiveScaling", 1)),
168  filmDeltaDry_(dict.get<scalar>("filmDeltaDry")),
169  filmDeltaWet_(dict.get<scalar>("filmDeltaWet"))
170 {
171  if (!isA<mappedPatchBase>(this->patch().patch()))
172  {
174  << "' not type '" << mappedPatchBase::typeName << "'"
175  << "\n for patch " << p.name()
176  << " of field " << internalField().name()
177  << " in file " << internalField().objectPath()
178  << exit(FatalError);
179  }
180 
181  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
182 
183  if (dict.found("refValue"))
184  {
185  // Full restart
186  refValue() = scalarField("refValue", dict, p.size());
187  refGrad() = scalarField("refGradient", dict, p.size());
188  valueFraction() = scalarField("valueFraction", dict, p.size());
189  }
190  else
191  {
192  // Start from user entered data. Assume fixedValue.
193  refValue() = *this;
194  refGrad() = 0.0;
195  valueFraction() = 1.0;
196  }
197 }
198 
199 
202 (
203  const filmPyrolysisRadiativeCoupledMixedFvPatchScalarField& psf,
204  const DimensionedField<scalar, volMesh>& iF
205 )
206 :
207  mixedFvPatchScalarField(psf, iF),
208  temperatureCoupledBase(patch(), psf),
209  filmRegionName_(psf.filmRegionName_),
210  pyrolysisRegionName_(psf.pyrolysisRegionName_),
211  TnbrName_(psf.TnbrName_),
212  qrName_(psf.qrName_),
213  convectiveScaling_(psf.convectiveScaling_),
214  filmDeltaDry_(psf.filmDeltaDry_),
215  filmDeltaWet_(psf.filmDeltaWet_)
216 {}
217 
218 
219 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
220 
222 (
223  const fvPatchFieldMapper& mapper
224 )
225 {
226  mixedFvPatchScalarField::autoMap(mapper);
228 }
229 
230 
232 (
233  const fvPatchField<scalar>& ptf,
234  const labelList& addr
235 )
236 {
237  mixedFvPatchScalarField::rmap(ptf, addr);
238 
240  refCast
241  <
243  >(ptf);
244 
245  temperatureCoupledBase::rmap(tiptf, addr);
246 }
247 
248 
250 {
251  if (updated())
252  {
253  return;
254  }
255 
256  // Get the coupling information from the mappedPatchBase
257  const mappedPatchBase& mpp =
258  refCast<const mappedPatchBase>(patch().patch());
259 
260  const label patchi = patch().index();
261  const label nbrPatchi = mpp.samplePolyPatch().index();
262  const polyMesh& mesh = patch().boundaryMesh().mesh();
263  const polyMesh& nbrMesh = mpp.sampleMesh();
264  const fvPatch& nbrPatch =
265  refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchi];
266 
267  scalarField intFld(patchInternalField());
268 
270  nbrField =
271  refCast
272  <
274  >
275  (
276  nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
277  );
278 
279  // Swap to obtain full local values of neighbour internal field
280  scalarField nbrIntFld(nbrField.patchInternalField());
281  mpp.distribute(nbrIntFld);
282 
283  scalarField& Tp = *this;
284 
285  const scalarField K(this->kappa(*this));
286  const scalarField nbrK(nbrField.kappa(*this));
287 
288  // Swap to obtain full local values of neighbour K*delta
289  scalarField KDeltaNbr(nbrK*nbrPatch.deltaCoeffs());
290  mpp.distribute(KDeltaNbr);
291 
292  scalarField myKDelta(K*patch().deltaCoeffs());
293 
294  scalarList Tfilm(patch().size(), Zero);
295  scalarList htcwfilm(patch().size(), Zero);
296  scalarList filmDelta(patch().size(), Zero);
297 
298  const pyrolysisModelType& pyrolysis = pyrModel();
299  const filmModelType& film = filmModel();
300 
301  // Obtain Rad heat (qr)
302  scalarField qr(patch().size(), Zero);
303 
304  label coupledPatchi = -1;
305  if (pyrolysisRegionName_ == mesh.name())
306  {
307  coupledPatchi = patchi;
308  if (qrName_ != "none")
309  {
310  qr = nbrPatch.lookupPatchField<volScalarField, scalar>(qrName_);
311  mpp.distribute(qr);
312  }
313  }
314  else if (pyrolysis.primaryMesh().name() == mesh.name())
315  {
316  coupledPatchi = nbrPatch.index();
317  if (qrName_ != "none")
318  {
319  qr = patch().lookupPatchField<volScalarField, scalar>(qrName_);
320  }
321  }
322  else
323  {
325  << type() << " condition is intended to be applied to either the "
326  << "primary or pyrolysis regions only"
327  << exit(FatalError);
328  }
329 
330  const label filmPatchi = pyrolysis.nbrCoupledPatchID(film, coupledPatchi);
331 
332  const scalarField htcw(film.htcw().h()().boundaryField()[filmPatchi]);
333 
334  // Obtain htcw
335  htcwfilm =
336  pyrolysis.mapRegionPatchField
337  (
338  film,
339  coupledPatchi,
340  filmPatchi,
341  htcw,
342  true
343  );
344 
345 
346  // Obtain Tfilm at the boundary through Ts.
347  // NOTE: Tf is not good as at the boundary it will retrieve Tp
348  Tfilm = film.Ts().boundaryField()[filmPatchi];
349  film.toPrimary(filmPatchi, Tfilm);
350 
351  // Obtain delta
352  filmDelta =
353  pyrolysis.mapRegionPatchField<scalar>
354  (
355  film,
356  "deltaf",
357  coupledPatchi,
358  true
359  );
360 
361  // Estimate wetness of the film (1: wet , 0: dry)
362  scalarField ratio
363  (
364  min
365  (
366  max
367  (
368  (filmDelta - filmDeltaDry_)/(filmDeltaWet_ - filmDeltaDry_),
369  scalar(0)
370  ),
371  scalar(1)
372  )
373  );
374 
375  scalarField qConv(ratio*htcwfilm*(Tfilm - Tp)*convectiveScaling_);
376 
377  scalarField qRad((1.0 - ratio)*qr);
378 
379  scalarField alpha(KDeltaNbr - (qRad + qConv)/Tp);
380 
381  valueFraction() = alpha/(alpha + (1.0 - ratio)*myKDelta);
382 
383  refValue() = ratio*Tfilm + (1.0 - ratio)*(KDeltaNbr*nbrIntFld)/alpha;
384 
385  mixedFvPatchScalarField::updateCoeffs();
386 
387  if (debug)
388  {
389  scalar Qc = gSum(qConv*patch().magSf());
390  scalar qr = gSum(qRad*patch().magSf());
391  scalar Qt = gSum((qConv + qRad)*patch().magSf());
392 
393  Info<< mesh.name() << ':'
394  << patch().name() << ':'
395  << this->internalField().name() << " <- "
396  << nbrMesh.name() << ':'
397  << nbrPatch.name() << ':'
398  << this->internalField().name() << " :" << nl
399  << " convective heat[W] : " << Qc << nl
400  << " radiative heat [W] : " << qr << nl
401  << " total heat [W] : " << Qt << nl
402  << " wall temperature "
403  << " min:" << gMin(*this)
404  << " max:" << gMax(*this)
405  << " avg:" << gAverage(*this)
406  << endl;
407  }
408 }
409 
410 
412 (
413  Ostream& os
414 ) const
415 {
417  os.writeEntryIfDifferent<word>
418  (
419  "filmRegion",
420  "surfaceFilmProperties",
421  filmRegionName_
422  );
423  os.writeEntryIfDifferent<word>
424  (
425  "pyrolysisRegion",
426  "pyrolysisProperties",
427  pyrolysisRegionName_
428  );
429  os.writeEntry("Tnbr", TnbrName_);
430  os.writeEntry("qr", qrName_);
431  os.writeEntry("convectiveScaling", convectiveScaling_);
432  os.writeEntry("filmDeltaDry", filmDeltaDry_);
433  os.writeEntry("filmDeltaWet", filmDeltaWet_);
435 }
436 
437 
438 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
439 
441 (
443  filmPyrolysisRadiativeCoupledMixedFvPatchScalarField
444 );
445 
446 
447 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
448 
449 } // End namespace Foam
450 
451 
452 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
faceListList boundary
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Mixed boundary condition for temperature, to be used in the flow and pyrolysis regions when a film re...
dictionary dict
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
Type gMin(const FieldField< Field, Type > &f)
tmp< scalarField > K() const
Get corresponding K field.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
Type & refCast(U &obj)
A dynamic_cast (for references) that generates FatalError on failed casts, uses the virtual type() me...
Definition: typeInfo.H:151
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
virtual tmp< scalarField > alpha(const scalarField &Tp) const
Given patch temperature calculate corresponding alphaEff field.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
void write(Ostream &os) const
Write.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
CGAL::Exact_predicates_exact_constructions_kernel K
filmPyrolysisRadiativeCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:752
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:61
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Common functions used in temperature coupled boundaries.
Type gAverage(const FieldField< Field, Type > &f)
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:352
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
HashTable< const Type * > lookupClass(const bool strict=false) const
Return all objects with a class satisfying isA<Type>
const Time & time() const noexcept
Return the reference to the time database.
Definition: regionModel.H:244
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157