Initial release of the new fireFoam solver and ancillary libraries.
[OpenFOAM-1.6.x.git] / src / thermophysicalModels / radiation / radiationModel / radiationModel / radiationModel.C
blobc56ef02488e9433cbb8d173b05cbe5c477e08185
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 "radiationModel.H"
28 #include "absorptionEmissionModel.H"
29 #include "scatterModel.H"
30 #include "fvm.H"
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 namespace Foam
36     namespace radiation
37     {
38         defineTypeNameAndDebug(radiationModel, 0);
39         defineRunTimeSelectionTable(radiationModel, dictionary);
40     }
44 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
46 Foam::radiation::radiationModel::radiationModel(const volScalarField& T)
48     IOdictionary
49     (
50         IOobject
51         (
52             "radiationProperties",
53             T.time().constant(),
54             T.mesh(),
55             IOobject::MUST_READ,
56             IOobject::NO_WRITE
57         )
58     ),
59     mesh_(T.mesh()),
60     time_(T.time()),
61     T_(T),
62     radiation_(false),
63     coeffs_(dictionary::null),
64     solverFreq_(0),
65     absorptionEmission_(NULL),
66     scatter_(NULL)
70 Foam::radiation::radiationModel::radiationModel
72     const word& type,
73     const volScalarField& T
76     IOdictionary
77     (
78         IOobject
79         (
80             "radiationProperties",
81             T.time().constant(),
82             T.mesh(),
83             IOobject::MUST_READ,
84             IOobject::NO_WRITE
85         )
86     ),
87     mesh_(T.mesh()),
88     time_(T.time()),
89     T_(T),
90     radiation_(lookup("radiation")),
91     coeffs_(subDict(type + "Coeffs")),
92     solverFreq_(readLabel(lookup("solverFreq"))),
93     absorptionEmission_(absorptionEmissionModel::New(*this, mesh_)),
94     scatter_(scatterModel::New(*this, mesh_))
96     solverFreq_ = max(1, solverFreq_);
100 // * * * * * * * * * * * * * * * * Destructor    * * * * * * * * * * * * * * //
102 Foam::radiation::radiationModel::~radiationModel()
106 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
108 bool Foam::radiation::radiationModel::read()
110     if (regIOobject::read())
111     {
112         lookup("radiation") >> radiation_;
113         coeffs_ = subDict(type() + "Coeffs");
115         return true;
116     }
117     else
118     {
119         return false;
120     }
124 void Foam::radiation::radiationModel::correct()
126     if (!radiation_)
127     {
128         return;
129     }
131     if (time_.timeIndex() % solverFreq_ == 0)
132     {
133         calculate();
134     }
138 Foam::tmp<Foam::fvScalarMatrix> Foam::radiation::radiationModel::Sh
140     basicThermo& thermo
141 ) const
143     volScalarField& h = thermo.h();
144     const volScalarField cp = thermo.Cp();
145     const volScalarField T3 = pow3(T_);
147     return
148     (
149         Ru()
150       - fvm::Sp(4.0*Rp()*T3/cp, h)
151       - Rp()*T3*(T_ - 4.0*h/cp)
152     );
155 Foam::tmp<Foam::fvScalarMatrix> Foam::radiation::radiationModel::Shs
157     basicThermo& thermo
158 ) const
160     volScalarField& hs = thermo.hs();
161     const volScalarField cp = thermo.Cp();
162     const volScalarField T3 = pow3(T_);
164     return
165     (
166         Ru()
167       - fvm::Sp(4.0*Rp()*T3/cp, hs)
168       - Rp()*T3*(T_ - 4.0*hs/cp)
169     );
172 // ************************************************************************* //