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-2020 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class BasePhaseSystem>
33 (
34  const fvMesh& mesh
35 )
36 :
37  BasePhaseSystem(mesh),
38 
39  populationBalances_
40  (
41  this->lookup("populationBalances"),
42  diameterModels::populationBalanceModel::iNew(*this, dmdtfs_)
43  )
44 {
45  forAll(populationBalances_, i)
46  {
48  populationBalances_[i];
49 
51  {
52  const phasePairKey& key = iter.key();
53 
54  if (!this->phasePairs_.found(key))
55  {
56  this->phasePairs_.insert
57  (
58  key,
59  autoPtr<phasePair>
60  (
61  new phasePair
62  (
63  this->phaseModels_[key.first()],
64  this->phaseModels_[key.second()]
65  )
66  )
67  );
68  }
69 
70  dmdtfs_.insert
71  (
72  key,
73  new volScalarField
74  (
75  IOobject
76  (
78  (
79  "populationBalance:dmdtf",
80  this->phasePairs_[key]->name()
81  ),
82  this->mesh().time().timeName(),
83  this->mesh(),
86  ),
87  this->mesh(),
89  )
90  );
91  }
92  }
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
98 template<class BasePhaseSystem>
101 {}
102 
103 
104 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
105 
106 template<class BasePhaseSystem>
109 (
110  const phasePairKey& key
111 ) const
112 {
113  tmp<volScalarField> tDmdt = BasePhaseSystem::dmdtf(key);
114 
115  if (!dmdtfs_.found(key))
116  {
117  const label dmdtSign(Pair<word>::compare(this->phasePairs_[key], key));
118 
119  tDmdt.ref() += dmdtSign**dmdtfs_[key];
120  }
121 
122  return tDmdt;
123 }
124 
125 
126 template<class BasePhaseSystem>
129 {
130  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
131 
132  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs_, dmdtfIter)
133  {
134  const phasePair& pair = this->phasePairs_[dmdtfIter.key()];
135  const volScalarField& pDmdt = *dmdtfIter();
136 
137  addField(pair.phase1(), "dmdt", pDmdt, dmdts);
138  addField(pair.phase2(), "dmdt", - pDmdt, dmdts);
139  }
140 
141  return dmdts;
142 }
143 
144 
145 template<class BasePhaseSystem>
148 {
149  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
151 
152  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
153 
154  this->addDmdtUfs(dmdtfs_, eqns);
155 
156  return eqnsPtr;
157 }
158 
159 
160 template<class BasePhaseSystem>
163 {
164  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
165  BasePhaseSystem::momentumTransferf();
166 
167  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
168 
169  this->addDmdtUfs(dmdtfs_, eqns);
170 
171  return eqnsPtr;
172 }
173 
174 
175 template<class BasePhaseSystem>
178 {
179  autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
180  BasePhaseSystem::heatTransfer();
181 
182  phaseSystem::heatTransferTable& eqns = eqnsPtr();
183 
184  this->addDmdtHefs(dmdtfs_, eqns);
185 
186  return eqnsPtr;
187 }
188 
189 
190 template<class BasePhaseSystem>
193 {
194  autoPtr<phaseSystem::specieTransferTable> eqnsPtr =
196 
197  phaseSystem::specieTransferTable& eqns = eqnsPtr();
198 
199  this->addDmdtYfs(dmdtfs_, eqns);
200 
201  return eqnsPtr;
202 }
203 
204 
205 template<class BasePhaseSystem>
207 {
208  if (BasePhaseSystem::read())
209  {
210  bool readOK = true;
211 
212  // Models ...
213 
214  return readOK;
215  }
216  else
217  {
218  return false;
219  }
220 }
221 
222 
223 template<class BasePhaseSystem>
225 (
226  const PtrList<volScalarField>& rAUs,
227  const PtrList<surfaceScalarField>& rAUfs
228 )
229 {
230  BasePhaseSystem::solve(rAUs, rAUfs);
231 
232  forAll(populationBalances_, i)
233  {
234  populationBalances_[i].solve();
235  }
236 }
237 
238 
239 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
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 ...
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:75
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:142
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.
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:77
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:79
rhoEqn solve()
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
virtual ~PopulationBalancePhaseSystem()
Destructor.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
virtual tmp< volScalarField > dmdtf(const phasePairKey &key) const
Return the mass transfer rate for an interface.
static word groupName(Name name, const word &group)
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdtfTable
Definition: phaseSystem.H:91
virtual void solve(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs)
Solve all population balance equations.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
word timeName
Definition: getTimeIndex.H:3
void addDmdtYfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from bulk mass transfers.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
phaseSystem::specieTransferTable & specieTransfer(specieTransferPtr())
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
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
phaseSystem::momentumTransferTable & momentumTransfer(momentumTransferPtr())
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 autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
virtual bool read()
Read base phaseProperties dictionary.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.