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-2021 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  this->template validateMassTransfer
71  <
72  diameterModels::populationBalanceModel
73  >(this->phasePairs_[key]);
74 
75  dmdtfs_.insert
76  (
77  key,
78  new volScalarField
79  (
80  IOobject
81  (
83  (
84  "populationBalance:dmdtf",
85  this->phasePairs_[key]->name()
86  ),
87  this->mesh().time().timeName(),
88  this->mesh()
89  ),
90  this->mesh(),
92  )
93  );
94  }
95  }
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
100 
101 template<class BasePhaseSystem>
104 {}
105 
106 
107 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
108 
109 template<class BasePhaseSystem>
112 (
113  const phasePairKey& key
114 ) const
115 {
116  tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
117 
118  if (dmdtfs_.found(key))
119  {
120  const label dmdtSign(Pair<word>::compare(this->phasePairs_[key], key));
121 
122  tDmdtf.ref() += dmdtSign**dmdtfs_[key];
123  }
124 
125  return tDmdtf;
126 }
127 
128 
129 template<class BasePhaseSystem>
132 {
133  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
134 
135  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs_, dmdtfIter)
136  {
137  const phasePair& pair = this->phasePairs_[dmdtfIter.key()];
138  const volScalarField& pDmdt = *dmdtfIter();
139 
140  addField(pair.phase1(), "dmdt", pDmdt, dmdts);
141  addField(pair.phase2(), "dmdt", - pDmdt, dmdts);
142  }
143 
144  return dmdts;
145 }
146 
147 
148 template<class BasePhaseSystem>
151 {
152  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
154 
155  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
156 
157  this->addDmdtUfs(dmdtfs_, eqns);
158 
159  return eqnsPtr;
160 }
161 
162 
163 template<class BasePhaseSystem>
166 {
167  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
168  BasePhaseSystem::momentumTransferf();
169 
170  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
171 
172  this->addDmdtUfs(dmdtfs_, eqns);
173 
174  return eqnsPtr;
175 }
176 
177 
178 template<class BasePhaseSystem>
181 {
182  autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
183  BasePhaseSystem::heatTransfer();
184 
185  phaseSystem::heatTransferTable& eqns = eqnsPtr();
186 
187  this->addDmdtHefs(dmdtfs_, eqns);
188 
189  return eqnsPtr;
190 }
191 
192 
193 template<class BasePhaseSystem>
196 {
197  autoPtr<phaseSystem::specieTransferTable> eqnsPtr =
199 
200  phaseSystem::specieTransferTable& eqns = eqnsPtr();
201 
202  this->addDmdtYfs(dmdtfs_, eqns);
203 
204  return eqnsPtr;
205 }
206 
207 
208 template<class BasePhaseSystem>
210 (
211  const PtrList<volScalarField>& rAUs,
212  const PtrList<surfaceScalarField>& rAUfs
213 )
214 {
215  BasePhaseSystem::solve(rAUs, rAUfs);
216 
217  forAll(populationBalances_, i)
218  {
219  populationBalances_[i].solve();
220  }
221 }
222 
223 
224 template<class BasePhaseSystem>
226 {
228 
229  forAll(populationBalances_, i)
230  {
231  populationBalances_[i].correct();
232  }
233 }
234 
235 
236 template<class BasePhaseSystem>
238 {
239  if (BasePhaseSystem::read())
240  {
241  bool readOK = true;
242 
243  // Read models ...
244 
245  return readOK;
246  }
247  else
248  {
249  return false;
250  }
251 }
252 
253 
254 // ************************************************************************* //
#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:90
Class that solves the univariate population balance equation by means of a class method (also called ...
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:78
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:153
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:80
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:82
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
const dimensionSet dimTime
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.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/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
static word groupName(Name name, const word &group)
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdtfTable
Definition: phaseSystem.H:94
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.
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
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.
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())
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
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.