Air Entrainment Model in interFoam
Air Entrainment Model in interFoam
In Proceedings of
CFD with OpenSource Software, 2018, Edited by Nilsson. H.,
[Link]
Implementation of an air-entrainment
model in interFoam
Author:
Silje Kreken Almeland Peer reviewed by:
Norwegian University of Ebrahim Ghahramani
Science and Technology Sarmakeeva Anastasi
[Link]@[Link]
Disclaimer: This is a student project work, done as part of a course where OpenFOAM and some
other OpenSource software are introduced to the students. Any reader should be aware that it
might not be free of errors. Still, it might be useful for someone who would like learn some details
similar to the ones presented in the report and in the accompanying files. The material has gone
through a review process. The role of the reviewer is to go through the tutorial and make sure that
it works, that it is possible to follow, and to some extent correct the writing. The reviewer has no
responsibility for the contents.
January 2, 2019
Learning outcomes
How it is implemented:
• About the design of the interFoam solver
• How the solution of the α-equation is implemented
How to modify it:
• How to modify the interFoam solver by implementing the described air entrainment model
• How to modify parts of the model that could improve its performance
1
Prerequisites
The reader is expected to know the following in order to get maximum benefit out of this report:
• Have basic knowledge of coputational fluid dynamics and the equations involved
2
Contents
1 Introduction 6
2 Theory 7
2.1 interFoam . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 7
2.1.1 Solution algorithm . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 8
2.1.2 interAirEntFoam . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 9
2.1.3 Mass Continuity and buoyancy . . . . . . . . . . . . . . . . . . . . . . . . . . 9
3 Design of interFoam 11
3.1 Calculation of αl - alphaEqn.H . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 13
3.1.1 Calculation scheme with MULES . . . . . . . . . . . . . . . . . . . . . . . . . 15
3.1.2 Calculation scheme without MULES . . . . . . . . . . . . . . . . . . . . . . . 18
4 Implementation of interAirEntFoam 20
4.1 Preparations . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 20
4.1.1 Make a local copy of interFoam . . . . . . . . . . . . . . . . . . . . . . . . . 20
4.1.2 Copy immiscibleIncompressibleTwoPhaseMixture . . . . . . . . . . . . . . 22
4.2 Implementations . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 25
4.2.1 Pt (k, ρ) . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 26
4.2.2 Pd (LT , k, g, ) . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 29
4.2.3 Cair and As . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 32
4.2.4 CalcSource . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 33
4.2.5 Adding the source term to the αl -equation . . . . . . . . . . . . . . . . . . . 35
4.2.6 Write additional fields for visualization . . . . . . . . . . . . . . . . . . . . . . 38
5 Verification cases 41
5.1 Vertical plunging jet . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 41
5.2 Possible improvements . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 43
3
Listings
4
LISTINGS LISTINGS
5
Chapter 1
Introduction
The aim of this tutorial is to show how a sub-grid model to predict the amount of air entrainment
in an air-water system can be implemented based on interFoam. In the framework of a volume of
fluid model, the intention of the model is to describe a time-averaging solution for air-concentration
values attained from self-aeration of the flow. Hirt (2003) modelled air-entrainment with the VOF
method, by including an additional air volume to cells in which some criteria for air entrainment
were satisfied. This code is implemented into Flow3D, which is not open source. Therefore, the
exact algorithm of this method is not known. Nevertheless, the expression of the source term is
given. The implementation shown in this tutorial will represent a try on reproducing this model.
Possible implementations of this term into interFoam will be discussed. In addition to the model of
Hirt (2003) other models have been published in the literature. An aim for this project will be to
make the code general so that it will be easy to use different air entrainment models on this base.
The implemented air entrainment model is developed based on interFoam. Tutorials describing
interFoam has been written by students in the OSCFD course a couple of times before. By Hassan
Hemida in 2008 (Hemida, 2008), describing general use of the solver, and giving a brief explanation
of equations that are solved. In 2017 Elin Olsson (Olsson, 2017) described the theory of the VoF
method. She also explained the concept of the MULES scheme. Nevertheless, the focus of her project
was on the solver interIsoFoam. Alessandro Manni described the theory and some implementations
of the MULES algorithm in his twoPhaseEulerFoam tutorial from 2014 (Manni, 2014). Also Surya
Kaundinya Oruganti (Oruganti, 2015) touched the theory and implementation of MULES when
describing the implementation of cavitation models in multiphaseEulerFoam.
In this tutorial the intention is to give an overview of the solution algorithm of interFoam, along
with a description of the design of the solver. Due to the fact that the implementation of the air
entrainment model will affect the solution of the α-field, a description of how the α-equation solved
and how the MULES routine is implemented are included. The theory of MULES is briefly described,
while the effect of the user variables on the calculation routine are emphasized. Then the theory of
the implemented air entrainment model is outlined, and its implementation are described. Finally,
some verification cases are given, before some comments are added regarding the performance of the
implemented solver.
6
Chapter 2
Theory
When the turbulence is high enough in free-surface flow, air will be entrained in the water phase
and transported with the water flow. The length scales involved in the air entrainment process are
in the range of millimeters. This means that a very fine grid is needed to capture these processes.
When the grid is not fine enough, the processes at the surface will not be captured, and the amount
of air within the water phase will be underestimated. Therfore, to be able to predict the amount of
air entrained into the water flow using a reasonable grid resolution, a subgrid model will be needed
to model the dynamics of the air entrainment process.
As described in Section 1 the aim of this tutorial is to describe how to implement the air en-
trainment model by Hirt (2003) into the framework of interFoam. In this tutorial the subgrid
model is implemented in a source term, which is added as an explicit source term to the αl equation.
This means that the solution routines of interFoam are used for the new solver as well. To get
these implementations into a the relevant context, a brief explanation of the interFoam solver is
given next, followed by the theory of the air entrainment model. To supplement the description of
interFoam given in the following section, the interested reader should consult the tutorials referred
to in Chapter 1.
2.1 interFoam
interFoam is a two-phase solver for incompressible, isothermal
and immiscible fluids, using the volume of fluid concept to cap-
ture the interface. It solves a single set of mass- and momentum
equations
∇· u = 0 (2.1)
7
2.1. INTERFOAM CHAPTER 2. THEORY
∂αl
+ ∇· (αl u) + ∇· [uc αl (1 − αl )] = 0 (2.4)
∂t
using the MULES scheme, which will be further outlined below. In Equation (2.4) αl is the volume
fraction of water, the corresponding volume fraction of air is
αair = 1 − αl (2.5)
uc is a relative velocity that acts as a velocity perpendicular to the interface. This term is meant to
estimate the relative velocity between the phases, in interFoam only one velocity field is available,
and the relative velocity is estimated according to
∇αl
uc = Cα |u| (2.6)
|∇αl |
∇αl
where Cα is a user variable, which determines the effect of the compression term. |∇α l|
is meant to
give the unit vector normal to the interface. The interface is localized at the cells where 0 < αl < 1.
An illustration of the calulation of the free surface is given in Figure 2.1. In the interface cells, ρ
and µ are weighted values given by
ρ = αl ρl + (1 − αl )ρair (2.7a)
µ = αl µl + (1 − αl )µair (2.7b)
where ρl , ρair and µl , µair are the density and viscosity in the water and air phase respectively.
The volumetric surface tension force is estimated by the Continuum Surface Force model (Brackbill,
Kothe, & Zemach, 1992).
f = σκ∇αl (2.8)
where σ is the surface tension and κ is the surface curvature calculated as
∇αl
κ= (2.9)
|∇αl |
interFoam can be run both in laminar and turbulent mode. And the turbulence can be solved using
the Reynolds average (RAS) framework, or based on large eddy simulations (LES).
• Initiate fields
8
2.1. INTERFOAM CHAPTER 2. THEORY
The VOF routine could in principle be used to solve all the scales in the air entrainment process in a
direct numerical simulation scheme (DNS). However, it would exceed the computational capacities
in the forseeable future. The problem with interface capturing methods is that they cannot solve
the dynamics with length scales smaller then the grid size, and this has lead to the development of
subgrid models for air entrainment. The solver interAirEntFoam, implemented in this tutorial, and
described in the following section, constitutes an attempt to implement one such subgrid model.
2.1.2 interAirEntFoam
One of the first attempts on calculating the mean flow characteristics and predict the air-entrainment
with a sub-grid model was introduced in FLOW-3D (Hirt, 2003). In this model air was entrained
into the water when turbulent energy per unit volume was larger than the stabilizing forces of gravity
and surface tension per unit volume. The expression for this source term is given as
s
δV Pt − Pd
= Cair As 2 (2.10)
δt ρ
where Cair is a calibration parameter, As the free surface area at each cell, ρ is the fluid density,
and Pt is expressed by
Pt = ρk (2.11)
and constitute the turbulent forces. Pd represents the stabilizing forces and is expressed by
σ
Pd = ρgn LT + (2.12)
LT
where gn is the gravity component normal to the free surface, σ the liquid-gas surface tension, and
the turbulence characteristic length scale LT is calculated according to
r
3 k 3/2
LT = Cµ (2.13)
2
where Cµ = 0.09 for the standard k turbulence model. The length scale LT should be taken as an
approximation to a length scale of perturbations. Pd is an expression for the energy density [J/m3 ]
associated with a disturbed fluid element raised over the free surface to a height LT . Pd represents a
surface stabilizing force, while the turbulent kinetic energy per unit volume, Pt , represents the per-
turbing component that makes the flow unstable. As given in Equation (2.10), when the perturbing
forces Pt exceed the stabilizing forces Pd , an air volume according to Equation (2.10) is added to
the cell in question.
9
2.1. INTERFOAM CHAPTER 2. THEORY
the calculation of the density and the fluxes. Nevertheless, in this tutorial the additional air term is
added to the αl equation as an explicit source term. Boundary effects and relative velocity between
the phases are given as features in the original model, these effects are not treated in the model
implemented in this tutorial.
10
Chapter 3
Design of interFoam
The structure of the solution algorithm in interFoam is as given in Section 2.1. Nevertheless, to be
able to implement a general solver according to the structure of OpenFOAM, an inspection into the
design of the solver is needed. The following tasks has to be done by the solver:
• Solution algorithm
• createFiles.H
• interFoam.C
• pEqn.H
• UEqn.H
• alphaEqn.H
• alphaEqnSubCycle.H
To start the solution procedure, the input parameters have to be read from the input files, and
stored in a proper manner. Also the solution variables, and output variables have to be handled.
The phase properties and turbulence properties needed in the calculation, are declared and initialized
createFields.H.
interFoam as an incompressible, immiscible, two-phase model, utilizes a physical phase model
named immiscibleIncompressibleTwoPhaseMixture to initialize the transport properties. This is
done in createFiles.H by code line 35
11
CHAPTER 3. DESIGN OF INTERFOAM
The model takes care of parsing and updating of the transport properties, like density and vis-
cosity, as well as the interface properties, during the calculation procedure of interFoam. immis-
cibleIncompressibleTwoPhaseMixture is a virtual subclass of incompressibleTwoPhaseMixture and
interfaceProperties. This means that the functions implemented in immiscibleIncompressibleT-
woPhaseMixture only calls functions in the parents classes with the same name, and apart from
that, no new function implementations are done for this class. For example the correct function
within immiscibleIncompressibleTwoPhaseMixture.H (located in $FOAM SRC/transportModels/
immiscibleIncompressibleTwoPhaseMixture) calls the respective correct-functions of incom-
pressibleTwoPhaseMixture and interfaceProperties as shown in the code listing
here it can be seen that the transport model is needed as an argument to create the turbulence
model. Towards the end of createFields.H the gravitational force is parsed from the input file,
and initiated as a variable. This is done within the file readGravitationalAcceleration.H (located
in $FOAM SRC/finiteVolume/cfdTools/general/include/).
12
3.1. CALCULATION OF αL - ALPHAEQN.H CHAPTER 3. DESIGN OF INTERFOAM
12 );
In addition to the phase- and mixture properties, the solver needs some input on the solver vari-
ables. The variables that are connected with the behavior of the solution algorithm are given in
system/fvSolution. These variables are parsed by different files depending on their role. For exam-
ple the parsing of the parameters releated to the PIMPLE routine is done in pimpleControl.H and
variables related solution of the αl -equation is parsed in alphaControl.H. These files are executed
from interFoam.C.
When solution variables, phases, and transport- and turbulence models are initiated, the solu-
tion algorithm can proceed. interFoam.C contains the solution algorithm as outlined in Section 2.1.
UEqn.H calculates the momentum equation, while the pressure field is solved in pEqn.H. The solution
of the αl equation is taken care of by alphaEqn.H and alphaEqnSubCycle.H. In the current imple-
mentation of the air entrainment model, a source term is added to the αl -equation, which makes its
solution routine relevant for the implementation of the air entrainment model, interAirEntFoam.
The αl -quation is solved using MULES compression. The solution procedure is described in the
following section.
The solution algorithm in alphaEqn.H solves Equation (3.1) successively using MULES compression.
MULES is an iterative implementation of the Flux Corrected Transport technique (FCT), used to
guarantee boundedness in the solution of hyperbolic problems. In this technique a corrected flux
FC , which is the weighted solution between the solution obtained with a high- and a low order
discretization scheme, is applied
FC = FL + λA = FL + λ(FH − FL ) (3.2)
here FL and FH denote the solution of the numerical equation, obtained with a low- and a high
order discretization scheme respectively. A = FH − FL is the anti-diffusive flux, and λ is a weighting
factor that limits the solution within a defined range. In the framework of interFoam, Equation
(3.2) is expressed as
lim,Corr
φα,l = φU D
α,l + φα,l = φU D Corr
α,l + θl · φα,l (3.3)
where φU D
α,l is the flux calculated by solving the alpha equation (3.1) without the compression term
∂αl
+ ∇ · (αl u) = 0 (3.4)
∂t
using the first order upwind discretization scheme for the convective term. And
φlim,Corr
α,l = θl · φCorr
α,l (3.5)
where θl is a weighting factor that limits the solutions within given maximum and minimum values.
Further the corrected flux
φCorr
α,l = φHO C UD
α,l + φα,l − φα,l (3.6)
where φHO
α,l is the higher order solution of Equation (3.4), and φCα,l is the correction flux calculated
as a compression term. The first two terms in Equation (3.6), φHO
α,l + φCα,l , can be recognized as FH
13
3.1. CALCULATION OF αL - ALPHAEQN.H CHAPTER 3. DESIGN OF INTERFOAM
φ linearInterpolate(U) & Sf uT Sf
φC = Cα · n̂ = Cα · n̂ = Cα · n̂ (3.8)
|Sf | |Sf | |Sf |
where Sf is the face area vector, φc is the compression flux, Cα is a user defined parameter, available
through cAlpha in the input file system/fvSolution, and n̂ is the interface normal flux
∇αl
n̂ = −8 · Sf (3.9)
|∇αl | + 10 1
V 3
were ∇αl is the gradient of αl and V is the cell volume. As evident from the notation in the listings
above, n̂ is calculated within the transport model. For the immiscibleIncompressibleTwoPhaseMix-
ture-model this is done within interfaceProperties.C. The relevant code lines are shown in the
listing below
14
3.1. CALCULATION OF αL - ALPHAEQN.H CHAPTER 3. DESIGN OF INTERFOAM
168 deltaN_
169 (
170 "deltaN",
171 1e-8/pow(average([Link]().V()), 1.0/3.0)
172 ),
The calculation of the αl fluxes in interFoam is, as mentioned, done successively resulting in
fluxes given by Equation (3.3), (3.5)-(3.6). Written in a more compact form
φα,l = φU D HO C UD
α,l + θl (φα,l + φα,l − φα,l ) (3.10)
Nevertheless, different user settings (supplied in system/fvSolutions) make Equation (3.10) look
different. The parameters that decides on the solutions procedure of the αl -equation, are found in the
directory solvers and [Link] in system/fvSolution. There is several parameters affecting
the solution of the αl field, the most important regarding the solution procedure, being MULESCorr
and cAlpha. Different combinations of values for these parameter will change the form of Equation
(3.10). Through these parameters it is possible to turn on and off the MULES compression, and
also to decide on the the extent of the effects of the compression term. How these parameters can
be set to achieve different behavior for the solution procedure will be given in the following sections.
This means that as long as one do not set MULESCorr in system/fvSolution, MULES-correction is
not applied. To apply MULES correction, the user parameter MULESCorr has to be set to true or
yes. In this section the possible variants within the MULES-correction path are investigated.
Using the MULES algorithm, the fluxes of αl is calculated according to Equation (3.10)
φα,l = φU D HO C UD
α,l + θl (φα,l + φα,l − φα,l )
99 if (MULESCorr)
100 {
101 #include "alphaSuSp.H"
102
15
3.1. CALCULATION OF αL - ALPHAEQN.H CHAPTER 3. DESIGN OF INTERFOAM
109 )
110 + fv::gaussConvectionScheme<scalar>
111 (
112 mesh,
113 phiCN,
114 upwind<scalar>(mesh, phiCN)
115 ).fvmDiv(phiCN, alpha1)
116 // - fvm::Sp(fvc::ddt(dimensionedScalar("1", dimless, 1), mesh)
117 // + fvc::div(phiCN), alpha1)
118 ==
119 Su + fvm::Sp(Sp + divU, alpha1)
120 );
121
122 [Link]();
After solving the αl field by the upwind scheme, the αg -field is updated
146 [Link]();
Next the higher order fluxes and the correction fluxes are stored in talphaPhi1Un
16
3.1. CALCULATION OF αL - ALPHAEQN.H CHAPTER 3. DESIGN OF INTERFOAM
168 alpharScheme
169 )
170 );
Then we have the following situation before entering the MULESCorr-loop at line 172
alphaPhi10 = φU
l
D
(3.11)
talphaPhi1Un = φHO
l + φC
l
talphaPhi1Corr = φHO
l + φC UD
l − φl (3.12)
which in the source code looks like
uT Sf
φC = C α · n̂
|Sf |
where Cα equals the user variable cAlpha given in system/fvSolution. According to this expres-
sion, the correction flux φC equals zero if cAlpha is set to zero by the user. In this situation φCorr
α,l
(Equation 3.6) will look like
φCorr
α,l = φHO UD
α,l − φα,l (3.13)
where φHO
α,l is the user chosen discretization scheme given for αl
div(phi,alpha)
φHO
α,l = φα,l
φdiv(phi,alpha)
α,l = φU D
α,l
φα,l = φU D
α,l (3.14)
To summarize, a simple upwind solution of the uncompressed αl -equation can be obtained by setting
MULESCorr to true and cAlpha to 0, while the upwind discretization scheme is chosen for αl . For
every other discretization scheme, the α-fluxes will be equal to
φα,l = φU D Corr
α,l + θ · φα,l
17
3.1. CALCULATION OF αL - ALPHAEQN.H CHAPTER 3. DESIGN OF INTERFOAM
α0l
∂αl αln − αl0 + Su − psiIf
+ ∇ · (αl u) = Sp αln + Su =⇒ + psiIf = Sp αln + Su =⇒ αln = 1
δt
δt − Sp
∂t δt
(3.15)
This calculation are found in CMULESTemplates.C. In Equation (3.15) Su = −Sp αl , which means
that in fact no source term is added to the equation. Probably this is a trick done to ease the
numerical solution of the equation. The terms Su and Sp in this equation are to be set in the
alphaSuSp.H-file, and is set to zeroFields by default.
The argument values which are changed in MULES::correct(...) are φ = talphaPhi1Corr and
α = αlU D . The solution is limited according to the FCT method. talphaPhi1Corr is limited in
(3.5). phiIf in the routine is intergrated from phiCorr (CMULESTemplates.C)
52 fvc::surfaceIntegrate(psiIf, phiCorr);
where phiCorr is the limited version of the input variable talphaPhi1Corr (ref. Equation (3.5)).
this means that the FCT procedure given in Equation (3.2) is not followed. Nevertheless, the
compression term φCα,l is still included, and the solution is limited. In the source code this is written
as (in alphaEqn.H)
204 MULES::explicitSolve
205 (
206 geometricOneField(),
207 alpha1,
208 phiCN,
209 alphaPhi10,
210 Sp,
211 (Su + divU*min(alpha1(), scalar(1)))(),
212 1,
213 0
214 );
talphaPhi1Un = φHO
l + φC
l
18
3.1. CALCULATION OF αL - ALPHAEQN.H CHAPTER 3. DESIGN OF INTERFOAM
If cAlpha is set to 0, φC
α,l will be 0, and φα,l equals
div(phi,alpha)
φα,l = θ · φHO
α,l = θ · φα,l (3.17)
which means that the fluxes applied for the αl -field, are the fluxes obtained from a higher order
(user chosen div(phi,alpha) scheme) solution of the uncompressed αl equation (3.4).
When MULESCorr is set to false the routine MULES::explicitSolve is used for solving the alpha
field. This routine is implemented in MULESTemplates.C and solves the αl equation according to
α0l
∂αl αln − αl0 + Su − psiIf
+ ∇ · (αl u) = Sp αln + Su =⇒ + psiIf = Sp αln + Su =⇒ αln = δt
1
δt − Sp
∂t δt
(3.18)
where Su = Su,input + divU · αl . Here Su,input is the Su -term which is supposed to be given in
alphaSuSp.H and is set to a zero-field by default. divU is also to be set in alphaSuSp.H and is a
zero field by default. psiIf is the net flux, which is integrated from phiPsi (MULESTemplates.C)
55 fvc::surfaceIntegrate(psiIf, phiPsi);
where phiPsi is the limited version of the input variable talphaPhi1Un (ref. Equation (3.5)). The
solution of αl is limited between 0 and 1.
19
Chapter 4
Implementation of interAirEntFoam
The implementation of the air entrainment model starts by copying the interFoam solver, change
its name to interAirEntFoam, and rename its files and functions accordingly. As described in
Chapter 3, interFoam uses a class named immiscibleIncompressibleTwoPhaseMixture to han-
dle the transport properties. This class takes care of the reading and calculation of the transport
properties and phase properties, like phase fraction, density and viscosity, as well as the interface
properties. These properties are also desired to be incorporated within the corresponding trans-
port model to be used for the air entraintrainment model. This means that we start by copying
immiscibleIncompressibleTwoPhaseMixture and change its name to airEntrainmentTwoPhaseMix-
ture. In addition to the features of immiscibleIncompressibleTwoPhaseMixture, the new class
should calculate a source term of air, and decide on whether air should be added to the water phase
or not. This calculation is based on Equation (2.10). According to Section 2.1.2, our air-entrainment
model requires the knowledge of:
– Cell volume needed to divide the source term by its cell volume (ref. Equation (2.10)).
– Information about cell size also needed to calculate the variable As in Equation (2.10).
This section will provide the steps in the implementation of the air entrainment model named,
interAirEntFoam. The starting point is to make a local copy of interFoam.
4.1 Preparations
4.1.1 Make a local copy of interFoam
Start by copying the interFoam solver from the OpenFOAM installation as
Listing 4.1: Copy the interFoam solver
20
4.1. PREPARATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
cd $WM_PROJECT_USER_DIR
mkdir -p applications/solvers/multiphase
cd !$
cp -r $FOAM_SOLVERS/multiphase/interFoam .
mv interFoam interAirEntFoam
cd !$
mv interFoam.C interAirEntFoam.C
rm -r interMixingFoam overInterDyMFoam
Then you have to change the name of the excutable path to the user bin. This is done by running
the following sed-command on the files as
Listing 4.5: Change the name of the executable path
21
4.1. PREPARATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
Then you have to update the Make/options file, to include the desired libraries. In the Make/options
file, copied from the interFoam solver, it is linked to the VoF-directory by a relative path. Change
this to the absolute path. Open the Make/options-file by your favorite editor as
emacs Make/options
Listing 4.6: The first two lines in Make/options for the new solver.
1 EXE_INC = \
2 -I$(FOAM_SOLVERS)/multiphase/VoF \
by this you link to the volume of fluid solver within the foam solvers path. At this point you have
a local version of the interFoam solver. Try to compile it to see that everything works as intended
wmake
and copy the output from this search into the terminal window as
Listing 4.8: Copy the transport model immiscibleIncompressibleTwoPhaseMixture
cp -r $FOAM_SRC/transportModels/immiscibleIncompressibleTwoPhaseMixture .
22
4.1. PREPARATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
mv immiscibleIncompressibleTwoPhaseMixture airEntrainmentTwoPhaseMixture
cd airEntrainmentTwoPhaseMixture
mv immiscibleIncompressibleTwoPhaseMixture.C airEntrainmentTwoPhaseMixture.C
mv immiscibleIncompressibleTwoPhaseMixture.H airEntrainmentTwoPhaseMixture.H
sed -i s/immiscibleIncompressible/airEntrainment/g Make/files
sed -i s/immiscibleIncompressible/airEntrainment/g airEntrainmentTwoPhaseMixture.*
cd ..
sed -i s/immiscibleIncompressible/airEntrainment/g createFields.H
sed -i s/immiscibleIncompressible/airEntrainment/g interAirEntFoam.C
emacs airEntrainmentTwoPhaseMixture/Make/options
23
4.1. PREPARATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
1 EXE_INC = \
2 -I$(LIB_SRC)/transportModels \
3 -I$(LIB_SRC)/transportModels/incompressible/lnInclude \
4 -I$(LIB_SRC)/transportModels/interfaceProperties/lnInclude \
5 -I$(LIB_SRC)/transportModels/twoPhaseMixture/lnInclude \
6 -I$(LIB_SRC)/finiteVolume/lnInclude
Then your new library, airEntrainmentTwoPhaseMixture, which at this point is a local version of
immiscibleIncompressibleTwoPhaseMixture, should be ready to use. Compile it to ensure that
the name changes and linking has been correctly done
wmake airEntrainmentTwoPhaseMixture
When this compiles with no error messages, you are ready to proceed. Now you want your new
solver interAirEntFoam to use the new library. Then you have to link to it in the Make/options
file within the solver. Open the file in your favorite editor as
emacs Make/options
and do this by adding the following lines at the end of your Make/options file
Listing 4.13: Link the solver to the new user library library airEntrainmentTwoPhaseMixture
25 -L$(FOAM_USER_LIBBIN) \
26 -lairEntrainmentTwoPhaseMixture
also remember to add a backslash to the line that used to be the last line, since it is no longer the
last line in the file. Since you no longer make use of immiscibleIncompressibleTwoPhaseMixture,
the linking to this library can be removed in the Make/options file as
Listing 4.14: Remove linking to immiscibleIncompressibleTwoPhaseMixture in Make/options
Then a link to the new directory has to be linked to under EXE INC. Open the Make/options file in
your favorite editor
24
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
emacs Make/options
Listing 4.15: Link the to the new directory, airEntrainmentTwoPhaseMixture, within Make/options
2 -IairEntrainmentTwoPhaseMixture \
Then the linking should be in place and you can compile both your library and your solver to
ensure that it works. Make sure that you are placed in the interAirEntFoam directory, and run the
compiling commands as
wmake airEntrainmentTwoPhaseMixture
wmake
When both the library and the solver compiles, you should copy the damBreak tutorial, and run it
with your new solver as
cp -r $FOAM_TUTORIALS/multiphase/interFoam/RAS/damBreak/damBreak $FOAM_RUN
run
cd damBreak
sed -i s/interFoam/interAirEntFoam/g system/controlDict
./Allrun
4.2 Implementations
The aim of the solver is to calculate and return the output of Equation (2.10). This term is a
function of several terms
δV
= f (Pt (k, ρ), Pd (LT , k, g, ), ρ, Cair , As ) (4.1)
δt
The air entrainment model will be implemented within airEntrainmentTwoPhaseMixture
(airEntrainmentTwoPhaseMixture.H, airEntrainmentTwoPhaseMixture.C), which at this point
is a local copy of immiscibleIncompressibleTwoPhaseMixture. We start by calculating the per-
turbing forces Pt (k, ρ) according to Equation (2.11).
25
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
4.2.1 Pt (k, ρ)
This variable is calculated as a function of k, and ρmix . The densities are available in our
airEntrainmentTwoPhaseMixture model through the incompressibleTwoPhaseMixture class.
Nevertheless, we want to calculate a density for the mixture, we do this by introducing a pri-
vate member variable in airEntrainmentTwoPhaseMixture. Go to the interAirEntFoam folder
and open airEntrainmentTwoPhaseMixture.H in your favorite editor as
ufoam
cd applications/solvers/multiphase/interAirEntFoam
emacs airEntrainmentTwoPhaseMixture/airEntrainmentTwoPhaseMixture.H
57 // Private data
58 volScalarField rhomix_;
39 interfaceProperties(alpha1(), U, *this),
40
41 rhomix_
42 (
43 min(max(alpha1(), scalar(0)), scalar(1))*rho1_ +
44 (scalar(1) - min(max(alpha1(), scalar(0)), scalar(1)))*rho2_
45 )
and remember to add a comma after the call to interfaceProperties, since this is no longer
the last command in the constructor. The mixture density should be updated each time we
call the source term from the solver. That means that an update function is needed for this
density. We implement this as a public member function in our class. Open the header file
airEntrainmentTwoPhaseMixture.H, and implement the update function as
26
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
59 void Foam::airEntrainmentTwoPhaseMixture::updateRhomix()
60 {
61 rhomix_ = min(max(alpha1(), scalar(0)), scalar(1))*rho1_ +
62 (scalar(1) - min(max(alpha1(), scalar(0)), scalar(1)))*rho2_;
63 }
From the solver directory, compile the library and the solver to check if it works
wmake airEntrainmentTwoPhaseMixture
wmake
To check that our new function works we try to call it from one of the files within the solver, we do
this in createFields.H. Open the file in your favorite editor
emacs createFields.H
37 [Link]();
wmake
When this works, we proceed to the implementation of Pt (ρ, k). As Pt is a function of the total kinetic
27
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
energy k, it needs the knowledge of this variable. Our class airEntrainmentTwoPhaseMixture does
not know about the turbulence variables. A simple (but perhaps not very elegant) way to attain the
knowledge of k, is to pass it as an argument to a function. We put the following in the header file
(airEntrainmentTwoPhaseMixture.H)
87 // Calculates the pertubing component that makes the flow unstable and
88 // returnes a volumScalarField
89 volScalarField Pt(const volScalarField k);
which means that we implement a public member function that returns a volScalarField. We
place the implementation of the function in the .C-file as
65 Foam::volScalarField Foam::airEntrainmentTwoPhaseMixture::Pt(const
66 volScalarField k)
67 {
68 updateRhomix();
69 volScalarField Pt = k*rhomix_;
70 return Pt;
71 }
Again, make a function call to your function from createFields.H to check whether the function
works. Since this function use k as an argument, the function call should be made after the initiation
of the turbulence model. Open the file in your favorite editor
emacs createFields.H
82 [Link](turbulence->k());
wmake airEntrainmentTwoPhaseMixture
wmake
28
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
4.2.2 Pd (LT , k, g, )
Since Pd is a function of the characteristic length, LT (Cµ , k, ), it is convenient to start by im-
plementing the function to calculate LT . We implement this as a member function in our class
(airEntrainmentTwoPhaseMixture), returning a volScalarField, similar to what was done for
Pt . Paste the following in your header file (airEntrainmentTwoPhaseMixture.H)
73 Foam::volScalarField Foam::airEntrainmentTwoPhaseMixture::Lt(
74 const volScalarField k,
75 const volScalarField epsilon)
76 {
77 dimensionedScalar dimCorr("dimCorr",dimLength/dimTime,1);
78 volScalarField newepsilon = epsilon +
79 dimensionedScalar("smallEpsilon",dimensionSet(0,2,-3,0,0,0,0),1e-10);
80 // dimensions did not follow the meth operations -> dimCorr
81 volScalarField Lt = 0.09*sqr(1.5)*pow(k,3/2)/newepsilon*dimCorr;
82 return Lt;
83 }
83 [Link](turbulence->k(), turbulence->epsilon());
Again, it should be noted that this function call has to be executed after the initialization of the
29
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
wmake airEntrainmentTwoPhaseMixture
wmake
interAirEntFoam -case $FOAM_RUN/damBreak
Then we can continue implementing the surface stabilizing term Pd (LT , σ, g). In addition to LT , this
variable is a function of surface tension σ and gravity g. These variables are not directly available in
out airEntrainmentTwoPhaseMixture model, but they can be read from the input files. These vari-
ables will be constant for a simulation, as long as they are not altered by the user. This means that
it make sense to add them as private member variables in our airEntrainmentTwoPhaseMixture
class. Add the following to airEntrainmentTwoPhaseMixture.H
Here we introduce a new variable type, which means that the following must be added to the header
file
40 #include "uniformDimensionedFields.H"
47 sigma_
48 (
49 "sigma", dimensionSet(1,0,-2,0,0,0,0), lookup("sigma")
50 ),
51
52 g_(
53 IOobject
54 (
55 "g",
56 [Link]().constant(),
57 [Link](),
58 IOobject::MUST_READ_IF_MODIFIED,
59 IOobject::NO_WRITE
60 )
61 )
30
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
and remember to add a comma before this implementation. As private member variables σ and g
will be available for the member functions of the solver. Compile the model to check whether it
works
wmake airEntrainmentTwoPhaseMixture
when this works, we should be ready to start on the implementation of Pd . We declare the function
accordingly in the header file as
At line 107 this function make use of the namespace fvc, which is not declared within this scope.
Include the fvCFD.H in airEntrainmentTwoPhaseMixture.H to fix this
41 #include "fvCFD.H"
31
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
6 -I$(LIB_SRC)/meshTools/lnInclude \
14 -lmeshTools \
84 [Link](turbulence->k(), turbulence->epsilon());
Compile the library, the solver, and run the damBreak case again
wmake airEntrainmentTwoPhaseMixture
wmake
interAirEntFoam -case $FOAM_RUN/damBreak
35 Cair 0.5;
32
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
51 Cair_
52 (
53 "Cair", dimensionSet(0,0,0,0,0,0,0), lookup("Cair")
54 ),
Compile the library, the solver, and run the damBreak case again
wmake airEntrainmentTwoPhaseMixture
wmake
interAirEntFoam -case $FOAM_RUN/damBreak
As should represent the free surface area in an interface cell. This could be calculated, or at least
estimated based on the gradients of the αl field. Nevertheless, in this model, this calculation was
simplified to be estimated based on the cell volume as
The solver interIsoFoam uses isosurfaces to capture the interface. Perhaps a method for a proper
calculation of As can be found from this solver.
4.2.4 CalcSource
Then we are ready to implement our final function that will return an estimate of the amount of
air that should be added to a cell. The source term should only be active at the interface. This
means that the model needs to know where the interface is located. A function nearInterface(),
which returns a volScalarField, is available in our model (inherited from interfaceProperties).
We try to make use of this function to notify the placement of the interface. Our end goal for
the airEntrainmentTwoPhaseMixture model is that it will provide a source term that can be
added as an explicit source term to the αl -equation. We declare the function calcSource() in
the airEntrainmentTwoPhaseMixture.H-file as
33
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
Since we are going to add the source term to the αl -equation, we have to add the volume as a
volume fraction relative to the relevant cell. Then this function needs the knowledge of the cell
volume, which is a reason for providing U as a function argument here. Nevertheless, since U is
already an input to the constructor, and the cell volume do not change for this solver, a better
solution might have been to add U as a member variable to the class. However, we proceed to
implement the function in the .C-file as
34
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
The formulation for the source term given in 2.10 gives the volume of air to be added to a cell. In
our solver, it is the water fraction that is calculated in the αl -equation. This means that the amount
calculated has to be subtracted from the αl -equation. This is taken care of in the calculation of the
source term above. Also the possible amount of air should be restricted according to the relevant cell
volume. This is done by line 148 in the listing above. The volume of the source term is restricted
to half of the cell volume, which is more of an arbitrary choice.
Test your model by introducing a function call to calcSource() in createFields.H, as you should
know how to do by now. Compile the model, solver, and run the damBreak case accordingly
wmake airEntrainmentTwoPhaseMixture
wmake
interAirEntFoam -case $FOAM_RUN/damBreak
when this works, the next step is to integrate our source term in the calculation of the α-field.
99 if (MULESCorr)
100 {
101 #include "alphaSuSp.H"
102
35
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
122 [Link]();
At line 101 above there is a call to the file alphaSuSp.H. This is the file where we are supposed to
add our source term. Nevertheless, we also want to visualize our source term, therefore we create
the field in createFields.H as
Then we update our source term in alphaSuSp.H by a function call to calcSource, by adding the
following
1 Su = [Link](U,turbulence->k(),
2 turbulence->epsilon());
3 Info << "Max(Su): " << max(Su).value()
4 << "\t Min(Su): " << min(Su).value() << endl;
Compile the solver, and run the damBreak case to test the implementation of the source term.
36
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
wmake
interAirEntFoam -case $FOAM_RUN/damBreak
Open the case in paraFoam and view the Su -field. If this field is not available in paraFoam something
is wrong.
In alphaSuSp.H the inplicit source term Sp and divU are created. In the interFoam solver these
fields are given as zeroFields, which is the correct way to do it when you don’t have any source
terms. In our solver we have already made the explicit source term Su unzero. We keep divU and
Sp as zero fields (in alphaSuSp.H)
6 zeroField Sp;
7 zeroField divU;
Then we have implemented our source term into the αl -equation. As described in Section 3.1, the
user can choose whether to apply MULES correction or not. By setting MULES to false/no in the
system/fvSolution file
solvers
{
"[Link].*"
{
nAlphaCorr 2;
nAlphaSubCycles 1;
cAlpha 1;
MULESCorr no;
The default behavior in the solver is false, which means that you actively have to turn MULES on
to apply it. If MULES is not activated, our source term Su is added to the solution algorithm in
the function argument to MULES::explicitSolve as
204 MULES::explicitSolve
205 (
206 geometricOneField(),
207 alpha1,
208 phiCN,
209 alphaPhi10,
210 Sp,
211 (Su + divU*min(alpha1(), scalar(1)))(),
212 1,
37
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
213 0
214 );
When MULES is activated, the source term Su is added to the solution of the upwind flux, as shown
some paragraphs above. The source term is then not added to the calculation of the correction flux
or the higher order flux term. Nevertheless, the calculation of the higher order fluxes are based on
the upwind solution of the αl -field. In MULES mode the αl -equation is solved through a function
call to MULES::Correct as
177 MULES::correct
178 (
179 geometricOneField(),
180 alpha1,
181 talphaPhi1Un(),
182 [Link](),
183 Sp,
184 (-Sp*alpha1)(),
185 1,
186 0
187 );
here, the source term Su is not explicitly included in this function call. Nevertheless, it is implicitly
included in the term [Link](), which equals φcorr α,l in Equation (3.5).
38
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
181 IOobject
182 (
183 "PtmPd",
184 [Link](),
185 mesh,
186 IOobject::READ_IF_PRESENT, //NO_READ,
187 IOobject::AUTO_WRITE //NO_WRITE
188 ),
189 mesh,
190 dimensionedScalar("PtmPd",dimMass/dimTime/dimTime/dimLength, 0)
191 );
We want to update these fields during the simulation. This seems convenient to do in alphaSuSp.H
after the function call to calcSource as
17 forAll(test,cellI)
18 {
19 if(test[cellI] < 0)
20 {
21 PtmPd[cellI] = pos(test[cellI]);
22 }else{
23 PtmPd[cellI] = test[cellI];
24 }
25 }
27 forAll([Link](),patchI)
28 {
29 forAll([Link]()[patchI], faceI)
30 {
31 [Link]()[patchI][faceI] = scalar(0.0);
32 }
33 }
34 [Link]();
Towards the end here, the boundary fields are set to 0, and corrected. This is perhaps not a decent
way to do it, but it seemed to be a problem if they were left without setting some values. Then
39
4.2. IMPLEMENTATIONS CHAPTER 4. IMPLEMENTATION OF INTERAIRENTFOAM
our intended implementations are done. Nevertheless, during the implementation work, we added
some functions in createFields.H to test our implementations. Remove these lines using the sed-
command, and compile and test the solver by typing the following in the command window
Ensure that the new fields are available in paraFoam. When this works, the implementation of the
solver is complete, and we can start to review its performance by a verification case. This is done
in the following chapter.
40
Chapter 5
Verification cases
cd $FOAM_RUN/jetintopool
./Allrun &
This will run the case for 10 seconds on a single core. Open the case in paraFoam, and ensure that
the new fields are available. Then copy the case, to a new folder, and run it with interFoam to
compare the behavior
cd ..
cp -r jetintopool jetintopool_interFoam
sed -i s/interAirEntFoam/interFoam/g system/controlDict
./Allrun &
Open the results side by side in paraFoam, using the buttoms in the upper rigth corner. Apply the
field [Link] and play the simulation. At the end time step, the results will be as shown
in Figure 5.1. From Figure (5.1) you can see that a larger region below the surface has a higher air
concentration compared to the results obtained by interFoam. This means that our air entrainment
model at least add some air to some regions of the water phase.
Let us have a closer look into the actual fields calculated by the air entrainment model. Go into the
interAirEntFoam version of the jetintopool case, and open a new paraFoam window. Color the
case by [Link], and apply a contour line showing [Link] = 0.5. Show the contour line
in Solid Color and Wirefram. This should give you the view of Figure 5.3 at the last time step.
41
5.1. VERTICAL PLUNGING JET CHAPTER 5. VERIFICATION CASES
Then you apply the positive source term Sg field, together with the countour line. Go back to the
first time step of the simulation and go through the simulation step by step. This will give you the
views shown in Figure 5.3 for some of the first time steps. This figure indicates that the solver at
least conceptually seems to work as intended for a vertical plunging jet on a coarse grid.
42
5.2. POSSIBLE IMPROVEMENTS CHAPTER 5. VERIFICATION CASES
43
Study questions
1. Which parameter(s) do you have to add or change in the input file to run a tutorial interFoam
case with interAirEntFoam?
2. What is the purpose of the implementations done to develop interAirEntFoam
3. Which method is used to compress the interface in interFoam
4. How can you manipulate the input variables to obtain the uncorrected higher order solution
for the αl -equation?
5. What is the criteria for getting air entrainment on a free water surface?
6. Why will it be convenient to model the air entrainment with a subgrid model?
7. How do you link to a user library in the Make/options file of a solver?
44
References
Brackbill, J., Kothe, D., & Zemach, C. (1992). A continuum method for modeling
surface tension. Journal of Computational Physics, 100 (2), 335 - 354. Retrieved
from [Link] doi:
[Link]
Hemida, H. (2008). Free surface tutorial using interfoam and rasinterfoam. In Proceedings of CFD
with OpenSource Software, 2008, Edited by Nilsson. H.. doi: [Link] -
CFD#YEAR 2008
Hirt, C. (2003). Modeling turbulent entrainment of air at a free surface. Flow Science, Inc.
Manni, A. (2014). An introduction to twophaseeulerfoam with addition of an heat exchange
model. In Proceedings of CFD with OpenSource Software, 2014, Edited by Nilsson. H.. doi:
[Link] CFD#YEAR 2014
Olsson, E. (2017). A description of isoadvector - a numerical method for improved surface sharpness
in two-phase flows. In Proceedings of CFD with OpenSource Software, 2017, Edited by Nilsson.
H.. doi: [Link] CFD#YEAR 2017
Oruganti, S. K. (2015). Coupled level-set with vof interfoam. In Proceedings of CFD with OpenSource
Software, 2015, Edited by Nilsson. H.. doi: [Link] CFD#YEAR -
2015
Wikipedia. (2018). Volume of fluid-metoden. [Link] of -
fluid-metoden. (Accessed: 2018.12.28)
45