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-2023 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,
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& sourceName,
54  const word& modelType,
55  const fvMesh& mesh,
56  const dictionary& dict
57 )
58 :
59  fvModel(sourceName, modelType, mesh, dict),
60  radiation_
61  (
63  (
64  mesh.lookupObject<basicThermo>(physicalProperties::typeName).T()
65  )
66  )
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
73 {
74  const basicThermo& thermo =
75  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
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>(physicalProperties::typeName);
90 
91  radiation_->correct();
92 
93  eqn += radiation_->Sh(thermo, eqn.psi());
94 }
95 
96 
98 {}
99 
100 
102 {}
103 
104 
106 {}
107 
108 
110 {
111  return true;
112 }
113 
114 
115 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
virtual void correct()=0
Update properties.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
Calculates and applies the radiation source to the energy equation.
Definition: radiation.H:65
virtual bool movePoints()
Update for mesh motion.
Definition: radiation.C:109
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
Definition: radiation.C:72
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 void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: radiation.C:97
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: radiation.C:105
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: radiation.C:101
radiation(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: radiation.C:52
An abstract base class for physical properties.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Top level model for radiation modelling.
A class for handling words, derived from string.
Definition: word.H:62
A special matrix type and solver, designed for finite volume solutions of scalar equations.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31