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-2021 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 {
46  defineTypeNameAndDebug(IATE, 0);
47  addToRunTimeSelectionTable(diameterModel, IATE, dictionary);
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().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_(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 
148  (
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 
173 bool Foam::diameterModels::IATE::read(const dictionary& phaseProperties)
174 {
175  diameterModel::read(phaseProperties);
176 
177  diameterProperties().lookup("dMax") >> dMax_;
178  diameterProperties().lookup("dMin") >> dMin_;
179 
180  // Re-create all the sources updating number, type and coefficients
181  PtrList<IATEsource>
182  (
183  diameterProperties().lookup("sources"),
184  IATEsource::iNew(*this)
185  ).transfer(sources_);
186 
187  return true;
188 }
189 
190 
191 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual tmp< volScalarField > d() const
Get the diameter field.
const dimensionSet dimless
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Fundamental dimensioned constants.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
Area-weighted average a surfaceField creating a volField.
virtual tmp< volScalarField > a() const
Get the surface area per unit volume field.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Calculate the first temporal derivative.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
stressControl lookup("compactNormalStress") >> compactNormalStress
const phaseModel & phase() const
Return the phase.
dynamicFvMesh & mesh
Foam::fvConstraints & fvConstraints
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
virtual ~IATE()
Destructor.
Finite volume models.
Definition: fvModels.H:60
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
phi
Definition: correctPhi.H:3
IATE(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
word timeName
Definition: getTimeIndex.H:3
Calculate the divergence of the given field.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
rho oldTime()
Foam::fvModels & fvModels
Calculate the matrix for the divergence of the given field and flux.
#define R(A, B, C, D, E, F, K, M)
Finite volume constraints.
Definition: fvConstraints.H:57
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
virtual void correct()
Correct the model.
A class for managing temporary objects.
Definition: PtrList.H:53
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844