initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / thermophysicalModels / radiation / derivedFvPatchFields / wideBandDiffusiveRadiation / wideBandDiffusiveRadiationMixedFvPatchScalarField.C
blob08dbb061ed6e9c8b63b9032c0b7e33ce8e6ee277
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-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 "wideBandDiffusiveRadiationMixedFvPatchScalarField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
32 #include "fvDOM.H"
33 #include "wideBandAbsorptionEmission.H"
34 #include "radiationConstants.H"
35 #include "mathematicalConstants.H"
37 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
39 Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::
40 wideBandDiffusiveRadiationMixedFvPatchScalarField
42     const fvPatch& p,
43     const DimensionedField<scalar, volMesh>& iF
46     mixedFvPatchScalarField(p, iF),
47     TName_("undefinedT"),
48     emissivity_(0.0)
50     refValue() = 0.0;
51     refGrad() = 0.0;
52     valueFraction() = 1.0;
56 Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::
57 wideBandDiffusiveRadiationMixedFvPatchScalarField
59     const wideBandDiffusiveRadiationMixedFvPatchScalarField& ptf,
60     const fvPatch& p,
61     const DimensionedField<scalar, volMesh>& iF,
62     const fvPatchFieldMapper& mapper
65     mixedFvPatchScalarField(ptf, p, iF, mapper),
66     TName_(ptf.TName_),
67     emissivity_(ptf.emissivity_)
71 Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::
72 wideBandDiffusiveRadiationMixedFvPatchScalarField
74     const fvPatch& p,
75     const DimensionedField<scalar, volMesh>& iF,
76     const dictionary& dict
79     mixedFvPatchScalarField(p, iF),
80     TName_(dict.lookup("T")),
81     emissivity_(readScalar(dict.lookup("emissivity")))
83     const scalarField& Tp =
84         patch().lookupPatchField<volScalarField, scalar>(TName_);
86     refValue() =
87         emissivity_*4.0*radiation::sigmaSB.value()*pow4(Tp)
88        /Foam::mathematicalConstant::pi;
89     refGrad() = 0.0;
91     if (dict.found("value"))
92     {
93         fvPatchScalarField::operator=
94         (
95             scalarField("value", dict, p.size())
96         );
97     }
98     else
99     {
100         fvPatchScalarField::operator=(refValue());
101     }
105 Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::
106 wideBandDiffusiveRadiationMixedFvPatchScalarField
108     const wideBandDiffusiveRadiationMixedFvPatchScalarField& ptf
111     mixedFvPatchScalarField(ptf),
112     TName_(ptf.TName_),
113     emissivity_(ptf.emissivity_)
117 Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::
118 wideBandDiffusiveRadiationMixedFvPatchScalarField
120     const wideBandDiffusiveRadiationMixedFvPatchScalarField& ptf,
121     const DimensionedField<scalar, volMesh>& iF
124     mixedFvPatchScalarField(ptf, iF),
125     TName_(ptf.TName_),
126     emissivity_(ptf.emissivity_)
130 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
132 void Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::
133 updateCoeffs()
135     if (this->updated())
136     {
137         return;
138     }
140     const radiationModel& radiation =
141         db().lookupObject<radiationModel>("radiationProperties");
143     const fvDOM& dom(refCast<const fvDOM>(radiation));
145     label rayId = -1;
146     label lambdaId = -1;
147     dom.setRayIdLambdaId(dimensionedInternalField().name(), rayId, lambdaId);
149     const label patchI = patch().index();
151     if (dom.nLambda() == 0)
152     {
153         FatalErrorIn
154         (
155             "Foam::radiation::"
156             "wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs"
157         )   << " a non-grey boundary condition is used with a grey "
158             << "absorption model" << nl << exit(FatalError);
159     }
161     scalarField& Iw = *this;
162     vectorField n = patch().Sf()/patch().magSf();
164     radiativeIntensityRay& ray =
165         const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
167     ray.Qr().boundaryField()[patchI] += Iw*(-n & ray.dAve());
169     const scalarField Eb =
170         dom.blackBody().bLambda(lambdaId).boundaryField()[patchI];
172     forAll(Iw, faceI)
173     {
174         scalar Ir = 0.0;
175         for (label rayI=0; rayI < dom.nRay(); rayI++)
176         {
177             const vector& d = dom.IRay(rayI).d();
179             const scalarField& IFace =
180                 dom.IRay(rayI).ILambda(lambdaId).boundaryField()[patchI];
182             if ((-n[faceI] & d) < 0.0) // qin into the wall
183             {
184                 const vector& dAve = dom.IRay(rayI).dAve();
185                 Ir = Ir + IFace[faceI]*mag(n[faceI] & dAve);
186             }
187         }
189         const vector& d = dom.IRay(rayId).d();
191         if ((-n[faceI] & d) > 0.0)
192         {
193             // direction out of the wall
194             refGrad()[faceI] = 0.0;
195             valueFraction()[faceI] = 1.0;
196             refValue()[faceI] =
197                 (
198                     Ir*(1.0 - emissivity_)
199                   + emissivity_*Eb[faceI]
200                 )
201                /mathematicalConstant::pi;
202         }
203         else
204         {
205             // direction into the wall
206             valueFraction()[faceI] = 0.0;
207             refGrad()[faceI] = 0.0;
208             refValue()[faceI] = 0.0; //not used
209         }
210     }
212     mixedFvPatchScalarField::updateCoeffs();
216 void Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::write
218     Ostream& os
219 ) const
221     fvPatchScalarField::write(os);
222     os.writeKeyword("T") << TName_ << token::END_STATEMENT << nl;
223     os.writeKeyword("emissivity") << emissivity_ << token::END_STATEMENT << nl;
224     writeEntry("value", os);
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 namespace Foam
232 namespace radiation
234     makePatchTypeField
235     (
236         fvPatchScalarField,
237         wideBandDiffusiveRadiationMixedFvPatchScalarField
238     );
243 // ************************************************************************* //