filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2018 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 "mappedPatchBase.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
38 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
39 filmModel() const
40 {
41  HashTable<const filmModelType*> models
42  = db().time().lookupClass<filmModelType>();
43 
44  forAllConstIter(HashTable<const filmModelType*>, models, iter)
45  {
46  if (iter()->regionMesh().name() == filmRegionName_)
47  {
48  return *iter();
49  }
50  }
51 
52  DynamicList<word> modelNames;
53  forAllConstIter(HashTable<const filmModelType*>, models, iter)
54  {
55  modelNames.append(iter()->regionMesh().name());
56  }
57 
59  << "Unable to locate film region " << filmRegionName_
60  << ". Available regions include: " << modelNames
61  << abort(FatalError);
62 
63  return **models.begin();
64 }
65 
66 
69 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
70 pyrModel() const
71 {
72  HashTable<const pyrolysisModelType*> models =
74 
75  forAllConstIter(HashTable<const pyrolysisModelType*>, models, iter)
76  {
77  if (iter()->regionMesh().name() == pyrolysisRegionName_)
78  {
79  return *iter();
80  }
81  }
82 
83  DynamicList<word> modelNames;
84  forAllConstIter(HashTable<const pyrolysisModelType*>, models, iter)
85  {
86  modelNames.append(iter()->regionMesh().name());
87  }
88 
89 
91  << "Unable to locate pyrolysis region " << pyrolysisRegionName_
92  << ". Available regions include: " << modelNames
93  << abort(FatalError);
94 
95  return **models.begin();
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100 
103 (
104  const fvPatch& p,
106 )
107 :
108  mixedFvPatchScalarField(p, iF),
109  temperatureCoupledBase(patch(), "undefined", "undefined", "undefined-K"),
110  filmRegionName_("surfaceFilmProperties"),
111  pyrolysisRegionName_("pyrolysisProperties"),
112  TnbrName_("undefined-Tnbr"),
113  qrName_("undefined-qr"),
114  convectiveScaling_(1.0),
115  filmDeltaDry_(0.0),
116  filmDeltaWet_(0.0)
117 {
118  this->refValue() = 0.0;
119  this->refGrad() = 0.0;
120  this->valueFraction() = 1.0;
121 }
122 
123 
126 (
128  const fvPatch& p,
130  const fvPatchFieldMapper& mapper
131 )
132 :
133  mixedFvPatchScalarField(psf, p, iF, mapper),
134  temperatureCoupledBase(patch(), psf),
135  filmRegionName_(psf.filmRegionName_),
136  pyrolysisRegionName_(psf.pyrolysisRegionName_),
137  TnbrName_(psf.TnbrName_),
138  qrName_(psf.qrName_),
139  convectiveScaling_(psf.convectiveScaling_),
140  filmDeltaDry_(psf.filmDeltaDry_),
141  filmDeltaWet_(psf.filmDeltaWet_)
142 {}
143 
144 
147 (
148  const fvPatch& p,
150  const dictionary& dict
151 )
152 :
153  mixedFvPatchScalarField(p, iF),
154  temperatureCoupledBase(patch(), dict),
155  filmRegionName_
156  (
157  dict.lookupOrDefault<word>("filmRegion", "surfaceFilmProperties")
158  ),
159  pyrolysisRegionName_
160  (
161  dict.lookupOrDefault<word>("pyrolysisRegion", "pyrolysisProperties")
162  ),
163  TnbrName_(dict.lookup("Tnbr")),
164  qrName_(dict.lookup("qr")),
165  convectiveScaling_(dict.lookupOrDefault<scalar>("convectiveScaling", 1.0)),
166  filmDeltaDry_(readScalar(dict.lookup("filmDeltaDry"))),
167  filmDeltaWet_(readScalar(dict.lookup("filmDeltaWet")))
168 {
169  if (!isA<mappedPatchBase>(this->patch().patch()))
170  {
172  << "' not type '" << mappedPatchBase::typeName << "'"
173  << "\n for patch " << p.name()
174  << " of field " << internalField().name()
175  << " in file " << internalField().objectPath()
176  << exit(FatalError);
177  }
178 
179  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
180 
181  if (dict.found("refValue"))
182  {
183  // Full restart
184  refValue() = scalarField("refValue", dict, p.size());
185  refGrad() = scalarField("refGradient", dict, p.size());
186  valueFraction() = scalarField("valueFraction", dict, p.size());
187  }
188  else
189  {
190  // Start from user entered data. Assume fixedValue.
191  refValue() = *this;
192  refGrad() = 0.0;
193  valueFraction() = 1.0;
194  }
195 }
196 
197 
200 (
203 )
204 :
205  mixedFvPatchScalarField(psf, iF),
206  temperatureCoupledBase(patch(), psf),
207  filmRegionName_(psf.filmRegionName_),
208  pyrolysisRegionName_(psf.pyrolysisRegionName_),
209  TnbrName_(psf.TnbrName_),
210  qrName_(psf.qrName_),
211  convectiveScaling_(psf.convectiveScaling_),
212  filmDeltaDry_(psf.filmDeltaDry_),
213  filmDeltaWet_(psf.filmDeltaWet_)
214 {}
215 
216 
217 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218 
220 {
221  if (updated())
222  {
223  return;
224  }
225 
226  // Get the coupling information from the mappedPatchBase
227  const mappedPatchBase& mpp =
228  refCast<const mappedPatchBase>(patch().patch());
229 
230  const label patchi = patch().index();
231  const label nbrPatchi = mpp.samplePolyPatch().index();
232  const polyMesh& mesh = patch().boundaryMesh().mesh();
233  const polyMesh& nbrMesh = mpp.sampleMesh();
234  const fvPatch& nbrPatch =
235  refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchi];
236 
237  scalarField intFld(patchInternalField());
238 
240  nbrField =
241  refCast
242  <
244  >
245  (
246  nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
247  );
248 
249  // Swap to obtain full local values of neighbour internal field
250  scalarField nbrIntFld(nbrField.patchInternalField());
251  mpp.distribute(nbrIntFld);
252 
253  scalarField& Tp = *this;
254 
255  const scalarField K(this->kappa(*this));
256  const scalarField nbrK(nbrField.kappa(*this));
257 
258  // Swap to obtain full local values of neighbour K*delta
259  scalarField KDeltaNbr(nbrK*nbrPatch.deltaCoeffs());
260  mpp.distribute(KDeltaNbr);
261 
262  scalarField myKDelta(K*patch().deltaCoeffs());
263 
264  scalarList Tfilm(patch().size(), 0.0);
265  scalarList htcwfilm(patch().size(), 0.0);
266  scalarList filmDelta(patch().size(), 0.0);
267 
268  const pyrolysisModelType& pyrolysis = pyrModel();
269  const filmModelType& film = filmModel();
270 
271  // Obtain Rad heat (qr)
272  scalarField qr(patch().size(), 0.0);
273 
274  label coupledPatchi = -1;
275  if (pyrolysisRegionName_ == mesh.name())
276  {
277  coupledPatchi = patchi;
278  if (qrName_ != "none")
279  {
280  qr = nbrPatch.lookupPatchField<volScalarField, scalar>(qrName_);
281  mpp.distribute(qr);
282  }
283  }
284  else if (pyrolysis.primaryMesh().name() == mesh.name())
285  {
286  coupledPatchi = nbrPatch.index();
287  if (qrName_ != "none")
288  {
289  qr = patch().lookupPatchField<volScalarField, scalar>(qrName_);
290  }
291  }
292  else
293  {
295  << type() << " condition is intended to be applied to either the "
296  << "primary or pyrolysis regions only"
297  << exit(FatalError);
298  }
299 
300  const label filmPatchi = pyrolysis.nbrCoupledPatchID(film, coupledPatchi);
301 
302  const scalarField htcw(film.htcw().h()().boundaryField()[filmPatchi]);
303 
304  // Obtain htcw
305  htcwfilm =
306  pyrolysis.mapRegionPatchField
307  (
308  film,
309  coupledPatchi,
310  filmPatchi,
311  htcw,
312  true
313  );
314 
315 
316  // Obtain Tfilm at the boundary through Ts.
317  // NOTE: Tf is not good as at the boundary it will retrieve Tp
318  Tfilm = film.Ts().boundaryField()[filmPatchi];
319  film.toPrimary(filmPatchi, Tfilm);
320 
321  // Obtain delta
322  filmDelta =
323  pyrolysis.mapRegionPatchField<scalar>
324  (
325  film,
326  "deltaf",
327  coupledPatchi,
328  true
329  );
330 
331  // Estimate wetness of the film (1: wet , 0: dry)
332  scalarField ratio
333  (
334  min
335  (
336  max
337  (
338  (filmDelta - filmDeltaDry_)/(filmDeltaWet_ - filmDeltaDry_),
339  scalar(0)
340  ),
341  scalar(1)
342  )
343  );
344 
345  scalarField qConv(ratio*htcwfilm*(Tfilm - Tp)*convectiveScaling_);
346 
347  scalarField qRad((1.0 - ratio)*qr);
348 
349  scalarField alpha(KDeltaNbr - (qRad + qConv)/Tp);
350 
351  valueFraction() = alpha/(alpha + (1.0 - ratio)*myKDelta);
352 
353  refValue() = ratio*Tfilm + (1.0 - ratio)*(KDeltaNbr*nbrIntFld)/alpha;
354 
355  mixedFvPatchScalarField::updateCoeffs();
356 
357  if (debug)
358  {
359  scalar Qc = gSum(qConv*patch().magSf());
360  scalar qr = gSum(qRad*patch().magSf());
361  scalar Qt = gSum((qConv + qRad)*patch().magSf());
362 
363  Info<< mesh.name() << ':'
364  << patch().name() << ':'
365  << this->internalField().name() << " <- "
366  << nbrMesh.name() << ':'
367  << nbrPatch.name() << ':'
368  << this->internalField().name() << " :" << nl
369  << " convective heat[W] : " << Qc << nl
370  << " radiative heat [W] : " << qr << nl
371  << " total heat [W] : " << Qt << nl
372  << " wall temperature "
373  << " min:" << gMin(*this)
374  << " max:" << gMax(*this)
375  << " avg:" << gAverage(*this)
376  << endl;
377  }
378 }
379 
380 
382 (
383  Ostream& os
384 ) const
385 {
387  writeEntryIfDifferent<word>
388  (
389  os,
390  "filmRegion",
391  "surfaceFilmProperties",
392  filmRegionName_
393  );
394  writeEntryIfDifferent<word>
395  (
396  os,
397  "pyrolysisRegion",
398  "pyrolysisProperties",
399  pyrolysisRegionName_
400  );
401  os.writeKeyword("Tnbr")<< TnbrName_ << token::END_STATEMENT << nl;
402  os.writeKeyword("qr")<< qrName_ << token::END_STATEMENT << nl;
403  os.writeKeyword("convectiveScaling") << convectiveScaling_
404  << token::END_STATEMENT << nl;
405  os.writeKeyword("filmDeltaDry") << filmDeltaDry_
406  << token::END_STATEMENT << nl;
407  os.writeKeyword("filmDeltaWet") << filmDeltaWet_
410 }
411 
412 
413 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
414 
416 (
419 );
420 
421 
422 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
423 
424 } // End namespace Foam
425 
426 
427 // ************************************************************************* //
Mixed boundary condition for temperature, to be used in the flow and pyrolysis regions when a film re...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:297
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Type gMin(const FieldField< Field, Type > &f)
const polyMesh & sampleMesh() const
Get the region mesh.
tmp< scalarField > K() const
Get corresponding K field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
virtual tmp< volScalarField > h() const =0
Return the heat transfer coefficient [W/m2/K].
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:256
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:66
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
filmPyrolysisRadiativeCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const word & name() const
Return name.
Definition: fvPatch.H:149
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:31
dynamicFvMesh & mesh
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.
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
Return the coupled patch ID paired with coupled patch.
Definition: regionModel.C:278
Foam::fvPatchFieldMapper.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
faceListList boundary(nPatches)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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:265
Type gMax(const FieldField< Field, Type > &f)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Lookup and return the patchField of the named field from the.
Type gAverage(const FieldField< Field, Type > &f)
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
Definition: fvPatch.C:164
label index() const
Return the index of this patch in the boundaryMesh.
tmp< Foam::Field< Type > > mapRegionPatchField(const regionModel &nbrRegion, const label regionPatchi, const label nbrPatchi, const Field< Type > &nbrField, const bool flip=false) const
Map patch field from another region model to local patch.
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...
messageStream Info
temperatureCoupledBase(const fvPatch &patch, const word &calculationMethod, const word &kappaName, const word &alphaAniName)
Construct from patch and K name.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
void write(Ostream &) const
Write.
Thermodynamic form of single-cell layer surface film model.
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