Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / lagrangian / intermediate / clouds / Templates / ThermoCloud / ThermoCloudTemplateI.H
blob2f28dde1c1147f6bf62bf1b1ff7ec117fbae724d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "radiationConstants.H"
28 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
30 template<class ParcelType>
31 inline const typename ParcelType::constantProperties&
32 Foam::ThermoCloud<ParcelType>::constProps() const
34     return constProps_;
38 template<class ParcelType>
39 inline const Foam::basicThermo&
40 Foam::ThermoCloud<ParcelType>::carrierThermo() const
42     return carrierThermo_;
46 template<class ParcelType>
47 inline Foam::basicThermo&
48 Foam::ThermoCloud<ParcelType>::carrierThermo()
50     return carrierThermo_;
54 template<class ParcelType>
55 inline const Foam::HeatTransferModel<Foam::ThermoCloud<ParcelType> >&
56 Foam::ThermoCloud<ParcelType>::heatTransfer() const
58     return heatTransferModel_;
62 template<class ParcelType>
63 inline const Foam::scalarIntegrationScheme&
64 Foam::ThermoCloud<ParcelType>::TIntegrator() const
66     return TIntegrator_;
70 template<class ParcelType>
71 inline bool Foam::ThermoCloud<ParcelType>::radiation() const
73     return radiation_;
77 template<class ParcelType>
78 inline Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
79 Foam::ThermoCloud<ParcelType>::hsTrans()
81     return hsTrans_;
85 template<class ParcelType>
86 inline Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
87 Foam::ThermoCloud<ParcelType>::Sh() const
89     tmp<DimensionedField<scalar, volMesh> > tSh
90     (
91         new DimensionedField<scalar, volMesh>
92         (
93             IOobject
94             (
95                 this->name() + "Sh",
96                 this->db().time().timeName(),
97                 this->mesh(),
98                 IOobject::NO_READ,
99                 IOobject::AUTO_WRITE,
100                 false
101             ),
102             hsTrans_/(this->mesh().V()*this->db().time().deltaT())
103         )
104     );
106     return tSh;
110 template<class ParcelType>
111 inline Foam::tmp<Foam::volScalarField>
112 Foam::ThermoCloud<ParcelType>::Ep() const
114     tmp<volScalarField> tEp
115     (
116         new volScalarField
117         (
118             IOobject
119             (
120                 this->name() + "radiation::Ep",
121                 this->db().time().timeName(),
122                 this->db(),
123                 IOobject::NO_READ,
124                 IOobject::NO_WRITE,
125                 false
126             ),
127             this->mesh(),
128             dimensionedScalar("zero", dimMass/dimLength/pow3(dimTime), 0.0)
129         )
130     );
132     // Need to check if coupled as field is created on-the-fly
133     if (radiation_ && this->coupled())
134     {
135         scalarField& Ep = tEp().internalField();
136         const scalarField& V = this->mesh().V();
137         const scalar epsilon = constProps_.epsilon0();
139         forAllConstIter(typename ThermoCloud<ParcelType>, *this, iter)
140         {
141             const ParcelType& p = iter();
142             const label cellI = p.cell();
143             Ep[cellI] += p.nParticle()*p.areaP()*pow4(p.T());
144         }
146         Ep *= epsilon*radiation::sigmaSB.value()/V;
147     }
149     return tEp;
153 template<class ParcelType>
154 inline Foam::tmp<Foam::volScalarField>
155 Foam::ThermoCloud<ParcelType>::ap() const
157     tmp<volScalarField> tap
158     (
159         new volScalarField
160         (
161             IOobject
162             (
163                 this->name() + "radiation::ap",
164                 this->db().time().timeName(),
165                 this->db(),
166                 IOobject::NO_READ,
167                 IOobject::NO_WRITE,
168                 false
169             ),
170             this->mesh(),
171             dimensionedScalar("zero", dimless/dimLength, 0.0)
172         )
173     );
175     // Need to check if coupled as field is created on-the-fly
176     if (radiation_ && this->coupled())
177     {
178         scalarField& ap = tap().internalField();
179         const scalarField& V = this->mesh().V();
180         const scalar epsilon = constProps_.epsilon0();
182         forAllConstIter(typename ThermoCloud<ParcelType>, *this, iter)
183         {
184             const ParcelType& p = iter();
185             const label cellI = p.cell();
186             ap[cellI] += p.nParticle()*p.areaP();
187         }
189         ap *= epsilon/V;
190     }
192     return tap;
196 template<class ParcelType>
197 inline Foam::tmp<Foam::volScalarField>
198 Foam::ThermoCloud<ParcelType>::sigmap() const
200     tmp<volScalarField> tsigmap
201     (
202         new volScalarField
203         (
204             IOobject
205             (
206                 this->name() + "radiation::sigmap",
207                 this->db().time().timeName(),
208                 this->db(),
209                 IOobject::NO_READ,
210                 IOobject::NO_WRITE,
211                 false
212             ),
213             this->mesh(),
214             dimensionedScalar("zero", dimless/dimLength, 0.0)
215         )
216     );
218     // Need to check if coupled as field is created on-the-fly
219     if (radiation_ && this->coupled())
220     {
221         scalarField& sigmap = tsigmap().internalField();
223         const scalarField& V = this->mesh().V();
224         const scalar epsilon = constProps_.epsilon0();
225         const scalar f = constProps_.f0();
227         forAllConstIter(typename ThermoCloud<ParcelType>, *this, iter)
228         {
229             const ParcelType& p = iter();
230             const label cellI = p.cell();
231             sigmap[cellI] += p.nParticle()*p.areaP();
232         }
234         sigmap *= (1.0 - f)*(1.0 - epsilon)/V;
235     }
237     return tsigmap;
241 // ************************************************************************* //