ReactingCloudI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
27 
28 template<class CloudType>
31 {
32  return cloudCopyPtr_();
33 }
34 
35 
36 template<class CloudType>
37 inline const typename CloudType::particleType::constantProperties&
39 {
40  return constProps_;
41 }
42 
43 
44 template<class CloudType>
45 inline typename CloudType::particleType::constantProperties&
47 {
48  return constProps_;
49 }
50 
51 
52 template<class CloudType>
55 {
56  return compositionModel_;
57 }
58 
59 
60 template<class CloudType>
63 {
64  return phaseChangeModel_;
65 }
66 
67 
68 template<class CloudType>
71 {
72  return phaseChangeModel_();
73 }
74 
75 
76 template<class CloudType>
79 {
80  return rhoTrans_[i];
81 }
82 
83 
84 template<class CloudType>
85 inline
88 {
89  return rhoTrans_;
90 }
91 
92 
93 template<class CloudType>
96 {
97  return rhoTrans_;
98 }
99 
100 
101 template<class CloudType>
103 (
104  const label i,
105  volScalarField& Yi
106 ) const
107 {
108  if (this->solution().coupled())
109  {
110  if (this->solution().semiImplicit("Yi"))
111  {
112  tmp<volScalarField> trhoTrans
113  (
115  (
116  this->name() + ":rhoTrans",
117  this->mesh(),
119  )
120  );
121 
122  volScalarField& sourceField = trhoTrans.ref();
123 
124  sourceField.primitiveFieldRef() =
125  rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
126 
127  const dimensionedScalar Yismall("Yismall", dimless, small);
128 
129  return
130  fvm::Sp(neg(sourceField)*sourceField/(Yi + Yismall), Yi)
131  + pos0(sourceField)*sourceField;
132  }
133  else
134  {
136  fvScalarMatrix& fvm = tfvm.ref();
137 
138  fvm.source() = -rhoTrans_[i]/this->db().time().deltaTValue();
139 
140  return tfvm;
141  }
142  }
143 
145 }
146 
147 
148 template<class CloudType>
151 {
153  (
155  (
156  this->name() + ":rhoTrans",
157  this->mesh(),
159  (
160  rhoTrans_[0].dimensions()/dimTime/dimVolume,
161  0
162  )
163  )
164  );
165 
166  if (this->solution().coupled())
167  {
168  scalarField& rhoi = tRhoi.ref();
169  rhoi = rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
170  }
171 
172  return tRhoi;
173 }
174 
175 
176 template<class CloudType>
179 {
181  (
183  (
184  this->name() + ":rhoTrans",
185  this->mesh(),
187  (
188  rhoTrans_[0].dimensions()/dimTime/dimVolume,
189  0
190  )
191  )
192  );
193 
194  if (this->solution().coupled())
195  {
196  scalarField& sourceField = trhoTrans.ref();
197  forAll(rhoTrans_, i)
198  {
199  sourceField += rhoTrans_[i];
200  }
201 
202  sourceField /= this->db().time().deltaTValue()*this->mesh().V();
203  }
204 
205  return trhoTrans;
206 }
207 
208 
209 template<class CloudType>
212 {
213  if (this->solution().coupled())
214  {
215  tmp<volScalarField> trhoTrans
216  (
218  (
219  this->name() + ":rhoTrans",
220  this->mesh(),
222  )
223  );
224 
225  scalarField& sourceField = trhoTrans.ref();
226 
227  if (this->solution().semiImplicit("rho"))
228  {
229 
230  forAll(rhoTrans_, i)
231  {
232  sourceField += rhoTrans_[i];
233  }
234  sourceField /= this->db().time().deltaTValue()*this->mesh().V();
235 
236  return fvm::SuSp(trhoTrans()/rho, rho);
237  }
238  else
239  {
241  fvScalarMatrix& fvm = tfvm.ref();
242 
243  forAll(rhoTrans_, i)
244  {
245  sourceField += rhoTrans_[i];
246  }
247 
248  fvm.source() = -trhoTrans()/this->db().time().deltaT();
249 
250  return tfvm;
251  }
252  }
253 
255 }
256 
257 
258 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< GeometricField< Type, fvPatchField, volMesh > > SuSp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:102
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Templated phase change model class.
Definition: ReactingCloud.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const PtrList< volScalarField::Internal > & rhoTrans() const
Return const access to mass source fields.
tmp< volScalarField::Internal > Srho() const
Return tmp total mass source for carrier phase.
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const ReactingCloud & cloudCopy() const
Return a reference to the cloud copy.
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
const parcelType::constantProperties & constProps() const
Return the constant properties.
dimensionedScalar pos0(const dimensionedScalar &ds)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
Field< Type > & source()
Definition: fvMatrix.H:292
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Templated base class for reacting cloud.
Definition: ReactingCloud.H:62
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const CompositionModel< ReactingCloud< CloudType > > & composition() const
Return const access to reacting composition model.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
Templated reacting parcel composition model class Consists of carrier species (via thermo package)...
Definition: ReactingCloud.H:52
tmp< fvScalarMatrix > SYi(const label i, volScalarField &Yi) const
Return mass source term for specie i - specie eqn.
const PhaseChangeModel< ReactingCloud< CloudType > > & phaseChange() const
Return const access to reacting phase change model.
zeroField Sp
Definition: alphaSuSp.H:2