IATE.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) 2013-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 "IATE.H"
30 #include "IATEsource.H"
31 #include "fvmDdt.H"
32 #include "fvmDiv.H"
33 #include "fvmSup.H"
34 #include "fvcDdt.H"
35 #include "fvcDiv.H"
36 #include "fvcAverage.H"
37 #include "fvOptions.H"
38 #include "mathematicalConstants.H"
39 #include "fundamentalConstants.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace diameterModels
47 {
48  defineTypeNameAndDebug(IATE, 0);
49 
51  (
52  diameterModel,
53  IATE,
54  dictionary
55  );
56 }
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().timeName(),
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_
85  (
86  IOobject
87  (
88  IOobject::groupName("d", phase.name()),
89  phase_.time().timeName(),
90  phase_.mesh(),
91  IOobject::NO_READ,
92  IOobject::AUTO_WRITE
93  ),
94  dsm()
95  ),
96  sources_
97  (
98  diameterProperties_.lookup("sources"),
99  IATEsource::iNew(*this)
100  )
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
105 
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
112 Foam::tmp<Foam::volScalarField> Foam::diameterModels::IATE::dsm() const
113 {
114  return max(6/max(kappai_, 6/dMax_), dMin_);
115 }
116 
118 {
119  volScalarField alphaAv
120  (
121  max
122  (
123  0.5*fvc::average(phase_ + phase_.oldTime()),
124  residualAlpha_
125  )
126  );
127 
128  // Initialise the accumulated source term to the dilatation effect
130  (
131  -fvm::SuSp
132  (
133  ((1.0/3.0)/alphaAv)
134  *(
135  fvc::ddt(phase_) + fvc::div(phase_.alphaPhi())
136  - phase_.continuityError()/phase_.rho()
137  ),
138  kappai_
139  )
140  );
141 
142  // Accumulate the run-time selectable sources
143  forAll(sources_, j)
144  {
145  R += sources_[j].R(alphaAv, kappai_);
146  }
147 
148  fv::options& fvOptions(fv::options::New(phase_.mesh()));
149 
150  // Construct the interfacial curvature equation
151  fvScalarMatrix kappaiEqn
152  (
153  fvm::ddt(kappai_) + fvm::div(phase_.phi(), kappai_)
154  - fvm::Sp(fvc::div(phase_.phi()), kappai_)
155  ==
156  R
157  + fvOptions(kappai_)
158  );
159 
160  kappaiEqn.relax();
161 
162  fvOptions.constrain(kappaiEqn);
163 
164  kappaiEqn.solve();
165 
166  // Update the Sauter-mean diameter
167  d_ = dsm();
168 }
169 
170 
172 {
174 
175  diameterProperties_.readEntry("dMax", dMax_);
176  diameterProperties_.readEntry("dMin", dMin_);
177 
178  // Re-create all the sources updating number, type and coefficients
180  (
181  diameterProperties_.lookup("sources"),
182  IATEsource::iNew(*this)
183  ).transfer(sources_);
184 
185  return true;
186 }
187 
188 
189 // ************************************************************************* //
addToRunTimeSelectionTable(diameterModel, constant, dictionary)
virtual bool read(const dictionary &phaseProperties)=0
Read phaseProperties dictionary.
Definition: diameterModel.C:89
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Helper class to manage multi-specie phase properties.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
IATE (Interfacial Area Transport Equation) bubble diameter model run-time selectable sources...
Definition: IATEsource.H:52
defineTypeNameAndDebug(constant, 0)
const dimensionSet dimless
Dimensionless.
Finite-volume options.
Definition: fvOptions.H:51
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Fundamental dimensioned constants.
Macros for easy insertion into run-time selection tables.
fv::options & fvOptions
Area-weighted average a surfaceField creating a volField.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:40
Calculate the first temporal derivative.
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual void correct()
Correct the diameter field.
Definition: IATE.C:110
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Definition: IATE.C:164
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:49
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
IATE(const dictionary &diameterProperties, const phaseModel &phase)
Construct from components.
Definition: IATE.C:56
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Calculate the divergence of the given field.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1101
Class used for the read-construction of.
Definition: IATEsource.H:93
Calculate the matrix for the divergence of the given field and flux.
virtual ~IATE()
Destructor.
Definition: IATE.C:99
#define R(A, B, C, D, E, F, K, M)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:50
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:41
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.