41 template<
class sol
idType>
50 mixedFvPatchScalarField(p, iF),
52 baffleActivated_(
true),
57 QrPrevious_(p.
size()),
59 QrName_(
"undefined-Qr")
63 template<
class sol
idType>
74 mixedFvPatchScalarField(ptf, p, iF, mapper),
76 baffleActivated_(ptf.baffleActivated_),
77 thickness_(ptf.thickness_, mapper),
79 solidDict_(ptf.solidDict_),
80 solidPtr_(ptf.solidPtr_),
81 QrPrevious_(ptf.QrPrevious_, mapper),
82 QrRelaxation_(ptf.QrRelaxation_),
87 template<
class sol
idType>
97 mixedFvPatchScalarField(p, iF),
104 QrPrevious_(p.
size(), 0.0),
110 if (dict.
found(
"thickness"))
115 if (dict.
found(
"Qs"))
120 if (dict.
found(
"QrPrevious"))
125 if (dict.
found(
"refValue") && baffleActivated_)
137 valueFraction() = 0.0;
143 template<
class sol
idType>
151 mixedFvPatchScalarField(ptf),
153 baffleActivated_(ptf.baffleActivated_),
154 thickness_(ptf.thickness_),
156 solidDict_(ptf.solidDict_),
157 solidPtr_(ptf.solidPtr_),
158 QrPrevious_(ptf.QrPrevious_),
159 QrRelaxation_(ptf.QrRelaxation_),
164 template<
class sol
idType>
173 mixedFvPatchScalarField(ptf, iF),
175 baffleActivated_(ptf.baffleActivated_),
176 thickness_(ptf.thickness_),
178 solidDict_(ptf.solidDict_),
179 solidPtr_(ptf.solidPtr_),
180 QrPrevious_(ptf.QrPrevious_),
181 QrRelaxation_(ptf.QrRelaxation_),
188 template<
class sol
idType>
193 const label nbrPatchi = samplePolyPatch().index();
195 return (patchi < nbrPatchi);
199 template<
class sol
idType>
204 if (solidPtr_.empty())
206 solidPtr_.reset(
new solidType(solidDict_));
216 refCast<const thermalBaffle1DFvPatchScalarField>
218 nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
221 return nbrField.solid();
226 template<
class sol
idType>
232 if (thickness_.size() != patch().size())
237 )<<
" Field thickness has not been specified " 238 <<
" for patch " << this->patch().name()
251 refCast<const thermalBaffle1DFvPatchScalarField>
253 nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
267 template<
class sol
idType>
282 refCast<const thermalBaffle1DFvPatchScalarField>
284 nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
295 template<
class sol
idType>
301 mixedFvPatchScalarField::autoMap(m);
305 thickness_.autoMap(m);
311 template<
class sol
idType>
318 mixedFvPatchScalarField::rmap(ptf, addr);
321 refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
325 thickness_.
rmap(tiptf.thickness_, addr);
326 Qs_.rmap(tiptf.Qs_, addr);
331 template<
class sol
idType>
347 const label nbrPatchi = samplePolyPatch().index();
349 if (baffleActivated_)
354 db().template lookupObject<compressible::turbulenceModel>
363 patch().template lookupPatchField<volScalarField, scalar>(TName_);
368 if (QrName_ !=
"none")
370 Qr = patch().template lookupPatchField<volScalarField, scalar>
373 Qr = QrRelaxation_*Qr + (1.0 - QrRelaxation_)*QrPrevious_;
379 scalarField myKDelta(patch().deltaCoeffs()*kappaw);
383 turbModel.transport().T().boundaryField()[nbrPatchi];
390 kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0);
397 valueFraction() = alpha/(alpha + myKDelta);
399 refValue() = (KDeltaSolid*nbrTp + Qs()/2.0)/
alpha;
404 Info<< patch().boundaryMesh().mesh().name() <<
':' 405 << patch().name() <<
':' 406 << this->internalField().name() <<
" <- " 407 << nbrPatch.
name() <<
':' 408 << this->internalField().name() <<
" :" 410 <<
" walltemperature " 411 <<
" min:" <<
gMin(*
this)
412 <<
" max:" <<
gMax(*
this)
421 mixedFvPatchScalarField::updateCoeffs();
424 template<
class sol
idType>
432 baffleThickness()().writeEntry(
"thickness", os);
433 Qs()().writeEntry(
"Qs", os);
437 QrPrevious_.writeEntry(
"QrPrevious", os);
virtual void write(Ostream &) const
Write as a dictionary.
thermalBaffle1DFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
#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.
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...
Type gMin(const FieldField< Field, Type > &f)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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...
static int & msgType()
Message tag of standard messages.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
const polyPatch & patch() const
Return the polyPatch.
const mapDistribute & map() const
Return reference to the parallel distribution map.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
virtual void write(Ostream &) const
Write.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Type gMax(const FieldField< Field, Type > &f)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Class containing processor-to-processor mapping information.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
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...
This BC solves a steady 1D thermal baffle.
virtual tmp< volScalarField > kappaEff() const
Return the effective turbulent thermal diffusivity for temperature.
A class for managing temporary objects.
volScalarField Qr(IOobject("Qr", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("Qr", dimMass/pow3(dimTime), 0.0))
T & ref() const
Return non-const reference or generate a fatal error.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.