Go to the documentation of this file.
44 wallAbsorptionEmissionModel,
54 const Foam::fvMesh& Foam::radiation::solidAbsorption::nbrRegion()
const
56 const mappedPatchBase& mpp = refCast<const mappedPatchBase>(
pp_);
57 return (refCast<const fvMesh>(mpp.sampleMesh()));
61 Foam::label Foam::radiation::solidAbsorption::nbrPatchIndex()
const
63 const mappedPatchBase& mpp = refCast<const mappedPatchBase>(pp_);
64 return (mpp.samplePolyPatch().index());
78 if (!isA<mappedPatchBase>(pp))
81 <<
"\n patch type '" << pp.type()
82 <<
"' not type '" << mappedPatchBase::typeName <<
"'"
83 <<
"\n for patch " << pp.
name()
109 const fvMesh& nbrMesh = nbrRegion();
114 "radiationProperties"
130 mpp.distribute(absorptivity);
147 return a(bandI,
nullptr,
nullptr)()[faceI];
163 const fvMesh& nbrMesh = nbrRegion();
168 "radiationProperties"
184 mpp.distribute(emissivity);
201 return e(bandI,
nullptr,
nullptr)()[faceI];
207 const fvMesh& nbrMesh = nbrRegion();
212 "radiationProperties"
215 return (
radiation.absorptionEmission().nBands());
A class for managing temporary objects.
solidAbsorption(const dictionary &dict, const polyPatch &pp)
Construct from components.
virtual ~solidAbsorption()
Destructor.
tmp< scalarField > e(const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Return emission coefficient.
const polyPatch & pp_
Reference to the polyPatch.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
A patch is a list of labels that address the faces in the global face list.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const Type & lookupObject(const word &name, const bool recursive=false) const
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
errorManip< error > abort(error &err)
Based class for wall absorption emission models.
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
tmp< scalarField > a(const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
absorptivity coefficient
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
label nBands() const
Number of bands.
static int & msgType() noexcept
Message tag of standard messages.
Top level model for radiation modelling.
const dimensionedScalar e
Elementary charge.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
const word & name() const noexcept
The patch name.
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)