fvDOM.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 "fvDOM.H"
28 #include "scatterModel.H"
29 #include "constants.H"
30 #include "fvm.H"
32 
33 using namespace Foam::constant;
34 using namespace Foam::constant::mathematical;
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  namespace radiation
41  {
42  defineTypeNameAndDebug(fvDOM, 0);
44  }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 {
52  // 3D
53  if (mesh_.nSolutionD() == 3)
54  {
55  nRay_ = 4*nPhi_*nTheta_;
56  IRay_.setSize(nRay_);
57  scalar deltaPhi = pi/(2.0*nPhi_);
58  scalar deltaTheta = pi/nTheta_;
59  label i = 0;
60  for (label n = 1; n <= nTheta_; n++)
61  {
62  for (label m = 1; m <= 4*nPhi_; m++)
63  {
64  scalar thetai = (2.0*n - 1.0)*deltaTheta/2.0;
65  scalar phii = (2.0*m - 1.0)*deltaPhi/2.0;
66  IRay_.set
67  (
68  i,
70  (
71  *this,
72  mesh_,
73  phii,
74  thetai,
75  deltaPhi,
76  deltaTheta,
77  nLambda_,
78  absorptionEmission_,
79  blackBody_,
80  i
81  )
82  );
83  i++;
84  }
85  }
86  }
87  // 2D
88  else if (mesh_.nSolutionD() == 2)
89  {
90  scalar thetai = piByTwo;
91  scalar deltaTheta = pi;
92  nRay_ = 4*nPhi_;
93  IRay_.setSize(nRay_);
94  scalar deltaPhi = pi/(2.0*nPhi_);
95  label i = 0;
96  for (label m = 1; m <= 4*nPhi_; m++)
97  {
98  scalar phii = (2.0*m - 1.0)*deltaPhi/2.0;
99  IRay_.set
100  (
101  i,
103  (
104  *this,
105  mesh_,
106  phii,
107  thetai,
108  deltaPhi,
109  deltaTheta,
110  nLambda_,
111  absorptionEmission_,
112  blackBody_,
113  i
114  )
115  );
116  i++;
117  }
118  }
119  // 1D
120  else
121  {
122  scalar thetai = piByTwo;
123  scalar deltaTheta = pi;
124  nRay_ = 2;
125  IRay_.setSize(nRay_);
126  scalar deltaPhi = pi;
127  label i = 0;
128  for (label m = 1; m <= 2; m++)
129  {
130  scalar phii = (2.0*m - 1.0)*deltaPhi/2.0;
131  IRay_.set
132  (
133  i,
135  (
136  *this,
137  mesh_,
138  phii,
139  thetai,
140  deltaPhi,
141  deltaTheta,
142  nLambda_,
143  absorptionEmission_,
144  blackBody_,
145  i
146  )
147  );
148  i++;
149  }
150  }
151 
152 
153  // Construct absorption field for each wavelength
154  forAll(aLambda_, lambdaI)
155  {
156  aLambda_.set
157  (
158  lambdaI,
159  new volScalarField
160  (
161  IOobject
162  (
163  "aLambda_" + Foam::name(lambdaI) ,
164  mesh_.time().timeName(),
165  mesh_,
168  ),
169  a_
170  )
171  );
172  }
173 
174  Info<< "fvDOM : Allocated " << IRay_.size()
175  << " rays with average orientation:" << nl;
176 
177  if (cacheDiv_)
178  {
179  Info<< "Caching div fvMatrix..."<< endl;
180  for (label lambdaI = 0; lambdaI < nLambda_; lambdaI++)
181  {
182  fvRayDiv_[lambdaI].setSize(nRay_);
183 
184  forAll(IRay_, rayId)
185  {
186  const surfaceScalarField Ji(IRay_[rayId].dAve() & mesh_.Sf());
187  const volScalarField& iRayLambdaI =
188  IRay_[rayId].ILambda(lambdaI);
189 
190  fvRayDiv_[lambdaI].set
191  (
192  rayId,
193  new fvScalarMatrix
194  (
195  fvm::div(Ji, iRayLambdaI, "div(Ji,Ii_h)")
196  )
197  );
198  }
199  }
200  }
201 
202  forAll(IRay_, rayId)
203  {
204  if (omegaMax_ < IRay_[rayId].omega())
205  {
206  omegaMax_ = IRay_[rayId].omega();
207  }
208  Info<< '\t' << IRay_[rayId].I().name() << " : " << "dAve : "
209  << '\t' << IRay_[rayId].dAve() << nl;
210  }
211 
212  Info<< endl;
213 
214  if (this->found("useSolarLoad"))
215  {
216  this->lookup("useSolarLoad") >> useSolarLoad_;
217  }
218 
219  if (useSolarLoad_)
220  {
221  const dictionary& solarDict = this->subDict("solarLoarCoeffs");
222  solarLoad_.reset
223  (
224  new solarLoad(solarDict, T_, externalRadHeatFieldName_)
225  );
226 
227  if (solarLoad_->nBands() > 1)
228  {
230  << "Requested solar radiation with fvDOM. Using "
231  << "more than one band for the solar load is not allowed"
232  << abort(FatalError);
233  }
234 
235  Info<< "Creating Solar Load Model " << nl;
236  }
237 }
238 
239 
240 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
241 
243 :
244  radiationModel(typeName, T),
245  G_
246  (
247  IOobject
248  (
249  "G",
250  mesh_.time().timeName(),
251  mesh_,
252  IOobject::NO_READ,
253  IOobject::AUTO_WRITE
254  ),
255  mesh_,
257  ),
258  Qr_
259  (
260  IOobject
261  (
262  "Qr",
263  mesh_.time().timeName(),
264  mesh_,
265  IOobject::READ_IF_PRESENT,
266  IOobject::AUTO_WRITE
267  ),
268  mesh_,
270  ),
271  Qem_
272  (
273  IOobject
274  (
275  "Qem",
276  mesh_.time().timeName(),
277  mesh_,
278  IOobject::NO_READ,
279  IOobject::NO_WRITE
280  ),
281  mesh_,
282  dimensionedScalar("Qem", dimMass/pow3(dimTime), 0.0)
283  ),
284  Qin_
285  (
286  IOobject
287  (
288  "Qin",
289  mesh_.time().timeName(),
290  mesh_,
291  IOobject::READ_IF_PRESENT,
292  IOobject::AUTO_WRITE
293  ),
294  mesh_,
295  dimensionedScalar("Qin", dimMass/pow3(dimTime), 0.0)
296  ),
297  a_
298  (
299  IOobject
300  (
301  "a",
302  mesh_.time().timeName(),
303  mesh_,
304  IOobject::NO_READ,
305  IOobject::NO_WRITE
306  ),
307  mesh_,
309  ),
310  nTheta_(readLabel(coeffs_.lookup("nTheta"))),
311  nPhi_(readLabel(coeffs_.lookup("nPhi"))),
312  nRay_(0),
313  nLambda_(absorptionEmission_->nBands()),
314  aLambda_(nLambda_),
315  blackBody_(nLambda_, T),
316  IRay_(0),
317  convergence_(coeffs_.lookupOrDefault<scalar>("convergence", 0.0)),
318  maxIter_(coeffs_.lookupOrDefault<label>("maxIter", 50)),
319  fvRayDiv_(nLambda_),
320  cacheDiv_(coeffs_.lookupOrDefault<bool>("cacheDiv", false)),
321  omegaMax_(0),
322  useSolarLoad_(false),
323  solarLoad_(),
324  meshOrientation_
325  (
326  coeffs_.lookupOrDefault<vector>("meshOrientation", vector::zero)
327  )
328 {
329  initialise();
330 }
331 
332 
334 (
335  const dictionary& dict,
336  const volScalarField& T
337 )
338 :
339  radiationModel(typeName, dict, T),
340  G_
341  (
342  IOobject
343  (
344  "G",
345  mesh_.time().timeName(),
346  mesh_,
349  ),
350  mesh_,
352  ),
353  Qr_
354  (
355  IOobject
356  (
357  "Qr",
358  mesh_.time().timeName(),
359  mesh_,
362  ),
363  mesh_,
365  ),
366  Qem_
367  (
368  IOobject
369  (
370  "Qem",
371  mesh_.time().timeName(),
372  mesh_,
375  ),
376  mesh_,
377  dimensionedScalar("Qem", dimMass/pow3(dimTime), 0.0)
378  ),
379  Qin_
380  (
381  IOobject
382  (
383  "Qin",
384  mesh_.time().timeName(),
385  mesh_,
388  ),
389  mesh_,
390  dimensionedScalar("Qin", dimMass/pow3(dimTime), 0.0)
391  ),
392  a_
393  (
394  IOobject
395  (
396  "a",
397  mesh_.time().timeName(),
398  mesh_,
401  ),
402  mesh_,
404  ),
405  nTheta_(readLabel(coeffs_.lookup("nTheta"))),
406  nPhi_(readLabel(coeffs_.lookup("nPhi"))),
407  nRay_(0),
408  nLambda_(absorptionEmission_->nBands()),
409  aLambda_(nLambda_),
410  blackBody_(nLambda_, T),
411  IRay_(0),
412  convergence_(coeffs_.lookupOrDefault<scalar>("convergence", 0.0)),
413  maxIter_(coeffs_.lookupOrDefault<label>("maxIter", 50)),
414  fvRayDiv_(nLambda_),
415  cacheDiv_(coeffs_.lookupOrDefault<bool>("cacheDiv", false)),
416  omegaMax_(0),
417  useSolarLoad_(false),
418  solarLoad_(),
419  meshOrientation_
420  (
421  coeffs_.lookupOrDefault<vector>("meshOrientation", vector::zero)
422  )
423 {
424  initialise();
425 }
426 
427 
428 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
429 
431 {}
432 
433 
434 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
435 
437 {
438  if (radiationModel::read())
439  {
440  // Only reading solution parameters - not changing ray geometry
441  coeffs_.readIfPresent("convergence", convergence_);
442  coeffs_.readIfPresent("maxIter", maxIter_);
443 
444  return true;
445  }
446  else
447  {
448  return false;
449  }
450 }
451 
452 
454 {
455  absorptionEmission_->correct(a_, aLambda_);
456 
457  updateBlackBodyEmission();
458 
459  if (useSolarLoad_)
460  {
461  solarLoad_->calculate();
462  }
463 
464  // Set rays convergence false
465  List<bool> rayIdConv(nRay_, false);
466 
467  scalar maxResidual = 0.0;
468  label radIter = 0;
469  do
470  {
471  Info<< "Radiation solver iter: " << radIter << endl;
472 
473  radIter++;
474  maxResidual = 0.0;
475  forAll(IRay_, rayI)
476  {
477  if (!rayIdConv[rayI])
478  {
479  scalar maxBandResidual = IRay_[rayI].correct();
480  maxResidual = max(maxBandResidual, maxResidual);
481 
482  if (maxBandResidual < convergence_)
483  {
484  rayIdConv[rayI] = true;
485  }
486  }
487  }
488 
489  } while (maxResidual > convergence_ && radIter < maxIter_);
490 
491  updateG();
492 }
493 
494 
496 {
497  return tmp<volScalarField>
498  (
499  new volScalarField
500  (
501  IOobject
502  (
503  "Rp",
504  mesh_.time().timeName(),
505  mesh_,
508  false
509  ),
510  4.0*a_*physicoChemical::sigma //absorptionEmission_->a()
511  )
512  );
513 }
514 
515 
518 {
519 
521  G_.dimensionedInternalField();
523  absorptionEmission_->ECont()().dimensionedInternalField();
525  a_.dimensionedInternalField();
526 
527  return a*G - E;
528 }
529 
530 
532 {
533  for (label j=0; j < nLambda_; j++)
534  {
535  blackBody_.correct(j, absorptionEmission_->bands(j));
536  }
537 }
538 
539 
541 {
542  G_ = dimensionedScalar("zero",dimMass/pow3(dimTime), 0.0);
543  Qr_ = dimensionedScalar("zero",dimMass/pow3(dimTime), 0.0);
544  Qem_ = dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0);
545  Qin_ = dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0);
546 
547  forAll(IRay_, rayI)
548  {
549  IRay_[rayI].addIntensity();
550  G_ += IRay_[rayI].I()*IRay_[rayI].omega();
551  Qr_.boundaryField() += IRay_[rayI].Qr().boundaryField();
552  Qem_.boundaryField() += IRay_[rayI].Qem().boundaryField();
553  Qin_.boundaryField() += IRay_[rayI].Qin().boundaryField();
554  }
555 }
556 
557 
559 (
560  const word& name,
561  label& rayId,
562  label& lambdaId
563 ) const
564 {
565  // assuming name is in the form: CHARS_rayId_lambdaId
566  size_type i1 = name.find_first_of("_");
567  size_type i2 = name.find_last_of("_");
568 
569  rayId = readLabel(IStringStream(name.substr(i1+1, i2-1))());
570  lambdaId = readLabel(IStringStream(name.substr(i2+1, name.size()-1))());
571 }
572 
573 
574 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::radiation::fvDOM::~fvDOM
virtual ~fvDOM()
Destructor.
Definition: fvDOM.C:430
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
scatterModel.H
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::constant
Collection of constants.
Definition: atomicConstants.C:37
Foam::radiation::fvDOM::fvDOM
fvDOM(const fvDOM &)
Disallow default bitwise copy construct.
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::radiation::fvDOM::updateG
void updateG()
Update G and calculate total heat flux on boundary.
Definition: fvDOM.C:540
Foam::radiation::radiativeIntensityRay
Radiation intensity for a ray in a given direction.
Definition: radiativeIntensityRay.H:55
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
fvDOM.H
Foam::radiation::fvDOM::initialise
void initialise()
Initialise.
Definition: fvDOM.C:50
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::constant::mathematical::piByTwo
const scalar piByTwo(0.5 *pi)
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::radiation::radiationModel::read
virtual bool read()=0
Read radiationProperties dictionary.
Definition: radiationModel.C:197
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::radiation::fvDOM::updateBlackBodyEmission
void updateBlackBodyEmission()
Update nlack body emission.
Definition: fvDOM.C:531
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::radiation::fvDOM::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: fvDOM.C:495
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
size_type
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
Foam::radiation::fvDOM::calculate
void calculate()
Solve radiation equation(s)
Definition: fvDOM.C:453
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
dict
dictionary dict
Definition: searchingEngine.H:14
fvm.H
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
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::IStringStream
Input from memory buffer stream.
Definition: IStringStream.H:49
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
constants.H
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::radiation::fvDOM::read
bool read()
Read radiation properties dictionary.
Definition: fvDOM.C:436
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::constant::mathematical
mathematical constants.
Foam::radiation::solarLoad
The solar load radiation model includes Sun primary hits, their reflective fluxes and diffusive sky r...
Definition: solarLoad.H:80
absorptionEmissionModel.H
Foam::Vector< scalar >
Foam::List< bool >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::constant::mathematical::pi
const scalar pi(M_PI)
addToRadiationRunTimeSelectionTables
#define addToRadiationRunTimeSelectionTables(model)
Definition: radiationModel.H:263
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::radiation::fvDOM::Ru
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: fvDOM.C:517
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::radiation::fvDOM::setRayIdLambdaId
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:559
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::zero
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:47