38 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
41 HashTable<const filmModelType*> models
46 if (iter()->regionMesh().
name() == filmRegionName_)
52 DynamicList<word> modelNames;
55 modelNames.append(iter()->regionMesh().
name());
59 <<
"Unable to locate film region " << filmRegionName_
60 <<
". Available regions include: " << modelNames
63 return **models.begin();
69 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
72 HashTable<const pyrolysisModelType*> models =
77 if (iter()->regionMesh().
name() == pyrolysisRegionName_)
83 DynamicList<word> modelNames;
86 modelNames.append(iter()->regionMesh().
name());
91 <<
"Unable to locate pyrolysis region " << pyrolysisRegionName_
92 <<
". Available regions include: " << modelNames
95 return **models.begin();
108 mixedFvPatchScalarField(p, iF),
110 filmRegionName_(
"surfaceFilmProperties"),
111 pyrolysisRegionName_(
"pyrolysisProperties"),
112 TnbrName_(
"undefined-Tnbr"),
113 QrName_(
"undefined-Qr"),
114 convectiveScaling_(1.0),
118 this->refValue() = 0.0;
119 this->refGrad() = 0.0;
120 this->valueFraction() = 1.0;
133 mixedFvPatchScalarField(psf, p, iF, mapper),
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_)
153 mixedFvPatchScalarField(p, iF),
163 TnbrName_(dict.
lookup(
"Tnbr")),
164 QrName_(dict.
lookup(
"Qr")),
165 convectiveScaling_(dict.
lookupOrDefault<scalar>(
"convectiveScaling", 1.0)),
169 if (!isA<mappedPatchBase>(this->patch().patch()))
172 <<
"' not type '" << mappedPatchBase::typeName <<
"'" 173 <<
"\n for patch " << p.
name()
174 <<
" of field " << internalField().name()
175 <<
" in file " << internalField().objectPath()
181 if (dict.
found(
"refValue"))
193 valueFraction() = 1.0;
205 mixedFvPatchScalarField(psf, iF),
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_)
228 refCast<const mappedPatchBase>(patch().patch());
235 refCast<const fvMesh>(nbrMesh).
boundary()[nbrPatchi];
250 scalarField nbrIntFld(nbrField.patchInternalField());
274 label coupledPatchi = -1;
275 if (pyrolysisRegionName_ == mesh.
name())
278 if (QrName_ !=
"none")
286 coupledPatchi = nbrPatch.
index();
287 if (QrName_ !=
"none")
295 <<
type() <<
" condition is intended to be applied to either the " 296 <<
"primary or pyrolysis regions only" 338 (filmDelta - filmDeltaDry_)/(filmDeltaWet_ - filmDeltaDry_),
345 scalarField qConv(ratio*htcwfilm*(Tfilm - Tp)*convectiveScaling_);
351 valueFraction() = alpha/(alpha + (1.0 - ratio)*myKDelta);
353 refValue() = ratio*Tfilm + (1.0 - ratio)*(KDeltaNbr*nbrIntFld)/
alpha;
355 mixedFvPatchScalarField::updateCoeffs();
359 scalar Qc =
gSum(qConv*patch().magSf());
360 scalar Qr =
gSum(qRad*patch().magSf());
361 scalar Qt =
gSum((qConv + qRad)*patch().magSf());
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)
387 writeEntryIfDifferent<word>
391 "surfaceFilmProperties",
394 writeEntryIfDifferent<word>
398 "pyrolysisProperties",
403 os.writeKeyword(
"convectiveScaling") << convectiveScaling_
405 os.writeKeyword(
"filmDeltaDry") << filmDeltaDry_
407 os.writeKeyword(
"filmDeltaWet") << filmDeltaWet_
Foam::regionModels::surfaceFilmModels::thermoSingleLayer filmModelType
Mixed boundary condition for temperature, to be used in the flow and pyrolysis regions when a film re...
Base class for pyrolysis models.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const Time & time() const
Return the reference to the time database.
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void write(Ostream &) const
Write.
A list of keyword definitions, which are a keyword followed by any number of values (e...
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Type gMin(const FieldField< Field, Type > &f)
To & refCast(From &r)
Reference type cast template function.
virtual tmp< volScalarField > h() const =0
Return the heat transfer coefficient [W/m2/K].
const word & name() const
Return name.
Ostream & endl(Ostream &os)
Add newline and flush stream.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const polyMesh & sampleMesh() const
Get the region mesh.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
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
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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.
Type gSum(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A class for handling words, derived from string.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
virtual label size() const
Return size.
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
Return the coupled patch ID paired with coupled patch.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
faceListList boundary(nPatches)
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=NULL, const Type *=NULL) const
Lookup and return the patchField of the named field from the.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
errorManip< error > abort(error &err)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Type gMax(const FieldField< Field, Type > &f)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
label index() const
Return the index of this patch in the fvBoundaryMesh.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
tmp< scalarField > K() const
Get corresponding K field.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
void toPrimary(const label regionPatchi, List< Type > ®ionField) const
Convert a local region field to the primary region.
Type gAverage(const FieldField< Field, Type > &f)
virtual void operator=(const UList< Type > &)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
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.
volScalarField Qr(IOobject("Qr", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("Qr", dimMass/pow3(dimTime), 0.0))
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Foam::regionModels::pyrolysisModels::pyrolysisModel pyrolysisModelType
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
label index() const
Return the index of this patch in the boundaryMesh.
const word & name() const
Return reference to name.
const word & name() const
Return name.
virtual void write(Ostream &) const
Write.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.