Loading...
Searching...
No Matches
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-2023 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;
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 (
125 Function1<scalarField>::New
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 (
279 IOobject
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
348Foam::radiation::fvDOM::fvDOM(const volScalarField& T)
349:
351 G_
352 (
354 (
355 "G",
356 mesh_.time().timeName(),
357 mesh_,
358 IOobject::NO_READ,
359 IOobject::AUTO_WRITE,
360 IOobject::REGISTER
361 ),
362 mesh_,
364 ),
365 qr_
366 (
368 (
369 "qr",
370 mesh_.time().timeName(),
371 mesh_,
372 IOobject::READ_IF_PRESENT,
373 IOobject::AUTO_WRITE,
374 IOobject::REGISTER
375 ),
376 mesh_,
378 ),
379 qem_
380 (
382 (
383 "qem",
384 mesh_.time().timeName(),
385 mesh_,
386 IOobject::NO_READ,
387 IOobject::NO_WRITE
388 ),
389 mesh_,
391 ),
392 qin_
393 (
395 (
396 "qin",
397 mesh_.time().timeName(),
398 mesh_,
399 IOobject::READ_IF_PRESENT,
400 IOobject::AUTO_WRITE,
401 IOobject::REGISTER
402 ),
403 mesh_,
405 ),
406 a_
407 (
409 (
410 "a",
411 mesh_.time().timeName(),
412 mesh_,
413 IOobject::NO_READ,
414 IOobject::NO_WRITE
415 ),
416 mesh_,
418 ),
419 nTheta_(coeffs_.get<label>("nTheta")),
420 nPhi_(coeffs_.get<label>("nPhi")),
421 nRay_(0),
422 nLambda_(absorptionEmission_->nBands()),
423 aLambda_(nLambda_),
424 blackBody_(nLambda_, T),
425 IRay_(0),
426 tolerance_
427 (
428 coeffs_.getOrDefaultCompat<scalar>
429 (
430 "tolerance",
431 {{"convergence", 1712}},
432 0
433 )
434 ),
435 maxIter_(coeffs_.getOrDefault<label>("maxIter", 50)),
436 omegaMax_(0),
437 useSolarLoad_(false),
438 solarLoad_(),
439 meshOrientation_
440 (
441 coeffs_.getOrDefault<vector>("meshOrientation", Zero)
442 ),
443 useExternalBeam_(false),
444 spectralDistribution_(),
445 spectralDistributions_(),
446 solarCalculator_(),
447 updateTimeIndex_(0)
448{
449 initialise();
450}
451
452
453Foam::radiation::fvDOM::fvDOM
454(
455 const dictionary& dict,
456 const volScalarField& T
457)
458:
460 G_
461 (
463 (
464 "G",
465 mesh_.time().timeName(),
466 mesh_,
467 IOobject::READ_IF_PRESENT,
468 IOobject::AUTO_WRITE,
469 IOobject::REGISTER
470 ),
471 mesh_,
473 ),
474 qr_
475 (
477 (
478 "qr",
479 mesh_.time().timeName(),
480 mesh_,
481 IOobject::READ_IF_PRESENT,
482 IOobject::AUTO_WRITE,
483 IOobject::REGISTER
484 ),
485 mesh_,
487 ),
488 qem_
489 (
491 (
492 "qem",
493 mesh_.time().timeName(),
494 mesh_,
495 IOobject::NO_READ,
496 IOobject::NO_WRITE
497 ),
498 mesh_,
500 ),
501 qin_
502 (
504 (
505 "qin",
506 mesh_.time().timeName(),
507 mesh_,
508 IOobject::READ_IF_PRESENT,
509 IOobject::AUTO_WRITE,
510 IOobject::REGISTER
511 ),
512 mesh_,
514 ),
515 a_
516 (
518 (
519 "a",
520 mesh_.time().timeName(),
521 mesh_,
522 IOobject::NO_READ,
523 IOobject::NO_WRITE
524 ),
525 mesh_,
527 ),
528 nTheta_(coeffs_.get<label>("nTheta")),
529 nPhi_(coeffs_.get<label>("nPhi")),
530 nRay_(0),
531 nLambda_(absorptionEmission_->nBands()),
532 aLambda_(nLambda_),
533 blackBody_(nLambda_, T),
534 IRay_(0),
535 tolerance_
536 (
537 coeffs_.getOrDefaultCompat<scalar>
538 (
539 "tolerance",
540 {{"convergence", 1712}},
541 0
542 )
543 ),
544 maxIter_(coeffs_.getOrDefault<label>("maxIter", 50)),
545 omegaMax_(0),
546 useSolarLoad_(false),
547 solarLoad_(),
548 meshOrientation_
549 (
550 coeffs_.getOrDefault<vector>("meshOrientation", Zero)
551 ),
552 useExternalBeam_(false),
553 spectralDistribution_(),
554 spectralDistributions_(),
555 solarCalculator_(),
556 updateTimeIndex_(0)
558 initialise();
559}
560
561
562// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
563
565{
567 {
568 // Only reading solution parameters - not changing ray geometry
569 coeffs_.readIfPresentCompat
570 (
571 "tolerance", {{"convergence", 1712}}, tolerance_
572 );
573 coeffs_.readIfPresent("maxIter", maxIter_);
574
575 return true;
576 }
577
578 return false;
579}
580
581
583{
584 absorptionEmission_->correct(a_, aLambda_);
585
586 updateBlackBodyEmission();
587
588 if (useSolarLoad_)
589 {
590 solarLoad_->calculate();
591 }
592
593 if (useExternalBeam_)
594 {
595 switch (solarCalculator_->sunDirectionModel())
596 {
598 {
599 break;
600 }
602 {
603 label updateIndex = label
604 (
605 mesh_.time().value()
606 /solarCalculator_->sunTrackingUpdateInterval()
607 );
608
609 if (updateIndex > updateTimeIndex_)
610 {
611 Info << "Updating Sun position..." << endl;
612 updateTimeIndex_ = updateIndex;
613 updateRaysDir();
614 }
615 break;
616 }
617 }
618 }
619
620 // Set rays convergence false
621 List<bool> rayIdConv(nRay_, false);
622
623 scalar maxResidual = 0;
624 label radIter = 0;
625 do
626 {
627 Info<< "Radiation solver iter: " << radIter << endl;
628
629 radIter++;
630 maxResidual = 0;
631 forAll(IRay_, rayI)
632 {
633 if (!rayIdConv[rayI])
634 {
635 scalar maxBandResidual = IRay_[rayI].correct();
636 maxResidual = max(maxBandResidual, maxResidual);
637
638 if (maxBandResidual < tolerance_)
639 {
640 rayIdConv[rayI] = true;
641 }
642 }
643 }
645 } while (maxResidual > tolerance_ && radIter < maxIter_);
646
647 updateG();
648}
649
650
652{
653 // Construct using contribution from first frequency band
654 auto tRp = volScalarField::New
655 (
656 "Rp",
658 (
659 4
661 * (aLambda_[0] - absorptionEmission_->aDisp(0)())
662 * blackBody_.deltaLambdaT(T_, absorptionEmission_->bands(0))
663 )
664 );
665 auto& Rp = tRp.ref();
666
667 // Add contributions over remaining frequency bands
668 for (label j=1; j < nLambda_; j++)
669 {
670 Rp +=
671 (
672 4
674 *(aLambda_[j] - absorptionEmission_->aDisp(j)())
675 *blackBody_.deltaLambdaT(T_, absorptionEmission_->bands(j))
676 );
678
679 return tRp;
680}
681
682
685{
687 (
688 "Ru",
690 mesh_,
691 dimensionedScalar(dimensionSet(1, -1, -3, 0, 0), Zero)
692 );
693 auto& Ru = tRu.ref();
694
695 // Sum contributions over all frequency bands
696 for (label j=0; j < nLambda_; j++)
697 {
698 // Compute total incident radiation within frequency band
700 (
701 IRay_[0].ILambda(j)()*IRay_[0].omega()
702 );
703
704 for (label rayI=1; rayI < nRay_; rayI++)
705 {
706 Gj.ref() += IRay_[rayI].ILambda(j)()*IRay_[rayI].omega();
707 }
708
709 Ru += (aLambda_[j]() - absorptionEmission_->aDisp(j)()())*Gj
710 - absorptionEmission_->ECont(j)()();
711 }
712
713 return tRu;
714}
715
716
717void Foam::radiation::fvDOM::updateBlackBodyEmission()
718{
719 for (label j=0; j < nLambda_; j++)
720 {
721 blackBody_.correct(j, absorptionEmission_->bands(j));
722 }
723}
724
725
727{
728 G_ = Zero;
729 qr_ = Zero;
730 qem_ = Zero;
731 qin_ = Zero;
732
733 forAll(IRay_, rayI)
734 {
735 IRay_[rayI].addIntensity();
736 G_ += IRay_[rayI].I()*IRay_[rayI].omega();
737 qr_.boundaryFieldRef() += IRay_[rayI].qr().boundaryField();
738 qem_.boundaryFieldRef() += IRay_[rayI].qem().boundaryField();
739 qin_.boundaryFieldRef() += IRay_[rayI].qin().boundaryField();
740 }
741}
742
743
745(
746 const word& name,
747 label& rayId,
748 label& lambdaId
749) const
750{
751 // Assuming name is in the form: CHARS_rayId_lambdaId
752 const auto i1 = name.find('_');
753 const auto i2 = name.find('_', i1+1);
754
755 rayId = readLabel(name.substr(i1+1, i2-i1-1));
756 lambdaId = readLabel(name.substr(i2+1));
757}
758
759
761{
762 return solarCalculator_();
763}
764
765
766// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field....
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
T getOrDefaultCompat(const word &keyword, std::initializer_list< std::pair< const char *, int > > compat, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value using any compatibility names if needed.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition fvDOM.H:117
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant).
Definition fvDOM.C:677
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4).
Definition fvDOM.C:644
const solarCalculator & solarCalc() const
Solar calculator.
Definition fvDOM.C:753
void alignClosestRayToSun(const vector &sunDir)
Align closest ray to sunDir.
Definition fvDOM.C:61
label nBands() const
Number of bands.
Definition fvDOM.H:65
void updateG()
Update G and calculate total heat flux on boundary.
Definition fvDOM.C:719
void rotateInitialRays(const vector &sunDir)
Rotate rays spheric equator to sunDir.
Definition fvDOM.C:47
void updateRaysDir()
Rotate rays according to Sun direction.
Definition fvDOM.C:92
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition fvDOM.C:738
bool read()
Read radiation properties dictionary.
Definition fvDOM.C:557
void calculate()
Solve radiation equation(s).
Definition fvDOM.C:575
Top level model for radiation modelling.
const fvMesh & mesh_
Reference to the mesh database.
virtual bool read()=0
Read radiationProperties dictionary.
dictionary coeffs_
Radiation model dictionary.
autoPtr< absorptionEmissionModel > absorptionEmission_
Absorption/emission model.
const volScalarField & T() const noexcept
Return access to the temperature field.
const volScalarField & T_
Reference to the temperature field.
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:163
A solar calculator model providing models for the solar direction and solar loads.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
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 radiation modelling.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition transform.H:47
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition label.H:63
messageStream Info
Information stream (stdout output on master, null elsewhere).
Tensor< scalar > tensor
Definition symmTensor.H:57
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
errorManip< error > abort(error &err)
Definition errorManip.H:139
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define addToRadiationRunTimeSelectionTables(model)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.