opaqueSolid.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 "opaqueSolid.H"
30 #include "fvMesh.H"
31 #include "Time.H"
32 #include "volFields.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  namespace radiation
40  {
41  defineTypeNameAndDebug(opaqueSolid, 0);
42 
44  }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 Foam::radiation::opaqueSolid::opaqueSolid(const volScalarField& T)
51 :
52  radiationModel(typeName, T)
53 {}
54 
55 
56 Foam::radiation::opaqueSolid::opaqueSolid
57 (
58  const dictionary& dict,
59  const volScalarField& T
60 )
61 :
62  radiationModel(typeName, dict, T)
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
75 {
76  return radiationModel::read();
77 }
78 
79 
81 {}
82 
83 
85 {
87  (
88  IOobject
89  (
90  "Rp",
91  mesh_.time().timeName(),
92  mesh_,
95  ),
96  mesh_,
98  (
100  )
101  );
102 }
103 
104 
107 {
109  (
110  IOobject
111  (
112  "Ru",
113  mesh_.time().timeName(),
114  mesh_,
117  ),
118  mesh_,
120  );
121 }
122 
123 
124 Foam::label Foam::radiation::opaqueSolid::nBands() const
125 {
126  return absorptionEmission_->nBands();
127 }
128 
129 // ************************************************************************* //
dictionary dict
Radiation for solid opaque solids - does nothing to energy equation source terms (returns zeros) but ...
Definition: opaqueSolid.H:53
addToRadiationRunTimeSelectionTables(fvDOM)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Ignore writing from objectRegistry::writeObject()
Macros for easy insertion into run-time selection tables.
virtual bool read()=0
Read radiationProperties dictionary.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Top level model for radiation modelling.
virtual ~opaqueSolid()
Destructor.
Definition: opaqueSolid.C:61
label nBands() const
Number of bands.
Definition: opaqueSolid.C:117
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const volScalarField & T
dimensionedScalar pow3(const dimensionedScalar &ds)
void calculate()
Solve radiation equation(s)
Definition: opaqueSolid.C:73
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
bool read()
Read radiationProperties dictionary.
Definition: opaqueSolid.C:67
tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: opaqueSolid.C:77
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: opaqueSolid.C:99
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127