powerLaw.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 -------------------------------------------------------------------------------
11 License
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 "powerLaw.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace diameterModels
37 {
38 namespace breakupModels
39 {
40  defineTypeNameAndDebug(powerLaw, 0);
41  addToRunTimeSelectionTable(breakupModel, powerLaw, dictionary);
42 }
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const populationBalanceModel& popBal,
52  const dictionary& dict
53 )
54 :
55  breakupModel(popBal, dict),
56  power_(dict.get<scalar>("power"))
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
63 (
64  volScalarField& breakupRate,
65  const label i
66 )
67 {
68  const sizeGroup& fi = popBal_.sizeGroups()[i];
69 
70  breakupRate.primitiveFieldRef() = pow(fi.x().value(), power_);
71 }
72 
73 
74 // ************************************************************************* //
powerLaw.H
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::diameterModels::breakupModels::powerLaw::powerLaw
powerLaw(const populationBalanceModel &popBal, const dictionary &dict)
Definition: powerLaw.C:50
Foam::diameterModels::breakupModels::powerLaw::setBreakupRate
virtual void setBreakupRate(volScalarField &breakupRate, const label i)
Set total breakupRate.
Definition: powerLaw.C:63
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::diameterModels::sizeGroup
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:96
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::diameterModels::populationBalanceModel
Class that solves the univariate population balance equation by means of a class method (also called ...
Definition: populationBalanceModel.H:179
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::diameterModels::breakupModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(breakupModel, exponential, dictionary)
Foam::diameterModels::breakupModel
Base class for breakup models which give a total breakup rate and a separate daughter size distributi...
Definition: breakupModel.H:55
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::diameterModels::breakupModels::defineTypeNameAndDebug
defineTypeNameAndDebug(exponential, 0)
Foam::diameterModels::sizeGroup::x
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:59