absorptionEmissionModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  namespace radiationModels
33  {
34  defineTypeNameAndDebug(absorptionEmissionModel, 0);
35  defineRunTimeSelectionTable(absorptionEmissionModel, dictionary);
36  }
37 }
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const dictionary& dict,
44  const fvMesh& mesh
45 )
46 :
47  dict_(dict),
48  mesh_(mesh)
49 {}
50 
51 
52 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
53 
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
59 
62 {
63  return aDisp(bandI) + aCont(bandI);
64 }
65 
66 
69 {
70  return volScalarField::New
71  (
72  "aCont",
73  mesh_,
75  );
76 }
77 
78 
81 {
82  return volScalarField::New
83  (
84  "aDisp",
85  mesh_,
87  );
88 }
89 
90 
93 {
94  return eDisp(bandI) + eCont(bandI);
95 }
96 
97 
100 {
101  return volScalarField::New
102  (
103  "eCont",
104  mesh_,
106  );
107 }
108 
109 
112 {
113  return volScalarField::New
114  (
115  "eDisp",
116  mesh_,
118  );
119 }
120 
121 
124 {
125  return EDisp(bandI) + ECont(bandI);
126 }
127 
128 
131 {
132  return volScalarField::New
133  (
134  "ECont",
135  mesh_,
137  );
138 }
139 
140 
143 {
144  return volScalarField::New
145  (
146  "EDisp",
147  mesh_,
149  );
150 }
151 
152 
154 {
155  return pTraits<label>::one;
156 }
157 
158 
161 {
162  return Vector2D<scalar>::one;
163 }
164 
165 
167 {
168  return false;
169 }
170 
171 
173 (
174  volScalarField& a,
176 ) const
177 {
178  a = this->a();
179  aj[0] = a;
180 }
181 
182 
183 // ************************************************************************* //
virtual tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual tmp< volScalarField > E(const label bandI=0) const
Emission contribution (net)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Traits class for primitives.
Definition: pTraits.H:50
virtual label nBands() const
Const access to the number of bands - defaults to 1 for grey.
const dimensionSet dimless
virtual const Vector2D< scalar > & bands(const label n) const
Const access to the bands - defaults to Vector2D::one for grey.
virtual tmp< volScalarField > e(const label bandI=0) const
Emission coefficient (net)
const dimensionSet dimLength
defineRunTimeSelectionTable(absorptionEmissionModel, dictionary)
virtual bool isGrey() const
Flag for whether the absorption/emission is for a grey gas.
virtual tmp< volScalarField > a(const label bandI=0) const
Absorption coefficient (net)
const dimensionSet dimTime
virtual tmp< volScalarField > eDisp(const label bandI=0) const
Return emission coefficient for dispersed phase.
virtual tmp< volScalarField > aDisp(const label bandI=0) const
Absorption coefficient for dispersed phase.
absorptionEmissionModel(const dictionary &dict, const fvMesh &mesh)
Construct from components.
const dimensionSet dimMass
virtual void correct(volScalarField &a, PtrList< volScalarField > &aj) const
Correct absorption coefficients.
defineTypeNameAndDebug(absorptionEmissionModel, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual tmp< volScalarField > EDisp(const label bandI=0) const
Emission contribution for dispersed phase.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual tmp< volScalarField > eCont(const label bandI=0) const
Return emission coefficient for continuous phase.
label n
A class for managing temporary objects.
Definition: PtrList.H:53
Templated 2D Vector derived from VectorSpace adding construction from 2 components, element access using x() and y() member functions and the inner-product (dot-product).
Definition: Vector2D.H:51
Namespace for OpenFOAM.