densityChange.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 
26 #include "densityChange.H"
28 #include "fvcDdt.H"
29 #include "fvcDiv.H"
30 #include "fvcSup.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace diameterModels
37 {
38 namespace driftModels
39 {
42 }
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const populationBalanceModel& popBal,
52  const dictionary& dict
53 )
54 :
55  driftModel(popBal, dict)
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60 
62 (
63  volScalarField& driftRate,
64  const label i
65 )
66 {
67  const sizeGroup& fi = popBal_.sizeGroups()[i];
68  const phaseModel& phase = fi.phase();
69  const volScalarField& alpha = phase;
70  const volScalarField& rho = phase.rho();
71 
72  driftRate -=
73  fi.x()
74  *(
76  - fvc::Sp(fvc::ddt(alpha) + fvc::div(phase.alphaPhi()), rho)
77  )/rho;
78 }
79 
80 
81 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Base class for drift models.
Definition: driftModel.H:54
Drift rate induced by changes in density.
Definition: densityChange.H:55
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.
Definition: densityChange.C:62
densityChangeDrift(const populationBalanceModel &popBal, const dictionary &dict)
Definition: densityChange.C:50
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:57
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:36
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
virtual tmp< surfaceScalarField > alphaRhoPhi() const =0
Return the mass flux of the phase.
virtual const volScalarField & rho() const =0
Return the density field.
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the field for explicit evaluation of implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(driftModel, constantDrift, dictionary)
defineTypeNameAndDebug(constantDrift, 0)
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
Namespace for OpenFOAM.
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
dictionary dict