FreeFOAM The Cross-Platform CFD Toolkit
radiationModel.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 | Copyright (C) 1991-2010 OpenCFD Ltd.
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 "radiationModel.H"
28 #include <radiation/scatterModel.H>
29 #include <finiteVolume/fvm.H>
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace radiation
36  {
37  defineTypeNameAndDebug(radiationModel, 0);
38  defineRunTimeSelectionTable(radiationModel, dictionary);
39  }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::radiation::radiationModel::radiationModel(const volScalarField& T)
46 :
48  (
49  IOobject
50  (
51  "radiationProperties",
52  T.time().constant(),
53  T.mesh(),
54  IOobject::MUST_READ,
55  IOobject::NO_WRITE
56  )
57  ),
58  mesh_(T.mesh()),
59  time_(T.time()),
60  T_(T),
61  radiation_(false),
62  coeffs_(dictionary::null),
63  solverFreq_(0),
64  absorptionEmission_(NULL),
65  scatter_(NULL)
66 {}
67 
68 
69 Foam::radiation::radiationModel::radiationModel
70 (
71  const word& type,
72  const volScalarField& T
73 )
74 :
76  (
77  IOobject
78  (
79  "radiationProperties",
80  T.time().constant(),
81  T.mesh(),
84  )
85  ),
86  mesh_(T.mesh()),
87  time_(T.time()),
88  T_(T),
89  radiation_(lookup("radiation")),
90  coeffs_(subDict(type + "Coeffs")),
91  solverFreq_(readLabel(lookup("solverFreq"))),
92  absorptionEmission_(absorptionEmissionModel::New(*this, mesh_)),
93  scatter_(scatterModel::New(*this, mesh_))
94 {
95  solverFreq_ = max(1, solverFreq_);
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
100 
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 {
109  if (regIOobject::read())
110  {
111  lookup("radiation") >> radiation_;
112  coeffs_ = subDict(type() + "Coeffs");
113 
114  return true;
115  }
116  else
117  {
118  return false;
119  }
120 }
121 
122 
124 {
125  if (!radiation_)
126  {
127  return;
128  }
129 
130  if (time_.timeIndex() % solverFreq_ == 0)
131  {
132  calculate();
133  }
134 }
135 
136 
138 (
140 ) const
141 {
142  volScalarField& h = thermo.h();
143  const volScalarField cp = thermo.Cp();
144  const volScalarField T3 = pow3(T_);
145 
146  return
147  (
148  Ru()
149  - fvm::Sp(4.0*Rp()*T3/cp, h)
150  - Rp()*T3*(T_ - 4.0*h/cp)
151  );
152 }
153 
155 (
156  basicThermo& thermo
157 ) const
158 {
159  volScalarField& hs = thermo.hs();
160  const volScalarField cp = thermo.Cp();
161  const volScalarField T3 = pow3(T_);
162 
163  return
164  (
165  Ru()
166  - fvm::Sp(4.0*Rp()*T3/cp, hs)
167  - Rp()*T3*(T_ - 4.0*hs/cp)
168  );
169 }
170 
171 // ************************ vim: set sw=4 sts=4 et: ************************ //