BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / lagrangian / intermediate / clouds / Templates / ReactingCloud / ReactingCloudI.H
blob68925178d39e62d96c89a7734bd48f8465f4791d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
28 template<class CloudType>
29 inline const Foam::ReactingCloud<CloudType>&
30 Foam::ReactingCloud<CloudType>::cloudCopy() const
32     return cloudCopyPtr_();
36 template<class CloudType>
37 inline const typename CloudType::particleType::constantProperties&
38 Foam::ReactingCloud<CloudType>::constProps() const
40     return constProps_;
44 template<class CloudType>
45 inline const Foam::CompositionModel<Foam::ReactingCloud<CloudType> >&
46 Foam::ReactingCloud<CloudType>::composition() const
48     return compositionModel_;
52 template<class CloudType>
53 inline const Foam::PhaseChangeModel<Foam::ReactingCloud<CloudType> >&
54 Foam::ReactingCloud<CloudType>::phaseChange() const
56     return phaseChangeModel_;
60 template<class CloudType>
61 inline Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
62 Foam::ReactingCloud<CloudType>::rhoTrans(const label i)
64     return rhoTrans_[i];
68 template<class CloudType>
69 inline
70 const Foam::PtrList<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >&
71 Foam::ReactingCloud<CloudType>::rhoTrans() const
73     return rhoTrans_;
77 template<class CloudType>
78 inline Foam::PtrList<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >&
79 Foam::ReactingCloud<CloudType>::rhoTrans()
81     return rhoTrans_;
85 template<class CloudType>
86 inline Foam::tmp<Foam::fvScalarMatrix> Foam::ReactingCloud<CloudType>::SYi
88     const label i,
89     volScalarField& Yi
90 ) const
92     if (this->solution().coupled())
93     {
94         if (this->solution().semiImplicit("Yi"))
95         {
96             tmp<volScalarField> trhoTrans
97             (
98                 new volScalarField
99                 (
100                     IOobject
101                     (
102                         this->name() + "rhoTrans",
103                         this->db().time().timeName(),
104                         this->db(),
105                         IOobject::NO_READ,
106                         IOobject::NO_WRITE,
107                         false
108                     ),
109                     this->mesh(),
110                     dimensionedScalar("zero", dimMass/dimTime/dimVolume, 0.0)
111                 )
112             );
114             volScalarField& sourceField = trhoTrans();
116             sourceField.internalField() =
117                 rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
119             const dimensionedScalar YiSMALL("YiSMALL", dimless, SMALL);
121             return
122                 fvm::Sp(neg(sourceField)*sourceField/(Yi + YiSMALL), Yi)
123               + pos(sourceField)*sourceField;
124         }
125         else
126         {
127             tmp<fvScalarMatrix> tfvm(new fvScalarMatrix(Yi, dimMass/dimTime));
128             fvScalarMatrix& fvm = tfvm();
130             fvm.source() = -rhoTrans_[i]/this->db().time().deltaTValue();
132             return tfvm;
133         }
134     }
136     return tmp<fvScalarMatrix>(new fvScalarMatrix(Yi, dimMass/dimTime));
140 template<class CloudType>
141 inline Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
142 Foam::ReactingCloud<CloudType>::Srho(const label i) const
144     tmp<DimensionedField<scalar, volMesh> > tRhoi
145     (
146         new DimensionedField<scalar, volMesh>
147         (
148             IOobject
149             (
150                 this->name() + "rhoTrans",
151                 this->db().time().timeName(),
152                 this->db(),
153                 IOobject::NO_READ,
154                 IOobject::NO_WRITE,
155                 false
156             ),
157             this->mesh(),
158             dimensionedScalar
159             (
160                 "zero",
161                 rhoTrans_[0].dimensions()/dimTime/dimVolume,
162                 0.0
163             )
164         )
165     );
167     if (this->solution().coupled())
168     {
169         scalarField& rhoi = tRhoi();
170         rhoi = rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
171     }
173     return tRhoi;
177 template<class CloudType>
178 inline Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
179 Foam::ReactingCloud<CloudType>::Srho() const
181     tmp<DimensionedField<scalar, volMesh> > trhoTrans
182     (
183         new DimensionedField<scalar, volMesh>
184         (
185             IOobject
186             (
187                 this->name() + "rhoTrans",
188                 this->db().time().timeName(),
189                 this->db(),
190                 IOobject::NO_READ,
191                 IOobject::NO_WRITE,
192                 false
193             ),
194             this->mesh(),
195             dimensionedScalar
196             (
197                 "zero",
198                 rhoTrans_[0].dimensions()/dimTime/dimVolume,
199                 0.0
200             )
201         )
202     );
204     if (this->solution().coupled())
205     {
206         scalarField& sourceField = trhoTrans();
207         forAll(rhoTrans_, i)
208         {
209             sourceField += rhoTrans_[i];
210         }
212         sourceField /= this->db().time().deltaTValue()*this->mesh().V();
213     }
215     return trhoTrans;
219 template<class CloudType>
220 inline Foam::tmp<Foam::fvScalarMatrix>
221 Foam::ReactingCloud<CloudType>::Srho(volScalarField& rho) const
223     if (this->solution().coupled())
224     {
225         tmp<volScalarField> trhoTrans
226         (
227             new volScalarField
228             (
229                 IOobject
230                 (
231                     this->name() + "rhoTrans",
232                     this->db().time().timeName(),
233                     this->db(),
234                     IOobject::NO_READ,
235                     IOobject::NO_WRITE,
236                     false
237                 ),
238                 this->mesh(),
239                 dimensionedScalar("zero", dimMass/dimTime/dimVolume, 0.0)
240             )
241         );
243         scalarField& sourceField = trhoTrans();
245         if (this->solution().semiImplicit("rho"))
246         {
248             forAll(rhoTrans_, i)
249             {
250                 sourceField += rhoTrans_[i];
251             }
252             sourceField /= this->db().time().deltaTValue()*this->mesh().V();
254             return fvm::SuSp(trhoTrans()/rho, rho);
255         }
256         else
257         {
258             tmp<fvScalarMatrix> tfvm(new fvScalarMatrix(rho, dimMass/dimTime));
259             fvScalarMatrix& fvm = tfvm();
261             forAll(rhoTrans_, i)
262             {
263                 sourceField += rhoTrans_[i];
264             }
266             fvm.source() = -trhoTrans()/this->db().time().deltaT();
268             return tfvm;
269         }
270     }
272     return tmp<fvScalarMatrix>(new fvScalarMatrix(rho, dimMass/dimTime));
276 // ************************************************************************* //