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-2022 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  populationBalances_
39  (
40  this->lookup("populationBalances"),
41  diameterModels::populationBalanceModel::iNew(*this, dmdtfs_)
42  )
43 {
44  forAll(populationBalances_, popBali)
45  {
46  const diameterModels::populationBalanceModel& popBal =
47  populationBalances_[popBali];
48 
50  (
51  HashTable<const diameterModels::velocityGroup*>,
52  popBal.velocityGroupPtrs(),
53  iter1
54  )
55  {
56  const diameterModels::velocityGroup& velGrp1 = *iter1();
57 
59  (
60  HashTable<const diameterModels::velocityGroup*>,
61  popBal.velocityGroupPtrs(),
62  iter2
63  )
64  {
65  const diameterModels::velocityGroup& velGrp2 = *iter2();
66 
67  const phaseInterface interface
68  (
69  velGrp1.phase(),
70  velGrp2.phase()
71  );
72 
73  if
74  (
75  &velGrp1 != &velGrp2
76  &&
77  !dmdtfs_.found(interface)
78  )
79  {
80  this->template validateMassTransfer
81  <diameterModels::populationBalanceModel>(interface);
82 
83  dmdtfs_.insert
84  (
85  interface,
86  new volScalarField
87  (
88  IOobject
89  (
91  (
92  "populationBalance:dmdtf",
93  interface.name()
94  ),
95  this->mesh().time().timeName(),
96  this->mesh()
97  ),
98  this->mesh(),
100  )
101  );
102  }
103  }
104  }
105  }
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
110 
111 template<class BasePhaseSystem>
114 {}
115 
116 
117 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
118 
119 template<class BasePhaseSystem>
122 (
123  const phaseInterfaceKey& key
124 ) const
125 {
126  tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
127 
128  if (dmdtfs_.found(key))
129  {
130  tDmdtf.ref() += *dmdtfs_[key];
131  }
132 
133  return tDmdtf;
134 }
135 
136 
137 template<class BasePhaseSystem>
140 {
141  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
142 
143  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs_, dmdtfIter)
144  {
145  const phaseInterface interface(*this, dmdtfIter.key());
146 
147  addField(interface.phase1(), "dmdt", *dmdtfIter(), dmdts);
148  addField(interface.phase2(), "dmdt", - *dmdtfIter(), dmdts);
149  }
150 
151  return dmdts;
152 }
153 
154 
155 template<class BasePhaseSystem>
158 {
159  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
161 
162  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
163 
164  this->addDmdtUfs(dmdtfs_, eqns);
165 
166  return eqnsPtr;
167 }
168 
169 
170 template<class BasePhaseSystem>
173 {
174  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
175  BasePhaseSystem::momentumTransferf();
176 
177  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
178 
179  this->addDmdtUfs(dmdtfs_, eqns);
180 
181  return eqnsPtr;
182 }
183 
184 
185 template<class BasePhaseSystem>
188 {
189  autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
190  BasePhaseSystem::heatTransfer();
191 
192  phaseSystem::heatTransferTable& eqns = eqnsPtr();
193 
194  this->addDmdtHefs(dmdtfs_, eqns);
195 
196  return eqnsPtr;
197 }
198 
199 
200 template<class BasePhaseSystem>
203 {
204  autoPtr<phaseSystem::specieTransferTable> eqnsPtr =
206 
207  phaseSystem::specieTransferTable& eqns = eqnsPtr();
208 
209  this->addDmdtYfs(dmdtfs_, eqns);
210 
211  return eqnsPtr;
212 }
213 
214 
215 template<class BasePhaseSystem>
217 (
218  const PtrList<volScalarField>& rAUs,
219  const PtrList<surfaceScalarField>& rAUfs
220 )
221 {
222  BasePhaseSystem::solve(rAUs, rAUfs);
223 
224  forAll(populationBalances_, i)
225  {
226  populationBalances_[i].solve();
227  }
228 }
229 
230 
231 template<class BasePhaseSystem>
233 {
235 
236  forAll(populationBalances_, i)
237  {
238  populationBalances_[i].correct();
239  }
240 }
241 
242 
243 template<class BasePhaseSystem>
245 {
246  if (BasePhaseSystem::read())
247  {
248  bool readOK = true;
249 
250  // Read models ...
251 
252  return readOK;
253  }
254  else
255  {
256  return false;
257  }
258 }
259 
260 
261 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:76
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:78
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:80
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
virtual ~PopulationBalancePhaseSystem()
Destructor.
fvMesh & mesh
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet dimTime
stressControl lookup("compactNormalStress") >> compactNormalStress
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
static word groupName(Name name, const word &group)
const dimensionSet dimDensity
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.
void addDmdtYfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from bulk mass transfers.
virtual void correct()
Correct derived properties.
phaseSystem::specieTransferTable & specieTransfer(specieTransferPtr())
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > dmdtfTable
Definition: phaseSystem.H:93
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
phaseSystem::momentumTransferTable & momentumTransfer(momentumTransferPtr())
PopulationBalancePhaseSystem(const fvMesh &)
Construct from fvMesh.
rhoEqn solve()
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
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
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
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.