1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2008-2011 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
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
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"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 defineTypeNameAndDebug(radiationModel, 0);
38 defineRunTimeSelectionTable(radiationModel, dictionary);
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 Foam::radiation::radiationModel::radiationModel(const volScalarField& T)
51 "radiationProperties",
54 IOobject::MUST_READ_IF_MODIFIED,
62 coeffs_(dictionary::null),
64 absorptionEmission_(NULL),
69 Foam::radiation::radiationModel::radiationModel
72 const volScalarField& T
79 "radiationProperties",
82 IOobject::MUST_READ_IF_MODIFIED,
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())
111 lookup("radiation") >> radiation_;
112 coeffs_ = subDict(type() + "Coeffs");
114 lookup("solverFreq") >> solverFreq_,
115 solverFreq_ = max(1, solverFreq_);
126 void Foam::radiation::radiationModel::correct()
133 if (time_.timeIndex() % solverFreq_ == 0)
140 Foam::tmp<Foam::fvScalarMatrix> Foam::radiation::radiationModel::Sh
145 volScalarField& h = thermo.h();
146 const volScalarField Cp(thermo.Cp());
147 const volScalarField T3(pow3(T_));
152 - fvm::Sp(4.0*Rp()*T3/Cp, h)
153 - Rp()*T3*(T_ - 4.0*h/Cp)
158 Foam::tmp<Foam::fvScalarMatrix> Foam::radiation::radiationModel::Shs
163 volScalarField& hs = thermo.hs();
164 const volScalarField Cp(thermo.Cp());
165 const volScalarField T3(pow3(T_));
170 - fvm::Sp(4.0*Rp()*T3/Cp, hs)
171 - Rp()*T3*(T_ - 4.0*hs/Cp)
176 // ************************************************************************* //