FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

radiativeIntensityRay.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 2008-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "radiativeIntensityRay.H"
00027 #include <finiteVolume/fvm.H>
00028 #include <radiation/fvDOM.H>
00029 
00030 const Foam::word
00031 Foam::radiation::radiativeIntensityRay::intensityPrefix("ILambda");
00032 
00033 
00034 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00035 
00036 Foam::radiation::radiativeIntensityRay::radiativeIntensityRay
00037 (
00038     const fvDOM& dom,
00039     const fvMesh& mesh,
00040     const scalar phi,
00041     const scalar theta,
00042     const scalar deltaPhi,
00043     const scalar deltaTheta,
00044     const label nLambda,
00045     const absorptionEmissionModel& absorptionEmission,
00046     const blackBodyEmission& blackBody,
00047     const label rayId
00048 )
00049 :
00050     dom_(dom),
00051     mesh_(mesh),
00052     absorptionEmission_(absorptionEmission),
00053     blackBody_(blackBody),
00054     I_
00055     (
00056         IOobject
00057         (
00058             "I" + name(rayId),
00059             mesh_.time().timeName(),
00060             mesh_,
00061             IOobject::NO_READ,
00062             IOobject::NO_WRITE
00063         ),
00064         mesh_,
00065         dimensionedScalar("I", dimMass/pow3(dimTime), 0.0)
00066     ),
00067     Qr_
00068     (
00069         IOobject
00070         (
00071             "Qr" + name(rayId),
00072             mesh_.time().timeName(),
00073             mesh_,
00074             IOobject::NO_READ,
00075             IOobject::NO_WRITE
00076         ),
00077         mesh_,
00078         dimensionedScalar("Qr", dimMass/pow3(dimTime), 0.0)
00079     ),
00080     Qin_
00081     (
00082         IOobject
00083         (
00084             "Qin" + name(rayId),
00085             mesh_.time().timeName(),
00086             mesh_,
00087             IOobject::NO_READ,
00088             IOobject::NO_WRITE
00089         ),
00090         mesh_,
00091         dimensionedScalar("Qin", dimMass/pow3(dimTime), 0.0)
00092     ),
00093     Qem_
00094     (
00095         IOobject
00096         (
00097             "Qem" + name(rayId),
00098             mesh_.time().timeName(),
00099             mesh_,
00100             IOobject::NO_READ,
00101             IOobject::NO_WRITE
00102         ),
00103         mesh_,
00104         dimensionedScalar("Qem", dimMass/pow3(dimTime), 0.0)
00105     ),
00106     d_(vector::zero),
00107     dAve_(vector::zero),
00108     theta_(theta),
00109     phi_(phi),
00110     omega_(0.0),
00111     nLambda_(nLambda),
00112     ILambda_(nLambda)
00113 {
00114     scalar sinTheta = Foam::sin(theta);
00115     scalar cosTheta = Foam::cos(theta);
00116     scalar sinPhi = Foam::sin(phi);
00117     scalar cosPhi = Foam::cos(phi);
00118 
00119     omega_ = 2.0*sinTheta*Foam::sin(deltaTheta/2.0)*deltaPhi;
00120     d_ = vector(sinTheta*sinPhi, sinTheta*cosPhi, cosTheta);
00121     dAve_ = vector
00122     (
00123         sinPhi
00124        *Foam::sin(0.5*deltaPhi)
00125        *(deltaTheta - Foam::cos(2.0*theta)
00126        *Foam::sin(deltaTheta)),
00127         cosPhi
00128        *Foam::sin(0.5*deltaPhi)
00129        *(deltaTheta - Foam::cos(2.0*theta)
00130        *Foam::sin(deltaTheta)),
00131         0.5*deltaPhi*Foam::sin(2.0*theta)*Foam::sin(deltaTheta)
00132     );
00133 
00134 
00135     autoPtr<volScalarField> IDefaultPtr;
00136 
00137     forAll(ILambda_, lambdaI)
00138     {
00139         IOobject IHeader
00140         (
00141             intensityPrefix + "_" + name(rayId) + "_" + name(lambdaI),
00142             mesh_.time().timeName(),
00143             mesh_,
00144             IOobject::MUST_READ,
00145             IOobject::AUTO_WRITE
00146         );
00147 
00148         // check if field exists and can be read
00149         if (IHeader.headerOk())
00150         {
00151             ILambda_.set
00152             (
00153                 lambdaI,
00154                 new volScalarField(IHeader, mesh_)
00155             );
00156         }
00157         else
00158         {
00159             // Demand driven load the IDefault field
00160             if (!IDefaultPtr.valid())
00161             {
00162                 IDefaultPtr.reset
00163                 (
00164                     new volScalarField
00165                     (
00166                         IOobject
00167                         (
00168                             "IDefault",
00169                             mesh_.time().timeName(),
00170                             mesh_,
00171                             IOobject::MUST_READ,
00172                             IOobject::NO_WRITE
00173                         ),
00174                         mesh_
00175                     )
00176                 );
00177             }
00178 
00179             // Reset the MUST_READ flag
00180             IOobject noReadHeader(IHeader);
00181             noReadHeader.readOpt() = IOobject::NO_READ;
00182 
00183             ILambda_.set
00184             (
00185                 lambdaI,
00186                 new volScalarField(noReadHeader, IDefaultPtr())
00187             );
00188         }
00189     }
00190 }
00191 
00192 
00193 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00194 
00195 Foam::radiation::radiativeIntensityRay::~radiativeIntensityRay()
00196 {}
00197 
00198 
00199 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00200 
00201 Foam::scalar Foam::radiation::radiativeIntensityRay::correct()
00202 {
00203     // reset boundary heat flux to zero
00204     Qr_.boundaryField() = 0.0;
00205 
00206     scalar maxResidual = -GREAT;
00207 
00208     forAll(ILambda_, lambdaI)
00209     {
00210         const volScalarField& k = dom_.aLambda(lambdaI);
00211 
00212         surfaceScalarField Ji = dAve_ & mesh_.Sf();
00213 
00214         fvScalarMatrix IiEq
00215         (
00216             fvm::div(Ji, ILambda_[lambdaI], "div(Ji,Ii_h)")
00217           + fvm::Sp(k*omega_, ILambda_[lambdaI])
00218          ==
00219             1.0/Foam::mathematicalConstant::pi*omega_
00220            *(
00221                 k*blackBody_.bLambda(lambdaI)
00222               + absorptionEmission_.ECont(lambdaI)/4
00223             )
00224         );
00225 
00226         IiEq.relax();
00227 
00228         scalar eqnResidual = solve
00229         (
00230             IiEq,
00231             mesh_.solver("Ii")
00232         ).initialResidual();
00233 
00234         maxResidual = max(eqnResidual, maxResidual);
00235     }
00236 
00237     return maxResidual;
00238 }
00239 
00240 
00241 void Foam::radiation::radiativeIntensityRay::addIntensity()
00242 {
00243     I_ = dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0);
00244 
00245     forAll(ILambda_, lambdaI)
00246     {
00247         I_ += absorptionEmission_.addIntensity(lambdaI, ILambda_[lambdaI]);
00248     }
00249 }
00250 
00251 
00252 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines