localDensityAbsorptionEmission.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) 2017 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace radiation
37  {
38  defineTypeNameAndDebug(localDensityAbsorptionEmission, 0);
39 
41  (
42  absorptionEmissionModel,
43  localDensityAbsorptionEmission,
44  dictionary
45  );
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 Foam::radiation::localDensityAbsorptionEmission::alpha(word alphaName) const
54 {
55  const volScalarField* ptr = mesh_.cfindObject<volScalarField>(alphaName);
56 
57  if (!ptr)
58  {
60  << "Unable to retrieve density field " << alphaName << " from "
61  << "database. Available objects:" << mesh_.sortedNames()
62  << exit(FatalError);
63  }
64 
65  return *ptr;
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
72 (
73  const dictionary& dict,
74  const fvMesh& mesh
75 )
76 :
77  absorptionEmissionModel(dict, mesh),
78  coeffsDict_(dict.subDict(typeName + "Coeffs")),
79  alphaNames_(coeffsDict_.lookup("alphaNames")),
80  aCoeff_(coeffsDict_.lookup("aCoeff")),
81  eCoeff_(coeffsDict_.lookup("eCoeff")),
82  ECoeff_(coeffsDict_.lookup("ECoeff"))
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
90 {
91  auto ta = volScalarField::New
92  (
93  "a",
95  mesh_,
97  );
98  auto& a = ta.ref();
99 
100  forAll(alphaNames_, i)
101  {
102  dimensionedScalar aPhase("a", dimless/dimLength, aCoeff_[i]);
103  a += max(alpha(alphaNames_[i]), scalar(0))*aPhase;
104  }
105 
106  return ta;
107 }
108 
109 
112 {
113  auto te = volScalarField::New
114  (
115  "e",
117  mesh_,
119  );
120  auto& e = te.ref();
121 
122  forAll(alphaNames_, i)
123  {
124  dimensionedScalar ePhase("e", dimless/dimLength, eCoeff_[i]);
125  e += max(alpha(alphaNames_[i]), scalar(0))*ePhase;
126  }
127 
128  return te;
129 }
130 
131 
134 {
135  auto tE = volScalarField::New
136  (
137  "E",
139  mesh_,
141  );
142 
143  scalarField& E = tE.ref().primitiveFieldRef();
144 
145  forAll(alphaNames_, i)
146  {
147  dimensionedScalar EPhase
148  (
149  "E",
151  ECoeff_[i]
152  );
153 
154  E += max(alpha(alphaNames_[i]), scalar(0))*EPhase;
155  }
156 
157  return tE;
158 }
159 
160 
161 // ************************************************************************* //
dictionary dict
localDensityAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
const dimensionSet dimless
Dimensionless.
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
const fvMesh & mesh_
Reference to the fvMesh.
Macros for easy insertion into run-time selection tables.
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
wordList sortedNames() const
The sorted names of all objects.
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Do not request registration (bool: false)
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127