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 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 #include "phasePairKey.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace driftModels
38 {
39  defineTypeNameAndDebug(phaseChange, 0);
40  addToRunTimeSelectionTable(driftModel, phaseChange, dictionary);
41 }
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const populationBalanceModel& popBal,
51  const dictionary& dict
52 )
53 :
54  driftModel(popBal, dict),
55  pairNames_(dict.lookup("pairNames")),
56  iDmdt_
57  (
58  IOobject
59  (
60  "iDmdt",
61  popBal.time().timeName(),
62  popBal.mesh()
63  ),
64  popBal.mesh(),
66  ),
67  N_
68  (
69  IOobject
70  (
71  "N",
72  popBal.mesh().time().timeName(),
73  popBal.mesh()
74  ),
75  popBal.mesh(),
77  )
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
83 
85 {
86  iDmdt_ *= 0.0;
87 
88  forAll(pairNames_, i)
89  {
90  const word& pairName = pairNames_[i];
91 
92  iDmdt_ +=
94  (
95  IOobject::groupName("iDmdt", pairName)
96  );
97  }
98 
99  N_ *= 0.0;
100 
102  {
103  const sizeGroup& fi = *popBal_.sizeGroups()[i];
104 
105  N_ += fi*max(fi.phase(), small)/fi.x();
106  }
107 }
108 
109 
111 (
112  volScalarField& driftRate,
113  const label i
114 )
115 {
116  const sizeGroup& fi = *popBal_.sizeGroups()[i];
117 
118  driftRate += iDmdt_/(N_*fi.phase().rho());
119 }
120 
121 
122 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
virtual void correct()
Correct diameter independent expressions.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: driftModel.H:57
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.
phaseChange(const populationBalanceModel &popBal, const dictionary &dict)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
static word groupName(Name name, const word &group)
const List< sizeGroup * > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const dimensionSet dimDensity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const fvMesh & mesh() const
Return reference to the mesh.
Namespace for OpenFOAM.