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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2018 OpenFOAM Foundation
9 Copyright (C) 2019-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "fvDOM.H"
31#include "scatterModel.H"
32#include "constants.H"
33#include "unitConversion.H"
34#include "fvm.H"
36
37using namespace Foam::constant;
38using namespace Foam::constant::mathematical;
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
44 namespace radiation
45 {
48 }
49}
50
51
52// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53
55{
56 // Rotate Y spherical cordinates to Sun direction.
57 // Solid angles on the equator are better fit for planar radiation
58 const tensor coordRot = rotationTensor(vector(0, 1, 0), sunDir);
59
60 forAll(IRay_, rayId)
61 {
62 IRay_[rayId].dAve() = coordRot & IRay_[rayId].dAve();
63 IRay_[rayId].d() = coordRot & IRay_[rayId].d();
64 }
65}
66
67
69{
70 label SunRayId(-1);
71 scalar maxSunRay = -GREAT;
72
73 // Looking for the ray closest to the Sun direction
74 forAll(IRay_, rayId)
75 {
76 const vector& iD = IRay_[rayId].d();
77 scalar dir = sunDir & iD;
78 if (dir > maxSunRay)
79 {
80 maxSunRay = dir;
81 SunRayId = rayId;
82 }
83 }
84
85 // Second rotation to align colimated radiation with the closest ray
86 const tensor coordRot = rotationTensor(IRay_[SunRayId].d(), sunDir);
87
88 forAll(IRay_, rayId)
89 {
90 IRay_[rayId].dAve() = coordRot & IRay_[rayId].dAve();
91 IRay_[rayId].d() = coordRot & IRay_[rayId].d();
92 }
93
94 Info << "Sun direction : " << sunDir << nl << endl;
95 Info << "Sun ray ID : " << SunRayId << nl << endl;
96}
97
98
100{
101 solarCalculator_->correctSunDirection();
102 const vector sunDir = solarCalculator_->direction();
103
104 // First iteration
105 if (updateTimeIndex_ == 0)
106 {
107 rotateInitialRays(sunDir);
108 alignClosestRayToSun(sunDir);
109 }
110 else if (updateTimeIndex_ > 0)
111 {
112 alignClosestRayToSun(sunDir);
113 }
114}
115
116
117void Foam::radiation::fvDOM::initialise()
118{
119 coeffs_.readIfPresent("useExternalBeam", useExternalBeam_);
120
121 if (useExternalBeam_)
122 {
123 spectralDistributions_.reset
124 (
126 (
127 "spectralDistribution",
128 coeffs_,
129 &mesh_
130 )
131 );
132
133 spectralDistribution_ =
134 spectralDistributions_->value(mesh_.time().timeOutputValue());
135
136 spectralDistribution_ =
137 spectralDistribution_/sum(spectralDistribution_);
138
139 const dictionary& solarDict = this->subDict("solarCalculatorCoeffs");
140 solarCalculator_.reset(new solarCalculator(solarDict, mesh_));
141
142 if (mesh_.nSolutionD() != 3)
143 {
145 << "External beam model only available in 3D meshes "
146 << abort(FatalError);
147 }
148
149 if (solarCalculator_->diffuseSolarRad() > 0)
150 {
152 << "External beam model does not support Diffuse "
153 << "Solar Radiation. Set diffuseSolarRad to zero"
154 << abort(FatalError);
155 }
156 if (spectralDistribution_.size() != nLambda_)
157 {
159 << "The epectral energy distribution has different bands "
160 << "than the absoprtivity model "
161 << abort(FatalError);
162 }
163 }
164
165 // 3D
166 if (mesh_.nSolutionD() == 3)
167 {
168 nRay_ = 4*nPhi_*nTheta_;
169
170 IRay_.setSize(nRay_);
171
172 const scalar deltaPhi = pi/(2*nPhi_);
173 const scalar deltaTheta = pi/nTheta_;
174
175 label i = 0;
176
177 for (label n = 1; n <= nTheta_; n++)
178 {
179 for (label m = 1; m <= 4*nPhi_; m++)
180 {
181 scalar thetai = (2*n - 1)*deltaTheta/2.0;
182 scalar phii = (2*m - 1)*deltaPhi/2.0;
183
184 IRay_.set
185 (
186 i,
188 (
189 *this,
190 mesh_,
191 phii,
192 thetai,
193 deltaPhi,
194 deltaTheta,
195 nLambda_,
196 *absorptionEmission_,
197 blackBody_,
198 i
199 )
200 );
201 i++;
202 }
203 }
204 }
205 // 2D
206 else if (mesh_.nSolutionD() == 2)
207 {
208 const scalar thetai = piByTwo;
209 const scalar deltaTheta = pi;
210 nRay_ = 4*nPhi_;
211 IRay_.setSize(nRay_);
212 const scalar deltaPhi = pi/(2.0*nPhi_);
213 label i = 0;
214 for (label m = 1; m <= 4*nPhi_; m++)
215 {
216 const scalar phii = (2*m - 1)*deltaPhi/2.0;
217 IRay_.set
218 (
219 i,
221 (
222 *this,
223 mesh_,
224 phii,
225 thetai,
226 deltaPhi,
227 deltaTheta,
228 nLambda_,
229 *absorptionEmission_,
230 blackBody_,
231 i
232 )
233 );
234 i++;
235 }
236 }
237 // 1D
238 else
239 {
240 const scalar thetai = piByTwo;
241 const scalar deltaTheta = pi;
242 nRay_ = 2;
243 IRay_.setSize(nRay_);
244 const scalar deltaPhi = pi;
245 label i = 0;
246 for (label m = 1; m <= 2; m++)
247 {
248 const scalar phii = (2*m - 1)*deltaPhi/2.0;
249 IRay_.set
250 (
251 i,
253 (
254 *this,
255 mesh_,
256 phii,
257 thetai,
258 deltaPhi,
259 deltaTheta,
260 nLambda_,
261 *absorptionEmission_,
262 blackBody_,
263 i
264 )
265 );
266 i++;
267 }
268 }
269
270
271 // Construct absorption field for each wavelength
272 forAll(aLambda_, lambdaI)
273 {
274 aLambda_.set
275 (
276 lambdaI,
278 (
280 (
281 "aLambda_" + Foam::name(lambdaI) ,
282 mesh_.time().timeName(),
283 mesh_,
286 ),
287 a_
288 )
289 );
290 }
291
292 Info<< "fvDOM : Allocated " << IRay_.size()
293 << " rays with average orientation:" << nl;
294
295 if (useExternalBeam_)
296 {
297 // Rotate rays for Sun direction
298 updateRaysDir();
299 }
300
301 scalar totalOmega = 0;
302 forAll(IRay_, rayId)
303 {
304 if (omegaMax_ < IRay_[rayId].omega())
305 {
306 omegaMax_ = IRay_[rayId].omega();
307 }
308 totalOmega += IRay_[rayId].omega();
309 Info<< '\t' << IRay_[rayId].I().name() << " : " << "dAve : "
310 << '\t' << IRay_[rayId].dAve() << " : " << "omega : "
311 << '\t' << IRay_[rayId].omega() << " : " << "d : "
312 << '\t' << IRay_[rayId].d() << nl;
313 }
314
315 Info << "Total omega : " << totalOmega << endl;
316
317 Info<< endl;
318
319 coeffs_.readIfPresent("useSolarLoad", useSolarLoad_);
320
321 if (useSolarLoad_)
322 {
323 if (useExternalBeam_)
324 {
326 << "External beam with fvDOM can not be used "
327 << "with the solar load model"
328 << abort(FatalError);
329 }
330 const dictionary& solarDict = this->subDict("solarLoadCoeffs");
331 solarLoad_.reset(new solarLoad(solarDict, T_));
332
333 if (solarLoad_->nBands() != this->nBands())
334 {
336 << "Requested solar radiation with fvDOM. Using "
337 << "different number of bands for the solar load is not allowed"
338 << abort(FatalError);
339 }
340
341 Info<< "Creating Solar Load Model " << nl;
342 }
343}
344
345
346// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
347
349:
350 radiationModel(typeName, T),
351 G_
352 (
354 (
355 "G",
356 mesh_.time().timeName(),
357 mesh_,
358 IOobject::NO_READ,
359 IOobject::AUTO_WRITE
360 ),
361 mesh_,
363 ),
364 qr_
365 (
367 (
368 "qr",
369 mesh_.time().timeName(),
370 mesh_,
371 IOobject::READ_IF_PRESENT,
372 IOobject::AUTO_WRITE
373 ),
374 mesh_,
376 ),
377 qem_
378 (
380 (
381 "qem",
382 mesh_.time().timeName(),
383 mesh_,
384 IOobject::NO_READ,
385 IOobject::NO_WRITE
386 ),
387 mesh_,
389 ),
390 qin_
391 (
393 (
394 "qin",
395 mesh_.time().timeName(),
396 mesh_,
397 IOobject::READ_IF_PRESENT,
398 IOobject::AUTO_WRITE
399 ),
400 mesh_,
402 ),
403 a_
404 (
406 (
407 "a",
408 mesh_.time().timeName(),
409 mesh_,
410 IOobject::NO_READ,
411 IOobject::NO_WRITE
412 ),
413 mesh_,
415 ),
416 nTheta_(coeffs_.get<label>("nTheta")),
417 nPhi_(coeffs_.get<label>("nPhi")),
418 nRay_(0),
419 nLambda_(absorptionEmission_->nBands()),
420 aLambda_(nLambda_),
421 blackBody_(nLambda_, T),
422 IRay_(0),
423 tolerance_
424 (
425 coeffs_.getOrDefaultCompat<scalar>
426 (
427 "tolerance",
428 {{"convergence", 1712}},
429 0
430 )
431 ),
432 maxIter_(coeffs_.getOrDefault<label>("maxIter", 50)),
433 omegaMax_(0),
434 useSolarLoad_(false),
435 solarLoad_(),
436 meshOrientation_
437 (
438 coeffs_.getOrDefault<vector>("meshOrientation", Zero)
439 ),
440 useExternalBeam_(false),
441 spectralDistribution_(),
442 spectralDistributions_(),
443 solarCalculator_(),
444 updateTimeIndex_(0)
445{
446 initialise();
447}
448
449
451(
452 const dictionary& dict,
453 const volScalarField& T
454)
455:
456 radiationModel(typeName, dict, T),
457 G_
458 (
460 (
461 "G",
462 mesh_.time().timeName(),
463 mesh_,
464 IOobject::READ_IF_PRESENT,
465 IOobject::AUTO_WRITE
466 ),
467 mesh_,
469 ),
470 qr_
471 (
473 (
474 "qr",
475 mesh_.time().timeName(),
476 mesh_,
477 IOobject::READ_IF_PRESENT,
478 IOobject::AUTO_WRITE
479 ),
480 mesh_,
482 ),
483 qem_
484 (
486 (
487 "qem",
488 mesh_.time().timeName(),
489 mesh_,
490 IOobject::NO_READ,
491 IOobject::NO_WRITE
492 ),
493 mesh_,
495 ),
496 qin_
497 (
499 (
500 "qin",
501 mesh_.time().timeName(),
502 mesh_,
503 IOobject::READ_IF_PRESENT,
504 IOobject::AUTO_WRITE
505 ),
506 mesh_,
508 ),
509 a_
510 (
512 (
513 "a",
514 mesh_.time().timeName(),
515 mesh_,
516 IOobject::NO_READ,
517 IOobject::NO_WRITE
518 ),
519 mesh_,
521 ),
522 nTheta_(coeffs_.get<label>("nTheta")),
523 nPhi_(coeffs_.get<label>("nPhi")),
524 nRay_(0),
525 nLambda_(absorptionEmission_->nBands()),
526 aLambda_(nLambda_),
527 blackBody_(nLambda_, T),
528 IRay_(0),
529 tolerance_
530 (
531 coeffs_.getOrDefaultCompat<scalar>
532 (
533 "tolerance",
534 {{"convergence", 1712}},
535 0
536 )
537 ),
538 maxIter_(coeffs_.getOrDefault<label>("maxIter", 50)),
539 omegaMax_(0),
540 useSolarLoad_(false),
541 solarLoad_(),
542 meshOrientation_
543 (
544 coeffs_.getOrDefault<vector>("meshOrientation", Zero)
545 ),
546 useExternalBeam_(false),
547 spectralDistribution_(),
548 spectralDistributions_(),
549 solarCalculator_(),
550 updateTimeIndex_(0)
551{
552 initialise();
553}
554
555
556// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
557
559{
561 {
562 // Only reading solution parameters - not changing ray geometry
563 coeffs_.readIfPresentCompat
564 (
565 "tolerance", {{"convergence", 1712}}, tolerance_
566 );
567 coeffs_.readIfPresent("maxIter", maxIter_);
568
569 return true;
570 }
571
572 return false;
573}
574
575
577{
578 absorptionEmission_->correct(a_, aLambda_);
579
580 updateBlackBodyEmission();
581
582 if (useSolarLoad_)
583 {
584 solarLoad_->calculate();
585 }
586
587 if (useExternalBeam_)
588 {
589 switch (solarCalculator_->sunDirectionModel())
590 {
592 {
593 break;
594 }
596 {
597 label updateIndex = label
598 (
599 mesh_.time().value()
600 /solarCalculator_->sunTrackingUpdateInterval()
601 );
602
603 if (updateIndex > updateTimeIndex_)
604 {
605 Info << "Updating Sun position..." << endl;
606 updateTimeIndex_ = updateIndex;
607 updateRaysDir();
608 }
609 break;
610 }
611 }
612 }
613
614 // Set rays convergence false
615 List<bool> rayIdConv(nRay_, false);
616
617 scalar maxResidual = 0;
618 label radIter = 0;
619 do
620 {
621 Info<< "Radiation solver iter: " << radIter << endl;
622
623 radIter++;
624 maxResidual = 0;
625 forAll(IRay_, rayI)
626 {
627 if (!rayIdConv[rayI])
628 {
629 scalar maxBandResidual = IRay_[rayI].correct();
630 maxResidual = max(maxBandResidual, maxResidual);
631
632 if (maxBandResidual < tolerance_)
633 {
634 rayIdConv[rayI] = true;
635 }
636 }
637 }
638
639 } while (maxResidual > tolerance_ && radIter < maxIter_);
640
641 updateG();
642}
643
644
646{
647 // Construct using contribution from first frequency band
649 (
651 (
653 (
654 "Rp",
655 mesh_.time().timeName(),
656 mesh_,
659 false
660 ),
661 (
662 4
664 *(aLambda_[0] - absorptionEmission_->aDisp(0)())
665 *blackBody_.deltaLambdaT(T_, absorptionEmission_->bands(0))
666 )
667 )
668 );
669
670 volScalarField& Rp=tRp.ref();
671
672 // Add contributions over remaining frequency bands
673 for (label j=1; j < nLambda_; j++)
674 {
675 Rp +=
676 (
677 4
679 *(aLambda_[j] - absorptionEmission_->aDisp(j)())
680 *blackBody_.deltaLambdaT(T_, absorptionEmission_->bands(j))
681 );
682 }
683
684 return tRp;
685}
686
687
690{
692 (
694 (
696 (
697 "Ru",
698 mesh_.time().timeName(),
699 mesh_,
702 false
703 ),
704 mesh_,
705 dimensionedScalar(dimensionSet(1, -1, -3, 0, 0), Zero)
706 )
707 );
708
710
711 // Sum contributions over all frequency bands
712 for (label j=0; j < nLambda_; j++)
713 {
714 // Compute total incident radiation within frequency band
716 (
717 IRay_[0].ILambda(j)()*IRay_[0].omega()
718 );
719
720 for (label rayI=1; rayI < nRay_; rayI++)
721 {
722 Gj.ref() += IRay_[rayI].ILambda(j)()*IRay_[rayI].omega();
723 }
724
725 Ru += (aLambda_[j]() - absorptionEmission_->aDisp(j)()())*Gj
726 - absorptionEmission_->ECont(j)()();
727 }
728
729 return tRu;
730}
731
732
733void Foam::radiation::fvDOM::updateBlackBodyEmission()
734{
735 for (label j=0; j < nLambda_; j++)
736 {
737 blackBody_.correct(j, absorptionEmission_->bands(j));
738 }
739}
740
741
743{
748
749 forAll(IRay_, rayI)
750 {
751 IRay_[rayI].addIntensity();
752 G_ += IRay_[rayI].I()*IRay_[rayI].omega();
753 qr_.boundaryFieldRef() += IRay_[rayI].qr().boundaryField();
754 qem_.boundaryFieldRef() += IRay_[rayI].qem().boundaryField();
755 qin_.boundaryFieldRef() += IRay_[rayI].qin().boundaryField();
756 }
757}
758
759
761(
762 const word& name,
763 label& rayId,
764 label& lambdaId
765) const
766{
767 // Assuming name is in the form: CHARS_rayId_lambdaId
768 const auto i1 = name.find('_');
769 const auto i2 = name.find('_', i1+1);
770
771 rayId = readLabel(name.substr(i1+1, i2-i1-1));
772 lambdaId = readLabel(name.substr(i2+1));
773}
774
775
777{
778 return solarCalculator_();
779}
780
781
782// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:121
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: fvDOM.C:689
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: fvDOM.C:645
const solarCalculator & solarCalc() const
Solar calculator.
Definition: fvDOM.C:776
void alignClosestRayToSun(const vector &sunDir)
Align closest ray to sunDir.
Definition: fvDOM.C:68
void updateG()
Update G and calculate total heat flux on boundary.
Definition: fvDOM.C:742
void rotateInitialRays(const vector &sunDir)
Rotate rays spheric equator to sunDir.
Definition: fvDOM.C:54
void updateRaysDir()
Rotate rays according to Sun direction.
Definition: fvDOM.C:99
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:761
bool read()
Read radiation properties dictionary.
Definition: fvDOM.C:558
void calculate()
Solve radiation equation(s)
Definition: fvDOM.C:576
Top level model for radiation modelling.
virtual bool read()=0
Read radiationProperties dictionary.
Radiation intensity for a ray in a given direction.
The solarLoad radiation model includes Sun primary hits, their reflective fluxes and diffusive sky ra...
Definition: solarLoad.H:165
A solar calculator model providing models for the solar direction and solar loads.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
word timeName
Definition: getTimeIndex.H:3
Mathematical constants.
constexpr scalar pi(M_PI)
constexpr scalar piByTwo(0.5 *M_PI)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:51
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition: label.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Vector< scalar > vector
Definition: vector.H:61
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define addToRadiationRunTimeSelectionTables(model)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Unit conversion functions.