cloudScatter.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-2016 OpenFOAM Foundation
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 
28 #include "cloudScatter.H"
30 #include "thermoCloud.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace radiation
37  {
38  defineTypeNameAndDebug(cloudScatter, 0);
39 
41  (
42  scatterModel,
43  cloudScatter,
44  dictionary
45  );
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const dictionary& dict,
55  const fvMesh& mesh
56 )
57 :
59  coeffsDict_(dict.subDict(typeName + "Coeffs")),
60  cloudNames_(coeffsDict_.lookup("cloudNames"))
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
74 {
75  auto tsigma = volScalarField::New
76  (
77  "sigma",
79  mesh_,
81  );
82 
83  for (const word& cloudName : cloudNames_)
84  {
85  const auto& tc =
86  mesh_.objectRegistry::lookupObject<thermoCloud>(cloudName);
87 
88  tsigma.ref() += tc.sigmap();
89  }
90 
91  return 3.0*tsigma;
92 }
93 
94 
95 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual ~cloudScatter()
Destructor.
Definition: cloudScatter.C:59
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
const dimensionSet dimless
Dimensionless.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Macros for easy insertion into run-time selection tables.
dynamicFvMesh & mesh
const word cloudName(propsDict.get< word >("cloud"))
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< volScalarField > sigmaEff() const
Return scatter coefficient.
Definition: cloudScatter.C:66
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...
Virtual abstract base class for templated ThermoCloud.
Definition: thermoCloud.H:46
cloudScatter(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: cloudScatter.C:46
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:46
virtual tmp< volScalarField > sigmap() const =0
Return tmp equivalent particulate scattering factor.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
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