Loading...
Searching...
No Matches
densityChange.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "densityChange.H"
30#include "phaseSystem.H"
31#include "fvcDdt.H"
32#include "fvcGrad.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36namespace Foam
37{
38namespace diameterModels
39{
40namespace driftModels
41{
45}
46}
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
52(
53 const populationBalanceModel& popBal,
54 const dictionary& dict
55)
56:
58{}
59
60
61// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62
63
65(
66 volScalarField& driftRate,
67 const label i
68)
69{
70 const sizeGroup& fi = popBal_.sizeGroups()[i];
71 volScalarField& rho = const_cast<volScalarField&>(fi.phase().rho()());
72
73 driftRate -= (fvc::ddt(rho) + (fvc::grad(rho)&popBal_.U()))
74 *popBal_.sizeGroups()[i].x()/rho;
75}
76
77
78// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Base class for drift models.
Definition driftModel.H:49
const populationBalanceModel & popBal() const
Return reference to the populationBalanceModel.
Definition driftModel.H:144
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition driftModel.H:57
driftModel(const populationBalanceModel &popBal, const dictionary &dict)
Definition driftModel.C:67
Drift rate induced by changes in density.
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.
densityChangeDrift(const populationBalanceModel &popBal, const dictionary &dict)
Class that solves the univariate population balance equation by means of a class method (also called ...
const volVectorField & U() const
Return average velocity.
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition sizeGroup.H:95
const phaseModel & phase() const
Return const-reference to the phase.
Definition sizeGroupI.H:31
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
Calculate the first temporal derivative.
Calculate the gradient of the given field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition fvcDdt.C:40
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dictionary dict