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-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 "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 "fvOptions.H"
35 #include "mathematicalConstants.H"
36 #include "fundamentalConstants.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace diameterModels
44 {
45  defineTypeNameAndDebug(IATE, 0);
46  addToRunTimeSelectionTable(diameterModel, IATE, dictionary);
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52 
54 {
55  return d_;
56 }
57 
58 
60 {
61  return phase()*kappai_;
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66 
68 (
69  const dictionary& diameterProperties,
70  const phaseModel& phase
71 )
72 :
73  diameterModel(diameterProperties, phase),
74  kappai_
75  (
76  IOobject
77  (
78  IOobject::groupName("kappai", phase.name()),
79  phase.time().timeName(),
80  phase.mesh(),
81  IOobject::MUST_READ,
82  IOobject::AUTO_WRITE
83  ),
84  phase.mesh()
85  ),
86  dMax_("dMax", dimLength, diameterProperties),
87  dMin_("dMin", dimLength, diameterProperties),
88  residualAlpha_("residualAlpha", dimless, diameterProperties),
89  d_(dRef()),
90  sources_(diameterProperties.lookup("sources"), IATEsource::iNew(*this))
91 {
92  d_ = dsm();
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
104 Foam::tmp<Foam::volScalarField> Foam::diameterModels::IATE::dsm() const
105 {
106  return max(6/max(kappai_, 6/dMax_), dMin_);
107 }
108 
109 
111 {
112  volScalarField alphaAv
113  (
114  max
115  (
116  0.5*fvc::average(phase() + phase().oldTime()),
117  residualAlpha_
118  )
119  );
120 
121  // Initialise the accumulated source term to the dilatation effect
123  (
124  -fvm::SuSp
125  (
126  ((1.0/3.0)/alphaAv)
127  *(
128  (
129  fvc::ddt(phase())
130  + fvc::div(phase().alphaPhi())
131  )
132  - (
133  fvc::ddt(phase(), phase().rho()())
134  + fvc::div(phase().alphaRhoPhi())
135  )/phase().rho()
136  ),
137  kappai_
138  )
139  );
140 
141  // Accumulate the run-time selectable sources
142  forAll(sources_, j)
143  {
144  R += sources_[j].R(alphaAv, kappai_);
145  }
146 
147  fv::options& fvOptions(fv::options::New(phase().mesh()));
148 
149  // Construct the interfacial curvature equation
150  fvScalarMatrix kappaiEqn
151  (
152  fvm::ddt(kappai_) + fvm::div(phase().phi(), kappai_)
153  - fvm::Sp(fvc::div(phase().phi()), kappai_)
154  ==
155  R
156  + fvOptions(kappai_)
157  );
158 
159  kappaiEqn.relax();
160 
161  fvOptions.constrain(kappaiEqn);
162 
163  kappaiEqn.solve();
164 
165  // Update the Sauter-mean diameter
166  d_ = dsm();
167 }
168 
169 
170 bool Foam::diameterModels::IATE::read(const dictionary& phaseProperties)
171 {
172  diameterModel::read(phaseProperties);
173 
174  diameterProperties().lookup("dMax") >> dMax_;
175  diameterProperties().lookup("dMin") >> dMin_;
176 
177  // Re-create all the sources updating number, type and coefficients
178  PtrList<IATEsource>
179  (
180  diameterProperties().lookup("sources"),
181  IATEsource::iNew(*this)
182  ).transfer(sources_);
183 
184  return true;
185 }
186 
187 
188 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
rho oldTime()
#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 > &)
fv::options & fvOptions
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
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Fundamental dimensioned constants.
Macros for easy insertion into run-time selection tables.
Area-weighted average a surfaceField creating a volField.
phi
Definition: pEqn.H:104
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
virtual tmp< volScalarField > calcA() const
Get the surface area per unit volume field.
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
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
virtual ~IATE()
Destructor.
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))
IATE(const dictionary &diameterProperties, const phaseModel &phase)
Construct from components.
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
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Calculate the matrix for the divergence of the given field and flux.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define R(A, B, C, D, E, F, K, M)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.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:46
virtual void correct()
Correct the diameter field.
virtual tmp< volScalarField > calcD() const
Get the diameter field.
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
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:812