Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / lagrangian / intermediate / clouds / Templates / ThermoCloud / ThermoCloudI.H
blobb4e54095a03813abdac9e52d88b5ed2380f341e8
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2011 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
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
19     for more details.
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
44     return constProps_;
48 template<class CloudType>
49 inline const Foam::SLGThermo& Foam::ThermoCloud<CloudType>::thermo() const
51     return thermo_;
55 template<class CloudType>
56 inline const Foam::volScalarField& Foam::ThermoCloud<CloudType>::T() const
58     return T_;
62 template<class CloudType>
63 inline const Foam::volScalarField& Foam::ThermoCloud<CloudType>::p() const
65     return p_;
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
81     return TIntegrator_;
85 template<class CloudType>
86 inline bool Foam::ThermoCloud<CloudType>::radiation() const
88     return radiation_;
92 template<class CloudType>
93 inline Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
94 Foam::ThermoCloud<CloudType>::hsTrans()
96     return hsTrans_();
100 template<class CloudType>
101 inline const Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
102 Foam::ThermoCloud<CloudType>::hsTrans() const
104     return hsTrans_();
108 template<class CloudType>
109 inline Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
110 Foam::ThermoCloud<CloudType>::hsCoeff()
112     return hsCoeff_();
116 template<class CloudType>
117 inline const Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
118 Foam::ThermoCloud<CloudType>::hsCoeff() const
120     return hsCoeff_();
124 template<class CloudType>
125 inline Foam::tmp<Foam::fvScalarMatrix>
126 Foam::ThermoCloud<CloudType>::Sh(volScalarField& hs) const
128     if (debug)
129     {
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;
134     }
136     if (this->solution().coupled())
137     {
138         if (this->solution().semiImplicit("hs"))
139         {
140             const volScalarField Cp(thermo_.thermo().Cp());
141             const DimensionedField<scalar, volMesh>
142                 Vdt(this->mesh().V()*this->db().time().deltaT());
144             return
145                 hsTrans()/Vdt
146               - fvm::SuSp(hsCoeff()/(Cp*Vdt), hs)
147               + hsCoeff()/(Cp*Vdt)*hs;
148         }
149         else
150         {
151             tmp<fvScalarMatrix> tfvm(new fvScalarMatrix(hs, dimEnergy/dimTime));
152             fvScalarMatrix& fvm = tfvm();
154             fvm.source() = -hsTrans()/(this->db().time().deltaT());
156             return tfvm;
157         }
158     }
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
168     (
169         new volScalarField
170         (
171             IOobject
172             (
173                 this->name() + "radiation::Ep",
174                 this->db().time().timeName(),
175                 this->db(),
176                 IOobject::NO_READ,
177                 IOobject::NO_WRITE,
178                 false
179             ),
180             this->mesh(),
181             dimensionedScalar("zero", dimMass/dimLength/pow3(dimTime), 0.0)
182         )
183     );
185     // Need to check if coupled as field is created on-the-fly
186     if (radiation_ && this->solution().coupled())
187     {
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)
193         {
194             const parcelType& p = iter();
195             const label cellI = p.cell();
196             Ep[cellI] += p.nParticle()*p.areaP()*pow4(p.T());
197         }
199         Ep *= epsilon*physicoChemical::sigma.value()/V;
200     }
202     return tEp;
206 template<class CloudType>
207 inline Foam::tmp<Foam::volScalarField> Foam::ThermoCloud<CloudType>::ap() const
209     tmp<volScalarField> tap
210     (
211         new volScalarField
212         (
213             IOobject
214             (
215                 this->name() + "radiation::ap",
216                 this->db().time().timeName(),
217                 this->db(),
218                 IOobject::NO_READ,
219                 IOobject::NO_WRITE,
220                 false
221             ),
222             this->mesh(),
223             dimensionedScalar("zero", dimless/dimLength, 0.0)
224         )
225     );
227     // Need to check if coupled as field is created on-the-fly
228     if (radiation_ && this->solution().coupled())
229     {
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)
235         {
236             const parcelType& p = iter();
237             const label cellI = p.cell();
238             ap[cellI] += p.nParticle()*p.areaP();
239         }
241         ap *= epsilon/V;
242     }
244     return tap;
248 template<class CloudType>
249 inline Foam::tmp<Foam::volScalarField>
250 Foam::ThermoCloud<CloudType>::sigmap() const
252     tmp<volScalarField> tsigmap
253     (
254         new volScalarField
255         (
256             IOobject
257             (
258                 this->name() + "radiation::sigmap",
259                 this->db().time().timeName(),
260                 this->db(),
261                 IOobject::NO_READ,
262                 IOobject::NO_WRITE,
263                 false
264             ),
265             this->mesh(),
266             dimensionedScalar("zero", dimless/dimLength, 0.0)
267         )
268     );
270     // Need to check if coupled as field is created on-the-fly
271     if (radiation_ && this->solution().coupled())
272     {
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)
280         {
281             const parcelType& p = iter();
282             const label cellI = p.cell();
283             sigmap[cellI] += p.nParticle()*p.areaP();
284         }
286         sigmap *= (1.0 - f)*(1.0 - epsilon)/V;
287     }
289     return tsigmap;
293 // ************************************************************************* //