constantNucleation.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) 2018-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 "constantNucleation.H"
27 #include "phaseSystem.H"
29 #include "mathematicalConstants.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace nucleationModels
38 {
39  defineTypeNameAndDebug(constantNucleation, 0);
41  (
42  nucleationModel,
43  constantNucleation,
44  dictionary
45  );
46 }
47 }
48 }
49 
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
57 (
58  const populationBalanceModel& popBal,
59  const dictionary& dict
60 )
61 :
62  nucleationModel(popBal, dict),
63  dNuc_("nucleationDiameter", dimLength, dict),
64  velGroup_
65  (
66  refCast<const velocityGroup>
67  (
68  popBal.mesh().lookupObject<phaseModel>
69  (
70  IOobject::groupName
71  (
72  "alpha",
73  dict.lookup("velocityGroup")
74  )
75  ).dPtr()()
76  )
77  )
78 {
79  if
80  (
81  dNuc_.value() < velGroup_.sizeGroups().first().dSph().value()
82  || dNuc_.value() > velGroup_.sizeGroups().last().dSph().value()
83  )
84  {
86  << "Nucleation diameter " << dNuc_.value() << "m outside of range ["
87  << velGroup_.sizeGroups().first().dSph().value() << ", "
88  << velGroup_.sizeGroups().last().dSph().value() << "]." << nl
89  << exit(FatalIOError);
90  }
91 }
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
96 void
98 (
99  volScalarField& nucleationRate,
100  const label i
101 )
102 {
103  const sizeGroup& fi = popBal_.sizeGroups()[i];
104  phaseModel& phase = const_cast<phaseModel&>(fi.phase());
105  volScalarField& rho = phase.thermoRef().rho();
106 
107  nucleationRate +=
108  popBal_.eta(i, pi/6.0*pow3(dNuc_))
109  *(popBal_.fluid().fvOptions()(phase, rho)&rho)/rho/fi.x();
110 }
111 
112 
113 // ************************************************************************* //
fv::options & fvOptions() const
Access the fvOptions.
Definition: phaseSystemI.H:164
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: driftModel.H:57
Macros for easy insertion into run-time selection tables.
const dimensionedScalar eta(const label i, const dimensionedScalar &v) const
Return allocation coefficient.
virtual void addToNucleationRate(volScalarField &nucleationRate, const label i)
Add to nucleationRate.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
const phaseSystem & fluid() const
Return reference to the phaseSystem.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
constantNucleation(const populationBalanceModel &popBal, const dictionary &dict)
Namespace for OpenFOAM.
IOerror FatalIOError