Loading...
Searching...
No Matches
absorptionEmissionModel.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) 2011-2015 OpenFOAM Foundation
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
30// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31
32namespace Foam
33{
34 namespace radiation
35 {
38 }
39}
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
44(
45 const dictionary& dict,
46 const fvMesh& mesh
47)
48:
51{}
52
53
54// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
55
60// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61
64{
65 return aDisp(bandI) + aCont(bandI);
66}
67
68
71{
73 (
74 "aCont",
86 (
87 "aDisp",
93
94
97{
98 return eDisp(bandI) + eCont(bandI);
99}
100
101
104{
106 (
107 "eCont",
119 (
120 "eDisp",
126
127
130{
131 return EDisp(bandI) + ECont(bandI);
132}
133
134
137{
139 (
140 "ECont",
152 (
153 "EDisp",
158}
159
160
165
166
171}
172
175{
176 return false;
177}
178
179
181(
184) const
185{
186 a = this->a();
187 aj[0] = a;
188}
189
190
191// ************************************************************************* //
label n
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Templated 2D Vector derived from VectorSpace adding construction from 2 components,...
Definition Vector2D.H:54
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
Model to supply absorption and emission coefficients for radiation modelling.
virtual tmp< volScalarField > e(const label bandI=0) const
Emission coefficient (net).
virtual bool isGrey() const
Flag for whether the absorption/emission is for a grey gas.
virtual const Vector2D< scalar > & bands(const label n) const
Const access to the bands - defaults to Vector2D::one for grey.
const fvMesh & mesh_
Reference to the fvMesh.
virtual tmp< volScalarField > aDisp(const label bandI=0) const
Absorption coefficient for dispersed phase.
virtual tmp< volScalarField > E(const label bandI=0) const
Emission contribution (net).
virtual tmp< volScalarField > EDisp(const label bandI=0) const
Emission contribution for dispersed phase.
virtual label nBands() const
Const access to the number of bands - defaults to 1 for grey.
virtual tmp< volScalarField > a(const label bandI=0) const
Absorption coefficient (net).
const dictionary dict_
Radiation model dictionary.
virtual tmp< volScalarField > eDisp(const label bandI=0) const
Return emission coefficient for dispersed phase.
const fvMesh & mesh() const
Reference to the mesh.
virtual tmp< volScalarField > eCont(const label bandI=0) const
Return emission coefficient for continuous phase.
absorptionEmissionModel(const dictionary &dict, const fvMesh &mesh)
Construct from components.
const dictionary & dict() const
Reference to the dictionary.
virtual tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
virtual void correct(volScalarField &a, PtrList< volScalarField > &aj) const
Correct absorption coefficients.
virtual tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
Namespace for radiation modelling.
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict