1 | program sfbasic |
---|
2 | c============================================================ |
---|
3 | c Combine evolution of w = h+z and h for shallow flows |
---|
4 | c============================================================ |
---|
5 | |
---|
6 | implicit none |
---|
7 | |
---|
8 | integer md,nxmax,nxd, mn |
---|
9 | |
---|
10 | parameter (md=3, nxmax=5000, nxd=nxmax+2*md, mn=2) |
---|
11 | real*8 h(nxd), uh(nxd), w(nxd) |
---|
12 | real*8 z(nxd),z2(nxd) |
---|
13 | real*8 x(nxd),x2(nxd) |
---|
14 | |
---|
15 | integer terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
16 | common /setup/ terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
17 | |
---|
18 | integer nlow,nhigh,nx,nxx,dn,version,isteps |
---|
19 | common /ivariables/ nlow,nhigh,nx,nxx,version,isteps |
---|
20 | |
---|
21 | real*8 toll,cfl,theta |
---|
22 | real*8 g,hf,dx,dt,dtmax |
---|
23 | real*8 Tout,tf,tc,uhf |
---|
24 | real*8 K,limith,length |
---|
25 | real*8 manning,hl,hr,zb |
---|
26 | common /rvariables/ toll,cfl,theta, |
---|
27 | . g,hf,dx,dt,dtmax, |
---|
28 | . Tout,tf,tc,uhf, |
---|
29 | . K,limith,length, |
---|
30 | . manning,hl,hr |
---|
31 | |
---|
32 | integer io, i |
---|
33 | |
---|
34 | call util_open_files |
---|
35 | |
---|
36 | open(1,file='MacDonald.bed') |
---|
37 | c |
---|
38 | c Parameters |
---|
39 | c |
---|
40 | nx = 101 |
---|
41 | dtmax = 0.5d0 |
---|
42 | Tout = 1.0d0 |
---|
43 | tf = 1500.0d0 |
---|
44 | length = 1000.0d0 |
---|
45 | |
---|
46 | nlow = md+2 |
---|
47 | nhigh = md+nx |
---|
48 | nxx = nx+2*md |
---|
49 | |
---|
50 | limith = 1 |
---|
51 | |
---|
52 | dx = 10000.d0/dble(nx-1) |
---|
53 | K = 0.0d0 |
---|
54 | |
---|
55 | if (nx.gt.nxmax) write(*,*)'nx too large' |
---|
56 | toll = 1.0d-10 |
---|
57 | cfl = 0.2d0 |
---|
58 | theta = 1.2d0 |
---|
59 | version = 10 |
---|
60 | g = 9.81d0 |
---|
61 | hl = 5.d0 |
---|
62 | hr = 0.d0 |
---|
63 | tc = 0.d0 |
---|
64 | |
---|
65 | c |
---|
66 | c ---- Setup Terrain, Boundary Values and Initial Values |
---|
67 | c |
---|
68 | |
---|
69 | c ---- Non-uniform terrain from MacDonald |
---|
70 | do i = 1,nxx |
---|
71 | read(1,*)x2(i),z2(i) |
---|
72 | z2(i) = z2(i)*5.d0 |
---|
73 | x(i) = x2(i) - dx/2.d0 |
---|
74 | end do |
---|
75 | do i = 2,nxx |
---|
76 | z(i) = (z2(i-1)+z2(i))/2.d0 |
---|
77 | end do |
---|
78 | z(1) = (3.d0*z2(1) - z2(2))/2.d0 |
---|
79 | z(nxx) = (3.d0*z2(nxx-1) - z2(nxx-2))/2.d0 |
---|
80 | |
---|
81 | call boundary_apply(x,h,uh,z) |
---|
82 | |
---|
83 | c ---- Rectangular pulse initial condition |
---|
84 | do i = md+1, nx+md |
---|
85 | if (x(i).lt.100.d0.or.x(i).gt.200.d0) then |
---|
86 | h(i) = hr |
---|
87 | else |
---|
88 | h(i) = hl |
---|
89 | endif |
---|
90 | uh(i) = 0.0d0 |
---|
91 | end do |
---|
92 | c |
---|
93 | c ---- Setup for Time looping |
---|
94 | c |
---|
95 | isteps = 0 |
---|
96 | tc = 0.0d0 |
---|
97 | call util_dump_solution(h,uh,z) |
---|
98 | c |
---|
99 | c ---- Evolve |
---|
100 | c |
---|
101 | call evolver_euler(x,h,uh,z,z2) |
---|
102 | c |
---|
103 | c ---- Close down |
---|
104 | c |
---|
105 | write(15,*)isteps,nxx,version |
---|
106 | do i=1,nxx |
---|
107 | write(16,*) x(i) |
---|
108 | enddo |
---|
109 | |
---|
110 | do i=5,nxx-3 |
---|
111 | write(25,*) x(i) |
---|
112 | enddo |
---|
113 | |
---|
114 | call util_close_files |
---|
115 | |
---|
116 | stop |
---|
117 | end |
---|
118 | c=========================================================== |
---|
119 | c Shallow Water Flux |
---|
120 | c=========================================================== |
---|
121 | subroutine flux_huh(f1,f2,h0,uh0,em_x,l1,l2,g,toll) |
---|
122 | implicit none |
---|
123 | |
---|
124 | real*8 f1,f2,h0,uh0,z,l1,l2 |
---|
125 | real*8 em_x,g,toll |
---|
126 | |
---|
127 | real*8 u,cvel,h,uh,w |
---|
128 | integer i |
---|
129 | |
---|
130 | c----------------------------------- |
---|
131 | |
---|
132 | if(em_x.lt.1.0d-15) then |
---|
133 | em_x = 1.d-15 |
---|
134 | endif |
---|
135 | |
---|
136 | uh = uh0 |
---|
137 | h = h0 |
---|
138 | |
---|
139 | if(h.le.toll)then |
---|
140 | u = 0.d0 |
---|
141 | h = 0.0d0 |
---|
142 | uh = 0.0d0 |
---|
143 | else |
---|
144 | u = uh/h |
---|
145 | endif |
---|
146 | |
---|
147 | cvel = dsqrt( g*h ) |
---|
148 | l1 = u - cvel |
---|
149 | l2 = u + cvel |
---|
150 | em_x = dmax1( em_x, dabs(u) + cvel ) |
---|
151 | f1 = uh |
---|
152 | f2 = uh*u + 0.5d0*g*h**2 |
---|
153 | |
---|
154 | return |
---|
155 | end |
---|
156 | c=========================================================== |
---|
157 | c Shallow Water Flux |
---|
158 | c=========================================================== |
---|
159 | subroutine flux_huh_old(f1,f2,h0,uh0,z,em_x,l1,l2,g,toll) |
---|
160 | implicit none |
---|
161 | |
---|
162 | real*8 f1,f2,h0,uh0,z,l1,l2 |
---|
163 | real*8 em_x,g,toll |
---|
164 | |
---|
165 | real*8 u,cvel,h,uh,w |
---|
166 | integer i |
---|
167 | |
---|
168 | c----------------------------------- |
---|
169 | |
---|
170 | if(em_x.lt.1.0d-15) then |
---|
171 | em_x = 1.d-15 |
---|
172 | endif |
---|
173 | |
---|
174 | w = h0+z |
---|
175 | uh = uh0 |
---|
176 | h = h0 |
---|
177 | |
---|
178 | if(h.le.toll)then |
---|
179 | u = 0.d0 |
---|
180 | h = 0.0d0 |
---|
181 | uh = 0.0d0 |
---|
182 | else |
---|
183 | u = uh / h |
---|
184 | endif |
---|
185 | if (abs(u).gt.1.0d1) then |
---|
186 | u = dsign(1.0d0,u)*1.0d3 |
---|
187 | uh = u*h |
---|
188 | end if |
---|
189 | cvel = dsqrt( g*h ) |
---|
190 | l1 = u - cvel |
---|
191 | l2 = u + cvel |
---|
192 | em_x = dmax1( em_x, dabs(u) + cvel ) |
---|
193 | f1 = uh |
---|
194 | f2 = uh*u + 0.5d0*g*h**2 |
---|
195 | |
---|
196 | return |
---|
197 | end |
---|
198 | c=========================================================== |
---|
199 | c Shallow Water Flux using h+z and uh variables |
---|
200 | c=========================================================== |
---|
201 | subroutine flux_wuh(f1,f2,w0,uh0,z,em_x,l1,l2,g,toll) |
---|
202 | implicit none |
---|
203 | |
---|
204 | real*8 f1,f2,w0,uh0,z |
---|
205 | real*8 g,toll |
---|
206 | |
---|
207 | real*8 h,u,cvel,uh,w,l1,l2,em_x |
---|
208 | integer i |
---|
209 | |
---|
210 | c----------------------------------- |
---|
211 | |
---|
212 | if(em_x.lt.1.0d-15) then |
---|
213 | em_x = 1.d-15 |
---|
214 | endif |
---|
215 | |
---|
216 | w = w0 |
---|
217 | uh = uh0 |
---|
218 | |
---|
219 | h = w - z |
---|
220 | if(h.le.toll)then |
---|
221 | u = 0.d0 |
---|
222 | h = 0.0d0 |
---|
223 | uh = 0.0d0 |
---|
224 | else |
---|
225 | u = uh / h |
---|
226 | endif |
---|
227 | if (abs(u).gt.1.0d3) then |
---|
228 | u = dsign(1.0d0,u)*1.0d3 |
---|
229 | uh = u*h |
---|
230 | end if |
---|
231 | cvel = dsqrt( g*h ) |
---|
232 | l1 = u - cvel |
---|
233 | l2 = u + cvel |
---|
234 | em_x = dmax1( em_x, dabs(u) + cvel ) |
---|
235 | f1 = uh |
---|
236 | f2 = uh*u + 0.5d0*g*h**2 |
---|
237 | |
---|
238 | return |
---|
239 | end |
---|
240 | c=============================================================== |
---|
241 | c array limiter |
---|
242 | c=============================================================== |
---|
243 | subroutine array_limit_h(h,uh,z2,z,uh_l,uh_r,h_l,h_r,h_2) |
---|
244 | |
---|
245 | implicit none |
---|
246 | integer md,nxmax,nxd, mn |
---|
247 | parameter (md=3, nxmax=5000, nxd=nxmax+2*md, mn=2) |
---|
248 | real*8 h(nxd), uh(nxd), w(nxd) |
---|
249 | real*8 z(nxd),z2(nxd) |
---|
250 | real*8 x(nxd),x2(nxd) |
---|
251 | |
---|
252 | integer terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
253 | common /setup/ terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
254 | |
---|
255 | integer nlow,nhigh,nx,nxx,dn,version,isteps |
---|
256 | common /ivariables/ nlow,nhigh,nx,nxx,version,isteps |
---|
257 | |
---|
258 | real*8 toll,cfl,theta |
---|
259 | real*8 g,hf,dx,dt,dtmax |
---|
260 | real*8 Tout,tf,tc,uhf |
---|
261 | real*8 K,limith,length |
---|
262 | real*8 manning,hl,hr |
---|
263 | common /rvariables/ toll,cfl,theta, |
---|
264 | . g,hf,dx,dt,dtmax, |
---|
265 | . Tout,tf,tc,uhf, |
---|
266 | . K,limith,length, |
---|
267 | . manning,hl,hr |
---|
268 | |
---|
269 | real*8 h_2(nxd) |
---|
270 | real*8 h_l(nxd),h_r(nxd) |
---|
271 | real*8 h1_l(nxd), h1_r(nxd) |
---|
272 | real*8 uh_l(nxd), uh_r(nxd) |
---|
273 | |
---|
274 | integer sslope(nxd) |
---|
275 | real*8 xmin,xmic,xminmod,a,b,c,t |
---|
276 | |
---|
277 | xmin(a,b) = 0.5d0*(dsign(1.d0,a)+dsign(1.d0,b))* |
---|
278 | . dmin1(dabs(a),dabs(b)) |
---|
279 | xmic(t,a,b) = xmin(t*xmin(a,b), 0.5d0*(a+b) ) |
---|
280 | xminmod(a,b,c) = xmin(a,xmin(b,c)) |
---|
281 | |
---|
282 | |
---|
283 | integer i |
---|
284 | real*8 dh, hmax, hmin , tt, cc |
---|
285 | real*8 w0,w1,w2,d0h,d0w,d1h,d1w,d2h,d2w,ddh,ddw |
---|
286 | real*8 h0,h1,h2,hz |
---|
287 | real*8 ss0,ss1,ss |
---|
288 | real*8 dz0,dz1,dz2 |
---|
289 | real*8 hx1,hx2,h_x |
---|
290 | real*8 froud,d1,d2 |
---|
291 | real*8 gh3,uh2,duh, w_x,dz |
---|
292 | |
---|
293 | c------------------------------------ |
---|
294 | g = 9.81d0 |
---|
295 | cc = 0.1d0 |
---|
296 | ss0 = 0.01d0 |
---|
297 | ss1 = 0.01d0 |
---|
298 | |
---|
299 | do i=md,nx+md+2 |
---|
300 | |
---|
301 | h_2(i) = h(i) |
---|
302 | h0 = h(i-1) |
---|
303 | h1 = h(i) |
---|
304 | h2 = h(i+1) |
---|
305 | |
---|
306 | dz0 = abs(z2(i-1)-z2(i-2)) |
---|
307 | dz1 = abs(z2(i) -z2(i-1)) |
---|
308 | dz2 = abs(z2(i+1)-z2(i) ) |
---|
309 | |
---|
310 | w0 = h0 + z(i-1) |
---|
311 | w1 = h1 + z(i) |
---|
312 | w2 = h2 + z(i+1) |
---|
313 | |
---|
314 | d1h = h1 - h0 |
---|
315 | d2h = h2 - h1 |
---|
316 | d1w = w1 - w0 |
---|
317 | d2w = w2 - w1 |
---|
318 | dz = (z2(i) - z2(i-1)) |
---|
319 | |
---|
320 | c ------- W = H+Z |
---|
321 | c h_x = xmin( d1w, d2w) - dz |
---|
322 | h_x = xmic(theta, d1h, d2h) |
---|
323 | |
---|
324 | c ------- Return the left and right values of h in an interval |
---|
325 | h_l(i) = h(i) - 0.5d0*h_x |
---|
326 | h_r(i) = h(i) + 0.5d0*h_x |
---|
327 | |
---|
328 | c ------- Test for h_l or h_r negative |
---|
329 | c$$$ hmin = dmin1(h0,h1,h2) |
---|
330 | c$$$ if ( h_l(i) .le. hmin ) then |
---|
331 | c$$$ h_l(i) = hmin |
---|
332 | c$$$ h_r(i) = h(i) + (h(i) - h_l(i)) |
---|
333 | c$$$ elseif ( h_r(i) .le. hmin ) then |
---|
334 | c$$$ h_r(i) = hmin |
---|
335 | c$$$ h_l(i) = h(i) + (h(i) - h_r(i)) |
---|
336 | c$$$ endif |
---|
337 | enddo |
---|
338 | c |
---|
339 | c ---- Enforce a minimum height |
---|
340 | c |
---|
341 | c$$$ do i=md,nx+md+2 |
---|
342 | c$$$ hmax = dmax1(h_l(i),h_r(i)) |
---|
343 | c$$$ if (hmax .lt. 1.0d-3 ) then |
---|
344 | c$$$c h_r(i) = 0.0d0 |
---|
345 | c$$$ uh_r(i) = 0.0d0 |
---|
346 | c$$$c h_l(i) = 0.0d0 |
---|
347 | c$$$ uh_l(i) = 0.0d0 |
---|
348 | c$$$c h(i) = 0.0d0 |
---|
349 | c$$$ uh(i) = 0.0d0 |
---|
350 | c$$$ endif |
---|
351 | c$$$ enddo |
---|
352 | c |
---|
353 | c ---- Try for a flat subinterval |
---|
354 | cc |
---|
355 | c do i=md,nx+md+2 |
---|
356 | c hmax = dmax1(h_l(i),h_r(i)) |
---|
357 | c dz = (z2(i) - z2(i-1)) |
---|
358 | c |
---|
359 | c if (h_l(i) .lt. 0.01d0*h_r(i) ) then |
---|
360 | c h_r(i) = sqrt(abs(2.0d0*dz*h(i))) |
---|
361 | c endif |
---|
362 | c |
---|
363 | c if (h_r(i) .lt. 0.01d0*h_l(i) ) then |
---|
364 | c h_l(i) = sqrt(abs(2.0d0*dz*h(i))) |
---|
365 | c endif |
---|
366 | c enddo |
---|
367 | |
---|
368 | return |
---|
369 | end |
---|
370 | c=============================================================== |
---|
371 | c array limiter |
---|
372 | c=============================================================== |
---|
373 | subroutine array_limit(u,u_l,u_r) |
---|
374 | |
---|
375 | implicit none |
---|
376 | integer md,nxmax,nxd, mn |
---|
377 | parameter (md=3, nxmax=5000, nxd=nxmax+2*md, mn=2) |
---|
378 | real*8 h(nxd), uh(nxd), w(nxd) |
---|
379 | real*8 z(nxd),z2(nxd) |
---|
380 | real*8 x(nxd),x2(nxd) |
---|
381 | |
---|
382 | integer terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
383 | common /setup/ terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
384 | |
---|
385 | integer nlow,nhigh,nx,nxx,dn,version,isteps |
---|
386 | common /ivariables/ nlow,nhigh,nx,nxx,version,isteps |
---|
387 | |
---|
388 | real*8 toll,cfl,theta |
---|
389 | real*8 g,hf,dx,dt,dtmax |
---|
390 | real*8 Tout,tf,tc,uhf |
---|
391 | real*8 K,limith,length |
---|
392 | real*8 manning,hl,hr |
---|
393 | common /rvariables/ toll,cfl,theta, |
---|
394 | . g,hf,dx,dt,dtmax, |
---|
395 | . Tout,tf,tc,uhf, |
---|
396 | . K,limith,length, |
---|
397 | . manning,hl,hr |
---|
398 | |
---|
399 | real*8 u(nxd), u_l(nxd),u_r(nxd) |
---|
400 | real*8 u_x,d1,d2,d0,ux1 |
---|
401 | real*8 u1_l(nxd),u1_r(nxd) |
---|
402 | integer i |
---|
403 | |
---|
404 | real*8 a,b,t,xmin,xmic |
---|
405 | |
---|
406 | xmin(a,b) = 0.5d0*(dsign(1.d0,a)+dsign(1.d0,b))* |
---|
407 | . dmin1(dabs(a),dabs(b)) |
---|
408 | xmic(t,a,b) = xmin(t*xmin(a,b), 0.5d0*(a+b) ) |
---|
409 | |
---|
410 | c------------------------------------ |
---|
411 | |
---|
412 | do i=md,nx+md+2 |
---|
413 | d1 = u(i) - u(i-1) |
---|
414 | d2 = u(i+1) - u(i) |
---|
415 | u_x = xmic( theta, d1, d2 ) |
---|
416 | c u_x = xmin( d1, d2 ) |
---|
417 | u_l(i) = u(i) - 0.5d0*u_x |
---|
418 | u_r(i) = u(i) + 0.5d0*u_x |
---|
419 | enddo |
---|
420 | |
---|
421 | return |
---|
422 | end |
---|
423 | c=========================================================== |
---|
424 | c Output Routine |
---|
425 | c=========================================================== |
---|
426 | subroutine util_dump_solution(h,uh,z) |
---|
427 | |
---|
428 | implicit none |
---|
429 | integer md,nxmax,nxd, mn |
---|
430 | parameter (md=3, nxmax=5000, nxd=nxmax+2*md, mn=2) |
---|
431 | real*8 h(nxd), uh(nxd), w(nxd) |
---|
432 | real*8 z(nxd),z2(nxd) |
---|
433 | real*8 x(nxd),x2(nxd) |
---|
434 | |
---|
435 | integer terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
436 | common /setup/ terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
437 | |
---|
438 | integer nlow,nhigh,nx,nxx,dn,version,isteps |
---|
439 | common /ivariables/ nlow,nhigh,nx,nxx,version,isteps |
---|
440 | |
---|
441 | real*8 toll,cfl,theta |
---|
442 | real*8 g,hf,dx,dt,dtmax |
---|
443 | real*8 Tout,tf,tc,uhf |
---|
444 | real*8 K,limith,length |
---|
445 | real*8 manning,hl,hr |
---|
446 | common /rvariables/ toll,cfl,theta, |
---|
447 | . g,hf,dx,dt,dtmax, |
---|
448 | . Tout,tf,tc,uhf, |
---|
449 | . K,limith,length, |
---|
450 | . manning,hl,hr |
---|
451 | |
---|
452 | integer ii,i |
---|
453 | |
---|
454 | write(12,*)tc |
---|
455 | |
---|
456 | isteps = isteps + 1 |
---|
457 | write(*,*)isteps, tc |
---|
458 | |
---|
459 | do i = 1, nx+2*md |
---|
460 | c------ Dump h |
---|
461 | write(13,200)sngl(h(i)) |
---|
462 | c------ Dump uh |
---|
463 | write(14,200)sngl(uh(i)) |
---|
464 | c------ Dump w |
---|
465 | write(11,200)sngl(h(i)+z(i)) |
---|
466 | 200 format(e20.8) |
---|
467 | end do |
---|
468 | |
---|
469 | return |
---|
470 | end |
---|
471 | c=============================================================== |
---|
472 | c The evolver using euler's method |
---|
473 | c=============================================================== |
---|
474 | subroutine evolver_euler(x,h,uh,z,z2) |
---|
475 | |
---|
476 | implicit none |
---|
477 | integer md,nxmax,nxd, mn |
---|
478 | parameter (md=3, nxmax=5000, nxd=nxmax+2*md, mn=2) |
---|
479 | real*8 h(nxd), uh(nxd), w(nxd) |
---|
480 | real*8 z(nxd),z2(nxd) |
---|
481 | real*8 x(nxd),x2(nxd) |
---|
482 | |
---|
483 | integer terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
484 | c common /setup/ terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
485 | |
---|
486 | integer nlow,nhigh,nx,nxx,dn,version,isteps |
---|
487 | common /ivariables/ nlow,nhigh,nx,nxx,version,isteps |
---|
488 | |
---|
489 | real*8 toll,cfl,theta |
---|
490 | real*8 g,hf,dx,dt,dtmax |
---|
491 | real*8 Tout,tf,tc,uhf |
---|
492 | real*8 K,limith,length |
---|
493 | real*8 manning,hl,hr |
---|
494 | common /rvariables/ toll,cfl,theta, |
---|
495 | . g,hf,dx,dt,dtmax, |
---|
496 | . Tout,tf,tc,uhf, |
---|
497 | . K,limith,length, |
---|
498 | . manning,hl,hr |
---|
499 | |
---|
500 | real*8 h_new(nxd),uh_new(nxd) |
---|
501 | real*8 h_h(nxd),uh_h(nxd) |
---|
502 | real*8 h_h_x(nxd),uh_h_x(nxd) |
---|
503 | real*8 h_x(nxd), uh_x(nxd) |
---|
504 | real*8 f_h(nxd), f_uh(nxd) |
---|
505 | real*8 f_h_l(nxd), f_uh_l(nxd) |
---|
506 | real*8 f_h_r(nxd), f_uh_r(nxd) |
---|
507 | real*8 z_l(nxd),z_r(nxd) |
---|
508 | real*8 h_2(nxd), uh_2(nxd) |
---|
509 | |
---|
510 | real*8 h_l(nxd), h_r(nxd), uh_l(nxd), uh_r(nxd) |
---|
511 | real*8 z2l, z2r, zl, zr |
---|
512 | |
---|
513 | integer io, i, jj, ipos,ip |
---|
514 | real*8 factor, Ki |
---|
515 | |
---|
516 | integer istop,iout,nt,ii,oi,i2,m |
---|
517 | real*8 a,b,t,Tleft,em_x,dtcdx,dtcdx2,dtcdx3 |
---|
518 | real*8 den, xmt,pre, Tnext, em_old, hh |
---|
519 | real*8 aterm1 ,term1, aterm2 ,term2, aterm3 ,term3 |
---|
520 | |
---|
521 | integer kk |
---|
522 | real*8 umax, umax2, vmin, q1,q2,q3,h1,h2 |
---|
523 | real*8 Kcal,k1,k2,za,zb,pi |
---|
524 | |
---|
525 | real*8 fh_l, fuh_l, fh_r, fuh_r |
---|
526 | real*8 l1_l, l2_l, l1_r, l2_r |
---|
527 | real*8 a_max, a_min, aa, axa |
---|
528 | |
---|
529 | c---------------------------------------------------------- |
---|
530 | c |
---|
531 | pi = 2.d0*dtan(1.d0) |
---|
532 | factor = 1.0d0 |
---|
533 | em_old = 1.0d0 |
---|
534 | tc = 0.d0 |
---|
535 | istop = 0 |
---|
536 | iout = 0 |
---|
537 | nt =0 |
---|
538 | Tleft = dmin1(Tout,tf) |
---|
539 | Tnext = dmin1(Tout,tf) |
---|
540 | dt = dmin1(dtmax,Tleft) |
---|
541 | |
---|
542 | do i=1,nxx |
---|
543 | f_h(i) = 0.0d0 |
---|
544 | f_uh(i) = 0.0d0 |
---|
545 | h_l(i) = h(i) |
---|
546 | h_r(i) = h(i) |
---|
547 | uh_l(i) = uh(i) |
---|
548 | uh_r(i) = uh(i) |
---|
549 | enddo |
---|
550 | |
---|
551 | call boundary_apply(x,h,uh,z) |
---|
552 | call array_limit (uh,uh_l,uh_r) |
---|
553 | call array_limit_h(h,uh,z2,z,uh_l,uh_r,h_l,h_r,h_2) |
---|
554 | |
---|
555 | do i=1,nxx |
---|
556 | if ( h(i) .lt. 0.0d0 ) then |
---|
557 | write(*,*)'i=',i,' h(i) = ',h(i) |
---|
558 | endif |
---|
559 | enddo |
---|
560 | |
---|
561 | c------------------------------------------------------ |
---|
562 | c |
---|
563 | do while (.true.) |
---|
564 | c write(*,*)nt |
---|
565 | nt = nt+1 |
---|
566 | ii = 1 |
---|
567 | |
---|
568 | call boundary_apply(x,h,uh,z) |
---|
569 | |
---|
570 | call array_limit (uh,uh_l,uh_r) |
---|
571 | call array_limit_h(h,uh,z2,z,uh_l,uh_r,h_l,h_r,h_2) |
---|
572 | |
---|
573 | do i=1,nxx |
---|
574 | if ( h(i) .lt. 0.0d0 ) then |
---|
575 | write(*,*)'i=',i,' h(i) = ',h(i),' nt = ',nt |
---|
576 | endif |
---|
577 | enddo |
---|
578 | |
---|
579 | c------------------------- |
---|
580 | c Flux Calculation |
---|
581 | c------------------------- |
---|
582 | call numerical_flux_central(h,uh,z2,z,uh_l,uh_r,h_l,h_r,h_2, |
---|
583 | . f_h,f_uh, em_x) |
---|
584 | |
---|
585 | c----------------------------------------------------- |
---|
586 | c Compute time step size according to the input CFL # |
---|
587 | c----------------------------------------------------- |
---|
588 | 1002 dt = dmin1(factor*cfl*dx/em_x,dtmax,Tout) |
---|
589 | c write(*,*)'dt = ',dt,' factor = ',factor |
---|
590 | if( ( tc + dt ) .ge. Tnext ) then |
---|
591 | dt = (Tnext - tc ) |
---|
592 | iout = 1 |
---|
593 | Tnext = Tnext+Tout |
---|
594 | end if |
---|
595 | if (Tnext.gt.tf ) then |
---|
596 | istop = 1 |
---|
597 | endif |
---|
598 | |
---|
599 | c------------------------------------------------- |
---|
600 | c Compute the values of 'u' at the next time level |
---|
601 | c------------------------------------------------- |
---|
602 | ipos = 0 |
---|
603 | |
---|
604 | do i = md+2, md+nx |
---|
605 | |
---|
606 | h_new(i) = h(i) - dt/dx*( f_h(i) - f_h(i-1) ) |
---|
607 | uh_new(i) = uh(i) - dt/dx*( f_uh(i) - f_uh(i-1)) |
---|
608 | . - dt/dx*g*(z2(i)-z2(i-1))*(h(i)) |
---|
609 | |
---|
610 | if (h_new(i).lt.-1.0d-5) then |
---|
611 | write(*,*)'h_new(',i,') = ',h_new(i) |
---|
612 | ipos=1 |
---|
613 | endif |
---|
614 | |
---|
615 | if (h_new(i).le.0.0d0) then |
---|
616 | h_new(i) = 0.0d0 |
---|
617 | c - my statement |
---|
618 | uh_new(i) = 0.d0 |
---|
619 | endif |
---|
620 | |
---|
621 | end do |
---|
622 | |
---|
623 | if(ipos.eq.1) then |
---|
624 | factor = factor/2.0d0 |
---|
625 | write(*,*)'factor = ',factor |
---|
626 | goto 1002 |
---|
627 | endif |
---|
628 | |
---|
629 | do i = md+2,nx+md |
---|
630 | h(i) = h_new(i) |
---|
631 | uh(i) = uh_new(i) |
---|
632 | enddo |
---|
633 | |
---|
634 | factor = 1.0d0 |
---|
635 | tc = tc + dt |
---|
636 | |
---|
637 | if(iout.eq.1) then |
---|
638 | call util_dump_solution(h,uh,z) |
---|
639 | |
---|
640 | call boundary_apply(x,h,uh,z) |
---|
641 | call array_limit (uh,uh_l,uh_r) |
---|
642 | call array_limit_h(h,uh,z2,z,uh_l,uh_r,h_l,h_r,h_2) |
---|
643 | iout = 0 |
---|
644 | endif |
---|
645 | if(istop.eq.1) go to 1001 |
---|
646 | c |
---|
647 | c End timestep |
---|
648 | c |
---|
649 | end do |
---|
650 | c |
---|
651 | 1001 write(*,*)isteps |
---|
652 | |
---|
653 | return |
---|
654 | |
---|
655 | end |
---|
656 | c=========================================================== |
---|
657 | c Boundary Conditions |
---|
658 | c |
---|
659 | c Naive Boundary Conditions |
---|
660 | c Inflow, or reflective |
---|
661 | c |
---|
662 | c=========================================================== |
---|
663 | subroutine boundary_apply(x,h,uh,z) |
---|
664 | |
---|
665 | implicit none |
---|
666 | integer md,nxmax,nxd, mn |
---|
667 | parameter (md=3, nxmax=5000, nxd=nxmax+2*md, mn=2) |
---|
668 | real*8 h(nxd), uh(nxd), w(nxd) |
---|
669 | real*8 z(nxd),z2(nxd) |
---|
670 | real*8 x(nxd),x2(nxd) |
---|
671 | |
---|
672 | integer terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
673 | common /setup/ terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
674 | |
---|
675 | integer nlow,nhigh,nx,nxx,dn,version,isteps |
---|
676 | common /ivariables/ nlow,nhigh,nx,nxx,version,isteps |
---|
677 | |
---|
678 | real*8 toll,cfl,theta |
---|
679 | real*8 g,hf,dx,dt,dtmax |
---|
680 | real*8 Tout,tf,tc,uhf |
---|
681 | real*8 K,limith,length |
---|
682 | real*8 manning,hl,hr |
---|
683 | common /rvariables/ toll,cfl,theta, |
---|
684 | . g,hf,dx,dt,dtmax, |
---|
685 | . Tout,tf,tc,uhf, |
---|
686 | . K,limith,length, |
---|
687 | . manning,hl,hr |
---|
688 | |
---|
689 | real*8 du,dh,uf,h0,u0,uh0 |
---|
690 | integer ii,i |
---|
691 | |
---|
692 | c-------------------------------------------------- |
---|
693 | |
---|
694 | do i = 1, md |
---|
695 | c================== |
---|
696 | c Left Boundary |
---|
697 | c================== |
---|
698 | |
---|
699 | c --- Zero |
---|
700 | h(i) = 0.d0 |
---|
701 | uh(i) = 0.d0 |
---|
702 | end do |
---|
703 | |
---|
704 | do i = 1, md |
---|
705 | c================== |
---|
706 | c Right Boundary |
---|
707 | c================== |
---|
708 | |
---|
709 | c --- Zero |
---|
710 | h(nx+md+i) = 0.d0 |
---|
711 | uh(nx+md+i) = 0.d0 |
---|
712 | c h(nx+md+1) = h(nx-1) |
---|
713 | c uh(nx+md+1) = uh(nx-1) |
---|
714 | end do |
---|
715 | |
---|
716 | return |
---|
717 | end |
---|
718 | c=========================================================== |
---|
719 | c Numerical Central Flux |
---|
720 | c=========================================================== |
---|
721 | subroutine numerical_flux_central(h,uh,z2,z,uh_l,uh_r,h_l,h_r, |
---|
722 | . h_2,f_h,f_uh, em_x) |
---|
723 | |
---|
724 | implicit none |
---|
725 | integer md,nxmax,nxd, mn |
---|
726 | parameter (md=3, nxmax=5000, nxd=nxmax+2*md, mn=2) |
---|
727 | real*8 h(nxd), uh(nxd), w(nxd) |
---|
728 | real*8 z(nxd),z2(nxd) |
---|
729 | real*8 x(nxd),x2(nxd) |
---|
730 | |
---|
731 | c integer terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
732 | c common /setup/ terrain,bc_l,bc_r,initial,s_h,s_uh,s_semi_uh |
---|
733 | |
---|
734 | integer nlow,nhigh,nx,nxx,dn,version,isteps |
---|
735 | common /ivariables/ nlow,nhigh,nx,nxx,version,isteps |
---|
736 | |
---|
737 | real*8 toll,cfl,theta |
---|
738 | real*8 g,hf,dx,dt,dtmax |
---|
739 | real*8 Tout,tf,tc,uhf |
---|
740 | real*8 K,limith,length |
---|
741 | real*8 manning,hl,hr |
---|
742 | common /rvariables/ toll,cfl,theta, |
---|
743 | . g,hf,dx,dt,dtmax, |
---|
744 | . Tout,tf,tc,uhf, |
---|
745 | . K,limith,length, |
---|
746 | . manning,hl,hr |
---|
747 | |
---|
748 | real*8 f_h(nxd), f_uh(nxd) |
---|
749 | real*8 h_2(nxd), uh_2(nxd) |
---|
750 | real*8 h_l(nxd), h_r(nxd), uh_l(nxd), uh_r(nxd) |
---|
751 | |
---|
752 | |
---|
753 | real*8 fh_l, fuh_l, fh_r, fuh_r |
---|
754 | real*8 l1_l, l2_l, l1_r, l2_r |
---|
755 | real*8 a_max, a_min, aa, axa |
---|
756 | |
---|
757 | integer istop,iout,nt,ii,oi,i2,m,i |
---|
758 | real*8 a,b,t,Tleft,em_x,dtcdx,dtcdx2,dtcdx3 |
---|
759 | c real*8 den, xmt,pre, Tnext, em_old, hh |
---|
760 | c real*8 aterm1 ,term1, aterm2 ,term2, aterm3 ,term3 |
---|
761 | |
---|
762 | c----------------------------------- |
---|
763 | |
---|
764 | em_x = 1.0d-15 |
---|
765 | |
---|
766 | do i = md,nxx-md+1 |
---|
767 | |
---|
768 | c |
---|
769 | c ---- Calculate flux |
---|
770 | c |
---|
771 | |
---|
772 | call flux_huh(fh_l,fuh_l,h_r(i),uh_r(i), |
---|
773 | . em_x,l1_l,l2_l,g,toll) |
---|
774 | call flux_huh(fh_r,fuh_r,h_l(i+1),uh_l(i+1), |
---|
775 | . em_x,l1_r,l2_r,g,toll) |
---|
776 | c |
---|
777 | c --- Calculate max wave speed |
---|
778 | c |
---|
779 | a_max = dmax1(l2_l,l2_r,0.0d0) |
---|
780 | a_min = dmin1(l1_l,l1_r,0.0d0) |
---|
781 | c |
---|
782 | c --- Compute the numerical flux |
---|
783 | c |
---|
784 | aa = a_max - a_min |
---|
785 | axa = a_max*a_min |
---|
786 | if ( aa .gt. dx*dt/1000.0d0) then |
---|
787 | h_2(i) =(a_max*h_l(i+1) - a_min*h_r(i) - |
---|
788 | . fh_r + fh_l )/aa |
---|
789 | f_h(i) =(a_max*fh_l - a_min*fh_r + |
---|
790 | . axa*(h_l(i+1)-h_r(i)))/aa |
---|
791 | f_uh(i) =(a_max*fuh_l - a_min*fuh_r + |
---|
792 | . axa*(uh_l(i+1)-uh_r(i)))/aa |
---|
793 | else |
---|
794 | h_2(i) = 0.0d0 |
---|
795 | f_h(i) = 0.0d0 |
---|
796 | f_uh(i) = 0.0d0 |
---|
797 | endif |
---|
798 | |
---|
799 | enddo |
---|
800 | |
---|
801 | return |
---|
802 | end |
---|
803 | c=========================================================== |
---|
804 | c Open Files |
---|
805 | c=========================================================== |
---|
806 | subroutine util_open_files |
---|
807 | |
---|
808 | open(11,file='data.w') |
---|
809 | open(12,file='data.t') |
---|
810 | open(13,file='data.h') |
---|
811 | open(14,file='data.uh') |
---|
812 | open(15,file='data.n') |
---|
813 | open(16,file='data.x') |
---|
814 | |
---|
815 | |
---|
816 | open(25,file='data.xx') |
---|
817 | |
---|
818 | return |
---|
819 | end |
---|
820 | c=========================================================== |
---|
821 | c Close Files |
---|
822 | c=========================================================== |
---|
823 | subroutine util_close_files |
---|
824 | |
---|
825 | close(11) |
---|
826 | close(12) |
---|
827 | close(13) |
---|
828 | close(14) |
---|
829 | close(15) |
---|
830 | close(16) |
---|
831 | |
---|
832 | close(25) |
---|
833 | |
---|
834 | return |
---|
835 | end |
---|