source: development/pyvolution-1d/dry_dam.py @ 3365

Last change on this file since 3365 was 3362, checked in by jakeman, 19 years ago

New files dam2.py and dry_dam.py further test numerical solution
against analytical solution of
Stoker 57. Limiting is still not working correctly.

File size: 3.0 KB
Line 
1import os
2from math import sqrt, pi
3from shallow_water_1d import *
4from Numeric import allclose, array, zeros, ones, Float, take, sqrt
5from config import g, epsilon
6
7def analytical_sol(C,t):
8   
9    #t  = 0.0     # time (s)
10    g  = 9.81    # gravity (m/s^2)
11    h1 = 10.0    # depth upstream (m)
12    h0 = 0.0     # depth downstream (m)
13    L = 2000.0   # length of stream/domain (m)
14    n = len(C)    # number of cells
15
16    u = zeros(n,Float)
17    h = zeros(n,Float)
18    x = C-L/2
19    #x = zeros(n,Float)
20    #for i in range(n):
21    #    x[i] = C[i]-1000.0
22
23    for i in range(n):
24        # Calculate Analytical Solution at time t > 0
25        u3 = 2/3*(sqrt(g*h1)+x[i]/t) 
26        h3 = 4/(9*g)*(sqrt(g*h1)-x[i]/(2*t))*(sqrt(g*h1)-x[i]/(2*t)) 
27
28        if ( x[i] <= -t*sqrt(g*h1) ):
29            u[i] = 0.0 
30            h[i] = h1
31        elif ( x[i] <= 2*t*sqrt(g*h1) ):
32            u[i] = u3
33            h[i] = h3
34        else:
35            u[i] = 0.0 
36            h[i] = h0
37
38    return h
39
40
41def newLinePlot(title='Simple Plot'):
42    import Gnuplot
43    g = Gnuplot.Gnuplot(persist=1)
44    g.title(title)
45    g('set data style linespoints') 
46    g.xlabel('x')
47    g.ylabel('y')
48    return g
49
50def linePlot(g,x1,y1,x2,y2):
51    import Gnuplot
52    plot1 = Gnuplot.PlotItems.Data(x1.flat,y1.flat,with="linespoints")
53    plot2 = Gnuplot.PlotItems.Data(x2.flat,y2.flat, with="lines 3")
54    g.plot(plot1,plot2)
55
56debug = False
57
58print "TEST 1D-SOLUTION III -- DRY BED"
59
60L = 2000.0     # Length of channel (m)
61N = 100        # Number of compuational cells
62cell_len = L/N # Origin = 0.0
63
64points = zeros(N+1,Float)
65for i in range(N+1):
66    points[i] = i*cell_len
67
68domain = Domain(points)
69
70def stage(x):
71    y = zeros(len(x),Float)
72    for i in range(len(x)):
73        if x[i]<=1000.0:
74            y[i] = 10.0
75        else:
76            y[i] = 0.0
77    return y
78
79domain.set_quantity('stage', stage)
80
81domain.order = 2
82domain.default_order = 2
83print "domain.order", domain.order
84
85if (debug == True):
86    area = domain.areas
87    for i in range(len(area)):
88        if area != 20:
89            print "Cell Areas are Wrong"
90           
91    L = domain.quantities['stage'].vertex_values
92    print "Initial Stage"
93    print L
94
95domain.set_boundary({'exterior': Reflective_boundary(domain)})
96
97X = domain.vertices
98C = domain.centroids
99plot1 = newLinePlot("Stage")
100
101import time
102t0 = time.time()
103yieldstep = 1.0
104finaltime = 30.0
105for t in domain.evolve(yieldstep = yieldstep, finaltime = finaltime):
106    domain.write_time()
107    if t > 0.0:
108        StageQ = domain.quantities['stage'].vertex_values
109        y = analytical_sol(X.flat,domain.time)
110        linePlot(plot1,X,StageQ,X,y)
111    raw_input('press_return')
112    #pass
113
114print 'That took %.2f seconds' %(time.time()-t0)
115
116C = domain.quantities['stage'].centroid_values
117
118if (debug == True):
119    L = domain.quantities['stage'].vertex_values
120    print "Final Stage Vertex Values"
121    print L
122    print "Final Stage Centroid Values"
123    print C
124
125#f = file('test_solution_I.out', 'w')
126#for i in range(N):
127#    f.write(str(C[i]))
128#    f.write("\n")
129#f.close
130
Note: See TracBrowser for help on using the repository browser.