The filmCloudTransfer fvModel now supports an optional ejection model which
provides transfer of film to cloud by dripping from an inverted surface or
curvature separation:
Class
Foam::filmEjectionModels::dripping
Description
Dripping film to cloud ejection transfer model
On an inverted surface if the film thickness is sufficient to generate a
valid parcel the equivalent mass is removed from the film and transfered to
the cloud as a parcel containing droplets with a diameter obtained from
the specified parcelDistribution.
Usage
Example usage:
\verbatim
filmCloudTransfer
{
type filmCloudTransfer;
libs ("libfilmCloudTransfer.so");
ejection
{
model dripping;
deltaStable 5e-4;
minParticlesPerParcel 10;
parcelDistribution
{
type RosinRammler;
Q 0;
min 1e-3;
max 2e-3;
d 7.5e-05;
n 0.5;
}
}
}
\endverbatim
Class
Foam::filmEjectionModels::BrunDripping
Description
Brun dripping film to cloud ejection transfer model
If the film thickness exceeds the critical value needed to generate one or
more drops, the equivalent mass is removed from the film. The critical film
thickness is calculated from the Rayleigh-Taylor stability analysis of film
flow on an inclined plane by Brun et.al.
Reference:
\verbatim
Brun, P. T., Damiano, A., Rieu, P., Balestra, G., & Gallaire, F. (2015).
Rayleigh-Taylor instability under an inclined plane.
Physics of Fluids (1994-present), 27(8), 084107.
\endverbatim
The diameter of the drops formed are obtained from the local capillary
length multiplied by the \c dCoeff coefficient which defaults to 3.3.
Reference:
\verbatim
Lefebvre, A. (1988).
Atomisation and sprays
(Vol. 1040, No. 2756). CRC press.
\endverbatim
Usage
Example usage:
\verbatim
filmCloudTransfer
{
type filmCloudTransfer;
libs ("libfilmCloudTransfer.so");
ejection
{
model BrunDripping;
deltaStable 5e-4;
}
}
\endverbatim
Class
Foam::filmEjectionModels::curvatureSeparation
Description
Curvature induced separation film to cloud ejection transfer model
Assesses film curvature via the mesh geometry and calculates a force
balance of the form:
F_sum = F_inertial + F_body + F_surface_tension
If F_sum < 0, the film separates and is transferred to the cloud
if F_sum >= 0 the film remains attached.
Reference:
\verbatim
Owen, I., & Ryley, D. J. (1985).
The flow of thin liquid films around corners.
International journal of multiphase flow, 11(1), 51-62.
\endverbatim
Usage
Example usage:
\verbatim
filmCloudTransfer
{
type filmCloudTransfer;
libs ("libfilmCloudTransfer.so");
ejection
{
model curvatureSeparation;
deltaStable 5e-4;
}
}
\endverbatim
The new tutorials/modules/multiRegion/film/cylinderDripping tutorial case
demonstrates a film dripping into the cloud. The standard cylinder case is
turned upside-down (by changing the orientation of gravity) with an initial
0.2mm film of water over the surface which drips when the thickness is greater
than 0.5mm. Settings for all three ejection models are provided in the
constant/film/fvModels dictionary with the standard dripping model selected.
153 lines
4.1 KiB
C++
153 lines
4.1 KiB
C++
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration | Website: https://openfoam.org
|
|
\\ / A nd | Copyright (C) 2023 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 "isothermalFilm.H"
|
|
#include "surfaceTensionModel.H"
|
|
#include "fvmDiv.H"
|
|
#include "fvcSnGrad.H"
|
|
#include "fvcLaplacian.H"
|
|
#include "fvcReconstruct.H"
|
|
|
|
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
|
|
|
|
Foam::tmp<Foam::volScalarField>
|
|
Foam::solvers::isothermalFilm::sigma() const
|
|
{
|
|
return constrainedField(surfaceTension->sigma());
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::surfaceScalarField>
|
|
Foam::solvers::isothermalFilm::pbByAlphaRhof() const
|
|
{
|
|
return fvc::interpolate
|
|
(
|
|
max(nHat & g, dimensionedScalar(g.dimensions(), 0))*VbyA
|
|
);
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::surfaceScalarField>
|
|
Foam::solvers::isothermalFilm::pbByAlphaf() const
|
|
{
|
|
return fvc::interpolate(rho)*pbByAlphaRhof();
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::surfaceScalarField>
|
|
Foam::solvers::isothermalFilm::pbByAlphaGradRhof() const
|
|
{
|
|
return pbByAlphaRhof()*fvc::snGrad(rho);
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField>
|
|
Foam::solvers::isothermalFilm::pc(const volScalarField& sigma) const
|
|
{
|
|
return -fvc::laplacian(sigma, delta);
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField>
|
|
Foam::solvers::isothermalFilm::pe() const
|
|
{
|
|
// Update the pressure, mapping from the fluid region as required
|
|
p.correctBoundaryConditions();
|
|
|
|
// Add the droplet impingement pressure
|
|
p.ref() += mesh.time().deltaT()*fvModels().source(p, "pi")().Su();
|
|
|
|
return p;
|
|
}
|
|
|
|
|
|
void Foam::solvers::isothermalFilm::momentumPredictor()
|
|
{
|
|
volVectorField& U = U_;
|
|
|
|
// Calculate the surface tension coefficient
|
|
const volScalarField sigma(this->sigma());
|
|
|
|
tUEqn =
|
|
(
|
|
fvm::ddt(alpha, rho, U) + fvm::div(alphaRhoPhi, U)
|
|
- fvm::Sp(contErr(), U)
|
|
+ momentumTransport->divDevTau(U)
|
|
==
|
|
contactForce(sigma)
|
|
+ fvModels().source(alpha, rho, U)
|
|
);
|
|
fvVectorMatrix& UEqn = tUEqn.ref();
|
|
|
|
// Thermocapillary force
|
|
if (thermocapillary)
|
|
{
|
|
UEqn -= fvc::grad(sigma)/VbyA;
|
|
}
|
|
|
|
UEqn.relax();
|
|
|
|
fvConstraints().constrain(UEqn);
|
|
|
|
if (pimple.momentumPredictor())
|
|
{
|
|
const surfaceScalarField alphaf(fvc::interpolate(alpha));
|
|
|
|
solve
|
|
(
|
|
UEqn
|
|
==
|
|
fvc::reconstruct
|
|
(
|
|
constrainedField
|
|
(
|
|
alphaf
|
|
*(
|
|
// Buoyancy force
|
|
fvc::interpolate(rho)*(g & mesh.Sf())
|
|
|
|
- (
|
|
// External and capillary pressure
|
|
fvc::snGrad(pe() + pc(sigma), "snGrad(p)")
|
|
|
|
// Buoyant pressure
|
|
+ pbByAlphaGradRhof()*alphaf
|
|
+ pbByAlphaf()*fvc::snGrad(alpha)
|
|
)*mesh.magSf()
|
|
)
|
|
)
|
|
)
|
|
);
|
|
|
|
// Remove film-normal component of velocity
|
|
U -= nHat*(nHat & U);
|
|
|
|
U.correctBoundaryConditions();
|
|
}
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|