46turbulentTemperatureRadCoupledMixedFvPatchScalarField::
47turbulentTemperatureRadCoupledMixedFvPatchScalarField
53 mixedFvPatchScalarField(
p, iF),
67 TnbrName_(
"undefined-Tnbr"),
68 qrNbrName_(
"undefined-qrNbr"),
69 qrName_(
"undefined-qr"),
70 thermalInertia_(false)
72 this->refValue() = 0.0;
73 this->refGrad() = 0.0;
74 this->valueFraction() = 1.0;
79turbulentTemperatureRadCoupledMixedFvPatchScalarField::
80turbulentTemperatureRadCoupledMixedFvPatchScalarField
88 mixedFvPatchScalarField(psf,
p, iF, mapper),
96 TnbrName_(psf.TnbrName_),
97 qrNbrName_(psf.qrNbrName_),
99 thicknessLayers_(psf.thicknessLayers_),
100 thicknessLayer_(psf.thicknessLayer_.clone(
p.patch())),
101 kappaLayers_(psf.kappaLayers_),
102 kappaLayer_(psf.kappaLayer_.clone(
p.patch())),
103 thermalInertia_(psf.thermalInertia_)
107turbulentTemperatureRadCoupledMixedFvPatchScalarField::
108turbulentTemperatureRadCoupledMixedFvPatchScalarField
115 mixedFvPatchScalarField(
p, iF),
123 TnbrName_(
dict.getOrDefault<
word>(
"Tnbr",
"T")),
124 qrNbrName_(
dict.getOrDefault<
word>(
"qrNbr",
"none")),
125 qrName_(
dict.getOrDefault<
word>(
"qr",
"none")),
126 thermalInertia_(
dict.getOrDefault<
Switch>(
"thermalInertia", false))
128 if (!isA<mappedPatchBase>(this->patch().patch()))
132 <<
"\n for patch " <<
p.
name()
133 <<
" of field " << internalField().name()
134 <<
" in file " << internalField().objectPath()
215 valueFraction() = 1.0;
221 this->useImplicit(boolVal);
234turbulentTemperatureRadCoupledMixedFvPatchScalarField::
235turbulentTemperatureRadCoupledMixedFvPatchScalarField
241 mixedFvPatchScalarField(psf, iF),
249 TnbrName_(psf.TnbrName_),
250 qrNbrName_(psf.qrNbrName_),
251 qrName_(psf.qrName_),
252 thicknessLayers_(psf.thicknessLayers_),
253 thicknessLayer_(psf.thicknessLayer_.clone(patch().patch())),
254 kappaLayers_(psf.kappaLayers_),
255 kappaLayer_(psf.kappaLayer_.clone(patch().patch())),
256 thermalInertia_(psf.thermalInertia_)
260turbulentTemperatureRadCoupledMixedFvPatchScalarField::
261turbulentTemperatureRadCoupledMixedFvPatchScalarField
266 mixedFvPatchScalarField(psf),
274 TnbrName_(psf.TnbrName_),
275 qrNbrName_(psf.qrNbrName_),
276 qrName_(psf.qrName_),
277 thicknessLayers_(psf.thicknessLayers_),
278 thicknessLayer_(psf.thicknessLayer_.clone(patch().patch())),
279 kappaLayers_(psf.kappaLayers_),
280 kappaLayer_(psf.kappaLayer_.clone(patch().patch())),
281 thermalInertia_(psf.thermalInertia_)
292 mixedFvPatchScalarField::autoMap(
mapper);
297 thicknessLayer_().autoMap(
mapper);
298 kappaLayer_().autoMap(
mapper);
309 mixedFvPatchScalarField::rmap(ptf, addr);
321 thicknessLayer_().rmap(tiptf.thicknessLayer_(), addr);
322 kappaLayer_().rmap(tiptf.kappaLayer_(), addr);
355 const label patchi = patch().index();
360 this->internalField()
367 const scalarField KDelta(kappaTp*patch().deltaCoeffs());
378 refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
387 TcNbr = nbrField.patchInternalField();
388 KDeltaNbr = nbrField.kappa(nbrField)*nbrPatch.
deltaCoeffs();
394 TcNbr = patchInternalField();
401 if (thicknessLayer_ || thicknessLayers_.
size())
409 const scalar t = db().time().timeOutputValue();
411 thicknessLayer_().value(t)
412 /kappaLayer_().value(t);
415 if (thicknessLayers_.
size())
417 forAll(thicknessLayers_, iLayer)
419 KDeltaC += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
422 KDeltaC = 1.0/(KDeltaC + SMALL);
429 if (qrName_ !=
"none")
435 if (qrNbrName_ !=
"none")
442 refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
458 <<
"thermalInertia not supported in combination with multi-world"
476 refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
478 thermo->p().boundaryField()[samplePatchi];
480 thermo->T().boundaryField()[samplePatchi];
484 thermo->Cp(ppn, Tpn, samplePatchi)
485 *
thermo->rho(samplePatchi)
510 thermo->Cp(pp, Tp, patchi)
512 / patch().deltaCoeffs()/dt
525 this->internalField().name()
533 scalarField c(KDeltaNbr*TcNbr + (mCpDt + mCpDtNbr)*TpOld);
534 refValue() = c/
alpha;
535 refGrad() = (qr + qrNbr)/kappaTp;
540 valueFraction() = KDeltaNbr/(KDeltaNbr +
alpha);
542 refGrad() = (qr + qrNbr)/kappaTp;
548 if (this->useImplicit())
553 valueFraction()*deltaH()
554 + (qr + qrNbr)/beta()
558 mixedFvPatchScalarField::updateCoeffs();
562 scalar Q =
gSum(kappaTp*patch().magSf()*snGrad());
564 Info<< patch().boundaryMesh().mesh().name() <<
':'
565 << patch().name() <<
':'
566 << this->internalField().name() <<
" <- "
569 << this->internalField().name() <<
" :"
570 <<
" heat transfer rate:" << Q
571 <<
" walltemperature "
572 <<
" min:" <<
gMin(Tp)
573 <<
" max:" <<
gMax(Tp)
591 <<
"This T BC does not support energy coupling "
592 <<
"It is implemented on he field "
605 <<
"This BC does not support energy coupling "
606 <<
"Use compressible::turbulentTemperatureRadCoupledMixed "
607 <<
"which has more functionalities and it can handle "
608 <<
"the assemble coupled option for energy. "
616turbulentTemperatureRadCoupledMixedFvPatchScalarField::alphaSfDelta()
const
618 return (
alpha(*
this)*patch().deltaCoeffs()*patch().magSf());
622tmp<scalarField> turbulentTemperatureRadCoupledMixedFvPatchScalarField::
625 const mappedPatchBase& mpp =
626 refCast<const mappedPatchBase>(patch().patch());
628 if (!mpp.sameWorld())
631 <<
"coupled energy not supported in combination with multi-world"
635 const label samplePatchi = mpp.samplePolyPatch().index();
636 const polyMesh& nbrMesh = mpp.sampleMesh();
638 const fvPatch& nbrPatch =
639 refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
650 mpp.distribute(TcNbr);
654 nbrField.alpha(TcNbr)*nbrPatch.deltaCoeffs()
656 mpp.distribute(alphaDeltaNbr);
663 return (alphaDeltaNbr + alphaDelta);
667tmp<scalarField> turbulentTemperatureRadCoupledMixedFvPatchScalarField::
670 const mappedPatchBase& mpp =
671 refCast<const mappedPatchBase>(
patch().
patch());
673 if (!mpp.sameWorld())
676 <<
"coupled energy not supported in combination with multi-world"
680 const polyMesh& nbrMesh = mpp.sampleMesh();
682 const basicThermo* nbrThermo =
685 const polyMesh&
mesh =
patch().boundaryMesh().mesh();
687 const basicThermo* localThermo =
691 if (nbrThermo && localThermo)
693 const label patchi =
patch().index();
694 const scalarField& pp = localThermo->p().boundaryField()[patchi];
697 const mappedPatchBase& mpp =
698 refCast<const mappedPatchBase>(
patch().
patch());
700 const label patchiNrb = mpp.samplePolyPatch().index();
702 const scalarField& ppNbr = nbrThermo->p().boundaryField()[patchiNrb];
712 - localThermo->he(pp, Tp, patchi)
713 + nbrThermo->he(ppNbr, Tp, patchiNrb)
719 <<
"Can't find thermos on mapped patch "
720 <<
" method, but thermo package not available"
733 mixedFvPatchScalarField::write(
os);
749 thicknessLayer_().writeData(
os);
750 kappaLayer_().writeData(
os);
752 if (thicknessLayers_.
size())
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
void setSize(const label n)
Alias for resize()
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
static autoPtr< PatchFunction1< Type > > NewIfPresent(const polyPatch &pp, const word &entryName, const dictionary &dict, const bool faceValues=true)
An optional selector.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
scalar deltaTValue() const noexcept
Return time step value.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
void writeEntry(Ostream &os) const
Write the UList with its compound type.
void size(const label n)
Older name for setAddressableSize.
static int & msgType() noexcept
Message tag of standard messages.
Abstract base-class for fluid and solid thermodynamic properties.
static const word dictName
const dictionary & coeffs() const
Return const dictionary of the model.
Mixed boundary condition for temperature and radiation heat transfer to be used for in multiregion ca...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual tmp< scalarField > kappa(const scalarField &Tp) const
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void manipulateMatrix(fvMatrix< scalar > &m, const label iMatrix, const direction cmpt)
Manipulate matrix.
turbulentTemperatureRadCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
virtual bool write()
Write the output fields.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const Time & time() const
Return the top-level database.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void operator=(const UList< scalar > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
const scalarField & deltaCoeffs() const
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
const polyMesh & sampleMesh() const
Get the region mesh.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE)
const word & sampleRegion() const
Region to sample.
bool sameWorld() const
Is sample world the local world?
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF)
Check that patch is of correct type.
void distribute(const word &fieldName, Field< T > &newValues) const
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
label index() const noexcept
The index of this patch in the boundaryMesh.
Mesh consisting of general polyhedral cells.
Common functions used in temperature coupled boundaries.
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const std::string patch
OpenFOAM patch number as a std::string.
Type gSum(const FieldField< Field, Type > &f)
To & refCast(From &r)
Reference type cast template function.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
errorManip< error > abort(error &err)
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Type gMax(const FieldField< Field, Type > &f)
#define forAll(list, i)
Loop across all elements in list.
static const char *const typeName
The type name used in ensight case files.