radiationModel.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) 2011-2018 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 "radiationModel.H"
28 #include "scatterModel.H"
29 #include "sootModel.H"
30 #include "fvmSup.H"
31 #include "basicThermo.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  namespace radiation
38  {
39  defineTypeNameAndDebug(radiationModel, 0);
40  defineRunTimeSelectionTable(radiationModel, T);
41  defineRunTimeSelectionTable(radiationModel, dictionary);
42  }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 Foam::IOobject Foam::radiation::radiationModel::createIOobject
49 (
50  const fvMesh& mesh
51 ) const
52 {
53  IOobject io
54  (
55  "radiationProperties",
56  mesh.time().constant(),
57  mesh,
60  );
61 
62  if (io.typeHeaderOk<IOdictionary>(true))
63  {
65  return io;
66  }
67  else
68  {
69  io.readOpt() = IOobject::NO_READ;
70  return io;
71  }
72 }
73 
74 
75 void Foam::radiation::radiationModel::initialise()
76 {
77  solverFreq_ = max(1, lookupOrDefault<label>("solverFreq", 1));
78 
79  absorptionEmission_.reset
80  (
81  absorptionEmissionModel::New(*this, mesh_).ptr()
82  );
83 
84  scatter_.reset(scatterModel::New(*this, mesh_).ptr());
85 
86  soot_.reset(sootModel::New(*this, mesh_).ptr());
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
92 Foam::radiation::radiationModel::radiationModel(const volScalarField& T)
93 :
95  (
96  IOobject
97  (
98  "radiationProperties",
99  T.time().constant(),
100  T.mesh(),
101  IOobject::NO_READ,
102  IOobject::NO_WRITE
103  )
104  ),
105  mesh_(T.mesh()),
106  time_(T.time()),
107  T_(T),
108  coeffs_(dictionary::null),
109  solverFreq_(0),
110  firstIter_(true),
111  absorptionEmission_(nullptr),
112  scatter_(nullptr),
113  soot_(nullptr)
114 {}
115 
116 
117 Foam::radiation::radiationModel::radiationModel
118 (
119  const word& type,
120  const volScalarField& T
121 )
122 :
123  IOdictionary(createIOobject(T.mesh())),
124  mesh_(T.mesh()),
125  time_(T.time()),
126  T_(T),
127  coeffs_(subOrEmptyDict(type + "Coeffs")),
128  solverFreq_(1),
129  firstIter_(true),
130  absorptionEmission_(nullptr),
131  scatter_(nullptr),
132  soot_(nullptr)
133 {
134  initialise();
135 }
136 
137 
138 Foam::radiation::radiationModel::radiationModel
139 (
140  const word& type,
141  const dictionary& dict,
142  const volScalarField& T
143 )
144 :
146  (
147  IOobject
148  (
149  "radiationProperties",
150  T.time().constant(),
151  T.mesh(),
154  ),
155  dict
156  ),
157  mesh_(T.mesh()),
158  time_(T.time()),
159  T_(T),
160  coeffs_(subOrEmptyDict(type + "Coeffs")),
161  solverFreq_(1),
162  firstIter_(true),
163  absorptionEmission_(nullptr),
164  scatter_(nullptr),
165  soot_(nullptr)
166 {
167  initialise();
168 }
169 
170 
171 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
172 
174 {}
175 
176 
177 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
178 
180 {
181  if (firstIter_ || (time_.timeIndex() % solverFreq_ == 0))
182  {
183  calculate();
184  firstIter_ = false;
185  }
186 
187  if (!soot_.empty())
188  {
189  soot_->correct();
190  }
191 }
192 
193 
195 {
196  if (regIOobject::read())
197  {
198  coeffs_ = subOrEmptyDict(type() + "Coeffs");
199 
200  solverFreq_ = lookupOrDefault<label>("solverFreq", 1);
201  solverFreq_ = max(1, solverFreq_);
202 
203  return true;
204  }
205  else
206  {
207  return false;
208  }
209 }
210 
211 
213 (
214  const basicThermo& thermo,
215  const volScalarField& he
216 ) const
217 {
218  const volScalarField Cpv(thermo.Cpv());
219  const volScalarField T3(pow3(T_));
220 
221  return
222  (
223  Ru()
224  - fvm::Sp(4.0*Rp()*T3/Cpv, he)
225  - Rp()*T3*(T_ - 4.0*he/Cpv)
226  );
227 }
228 
229 
231 (
232  const dimensionedScalar& rhoCp,
233  volScalarField& T
234 ) const
235 {
236  return
237  (
238  Ru()/rhoCp
239  - fvm::Sp(Rp()*pow3(T)/rhoCp, T)
240  );
241 }
242 
243 
246 {
247  if (!absorptionEmission_.valid())
248  {
250  << "Requested radiation absorptionEmission model, but model is "
251  << "not activate" << abort(FatalError);
252  }
253 
254  return absorptionEmission_();
255 }
256 
257 
260 {
261  if (!soot_.valid())
262  {
264  << "Requested radiation sootModel model, but model is "
265  << "not activate" << abort(FatalError);
266  }
267 
268  return soot_();
269 }
270 
271 
272 // ************************************************************************* //
static autoPtr< scatterModel > New(const dictionary &dict, const fvMesh &mesh)
dictionary coeffs_
Radiation model dictionary.
autoPtr< scatterModel > scatter_
Scatter model.
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
virtual bool read()
Read object.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const volScalarField & T_
Reference to the temperature field.
static autoPtr< sootModel > New(const dictionary &dict, const fvMesh &mesh)
Selector.
Definition: sootModelNew.C:33
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
const sootModel & soot() const
Access to soot Model.
virtual void calculate()=0
Solve radiation equation(s)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Model to supply absorption and emission coefficients for radiation modelling.
defineRunTimeSelectionTable(radiationModel, T)
virtual tmp< volScalarField > Cpv() const =0
Heat capacity at constant pressure/volume [J/kg/K].
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
virtual ~radiationModel()
Destructor.
dynamicFvMesh & mesh
Base class for soor models.
Definition: sootModel.H:53
virtual bool read()=0
Read radiationProperties dictionary.
autoPtr< sootModel > soot_
Soot model.
A class for handling words, derived from string.
Definition: word.H:59
static autoPtr< absorptionEmissionModel > New(const dictionary &dict, const fvMesh &mesh)
Selector.
bool firstIter_
Flag to enable radiation model to be evaluated on first iteration.
const fvMesh & mesh_
Reference to the mesh database.
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
virtual void correct()
Main update/correction routine.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const Mesh & mesh() const
Return mesh.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
virtual tmp< volScalarField::Internal > Ru() const =0
Source term component (constant)
virtual tmp< volScalarField > Rp() const =0
Source term component (for power of T^4)
const absorptionEmissionModel & absorptionEmission() const
Access to absorptionEmission model.
virtual tmp< fvScalarMatrix > Sh(const basicThermo &thermo, const volScalarField &he) const
Energy source term.
dimensionedScalar pow3(const dimensionedScalar &ds)
IOdictionary(const IOobject &)
Construct given an IOobject.
Definition: IOdictionary.C:30
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
const Time & time_
Reference to the time database.
virtual tmp< fvScalarMatrix > ST(const dimensionedScalar &rhoCp, volScalarField &T) const
Temperature source term.
const Time & time() const
Return time.
Definition: IOobject.C:367
A class for managing temporary objects.
Definition: PtrList.H:53
label solverFreq_
Radiation solver frequency - number flow solver iterations per.
readOption readOpt() const
Definition: IOobject.H:359
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
autoPtr< absorptionEmissionModel > absorptionEmission_
Absorption/emission model.
Calculate the matrix for implicit and explicit sources.
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:727
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.