IATE.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) 2013-2023 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 "IATE.H"
27 #include "IATEsource.H"
28 #include "fvmDdt.H"
29 #include "fvmDiv.H"
30 #include "fvmSup.H"
31 #include "fvcDdt.H"
32 #include "fvcDiv.H"
33 #include "fvcAverage.H"
34 #include "fvModels.H"
35 #include "fvConstraints.H"
36 #include "mathematicalConstants.H"
37 #include "fundamentalConstants.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace diameterModels
45 {
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
53 
54 Foam::tmp<Foam::volScalarField> Foam::diameterModels::IATE::dsm() const
55 {
56  return max(6/max(kappai_, 6/dMax_), dMin_);
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const dictionary& diameterProperties,
65  const phaseModel& phase
66 )
67 :
68  diameterModel(diameterProperties, phase),
69  kappai_
70  (
71  IOobject
72  (
73  IOobject::groupName("kappai", phase.name()),
74  phase.time().name(),
75  phase.mesh(),
76  IOobject::MUST_READ,
77  IOobject::AUTO_WRITE
78  ),
79  phase.mesh()
80  ),
81  dMax_("dMax", dimLength, diameterProperties),
82  dMin_("dMin", dimLength, diameterProperties),
83  residualAlpha_("residualAlpha", dimless, diameterProperties),
84  d_(IOobject::groupName("d", phase.name()), dsm()),
85  sources_(diameterProperties.lookup("sources"), IATEsource::iNew(*this))
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 {
99  return d_;
100 }
101 
102 
104 {
105  return phase()*kappai_;
106 }
107 
108 
110 {
111  volScalarField alphaAv
112  (
113  max
114  (
115  0.5*fvc::average(phase() + phase().oldTime()),
116  residualAlpha_
117  )
118  );
119 
120  // Initialise the accumulated source term to the dilatation effect
122  (
123  -fvm::SuSp
124  (
125  ((1.0/3.0)/alphaAv)
126  *(
127  (
128  fvc::ddt(phase())
129  + fvc::div(phase().alphaPhi())
130  )
131  - (
132  fvc::ddt(phase(), phase().rho())
133  + fvc::div(phase().alphaRhoPhi())
134  )/phase().rho()
135  ),
136  kappai_
137  )
138  );
139 
140  // Accumulate the run-time selectable sources
141  forAll(sources_, j)
142  {
143  R += sources_[j].R(alphaAv, kappai_);
144  }
145 
146  const Foam::fvModels& fvModels(Foam::fvModels::New(phase().mesh()));
148  (
149  Foam::fvConstraints::New(phase().mesh())
150  );
151 
152  // Construct the interfacial curvature equation
153  fvScalarMatrix kappaiEqn
154  (
155  fvm::ddt(kappai_) + fvm::div(phase().phi(), kappai_)
156  - fvm::Sp(fvc::div(phase().phi()), kappai_)
157  ==
158  R
159  + fvModels.source(kappai_)
160  );
161 
162  kappaiEqn.relax();
163 
164  fvConstraints.constrain(kappaiEqn);
165 
166  kappaiEqn.solve();
167 
168  // Update the Sauter-mean diameter
169  d_ = dsm();
170 }
171 
172 
174 {
176 
177  diameterProperties().lookup("dMax") >> dMax_;
178  diameterProperties().lookup("dMin") >> dMin_;
179 
180  // Re-create all the sources updating number, type and coefficients
182  (
183  diameterProperties().lookup("sources"),
184  IATEsource::iNew(*this)
185  ).transfer(sources_);
186 
187  return true;
188 }
189 
190 
191 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
void transfer(PtrList< T > &)
Transfer the contents of the argument PtrList into this PtrList.
Definition: PtrList.C:189
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:52
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Definition: diameterModel.C:62
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition: IATE.H:69
virtual ~IATE()
Destructor.
Definition: IATE.C:91
virtual void correct()
Correct the model.
Definition: IATE.C:109
virtual tmp< volScalarField > d() const
Get the diameter field.
Definition: IATE.C:97
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Definition: IATE.C:173
virtual tmp< volScalarField > Av() const
Get the surface area per unit volume field.
Definition: IATE.C:103
IATE(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
Definition: IATE.C:63
Class used for the read-construction of.
Definition: IATEsource.H:88
IATE (Interfacial Area Transport Equation) bubble diameter model run-time selectable sources.
Definition: IATEsource.H:54
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Finite volume constraints.
Definition: fvConstraints.H:62
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:604
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
Finite volume models.
Definition: fvModels.H:65
Helper class to manage multi-specie phase properties.
A class for managing temporary objects.
Definition: tmp.H:55
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
Fundamental dimensioned constants.
Area-weighted average a surfaceField creating a volField.
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
addToRunTimeSelectionTable(diameterModel, constant, dictionary)
defineTypeNameAndDebug(constant, 0)
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > average(const SurfaceField< Type > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
const dimensionSet dimless
const dimensionSet dimLength
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47