All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
phaseChange.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) 2018-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 
26 #include "phaseChange.H"
28 #include "phaseSystem.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace driftModels
37 {
39  addToRunTimeSelectionTable(driftModel, phaseChange, dictionary);
40 }
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const populationBalanceModel& popBal,
50  const dictionary& dict
51 )
52 :
53  driftModel(popBal, dict),
54  interfaces_
55  (
56  dict.lookup("interfaces"),
57  phaseInterface::iNew(popBal_.fluid())
58  ),
59  numberWeighted_(dict.lookupOrDefault<Switch>("numberWeighted", false)),
60  W_(interfaces_.size()),
61  dmdtfName_(dict.lookup("dmdtf")),
62  specieName_(dict.lookupOrDefault("specie", word()))
63 {
64  forAll(interfaces_, i)
65  {
66  W_.set
67  (
68  i,
69  new volScalarField
70  (
71  IOobject
72  (
73  IOobject::groupName(type() + ":W", interfaces_[i].name()),
74  popBal_.mesh().time().timeName(),
75  popBal_.mesh()
76  ),
77  popBal_.mesh(),
79  (
80  inv(numberWeighted_ ? dimVolume : dimLength),
81  Zero
82  )
83  )
84  );
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
92 {
93  forAll(interfaces_, i)
94  {
95  W_[i] = Zero;
96  }
97 
98  forAll(interfaces_, k)
99  {
101  (
102  HashTable<const velocityGroup*>,
104  iter
105  )
106  {
107  const velocityGroup& velGrp = *iter();
108 
109  if (interfaces_[k].contains(velGrp.phase()))
110  {
111  forAll(velGrp.sizeGroups(), i)
112  {
113  const sizeGroup& fi = velGrp.sizeGroups()[i];
114 
115  if (numberWeighted_)
116  {
117  W_[k] += fi*max(fi.phase(), small)/fi.x();
118  }
119  else
120  {
121  W_[k] += fi*max(fi.phase(), small)/fi.x()*fi.a();
122  }
123  }
124  }
125  }
126  }
127 }
128 
129 
131 (
132  volScalarField& driftRate,
133  const label i
134 )
135 {
136  const velocityGroup& velGrp = popBal_.sizeGroups()[i].VelocityGroup();
137 
138  forAll(interfaces_, k)
139  {
140  if (interfaces_[k].contains(velGrp.phase()))
141  {
142  const volScalarField& dmidtf =
144  (
146  (
148  (
149  dmdtfName_,
150  specieName_
151  ),
152  interfaces_[k].name()
153  )
154  );
155 
156  const scalar dmidtfSign =
157  interfaces_[k].index(velGrp.phase()) == 0 ? +1 : -1;
158 
159  const sizeGroup& fi = popBal_.sizeGroups()[i];
160 
161  tmp<volScalarField> dDriftRate
162  (
163  dmidtfSign*dmidtf/(fi.phase().rho()*W_[k])
164  );
165 
166  if (!numberWeighted_)
167  {
168  dDriftRate.ref() *= fi.a();
169  }
170 
171  driftRate += velGrp.phase()*dDriftRate;
172  }
173  }
174 }
175 
176 
177 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
virtual void precompute()
Precompute diameter independent expressions.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: driftModel.H:60
twoPhaseChangeModel & phaseChange
label k
Boltzmann constant.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
phaseChange(const populationBalanceModel &popBal, const dictionary &dict)
Construct from a population balance model and a dictionary.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
static word groupName(Name name, const word &group)
phaseSystem & fluid
Definition: createFields.H:11
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const dimensionSet dimVolume
const fvMesh & mesh() const
Return reference to the mesh.
Namespace for OpenFOAM.
const HashTable< const velocityGroup * > & velocityGroupPtrs() const
Return the velocity groups belonging to this populationBalance.