Loading...
Searching...
No Matches
exponentialBreakup.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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "exponentialBreakup.H"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34namespace Foam
35{
36namespace diameterModels
37{
38namespace breakupModels
39{
43}
44}
45
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
50(
51 const populationBalanceModel& popBal,
52 const dictionary& dict
53)
54:
55 breakupModel(popBal, dict),
56 exponent_(dict.get<scalar>("exponent")),
57 C_(dict.get<scalar>("C"))
58{}
59
60
61// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62
64(
65 volScalarField& breakupRate,
66 const label i
67)
68{
69 const sizeGroup& fi = popBal_.sizeGroups()[i];
70
71 breakupRate.primitiveFieldRef() =
72 C_*exp(exponent_*fi.x().value());
73}
74
75
76// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Base class for breakup models which give a total breakup rate and a separate daughter size distributi...
const populationBalanceModel & popBal() const
Return reference to the populationBalanceModel.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
breakupModel(const populationBalanceModel &popBal, const dictionary &dict)
Exponential kernel. Used for verification and validation of the breakup formulation implemented in th...
virtual void setBreakupRate(volScalarField &breakupRate, const label i)
Set total breakupRate.
exponential(const populationBalanceModel &popBal, const dictionary &dict)
Class that solves the univariate population balance equation by means of a class method (also called ...
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 dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition sizeGroupI.H:52
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Type & value() const noexcept
Return const reference to value.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dictionary dict