36 namespace compressible
48 mixedFvPatchScalarField(p, iF),
50 TnbrName_(
"undefined-Tnbr"),
55 this->refValue() = 0.0;
56 this->refGrad() = 0.0;
57 this->valueFraction() = 1.0;
69 mixedFvPatchScalarField(p, iF),
71 TnbrName_(dict.
lookup(
"Tnbr")),
76 if (!isA<mappedPatchBase>(this->patch().patch()))
79 <<
"' not type '" << mappedPatchBase::typeName <<
"'" 80 <<
"\n for patch " << p.
name()
81 <<
" of field " << internalField().name()
82 <<
" in file " << internalField().objectPath()
86 if (dict.
found(
"thicknessLayers"))
88 dict.
lookup(
"thicknessLayers") >> thicknessLayers_;
89 dict.
lookup(
"kappaLayers") >> kappaLayers_;
91 if (thicknessLayers_.size() > 0)
94 forAll(thicknessLayers_, iLayer)
96 contactRes_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
98 contactRes_ = 1.0/contactRes_;
104 if (dict.
found(
"refValue"))
116 valueFraction() = 1.0;
130 mixedFvPatchScalarField(ptf, p, iF, mapper),
132 TnbrName_(ptf.TnbrName_),
133 thicknessLayers_(ptf.thicknessLayers_),
134 kappaLayers_(ptf.kappaLayers_),
135 contactRes_(ptf.contactRes_)
146 mixedFvPatchScalarField(wtcsf, iF),
148 TnbrName_(wtcsf.TnbrName_),
149 thicknessLayers_(wtcsf.thicknessLayers_),
150 kappaLayers_(wtcsf.kappaLayers_),
151 contactRes_(wtcsf.contactRes_)
171 refCast<const mappedPatchBase>(patch().patch());
175 refCast<const fvMesh>(nbrMesh).
boundary()[samplePatchi];
185 if (!isA<thisType>(nbrTp))
188 <<
"Patch field for " << internalField().name() <<
" on " 189 << patch().name() <<
" is of type " << thisType::typeName
190 <<
endl <<
"The neighbouring patch field " << TnbrName_ <<
" on " 191 << nbrPatch.
name() <<
" is required to be the same, but is " 195 const thisType& nbrField = refCast<const thisType>(nbrTp);
201 if (contactRes_ == 0.0)
203 nbrIntFld.
ref() = nbrField.patchInternalField();
204 nbrKDelta.ref() = nbrField.kappa(nbrField)*nbrPatch.
deltaCoeffs();
208 nbrIntFld.ref() = nbrField;
209 nbrKDelta.ref() = contactRes_;
233 this->refValue() = nbrIntFld();
234 this->refGrad() = 0.0;
235 this->valueFraction() = nbrKDelta()/(nbrKDelta() + myKDelta());
237 mixedFvPatchScalarField::updateCoeffs();
243 Info<< patch().boundaryMesh().mesh().name() <<
':' 244 << patch().name() <<
':' 245 << this->internalField().name() <<
" <- " 246 << nbrMesh.name() <<
':' 247 << nbrPatch.
name() <<
':' 248 << this->internalField().name() <<
" :" 249 <<
" heat transfer rate:" << Q
250 <<
" walltemperature " 251 <<
" min:" <<
gMin(*
this)
252 <<
" max:" <<
gMax(*
this)
269 writeEntry(os,
"thicknessLayers", thicknessLayers_);
Mixed boundary condition for temperature, to be used for heat-transfer on back-to-back baffles...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
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)
const polyMesh & sampleMesh() const
Get the region mesh.
T & ref() const
Return non-const reference or generate a fatal error.
virtual void write(Ostream &) const
Write.
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...
static int & msgType()
Message tag of standard messages.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Type gSum(const FieldField< Field, Type > &f)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
faceListList boundary(nPatches)
virtual label size() const
Return size.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Type gMax(const FieldField< Field, Type > &f)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Common functions used in temperature coupled boundaries.
tmp< scalarField > kappa(const fvPatchScalarField &Tp) const
Given patch temperature calculate corresponding K field.
Type gAverage(const FieldField< Field, Type > &f)
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
label index() const
Return the index of this patch in the boundaryMesh.
virtual void operator=(const UList< Type > &)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Mesh consisting of general polyhedral cells.
virtual const word & name() const
Return name.
A class for managing temporary objects.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
void write(Ostream &) const
Write.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.