Note
Click here to download the full example code
Riemann problem 2d¶
This case is a simple Riemann problem in 2D. Geometry is a square which is split in four. On each part initial height is set to a constant value. Initial velocity is null and all sides are walls.
import os, sys
import argparse
import pylab as plt
import freshkiss3d as fk
import freshkiss3d.extra.plots as fk_plt
parser = argparse.ArgumentParser()
parser.add_argument('--nographics', action='store_true')
args = parser.parse_args()
#sphinx_gallery_thumbnail_number = 2
Time loop:¶
simutime = fk.SimuTime(final_time=5., time_iteration_max=2000, second_order=True)
create_figure_scheduler = fk.schedules(times=[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/square.msh')
if not args.nographics:
fk_plt.plot_mesh(triangular_mesh, plot_labels=False)
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 = 1.
if x < 2.5 and y < 2.5: h = 4.
if x < 2.5 and y > 2.5: h = 3.
if x > 2.5 and y < 2.5: h = 1.
if x > 2.5 and y > 2.5: h = 2.
return h
if not args.nographics:
fk_plt.plot_init_3d(triangular_mesh, H_0)
primitives = fk.Primitive(triangular_mesh, layer, height_funct=H_0)
Boundary conditions:¶
slides = [fk.Slide(ref=r) for r in [1, 2, 3, 4]]
Problem definition:¶
problem = fk.Problem(simutime, triangular_mesh, layer, primitives,
slides=slides,
numerical_parameters={'ipres':True, 'space_second_order':True},
custom_funct={'plot':fk_plt.plot_height_2d_3d},
custom_funct_scheduler=create_figure_scheduler)
===================================================================
| INITIALIZATION |
===================================================================
Problem size: Nodes=2601, Layers=1, Triangles=5000,
Iter = 0 ; Dt = 0.0000s ; Time = 0.00s ; ETA = 0.00s
Problem solving:¶
problem.solve()
if not args.nographics:
plt.show()
===================================================================
| TIME LOOP |
===================================================================
Iter = 40 ; Dt = 0.0026s ; Time = 0.10s ; ETA = 16.77s
Iter = 79 ; Dt = 0.0026s ; Time = 0.20s ; ETA = 15.05s
Iter = 119 ; Dt = 0.0026s ; Time = 0.31s ; ETA = 16.20s
Iter = 158 ; Dt = 0.0026s ; Time = 0.41s ; ETA = 14.10s
Iter = 196 ; Dt = 0.0028s ; Time = 0.51s ; ETA = 13.11s
Iter = 233 ; Dt = 0.0027s ; Time = 0.61s ; ETA = 13.79s
Iter = 270 ; Dt = 0.0028s ; Time = 0.71s ; ETA = 12.95s
Iter = 305 ; Dt = 0.0031s ; Time = 0.82s ; ETA = 11.48s
Iter = 337 ; Dt = 0.0032s ; Time = 0.92s ; ETA = 10.72s
Iter = 370 ; Dt = 0.0032s ; Time = 1.02s ; ETA = 10.09s
Iter = 402 ; Dt = 0.0032s ; Time = 1.12s ; ETA = 10.07s
Iter = 437 ; Dt = 0.0025s ; Time = 1.23s ; ETA = 11.86s
Iter = 480 ; Dt = 0.0023s ; Time = 1.33s ; ETA = 12.86s
Iter = 523 ; Dt = 0.0024s ; Time = 1.43s ; ETA = 11.96s
Iter = 566 ; Dt = 0.0024s ; Time = 1.53s ; ETA = 11.41s
Iter = 608 ; Dt = 0.0024s ; Time = 1.63s ; ETA = 11.19s
Iter = 650 ; Dt = 0.0024s ; Time = 1.74s ; ETA = 10.81s
Iter = 691 ; Dt = 0.0025s ; Time = 1.84s ; ETA = 10.27s
Iter = 732 ; Dt = 0.0025s ; Time = 1.94s ; ETA = 9.90s
Iter = 772 ; Dt = 0.0026s ; Time = 2.04s ; ETA = 9.28s
Iter = 811 ; Dt = 0.0026s ; Time = 2.14s ; ETA = 8.80s
Iter = 849 ; Dt = 0.0027s ; Time = 2.25s ; ETA = 8.27s
Iter = 887 ; Dt = 0.0028s ; Time = 2.35s ; ETA = 7.76s
Iter = 923 ; Dt = 0.0029s ; Time = 2.45s ; ETA = 7.22s
Iter = 958 ; Dt = 0.0029s ; Time = 2.55s ; ETA = 6.72s
Iter = 993 ; Dt = 0.0030s ; Time = 2.66s ; ETA = 6.41s
Iter = 1026 ; Dt = 0.0031s ; Time = 2.76s ; ETA = 5.90s
Iter = 1059 ; Dt = 0.0031s ; Time = 2.86s ; ETA = 5.55s
Iter = 1092 ; Dt = 0.0032s ; Time = 2.96s ; ETA = 5.11s
Iter = 1123 ; Dt = 0.0032s ; Time = 3.06s ; ETA = 4.91s
Iter = 1155 ; Dt = 0.0032s ; Time = 3.17s ; ETA = 7.49s
Iter = 1187 ; Dt = 0.0032s ; Time = 3.27s ; ETA = 4.34s
Iter = 1219 ; Dt = 0.0032s ; Time = 3.37s ; ETA = 4.08s
Iter = 1251 ; Dt = 0.0032s ; Time = 3.47s ; ETA = 3.86s
Iter = 1283 ; Dt = 0.0029s ; Time = 3.57s ; ETA = 3.95s
Iter = 1329 ; Dt = 0.0020s ; Time = 3.67s ; ETA = 5.28s
Iter = 1377 ; Dt = 0.0022s ; Time = 3.78s ; ETA = 4.33s
Iter = 1421 ; Dt = 0.0024s ; Time = 3.88s ; ETA = 3.81s
Iter = 1464 ; Dt = 0.0024s ; Time = 3.98s ; ETA = 3.31s
Iter = 1505 ; Dt = 0.0025s ; Time = 4.08s ; ETA = 2.92s
Iter = 1545 ; Dt = 0.0026s ; Time = 4.19s ; ETA = 2.53s
Iter = 1584 ; Dt = 0.0027s ; Time = 4.29s ; ETA = 2.14s
Iter = 1622 ; Dt = 0.0027s ; Time = 4.39s ; ETA = 1.81s
Iter = 1659 ; Dt = 0.0028s ; Time = 4.49s ; ETA = 1.51s
Iter = 1695 ; Dt = 0.0028s ; Time = 4.59s ; ETA = 1.18s
Iter = 1731 ; Dt = 0.0029s ; Time = 4.70s ; ETA = 0.84s
Iter = 1766 ; Dt = 0.0029s ; Time = 4.80s ; ETA = 0.55s
Iter = 1801 ; Dt = 0.0030s ; Time = 4.90s ; ETA = 0.28s
Iter = 1835 ; Dt = 0.0020s ; Time = 5.00s ; ETA = 0.00s
===================================================================
| END |
===================================================================
Problem.solve() completed in 21.71368384361267s (wall time)
Total running time of the script: ( 0 minutes 24.514 seconds)