cloudScatter.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-2018 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 
26 #include "cloudScatter.H"
28 #include "thermoCloud.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  namespace radiation
35  {
36  defineTypeNameAndDebug(cloudScatter, 0);
37 
39  (
40  scatterModel,
41  cloudScatter,
42  dictionary
43  );
44  }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const dictionary& dict,
53  const fvMesh& mesh
54 )
55 :
56  scatterModel(dict, mesh),
57  coeffsDict_(dict.subDict(typeName + "Coeffs")),
58  cloudNames_(coeffsDict_.lookup("cloudNames"))
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
72 {
73  tmp<volScalarField> tsigma
74  (
75  new volScalarField
76  (
77  IOobject
78  (
79  "sigma",
80  mesh_.time().timeName(),
81  mesh_,
84  false
85  ),
86  mesh_,
87  dimensionedScalar("zero", dimless/dimLength, 0.0)
88  )
89  );
90 
91  forAll(cloudNames_, i)
92  {
93  const thermoCloud& tc
94  (
95  mesh_.objectRegistry::lookupObject<thermoCloud>(cloudNames_[i])
96  );
97 
98  tsigma.ref() += tc.sigmap();
99  }
100 
101  return 3.0*tsigma;
102 }
103 
104 
105 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual ~cloudScatter()
Destructor.
Definition: cloudScatter.C:64
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Macros for easy insertion into run-time selection tables.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
tmp< volScalarField > sigmaEff() const
Return scatter coefficient.
Definition: cloudScatter.C:71
Virtual abstract base class for templated ThermoCloud.
Definition: thermoCloud.H:48
cloudScatter(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: cloudScatter.C:51
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
Base class for radiation scattering.
Definition: scatterModel.H:50
virtual tmp< volScalarField > sigmap() const =0
Return tmp equivalent particulate scattering factor.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.