PopulationBalancePhaseSystem.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29
30
31// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
32
33template<class BasePhaseSystem>
36(
37 const phasePairKey& key
38) const
39{
40 if (!pDmdt_.found(key))
41 {
42 return phaseSystem::dmdt(key);
43 }
44
45 const scalar pDmdtSign(Pair<word>::compare(pDmdt_.find(key).key(), key));
46
47 return pDmdtSign**pDmdt_[key];
48}
49
50
51// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52
53template<class BasePhaseSystem>
56(
57 const fvMesh& mesh
58)
59:
60 BasePhaseSystem(mesh),
61
62 populationBalances_
63 (
64 this->lookup("populationBalances"),
65 diameterModels::populationBalanceModel::iNew(*this, pDmdt_)
66 )
67{
69 {
72
74 {
75 const phasePairKey& key = iter.key();
76
77 if (!this->phasePairs_.found(key))
78 {
79 this->phasePairs_.insert
80 (
81 key,
83 (
84 new phasePair
85 (
86 this->phaseModels_[key.first()],
87 this->phaseModels_[key.second()]
88 )
89 )
90 );
91 }
92 }
93 }
94
96 (
98 this->phasePairs_,
99 phasePairIter
100 )
101 {
102 const phasePair& pair(phasePairIter());
103
104 if (pair.ordered())
105 {
106 continue;
107 }
108
109 // Initially assume no mass transfer
110 pDmdt_.set
111 (
112 pair,
114 (
116 (
117 IOobject::groupName("pDmdt", pair.name()),
118 this->mesh().time().timeName(),
119 this->mesh(),
122 ),
123 this->mesh(),
125 )
126 );
127 }
128}
129
130
131// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
132
133template<class BasePhaseSystem>
136{}
137
138
139// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
140
141template<class BasePhaseSystem>
144(
145 const phasePairKey& key
146) const
147{
148 return BasePhaseSystem::dmdt(key) + this->pDmdt(key);
149}
150
151
152template<class BasePhaseSystem>
155{
156 PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
157
158 forAllConstIter(pDmdtTable, pDmdt_, pDmdtIter)
159 {
160 const phasePair& pair = this->phasePairs_[pDmdtIter.key()];
161 const volScalarField& pDmdt = *pDmdtIter();
162
163 this->addField(pair.phase1(), "dmdt", pDmdt, dmdts);
164 this->addField(pair.phase2(), "dmdt", - pDmdt, dmdts);
165 }
166
167 return dmdts;
168}
169
170
171template<class BasePhaseSystem>
174{
176 BasePhaseSystem::massTransfer();
177
178 phaseSystem::massTransferTable& eqns = eqnsPtr();
179
181 (
183 this->phasePairs_,
184 phasePairIter
185 )
186 {
187 const phasePair& pair(phasePairIter());
188
189 if (pair.ordered())
190 {
191 continue;
192 }
193
194 const phaseModel& phase = pair.phase1();
195 const phaseModel& otherPhase = pair.phase2();
196
197 // Note that the phase YiEqn does not contain a continuity error term,
198 // so these additions represent the entire mass transfer
199
200 const volScalarField dmdt(this->pDmdt(pair));
201 const volScalarField dmdt12(negPart(dmdt));
202 const volScalarField dmdt21(posPart(dmdt));
203
204 const PtrList<volScalarField>& Yi = phase.Y();
205
206 forAll(Yi, i)
207 {
208 const word name
209 (
210 IOobject::groupName(Yi[i].member(), phase.name())
211 );
212
213 const word otherName
214 (
215 IOobject::groupName(Yi[i].member(), otherPhase.name())
216 );
217
218 *eqns[name] +=
219 dmdt21*eqns[otherName]->psi()
220 + fvm::Sp(dmdt12, eqns[name]->psi());
221
222 *eqns[otherName] -=
223 dmdt12*eqns[name]->psi()
224 + fvm::Sp(dmdt21, eqns[otherName]->psi());
225 }
226 }
227
228 return eqnsPtr;
229}
230
231
232template<class BasePhaseSystem>
234{
235 if (BasePhaseSystem::read())
236 {
237 bool readOK = true;
238
239 // Models ...
240
241 return readOK;
242 }
243 else
244 {
245 return false;
246 }
247}
248
249
250template<class BasePhaseSystem>
252{
253 BasePhaseSystem::solve();
254
255 forAll(populationBalances_, i)
256 {
257 populationBalances_[i].solve();
258 }
259}
260
261
262// ************************************************************************* //
tmp< volScalarField > dmdt() const
Return the blended mass transfer rate.
bool set(const Key &key, T *ptr)
Assign a new entry, overwriting existing entries.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
Class which provides population balance functionality.
PtrList< diameterModels::populationBalanceModel > populationBalances_
populationBalanceModels
pDmdtTable pDmdt_
Interfacial Mass transfer rate.
virtual tmp< volScalarField > pDmdt(const phasePairKey &key) const
Return the population balance mass transfer rate.
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
virtual void solve()
Solve all population balance equations.
virtual bool read()
Read base phaseProperties dictionary.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Class that solves the univariate population balance equation by means of a class method (also called ...
const phasePairTable & phasePairs() const
Return list of unordered phasePairs in this populationBalance.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
const word & name() const
Definition: phaseModel.H:144
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:68
bool ordered() const noexcept
Return the ordered flag.
Definition: phasePairKey.H:98
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
virtual word name() const
Pair name.
Definition: phasePair.C:69
const multiphaseInter::phaseModel & phase1() const
Definition: phasePairI.H:30
const multiphaseInter::phaseModel & phase2() const
Definition: phasePairI.H:36
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
const word & name() const
Definition: phase.H:111
Lookup type of boundary radiation properties.
Definition: lookup.H:66
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
const volScalarField & psi
dynamicFvMesh & mesh
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensionedScalar negPart(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimDensity
dimensionedScalar posPart(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:381