reactionDriven.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) 2019-2024 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 "reactionDriven.H"
28 #include "phaseSystem.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace nucleationModels
37 {
40  (
44  );
45 }
46 }
47 }
48 
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
56 (
57  const populationBalanceModel& popBal,
58  const dictionary& dict
59 )
60 :
61  nucleationModel(popBal, dict),
62  dNuc_("nucleationDiameter", dimLength, dict),
63  reactingPhase_
64  (
65  popBal_.fluid().phases()[dict.lookup<word>("reactingPhase")]
66  ),
67  dmdtfName_(dict.lookup("dmdtf")),
68  specieName_(dict.lookup("specie"))
69 {
70  if
71  (
72  dNuc_.value() < velGroup_.sizeGroups().first().dSph().value()
73  || dNuc_.value() > velGroup_.sizeGroups().last().dSph().value()
74  )
75  {
77  << "Nucleation diameter " << dNuc_.value() << "m outside of range ["
78  << velGroup_.sizeGroups().first().dSph().value() << ", "
79  << velGroup_.sizeGroups().last().dSph().value() << "]." << nl
80  << exit(FatalIOError);
81  }
82 }
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 void
89 (
90  volScalarField& nucleationRate,
91  const label i
92 )
93 {
94  const sizeGroup& fi = popBal_.sizeGroups()[i];
95  const volScalarField& rho = fi.phase().rho();
96 
97  const phaseInterface interface(velGroup_.phase(), reactingPhase_);
98 
99  const volScalarField& dmidtf =
100  popBal_.mesh().lookupObject<volScalarField>
101  (
103  (
104  IOobject::groupName(dmdtfName_, specieName_),
105  interface.name()
106  )
107  );
108 
109  const scalar dmidtfSign =
110  interface.index(velGroup_.phase()) == 0 ? +1 : -1;
111 
112  nucleationRate +=
113  popBal_.eta(i, pi/6*pow3(dNuc_))*dmidtfSign*dmidtf/rho/fi.x();
114 }
115 
116 
117 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
static word groupName(Name name, const word &group)
Base class for nucleation models.
const dictionary & dict() const
Return reference to model dictionary.
const velocityGroup & velGroup_
Velocity group in which the nucleation occurs.
virtual void addToNucleationRate(volScalarField &nucleationRate, const label i)
Add to nucleationRate.
reactionDriven(const populationBalanceModel &popBal, const dictionary &dict)
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:58
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
virtual word name() const
Name.
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
virtual const volScalarField & rho() const =0
Return the density field.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
addToRunTimeSelectionTable(nucleationModel, reactionDriven, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength
IOerror FatalIOError
static const char nl
Definition: Ostream.H:266
dictionary dict