1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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 // * * * * * * * * * * * * * * * 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
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)
68 template<class CloudType>
70 const Foam::PtrList<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >&
71 Foam::ReactingCloud<CloudType>::rhoTrans() const
77 template<class CloudType>
78 inline Foam::PtrList<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >&
79 Foam::ReactingCloud<CloudType>::rhoTrans()
85 template<class CloudType>
86 inline Foam::tmp<Foam::fvScalarMatrix> Foam::ReactingCloud<CloudType>::SYi
92 if (this->solution().coupled())
94 if (this->solution().semiImplicit("Yi"))
96 tmp<volScalarField> trhoTrans
102 this->name() + "rhoTrans",
103 this->db().time().timeName(),
110 dimensionedScalar("zero", dimMass/dimTime/dimVolume, 0.0)
114 volScalarField& sourceField = trhoTrans();
116 sourceField.internalField() =
117 rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
119 const dimensionedScalar YiSMALL("YiSMALL", dimless, SMALL);
122 fvm::Sp(neg(sourceField)*sourceField/(Yi + YiSMALL), Yi)
123 + pos(sourceField)*sourceField;
127 tmp<fvScalarMatrix> tfvm(new fvScalarMatrix(Yi, dimMass/dimTime));
128 fvScalarMatrix& fvm = tfvm();
130 fvm.source() = -rhoTrans_[i]/this->db().time().deltaTValue();
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
146 new DimensionedField<scalar, volMesh>
150 this->name() + "rhoTrans",
151 this->db().time().timeName(),
161 rhoTrans_[0].dimensions()/dimTime/dimVolume,
167 if (this->solution().coupled())
169 scalarField& rhoi = tRhoi();
170 rhoi = rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
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
183 new DimensionedField<scalar, volMesh>
187 this->name() + "rhoTrans",
188 this->db().time().timeName(),
198 rhoTrans_[0].dimensions()/dimTime/dimVolume,
204 if (this->solution().coupled())
206 scalarField& sourceField = trhoTrans();
209 sourceField += rhoTrans_[i];
212 sourceField /= this->db().time().deltaTValue()*this->mesh().V();
219 template<class CloudType>
220 inline Foam::tmp<Foam::fvScalarMatrix>
221 Foam::ReactingCloud<CloudType>::Srho(volScalarField& rho) const
223 if (this->solution().coupled())
225 tmp<volScalarField> trhoTrans
231 this->name() + "rhoTrans",
232 this->db().time().timeName(),
239 dimensionedScalar("zero", dimMass/dimTime/dimVolume, 0.0)
243 scalarField& sourceField = trhoTrans();
245 if (this->solution().semiImplicit("rho"))
250 sourceField += rhoTrans_[i];
252 sourceField /= this->db().time().deltaTValue()*this->mesh().V();
254 return fvm::SuSp(trhoTrans()/rho, rho);
258 tmp<fvScalarMatrix> tfvm(new fvScalarMatrix(rho, dimMass/dimTime));
259 fvScalarMatrix& fvm = tfvm();
263 sourceField += rhoTrans_[i];
266 fvm.source() = -trhoTrans()/this->db().time().deltaT();
272 return tmp<fvScalarMatrix>(new fvScalarMatrix(rho, dimMass/dimTime));
276 // ************************************************************************* //