initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / thermophysicalModels / radiation / radiationModel / P1 / P1.C
blobe11c73be0a5ab69f66e279229a4bf1b3a79e703b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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 "P1.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvm.H"
31 #include "absorptionEmissionModel.H"
32 #include "scatterModel.H"
33 #include "mathematicalConstants.H"
34 #include "radiationConstants.H"
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 namespace Foam
40     namespace radiation
41     {
42         defineTypeNameAndDebug(P1, 0);
44         addToRunTimeSelectionTable
45         (
46             radiationModel,
47             P1,
48             dictionary
49         );
50     }
53 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
55 // Construct from components
56 Foam::radiation::P1::P1(const volScalarField& T)
58     radiationModel(typeName, T),
59     G_
60     (
61         IOobject
62         (
63             "G",
64             mesh_.time().timeName(),
65             mesh_,
66             IOobject::MUST_READ,
67             IOobject::AUTO_WRITE
68         ),
69         mesh_
70     ),
71     a_
72     (
73         IOobject
74         (
75             "a",
76             mesh_.time().timeName(),
77             mesh_,
78             IOobject::NO_READ,
79             IOobject::NO_WRITE
80         ),
81         mesh_,
82         dimensionedScalar("a", dimless/dimLength, 0.0)
83     ),
84     e_
85     (
86         IOobject
87         (
88             "e",
89             mesh_.time().timeName(),
90             mesh_,
91             IOobject::NO_READ,
92             IOobject::NO_WRITE
93         ),
94         mesh_,
95         dimensionedScalar("a", dimless/dimLength, 0.0)
96     ),
97     E_
98     (
99         IOobject
100         (
101             "E",
102             mesh_.time().timeName(),
103             mesh_,
104             IOobject::NO_READ,
105             IOobject::NO_WRITE
106         ),
107         mesh_,
108         dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
109     )
113 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
115 Foam::radiation::P1::~P1()
119 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
121 bool Foam::radiation::P1::read()
123     if (radiationModel::read())
124     {
125         // nothing to read
127         return true;
128     }
129     else
130     {
131         return false;
132     }
136 void Foam::radiation::P1::correct()
138     if (!radiation_)
139     {
140         return;
141     }
142     a_ = absorptionEmission_->a();
143     e_ = absorptionEmission_->e();
144     E_ = absorptionEmission_->E();
145     const volScalarField sigmaEff = scatter_->sigmaEff();
147     // Construct diffusion
148     const volScalarField gamma
149     (
150         IOobject
151         (
152             "gammaRad",
153             G_.mesh().time().timeName(),
154             G_.mesh(),
155             IOobject::NO_READ,
156             IOobject::NO_WRITE
157         ),
158         1.0/(3.0*a_ + sigmaEff)
159     );
161     // Solve G transport equation
162     solve
163     (
164         fvm::laplacian(gamma, G_)
165       - fvm::Sp(a_, G_)
166      ==
167       - 4.0*(e_*radiation::sigmaSB*pow4(T_) + E_)
168     );
172 Foam::tmp<Foam::volScalarField> Foam::radiation::P1::Rp() const
174     return tmp<volScalarField>
175     (
176         new volScalarField
177         (
178             IOobject
179             (
180                 "Rp",
181                 mesh_.time().timeName(),
182                 mesh_,
183                 IOobject::NO_READ,
184                 IOobject::NO_WRITE,
185                 false
186             ),
187             4.0*absorptionEmission_->eCont()*radiation::sigmaSB
188         )
189     );
193 Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
194 Foam::radiation::P1::Ru() const
196     const DimensionedField<scalar, volMesh>& G =
197         G_.dimensionedInternalField();
198     const DimensionedField<scalar, volMesh> E =
199         absorptionEmission_->ECont()().dimensionedInternalField();
200     const DimensionedField<scalar, volMesh> a =
201         absorptionEmission_->aCont()().dimensionedInternalField();
203     return  a*G - 4.0*E;
207 // ************************************************************************* //