Go to the documentation of this file.
48 mixedFvPatchScalarField(
p, iF),
50 TnbrName_(
"undefined-Tnbr"),
51 QrNbrName_(
"undefined-QrNbr"),
52 QrName_(
"undefined-Qr"),
57 this->refValue() = 0.0;
58 this->refGrad() = 0.0;
59 this->valueFraction() = 1.0;
72 mixedFvPatchScalarField(psf,
p, iF, mapper),
91 mixedFvPatchScalarField(
p, iF),
93 TnbrName_(
dict.lookupOrDefault<
word>(
"Tnbr",
"T")),
94 QrNbrName_(
dict.lookupOrDefault<
word>(
"QrNbr",
"none")),
95 QrName_(
dict.lookupOrDefault<
word>(
"Qr",
"none")),
100 if (!isA<mappedPatchBase>(this->patch().patch()))
103 <<
"' not type '" << mappedPatchBase::typeName <<
"'"
104 <<
"\n for patch " <<
p.name()
110 if (
dict.found(
"thicknessLayers"))
112 dict.lookup(
"thicknessLayers") >> thicknessLayers_;
113 dict.lookup(
"kappaLayers") >> kappaLayers_;
115 if (thicknessLayers_.size() > 0)
118 forAll (thicknessLayers_, iLayer)
120 contactRes_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
122 contactRes_ = 1.0/contactRes_;
128 if (
dict.found(
"refValue"))
140 valueFraction() = 1.0;
152 mixedFvPatchScalarField(psf, iF),
179 refCast<const mappedPatchBase>(patch().patch());
183 refCast<const fvMesh>(nbrMesh).boundary()[samplePatchI];
228 valueFraction() = KDeltaNbr/(KDeltaNbr + KDelta);
230 refGrad() = (
Qr + QrNbr)/
kappa(Tp);
232 mixedFvPatchScalarField::updateCoeffs();
238 Info<< patch().boundaryMesh().mesh().name() <<
':'
239 << patch().name() <<
':'
241 << nbrMesh.
name() <<
':'
242 << nbrPatch.
name() <<
':'
244 <<
" heat transfer rate:" << Q
245 <<
" walltemperature "
246 <<
" min:" <<
gMin(Tp)
247 <<
" max:" <<
gMax(Tp)
266 thicknessLayers_.writeEntry(
"thicknessLayers", os);
267 kappaLayers_.writeEntry(
"kappaLayers", os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void operator=(const UList< Type > &)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
A class for handling words, derived from string.
scalarList kappaLayers_
Conductivity of layers.
#define forAll(list, i)
Loop across all elements in list.
Type gAverage(const FieldField< Field, Type > &f)
rDeltaT dimensionedInternalField()
const polyPatch & samplePolyPatch() const
Get the patch on the region.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Common functions for use in temperature coupled boundaries.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Type gSum(const FieldField< Field, Type > &f)
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Mesh consisting of general polyhedral cells.
const word & name() const
Return name.
This function object calculates and outputs the second invariant of the velocity gradient tensor [1/s...
scalar contactRes_
Total contact resistance.
const polyMesh & sampleMesh() const
Get the region mesh.
scalarList thicknessLayers_
Thickness of layers.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
volScalarField Qr(IOobject("Qr", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("Qr", dimMass/pow3(dimTime), 0.0))
const word TnbrName_
Name of field on the neighbour region.
const word & name() const
Return name.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word QrName_
Name of the radiative heat flux in local region.
const word QrNbrName_
Name of the radiative heat flux in the neighbout region.
turbulentTemperatureRadCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by any number of values (e....
virtual void write(Ostream &) const
Write.
Macros for easy insertion into run-time selection tables.
To & refCast(From &r)
Reference type cast template function.
errorManipArg< error, int > exit(error &err, const int errNo=1)
static int & msgType()
Message tag of standard messages.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Mixed boundary condition for temperature and radiation heat transfer to be used for in multiregion ca...
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
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.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Type gMin(const FieldField< Field, Type > &f)
void write(Ostream &) const
Write.
label index() const
Return the index of this patch in the boundaryMesh.
Type gMax(const FieldField< Field, Type > &f)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...