1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 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 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
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 "radiationConstants.H"
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 template<class ParcelType>
32 inline const typename ParcelType::constantProperties&
33 Foam::ThermoCloud<ParcelType>::constProps() const
39 template<class ParcelType>
40 inline const Foam::basicThermo&
41 Foam::ThermoCloud<ParcelType>::carrierThermo() const
43 return carrierThermo_;
47 template<class ParcelType>
48 inline Foam::basicThermo&
49 Foam::ThermoCloud<ParcelType>::carrierThermo()
51 return carrierThermo_;
55 template<class ParcelType>
56 inline const Foam::HeatTransferModel<Foam::ThermoCloud<ParcelType> >&
57 Foam::ThermoCloud<ParcelType>::heatTransfer() const
59 return heatTransferModel_;
63 template<class ParcelType>
64 inline const Foam::scalarIntegrationScheme&
65 Foam::ThermoCloud<ParcelType>::TIntegrator() const
71 template<class ParcelType>
72 inline bool Foam::ThermoCloud<ParcelType>::radiation() const
78 template<class ParcelType>
79 inline Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
80 Foam::ThermoCloud<ParcelType>::hsTrans()
86 template<class ParcelType>
87 inline Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
88 Foam::ThermoCloud<ParcelType>::Shs() const
90 tmp<DimensionedField<scalar, volMesh> > tShs
92 new DimensionedField<scalar, volMesh>
97 this->db().time().timeName(),
100 IOobject::AUTO_WRITE,
107 dimMass/dimLength/pow3(dimTime),
113 scalarField& Shs = tShs().field();
114 Shs = hsTrans_/(this->mesh().V()*this->db().time().deltaT());
120 template<class ParcelType>
121 inline Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
122 Foam::ThermoCloud<ParcelType>::hcTrans()
128 template<class ParcelType>
129 inline Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
130 Foam::ThermoCloud<ParcelType>::Shc() const
132 tmp<DimensionedField<scalar, volMesh> > tShc
134 new DimensionedField<scalar, volMesh>
138 this->name() + "Shc",
139 this->db().time().timeName(),
142 IOobject::AUTO_WRITE,
149 dimMass/dimLength/pow3(dimTime),
155 scalarField& Shc = tShc().field();
156 Shc = hcTrans_/(this->mesh().V()*this->db().time().deltaT());
162 template<class ParcelType>
163 inline Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
164 Foam::ThermoCloud<ParcelType>::Sh() const
166 tmp<DimensionedField<scalar, volMesh> > tSh
168 new DimensionedField<scalar, volMesh>
173 this->db().time().timeName(),
176 IOobject::AUTO_WRITE,
183 dimMass/dimLength/pow3(dimTime),
189 scalarField& Sh = tSh().field();
190 Sh = (hsTrans_ + hcTrans_)/(this->mesh().V()*this->db().time().deltaT());
196 template<class ParcelType>
197 inline Foam::tmp<Foam::volScalarField>
198 Foam::ThermoCloud<ParcelType>::Ep() const
200 tmp<volScalarField> tEp
206 this->name() + "radiationEp",
207 this->db().time().timeName(),
214 dimensionedScalar("zero", dimMass/dimLength/pow3(dimTime), 0.0)
218 // Need to check if coupled as field is created on-the-fly
219 if (radiation_ && this->coupled())
221 scalarField& Ep = tEp().internalField();
222 const scalarField& V = this->mesh().V();
223 const scalar epsilon = constProps_.epsilon0();
225 forAllConstIter(typename ThermoCloud<ParcelType>, *this, iter)
227 const ParcelType& p = iter();
228 const label cellI = p.cell();
229 Ep[cellI] += p.nParticle()*p.areaP()*pow4(p.T());
232 Ep *= epsilon*radiation::sigmaSB.value()/V;
239 template<class ParcelType>
240 inline Foam::tmp<Foam::volScalarField>
241 Foam::ThermoCloud<ParcelType>::ap() const
243 tmp<volScalarField> tap
249 this->name() + "radiationAp",
250 this->db().time().timeName(),
257 dimensionedScalar("zero", dimless/dimLength, 0.0)
261 // Need to check if coupled as field is created on-the-fly
262 if (radiation_ && this->coupled())
264 scalarField& ap = tap().internalField();
265 const scalarField& V = this->mesh().V();
266 const scalar epsilon = constProps_.epsilon0();
268 forAllConstIter(typename ThermoCloud<ParcelType>, *this, iter)
270 const ParcelType& p = iter();
271 const label cellI = p.cell();
272 ap[cellI] += p.nParticle()*p.areaP();
282 template<class ParcelType>
283 inline Foam::tmp<Foam::volScalarField>
284 Foam::ThermoCloud<ParcelType>::sigmap() const
286 tmp<volScalarField> tsigmap
292 this->name() + "radiationSigmap",
293 this->db().time().timeName(),
300 dimensionedScalar("zero", dimless/dimLength, 0.0)
304 // Need to check if coupled as field is created on-the-fly
305 if (radiation_ && this->coupled())
307 scalarField& sigmap = tsigmap().internalField();
309 const scalarField& V = this->mesh().V();
310 const scalar epsilon = constProps_.epsilon0();
311 const scalar f = constProps_.f0();
313 forAllConstIter(typename ThermoCloud<ParcelType>, *this, iter)
315 const ParcelType& p = iter();
316 const label cellI = p.cell();
317 sigmap[cellI] += p.nParticle()*p.areaP();
320 sigmap *= (1.0 - f)*(1.0 - epsilon)/V;
327 // ************************************************************************* //