1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2008-2011 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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "physicoChemicalConstants.H"
28 using namespace Foam::constant;
30 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 template<class CloudType>
33 inline const Foam::ThermoCloud<CloudType>&
34 Foam::ThermoCloud<CloudType>::cloudCopy() const
36 return cloudCopyPtr_();
40 template<class CloudType>
41 inline const typename CloudType::particleType::constantProperties&
42 Foam::ThermoCloud<CloudType>::constProps() const
48 template<class CloudType>
49 inline const Foam::SLGThermo& Foam::ThermoCloud<CloudType>::thermo() const
55 template<class CloudType>
56 inline const Foam::volScalarField& Foam::ThermoCloud<CloudType>::T() const
62 template<class CloudType>
63 inline const Foam::volScalarField& Foam::ThermoCloud<CloudType>::p() const
69 template<class CloudType>
70 inline const Foam::HeatTransferModel<Foam::ThermoCloud<CloudType> >&
71 Foam::ThermoCloud<CloudType>::heatTransfer() const
73 return heatTransferModel_;
77 template<class CloudType>
78 inline const Foam::scalarIntegrationScheme&
79 Foam::ThermoCloud<CloudType>::TIntegrator() const
85 template<class CloudType>
86 inline bool Foam::ThermoCloud<CloudType>::radiation() const
92 template<class CloudType>
93 inline Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
94 Foam::ThermoCloud<CloudType>::hsTrans()
100 template<class CloudType>
101 inline const Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
102 Foam::ThermoCloud<CloudType>::hsTrans() const
108 template<class CloudType>
109 inline Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
110 Foam::ThermoCloud<CloudType>::hsCoeff()
116 template<class CloudType>
117 inline const Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
118 Foam::ThermoCloud<CloudType>::hsCoeff() const
124 template<class CloudType>
125 inline Foam::tmp<Foam::fvScalarMatrix>
126 Foam::ThermoCloud<CloudType>::Sh(volScalarField& hs) const
130 Info<< "hsTrans min/max = " << min(hsTrans()).value() << ", "
131 << max(hsTrans()).value() << nl
132 << "hsCoeff min/max = " << min(hsCoeff()).value() << ", "
133 << max(hsCoeff()).value() << endl;
136 if (this->solution().coupled())
138 if (this->solution().semiImplicit("hs"))
140 const volScalarField Cp(thermo_.thermo().Cp());
141 const DimensionedField<scalar, volMesh>
142 Vdt(this->mesh().V()*this->db().time().deltaT());
146 - fvm::SuSp(hsCoeff()/(Cp*Vdt), hs)
147 + hsCoeff()/(Cp*Vdt)*hs;
151 tmp<fvScalarMatrix> tfvm(new fvScalarMatrix(hs, dimEnergy/dimTime));
152 fvScalarMatrix& fvm = tfvm();
154 fvm.source() = -hsTrans()/(this->db().time().deltaT());
160 return tmp<fvScalarMatrix>(new fvScalarMatrix(hs, dimEnergy/dimTime));
164 template<class CloudType>
165 inline Foam::tmp<Foam::volScalarField> Foam::ThermoCloud<CloudType>::Ep() const
167 tmp<volScalarField> tEp
173 this->name() + "radiation::Ep",
174 this->db().time().timeName(),
181 dimensionedScalar("zero", dimMass/dimLength/pow3(dimTime), 0.0)
185 // Need to check if coupled as field is created on-the-fly
186 if (radiation_ && this->solution().coupled())
188 scalarField& Ep = tEp().internalField();
189 const scalarField& V = this->mesh().V();
190 const scalar epsilon = constProps_.epsilon0();
192 forAllConstIter(typename ThermoCloud<CloudType>, *this, iter)
194 const parcelType& p = iter();
195 const label cellI = p.cell();
196 Ep[cellI] += p.nParticle()*p.areaP()*pow4(p.T());
199 Ep *= epsilon*physicoChemical::sigma.value()/V;
206 template<class CloudType>
207 inline Foam::tmp<Foam::volScalarField> Foam::ThermoCloud<CloudType>::ap() const
209 tmp<volScalarField> tap
215 this->name() + "radiation::ap",
216 this->db().time().timeName(),
223 dimensionedScalar("zero", dimless/dimLength, 0.0)
227 // Need to check if coupled as field is created on-the-fly
228 if (radiation_ && this->solution().coupled())
230 scalarField& ap = tap().internalField();
231 const scalarField& V = this->mesh().V();
232 const scalar epsilon = constProps_.epsilon0();
234 forAllConstIter(typename ThermoCloud<CloudType>, *this, iter)
236 const parcelType& p = iter();
237 const label cellI = p.cell();
238 ap[cellI] += p.nParticle()*p.areaP();
248 template<class CloudType>
249 inline Foam::tmp<Foam::volScalarField>
250 Foam::ThermoCloud<CloudType>::sigmap() const
252 tmp<volScalarField> tsigmap
258 this->name() + "radiation::sigmap",
259 this->db().time().timeName(),
266 dimensionedScalar("zero", dimless/dimLength, 0.0)
270 // Need to check if coupled as field is created on-the-fly
271 if (radiation_ && this->solution().coupled())
273 scalarField& sigmap = tsigmap().internalField();
275 const scalarField& V = this->mesh().V();
276 const scalar epsilon = constProps_.epsilon0();
277 const scalar f = constProps_.f0();
279 forAllConstIter(typename ThermoCloud<CloudType>, *this, iter)
281 const parcelType& p = iter();
282 const label cellI = p.cell();
283 sigmap[cellI] += p.nParticle()*p.areaP();
286 sigmap *= (1.0 - f)*(1.0 - epsilon)/V;
293 // ************************************************************************* //