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-2021 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  fvModel,
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  fvModel(sourceName, modelType, dict, mesh),
60  radiation_
61  (
63  (
65  )
66  )
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
73 {
74  const basicThermo& thermo =
75  mesh().lookupObject<basicThermo>(basicThermo::dictName);
76 
77  return wordList(1, thermo.he().name());
78 }
79 
80 
82 (
83  const volScalarField& rho,
84  fvMatrix<scalar>& eqn,
85  const word& fieldName
86 ) const
87 {
88  const basicThermo& thermo =
89  mesh().lookupObject<basicThermo>(basicThermo::dictName);
90 
91  radiation_->correct();
92 
93  eqn += radiation_->Sh(thermo, eqn.psi());
94 }
95 
96 
97 // ************************************************************************* //
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
const word & name() const
Return name.
Definition: IOobject.H:303
fluidReactionThermo & thermo
Definition: createFields.H:28
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:77
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
Finite volume model abstract base class.
Definition: fvModel.H:55
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
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].
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to compressible energy equation.
Definition: radiation.C:82
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
Definition: radiation.C:72
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
static const word dictName
Name of the thermophysical properties dictionary.
Definition: basicThermo.H:129
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
virtual const volScalarField & T() const =0
Temperature [K].
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 void correct()=0
Update properties.
Namespace for OpenFOAM.