All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
radiation.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) 2017-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 "radiation.H"
27 #include "fluidThermo.H"
28 #include "fvMatrices.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
38 
40  (
41  option,
42  radiation,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& sourceName,
54  const word& modelType,
55  const dictionary& dict,
56  const fvMesh& mesh
57 )
58 :
59  option(sourceName, modelType, dict, mesh)
60 {
61  const basicThermo& thermo =
62  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
63 
64  fieldNames_.setSize(1);
65  fieldNames_[0] = thermo.he().name();
66  applied_.setSize(fieldNames_.size(), false);
67 
68  radiation_ = radiationModel::New(thermo.T());
69 }
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
75 {
76  return option::read(dict);
77 }
78 
79 
81 (
82  const volScalarField& rho,
83  fvMatrix<scalar>& eqn,
84  const label fieldi
85 )
86 {
87  const basicThermo& thermo =
88  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
89 
90  radiation_->correct();
91 
92  eqn += radiation_->Sh(thermo, eqn.psi());
93 }
94 
95 
96 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:295
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:489
rhoReactionThermo & thermo
Definition: createFields.H:28
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
Macros for easy insertion into run-time selection tables.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:53
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
radiation(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Definition: radiation.C:52
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
autoPtr< radiationModel > radiation(radiationModel::New(T))
A special matrix type and solver, designed for finite volume solutions of scalar equations.
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: radiation.C:74
virtual void correct()=0
Update properties.
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible momentum equation.
Definition: radiation.C:81
Namespace for OpenFOAM.
Finite volume options abstract base class. Provides a base set of controls, e.g.: ...
Definition: fvOption.H:66