1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2008-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software; you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "blackBodyEmission.H"
28 #include "dimensionedConstants.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 const Foam::List<Foam::Tuple2<Foam::scalar, Foam::scalar> >
33 Foam::radiation::blackBodyEmission::emissivePowerTable
140 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 Foam::radiation::blackBodyEmission::blackBodyEmission
145 const volScalarField& T
151 interpolationTable<scalar>::CLAMP,
152 "blackBodyEmissivePower"
154 C1_("C1", dimensionSet(1, 4, 3, 0, 0, 0, 0), 3.7419e-16),
155 C2_("C2", dimensionSet(0, 1, 0, 1, 0, 0, 0), 14.388e-6),
159 forAll(bLambda_, lambdaI)
168 "bLambda_" + Foam::name(lambdaI) ,
169 T.mesh().time().timeName(),
174 radiation::sigmaSB*pow4(T)
182 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
184 Foam::radiation::blackBodyEmission::~blackBodyEmission()
188 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
190 Foam::scalar Foam::radiation::blackBodyEmission::fLambdaT
195 return table_(lambdaT*1.0e6);
199 Foam::tmp<Foam::volScalarField>
200 Foam::radiation::blackBodyEmission::EbDeltaLambdaT
202 const volScalarField& T,
203 const Vector2D<scalar>& band
206 tmp<volScalarField> Eb
213 T.mesh().time().timeName(),
218 radiation::sigmaSB*pow4(T)
223 if (band == Vector2D<scalar>::one)
231 scalar T1 = fLambdaT(band[1]*T[i]);
232 scalar T2 = fLambdaT(band[0]*T[i]);
233 dimensionedScalar fLambdaDelta
239 Eb()[i] = Eb()[i]*fLambdaDelta.value();
246 void Foam::radiation::blackBodyEmission::correct
249 const Vector2D<scalar>& band
252 bLambda_[lambdaI] = EbDeltaLambdaT(T_, band);
256 // ************************************************************************* //