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-2015 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "IATE.H"
29#include "IATEsource.H"
30#include "twoPhaseSystem.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"
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44namespace Foam
45{
46namespace diameterModels
47{
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_.U().time().timeName(),
75 phase_.U().mesh(),
76 IOobject::MUST_READ,
77 IOobject::AUTO_WRITE
78 ),
79 phase_.U().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_.U().time().timeName(),
90 phase_.U().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
112Foam::tmp<Foam::volScalarField> Foam::diameterModels::IATE::dsm() const
113{
114 return max(6/max(kappai_, 6/dMax_), dMin_);
115}
116
117// Placeholder for the nucleation/condensation model
118// Foam::tmp<Foam::volScalarField> Foam::diameterModels::IATE::Rph() const
119// {
120// const volScalarField& T = phase_.thermo().T();
121// const volScalarField& p = phase_.thermo().p();
122//
123// scalar A, B, C, sigma, vm, Rph;
124//
125// volScalarField ps(1e5*pow(10, A - B/(T + C)));
126// volScalarField Dbc
127// (
128// 4*sigma*vm/(constant::physicoChemical::k*T*log(p/ps))
129// );
130//
131// return constant::mathematical::pi*sqr(Dbc)*Rph;
132// }
133
135{
136 // Initialise the accumulated source term to the dilatation effect
138 (
139 (
140 (1.0/3.0)
141 /max
142 (
143 fvc::average(phase_ + phase_.oldTime()),
144 residualAlpha_
145 )
146 )
147 *(fvc::ddt(phase_) + fvc::div(phase_.alphaPhi()))
148 );
149
150 // Accumulate the run-time selectable sources
151 forAll(sources_, j)
152 {
153 R -= sources_[j].R();
154 }
155
156 fv::options& fvOptions(fv::options::New(phase_.mesh()));
157
158 // Construct the interfacial curvature equation
159 fvScalarMatrix kappaiEqn
160 (
161 fvm::ddt(kappai_) + fvm::div(phase_.phi(), kappai_)
162 - fvm::Sp(fvc::div(phase_.phi()), kappai_)
163 ==
164 - fvm::SuSp(R, kappai_)
165 //+ Rph() // Omit the nucleation/condensation term
166 + fvOptions(kappai_)
167 );
168
169 kappaiEqn.relax();
170
171 fvOptions.constrain(kappaiEqn);
172
173 kappaiEqn.solve();
174
175 // Update the Sauter-mean diameter
176 d_ = dsm();
177}
178
179
180bool Foam::diameterModels::IATE::read(const dictionary& phaseProperties)
181{
182 diameterModel::read(phaseProperties);
183
184 diameterProperties_.readEntry("dMax", dMax_);
185 diameterProperties_.readEntry("dMin", dMin_);
186
187 // Re-create all the sources updating number, type and coefficients
188 PtrList<IATEsource>
189 (
190 diameterProperties_.lookup("sources"),
191 IATEsource::iNew(*this)
192 ).transfer(sources_);
193
194 return true;
195}
196
197
198// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
fv::options & fvOptions
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition: IATE.H:71
virtual ~IATE()
Destructor.
Definition: IATE.C:106
virtual void correct()
Correct the diameter field.
Definition: IATE.C:117
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
dynamicFvMesh & 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 finiteVolume matrix for implicit and explicit sources.
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:48
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
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.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333