1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
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>::Sh() const
90 tmp<DimensionedField<scalar, volMesh> > tSh
92 new DimensionedField<scalar, volMesh>
97 this->db().time().timeName(),
100 IOobject::AUTO_WRITE,
103 hsTrans_/(this->mesh().V()*this->db().time().deltaT())
111 template<class ParcelType>
112 inline Foam::tmp<Foam::volScalarField>
113 Foam::ThermoCloud<ParcelType>::Ep() const
115 tmp<volScalarField> tEp
121 this->name() + "radiation::Ep",
122 this->db().time().timeName(),
129 dimensionedScalar("zero", dimMass/dimLength/pow3(dimTime), 0.0)
133 // Need to check if coupled as field is created on-the-fly
134 if (radiation_ && this->coupled())
136 scalarField& Ep = tEp().internalField();
137 const scalarField& V = this->mesh().V();
138 const scalar epsilon = constProps_.epsilon0();
140 forAllConstIter(typename ThermoCloud<ParcelType>, *this, iter)
142 const ParcelType& p = iter();
143 const label cellI = p.cell();
144 Ep[cellI] += p.nParticle()*p.areaP()*pow4(p.T());
147 Ep *= epsilon*radiation::sigmaSB.value()/V;
154 template<class ParcelType>
155 inline Foam::tmp<Foam::volScalarField>
156 Foam::ThermoCloud<ParcelType>::ap() const
158 tmp<volScalarField> tap
164 this->name() + "radiation::ap",
165 this->db().time().timeName(),
172 dimensionedScalar("zero", dimless/dimLength, 0.0)
176 // Need to check if coupled as field is created on-the-fly
177 if (radiation_ && this->coupled())
179 scalarField& ap = tap().internalField();
180 const scalarField& V = this->mesh().V();
181 const scalar epsilon = constProps_.epsilon0();
183 forAllConstIter(typename ThermoCloud<ParcelType>, *this, iter)
185 const ParcelType& p = iter();
186 const label cellI = p.cell();
187 ap[cellI] += p.nParticle()*p.areaP();
197 template<class ParcelType>
198 inline Foam::tmp<Foam::volScalarField>
199 Foam::ThermoCloud<ParcelType>::sigmap() const
201 tmp<volScalarField> tsigmap
207 this->name() + "radiation::sigmap",
208 this->db().time().timeName(),
215 dimensionedScalar("zero", dimless/dimLength, 0.0)
219 // Need to check if coupled as field is created on-the-fly
220 if (radiation_ && this->coupled())
222 scalarField& sigmap = tsigmap().internalField();
224 const scalarField& V = this->mesh().V();
225 const scalar epsilon = constProps_.epsilon0();
226 const scalar f = constProps_.f0();
228 forAllConstIter(typename ThermoCloud<ParcelType>, *this, iter)
230 const ParcelType& p = iter();
231 const label cellI = p.cell();
232 sigmap[cellI] += p.nParticle()*p.areaP();
235 sigmap *= (1.0 - f)*(1.0 - epsilon)/V;
242 // ************************************************************************* //