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-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 "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 {
40  defineTypeNameAndDebug(densityChangeDrift, 0);
41  addToRunTimeSelectionTable(driftModel, densityChangeDrift, dictionary);
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.thermo().rho();
71 
72  driftRate -=
73  fi.x()
74  *(
75  fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi())
76  - fvc::Sp(fvc::ddt(alpha) + fvc::div(phase.alphaPhi()), rho)
77  )/rho;
78 }
79 
80 
81 // ************************************************************************* //
dictionary dict
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: driftModel.H:60
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
Macros for easy insertion into run-time selection tables.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Calculate the first temporal derivative.
Calculate the field for explicit evaluation of implicit and explicit sources.
Calculate the divergence of the given field.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
densityChangeDrift(const populationBalanceModel &popBal, const dictionary &dict)
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.
Namespace for OpenFOAM.