greyMeanSolidAbsorptionEmission.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-2017 OpenFOAM Foundation
9  Copyright (C) 2015 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "unitConversion.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  namespace radiation
39  {
40  defineTypeNameAndDebug(greyMeanSolidAbsorptionEmission, 0);
41 
43  (
44  absorptionEmissionModel,
45  greyMeanSolidAbsorptionEmission,
46  dictionary
47  );
48  }
49 }
50 
51 // * * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * //
52 
53 Foam::tmp<Foam::scalarField> Foam::radiation::
54 greyMeanSolidAbsorptionEmission::X(const word specie) const
55 {
56  const volScalarField& T = thermo_.T();
57  const volScalarField& p = thermo_.p();
58 
59  tmp<scalarField> tXj(new scalarField(T.primitiveField().size(), Zero));
60  scalarField& Xj = tXj.ref();
61 
62  tmp<scalarField> tRhoInv(new scalarField(T.primitiveField().size(), Zero));
63  scalarField& rhoInv = tRhoInv.ref();
64 
65  forAll(mixture_.Y(), specieI)
66  {
67  const scalarField& Yi = mixture_.Y()[specieI];
68 
69  forAll(rhoInv, iCell)
70  {
71  rhoInv[iCell] +=
72  Yi[iCell]/mixture_.rho(specieI, p[iCell], T[iCell]);
73  }
74  }
75  const scalarField& Yj = mixture_.Y(specie);
76  const label mySpecieI = mixture_.species()[specie];
77  forAll(Xj, iCell)
78  {
79  Xj[iCell] = Yj[iCell]/mixture_.rho(mySpecieI, p[iCell], T[iCell]);
80  }
81 
82  return (Xj/rhoInv);
83 }
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
89 (
90  const dictionary& dict,
91  const fvMesh& mesh
92 )
93 :
95  coeffsDict_((dict.optionalSubDict(typeName + "Coeffs"))),
96  thermo_(mesh.lookupObject<solidThermo>(basicThermo::dictName)),
97  speciesNames_(0),
98  mixture_(dynamic_cast<const basicSpecieMixture&>(thermo_)),
99  solidData_(mixture_.Y().size())
100 {
101  if (!isA<basicSpecieMixture>(thermo_))
102  {
104  << "Model requires a multi-component thermo package"
105  << abort(FatalError);
106  }
107 
108  label nFunc = 0;
109  const dictionary& functionDicts = dict.optionalSubDict(typeName + "Coeffs");
110 
111  for (const entry& dEntry : functionDicts)
112  {
113  if (!dEntry.isDict()) // safety
114  {
115  continue;
116  }
117 
118  const word& key = dEntry.keyword();
119  const dictionary& dict = dEntry.dict();
120 
121  if (!mixture_.contains(key))
122  {
124  << " specie: " << key << " is not found in the solid mixture"
125  << nl
126  << " specie is the mixture are:" << mixture_.species() << nl
127  << nl << endl;
128  }
129  speciesNames_.insert(key, nFunc);
130 
131  dict.readEntry("absorptivity", solidData_[nFunc][absorptivity]);
132  dict.readEntry("emissivity", solidData_[nFunc][emissivity]);
133 
134  nFunc++;
135  }
136 }
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
142 Foam::radiation::greyMeanSolidAbsorptionEmission::
143 calc(const label propertyId) const
144 {
146  (
147  new volScalarField
148  (
149  IOobject
150  (
151  "a",
152  mesh().time().timeName(),
153  mesh(),
156  ),
157  mesh(),
159  extrapolatedCalculatedFvPatchVectorField::typeName
160  )
161  );
162 
163  scalarField& a = ta.ref().primitiveFieldRef();
164 
165  forAllConstIters(speciesNames_, iter)
166  {
167  if (mixture_.contains(iter.key()))
168  {
169  a += solidData_[iter()][propertyId]*X(iter.key());
170  }
171  }
172 
173  ta.ref().correctBoundaryConditions();
174  return ta;
175 }
176 
177 
180 (
181  const label bandI
182 ) const
183 {
184  return calc(emissivity);
185 }
186 
187 
190 (
191  const label bandI
192 ) const
193 {
194  return calc(absorptivity);
195 }
196 
197 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::basicSpecieMixture::rho
virtual scalar rho(const label speciei, const scalar p, const scalar T) const =0
Density [kg/m3].
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::basicSpecieMixture
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Definition: basicSpecieMixture.H:58
Foam::basicThermo::p
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:602
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::basicMultiComponentMixture::species
const speciesTable & species() const
Return the table of species.
Definition: basicMultiComponentMixtureI.H:29
Foam::glTF::key
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:108
Foam::radiation::greyMeanSolidAbsorptionEmission::eCont
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
Definition: greyMeanSolidAbsorptionEmission.C:180
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::basicMultiComponentMixture::Y
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Definition: basicMultiComponentMixtureI.H:69
greyMeanSolidAbsorptionEmission.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::solidThermo
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:52
Foam::Field< scalar >
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::radiation::greyMeanSolidAbsorptionEmission::aCont
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
Definition: greyMeanSolidAbsorptionEmission.C:190
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::radiation::greyMeanSolidAbsorptionEmission::greyMeanSolidAbsorptionEmission
greyMeanSolidAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: greyMeanSolidAbsorptionEmission.C:89
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::radiation::defineTypeNameAndDebug
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::basicThermo::T
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:614
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::radiation::absorptionEmissionModel
Model to supply absorption and emission coefficients for radiation modelling.
Definition: absorptionEmissionModel.H:54
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::radiation::addToRunTimeSelectionTable
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
extrapolatedCalculatedFvPatchFields.H
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60