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 | Copyright (C) 2013-2014 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  (
60  "const filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::"
61  "filmModelType& "
62  "filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::"
63  "filmModel() const"
64  )
65  << "Unable to locate film region " << filmRegionName_
66  << ". Available regions include: " << modelNames
67  << abort(FatalError);
68 
69  return **models.begin();
70 }
71 
72 
75 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
76 pyrModel() const
77 {
78  HashTable<const pyrolysisModelType*> models =
80 
81  forAllConstIter(HashTable<const pyrolysisModelType*>, models, iter)
82  {
83  if (iter()->regionMesh().name() == pyrolysisRegionName_)
84  {
85  return *iter();
86  }
87  }
88 
89  DynamicList<word> modelNames;
90  forAllConstIter(HashTable<const pyrolysisModelType*>, models, iter)
91  {
92  modelNames.append(iter()->regionMesh().name());
93  }
94 
95 
97  (
98  "const filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::"
99  "pyrolysisModelType& "
100  "filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::"
101  "pyrModel() const"
102  )
103  << "Unable to locate pyrolysis region " << pyrolysisRegionName_
104  << ". Available regions include: " << modelNames
105  << abort(FatalError);
106 
107  return **models.begin();
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
112 
115 (
116  const fvPatch& p,
118 )
119 :
120  mixedFvPatchScalarField(p, iF),
121  temperatureCoupledBase(patch(), "undefined", "undefined", "undefined-K"),
122  filmRegionName_("surfaceFilmProperties"),
123  pyrolysisRegionName_("pyrolysisProperties"),
124  TnbrName_("undefined-Tnbr"),
125  QrName_("undefined-Qr"),
126  convectiveScaling_(1.0),
127  filmDeltaDry_(0.0),
128  filmDeltaWet_(0.0)
129 {
130  this->refValue() = 0.0;
131  this->refGrad() = 0.0;
132  this->valueFraction() = 1.0;
133 }
134 
135 
138 (
140  const fvPatch& p,
142  const fvPatchFieldMapper& mapper
143 )
144 :
145  mixedFvPatchScalarField(psf, p, iF, mapper),
146  temperatureCoupledBase(patch(), psf),
147  filmRegionName_(psf.filmRegionName_),
148  pyrolysisRegionName_(psf.pyrolysisRegionName_),
149  TnbrName_(psf.TnbrName_),
150  QrName_(psf.QrName_),
151  convectiveScaling_(psf.convectiveScaling_),
152  filmDeltaDry_(psf.filmDeltaDry_),
153  filmDeltaWet_(psf.filmDeltaWet_)
154 {}
155 
156 
159 (
160  const fvPatch& p,
162  const dictionary& dict
163 )
164 :
165  mixedFvPatchScalarField(p, iF),
166  temperatureCoupledBase(patch(), dict),
167  filmRegionName_
168  (
169  dict.lookupOrDefault<word>("filmRegion", "surfaceFilmProperties")
170  ),
171  pyrolysisRegionName_
172  (
173  dict.lookupOrDefault<word>("pyrolysisRegion", "pyrolysisProperties")
174  ),
175  TnbrName_(dict.lookup("Tnbr")),
176  QrName_(dict.lookup("Qr")),
177  convectiveScaling_(dict.lookupOrDefault<scalar>("convectiveScaling", 1.0)),
178  filmDeltaDry_(readScalar(dict.lookup("filmDeltaDry"))),
179  filmDeltaWet_(readScalar(dict.lookup("filmDeltaWet")))
180 {
181  if (!isA<mappedPatchBase>(this->patch().patch()))
182  {
184  (
185  "filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::"
186  "filmPyrolysisRadiativeCoupledMixedFvPatchScalarField\n"
187  "(\n"
188  " const fvPatch& p,\n"
189  " const DimensionedField<scalar, volMesh>& iF,\n"
190  " const dictionary& dict\n"
191  ")\n"
192  ) << "\n patch type '" << p.type()
193  << "' not type '" << mappedPatchBase::typeName << "'"
194  << "\n for patch " << p.name()
195  << " of field " << dimensionedInternalField().name()
196  << " in file " << dimensionedInternalField().objectPath()
197  << exit(FatalError);
198  }
199 
200  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
201 
202  if (dict.found("refValue"))
203  {
204  // Full restart
205  refValue() = scalarField("refValue", dict, p.size());
206  refGrad() = scalarField("refGradient", dict, p.size());
207  valueFraction() = scalarField("valueFraction", dict, p.size());
208  }
209  else
210  {
211  // Start from user entered data. Assume fixedValue.
212  refValue() = *this;
213  refGrad() = 0.0;
214  valueFraction() = 1.0;
215  }
216 }
217 
218 
221 (
224 )
225 :
226  mixedFvPatchScalarField(psf, iF),
227  temperatureCoupledBase(patch(), psf),
228  filmRegionName_(psf.filmRegionName_),
229  pyrolysisRegionName_(psf.pyrolysisRegionName_),
230  TnbrName_(psf.TnbrName_),
231  QrName_(psf.QrName_),
232  convectiveScaling_(psf.convectiveScaling_),
233  filmDeltaDry_(psf.filmDeltaDry_),
234  filmDeltaWet_(psf.filmDeltaWet_)
235 {}
236 
237 
238 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
239 
241 {
242  if (updated())
243  {
244  return;
245  }
246 
247  // Get the coupling information from the mappedPatchBase
248  const mappedPatchBase& mpp =
249  refCast<const mappedPatchBase>(patch().patch());
250 
251  const label patchI = patch().index();
252  const label nbrPatchI = mpp.samplePolyPatch().index();
253  const polyMesh& mesh = patch().boundaryMesh().mesh();
254  const polyMesh& nbrMesh = mpp.sampleMesh();
255  const fvPatch& nbrPatch =
256  refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchI];
257 
258  scalarField intFld(patchInternalField());
259 
261  nbrField =
262  refCast
263  <
265  >
266  (
267  nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
268  );
269 
270  // Swap to obtain full local values of neighbour internal field
271  scalarField nbrIntFld(nbrField.patchInternalField());
272  mpp.distribute(nbrIntFld);
273 
274  scalarField& Tp = *this;
275 
276  const scalarField K(this->kappa(*this));
277  const scalarField nbrK(nbrField.kappa(*this));
278 
279  // Swap to obtain full local values of neighbour K*delta
280  scalarField KDeltaNbr(nbrK*nbrPatch.deltaCoeffs());
281  mpp.distribute(KDeltaNbr);
282 
283  scalarField myKDelta(K*patch().deltaCoeffs());
284 
285  scalarList Tfilm(patch().size(), 0.0);
286  scalarList htcwfilm(patch().size(), 0.0);
287  scalarList filmDelta(patch().size(), 0.0);
288 
289  const pyrolysisModelType& pyrolysis = pyrModel();
290  const filmModelType& film = filmModel();
291 
292  // Obtain Rad heat (Qr)
293  scalarField Qr(patch().size(), 0.0);
294 
295  label coupledPatchI = -1;
296  if (pyrolysisRegionName_ == mesh.name())
297  {
298  coupledPatchI = patchI;
299  if (QrName_ != "none")
300  {
301  Qr = nbrPatch.lookupPatchField<volScalarField, scalar>(QrName_);
302  mpp.distribute(Qr);
303  }
304  }
305  else if (pyrolysis.primaryMesh().name() == mesh.name())
306  {
307  coupledPatchI = nbrPatch.index();
308  if (QrName_ != "none")
309  {
310  Qr = patch().lookupPatchField<volScalarField, scalar>(QrName_);
311  }
312  }
313  else
314  {
316  (
317  "void filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::"
318  "updateCoeffs()"
319  )
320  << type() << " condition is intended to be applied to either the "
321  << "primary or pyrolysis regions only"
322  << exit(FatalError);
323  }
324 
325  const label filmPatchI = pyrolysis.nbrCoupledPatchID(film, coupledPatchI);
326 
327  const scalarField htcw(film.htcw().h()().boundaryField()[filmPatchI]);
328 
329  // Obtain htcw
330  htcwfilm =
331  pyrolysis.mapRegionPatchField
332  (
333  film,
334  coupledPatchI,
335  filmPatchI,
336  htcw,
337  true
338  );
339 
340 
341  // Obtain Tfilm at the boundary through Ts.
342  // NOTE: Tf is not good as at the boundary it will retrieve Tp
343  Tfilm = film.Ts().boundaryField()[filmPatchI];
344  film.toPrimary(filmPatchI, Tfilm);
345 
346  // Obtain delta
347  filmDelta =
348  pyrolysis.mapRegionPatchField<scalar>
349  (
350  film,
351  "deltaf",
352  coupledPatchI,
353  true
354  );
355 
356  // Estimate wetness of the film (1: wet , 0: dry)
357  scalarField ratio
358  (
359  min
360  (
361  max
362  (
363  (filmDelta - filmDeltaDry_)/(filmDeltaWet_ - filmDeltaDry_),
364  scalar(0.0)
365  ),
366  scalar(1.0)
367  )
368  );
369 
370  scalarField qConv(ratio*htcwfilm*(Tfilm - Tp)*convectiveScaling_);
371 
372  scalarField qRad((1.0 - ratio)*Qr);
373 
374  scalarField alpha(KDeltaNbr - (qRad + qConv)/Tp);
375 
376  valueFraction() = alpha/(alpha + (1.0 - ratio)*myKDelta);
377 
378  refValue() = ratio*Tfilm + (1.0 - ratio)*(KDeltaNbr*nbrIntFld)/alpha;
379 
380  mixedFvPatchScalarField::updateCoeffs();
381 
382  if (debug)
383  {
384  scalar Qc = gSum(qConv*patch().magSf());
385  scalar Qr = gSum(qRad*patch().magSf());
386  scalar Qt = gSum((qConv + qRad)*patch().magSf());
387 
388  Info<< mesh.name() << ':'
389  << patch().name() << ':'
390  << this->dimensionedInternalField().name() << " <- "
391  << nbrMesh.name() << ':'
392  << nbrPatch.name() << ':'
393  << this->dimensionedInternalField().name() << " :" << nl
394  << " convective heat[W] : " << Qc << nl
395  << " radiative heat [W] : " << Qr << nl
396  << " total heat [W] : " << Qt << nl
397  << " wall temperature "
398  << " min:" << gMin(*this)
399  << " max:" << gMax(*this)
400  << " avg:" << gAverage(*this)
401  << endl;
402  }
403 }
404 
405 
407 (
408  Ostream& os
409 ) const
410 {
412  writeEntryIfDifferent<word>
413  (
414  os,
415  "filmRegion",
416  "surfaceFilmProperties",
417  filmRegionName_
418  );
419  writeEntryIfDifferent<word>
420  (
421  os,
422  "pyrolysisRegion",
423  "pyrolysisProperties",
424  pyrolysisRegionName_
425  );
426  os.writeKeyword("Tnbr")<< TnbrName_ << token::END_STATEMENT << nl;
427  os.writeKeyword("Qr")<< QrName_ << token::END_STATEMENT << nl;
428  os.writeKeyword("convectiveScaling") << convectiveScaling_
429  << token::END_STATEMENT << nl;
430  os.writeKeyword("filmDeltaDry") << filmDeltaDry_
431  << token::END_STATEMENT << nl;
432  os.writeKeyword("filmDeltaWet") << filmDeltaWet_
435 }
436 
437 
438 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
439 
441 (
444 );
445 
446 
447 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
448 
449 } // End namespace Foam
450 
451 
452 // ************************************************************************* //
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
faceListList boundary(nPatches)
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
filmPyrolysisRadiativeCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const polyMesh & sampleMesh() const
Get the region mesh.
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
Definition: fvPatch.C:164
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
tmp< scalarField > K() const
Get corresponding K field.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
label index() const
Return the index of this patch in the boundaryMesh.
A class for handling words, derived from string.
Definition: word.H:59
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
void write(Ostream &) const
Write.
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchI) const
Return the coupled patch ID paired with coupled patch.
Definition: regionModel.C:303
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::fvPatchFieldMapper.
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
messageStream Info
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
Namespace for OpenFOAM.
const word & name() const
Return name.
Definition: fvPatch.H:149
runTime write()
const GeometricField::PatchFieldType & lookupPatchField(const word &name, const GeometricField *=NULL, const Type *=NULL) const
Lookup and return the patchField of the named field from the.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Type gSum(const FieldField< Field, Type > &f)
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const polyPatch & samplePolyPatch() const
Get the patch on the region.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
virtual tmp< volScalarField > h() const =0
Return the heat transfer coefficient [W/m2/K].
Macros for easy insertion into run-time selection tables.
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
const word & name() const
Return name.
Definition: IOobject.H:260
Type gMin(const FieldField< Field, Type > &f)
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
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.
rDeltaT dimensionedInternalField()
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:31
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,.
temperatureCoupledBase(const fvPatch &patch, const word &calculationMethod, const word &kappaName, const word &alphaAniName)
Construct from patch and K name.
tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual label size() const
Return size.
Definition: fvPatch.H:161
Type gMax(const FieldField< Field, Type > &f)
void toPrimary(const label regionPatchI, List< Type > &regionField) const
Convert a local region field to the primary region.