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

greyDiffusiveRadiationMixedFvPatchScalarField.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) 2008-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 "greyDiffusiveRadiationMixedFvPatchScalarField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/volFields.H>
00030 
00031 #include <radiation/fvDOM.H>
00032 #include <radiation/radiationConstants.H>
00033 #include <OpenFOAM/mathematicalConstants.H>
00034 
00035 
00036 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00037 
00038 Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
00039 greyDiffusiveRadiationMixedFvPatchScalarField
00040 (
00041     const fvPatch& p,
00042     const DimensionedField<scalar, volMesh>& iF
00043 )
00044 :
00045     mixedFvPatchScalarField(p, iF),
00046     TName_("undefinedT"),
00047     emissivity_(0.0)
00048 {
00049     refValue() = 0.0;
00050     refGrad() = 0.0;
00051     valueFraction() = 1.0;
00052 }
00053 
00054 
00055 Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
00056 greyDiffusiveRadiationMixedFvPatchScalarField
00057 (
00058     const greyDiffusiveRadiationMixedFvPatchScalarField& ptf,
00059     const fvPatch& p,
00060     const DimensionedField<scalar, volMesh>& iF,
00061     const fvPatchFieldMapper& mapper
00062 )
00063 :
00064     mixedFvPatchScalarField(ptf, p, iF, mapper),
00065     TName_(ptf.TName_),
00066     emissivity_(ptf.emissivity_)
00067 {}
00068 
00069 
00070 Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
00071 greyDiffusiveRadiationMixedFvPatchScalarField
00072 (
00073     const fvPatch& p,
00074     const DimensionedField<scalar, volMesh>& iF,
00075     const dictionary& dict
00076 )
00077 :
00078     mixedFvPatchScalarField(p, iF),
00079     TName_(dict.lookup("T")),
00080     emissivity_(readScalar(dict.lookup("emissivity")))
00081 {
00082     if (dict.found("refValue"))
00083     {
00084         fvPatchScalarField::operator=
00085         (
00086             scalarField("value", dict, p.size())
00087         );
00088         refValue() = scalarField("refValue", dict, p.size());
00089         refGrad() = scalarField("refGradient", dict, p.size());
00090         valueFraction() = scalarField("valueFraction", dict, p.size());
00091     }
00092     else
00093     {
00094         // No value given. Restart as fixedValue b.c.
00095 
00096         const scalarField& Tp =
00097             patch().lookupPatchField<volScalarField, scalar>(TName_);
00098 
00099         refValue() =
00100             emissivity_*4.0*radiation::sigmaSB.value()*pow4(Tp)
00101            /Foam::mathematicalConstant::pi;
00102 
00103         refGrad() = 0.0;
00104         valueFraction() = 1.0;
00105 
00106         fvPatchScalarField::operator=(refValue());
00107     }
00108 }
00109 
00110 
00111 Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
00112 greyDiffusiveRadiationMixedFvPatchScalarField
00113 (
00114     const greyDiffusiveRadiationMixedFvPatchScalarField& ptf
00115 )
00116 :
00117     mixedFvPatchScalarField(ptf),
00118     TName_(ptf.TName_),
00119     emissivity_(ptf.emissivity_)
00120 {}
00121 
00122 
00123 Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
00124 greyDiffusiveRadiationMixedFvPatchScalarField
00125 (
00126     const greyDiffusiveRadiationMixedFvPatchScalarField& ptf,
00127     const DimensionedField<scalar, volMesh>& iF
00128 )
00129 :
00130     mixedFvPatchScalarField(ptf, iF),
00131     TName_(ptf.TName_),
00132     emissivity_(ptf.emissivity_)
00133 {}
00134 
00135 
00136 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00137 
00138 void Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
00139 updateCoeffs()
00140 {
00141     if (this->updated())
00142     {
00143         return;
00144     }
00145 
00146     const scalarField& Tp =
00147         patch().lookupPatchField<volScalarField, scalar>(TName_);
00148 
00149     const radiationModel& radiation =
00150         db().lookupObject<radiationModel>("radiationProperties");
00151 
00152     const fvDOM& dom(refCast<const fvDOM>(radiation));
00153 
00154     label rayId = -1;
00155     label lambdaId = -1;
00156     dom.setRayIdLambdaId(dimensionedInternalField().name(), rayId, lambdaId);
00157 
00158     const label patchI = patch().index();
00159 
00160     if (dom.nLambda() != 1)
00161     {
00162         FatalErrorIn
00163         (
00164             "Foam::radiation::"
00165             "greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs"
00166         )   << " a grey boundary condition is used with a non-grey "
00167             << "absorption model" << nl << exit(FatalError);
00168     }
00169 
00170     scalarField& Iw = *this;
00171     vectorField n = patch().Sf()/patch().magSf();
00172 
00173     radiativeIntensityRay& ray =
00174         const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
00175 
00176     ray.Qr().boundaryField()[patchI] += Iw*(n & ray.dAve());
00177 
00178     forAll(Iw, faceI)
00179     {
00180         scalar Ir = 0.0;
00181 
00182         for (label rayI=0; rayI < dom.nRay(); rayI++)
00183         {
00184             const vector& d = dom.IRay(rayI).d();
00185 
00186             const scalarField& IFace =
00187                 dom.IRay(rayI).ILambda(lambdaId).boundaryField()[patchI];
00188 
00189             if ((-n[faceI] & d) < 0.0)
00190             {
00191                 // q into the wall
00192                 const vector& dAve = dom.IRay(rayI).dAve();
00193                 Ir += IFace[faceI]*mag(n[faceI] & dAve);
00194             }
00195         }
00196 
00197         const vector& d = dom.IRay(rayId).d();
00198 
00199         if ((-n[faceI] & d) > 0.0)
00200         {
00201             // direction out of the wall
00202             refGrad()[faceI] = 0.0;
00203             valueFraction()[faceI] = 1.0;
00204             refValue()[faceI] =
00205                 (
00206                     Ir*(1.0 - emissivity_)
00207                   + emissivity_*radiation::sigmaSB.value()*pow4(Tp[faceI])
00208                 )
00209                /mathematicalConstant::pi;
00210 
00211             // Emmited heat flux from this ray direction
00212             ray.Qem().boundaryField()[patchI][faceI] =
00213                 refValue()[faceI]*(n[faceI] & ray.dAve());
00214 
00215         }
00216         else
00217         {
00218             // direction into the wall
00219             valueFraction()[faceI] = 0.0;
00220             refGrad()[faceI] = 0.0;
00221             refValue()[faceI] = 0.0; //not used
00222 
00223             // Incident heat flux on this ray direction
00224             ray.Qin().boundaryField()[patchI][faceI] =
00225                 Iw[faceI]*(n[faceI] & ray.dAve());
00226         }
00227     }
00228 
00229     mixedFvPatchScalarField::updateCoeffs();
00230 }
00231 
00232 
00233 void Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::write
00234 (
00235     Ostream& os
00236 ) const
00237 {
00238     mixedFvPatchScalarField::write(os);
00239     os.writeKeyword("T") << TName_ << token::END_STATEMENT << nl;
00240     os.writeKeyword("emissivity") << emissivity_ << token::END_STATEMENT << nl;
00241 }
00242 
00243 
00244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00245 
00246 namespace Foam
00247 {
00248 namespace radiation
00249 {
00250     makePatchTypeField
00251     (
00252         fvPatchScalarField,
00253         greyDiffusiveRadiationMixedFvPatchScalarField
00254     );
00255 }
00256 }
00257 
00258 
00259 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines