36 namespace compressible
41 template<
class sol
idType>
50 mixedFvPatchScalarField(
p, iF,
dict, false),
52 baffleActivated_(
dict.lookupOrDefault<bool>(
"baffleActivated", true)),
57 qrPrevious_(
p.size(), 0.0),
58 qrRelaxation_(
dict.lookupOrDefault<scalar>(
"qrRelaxation",
dimless, 1)),
59 qrName_(
dict.lookupOrDefault<
word>(
"qr",
"none"))
75 if (
dict.found(
"thickness"))
85 if (
dict.found(
"qrPrevious"))
91 if (
dict.found(
"refValue") && baffleActivated_)
111 valueFraction() = 0.0;
116 template<
class sol
idType>
126 mixedFvPatchScalarField(ptf,
p, iF, mapper),
128 baffleActivated_(ptf.baffleActivated_),
129 thickness_(mapper(ptf.thickness_)),
130 qs_(mapper(ptf.qs_)),
131 solidDict_(ptf.solidDict_),
132 solidPtr_(ptf.solidPtr_),
133 qrPrevious_(mapper(ptf.qrPrevious_)),
134 qrRelaxation_(ptf.qrRelaxation_),
139 template<
class sol
idType>
147 mixedFvPatchScalarField(ptf, iF),
149 baffleActivated_(ptf.baffleActivated_),
150 thickness_(ptf.thickness_),
152 solidDict_(ptf.solidDict_),
153 solidPtr_(ptf.solidPtr_),
154 qrPrevious_(ptf.qrPrevious_),
155 qrRelaxation_(ptf.qrRelaxation_),
162 template<
class sol
idType>
171 template<
class sol
idType>
172 const thermalBaffle1DFvPatchScalarField<solidType>&
173 thermalBaffle1DFvPatchScalarField<solidType>::nbrField()
const
180 refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchi];
183 refCast<const thermalBaffle1DFvPatchScalarField>
185 nbrPatch.template lookupPatchField<volScalarField, scalar>
193 template<
class sol
idType>
194 const solidType& thermalBaffle1DFvPatchScalarField<solidType>::solid()
const
198 if (solidPtr_.empty())
200 solidPtr_.reset(
new solidType(
"solid", solidDict_));
207 return nbrField().solid();
212 template<
class sol
idType>
213 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::
214 baffleThickness()
const
218 if (thickness_.size() != patch().size())
221 <<
" Field thickness has not been specified "
222 <<
" for patch " << this->patch().name()
230 const mappedFvPatchBaseBase& mapper =
232 return mapper.fromNeighbour(nbrField().baffleThickness());
237 template<
class sol
idType>
238 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::qs()
const
246 const mappedFvPatchBaseBase& mapper =
248 return mapper.fromNeighbour(nbrField().qs());
253 template<
class sol
idType>
260 mixedFvPatchScalarField::map(ptf, mapper);
263 refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
267 mapper(thickness_, tiptf.thickness_);
268 mapper(qs_, tiptf.qs_);
273 template<
class sol
idType>
279 mixedFvPatchScalarField::reset(ptf);
282 refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
286 thickness_.reset(tiptf.thickness_);
287 qs_.reset(tiptf.qs_);
292 template<
class sol
idType>
308 if (baffleActivated_)
311 db().objectRegistry::template lookupObject
318 thermophysicalTransportModel::typeName,
319 internalField().
group()
325 patch().template lookupPatchField<volScalarField, scalar>(TName_);
331 if (qrName_ !=
"none")
333 qr = patch().template lookupPatchField<volScalarField, scalar>
336 qr = qrRelaxation_*qr + (1.0 - qrRelaxation_)*qrPrevious_;
340 const scalarField kappaDelta(kappap*patch().deltaCoeffs());
343 const scalarField nbrTp(mapper.fromNeighbour(nbrField()));
349 kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0);
352 const scalarField kappaByDeltaSolid(kappas/baffleThickness());
358 refValue() = (kappaByDeltaSolid*nbrTp + qs()/2.0)/
alpha;
363 Info<< patch().boundaryMesh().mesh().name() <<
':'
364 << patch().name() <<
':'
365 << this->internalField().name() <<
" <- "
366 << nbrField().patch().name() <<
':'
367 << this->internalField().name() <<
" :"
369 <<
" walltemperature "
370 <<
" min:" <<
gMin(*
this)
371 <<
" max:" <<
gMax(*
this)
380 mixedFvPatchScalarField::updateCoeffs();
384 template<
class sol
idType>
391 writeEntry(os,
"thickness", baffleThickness()());
398 writeEntry(os,
"qrRelaxation", qrRelaxation_);
#define forAll(list, i)
Loop across all elements in list.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
static word groupName(Name name, const word &group)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
static int & msgType()
Message tag of standard messages.
This BC solves a steady 1D thermal baffle.
virtual void write(Ostream &) const
Write.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
thermalBaffle1DFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Abstract base class for field mapping.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
label index() const
Return the index of this patch in the fvBoundaryMesh.
Base class for fv patches that provide mapping between two fv patches.
const fvPatch & nbrFvPatch() const
Get the patch to map from.
static const mappedFvPatchBaseBase & getMap(const fvPatch &patch)
Cast the given fvPatch to a mappedFvPatchBaseBase. Handle errors.
const fvMesh & nbrMesh() const
Get the mesh for the region to map from.
static void validateMapForField(const PatchField &field, const FieldType &iF, const dictionary &context, const label froms=from::any)
Validate that the map exists and is appropriate for the given.
Mesh consisting of general polyhedral cells.
Abstract base class for all fluid and solid thermophysical transport models.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
A class for handling words, derived from string.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const char *const group
Group name for atomic constants.
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
errorManipArg< error, int > exit(error &err, const int errNo=1)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const dimensionSet dimPower
const dimensionSet dimless
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Type gAverage(const FieldField< Field, Type > &f)
const dimensionSet dimArea
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
const unitConversion unitFraction
static const label differentPatch
static const label sameRegion