fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / lagrangian / intermediate / clouds / Templates / ThermoCloud / ThermoCloudI.H
blob77b619183da7e5ccb3e274f58b02d07b4e5e29a4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 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
19     for more details.
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
35     return constProps_;
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
67     return TIntegrator_;
71 template<class ParcelType>
72 inline bool Foam::ThermoCloud<ParcelType>::radiation() const
74     return radiation_;
78 template<class ParcelType>
79 inline Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
80 Foam::ThermoCloud<ParcelType>::hsTrans()
82     return 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
91     (
92         new DimensionedField<scalar, volMesh>
93         (
94             IOobject
95             (
96                 this->name() + "Sh",
97                 this->db().time().timeName(),
98                 this->mesh(),
99                 IOobject::NO_READ,
100                 IOobject::AUTO_WRITE,
101                 false
102             ),
103             hsTrans_/(this->mesh().V()*this->db().time().deltaT())
104         )
105     );
107     return tSh;
111 template<class ParcelType>
112 inline Foam::tmp<Foam::volScalarField>
113 Foam::ThermoCloud<ParcelType>::Ep() const
115     tmp<volScalarField> tEp
116     (
117         new volScalarField
118         (
119             IOobject
120             (
121                 this->name() + "radiation::Ep",
122                 this->db().time().timeName(),
123                 this->db(),
124                 IOobject::NO_READ,
125                 IOobject::NO_WRITE,
126                 false
127             ),
128             this->mesh(),
129             dimensionedScalar("zero", dimMass/dimLength/pow3(dimTime), 0.0)
130         )
131     );
133     // Need to check if coupled as field is created on-the-fly
134     if (radiation_ && this->coupled())
135     {
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)
141         {
142             const ParcelType& p = iter();
143             const label cellI = p.cell();
144             Ep[cellI] += p.nParticle()*p.areaP()*pow4(p.T());
145         }
147         Ep *= epsilon*radiation::sigmaSB.value()/V;
148     }
150     return tEp;
154 template<class ParcelType>
155 inline Foam::tmp<Foam::volScalarField>
156 Foam::ThermoCloud<ParcelType>::ap() const
158     tmp<volScalarField> tap
159     (
160         new volScalarField
161         (
162             IOobject
163             (
164                 this->name() + "radiation::ap",
165                 this->db().time().timeName(),
166                 this->db(),
167                 IOobject::NO_READ,
168                 IOobject::NO_WRITE,
169                 false
170             ),
171             this->mesh(),
172             dimensionedScalar("zero", dimless/dimLength, 0.0)
173         )
174     );
176     // Need to check if coupled as field is created on-the-fly
177     if (radiation_ && this->coupled())
178     {
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)
184         {
185             const ParcelType& p = iter();
186             const label cellI = p.cell();
187             ap[cellI] += p.nParticle()*p.areaP();
188         }
190         ap *= epsilon/V;
191     }
193     return tap;
197 template<class ParcelType>
198 inline Foam::tmp<Foam::volScalarField>
199 Foam::ThermoCloud<ParcelType>::sigmap() const
201     tmp<volScalarField> tsigmap
202     (
203         new volScalarField
204         (
205             IOobject
206             (
207                 this->name() + "radiation::sigmap",
208                 this->db().time().timeName(),
209                 this->db(),
210                 IOobject::NO_READ,
211                 IOobject::NO_WRITE,
212                 false
213             ),
214             this->mesh(),
215             dimensionedScalar("zero", dimless/dimLength, 0.0)
216         )
217     );
219     // Need to check if coupled as field is created on-the-fly
220     if (radiation_ && this->coupled())
221     {
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)
229         {
230             const ParcelType& p = iter();
231             const label cellI = p.cell();
232             sigmap[cellI] += p.nParticle()*p.areaP();
233         }
235         sigmap *= (1.0 - f)*(1.0 - epsilon)/V;
236     }
238     return tsigmap;
242 // ************************************************************************* //