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-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 "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  defineTypeNameAndDebug(radiationModel, 0);
38  defineRunTimeSelectionTable(radiationModel, T);
39  defineRunTimeSelectionTable(radiationModel, dictionary);
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 Foam::IOobject Foam::radiationModel::createIOobject(const fvMesh& mesh) const
46 {
47  typeIOobject<IOdictionary> io
48  (
49  "radiationProperties",
50  mesh.time().constant(),
51  mesh,
54  );
55 
56  if (io.headerOk())
57  {
59  return io;
60  }
61  else
62  {
63  io.readOpt() = IOobject::NO_READ;
64  return io;
65  }
66 }
67 
68 
69 void Foam::radiationModel::initialise()
70 {
71  solverFreq_ = max(1, lookupOrDefault<label>("solverFreq", 1));
72 
73  absorptionEmission_.reset
74  (
76  );
77 
78  scatter_.reset(radiationModels::scatterModel::New(*this, mesh_).ptr());
79 
80  soot_.reset(radiationModels::sootModel::New(*this, mesh_).ptr());
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85 
87 :
89  (
90  IOobject
91  (
92  "radiationProperties",
93  T.time().constant(),
94  T.mesh(),
95  IOobject::NO_READ,
96  IOobject::NO_WRITE
97  )
98  ),
99  mesh_(T.mesh()),
100  time_(T.time()),
101  T_(T),
102  coeffs_(dictionary::null),
103  solverFreq_(0),
104  firstIter_(true),
105  absorptionEmission_(nullptr),
106  scatter_(nullptr),
107  soot_(nullptr)
108 {}
109 
110 
112 :
113  IOdictionary(createIOobject(T.mesh())),
114  mesh_(T.mesh()),
115  time_(T.time()),
116  T_(T),
117  coeffs_(subOrEmptyDict(type + "Coeffs")),
118  solverFreq_(1),
119  firstIter_(true),
120  absorptionEmission_(nullptr),
121  scatter_(nullptr),
122  soot_(nullptr)
123 {
124  initialise();
125 }
126 
127 
129 (
130  const word& type,
131  const dictionary& dict,
132  const volScalarField& T
133 )
134 :
136  (
137  IOobject
138  (
139  "radiationProperties",
140  T.time().constant(),
141  T.mesh(),
144  ),
145  dict
146  ),
147  mesh_(T.mesh()),
148  time_(T.time()),
149  T_(T),
150  coeffs_(subOrEmptyDict(type + "Coeffs")),
151  solverFreq_(1),
152  firstIter_(true),
153  absorptionEmission_(nullptr),
154  scatter_(nullptr),
155  soot_(nullptr)
156 {
157  initialise();
158 }
159 
160 
161 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
162 
164 {}
165 
166 
167 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
168 
170 {
171  if (firstIter_ || (time_.timeIndex() % solverFreq_ == 0))
172  {
173  calculate();
174  firstIter_ = false;
175  }
176 
177  if (!soot_.empty())
178  {
179  soot_->correct();
180  }
181 }
182 
183 
185 {
186  if (regIOobject::read())
187  {
188  coeffs_ = subOrEmptyDict(type() + "Coeffs");
189 
190  solverFreq_ = lookupOrDefault<label>("solverFreq", 1);
191  solverFreq_ = max(1, solverFreq_);
192 
193  return true;
194  }
195  else
196  {
197  return false;
198  }
199 }
200 
201 
203 (
204  const basicThermo& thermo,
205  const volScalarField& he
206 ) const
207 {
208  const volScalarField Cpv(thermo.Cpv());
209  const volScalarField T3(pow3(T_));
210 
211  return
212  (
213  Ru()
214  - fvm::Sp(4.0*Rp()*T3/Cpv, he)
215  - Rp()*T3*(T_ - 4.0*he/Cpv)
216  );
217 }
218 
219 
221 (
222  const dimensionedScalar& rhoCp,
223  volScalarField& T
224 ) const
225 {
226  return
227  (
228  Ru()/rhoCp
229  - fvm::Sp(Rp()*pow3(T)/rhoCp, T)
230  );
231 }
232 
233 
236 {
237  if (!absorptionEmission_.valid())
238  {
240  << "Requested radiation absorptionEmission model, but model is "
241  << "not active" << abort(FatalError);
242  }
243 
244  return absorptionEmission_();
245 }
246 
247 
249 {
250  if (!soot_.valid())
251  {
253  << "Requested radiation sootModel model, but model is "
254  << "not active" << abort(FatalError);
255  }
256 
257  return soot_();
258 }
259 
260 
261 // ************************************************************************* //
autoPtr< radiationModels::scatterModel > scatter_
Scatter model.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static autoPtr< scatterModel > New(const dictionary &dict, const fvMesh &mesh)
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:77
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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
bool firstIter_
Flag to enable radiation model to be evaluated on first iteration.
virtual bool read()=0
Read radiationProperties dictionary.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvMesh & mesh
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:53
static autoPtr< absorptionEmissionModel > New(const dictionary &dict, const fvMesh &mesh)
Selector.
const volScalarField & T_
Reference to the temperature field.
const radiationModels::absorptionEmissionModel & absorptionEmission() const
Access to absorptionEmission model.
virtual void correct()
Main update/correction routine.
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
virtual tmp< volScalarField::Internal > Ru() const =0
Source term component (constant)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
virtual ~radiationModel()
Destructor.
IOdictionary(const IOobject &io, const word &wantedType)
Construct given an IOobject, supply wanted typeName.
Definition: IOdictionary.C:47
virtual tmp< volScalarField > Rp() const =0
Source term component (for power of T^4)
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
const radiationModels::sootModel & soot() const
Access to soot Model.
label solverFreq_
Radiation solver frequency - number flow solver iterations per.
dimensionedScalar pow3(const dimensionedScalar &ds)
autoPtr< radiationModels::sootModel > soot_
Soot model.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
virtual void calculate()=0
Solve radiation equation(s)
static autoPtr< sootModel > New(const dictionary &dict, const fvMesh &mesh)
Selector.
Definition: sootModelNew.C:35
const Time & time() const
Return time.
Definition: IOobject.C:318
radiationModel(const volScalarField &T)
Null constructor.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Base class for soot models.
Definition: sootModel.H:50
const fvMesh & mesh_
Reference to the mesh database.
Model to supply absorption and emission coefficients for radiation modelling.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< fvScalarMatrix > ST(const dimensionedScalar &rhoCp, volScalarField &T) const
Temperature source term.
readOption readOpt() const
Definition: IOobject.H:365
virtual tmp< fvScalarMatrix > Sh(const basicThermo &thermo, const volScalarField &he) const
Energy source term.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
autoPtr< radiationModels::absorptionEmissionModel > absorptionEmission_
Absorption/emission model.
const Time & time_
Reference to the time database.
dictionary coeffs_
Radiation model dictionary.
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:1033
Namespace for OpenFOAM.