FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

radiationModel.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "radiationModel.H"
00027 #include <radiation/absorptionEmissionModel.H>
00028 #include <radiation/scatterModel.H>
00029 #include <finiteVolume/fvm.H>
00030 
00031 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00032 
00033 namespace Foam
00034 {
00035     namespace radiation
00036     {
00037         defineTypeNameAndDebug(radiationModel, 0);
00038         defineRunTimeSelectionTable(radiationModel, dictionary);
00039     }
00040 }
00041 
00042 
00043 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00044 
00045 Foam::radiation::radiationModel::radiationModel(const volScalarField& T)
00046 :
00047     IOdictionary
00048     (
00049         IOobject
00050         (
00051             "radiationProperties",
00052             T.time().constant(),
00053             T.mesh(),
00054             IOobject::MUST_READ,
00055             IOobject::NO_WRITE
00056         )
00057     ),
00058     mesh_(T.mesh()),
00059     time_(T.time()),
00060     T_(T),
00061     radiation_(false),
00062     coeffs_(dictionary::null),
00063     solverFreq_(0),
00064     absorptionEmission_(NULL),
00065     scatter_(NULL)
00066 {}
00067 
00068 
00069 Foam::radiation::radiationModel::radiationModel
00070 (
00071     const word& type,
00072     const volScalarField& T
00073 )
00074 :
00075     IOdictionary
00076     (
00077         IOobject
00078         (
00079             "radiationProperties",
00080             T.time().constant(),
00081             T.mesh(),
00082             IOobject::MUST_READ,
00083             IOobject::NO_WRITE
00084         )
00085     ),
00086     mesh_(T.mesh()),
00087     time_(T.time()),
00088     T_(T),
00089     radiation_(lookup("radiation")),
00090     coeffs_(subDict(type + "Coeffs")),
00091     solverFreq_(readLabel(lookup("solverFreq"))),
00092     absorptionEmission_(absorptionEmissionModel::New(*this, mesh_)),
00093     scatter_(scatterModel::New(*this, mesh_))
00094 {
00095     solverFreq_ = max(1, solverFreq_);
00096 }
00097 
00098 
00099 // * * * * * * * * * * * * * * * * Destructor    * * * * * * * * * * * * * * //
00100 
00101 Foam::radiation::radiationModel::~radiationModel()
00102 {}
00103 
00104 
00105 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00106 
00107 bool Foam::radiation::radiationModel::read()
00108 {
00109     if (regIOobject::read())
00110     {
00111         lookup("radiation") >> radiation_;
00112         coeffs_ = subDict(type() + "Coeffs");
00113 
00114         return true;
00115     }
00116     else
00117     {
00118         return false;
00119     }
00120 }
00121 
00122 
00123 void Foam::radiation::radiationModel::correct()
00124 {
00125     if (!radiation_)
00126     {
00127         return;
00128     }
00129 
00130     if (time_.timeIndex() % solverFreq_ == 0)
00131     {
00132         calculate();
00133     }
00134 }
00135 
00136 
00137 Foam::tmp<Foam::fvScalarMatrix> Foam::radiation::radiationModel::Sh
00138 (
00139     basicThermo& thermo
00140 ) const
00141 {
00142     volScalarField& h = thermo.h();
00143     const volScalarField cp = thermo.Cp();
00144     const volScalarField T3 = pow3(T_);
00145 
00146     return
00147     (
00148         Ru()
00149       - fvm::Sp(4.0*Rp()*T3/cp, h)
00150       - Rp()*T3*(T_ - 4.0*h/cp)
00151     );
00152 }
00153 
00154 Foam::tmp<Foam::fvScalarMatrix> Foam::radiation::radiationModel::Shs
00155 (
00156     basicThermo& thermo
00157 ) const
00158 {
00159     volScalarField& hs = thermo.hs();
00160     const volScalarField cp = thermo.Cp();
00161     const volScalarField T3 = pow3(T_);
00162 
00163     return
00164     (
00165         Ru()
00166       - fvm::Sp(4.0*Rp()*T3/cp, hs)
00167       - Rp()*T3*(T_ - 4.0*hs/cp)
00168     );
00169 }
00170 
00171 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines