PopulationBalancePhaseSystem.C
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) 2017-2019 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 
27 
28 
29 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
30 
31 template<class BasePhaseSystem>
34 (
35  const phasePairKey& key
36 ) const
37 {
38  if (!pDmdt_.found(key))
39  {
40  return phaseSystem::dmdt(key);
41  }
42 
43  const scalar pDmdtSign(Pair<word>::compare(pDmdt_.find(key).key(), key));
44 
45  return pDmdtSign**pDmdt_[key];
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 template<class BasePhaseSystem>
54 (
55  const fvMesh& mesh
56 )
57 :
58  BasePhaseSystem(mesh),
59 
60  populationBalances_
61  (
62  this->lookup("populationBalances"),
63  diameterModels::populationBalanceModel::iNew(*this, pDmdt_)
64  )
65 {
66  forAll(populationBalances_, i)
67  {
69  populationBalances_[i];
70 
72  {
73  const phasePairKey& key = iter.key();
74 
75  if (!this->phasePairs_.found(key))
76  {
77  this->phasePairs_.insert
78  (
79  key,
80  autoPtr<phasePair>
81  (
82  new phasePair
83  (
84  this->phaseModels_[key.first()],
85  this->phaseModels_[key.second()]
86  )
87  )
88  );
89  }
90  }
91  }
92 
94  (
96  this->phasePairs_,
97  phasePairIter
98  )
99  {
100  const phasePair& pair(phasePairIter());
101 
102  if (pair.ordered())
103  {
104  continue;
105  }
106 
107  // Initially assume no mass transfer
108  pDmdt_.insert
109  (
110  pair,
111  new volScalarField
112  (
113  IOobject
114  (
115  IOobject::groupName("pDmdt", pair.name()),
116  this->mesh().time().timeName(),
117  this->mesh(),
120  ),
121  this->mesh(),
123  )
124  );
125  }
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
130 
131 template<class BasePhaseSystem>
134 {}
135 
136 
137 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
138 
139 template<class BasePhaseSystem>
142 (
143  const phasePairKey& key
144 ) const
145 {
146  return BasePhaseSystem::dmdt(key) + this->pDmdt(key);
147 }
148 
149 
150 template<class BasePhaseSystem>
153 {
154  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
155 
156  forAllConstIter(pDmdtTable, pDmdt_, pDmdtIter)
157  {
158  const phasePair& pair = this->phasePairs_[pDmdtIter.key()];
159  const volScalarField& pDmdt = *pDmdtIter();
160 
161  this->addField(pair.phase1(), "dmdt", pDmdt, dmdts);
162  this->addField(pair.phase2(), "dmdt", - pDmdt, dmdts);
163  }
164 
165  return dmdts;
166 }
167 
168 
169 template<class BasePhaseSystem>
172 {
173  autoPtr<phaseSystem::massTransferTable> eqnsPtr =
175 
176  phaseSystem::massTransferTable& eqns = eqnsPtr();
177 
179  (
181  this->phasePairs_,
182  phasePairIter
183  )
184  {
185  const phasePair& pair(phasePairIter());
186 
187  if (pair.ordered())
188  {
189  continue;
190  }
191 
192  const phaseModel& phase = pair.phase1();
193  const phaseModel& otherPhase = pair.phase2();
194 
195  // Note that the phase YiEqn does not contain a continuity error term,
196  // so these additions represent the entire mass transfer
197 
198  const volScalarField dmdt(this->pDmdt(pair));
199  const volScalarField dmdt12(negPart(dmdt));
200  const volScalarField dmdt21(posPart(dmdt));
201 
202  const PtrList<volScalarField>& Yi = phase.Y();
203 
204  forAll(Yi, i)
205  {
206  const word name
207  (
208  IOobject::groupName(Yi[i].member(), phase.name())
209  );
210 
211  const word otherName
212  (
213  IOobject::groupName(Yi[i].member(), otherPhase.name())
214  );
215 
216  *eqns[name] +=
217  dmdt21*eqns[otherName]->psi()
218  + fvm::Sp(dmdt12, eqns[name]->psi());
219 
220  *eqns[otherName] -=
221  dmdt12*eqns[name]->psi()
222  + fvm::Sp(dmdt21, eqns[otherName]->psi());
223  }
224  }
225 
226  return eqnsPtr;
227 }
228 
229 
230 template<class BasePhaseSystem>
232 {
233  if (BasePhaseSystem::read())
234  {
235  bool readOK = true;
236 
237  // Models ...
238 
239  return readOK;
240  }
241  else
242  {
243  return false;
244  }
245 }
246 
247 
248 template<class BasePhaseSystem>
250 {
252 
253  forAll(populationBalances_, i)
254  {
255  populationBalances_[i].solve();
256  }
257 }
258 
259 
260 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual void solve()
Solve all population balance equations.
HashPtrTable< fvScalarMatrix > massTransferTable
Definition: phaseSystem.H:79
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition: phaseSystem.H:87
Class that solves the univariate population balance equation by means of a class method (also called ...
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:142
phaseSystem::massTransferTable & massTransfer(massTransferPtr())
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
const phasePairTable & phasePairs() const
Return list of unordered phasePairs in this populationBalance.
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual ~PopulationBalancePhaseSystem()
Destructor.
dimensionedScalar posPart(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
static word groupName(Name name, const word &group)
rhoEqn solve()
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimDensity
PopulationBalancePhaseSystem(const fvMesh &)
Construct from fvMesh.
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
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const volScalarField & psi
A class for managing temporary objects.
Definition: PtrList.H:53
virtual bool read()
Read base phaseProperties dictionary.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
virtual tmp< volScalarField > pDmdt(const phasePairKey &key) const
Return the population balance mass transfer rate.
dimensionedScalar negPart(const dimensionedScalar &ds)