Commit 60fe4f9f authored by sergio's avatar sergio Committed by Andrew Heather
Browse files

ENH: Adding collimate beam source to fvDOM and multiBandZoneAbsorption

parent 35d77aae
......@@ -41,6 +41,7 @@ submodels/absorptionEmissionModel/greyMeanAbsorptionEmission/greyMeanAbsorptionE
submodels/absorptionEmissionModel/wideBandAbsorptionEmission/wideBandAbsorptionEmission.C
submodels/absorptionEmissionModel/greyMeanSolidAbsorptionEmission/greyMeanSolidAbsorptionEmission.C
submodels/absorptionEmissionModel/multiBandAbsorptionEmission/multiBandAbsorptionEmission.C
submodels/absorptionEmissionModel/multiBandZoneAbsorptionEmission/multiBandZoneAbsorptionEmission.C
submodels/boundaryRadiationProperties/boundaryRadiationProperties.C
submodels/boundaryRadiationProperties/boundaryRadiationPropertiesPatch.C
......
......@@ -33,6 +33,7 @@ License
#include "fvDOM.H"
#include "constants.H"
#include "unitConversion.H"
using namespace Foam::constant;
using namespace Foam::constant::mathematical;
......@@ -98,6 +99,7 @@ greyDiffusiveRadiationMixedFvPatchScalarField
fvPatchScalarField::operator=(refValue());
}
}
......@@ -178,6 +180,13 @@ updateCoeffs()
const scalarField& emissivity = temissivity();
const tmp<scalarField> ttransmissivity
(
boundaryRadiation.transmissivity(patch().index())
);
const scalarField& transmissivity = ttransmissivity();
scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
......@@ -211,6 +220,39 @@ updateCoeffs()
}
}
scalarField Iexternal(this->size(), 0.0);
if (dom.useExternalBeam())
{
const vector sunDir = dom.solarCalc().direction();
const scalar directSolarRad = dom.solarCalc().directSolarRad();
//label nRaysBeam = dom.nRaysBeam();
label SunRayId(-1);
scalar maxSunRay = -GREAT;
// Looking for the ray closest to the Sun direction
for (label rayI=0; rayI < dom.nRay(); rayI++)
{
const vector& iD = dom.IRay(rayI).d();
scalar dir = sunDir & iD;
if (dir > maxSunRay)
{
maxSunRay = dir;
SunRayId = rayI;
}
}
if (rayId == SunRayId)
{
const scalarField nAve(n & dom.IRay(rayId).dAve());
forAll(Iexternal, faceI)
{
Iexternal[faceI] = directSolarRad/mag(dom.IRay(rayId).dAve());
}
}
}
forAll(Iw, faceI)
{
if ((-n[faceI] & myRayId) > 0.0)
......@@ -219,7 +261,8 @@ updateCoeffs()
refGrad()[faceI] = 0.0;
valueFraction()[faceI] = 1.0;
refValue()[faceI] =
(
Iexternal[faceI]*transmissivity[faceI]
+ (
Ir[faceI]*(scalar(1) - emissivity[faceI])
+ emissivity[faceI]*physicoChemical::sigma.value()
* pow4(Tp[faceI])
......
......@@ -78,7 +78,7 @@ wideBandDiffusiveRadiationMixedFvPatchScalarField
:
mixedFvPatchScalarField(p, iF)
{
if (dict.found("value"))
if (dict.found("refValue"))
{
fvPatchScalarField::operator=
(
......@@ -176,9 +176,14 @@ updateCoeffs()
(
boundaryRadiation.emissivity(patch().index(), lambdaId)
);
const scalarField& emissivity = temissivity();
const tmp<scalarField> ttransmissivity
(
boundaryRadiation.transmissivity(patch().index(), lambdaId)
);
const scalarField& transmissivity = ttransmissivity();
scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
......@@ -232,6 +237,41 @@ updateCoeffs()
}
}
scalarField Iexternal(this->size(), 0.0);
if (dom.useExternalBeam())
{
const vector sunDir = dom.solarCalc().direction();
const scalar directSolarRad =
dom.solarCalc().directSolarRad()
*dom.spectralDistribution()[lambdaId];
//label nRaysBeam = dom.nRaysBeam();
label SunRayId(-1);
scalar maxSunRay = -GREAT;
// Looking for the ray closest to the Sun direction
for (label rayI=0; rayI < dom.nRay(); rayI++)
{
const vector& iD = dom.IRay(rayI).d();
scalar dir = sunDir & iD;
if (dir > maxSunRay)
{
maxSunRay = dir;
SunRayId = rayI;
}
}
if (rayId == SunRayId)
{
const scalarField nAve(n & dom.IRay(rayId).dAve());
forAll(Iexternal, faceI)
{
Iexternal[faceI] = directSolarRad/mag(dom.IRay(rayId).dAve());
}
}
}
forAll(Iw, facei)
{
const vector& d = dom.IRay(rayId).d();
......@@ -242,7 +282,8 @@ updateCoeffs()
refGrad()[facei] = 0.0;
valueFraction()[facei] = 1.0;
refValue()[facei] =
(
Iexternal[facei]*transmissivity[facei]
+ (
Ir[facei]*(1.0 - emissivity[facei])
+ emissivity[facei]*Eb[facei]
)/pi;
......
......@@ -29,6 +29,7 @@ License
#include "absorptionEmissionModel.H"
#include "scatterModel.H"
#include "constants.H"
#include "unitConversion.H"
#include "fvm.H"
#include "addToRunTimeSelectionTable.H"
......@@ -49,22 +50,125 @@ namespace Foam
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::radiation::fvDOM::rotateInitialRays(const vector& sunDir)
{
// Rotate Y spherical cordinates to Sun direction.
// Solid angles on the equator are better fit for planar radiation
const tensor coordRot = rotationTensor(vector(0, 1, 0), sunDir);
forAll(IRay_, rayId)
{
IRay_[rayId].dAve() = coordRot & IRay_[rayId].dAve();
IRay_[rayId].d() = coordRot & IRay_[rayId].d();
}
}
void Foam::radiation::fvDOM:: alignClosestRayToSun(const vector& sunDir)
{
label SunRayId(-1);
scalar maxSunRay = -GREAT;
// Looking for the ray closest to the Sun direction
forAll(IRay_, rayId)
{
const vector& iD = IRay_[rayId].d();
scalar dir = sunDir & iD;
if (dir > maxSunRay)
{
maxSunRay = dir;
SunRayId = rayId;
}
}
// Second rotation to align colimated radiation with the closest ray
const tensor coordRot = rotationTensor(IRay_[SunRayId].d(), sunDir);
forAll(IRay_, rayId)
{
IRay_[rayId].dAve() = coordRot & IRay_[rayId].dAve();
IRay_[rayId].d() = coordRot & IRay_[rayId].d();
}
Info << "Sun direction : " << sunDir << nl << endl;
Info << "Sun ray ID : " << SunRayId << nl << endl;
}
void Foam::radiation::fvDOM::updateRaysDir()
{
solarCalculator_->correctSunDirection();
const vector sunDir = solarCalculator_->direction();
// First iteration
if (updateTimeIndex_ == 0)
{
rotateInitialRays(sunDir);
alignClosestRayToSun(sunDir);
}
else if (updateTimeIndex_ > 0)
{
alignClosestRayToSun(sunDir);
}
}
void Foam::radiation::fvDOM::initialise()
{
coeffs_.readIfPresent("useExternalBeam", useExternalBeam_);
if (useExternalBeam_)
{
coeffs_.readEntry("spectralDistribution", spectralDistribution_);
spectralDistribution_ =
spectralDistribution_/sum(spectralDistribution_);
const dictionary& solarDict = this->subDict("solarCalculatorCoeffs");
solarCalculator_.reset(new solarCalculator(solarDict, mesh_));
if (mesh_.nSolutionD() != 3)
{
FatalErrorInFunction
<< "External beam model only available in 3D meshes "
<< abort(FatalError);
}
if (solarCalculator_->diffuseSolarRad() > 0)
{
FatalErrorInFunction
<< "External beam model does not support Diffuse "
<< "Solar Radiation. Set diffuseSolarRad to zero"
<< abort(FatalError);
}
if (spectralDistribution_.size() != nLambda_)
{
FatalErrorInFunction
<< "The epectral energy distribution has different bands "
<< "than the absoprtivity model "
<< abort(FatalError);
}
}
// 3D
if (mesh_.nSolutionD() == 3)
{
nRay_ = 4*nPhi_*nTheta_;
IRay_.setSize(nRay_);
const scalar deltaPhi = pi/(2.0*nPhi_);
const scalar deltaPhi = pi/(2*nPhi_);
const scalar deltaTheta = pi/nTheta_;
label i = 0;
for (label n = 1; n <= nTheta_; n++)
{
for (label m = 1; m <= 4*nPhi_; m++)
{
const scalar thetai = (2*n - 1)*deltaTheta/2.0;
const scalar phii = (2*m - 1)*deltaPhi/2.0;
scalar thetai = (2*n - 1)*deltaTheta/2.0;
scalar phii = (2*m - 1)*deltaPhi/2.0;
IRay_.set
(
i,
......@@ -176,22 +280,41 @@ void Foam::radiation::fvDOM::initialise()
Info<< "fvDOM : Allocated " << IRay_.size()
<< " rays with average orientation:" << nl;
if (useExternalBeam_)
{
// Rotate rays for Sun direction
updateRaysDir();
}
scalar totalOmega = 0;
forAll(IRay_, rayId)
{
if (omegaMax_ < IRay_[rayId].omega())
{
omegaMax_ = IRay_[rayId].omega();
}
totalOmega += IRay_[rayId].omega();
Info<< '\t' << IRay_[rayId].I().name() << " : " << "dAve : "
<< '\t' << IRay_[rayId].dAve() << nl;
<< '\t' << IRay_[rayId].dAve() << " : " << "omega : "
<< '\t' << IRay_[rayId].omega() << " : " << "d : "
<< '\t' << IRay_[rayId].d() << nl;
}
Info << "Total omega : " << totalOmega << endl;
Info<< endl;
coeffs_.readIfPresent("useSolarLoad", useSolarLoad_);
if (useSolarLoad_)
{
if (useExternalBeam_)
{
FatalErrorInFunction
<< "External beam with fvDOM can not be used "
<< "with the solar load model"
<< abort(FatalError);
}
const dictionary& solarDict = this->subDict("solarLoadCoeffs");
solarLoad_.reset(new solarLoad(solarDict, T_));
......@@ -296,7 +419,11 @@ Foam::radiation::fvDOM::fvDOM(const volScalarField& T)
meshOrientation_
(
coeffs_.lookupOrDefault<vector>("meshOrientation", Zero)
)
),
useExternalBeam_(false),
spectralDistribution_(),
solarCalculator_(),
updateTimeIndex_(0)
{
initialise();
}
......@@ -392,7 +519,11 @@ Foam::radiation::fvDOM::fvDOM
meshOrientation_
(
coeffs_.lookupOrDefault<vector>("meshOrientation", Zero)
)
),
useExternalBeam_(false),
spectralDistribution_(),
solarCalculator_(),
updateTimeIndex_(0)
{
initialise();
}
......@@ -435,6 +566,33 @@ void Foam::radiation::fvDOM::calculate()
solarLoad_->calculate();
}
if (useExternalBeam_)
{
switch (solarCalculator_->sunDirectionModel())
{
case solarCalculator::mSunDirConstant:
{
break;
}
case solarCalculator::mSunDirTracking:
{
label updateIndex = label
(
mesh_.time().value()
/solarCalculator_->sunTrackingUpdateInterval()
);
if (updateIndex > updateTimeIndex_)
{
Info << "Updating Sun position..." << endl;
updateTimeIndex_ = updateIndex;
updateRaysDir();
}
break;
}
}
}
// Set rays convergence false
List<bool> rayIdConv(nRay_, false);
......@@ -597,4 +755,10 @@ void Foam::radiation::fvDOM::setRayIdLambdaId
}
const Foam::solarCalculator& Foam::radiation::fvDOM::solarCalc() const
{
return solarCalculator_();
}
// ************************************************************************* //
......@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016 OpenCFD Ltd.
\\ / A nd | Copyright (C) 2019 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
| Copyright (C) 2011-2017 OpenFOAM Foundation
......@@ -32,12 +32,19 @@ Group
Description
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n
directions in a participating media, not including scatter.
directions in a participating media, not including scatter and reflective
walls.
Available absorption models:
- constantAbsorptionEmission
- greyMeanAbsoprtionEmission
- wideBandAbsorptionEmission
- multiBandAbsorptionEmission
This model can handle non-grey participating media using
multiBandAbsorptionEmission model. Accordingly the BC for rays should
be wideBandDiffussive type
Usage
\verbatim
......@@ -50,6 +57,10 @@ Usage
// iteration
maxIter 4; // maximum number of iterations
meshOrientation (1 1 1); //Mesh orientation used for 2D and 1D
useSolarLoad false;
useExternalBeam true;
spectralDistribution (2 1);
}
solverFreq 1; // Number of flow iterations per radiation iteration
......@@ -71,6 +82,16 @@ Usage
- rays geberated in 3-D using the \c nPhi and \c nTheta entries
- \c meshOrientation vector is not applicable.
useSolarLoad calculates the primary and diffusive Sun fluxes on walls in
addition to the RTE equations
useExternalBeam add an external collimated beam to the domain. This option
is not available if useSolarLoad is true.
spectralDistribution is the energy spectral distribution of the collimated
external beam.
SourceFiles
fvDOM.C
......@@ -156,6 +177,18 @@ class fvDOM
//- Mesh orientation vector
vector meshOrientation_;
//- Use external parallel irradiation beam
bool useExternalBeam_;
//- Spectral energy distribution for the external beam
scalarList spectralDistribution_;
//- Solar calculator
autoPtr<solarCalculator> solarCalculator_;
//- Update Sun position index
label updateTimeIndex_;
// Private Member Functions
......@@ -213,6 +246,16 @@ public:
label& lambdaId
) const;
//- Rotate rays according to Sun direction
void updateRaysDir();
//- Rotate rays spheric equator to sunDir
void rotateInitialRays(const vector& sunDir);
//- Align closest ray to sunDir
void alignClosestRayToSun(const vector& sunDir);
//- Source term component (for power of T^4)
virtual tmp<volScalarField> Rp() const;
......@@ -222,6 +265,9 @@ public:
// Access
//- Solar calculator
const solarCalculator& solarCalc() const;
//- Ray intensity for rayI
inline const radiativeIntensityRay& IRay(const label rayI) const;
......@@ -276,6 +322,12 @@ public:
//- Use solar load
inline bool useSolarLoad() const;
//- Use external beam
inline bool useExternalBeam() const;
//- Energy spectral distribution for external beam
inline const scalarList& spectralDistribution() const;
};
......
......@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd |
\\ / A nd | Copyright (C) 2008-2011 2019 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
| Copyright (C) 2011-2017 OpenFOAM Foundation
......@@ -134,4 +134,18 @@ inline bool Foam::radiation::fvDOM::useSolarLoad() const
return useSolarLoad_;
}
inline bool Foam::radiation::fvDOM::useExternalBeam() const
{
return useExternalBeam_;
}
inline const Foam::scalarList& Foam::radiation::fvDOM::
spectralDistribution() const
{
return spectralDistribution_;
}
// ************************************************************************* //
......@@ -201,6 +201,12 @@ public:
//- Return the average vector inside the solid angle
inline const vector& dAve() const;
//- Return direction
inline vector& d();
//- Return the average vector inside the solid angle
inline vector& dAve();
//- Return the number of bands
inline scalar nLambda() const;
......
......@@ -82,6 +82,18 @@ inline const Foam::vector& Foam::radiation::radiativeIntensityRay::dAve() const
}
inline Foam::vector& Foam::radiation::radiativeIntensityRay::d()
{
return d_;
}
inline Foam::vector& Foam::radiation::radiativeIntensityRay::dAve()
{
return dAve_;
}
inline Foam::scalar Foam<