Loading...
Searching...
No Matches
constantNucleation.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) 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 "constantNucleation.H"
29#include "phaseSystem.H"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34namespace Foam
35{
36namespace diameterModels
37{
38namespace nucleationModels
39{
42 (
46 );
48}
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
56(
57 const populationBalanceModel& popBal,
58 const dictionary& dict
59)
60:
61 nucleationModel(popBal, dict),
62 d_("departureDiameter", dimLength, dict),
63 velGroup_
64 (
65 refCast<const velocityGroup>
66 (
67 popBal.mesh().lookupObject<phaseModel>
68 (
69 IOobject::groupName
70 (
71 "alpha",
72 dict.get<word>("velocityGroup")
73 )
74 ).dPtr()()
75 )
76 )
77{}
78
79
80// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81
83{
84 if
85 (
86 d_.value() < velGroup_.sizeGroups().first().d().value()
87 || d_.value() > velGroup_.sizeGroups().last().d().value()
88 )
89 {
91 << "Departure diameter " << d_.value() << " m outside of range ["
92 << velGroup_.sizeGroups().first().d().value() << ", "
93 << velGroup_.sizeGroups().last().d().value() << "] m" << endl
94 << " The nucleation rate is set to zero." << endl
95 << " Adjust discretization over property space to suppress this"
96 << " warning."
97 << endl;
98 }
99}
100
101
102void
104(
105 volScalarField& nucleationRate,
106 const label i
107)
108{
109 const sizeGroup& fi = popBal_.sizeGroups()[i];
110 phaseModel& phase = const_cast<phaseModel&>(fi.phase());
111 volScalarField& rho = phase.thermoRef().rho();
112
113 nucleationRate +=
114 popBal_.gamma(i, velGroup_.formFactor()*pow3(d_))
115 *(popBal_.fluid().fvOptions()(phase, rho)&rho)/rho/fi.x();
116}
117
118
119// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
fv::options & fvOptions
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Base class for nucleation models.
const populationBalanceModel & popBal() const
Return reference to the populationBalanceModel.
nucleationModel(const populationBalanceModel &popBal, const dictionary &dict)
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Constant nucleation rate within all classes. Used for verification and validation of the nucleation f...
virtual void correct()
Correct diameter independent expressions.
virtual void addToNucleationRate(volScalarField &nucleationRate, const label i)
Add to nucleationRate.
constantNucleation(const populationBalanceModel &popBal, const dictionary &dict)
Class that solves the univariate population balance equation by means of a class method (also called ...
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
const phaseModel & phase() const
Return const-reference to the phase.
Definition sizeGroupI.H:31
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition phase.H:151
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dictionary dict