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-2020 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 "phaseSystem.H"
29 #include "fvcDdt.H"
30 #include "fvcDiv.H"
31 #include "fvcSup.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace diameterModels
38 {
39 namespace driftModels
40 {
41  defineTypeNameAndDebug(densityChangeDrift, 0);
42  addToRunTimeSelectionTable(driftModel, densityChangeDrift, dictionary);
43 }
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const populationBalanceModel& popBal,
53  const dictionary& dict
54 )
55 :
56  driftModel(popBal, dict)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
63 (
64  volScalarField& driftRate,
65  const label i
66 )
67 {
68  const sizeGroup& fi = popBal_.sizeGroups()[i];
69  const phaseModel& phase = fi.phase();
70  const volScalarField& alpha = phase;
71  const volScalarField& rho = phase.thermo().rho();
72 
73  driftRate -=
74  fi.x()/(rho*max(alpha, phase.residualAlpha()))
75  *(
76  fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi())
77  - fvc::Sp(fvc::ddt(alpha) + fvc::div(phase.alphaPhi()), rho)
78  );
79 }
80 
81 
82 // ************************************************************************* //
dictionary dict
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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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:57
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:57
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 sizeGroups belonging to this populationBalance.
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.
Namespace for OpenFOAM.