Water drop

This case is a simple water drop showing the ability of freshkiss3d to simulate wave propagation and reflections on solid boundaries.

import os, sys
import matplotlib.pylab as plt
import numpy as np
import freshkiss3d as fk
import freshkiss3d.extra.plots as fk_plt
import argparse

os.system('rm -r outputs')

#sphinx_gallery_thumbnail_number = 2

parser = argparse.ArgumentParser()
parser.add_argument('--nographics', action='store_true')
args = parser.parse_args()

Parameters:

NUM_PARAMS={'ipres':True,
            'space_second_order':False,
            'flux_type':2,
            'implicit_exchanges':True,
            'implicit_vertical_viscosity':True}

Time loop:

simutime = fk.SimuTime(final_time=4., time_iteration_max=20000, second_order=True)

create_figure_scheduler = fk.schedules(times=[0., 1., 2., 3., 4., 5.])

Mesh:

dir_path = os.path.abspath(os.path.dirname(sys.argv[0]))
triangular_mesh = fk.TriangularMesh.from_msh_file(dir_path + '/inputs/square2.mesh')

if not args.nographics:
    fk_plt.plot_mesh(triangular_mesh)
Mesh

Layers:

NL = 1
layer = fk.Layer(NL, triangular_mesh, topography=0.)

Primitives:

Initial height is set with a user defined function.

def H_0(x, y):
    h = 2.4*(1.0 + np.exp(-0.25*((x-10.05)**2+(y-10.05)**2)))
    return h

primitives = fk.Primitive(triangular_mesh, layer, height_funct=H_0)

Boundary conditions:

case = 1

if case == 1:
    slides = [fk.Slide(ref=r) for r in [1, 2, 3, 4]]
    damping_outs = []

if case == 2:
    slides = []
    damping_outs = [fk.DampingOut(ref=r, Kappa=1000., Theta=1000.) \
    for r in [1, 2, 3, 4]]

if case == 3:
    slides = [fk.Slide(ref=r) for r in [1, 3, 4]]
    damping_outs = [fk.DampingOut(ref=2, Kappa=1000., Theta=1000.)]

Writter:

vtk_writer = fk.VTKWriter(triangular_mesh,
                          scheduler=fk.schedules(count=10), scale_h=1.)

Problem definition:

problem = fk.Problem(simutime, triangular_mesh, layer, primitives,
                     slides=slides,
                     damping_outs=damping_outs,
                     numerical_parameters=NUM_PARAMS,
                     vtk_writer=vtk_writer,
                     custom_funct={'plot':fk_plt.plot_freesurface_3d},
                     custom_funct_scheduler=create_figure_scheduler)
===================================================================
|                          INITIALIZATION                         |
===================================================================
Problem size: Nodes=3303, Layers=1, Triangles=6408,
Iter =        0 ; Dt = 0.0000s ; Time =     0.00s ; ETA =     0.00s

Problem solving:

problem.solve()
if not args.nographics:
    plt.show()
  • Height at time = 0.0s
  • Height at time = 1.0s
  • Height at time = 2.0s
  • Height at time = 3.0s
===================================================================
|                            TIME LOOP                            |
===================================================================
Iter =       13 ; Dt = 0.0064s ; Time =     0.09s ; ETA =     4.19s
Iter =       26 ; Dt = 0.0060s ; Time =     0.17s ; ETA =     4.47s
Iter =       40 ; Dt = 0.0057s ; Time =     0.25s ; ETA =     4.57s
Iter =       54 ; Dt = 0.0057s ; Time =     0.33s ; ETA =     4.41s
Iter =       69 ; Dt = 0.0057s ; Time =     0.41s ; ETA =     4.30s
Iter =       83 ; Dt = 0.0057s ; Time =     0.49s ; ETA =     4.23s
Iter =       97 ; Dt = 0.0058s ; Time =     0.57s ; ETA =     4.05s
Iter =      111 ; Dt = 0.0060s ; Time =     0.66s ; ETA =     3.80s
Iter =      124 ; Dt = 0.0062s ; Time =     0.74s ; ETA =     3.65s
Iter =      137 ; Dt = 0.0064s ; Time =     0.82s ; ETA =     3.44s
Iter =      150 ; Dt = 0.0063s ; Time =     0.90s ; ETA =    93.62s
Iter =      163 ; Dt = 0.0063s ; Time =     0.98s ; ETA =     3.29s
Iter =      176 ; Dt = 0.0063s ; Time =     1.06s ; ETA =     3.06s
Iter =      189 ; Dt = 0.0064s ; Time =     1.15s ; ETA =     2.99s
Iter =      202 ; Dt = 0.0065s ; Time =     1.23s ; ETA =     2.85s
Iter =      214 ; Dt = 0.0065s ; Time =     1.31s ; ETA =     2.78s
Iter =      227 ; Dt = 0.0064s ; Time =     1.39s ; ETA =     2.68s
Iter =      240 ; Dt = 0.0064s ; Time =     1.47s ; ETA =     2.62s
Iter =      252 ; Dt = 0.0064s ; Time =     1.55s ; ETA =     2.48s
Iter =      265 ; Dt = 0.0064s ; Time =     1.64s ; ETA =     2.43s
Iter =      278 ; Dt = 0.0065s ; Time =     1.72s ; ETA =     2.35s
Iter =      290 ; Dt = 0.0065s ; Time =     1.80s ; ETA =     2.34s
Iter =      303 ; Dt = 0.0065s ; Time =     1.88s ; ETA =     2.29s
Iter =      315 ; Dt = 0.0065s ; Time =     1.96s ; ETA =     2.23s
Iter =      328 ; Dt = 0.0065s ; Time =     2.04s ; ETA =     2.08s
Iter =      341 ; Dt = 0.0065s ; Time =     2.13s ; ETA =     2.00s
Iter =      353 ; Dt = 0.0066s ; Time =     2.21s ; ETA =     1.89s
Iter =      365 ; Dt = 0.0067s ; Time =     2.29s ; ETA =     1.78s
Iter =      378 ; Dt = 0.0066s ; Time =     2.37s ; ETA =     1.68s
Iter =      390 ; Dt = 0.0066s ; Time =     2.45s ; ETA =     1.67s
Iter =      402 ; Dt = 0.0065s ; Time =     2.53s ; ETA =     1.57s
Iter =      415 ; Dt = 0.0065s ; Time =     2.62s ; ETA =     1.46s
Iter =      427 ; Dt = 0.0065s ; Time =     2.69s ; ETA =     1.39s
Iter =      440 ; Dt = 0.0065s ; Time =     2.78s ; ETA =     1.32s
Iter =      453 ; Dt = 0.0065s ; Time =     2.86s ; ETA =     1.23s
Iter =      465 ; Dt = 0.0065s ; Time =     2.94s ; ETA =     1.14s
Iter =      478 ; Dt = 0.0065s ; Time =     3.03s ; ETA =     1.04s
Iter =      490 ; Dt = 0.0066s ; Time =     3.10s ; ETA =     0.95s
Iter =      502 ; Dt = 0.0066s ; Time =     3.18s ; ETA =     0.86s
Iter =      515 ; Dt = 0.0067s ; Time =     3.27s ; ETA =     0.77s
Iter =      527 ; Dt = 0.0067s ; Time =     3.35s ; ETA =     0.68s
Iter =      539 ; Dt = 0.0068s ; Time =     3.43s ; ETA =     0.59s
Iter =      551 ; Dt = 0.0068s ; Time =     3.51s ; ETA =     0.50s
Iter =      563 ; Dt = 0.0068s ; Time =     3.59s ; ETA =     0.42s
Iter =      575 ; Dt = 0.0068s ; Time =     3.68s ; ETA =     0.33s
Iter =      587 ; Dt = 0.0068s ; Time =     3.76s ; ETA =     0.25s
Iter =      599 ; Dt = 0.0067s ; Time =     3.84s ; ETA =     0.17s
Iter =      612 ; Dt = 0.0067s ; Time =     3.92s ; ETA =     0.08s
Iter =      624 ; Dt = 0.0014s ; Time =     4.00s ; ETA =     0.00s
===================================================================
|                               END                               |
===================================================================
Problem.solve() completed in 6.52889609336853s (wall time)

Total running time of the script: ( 0 minutes 8.855 seconds)

Gallery generated by Sphinx-Gallery