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) 2011-2015 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 "fluidThermo.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 
46  "QrExt";
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 (
52  const fvMesh& mesh
53 ) const
54 {
55  IOobject io
56  (
57  "radiationProperties",
58  mesh.time().constant(),
59  mesh,
60  IOobject::MUST_READ,
61  IOobject::NO_WRITE
62  );
63 
64  if (io.headerOk())
65  {
66  io.readOpt() = IOobject::MUST_READ_IF_MODIFIED;
67  return io;
68  }
69  else
70  {
71  io.readOpt() = IOobject::NO_READ;
72  return io;
73  }
74 }
75 
76 
78 {
79  if (radiation_)
80  {
81  solverFreq_ = max(1, lookupOrDefault<label>("solverFreq", 1));
82 
84  (
86  );
87 
88  scatter_.reset(scatterModel::New(*this, mesh_).ptr());
89 
90  soot_.reset(sootModel::New(*this, mesh_).ptr());
91 
92  transmissivity_.reset(transmissivityModel::New(*this, mesh_).ptr());
93  }
94 }
95 
96 
97 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
98 
100 :
102  (
103  IOobject
104  (
105  "radiationProperties",
106  T.time().constant(),
107  T.mesh(),
108  IOobject::NO_READ,
109  IOobject::NO_WRITE
110  )
111  ),
112  mesh_(T.mesh()),
113  time_(T.time()),
114  T_(T),
115  radiation_(false),
116  coeffs_(dictionary::null),
117  solverFreq_(0),
118  firstIter_(true),
119  absorptionEmission_(NULL),
120  scatter_(NULL),
121  soot_(NULL),
122  transmissivity_(NULL)
123 {}
124 
125 
127 (
128  const word& type,
129  const volScalarField& T
130 )
131 :
132  IOdictionary(createIOobject(T.mesh())),
133  mesh_(T.mesh()),
134  time_(T.time()),
135  T_(T),
136  radiation_(lookupOrDefault("radiation", true)),
137  coeffs_(subOrEmptyDict(type + "Coeffs")),
138  solverFreq_(1),
139  firstIter_(true),
140  absorptionEmission_(NULL),
141  scatter_(NULL),
142  soot_(NULL),
143  transmissivity_(NULL)
144 {
145  if (readOpt() == IOobject::NO_READ)
146  {
147  radiation_ = false;
148  }
149 
150  initialise();
151 }
152 
153 
155 (
156  const word& type,
157  const dictionary& dict,
158  const volScalarField& T
159 )
160 :
162  (
163  IOobject
164  (
165  "radiationProperties",
166  T.time().constant(),
167  T.mesh(),
170  ),
171  dict
172  ),
173  mesh_(T.mesh()),
174  time_(T.time()),
175  T_(T),
176  radiation_(lookupOrDefault("radiation", true)),
177  coeffs_(subOrEmptyDict(type + "Coeffs")),
178  solverFreq_(1),
179  firstIter_(true),
180  absorptionEmission_(NULL),
181  scatter_(NULL),
182  soot_(NULL),
183  transmissivity_(NULL)
184 {
185  initialise();
186 }
187 
188 
189 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
190 
192 {}
193 
194 
195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 
198 {
199  if (regIOobject::read())
200  {
201  lookup("radiation") >> radiation_;
202  coeffs_ = subOrEmptyDict(type() + "Coeffs");
203 
204  solverFreq_ = lookupOrDefault<label>("solverFreq", 1);
205  solverFreq_ = max(1, solverFreq_);
206 
207  return true;
208  }
209  else
210  {
211  return false;
212  }
213 }
214 
215 
217 {
218  if (!radiation_)
219  {
220  return;
221  }
222 
223  if (firstIter_ || (time_.timeIndex() % solverFreq_ == 0))
224  {
225  calculate();
226  firstIter_ = false;
227  }
228 
229  if (!soot_.empty())
230  {
231  soot_->correct();
232  }
233 }
234 
235 
237 (
239 ) const
240 {
241  volScalarField& he = thermo.he();
242  const volScalarField Cpv(thermo.Cpv());
243  const volScalarField T3(pow3(T_));
244 
245  return
246  (
247  Ru()
248  - fvm::Sp(4.0*Rp()*T3/Cpv, he)
249  - Rp()*T3*(T_ - 4.0*he/Cpv)
250  );
251 }
252 
253 
255 (
256  const dimensionedScalar& rhoCp,
258 ) const
259 {
260  return
261  (
262  Ru()/rhoCp
263  - fvm::Sp(Rp()*pow3(T)/rhoCp, T)
264  );
265 }
266 
267 
270 {
271  if (!absorptionEmission_.valid())
272  {
274  << "Requested radiation absorptionEmission model, but model is "
275  << "not activate" << abort(FatalError);
276  }
277 
278  return absorptionEmission_();
279 }
280 
281 
284 {
285  if (!soot_.valid())
286  {
288  << "Requested radiation sootModel model, but model is "
289  << "not activate" << abort(FatalError);
290  }
291 
292  return soot_();
293 }
294 
295 
298 {
299  if (!transmissivity_.valid())
300  {
302  << "Requested radiation sootModel model, but model is "
303  << "not activate" << abort(FatalError);
304  }
305 
306  return transmissivity_();
307 }
308 
309 
310 // ************************************************************************* //
Foam::radiation::radiationModel::scatter_
autoPtr< scatterModel > scatter_
Scatter model.
Definition: radiationModel.H:119
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::radiation::radiationModel::solverFreq_
label solverFreq_
Radiation solver frequency - number flow solver iterations per.
Definition: radiationModel.H:107
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::radiation::radiationModel::transmissivity
const transmissivityModel & transmissivity() const
Access to transmissivity Model.
Definition: radiationModel.C:297
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::radiation::absorptionEmissionModel::New
static autoPtr< absorptionEmissionModel > New(const dictionary &dict, const fvMesh &mesh)
Selector.
Definition: absorptionEmissionModelNew.C:33
scatterModel.H
Foam::radiation::radiationModel::radiationModel
radiationModel(const radiationModel &)
Disallow default bitwise copy construct.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
fluidThermo.H
Foam::regIOobject::read
virtual bool read()
Read object.
Definition: regIOobjectRead.C:171
Foam::fluidThermo
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:49
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::radiation::transmissivityModel::New
static autoPtr< transmissivityModel > New(const dictionary &dict, const fvMesh &mesh)
Definition: transmissivityModelNew.C:33
Foam::radiation::transmissivityModel
Base class for radiation scattering.
Definition: transmissivityModel.H:50
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
constant
Constant dispersed-phase particle diameter model.
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::IOobject::headerOk
bool headerOk()
Read and check header info.
Definition: IOobject.C:439
Foam::radiation::scatterModel::New
static autoPtr< scatterModel > New(const dictionary &dict, const fvMesh &mesh)
Definition: scatterModelNew.C:32
Foam::radiation::radiationModel::read
virtual bool read()=0
Read radiationProperties dictionary.
Definition: radiationModel.C:197
Foam::radiation::sootModel::New
static autoPtr< sootModel > New(const dictionary &dict, const fvMesh &mesh)
Selector.
Definition: sootModelNew.C:33
Foam::radiation::radiationModel::soot
const sootModel & soot() const
Access to soot Model.
Definition: radiationModel.C:283
Foam::radiation::radiationModel::soot_
autoPtr< sootModel > soot_
Soot model.
Definition: radiationModel.H:122
Foam::IOobject::readOpt
readOption readOpt() const
Definition: IOobject.H:317
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::radiation::radiationModel::~radiationModel
virtual ~radiationModel()
Destructor.
Definition: radiationModel.C:191
Foam::radiation::radiationModel::createIOobject
IOobject createIOobject(const fvMesh &mesh) const
Create IO object if dictionary is present.
Definition: radiationModel.C:51
sootModel.H
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::radiation::radiationModel::ST
virtual tmp< fvScalarMatrix > ST(const dimensionedScalar &rhoCp, volScalarField &T) const
Temperature source term.
Definition: radiationModel.C:255
Foam::radiation::radiationModel::correct
virtual void correct()
Main update/correction routine.
Definition: radiationModel.C:216
Foam::radiation::radiationModel::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: radiationModel.H:91
Foam::radiation::radiationModel::absorptionEmission
const absorptionEmissionModel & absorptionEmission() const
Access to absorptionEmission model.
Definition: radiationModel.C:269
Foam::radiation::defineTypeNameAndDebug
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::radiation::defineRunTimeSelectionTable
defineRunTimeSelectionTable(radiationModel, T)
Foam::radiation::radiationModel::Sh
virtual tmp< fvScalarMatrix > Sh(fluidThermo &thermo) const
Energy source term.
Definition: radiationModel.C:237
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
he
volScalarField & he
Definition: YEEqn.H:56
absorptionEmissionModel.H
Foam::radiation::radiationModel::initialise
void initialise()
Initialise.
Definition: radiationModel.C:77
Foam::radiation::sootModel
Base class for soor models.
Definition: sootModel.H:53
Foam::radiation::absorptionEmissionModel
Model to supply absorption and emission coefficients for radiation modelling.
Definition: absorptionEmissionModel.H:52
Foam::radiation::radiationModel::absorptionEmission_
autoPtr< absorptionEmissionModel > absorptionEmission_
Absorption/emission model.
Definition: radiationModel.H:116
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::radiation::radiationModel::radiation_
Switch radiation_
Radiation model on/off flag.
Definition: radiationModel.H:100
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::radiation::radiationModel::transmissivity_
autoPtr< transmissivityModel > transmissivity_
Transmissivity model.
Definition: radiationModel.H:125
radiationModel.H
Foam::radiation::radiationModel::externalRadHeatFieldName_
static const word externalRadHeatFieldName_
Static name external radiative fluxes.
Definition: radiationModel.H:82
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress