Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37namespace Foam
38{
39 namespace radiation
40 {
42
44 }
45}
46
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50Foam::radiation::opaqueSolid::opaqueSolid(const volScalarField& T)
51:
52 radiationModel(typeName, T)
53{}
54
55
56Foam::radiation::opaqueSolid::opaqueSolid
57(
58 const dictionary& dict,
59 const volScalarField& T
60)
62 radiationModel(typeName, dict, T)
63{}
64
65
66// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
71
72// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75{
76 return radiationModel::read();
90 mesh_,
92 (
103 (
104 "Ru",
109}
110
111
112Foam::label Foam::radiation::opaqueSolid::nBands() const
113{
114 return absorptionEmission_->nBands();
115}
116
117
118// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
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())
@ NO_REGISTER
Do not request registration (bool: false).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Radiation for solid opaque solids - does nothing to energy equation source terms (returns zeros) but ...
Definition opaqueSolid.H:56
tmp< volScalarField::Internal > Ru() const
Source term component (constant).
Definition opaqueSolid.C:93
tmp< volScalarField > Rp() const
Source term component (for power of T^4).
Definition opaqueSolid.C:77
label nBands() const
Number of bands.
virtual ~opaqueSolid()
Destructor.
Definition opaqueSolid.C:61
bool read()
Read radiationProperties dictionary.
Definition opaqueSolid.C:67
void calculate()
Solve radiation equation(s).
Definition opaqueSolid.C:73
const fvMesh & mesh_
Reference to the mesh database.
virtual bool read()=0
Read radiationProperties dictionary.
autoPtr< absorptionEmissionModel > absorptionEmission_
Absorption/emission model.
const volScalarField & T() const noexcept
Return access to the temperature field.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Namespace for radiation modelling.
Namespace for OpenFOAM.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
#define addToRadiationRunTimeSelectionTables(model)
dictionary dict