Go to the documentation of this file.
41 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
44 HashTable<const filmModelType*> models
49 if (iter()->regionMesh().
name() == filmRegionName_)
55 DynamicList<word> modelNames;
58 modelNames.append(iter()->regionMesh().
name());
62 <<
"Unable to locate film region " << filmRegionName_
63 <<
". Available regions include: " << modelNames
66 return **models.begin();
72 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
75 HashTable<const pyrolysisModelType*> models =
80 if (iter()->regionMesh().
name() == pyrolysisRegionName_)
86 DynamicList<word> modelNames;
89 modelNames.append(iter()->regionMesh().
name());
93 <<
"Unable to locate pyrolysis region " << pyrolysisRegionName_
94 <<
". Available regions include: " << modelNames
97 return **models.begin();
110 mixedFvPatchScalarField(
p, iF),
119 filmRegionName_(
"surfaceFilmProperties"),
120 pyrolysisRegionName_(
"pyrolysisProperties"),
121 TnbrName_(
"undefined-Tnbr"),
122 qrName_(
"undefined-qr"),
123 convectiveScaling_(1.0),
127 this->refValue() = 0.0;
128 this->refGrad() = 0.0;
129 this->valueFraction() = 1.0;
142 mixedFvPatchScalarField(psf,
p, iF, mapper),
144 filmRegionName_(psf.filmRegionName_),
145 pyrolysisRegionName_(psf.pyrolysisRegionName_),
146 TnbrName_(psf.TnbrName_),
147 qrName_(psf.qrName_),
148 convectiveScaling_(psf.convectiveScaling_),
149 filmDeltaDry_(psf.filmDeltaDry_),
150 filmDeltaWet_(psf.filmDeltaWet_)
162 mixedFvPatchScalarField(
p, iF),
166 dict.getOrDefault<
word>(
"filmRegion",
"surfaceFilmProperties")
170 dict.getOrDefault<
word>(
"pyrolysisRegion",
"pyrolysisProperties")
174 convectiveScaling_(
dict.getOrDefault<scalar>(
"convectiveScaling", 1)),
175 filmDeltaDry_(
dict.
get<scalar>(
"filmDeltaDry")),
176 filmDeltaWet_(
dict.
get<scalar>(
"filmDeltaWet"))
178 if (!isA<mappedPatchBase>(this->
patch().
patch()))
181 <<
"' not type '" << mappedPatchBase::typeName <<
"'"
182 <<
"\n for patch " <<
p.name()
183 <<
" of field " << internalField().name()
184 <<
" in file " << internalField().objectPath()
190 if (
dict.found(
"refValue"))
202 valueFraction() = 1.0;
210 const filmPyrolysisRadiativeCoupledMixedFvPatchScalarField& psf,
211 const DimensionedField<scalar, volMesh>& iF
214 mixedFvPatchScalarField(psf, iF),
215 temperatureCoupledBase(
patch(), psf),
216 filmRegionName_(psf.filmRegionName_),
217 pyrolysisRegionName_(psf.pyrolysisRegionName_),
218 TnbrName_(psf.TnbrName_),
219 qrName_(psf.qrName_),
220 convectiveScaling_(psf.convectiveScaling_),
221 filmDeltaDry_(psf.filmDeltaDry_),
222 filmDeltaWet_(psf.filmDeltaWet_)
233 mixedFvPatchScalarField::autoMap(mapper);
244 mixedFvPatchScalarField::rmap(ptf, addr);
264 const mappedPatchBase& mpp =
265 refCast<const mappedPatchBase>(
patch().
patch());
267 const label patchi =
patch().index();
268 const label nbrPatchi = mpp.samplePolyPatch().index();
269 const polyMesh&
mesh =
patch().boundaryMesh().mesh();
270 const polyMesh& nbrMesh = mpp.sampleMesh();
271 const fvPatch& nbrPatch =
272 refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchi];
287 scalarField nbrIntFld(nbrField.patchInternalField());
288 mpp.distribute(nbrIntFld);
296 scalarField KDeltaNbr(nbrK*nbrPatch.deltaCoeffs());
297 mpp.distribute(KDeltaNbr);
311 label coupledPatchi = -1;
312 if (pyrolysisRegionName_ ==
mesh.name())
314 coupledPatchi = patchi;
315 if (qrName_ !=
"none")
321 else if (pyrolysis.primaryMesh().name() ==
mesh.name())
323 coupledPatchi = nbrPatch.index();
324 if (qrName_ !=
"none")
332 <<
type() <<
" condition is intended to be applied to either the "
333 <<
"primary or pyrolysis regions only"
337 const label filmPatchi = pyrolysis.nbrCoupledPatchID(film, coupledPatchi);
339 const scalarField htcw(film.htcw().h()().boundaryField()[filmPatchi]);
343 pyrolysis.mapRegionPatchField
355 Tfilm = film.Ts().boundaryField()[filmPatchi];
356 film.toPrimary(filmPatchi, Tfilm);
360 pyrolysis.mapRegionPatchField<scalar>
375 (filmDelta - filmDeltaDry_)/(filmDeltaWet_ - filmDeltaDry_),
382 scalarField qConv(ratio*htcwfilm*(Tfilm - Tp)*convectiveScaling_);
388 valueFraction() =
alpha/(
alpha + (1.0 - ratio)*myKDelta);
390 refValue() = ratio*Tfilm + (1.0 - ratio)*(KDeltaNbr*nbrIntFld)/
alpha;
392 mixedFvPatchScalarField::updateCoeffs();
398 scalar Qt =
gSum((qConv + qRad)*
patch().magSf());
401 <<
patch().name() <<
':'
402 << this->internalField().name() <<
" <- "
403 << nbrMesh.name() <<
':'
404 << nbrPatch.name() <<
':'
405 << this->internalField().name() <<
" :" <<
nl
406 <<
" convective heat[W] : " << Qc <<
nl
407 <<
" radiative heat [W] : " << qr <<
nl
408 <<
" total heat [W] : " << Qt <<
nl
409 <<
" wall temperature "
410 <<
" min:" <<
gMin(*
this)
411 <<
" max:" <<
gMax(*
this)
424 os.writeEntryIfDifferent<word>
427 "surfaceFilmProperties",
430 os.writeEntryIfDifferent<word>
433 "pyrolysisProperties",
436 os.writeEntry(
"Tnbr", TnbrName_);
437 os.writeEntry(
"qr", qrName_);
438 os.writeEntry(
"convectiveScaling", convectiveScaling_);
439 os.writeEntry(
"filmDeltaDry", filmDeltaDry_);
440 os.writeEntry(
"filmDeltaWet", filmDeltaWet_);
450 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField
Base class for pyrolysis models.
fvPatchField< scalar > fvPatchScalarField
virtual void operator=(const UList< Type > &)
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< scalar > scalarList
A List of scalars.
A class for handling words, derived from Foam::string.
static constexpr const zero Zero
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Type gAverage(const FieldField< Field, Type > &f)
Mixed boundary condition for temperature, to be used in the flow and pyrolysis regions when a film re...
virtual void write(Ostream &) const
const polyPatch & samplePolyPatch() const
Common functions used in temperature coupled boundaries.
const Time & time() const
Ostream & endl(Ostream &os)
Type gSum(const FieldField< Field, Type > &f)
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
virtual const word & name() const
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
virtual void autoMap(const fvPatchFieldMapper &)=0
label min(const labelHashSet &set, label minValue=labelMax)
Mesh consisting of general polyhedral cells.
virtual void autoMap(const fvPatchFieldMapper &)
Foam::regionModels::pyrolysisModels::pyrolysisModel pyrolysisModelType
CGAL::Exact_predicates_exact_constructions_kernel K
const fvMesh & primaryMesh() const
void toPrimary(const label regionPatchi, List< Type > ®ionField) const
HashTable< const Type * > lookupClass(const bool strict=false) const
const polyMesh & sampleMesh() const
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
tmp< scalarField > K() const
Thermodynamic form of single-cell layer surface film model.
Generic templated field type.
void distribute(List< Type > &lst) const
virtual tmp< volScalarField > h() const =0
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual tmp< scalarField > alpha(const scalarField &Tp) const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
label max(const labelHashSet &set, label maxValue=labelMin)
Lookup type of boundary radiation properties.
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
label index() const noexcept
errorManip< error > abort(error &err)
errorManipArg< error, int > exit(error &err, const int errNo=1)
filmPyrolysisRadiativeCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
const word & name() const noexcept
#define FatalErrorInFunction
forAllConstIters(mixture.phases(), phase)
const heatTransferModel & htcw() const
tmp< Foam::Field< Type > > mapRegionPatchField(const regionModel &nbrRegion, const label regionPatchi, const label nbrPatchi, const Field< Type > &nbrField, const bool flip=false) const
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
const scalarField & deltaCoeffs() const
fileName::Type type(const fileName &name, const bool followLink=true)
virtual tmp< scalarField > kappa(const scalarField &Tp) const
virtual void updateCoeffs()
Ostream & writeEntry(const keyType &key, const T &value)
void write(vtk::formatter &fmt, const Type &val, const label n=1)
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
word name(const expressions::valueTypeCode typeCode)
Foam::regionModels::surfaceFilmModels::thermoSingleLayer filmModelType
Foam::fvPatchFieldMapper.
void write(Ostream &os) const
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Type gMin(const FieldField< Field, Type > &f)
virtual const volScalarField & Ts() const
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Type gMax(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const word & name() const