2.3 Breaking of a dam

In this tutorial we shall solve a problem of a simplified dam break in 2 dimensions using the interFoam solver. The feature of the problem is a transient flow of two fluids separated by a sharp interface, or free surface. The two-phase algorithm in interFoam is based on the volume of fluid (VOF) method in which a phase transport equation is used to determine the relative volume fraction of the two phases, or phase fraction eqn, in each computational cell. Physical properties are calculated as weighted averages based on this fraction. The nature of the VOF method means that an interface between the phases is not explicitly computed, but rather emerges as a property of the phase fraction field. Since the phase fraction can have any value between 0 and 1, the interface is never sharply defined, but occupies a volume around the region where a sharp interface should exist.

The test setup consists of a column of water at rest located behind a membrane on the left side of a tank. At time eqn, the membrane is removed and the column of water collapses. During the collapse, the water impacts an obstacle at the bottom of the tank and creates a complicated flow structure, including several captured pockets of air. The geometry and the initial setup is shown in Figure 2.21 .


 0.584 m water column 0.584 m 0.292m 0.048 m 0.1461 m 0.1459 m 0.024 m \relax \special {t4ht=


Figure 2.21: Geometry of the dam break.


2.3.1 Mesh generation

The user should go to their run directory and copy the damBreak case from the $FOAM_TUTORIALS/multiphase/interFoam/laminar/damBreak directory, i.e.


    run
    cp -r $FOAM_TUTORIALS/multiphase/interFoam/laminar/damBreak/damBreak .
Go into the damBreak case directory and generate the mesh running blockMesh as described previously. The damBreak mesh consist of 5 blocks; the blockMeshDict entries are given below.
16convertToMeters 0.146;
17
18vertices
19(
20    (0 0 0)
21    (2 0 0)
22    (2.16438 0 0)
23    (4 0 0)
24    (0 0.32876 0)
25    (2 0.32876 0)
26    (2.16438 0.32876 0)
27    (4 0.32876 0)
28    (0 4 0)
29    (2 4 0)
30    (2.16438 4 0)
31    (4 4 0)
32    (0 0 0.1)
33    (2 0 0.1)
34    (2.16438 0 0.1)
35    (4 0 0.1)
36    (0 0.32876 0.1)
37    (2 0.32876 0.1)
38    (2.16438 0.32876 0.1)
39    (4 0.32876 0.1)
40    (0 4 0.1)
41    (2 4 0.1)
42    (2.16438 4 0.1)
43    (4 4 0.1)
44);
45
46blocks
47(
48    hex (0 1 5 4 12 13 17 16) (23 8 1) simpleGrading (1 1 1)
49    hex (2 3 7 6 14 15 19 18) (19 8 1) simpleGrading (1 1 1)
50    hex (4 5 9 8 16 17 21 20) (23 42 1) simpleGrading (1 1 1)
51    hex (5 6 10 9 17 18 22 21) (4 42 1) simpleGrading (1 1 1)
52    hex (6 7 11 10 18 19 23 22) (19 42 1) simpleGrading (1 1 1)
53);
54
55defaultPatch
56{
57    type empty;
58}
59
60boundary
61(
62    leftWall
63    {
64        type wall;
65        faces
66        (
67            (0 12 16 4)
68            (4 16 20 8)
69        );
70    }
71    rightWall
72    {
73        type wall;
74        faces
75        (
76            (7 19 15 3)
77            (11 23 19 7)
78        );
79    }
80    lowerWall
81    {
82        type wall;
83        faces
84        (
85            (0 1 13 12)
86            (1 5 17 13)
87            (5 6 18 17)
88            (2 14 18 6)
89            (2 3 15 14)
90        );
91    }
92    atmosphere
93    {
94        type patch;
95        faces
96        (
97            (8 20 21 9)
98            (9 21 22 10)
99            (10 22 23 11)
100        );
101    }
102);
103
104
105// ************************************************************************* //

2.3.2 Boundary conditions

The user can examine the boundary geometry generated by blockMesh by viewing the boundary file in the constant/polyMesh directory. The file contains a list of 5 boundary patches: leftWall, rightWall, lowerWall, atmosphere and defaultFaces.

The user should notice the type of the patches. Firstly, the atmosphere is a standard patch, i.e. has no special attributes, merely an entity on which boundary conditions can be specified. Then, the defaultFaces patch is formed of block faces that are omitted from the boundary sub-dictionary in the blockMeshDict file. Those block faces form a patch whose properties are specified in a defaultPatch sub-dictionary in in the blockMeshDict file. In this case, the default type is set to empty since the patch normal is in the direction we will not solve in this 2D case.

The leftWall, rightWall and lowerWall patches are each a wall. Like the generic patch, the wall type contains no geometric or topological information about the mesh and only differs from the plain patch in that it identifies the patch as a wall, should an application need to know, e.g. to apply special wall surface modelling. For example, the interFoam solver includes modelling of surface tension and can include wall adhesion at the contact point between the interface and wall surface. Wall adhesion models can be applied through a special boundary condition on the alpha (eqn) field, e.g. the constantAlphaContactAngle boundary condition, which requires the user to specify a static contact angle, theta0.

In this tutorial we would like to ignore surface tension effects between the wall and interface. We can do this by setting the static contact angle, eqn. However, rather than using the constantAlphaContactAngle boundary condition, the simpler zeroGradient can be applied to alpha on the walls.

The top boundary is free to the atmosphere so needs to permit both outflow and inflow according to the internal flow. We therefore use a combination of boundary conditions for pressure and velocity that does this while maintaining stability. They are:

  • prghTotalPressure, applied to the pressure field, without the hydrostatic component, eqn, which is a fixedValue condition calculated from a specified total pressure p0 and local velocity U;
  • pressureInletOutletVelocity, applied to velocity eqn, which sets zeroGradient on all components of eqn, except where there is inflow, in which case a fixedValue condition is applied to the tangential component;
  • inletOutlet applied to other fields, which is a zeroGradient condition when flow outwards, fixedValue when flow is inwards.

At all wall boundaries, the fixedFluxPressure boundary condition is applied to the pressure field, which adjusts the pressure gradient so that the boundary flux matches the velocity boundary condition for solvers that include body forces such as gravity and surface tension.

The defaultFaces patch representing the front and back planes of the 2D problem, is, as usual, an empty type.

2.3.3 Phases

The fluid phases are specified in the phaseProperties file in the constant directory as follows:

16
17phases          (water air);
18
19sigma           0.07;
20
21
22// ************************************************************************* //

It lists two phases, water and air. Equations for phase fraction are solved for the phases in the list, except the last phase listed, i.e. air in this case. Since there are only two phases, only one phase fraction equation is solved in this case, for the water phase fraction eqn, specified in the file alpha.water in the 0 directory.

The phaseProperties file also contains an entry for the surface tension between the two phases, specified by the keyword sigma in units eqn.

2.3.4 Setting initial fields

Unlike the previous cases, we shall now specify a non-uniform initial condition for the phase fraction eqn where

 { α = 1 for the water phase water 0 for the air phase \relax \special {t4ht=
(2.15)
This will be done by running the setFields utility. It requires a setFieldsDict dictionary, located in the system directory, whose entries for this case are shown below.
16
17defaultFieldValues
18(
19    volScalarFieldValue alpha.water 0
20);
21
22regions
23(
24    boxToCell
25    {
26        box (0 0 -1) (0.1461 0.292 1);
27        fieldValues
28        (
29            volScalarFieldValue alpha.water 1
30        );
31    }
32);
33
34
35// ************************************************************************* //

The defaultFieldValues sets the default value of the fields, i.e. the value the field takes unless specified otherwise in the regions sub-dictionary. That sub-dictionary contains a list of subdictionaries containing fieldValues that override the defaults in a specified region. The region creates a set of points, cells or faces based on some topological constraint. Here, boxToCell creates a bounding box within a vector minimum and maximum to define the set of cells of the water region. The phase fraction eqn is defined as 1 in this region.

The setFields utility reads fields from file and, after re-calculating those fields, will write them back to file. In the damBreak tutorial, the alpha.water field is initially stored as a backup named alpha.water.orig. A field file with the .orig extension is read in when the actual file does not exist, so setFields will read alpha.water.orig but write the resulting output to alpha.water (or alpha.water.gz if compression is switched on). This way the original file is not overwritten, so can be reused.

The user should therefore execute setFields like any other utility by:


    setFields
Using paraFoam, check that the initial alpha.water field corresponds to the desired distribution as in Figure 2.22.

Phase fraction, α1 1.0 0.9 0.8 0.7 0.6 0.5 0.4 0.3 0.2 0.1 0.0 \relax \special {t4ht=


Figure 2.22: Initial conditions for phase fraction alpha.water.


2.3.5 Fluid properties

The physical properties for the air and water phases are specified in physicalProperties.air and physicalProperties.water files, respectively, in the constant directory. In each file, the viscosity model is selected by the viscosityModel keyword. The constant model is selected for both fluids, i.e. representing them as Newtonian fluids with a kinematic viscosity specified by the nu keyword in units eqn. Other non-Newtonian viscosity models are described in section 7.3 . The density of each fluid is specified by the keyword rho in units eqn. The physicalProperties.air file is shown below as an example:

16
17viscosityModel  constant;
18
19nu              1.48e-05;
20
21rho             1;
22
23
24// ************************************************************************* //

Gravitational acceleration is uniform across the domain and is specified in a file named g in the constant directory. Unlike a normal field file, e.g. U and p, g is a uniformDimensionedVectorField and so simply contains a set of dimensions and a value that represents eqn for this tutorial:

16
17dimensions      [0 1 -2 0 0 0 0];
18value           (0 -9.81 0);
19
20
21// ************************************************************************* //

2.3.6 Turbulence modelling

As in the cavity example, the choice of turbulence modelling method is selectable at run-time through the simulationType keyword in momentumTransport dictionary. In this example, we wish to run without turbulence modelling so we set laminar:

16
17simulationType  laminar;
18
19
20// ************************************************************************* //

2.3.7 Time step control

Time step control is an important issue in transient simulation and the surface-tracking algorithm in interface-capturing solvers. The interFoam solver uses the multidimensional universal limiter for explicit solution (MULES) method, created by Henry Weller, to maintain boundedness of the phase fraction. The Courant number eqn needs to be limited depending on the choice of algorithm: with the explicit MULES algorithm, an upper limit of eqn for stability is typical in the region of the interface; but with semi-implicit MULES, specified by the MULESCorr keyword in the fvSolution file, there is really no upper limit in eqn for stability, but instead the level is determined by requirements of temporal accuracy.

In general it is difficult to specify a fixed time-step to satisfy the eqn criterion, so interFoam offers automatic adjustment of the time step as standard in the controlDict. The user should specify adjustTimeStep to be on and the the maximum eqn for the phase fields, maxAlphaCo, and other fields, maxCo, to be 1.0. The upper limit on time step maxDeltaT can be set to a value that will not be exceeded in this simulation, e.g. 1.0.

By using automatic time step control, the steps themselves are never rounded to a convenient value. Consequently if we request that OpenFOAM saves results at a fixed number of time step intervals, the times at which results are saved are somewhat arbitrary. However even with automatic time step adjustment, OpenFOAM allows the user to specify that results are written at fixed times; in this case OpenFOAM forces the automatic time stepping procedure to adjust time steps so that it ‘hits’ on the exact times specified for write output. The user selects this with the adjustableRunTime option for writeControl in the controlDict dictionary. The controlDict dictionary entries should be:

16
17application     interFoam;
18
19startFrom       startTime;
20
21startTime       0;
22
23stopAt          endTime;
24
25endTime         1;
26
27deltaT          0.001;
28
29writeControl    adjustableRunTime;
30
31writeInterval   0.05;
32
33purgeWrite      0;
34
35writeFormat     binary;
36
37writePrecision  6;
38
39writeCompression off;
40
41timeFormat      general;
42
43timePrecision   6;
44
45runTimeModifiable yes;
46
47adjustTimeStep  yes;
48
49maxCo           1;
50maxAlphaCo      1;
51
52maxDeltaT       1;
53
54
55// ************************************************************************* //

2.3.8 Discretisation schemes

The MULES method, used by the interFoam solver, maintains boundedness of the phase fraction independently of the underlying numerical scheme, mesh structure, etc. The choice of schemes for convection are therfore not restricted to those that are strongly stable or bounded, e.g. upwind differencing.

The convection schemes settings are made in the divSchemes sub-dictionary of the fvSchemes dictionary. In this example, the convection term in the momentum equation (eqn), denoted by the div(rhoPhi,U) keyword, uses Gauss linearUpwind grad(U) to produce good accuracy. Here, we have opted for best stability with eqn. The eqn term, represented by the div(phi,alpha) keyword uses a bespoke interfaceCompression scheme where the specified coefficient is a factor that controls the compression of the interface where: 0 corresponds to no compression; 1 corresponds to conservative compression; and, anything larger than 1, relates to enhanced compression of the interface. We generally adopt a value of 1.0 which is employed in this example.

The other discretised terms use commonly employed schemes so that the fvSchemes dictionary entries should therefore be:

16
17ddtSchemes
18{
19    default         Euler;
20}
21
22gradSchemes
23{
24    default         Gauss linear;
25}
26
27divSchemes
28{
29    div(rhoPhi,U)  Gauss linearUpwind grad(U);
30    div(phi,alpha)  Gauss interfaceCompression vanLeer 1;
31    div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
32}
33
34laplacianSchemes
35{
36    default         Gauss linear corrected;
37}
38
39interpolationSchemes
40{
41    default         linear;
42}
43
44snGradSchemes
45{
46    default         corrected;
47}
48
49
50// ************************************************************************* //

2.3.9 Linear-solver control

In the fvSolution file, the alpha.water sub-dictionary in solvers contains elements that are specific to interFoam. Of particular interest is the nAlphaCorr keyword which controls the number of iterations of the phase fraction equation within a solution step. The iteration is used to overcome nonlinearities in the advection which are present in this case due to the interfaceCompression scheme.

2.3.10 Running the code

Running of the code has been described in detail in previous tutorials. Try the following, that uses tee, a command that enables output to be written to both standard output and files:


    cd $FOAM_RUN/damBreak
    interFoam | tee log
The code will now be run interactively, with a copy of output stored in the log file.

2.3.11 Post-processing

Post-processing of the results can now be done in the usual way. The user can monitor the development of the phase fraction alpha.water in time, e.g. see Figure 2.23 .


Phase fraction, α1 1.0 0.9 0.8 0.7 0.6 0.5 0.4 0.3 0.2 0.1 0.0 \relax \special {t4ht=
(a) At eqn.
Phase fraction, α1 1.0 0.9 0.8 0.7 0.6 0.5 0.4 0.3 0.2 0.1 0.0 \relax \special {t4ht=
(b) At eqn.

Figure 2.23: Snapshots of phase eqn.


2.3.12 Running in parallel

The results from the previous example are generated using a fairly coarse mesh. We now wish to increase the mesh resolution and re-run the case. The new case will typically take a few hours to run with a single processor so, should the user have access to multiple processors, we can demonstrate the parallel processing capability of OpenFOAM.

The user should first clone the damBreak case, e.g. by


    run
    foamCloneCase damBreak damBreakFine
Enter the new case directory and change the blocks description in the blockMeshDict dictionary to


    blocks
    (
        hex (0 1 5 4 12 13 17 16) (46 10 1) simpleGrading (1 1 1)
        hex (2 3 7 6 14 15 19 18) (40 10 1) simpleGrading (1 1 1)
        hex (4 5 9 8 16 17 21 20) (46 76 1) simpleGrading (1 2 1)
        hex (5 6 10 9 17 18 22 21) (4 76 1) simpleGrading (1 2 1)
        hex (6 7 11 10 18 19 23 22) (40 76 1) simpleGrading (1 2 1)
    );
Here, the entry is presented as printed from the blockMeshDict file; in short the user must change the mesh densities, e.g. the 46 10 1 entry, and some of the mesh grading entries to 1 2 1. Once the dictionary is correct, generate the mesh by running blockMesh.

As the mesh has now changed from the damBreak example, the user must re-initialise the phase field alpha.water in the 0 time directory since it contains a number of elements that is inconsistent with the new mesh. Note that there is no need to change the U and p_rgh fields since they are specified as uniform which is independent of the number of elements in the field. We wish to initialise the field with a sharp interface, i.e. it elements would have eqn or eqn. Updating the field with mapFields may produce interpolated values eqn at the interface, so it is better to rerun the setFields utility.

The mesh size is now inconsistent with the number of elements in the alpha.water file in the 0 directory, so the user must delete that file so that the original alpha.water.orig file is used instead.


    rm 0/alpha.water
    setFields

The method of parallel computing used by OpenFOAM is known as domain decomposition, in which the geometry and associated fields are broken into pieces and allocated to separate processors for solution. The first step required to run a parallel case is therefore to decompose the domain using the decomposePar utility. There is a dictionary associated with decomposePar named decomposeParDict which is located in the system directory of the tutorial case. Also, sample dictionaries can be found within the etc directory in the OpenFOAM installation, which can be copied to the case directory by running the foamGet script, e.g. (you do not need to do this):


    foamGet decomposeParDict

The first entry is numberOfSubdomains which specifies the number of subdomains into which the case will be decomposed, usually corresponding to the number of processors available for the case.

In this tutorial, the method of decomposition should be simple and the corresponding simpleCoeffs should be edited according to the following criteria. The domain is split into pieces, or subdomains, in the eqn, eqn and eqn directions, the number of subdomains in each direction being given by the vector eqn. As this geometry is 2 dimensional, the 3rd direction, eqn, cannot be split, hence eqn must equal 1. The eqn and eqn components of eqn split the domain in the eqn and eqn directions and must be specified so that the number of subdomains specified by eqn and eqn equals the specified numberOfSubdomains, i.e. eqn numberOfSubdomains. It is beneficial to keep the number of cell faces adjoining the subdomains to a minimum so, for a square geometry, it is best to keep the split between the eqn and eqn directions should be fairly even.

For example, let us assume we wish to run on 4 processors. We would set numberOfSubdomains to 4 and eqn. The user should run decomposePar with:


    decomposePar
The terminal output shows that the decomposition is distributed fairly even between the processors.

The user should consult section 3.4 for details of how to run a case in parallel; in this tutorial we merely present an example of running in parallel. We use the openMPI implementation of the standard message-passing interface (MPI). As a test here, the user can run in parallel on a single node, the local host only, by typing:


    mpirun -np 4 interFoam -parallel > log &

The user may run on more nodes over a network by creating a file that lists the host names of the machines on which the case is to be run as described in section 3.4.3 . The case should run in the background and the user can follow its progress by monitoring the log file as usual.


 \relax \special {t4ht=


Figure 2.24: Mesh of processor 2 in parallel processed case.



Phase fraction, α1 1.0 0.9 0.8 0.7 0.6 0.5 0.4 0.3 0.2 0.1 0.0 \relax \special {t4ht=
(a) At eqn.
Phase fraction, α1 1.0 0.9 0.8 0.7 0.6 0.5 0.4 0.3 0.2 0.1 0.0 \relax \special {t4ht=
(b) At eqn.

Figure 2.25: Snapshots of phase eqn with refined mesh.


2.3.13 Post-processing a case run in parallel

Once the case has completed running, the decomposed fields and mesh can be reassembled for post-processing using the reconstructPar utility. Simply execute it from the command line. The results from the fine mesh are shown in Figure 2.25 . The user can see that the resolution of interface has improved significantly compared to the coarse mesh.

The user may also post-process an individual region of the decomposed domain individually by simply treating the individual processor directory as a case in its own right. For example if the user starts paraFoam by


    paraFoam -case processor1
then processor1 will appear as a case module in ParaView. Figure 2.24 shows the mesh from processor 1 following the decomposition of the domain using the simple method.
OpenFOAM v10 User Guide - 2.3 Breaking of a dam