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 | Copyright (C) 2011-2015 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 "opaqueSolid.H"
28 #include "fvMesh.H"
29 #include "Time.H"
30 #include "volFields.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  namespace radiation
38  {
39  defineTypeNameAndDebug(opaqueSolid, 0);
40 
42  }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 :
50  radiationModel(typeName, T)
51 {}
52 
53 
55 (
56  const dictionary& dict,
57  const volScalarField& T
58 )
59 :
60  radiationModel(typeName, dict, T)
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
73 {
74  return radiationModel::read();
75 }
76 
77 
79 {
80  // Do nothing
81 }
82 
83 
85 {
86  return tmp<volScalarField>
87  (
88  new volScalarField
89  (
90  IOobject
91  (
92  "Rp",
93  mesh_.time().timeName(),
94  mesh_,
97  ),
98  mesh_,
100  (
101  "Rp",
103  0.0
104  )
105  )
106  );
107 }
108 
109 
112 {
114  (
116  (
117  IOobject
118  (
119  "Ru",
120  mesh_.time().timeName(),
121  mesh_,
124  ),
125  mesh_,
127  (
128  "Ru", dimMass/dimLength/pow3(dimTime), 0.0
129  )
130  )
131  );
132 }
133 
134 
135 // ************************************************************************* //
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::radiation::opaqueSolid::Ru
tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: opaqueSolid.C:111
Foam::radiation::opaqueSolid::calculate
void calculate()
Solve radiation equation(s)
Definition: opaqueSolid.C:78
Foam::radiation::opaqueSolid::~opaqueSolid
virtual ~opaqueSolid()
Destructor.
Definition: opaqueSolid.C:66
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
opaqueSolid.H
Foam::radiation::opaqueSolid::read
bool read()
Read radiationProperties dictionary.
Definition: opaqueSolid.C:72
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::radiation::radiationModel::read
virtual bool read()=0
Read radiationProperties dictionary.
Definition: radiationModel.C:197
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::radiation::defineTypeNameAndDebug
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Foam::radiation::opaqueSolid::Rp
tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: opaqueSolid.C:84
Foam::radiation::addToRadiationRunTimeSelectionTables
addToRadiationRunTimeSelectionTables(fvDOM)
T
const volScalarField & T
Definition: createFields.H:25
physicoChemicalConstants.H
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
Foam::radiation::opaqueSolid::opaqueSolid
opaqueSolid(const opaqueSolid &)
Disallow default bitwise copy construct.
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51