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-2023 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)
42  )
43 {}
44 
45 
46 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
47 
48 template<class BasePhaseSystem>
51 {}
52 
53 
54 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
55 
56 template<class BasePhaseSystem>
59 (
60  const phaseInterfaceKey& key
61 ) const
62 {
63  tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
64 
65  forAll(populationBalances_, popBali)
66  {
67  if (populationBalances_[popBali].dmdtfs().found(key))
68  {
69  tDmdtf.ref() += *populationBalances_[popBali].dmdtfs()[key];
70  }
71  }
72 
73  return tDmdtf;
74 }
75 
76 
77 template<class BasePhaseSystem>
80 {
81  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
82 
83  forAll(populationBalances_, popBali)
84  {
86  (
88  populationBalances_[popBali].dmdtfs(),
89  dmdtfIter
90  )
91  {
92  const phaseInterface interface(*this, dmdtfIter.key());
93 
94  addField(interface.phase1(), "dmdt", *dmdtfIter(), dmdts);
95  addField(interface.phase2(), "dmdt", - *dmdtfIter(), dmdts);
96  }
97  }
98 
99  return dmdts;
100 }
101 
102 
103 template<class BasePhaseSystem>
106 {
108  BasePhaseSystem::momentumTransfer();
109 
110  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
111 
112  forAll(populationBalances_, popBali)
113  {
114  this->addDmdtUfs(populationBalances_[popBali].dmdtfs(), eqns);
115  }
116 
117  return eqnsPtr;
118 }
119 
120 
121 template<class BasePhaseSystem>
124 {
126  BasePhaseSystem::momentumTransferf();
127 
128  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
129 
130  forAll(populationBalances_, popBali)
131  {
132  this->addDmdtUfs(populationBalances_[popBali].dmdtfs(), eqns);
133  }
134 
135  return eqnsPtr;
136 }
137 
138 
139 template<class BasePhaseSystem>
142 {
144  BasePhaseSystem::heatTransfer();
145 
146  phaseSystem::heatTransferTable& eqns = eqnsPtr();
147 
148  forAll(populationBalances_, popBali)
149  {
150  this->addDmdtHefs(populationBalances_[popBali].dmdtfs(), eqns);
151  }
152 
153  return eqnsPtr;
154 }
155 
156 
157 template<class BasePhaseSystem>
160 {
162  BasePhaseSystem::specieTransfer();
163 
164  phaseSystem::specieTransferTable& eqns = eqnsPtr();
165 
166  forAll(populationBalances_, popBali)
167  {
168  this->addDmdtYfs(populationBalances_[popBali].dmdtfs(), eqns);
169  }
170 
171  return eqnsPtr;
172 }
173 
174 
175 template<class BasePhaseSystem>
177 (
178  const PtrList<volScalarField>& rAs
179 )
180 {
182 
183  forAll(populationBalances_, i)
184  {
185  populationBalances_[i].solve();
186  }
187 }
188 
189 
190 template<class BasePhaseSystem>
192 {
194 
195  forAll(populationBalances_, i)
196  {
197  populationBalances_[i].correct();
198  }
199 }
200 
201 
202 template<class BasePhaseSystem>
204 {
205  if (BasePhaseSystem::read())
206  {
207  bool readOK = true;
208 
209  // Read models ...
210 
211  return readOK;
212  }
213  else
214  {
215  return false;
216  }
217 }
218 
219 
220 // ************************************************************************* //
bool found
#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
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
virtual void correct()
Correct derived properties.
virtual void solve(const PtrList< volScalarField > &rAs)
Solve all population balance equations.
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
PopulationBalancePhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
virtual bool read()
Read base phaseProperties dictionary.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Word-pair based class used for keying interface models in hash tables.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const phaseModel & phase1() const
Return phase 1.
const phaseModel & phase2() const
Return phase 2.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
void addField(const label phasei, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
Definition: phaseSystemI.H:263
rhoEqn solve()