Loading...
Searching...
No Matches
multiBandAbsorption.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2015-2018 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "multiBandAbsorption.H"
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33namespace Foam
34{
35 namespace radiation
36 {
38
40 (
44 );
45 }
46}
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
52(
53 const dictionary& dict,
54 const polyPatch& pp
55)
56:
58 coeffsDict_(dict),
59 aCoeffs_(),
60 eCoeffs_(),
61 nBands_(0)
62{
63 coeffsDict_.readEntry("absorptivity", aCoeffs_);
64 coeffsDict_.readEntry("emissivity", eCoeffs_);
65 nBands_ = aCoeffs_.size();
66}
67
68
69// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
70
72{}
73
74
75// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76
79(
80 const label bandI,
81 const vectorField* incomingDirection,
82 const scalarField* T
83) const
84{
85 return tmp<scalarField>::New(pp_.size(), aCoeffs_[bandI]);
86}
87
89(
90 const label faceI,
91 const label bandI,
92 const vector dir,
93 const scalar T
94) const
95{
96 return aCoeffs_[bandI];
97}
98
99
101(
102 const label bandI,
103 const vectorField* incomingDirection,
105) const
106{
107 return tmp<scalarField>::New(pp_.size(), eCoeffs_[bandI]);
108}
109
110
112(
113 const label faceI,
114 const label bandI,
115 const vector dir,
116 const scalar T
117) const
118{
119 return eCoeffs_[bandI];
120}
121
122
123// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
multiBandAbsorption radiation transmissivity for solids.
tmp< scalarField > a(const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
absorptivity coefficient
tmp< scalarField > e(const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Return emission coefficient.
multiBandAbsorption(const dictionary &dict, const polyPatch &pp)
Construct from components.
Based class for wall absorption emission models.
wallAbsorptionEmissionModel(const dictionary &dict, const polyPatch &pp)
Construct from components.
const polyPatch & pp_
Reference to the polyPatch.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
Namespace for radiation modelling.
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Vector< scalar > vector
Definition vector.H:57
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict