Go to the documentation of this file.
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),
87 template<
class sol
idType>
97 mixedFvPatchScalarField(
p, iF),
99 baffleActivated_(
dict.lookupOrDefault<
bool>(
"baffleActivated",
true)),
104 QrPrevious_(
p.size(), 0.0),
105 QrRelaxation_(
dict.lookupOrDefault<scalar>(
"relaxation", 1)),
106 QrName_(
dict.lookupOrDefault<
word>(
"Qr",
"none"))
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),
164 template<
class sol
idType>
173 mixedFvPatchScalarField(ptf, iF),
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);
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);
399 refValue() = (KDeltaSolid*nbrTp + Qs()/2.0)/
alpha;
404 Info<< patch().boundaryMesh().mesh().name() <<
':'
405 << patch().name() <<
':'
407 << nbrPatch.
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);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void operator=(const UList< Type > &)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual void write(Ostream &) const
Write.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void write(Ostream &) const
Write as a dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from string.
word TName_
Name of the temperature field.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const mapDistribute & map() const
Return reference to the parallel distribution map.
Type gAverage(const FieldField< Field, Type > &f)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
rDeltaT dimensionedInternalField()
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
dictionary solidDict_
Solid dictionary.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
scalar QrRelaxation_
Relaxation for Qr.
const word & name() const
Return name.
This function object calculates and outputs the second invariant of the velocity gradient tensor [1/s...
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
scalarField QrPrevious_
Cache Qr for relaxation.
Class containing processor-to-processor mapping information.
tmp< scalarField > Qs() const
Return Qs from master.
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))
bool baffleActivated_
Baffle is activated.
autoPtr< solidType > solidPtr_
Solid thermo.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
bool owner() const
Is Owner.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const solidType & solid() const
Return const solid thermo.
tmp< scalarField > baffleThickness() const
Return thickness from master.
errorManipArg< error, int > exit(error &err, const int errNo=1)
virtual tmp< volScalarField > kappaEff() const
Return the effective turbulent thermal diffusivity for temperature.
static int & msgType()
Message tag of standard messages.
const word QrName_
Name of the radiative heat flux in local region.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
thermalBaffle1DFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
scalarField Qs_
Superficial heat source [W/m2].
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
scalarField thickness_
Baffle thickness [m].
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.