Loading...
Searching...
No Matches
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-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
35namespace Foam
36{
37
38// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39
41filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::filmModel() const
42{
43 typedef filmModelType ModelType;
44
45 const UPtrList<const ModelType> models
46 (
47 db().time().csorted<ModelType>()
48 );
49
50 for (const ModelType& mdl : models)
51 {
52 if (mdl.regionMesh().name() == filmRegionName_)
53 {
54 return mdl;
55 }
56 }
57
58 DynamicList<word> modelNames(models.size());
59 for (const ModelType& mdl : models)
60 {
61 modelNames.push_back(mdl.regionMesh().name());
62 }
63
65 << "Unable to locate film region " << filmRegionName_
66 << ". Available regions include: " << modelNames
67 << abort(FatalError);
68
69 return models.front(); // Failure
70}
71
72
75filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::pyrModel() const
76{
77 typedef pyrolysisModelType ModelType;
78
79 const UPtrList<const ModelType> models
80 (
81 db().time().csorted<ModelType>()
82 );
83
84 for (const ModelType& mdl : models)
85 {
86 if (mdl.regionMesh().name() == pyrolysisRegionName_)
87 {
88 return mdl;
89 }
90 }
91
92 DynamicList<word> modelNames(models.size());
93 for (const ModelType& mdl : models)
94 {
95 modelNames.push_back(mdl.regionMesh().name());
96 }
97
99 << "Unable to locate pyrolysis region " << pyrolysisRegionName_
100 << ". Available regions include: " << modelNames
101 << abort(FatalError);
103 return models.front(); // Failure
104}
105
106
107// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108
111(
112 const fvPatch& p,
114)
115:
116 mixedFvPatchScalarField(p, iF),
117 temperatureCoupledBase(patch()), // default method (fluidThermo)
118 filmRegionName_("surfaceFilmProperties"),
119 pyrolysisRegionName_("pyrolysisProperties"),
120 TnbrName_("undefined-Tnbr"),
121 qrName_("undefined-qr"),
122 convectiveScaling_(1.0),
123 filmDeltaDry_(0.0),
124 filmDeltaWet_(0.0)
126 this->refValue() = 0.0;
127 this->refGrad() = 0.0;
128 this->valueFraction() = 1.0;
129}
130
131
134(
136 const fvPatch& p,
138 const fvPatchFieldMapper& mapper
139)
140:
141 mixedFvPatchScalarField(psf, p, iF, mapper),
142 temperatureCoupledBase(patch(), psf),
143 filmRegionName_(psf.filmRegionName_),
144 pyrolysisRegionName_(psf.pyrolysisRegionName_),
145 TnbrName_(psf.TnbrName_),
146 qrName_(psf.qrName_),
147 convectiveScaling_(psf.convectiveScaling_),
148 filmDeltaDry_(psf.filmDeltaDry_),
149 filmDeltaWet_(psf.filmDeltaWet_)
150{}
151
152
155(
156 const fvPatch& p,
158 const dictionary& dict
159)
160:
161 mixedFvPatchScalarField(p, iF),
163 filmRegionName_
164 (
165 dict.getOrDefault<word>("filmRegion", "surfaceFilmProperties")
166 ),
167 pyrolysisRegionName_
168 (
169 dict.getOrDefault<word>("pyrolysisRegion", "pyrolysisProperties")
170 ),
171 TnbrName_(dict.lookup("Tnbr")),
172 qrName_(dict.lookup("qr")),
173 convectiveScaling_(dict.getOrDefault<scalar>("convectiveScaling", 1)),
174 filmDeltaDry_(dict.get<scalar>("filmDeltaDry")),
175 filmDeltaWet_(dict.get<scalar>("filmDeltaWet"))
176{
177 if (!isA<mappedPatchBase>(this->patch().patch()))
178 {
180 << "' not type '" << mappedPatchBase::typeName << "'"
181 << "\n for patch " << p.name()
182 << " of field " << internalField().name()
183 << " in file " << internalField().objectPath()
184 << exit(FatalError);
185 }
186
187 this->readValueEntry(dict, IOobjectOption::MUST_READ);
188
189 if (this->readMixedEntries(dict))
190 {
191 // Full restart
192 }
193 else
194 {
195 // Start from user entered data. Assume fixedValue.
196 refValue() = *this;
197 refGrad() = 0.0;
198 valueFraction() = 1.0;
199 }
200}
201
202
205(
206 const filmPyrolysisRadiativeCoupledMixedFvPatchScalarField& psf,
207 const DimensionedField<scalar, volMesh>& iF
208)
209:
210 mixedFvPatchScalarField(psf, iF),
211 temperatureCoupledBase(patch(), psf),
212 filmRegionName_(psf.filmRegionName_),
213 pyrolysisRegionName_(psf.pyrolysisRegionName_),
214 TnbrName_(psf.TnbrName_),
215 qrName_(psf.qrName_),
216 convectiveScaling_(psf.convectiveScaling_),
217 filmDeltaDry_(psf.filmDeltaDry_),
218 filmDeltaWet_(psf.filmDeltaWet_)
219{}
220
221
222// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
223
225(
226 const fvPatchFieldMapper& mapper
228{
229 mixedFvPatchScalarField::autoMap(mapper);
231}
232
233
235(
236 const fvPatchField<scalar>& ptf,
237 const labelList& addr
238)
239{
240 mixedFvPatchScalarField::rmap(ptf, addr);
242 const auto& fpptf = refCast<const myType>(ptf);
243
244 temperatureCoupledBase::rmap(fpptf, addr);
245}
246
247
249{
250 if (updated())
251 {
252 return;
253 }
254
255 // Get the coupling information from the mappedPatchBase
256 const auto& mpp = refCast<const mappedPatchBase>(patch().patch());
257
258 const label patchi = patch().index();
259 const label nbrPatchi = mpp.samplePolyPatch().index();
260 const polyMesh& mesh = patch().boundaryMesh().mesh();
261 const polyMesh& nbrMesh = mpp.sampleMesh();
262 const fvPatch& nbrPatch =
263 refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchi];
264
265 scalarField intFld(patchInternalField());
266
267 const auto& nbrField =
269 (
270 nbrPatch.lookupPatchField<volScalarField>(TnbrName_)
271 );
272
273 // Swap to obtain full local values of neighbour internal field
274 scalarField nbrIntFld(nbrField.patchInternalField());
275 mpp.distribute(nbrIntFld);
276
277 scalarField& Tp = *this;
278
279 const scalarField K(this->kappa(*this));
280 const scalarField nbrK(nbrField.kappa(*this));
281
282 // Swap to obtain full local values of neighbour K*delta
283 scalarField KDeltaNbr(nbrK*nbrPatch.deltaCoeffs());
284 mpp.distribute(KDeltaNbr);
285
286 const scalarField myKDelta(K*patch().deltaCoeffs());
287
288 scalarList Tfilm(patch().size(), Zero);
289 scalarList htcwfilm(patch().size(), Zero);
290 scalarList filmDelta(patch().size(), Zero);
291
292 const pyrolysisModelType& pyrolysis = pyrModel();
293 const filmModelType& film = filmModel();
294
295 // Obtain Rad heat (qr)
296 scalarField qr(patch().size(), Zero);
297
298 label coupledPatchi = -1;
299 label filmPatchi = -1;
300 if (pyrolysisRegionName_ == mesh.name())
301 {
302 // Working on the pyrolysis mesh
303
304 coupledPatchi = patchi;
305 if (qrName_ != "none")
306 {
307 qr = nbrPatch.lookupPatchField<volScalarField>(qrName_);
308 mpp.distribute(qr);
309 }
310
311 filmPatchi = pyrolysis.nbrCoupledPatchID(film, coupledPatchi);
312
313 const scalarField htcw(film.htcw().h()().boundaryField()[filmPatchi]);
314
315 // Obtain htcw
316 htcwfilm =
317 pyrolysis.mapRegionPatchField
318 (
319 film,
320 coupledPatchi,
321 filmPatchi,
322 htcw,
323 true
324 );
325
326 // Obtain Tfilm at the boundary through Ts.
327 // NOTE: Tf is not good as at the boundary it will retrieve Tp
328 const scalarField Ts(film.Ts().boundaryField()[filmPatchi]);
329 Tfilm =
330 pyrolysis.mapRegionPatchField
331 (
332 film,
333 coupledPatchi,
334 filmPatchi,
335 Ts,
336 true
337 );
338
339 // Obtain delta
340 filmDelta =
341 pyrolysis.mapRegionPatchField<scalar>
342 (
343 film,
344 "deltaf",
345 coupledPatchi,
346 true
347 );
348 }
349 else if (pyrolysis.primaryMesh().name() == mesh.name())
350 {
351 // Working on the primary mesh
352
353 coupledPatchi = nbrPatch.index();
354 if (qrName_ != "none")
355 {
356 qr = patch().lookupPatchField<volScalarField>(qrName_);
357 }
358
359 filmPatchi = pyrolysis.nbrCoupledPatchID(film, coupledPatchi);
360
361 htcwfilm = film.htcw().h()().boundaryField()[filmPatchi];
362 film.toPrimary(filmPatchi, htcwfilm);
363
364 // Obtain Tfilm at the boundary through Ts.
365 // NOTE: Tf is not good as at the boundary it will retrieve Tp
366 Tfilm = film.Ts().boundaryField()[filmPatchi];
367 film.toPrimary(filmPatchi, Tfilm);
368
369 filmDelta = film.delta().boundaryField()[filmPatchi];
370 film.toPrimary(filmPatchi, filmDelta);
371 }
372 else
373 {
375 << type() << " condition is intended to be applied to either the "
376 << "primary or pyrolysis regions only"
377 << exit(FatalError);
378 }
379
380 // Estimate wetness of the film (1: wet , 0: dry)
381 const scalarField ratio
382 (
383 min
384 (
385 max
386 (
387 (filmDelta - filmDeltaDry_)/(filmDeltaWet_ - filmDeltaDry_),
388 scalar(0)
389 ),
390 scalar(1)
391 )
392 );
393
394 const scalarField qConv(ratio*htcwfilm*(Tfilm - Tp)*convectiveScaling_);
395 const scalarField qRad((1.0 - ratio)*qr);
396 const scalarField alpha(KDeltaNbr - (qRad + qConv)/Tp);
397
398 valueFraction() = alpha/(alpha + (1.0 - ratio)*myKDelta);
399 refValue() = ratio*Tfilm + (1.0 - ratio)*(KDeltaNbr*nbrIntFld)/alpha;
400
401 mixedFvPatchScalarField::updateCoeffs();
402
403 if (debug)
404 {
405 scalar Qc = gWeightedSum(patch().magSf(), qConv);
406 scalar qr = gWeightedSum(patch().magSf(), qRad);
407 scalar Qt = gSum((qConv + qRad)*patch().magSf());
408
409 MinMax<scalar> limits = gMinMax(*this);
410
411 Info<< mesh.name() << ':'
412 << patch().name() << ':'
413 << this->internalField().name() << " <- "
414 << nbrMesh.name() << ':'
415 << nbrPatch.name() << ':'
416 << this->internalField().name() << " :" << nl
417 << " convective heat[W] : " << Qc << nl
418 << " radiative heat [W] : " << qr << nl
419 << " total heat [W] : " << Qt << nl
420 << " wall temperature "
421 << " min:" << limits.min()
422 << " max:" << limits.max()
423 << " avg:" << gAverage(*this)
424 << endl;
425 }
426}
427
428
430(
431 Ostream& os
432) const
433{
435 os.writeEntryIfDifferent<word>
436 (
437 "filmRegion",
438 "surfaceFilmProperties",
439 filmRegionName_
440 );
441 os.writeEntryIfDifferent<word>
442 (
443 "pyrolysisRegion",
444 "pyrolysisProperties",
445 pyrolysisRegionName_
446 );
447 os.writeEntry("Tnbr", TnbrName_);
448 os.writeEntry("qr", qrName_);
449 os.writeEntry("convectiveScaling", convectiveScaling_);
450 os.writeEntry("filmDeltaDry", filmDeltaDry_);
451 os.writeEntry("filmDeltaWet", filmDeltaWet_);
453}
454
455
456// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
457
459(
461 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField
462);
463
464
465// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
466
467} // End namespace Foam
468
469
470// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ MUST_READ
Reading required.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition MinMax.H:106
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mixed boundary condition for temperature, to be used in the flow and pyrolysis regions when a film re...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
filmPyrolysisRadiativeCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
tmp< scalarField > K() const
Get corresponding K field.
const word & name() const
Return reference to name.
Definition fvMesh.H:387
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual const word & name() const
Return name.
Definition fvPatch.H:210
label index() const noexcept
The index of this patch in the boundary mesh.
Definition fvPatch.H:218
const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient except for coupled patches for which the cell-centre to c...
Definition fvPatch.C:183
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup the named field from the local registry and return the patch field corresponding to this patch...
virtual void write(Ostream &) const
Write.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Lookup type of boundary radiation properties.
Definition lookup.H:60
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
Return the coupled patch ID paired with coupled patch.
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
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 tmp< volScalarField > h() const =0
Return the heat transfer coefficient [W/m2/K].
const volScalarField & delta() const
Return const access to the film thickness [m].
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
Common functions used in temperature coupled boundaries.
virtual tmp< scalarField > alpha(const scalarField &Tp) const
Given patch temperature calculate corresponding alphaEff field.
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
void write(Ostream &os) const
Write.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
temperatureCoupledBase(const fvPatch &patch, const KMethodType method=KMethodType::mtFluidThermo)
Default construct from patch, using fluidThermo (default) or specified method.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
Type gSum(const FieldField< Field, Type > &f)
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
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
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict