opaqueSolid.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-2019 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 radiationModels
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 
81 
83 {
84  return volScalarField::New
85  (
86  "Rp",
87  mesh_,
89  (
91  0
92  )
93  );
94 }
95 
96 
99 {
101  (
102  "Ru",
103  mesh_,
105  );
106 }
107 
108 
109 // ************************************************************************* //
virtual ~opaqueSolid()
Destructor.
Definition: opaqueSolid.C:66
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: opaqueSolid.C:98
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
virtual bool read()=0
Read radiationProperties dictionary.
opaqueSolid(const volScalarField &T)
Construct from components.
Definition: opaqueSolid.C:48
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
void calculate()
Solve radiation equation(s)
Definition: opaqueSolid.C:78
Top level model for radiation modelling.
bool read()
Read radiationProperties dictionary.
Definition: opaqueSolid.C:72
defineTypeNameAndDebug(absorptionEmissionModel, 0)
addToRadiationRunTimeSelectionTables(fvDOM)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
radiationModel(const volScalarField &T)
Null constructor.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const fvMesh & mesh_
Reference to the mesh database.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: opaqueSolid.C:82
Namespace for OpenFOAM.