1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-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 \*---------------------------------------------------------------------------*/
28 #include "addToRunTimeSelectionTable.H"
31 #include "absorptionEmissionModel.H"
32 #include "scatterModel.H"
33 #include "mathematicalConstants.H"
34 #include "radiationConstants.H"
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 defineTypeNameAndDebug(P1, 0);
44 addToRunTimeSelectionTable
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 Foam::radiation::P1::P1(const volScalarField& T)
58 radiationModel(typeName, T),
64 mesh_.time().timeName(),
76 mesh_.time().timeName(),
82 dimensionedScalar("a", dimless/dimLength, 0.0)
89 mesh_.time().timeName(),
95 dimensionedScalar("a", dimless/dimLength, 0.0)
102 mesh_.time().timeName(),
108 dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
113 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
115 Foam::radiation::P1::~P1()
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 bool Foam::radiation::P1::read()
123 if (radiationModel::read())
136 void Foam::radiation::P1::calculate()
138 a_ = absorptionEmission_->a();
139 e_ = absorptionEmission_->e();
140 E_ = absorptionEmission_->E();
141 const volScalarField sigmaEff = scatter_->sigmaEff();
143 // Construct diffusion
144 const volScalarField gamma
149 G_.mesh().time().timeName(),
154 1.0/(3.0*a_ + sigmaEff)
157 // Solve G transport equation
160 fvm::laplacian(gamma, G_)
163 - 4.0*(e_*radiation::sigmaSB*pow4(T_) + E_)
168 Foam::tmp<Foam::volScalarField> Foam::radiation::P1::Rp() const
170 return tmp<volScalarField>
177 mesh_.time().timeName(),
183 4.0*absorptionEmission_->eCont()*radiation::sigmaSB
189 Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
190 Foam::radiation::P1::Ru() const
192 const DimensionedField<scalar, volMesh>& G =
193 G_.dimensionedInternalField();
194 const DimensionedField<scalar, volMesh> E =
195 absorptionEmission_->ECont()().dimensionedInternalField();
196 const DimensionedField<scalar, volMesh> a =
197 absorptionEmission_->aCont()().dimensionedInternalField();
203 // ************************************************************************* //