Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / thermophysicalModels / radiationModels / radiationModel / radiationModel / radiationModel.C
blobceeb9aeed498e6c2ee35a3921a40adcb151e4a6f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2011 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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "radiationModel.H"
27 #include "absorptionEmissionModel.H"
28 #include "scatterModel.H"
29 #include "fvmSup.H"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 namespace Foam
35     namespace radiation
36     {
37         defineTypeNameAndDebug(radiationModel, 0);
38         defineRunTimeSelectionTable(radiationModel, dictionary);
39     }
43 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
45 Foam::radiation::radiationModel::radiationModel(const volScalarField& T)
47     IOdictionary
48     (
49         IOobject
50         (
51             "radiationProperties",
52             T.time().constant(),
53             T.mesh(),
54             IOobject::MUST_READ_IF_MODIFIED,
55             IOobject::NO_WRITE
56         )
57     ),
58     mesh_(T.mesh()),
59     time_(T.time()),
60     T_(T),
61     radiation_(false),
62     coeffs_(dictionary::null),
63     solverFreq_(0),
64     absorptionEmission_(NULL),
65     scatter_(NULL)
69 Foam::radiation::radiationModel::radiationModel
71     const word& type,
72     const volScalarField& T
75     IOdictionary
76     (
77         IOobject
78         (
79             "radiationProperties",
80             T.time().constant(),
81             T.mesh(),
82             IOobject::MUST_READ_IF_MODIFIED,
83             IOobject::NO_WRITE
84         )
85     ),
86     mesh_(T.mesh()),
87     time_(T.time()),
88     T_(T),
89     radiation_(lookup("radiation")),
90     coeffs_(subDict(type + "Coeffs")),
91     solverFreq_(readLabel(lookup("solverFreq"))),
92     absorptionEmission_(absorptionEmissionModel::New(*this, mesh_)),
93     scatter_(scatterModel::New(*this, mesh_))
95     solverFreq_ = max(1, solverFreq_);
99 // * * * * * * * * * * * * * * * * Destructor    * * * * * * * * * * * * * * //
101 Foam::radiation::radiationModel::~radiationModel()
105 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
107 bool Foam::radiation::radiationModel::read()
109     if (regIOobject::read())
110     {
111         lookup("radiation") >> radiation_;
112         coeffs_ = subDict(type() + "Coeffs");
114         lookup("solverFreq") >> solverFreq_,
115         solverFreq_ = max(1, solverFreq_);
117         return true;
118     }
119     else
120     {
121         return false;
122     }
126 void Foam::radiation::radiationModel::correct()
128     if (!radiation_)
129     {
130         return;
131     }
133     if (time_.timeIndex() % solverFreq_ == 0)
134     {
135         calculate();
136     }
140 Foam::tmp<Foam::fvScalarMatrix> Foam::radiation::radiationModel::Sh
142     basicThermo& thermo
143 ) const
145     volScalarField& h = thermo.h();
146     const volScalarField Cp(thermo.Cp());
147     const volScalarField T3(pow3(T_));
149     return
150     (
151         Ru()
152       - fvm::Sp(4.0*Rp()*T3/Cp, h)
153       - Rp()*T3*(T_ - 4.0*h/Cp)
154     );
158 Foam::tmp<Foam::fvScalarMatrix> Foam::radiation::radiationModel::Shs
160     basicThermo& thermo
161 ) const
163     volScalarField& hs = thermo.hs();
164     const volScalarField Cp(thermo.Cp());
165     const volScalarField T3(pow3(T_));
167     return
168     (
169         Ru()
170       - fvm::Sp(4.0*Rp()*T3/Cp, hs)
171       - Rp()*T3*(T_ - 4.0*hs/Cp)
172     );
176 // ************************************************************************* //