diff --git a/src/thermophysicalModels/radiationModels/radiationModel/fvDOM/fvDOM/fvDOM.C b/src/thermophysicalModels/radiationModels/radiationModel/fvDOM/fvDOM/fvDOM.C index 4c0ad425570e421c7a92b021595f573fe5b715f6..7ff7c41c6f8e4b74d33a5df087553238efdabd69 100644 --- a/src/thermophysicalModels/radiationModels/radiationModel/fvDOM/fvDOM/fvDOM.C +++ b/src/thermophysicalModels/radiationModels/radiationModel/fvDOM/fvDOM/fvDOM.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -49,7 +49,8 @@ namespace Foam void Foam::radiation::fvDOM::initialise() { - if (mesh_.nSolutionD() == 3) //3D + // 3D + if (mesh_.nSolutionD() == 3) { nRay_ = 4*nPhi_*nTheta_; IRay_.setSize(nRay_); @@ -83,69 +84,84 @@ void Foam::radiation::fvDOM::initialise() } } } - else + // 2D + else if (mesh_.nSolutionD() == 2) { - if (mesh_.nSolutionD() == 2) //2D (X & Y) + // Currently 2D solution is limited to the x-y plane + if (mesh_.solutionD()[vector::Z] != -1) { - scalar thetai = piByTwo; - scalar deltaTheta = pi; - nRay_ = 4*nPhi_; - IRay_.setSize(nRay_); - scalar deltaPhi = pi/(2.0*nPhi_); - label i = 0; - for (label m = 1; m <= 4*nPhi_; m++) - { - scalar phii = (2.0*m - 1.0)*deltaPhi/2.0; - IRay_.set + FatalErrorIn("fvDOM::initialise()") + << "Currently 2D solution is limited to the x-y plane" + << exit(FatalError); + } + + scalar thetai = piByTwo; + scalar deltaTheta = pi; + nRay_ = 4*nPhi_; + IRay_.setSize(nRay_); + scalar deltaPhi = pi/(2.0*nPhi_); + label i = 0; + for (label m = 1; m <= 4*nPhi_; m++) + { + scalar phii = (2.0*m - 1.0)*deltaPhi/2.0; + IRay_.set + ( + i, + new radiativeIntensityRay ( - i, - new radiativeIntensityRay - ( - *this, - mesh_, - phii, - thetai, - deltaPhi, - deltaTheta, - nLambda_, - absorptionEmission_, - blackBody_, - i - ) - ); - i++; - } + *this, + mesh_, + phii, + thetai, + deltaPhi, + deltaTheta, + nLambda_, + absorptionEmission_, + blackBody_, + i + ) + ); + i++; } - else //1D (X) + } + // 1D + else + { + // Currently 1D solution is limited to the x-direction + if (mesh_.solutionD()[vector::X] != 1) { - scalar thetai = piByTwo; - scalar deltaTheta = pi; - nRay_ = 2; - IRay_.setSize(nRay_); - scalar deltaPhi = pi; - label i = 0; - for (label m = 1; m <= 2; m++) - { - scalar phii = (2.0*m - 1.0)*deltaPhi/2.0; - IRay_.set + FatalErrorIn("fvDOM::initialise()") + << "Currently 1D solution is limited to the x-direction" + << exit(FatalError); + } + + scalar thetai = piByTwo; + scalar deltaTheta = pi; + nRay_ = 2; + IRay_.setSize(nRay_); + scalar deltaPhi = pi; + label i = 0; + for (label m = 1; m <= 2; m++) + { + scalar phii = (2.0*m - 1.0)*deltaPhi/2.0; + IRay_.set + ( + i, + new radiativeIntensityRay ( - i, - new radiativeIntensityRay - ( - *this, - mesh_, - phii, - thetai, - deltaPhi, - deltaTheta, - nLambda_, - absorptionEmission_, - blackBody_, - i - ) - ); - i++; - } + *this, + mesh_, + phii, + thetai, + deltaPhi, + deltaTheta, + nLambda_, + absorptionEmission_, + blackBody_, + i + ) + ); + i++; } } @@ -388,6 +404,7 @@ Foam::radiation::fvDOM::fvDOM initialise(); } + // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // Foam::radiation::fvDOM::~fvDOM()