1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2013-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "IATE.H"
#include "IATEsource.H"
#include "twoPhaseSystem.H"
#include "fvmDdt.H"
#include "fvmDiv.H"
#include "fvmSup.H"
#include "fvcDdt.H"
#include "fvcDiv.H"
#include "fvcAverage.H"
#include "fvOptions.H"
#include "mathematicalConstants.H"
#include "fundamentalConstants.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace diameterModels
{
defineTypeNameAndDebug(IATE, 0);
addToRunTimeSelectionTable
(
diameterModel,
IATE,
dictionary
);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::diameterModels::IATE::IATE
(
const dictionary& diameterProperties,
const phaseModel& phase
)
:
diameterModel(diameterProperties, phase),
kappai_
(
IOobject
(
IOobject::groupName("kappai", phase.name()),
phase_.U().time().timeName(),
phase_.U().mesh(),
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
phase_.U().mesh()
),
dMax_("dMax", dimLength, diameterProperties_),
dMin_("dMin", dimLength, diameterProperties_),
residualAlpha_
(
"residualAlpha",
dimless,
diameterProperties_
),
d_
(
IOobject
(
IOobject::groupName("d", phase.name()),
phase_.U().time().timeName(),
phase_.U().mesh(),
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
dsm()
),
sources_
(
diameterProperties_.lookup("sources"),
IATEsource::iNew(*this)
)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::diameterModels::IATE::~IATE()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp<Foam::volScalarField> Foam::diameterModels::IATE::dsm() const
{
return max(6/max(kappai_, 6/dMax_), dMin_);
}
// Placeholder for the nucleation/condensation model
// Foam::tmp<Foam::volScalarField> Foam::diameterModels::IATE::Rph() const
// {
// const volScalarField& T = phase_.thermo().T();
// const volScalarField& p = phase_.thermo().p();
//
// scalar A, B, C, sigma, vm, Rph;
//
// volScalarField ps(1e5*pow(10, A - B/(T + C)));
// volScalarField Dbc
// (
// 4*sigma*vm/(constant::physicoChemical::k*T*log(p/ps))
// );
//
// return constant::mathematical::pi*sqr(Dbc)*Rph;
// }
void Foam::diameterModels::IATE::correct()
{
// Initialise the accumulated source term to the dilatation effect
volScalarField R
(
(
(1.0/3.0)
/max
(
fvc::average(phase_ + phase_.oldTime()),
residualAlpha_
)
)
*(fvc::ddt(phase_) + fvc::div(phase_.alphaPhi()))
);
// Accumulate the run-time selectable sources
forAll(sources_, j)
{
R -= sources_[j].R();
}
fv::options& fvOptions(fv::options::New(phase_.mesh()));
// Construct the interfacial curvature equation
fvScalarMatrix kappaiEqn
(
fvm::ddt(kappai_) + fvm::div(phase_.phi(), kappai_)
- fvm::Sp(fvc::div(phase_.phi()), kappai_)
==
- fvm::SuSp(R, kappai_)
//+ Rph() // Omit the nucleation/condensation term
+ fvOptions(kappai_)
);
kappaiEqn.relax();
fvOptions.constrain(kappaiEqn);
kappaiEqn.solve();
// Update the Sauter-mean diameter
d_ = dsm();
}
bool Foam::diameterModels::IATE::read(const dictionary& phaseProperties)
{
diameterModel::read(phaseProperties);
diameterProperties_.lookup("dMax") >> dMax_;
diameterProperties_.lookup("dMin") >> dMin_;
// Re-create all the sources updating number, type and coefficients
PtrList<IATEsource>
(
diameterProperties_.lookup("sources"),
IATEsource::iNew(*this)
).transfer(sources_);
return true;
}
// ************************************************************************* //
|