Go to the documentation of this file.
45 "fixed_heat_transfer_coefficient",
67 mixedFvPatchScalarField(
p, iF),
73 QrPrevious_(
p.size()),
75 QrName_(
"undefined-Qr"),
81 valueFraction() = 1.0;
94 mixedFvPatchScalarField(ptf,
p, iF, mapper),
116 mixedFvPatchScalarField(
p, iF),
122 QrPrevious_(
p.size(), 0.0),
123 QrRelaxation_(
dict.lookupOrDefault<scalar>(
"relaxation", 1)),
124 QrName_(
dict.lookupOrDefault<
word>(
"Qr",
"none")),
128 if (
dict.found(
"q") && !
dict.found(
"h") && !
dict.found(
"Ta"))
130 mode_ = fixedHeatFlux;
133 else if (
dict.found(
"h") &&
dict.found(
"Ta") && !
dict.found(
"q"))
135 mode_ = fixedHeatTransferCoeff;
138 if (
dict.found(
"thicknessLayers"))
140 if (
dict.readIfPresent(
"thicknessLayers", thicknessLayers_))
142 dict.lookup(
"kappaLayers") >> kappaLayers_;
144 if (thicknessLayers_.size() != kappaLayers_.size())
147 <<
"\n number of layers for thicknessLayers and "
148 <<
"kappaLayers must be the same"
149 <<
"\n for patch " <<
p.name()
159 <<
"\n patch type '" <<
p.type()
160 <<
"' either q or h and Ta were not found '"
161 <<
"\n for patch " <<
p.name()
169 if (
dict.found(
"QrPrevious"))
174 if (
dict.found(
"refValue"))
186 valueFraction() = 1.0;
197 mixedFvPatchScalarField(tppsf),
218 mixedFvPatchScalarField(tppsf, iF),
239 mixedFvPatchScalarField::autoMap(m);
252 mixedFvPatchScalarField::rmap(ptf, addr);
255 refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
257 q_.rmap(tiptf.
q_, addr);
258 h_.rmap(tiptf.
h_, addr);
259 Ta_.rmap(tiptf.
Ta_, addr);
288 valueFraction() = 0.0;
294 scalar totalSolidRes = 0.0;
306 hp = 1.0/(1.0/
h_ + totalSolidRes);
310 refValue() = hp*
Ta_/(hp -
Qr);
312 (hp -
Qr)/((hp -
Qr) +
kappa(Tp)*patch().deltaCoeffs());
324 mixedFvPatchScalarField::updateCoeffs();
330 Info<< patch().boundaryMesh().mesh().name() <<
':'
331 << patch().name() <<
':'
333 <<
" heat transfer rate:" << Q
334 <<
" wall temperature "
335 <<
" min:" <<
gMin(*
this)
336 <<
" max:" <<
gMax(*
this)
351 QrPrevious_.writeEntry(
"QrPrevious", os);
361 q_.writeEntry(
"q", os);
364 case fixedHeatTransferCoeff:
366 h_.writeEntry(
"h", os);
367 Ta_.writeEntry(
"Ta", os);
368 thicknessLayers_.writeEntry(
"thicknessLayers", os);
369 kappaLayers_.writeEntry(
"kappaLayers", os);
375 <<
"Illegal heat flux mode " << operationModeNames[mode_]
scalarList kappaLayers_
Conductivity of layers.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
scalarField q_
Heat flux / [W/m2].
operationMode mode_
Operation mode.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
A class for handling words, derived from string.
#define forAll(list, i)
Loop across all elements in list.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
scalarField h_
Heat transfer coefficient / [W/m2K].
Type gAverage(const FieldField< Field, Type > &f)
rDeltaT dimensionedInternalField()
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Common functions for use in temperature coupled boundaries.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Type gSum(const FieldField< Field, Type > &f)
scalarField Ta_
Ambient temperature / [K].
This function object calculates and outputs the second invariant of the velocity gradient tensor [1/s...
Pre-declare SubField and related Field type.
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))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
operationMode
Operation mode enumeration.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Macros for easy insertion into run-time selection tables.
errorManip< error > abort(error &err)
scalarList thicknessLayers_
Thickness of layers.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
errorManipArg< error, int > exit(error &err, const int errNo=1)
externalWallHeatFluxTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
void write(Ostream &) const
Write.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
scalarField QrPrevious_
Chache Qr for relaxation.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
volScalarField scalarField(fieldObject, mesh)
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.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Foam::fvPatchFieldMapper.
void size(const label)
Override size to be inconsistent with allocated storage.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Type gMin(const FieldField< Field, Type > &f)
This boundary condition supplies a heat flux condition for temperature on an external wall....
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void write(Ostream &) const
Write.
static const NamedEnum< operationMode, 3 > operationModeNames
scalar QrRelaxation_
Relaxation for Qr.
const word QrName_
Name of the radiative heat flux.
Type gMax(const FieldField< Field, Type > &f)
Initialise the NamedEnum HashTable from the static list of names.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...