45 "fixed_heat_transfer_coefficient",
67 mixedFvPatchScalarField(p, iF),
73 QrPrevious_(p.
size(), 0.0),
75 QrName_(
"undefined-Qr"),
81 valueFraction() = 1.0;
94 mixedFvPatchScalarField(ptf, p, iF, mapper),
100 QrPrevious_(ptf.QrPrevious_, mapper),
101 QrRelaxation_(ptf.QrRelaxation_),
102 QrName_(ptf.QrName_),
103 thicknessLayers_(ptf.thicknessLayers_),
104 kappaLayers_(ptf.kappaLayers_)
116 mixedFvPatchScalarField(p, iF),
122 QrPrevious_(p.
size(), 0.0),
130 mode_ = fixedHeatFlux;
135 mode_ = fixedHeatTransferCoeff;
138 if (dict.
found(
"thicknessLayers"))
140 dict.
lookup(
"thicknessLayers") >> thicknessLayers_;
141 dict.
lookup(
"kappaLayers") >> kappaLayers_;
147 <<
"\n patch type '" << p.type()
148 <<
"' either q or h and Ta were not found '" 149 <<
"\n for patch " << p.
name()
150 <<
" of field " << internalField().name()
151 <<
" in file " << internalField().objectPath()
157 if (dict.
found(
"QrPrevious"))
162 if (dict.
found(
"refValue"))
174 valueFraction() = 1.0;
185 mixedFvPatchScalarField(tppsf),
191 QrPrevious_(tppsf.QrPrevious_),
192 QrRelaxation_(tppsf.QrRelaxation_),
193 QrName_(tppsf.QrName_),
194 thicknessLayers_(tppsf.thicknessLayers_),
195 kappaLayers_(tppsf.kappaLayers_)
206 mixedFvPatchScalarField(tppsf, iF),
212 QrPrevious_(tppsf.QrPrevious_),
213 QrRelaxation_(tppsf.QrRelaxation_),
214 QrName_(tppsf.QrName_),
215 thicknessLayers_(tppsf.thicknessLayers_),
216 kappaLayers_(tppsf.kappaLayers_)
227 mixedFvPatchScalarField::autoMap(m);
231 QrPrevious_.autoMap(m);
241 mixedFvPatchScalarField::rmap(ptf, addr);
244 refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
246 q_.
rmap(tiptf.q_, addr);
247 h_.rmap(tiptf.h_, addr);
248 Ta_.rmap(tiptf.Ta_, addr);
249 QrPrevious_.rmap(tiptf.QrPrevious_, addr);
264 if (QrName_ !=
"none")
268 Qr = QrRelaxation_*Qr + (1.0 - QrRelaxation_)*QrPrevious_;
276 refGrad() = (q_ +
Qr)/
kappa(Tp);
278 valueFraction() = 0.0;
282 case fixedHeatTransferCoeff:
284 scalar totalSolidRes = 0.0;
285 if (thicknessLayers_.size() > 0)
287 forAll(thicknessLayers_, iLayer)
289 const scalar l = thicknessLayers_[iLayer];
290 if (kappaLayers_[iLayer] > 0.0)
292 totalSolidRes += l/kappaLayers_[iLayer];
296 hp = 1.0/(1.0/h_ + totalSolidRes);
300 refValue() = hp*Ta_/(hp -
Qr);
302 (hp -
Qr)/((hp - Qr) +
kappa(Tp)*patch().deltaCoeffs());
309 <<
"Illegal heat flux mode " << operationModeNames[mode_]
314 mixedFvPatchScalarField::updateCoeffs();
320 Info<< patch().boundaryMesh().mesh().name() <<
':' 321 << patch().name() <<
':' 322 << this->internalField().name() <<
" :" 323 <<
" heat transfer rate:" << Q
324 <<
" walltemperature " 325 <<
" min:" <<
gMin(*
this)
326 <<
" max:" <<
gMax(*
this)
341 QrPrevious_.writeEntry(
"QrPrevious", os);
351 q_.writeEntry(
"q", os);
354 case fixedHeatTransferCoeff:
356 h_.writeEntry(
"h", os);
357 Ta_.writeEntry(
"Ta", os);
358 thicknessLayers_.writeEntry(
"thicknessLayers", os);
359 kappaLayers_.writeEntry(
"kappaLayers", os);
365 <<
"Illegal heat flux mode " << operationModeNames[mode_]
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
#define forAll(list, i)
Loop across all elements in list.
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...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Type gMin(const FieldField< Field, Type > &f)
void size(const label)
Override size to be inconsistent with allocated storage.
const word & name() const
Return name.
Ostream & endl(Ostream &os)
Add newline and flush stream.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Initialise the NamedEnum HashTable from the static list of names.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from string.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual label size() const
Return size.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
errorManip< error > abort(error &err)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
void write(Ostream &) const
Write.
Type gMax(const FieldField< Field, Type > &f)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
static const NamedEnum< operationMode, 3 > operationModeNames
Common functions used in temperature coupled boundaries.
operationMode
Operation mode enumeration.
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 void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
externalWallHeatFluxTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
volScalarField Qr(IOobject("Qr", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("Qr", dimMass/pow3(dimTime), 0.0))
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
This boundary condition supplies a heat flux condition for temperature on an external wall...
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.