1 | // Python - C extension module for shallow_water.py |
---|
2 | // |
---|
3 | // To compile (Python2.3): |
---|
4 | // gcc -c domain_ext.c -I/usr/include/python2.3 -o domain_ext.o -Wall -O |
---|
5 | // gcc -shared domain_ext.o -o domain_ext.so |
---|
6 | // |
---|
7 | // or use python compile.py |
---|
8 | // |
---|
9 | // See the module shallow_water.py |
---|
10 | // |
---|
11 | // |
---|
12 | // Ole Nielsen, GA 2004 |
---|
13 | |
---|
14 | |
---|
15 | #include "Python.h" |
---|
16 | #include "Numeric/arrayobject.h" |
---|
17 | #include "math.h" |
---|
18 | |
---|
19 | //Shared code snippets |
---|
20 | #include "util_ext.h" |
---|
21 | |
---|
22 | const double pi = 3.14159265358979; |
---|
23 | |
---|
24 | // Computational function for rotation |
---|
25 | int _rotate(double *q, double n1, double n2) { |
---|
26 | /*Rotate the momentum component q (q[1], q[2]) |
---|
27 | from x,y coordinates to coordinates based on normal vector (n1, n2). |
---|
28 | |
---|
29 | Result is returned in array 3x1 r |
---|
30 | To rotate in opposite direction, call rotate with (q, n1, -n2) |
---|
31 | |
---|
32 | Contents of q are changed by this function */ |
---|
33 | |
---|
34 | |
---|
35 | double q1, q2; |
---|
36 | |
---|
37 | //Shorthands |
---|
38 | q1 = q[1]; //uh momentum |
---|
39 | q2 = q[2]; //vh momentum |
---|
40 | |
---|
41 | //Rotate |
---|
42 | q[1] = n1*q1 + n2*q2; |
---|
43 | q[2] = -n2*q1 + n1*q2; |
---|
44 | |
---|
45 | return 0; |
---|
46 | } |
---|
47 | |
---|
48 | |
---|
49 | |
---|
50 | // Computational function for flux computation (using stage w=z+h) |
---|
51 | int flux_function(double *q_left, double *q_right, |
---|
52 | double z_left, double z_right, |
---|
53 | double n1, double n2, |
---|
54 | double epsilon, double g, |
---|
55 | double *edgeflux, double *max_speed) { |
---|
56 | |
---|
57 | /*Compute fluxes between volumes for the shallow water wave equation |
---|
58 | cast in terms of the 'stage', w = h+z using |
---|
59 | the 'central scheme' as described in |
---|
60 | |
---|
61 | Kurganov, Noelle, Petrova. 'Semidiscrete Central-Upwind Schemes For |
---|
62 | Hyperbolic Conservation Laws and Hamilton-Jacobi Equations'. |
---|
63 | Siam J. Sci. Comput. Vol. 23, No. 3, pp. 707-740. |
---|
64 | |
---|
65 | The implemented formula is given in equation (3.15) on page 714 |
---|
66 | */ |
---|
67 | |
---|
68 | int i; |
---|
69 | |
---|
70 | double w_left, h_left, uh_left, vh_left, u_left; |
---|
71 | double w_right, h_right, uh_right, vh_right, u_right; |
---|
72 | double s_min, s_max, soundspeed_left, soundspeed_right; |
---|
73 | double denom, z; |
---|
74 | double q_left_copy[3], q_right_copy[3]; |
---|
75 | double flux_right[3], flux_left[3]; |
---|
76 | |
---|
77 | //Copy conserved quantities to protect from modification |
---|
78 | for (i=0; i<3; i++) { |
---|
79 | q_left_copy[i] = q_left[i]; |
---|
80 | q_right_copy[i] = q_right[i]; |
---|
81 | } |
---|
82 | |
---|
83 | //Align x- and y-momentum with x-axis |
---|
84 | _rotate(q_left_copy, n1, n2); |
---|
85 | _rotate(q_right_copy, n1, n2); |
---|
86 | |
---|
87 | z = (z_left+z_right)/2; //Take average of field values |
---|
88 | |
---|
89 | //Compute speeds in x-direction |
---|
90 | w_left = q_left_copy[0]; // h+z |
---|
91 | h_left = w_left-z; |
---|
92 | uh_left = q_left_copy[1]; |
---|
93 | |
---|
94 | if (h_left < epsilon) { |
---|
95 | h_left = 0.0; //Could have been negative |
---|
96 | u_left = 0.0; |
---|
97 | } else { |
---|
98 | u_left = uh_left/h_left; |
---|
99 | } |
---|
100 | |
---|
101 | w_right = q_right_copy[0]; |
---|
102 | h_right = w_right-z; |
---|
103 | uh_right = q_right_copy[1]; |
---|
104 | |
---|
105 | if (h_right < epsilon) { |
---|
106 | h_right = 0.0; //Could have been negative |
---|
107 | u_right = 0.0; |
---|
108 | } else { |
---|
109 | u_right = uh_right/h_right; |
---|
110 | } |
---|
111 | |
---|
112 | //Momentum in y-direction |
---|
113 | vh_left = q_left_copy[2]; |
---|
114 | vh_right = q_right_copy[2]; |
---|
115 | |
---|
116 | |
---|
117 | //Maximal and minimal wave speeds |
---|
118 | soundspeed_left = sqrt(g*h_left); |
---|
119 | soundspeed_right = sqrt(g*h_right); |
---|
120 | |
---|
121 | s_max = max(u_left+soundspeed_left, u_right+soundspeed_right); |
---|
122 | if (s_max < 0.0) s_max = 0.0; |
---|
123 | |
---|
124 | s_min = min(u_left-soundspeed_left, u_right-soundspeed_right); |
---|
125 | if (s_min > 0.0) s_min = 0.0; |
---|
126 | |
---|
127 | //Flux formulas |
---|
128 | flux_left[0] = u_left*h_left; |
---|
129 | flux_left[1] = u_left*uh_left + 0.5*g*h_left*h_left; |
---|
130 | flux_left[2] = u_left*vh_left; |
---|
131 | |
---|
132 | flux_right[0] = u_right*h_right; |
---|
133 | flux_right[1] = u_right*uh_right + 0.5*g*h_right*h_right; |
---|
134 | flux_right[2] = u_right*vh_right; |
---|
135 | |
---|
136 | |
---|
137 | //Flux computation |
---|
138 | denom = s_max-s_min; |
---|
139 | if (denom == 0.0) { |
---|
140 | for (i=0; i<3; i++) edgeflux[i] = 0.0; |
---|
141 | *max_speed = 0.0; |
---|
142 | } else { |
---|
143 | for (i=0; i<3; i++) { |
---|
144 | edgeflux[i] = s_max*flux_left[i] - s_min*flux_right[i]; |
---|
145 | edgeflux[i] += s_max*s_min*(q_right_copy[i]-q_left_copy[i]); |
---|
146 | edgeflux[i] /= denom; |
---|
147 | } |
---|
148 | |
---|
149 | //Maximal wavespeed |
---|
150 | *max_speed = max(fabs(s_max), fabs(s_min)); |
---|
151 | |
---|
152 | //Rotate back |
---|
153 | _rotate(edgeflux, n1, -n2); |
---|
154 | } |
---|
155 | return 0; |
---|
156 | } |
---|
157 | |
---|
158 | void _manning_friction(double g, double eps, int N, |
---|
159 | double* w, double* z, |
---|
160 | double* uh, double* vh, |
---|
161 | double* eta, double* xmom, double* ymom) { |
---|
162 | |
---|
163 | int k; |
---|
164 | double S, h; |
---|
165 | |
---|
166 | for (k=0; k<N; k++) { |
---|
167 | if (eta[k] > eps) { |
---|
168 | h = w[k]-z[k]; |
---|
169 | if (h >= eps) { |
---|
170 | S = -g * eta[k]*eta[k] * sqrt((uh[k]*uh[k] + vh[k]*vh[k])); |
---|
171 | S /= pow(h, 7.0/3); //Expensive (on Ole's home computer) |
---|
172 | //S /= h*h*(1 + h/3.0 - h*h/9.0); //FIXME: Could use a Taylor expansion |
---|
173 | |
---|
174 | |
---|
175 | //Update momentum |
---|
176 | xmom[k] += S*uh[k]; |
---|
177 | ymom[k] += S*vh[k]; |
---|
178 | } |
---|
179 | } |
---|
180 | } |
---|
181 | } |
---|
182 | |
---|
183 | |
---|
184 | |
---|
185 | int _balance_deep_and_shallow(int N, |
---|
186 | double* wc, |
---|
187 | double* zc, |
---|
188 | double* hc, |
---|
189 | double* wv, |
---|
190 | double* zv, |
---|
191 | double* hv, |
---|
192 | double* xmomc, |
---|
193 | double* ymomc, |
---|
194 | double* xmomv, |
---|
195 | double* ymomv) { |
---|
196 | |
---|
197 | int k, k3, i; |
---|
198 | double dz, hmin, alpha; |
---|
199 | |
---|
200 | //Compute linear combination between constant levels and and |
---|
201 | //levels parallel to the bed elevation. |
---|
202 | |
---|
203 | for (k=0; k<N; k++) { |
---|
204 | // Compute maximal variation in bed elevation |
---|
205 | // This quantitiy is |
---|
206 | // dz = max_i abs(z_i - z_c) |
---|
207 | // and it is independent of dimension |
---|
208 | // In the 1d case zc = (z0+z1)/2 |
---|
209 | // In the 2d case zc = (z0+z1+z2)/3 |
---|
210 | |
---|
211 | k3 = 3*k; |
---|
212 | |
---|
213 | //FIXME: Try with this one precomputed |
---|
214 | dz = 0.0; |
---|
215 | hmin = hv[k3]; |
---|
216 | for (i=0; i<3; i++) { |
---|
217 | dz = max(dz, fabs(zv[k3+i]-zc[k])); |
---|
218 | hmin = min(hmin, hv[k3+i]); |
---|
219 | } |
---|
220 | |
---|
221 | |
---|
222 | //Create alpha in [0,1], where alpha==0 means using shallow |
---|
223 | //first order scheme and alpha==1 means using the stage w as |
---|
224 | //computed by the gradient limiter (1st or 2nd order) |
---|
225 | // |
---|
226 | //If hmin > dz/2 then alpha = 1 and the bed will have no effect |
---|
227 | //If hmin < 0 then alpha = 0 reverting to constant height above bed. |
---|
228 | |
---|
229 | if (dz > 0.0) |
---|
230 | alpha = max( min( 2*hmin/dz, 1.0), 0.0 ); |
---|
231 | else |
---|
232 | alpha = 1.0; //Flat bed |
---|
233 | |
---|
234 | |
---|
235 | //Weighted balance between stage parallel to bed elevation |
---|
236 | //(wvi = zvi + hc) and stage as computed by 1st or 2nd |
---|
237 | //order gradient limiter |
---|
238 | //(wvi = zvi + hvi) where i=0,1,2 denotes the vertex ids |
---|
239 | // |
---|
240 | //It follows that the updated wvi is |
---|
241 | // wvi := (1-alpha)*(zvi+hc) + alpha*(zvi+hvi) = |
---|
242 | // zvi + hc + alpha*(hvi - hc) |
---|
243 | // |
---|
244 | //Note that hvi = zc+hc-zvi in the first order case (constant). |
---|
245 | |
---|
246 | if (alpha < 1) { |
---|
247 | for (i=0; i<3; i++) { |
---|
248 | wv[k3+i] = zv[k3+i] + hc[k] + alpha*(hv[k3+i]-hc[k]); |
---|
249 | |
---|
250 | |
---|
251 | //Update momentum as a linear combination of |
---|
252 | //xmomc and ymomc (shallow) and momentum |
---|
253 | //from extrapolator xmomv and ymomv (deep). |
---|
254 | xmomv[k3+i] = (1-alpha)*xmomc[k] + alpha*xmomv[k3+i]; |
---|
255 | ymomv[k3+i] = (1-alpha)*ymomc[k] + alpha*ymomv[k3+i]; |
---|
256 | } |
---|
257 | } |
---|
258 | } |
---|
259 | return 0; |
---|
260 | } |
---|
261 | |
---|
262 | |
---|
263 | |
---|
264 | int _protect(int N, |
---|
265 | double minimum_allowed_height, |
---|
266 | double* wc, |
---|
267 | double* zc, |
---|
268 | double* xmomc, |
---|
269 | double* ymomc) { |
---|
270 | |
---|
271 | int k; |
---|
272 | double hc; |
---|
273 | |
---|
274 | //Protect against initesimal and negative heights |
---|
275 | |
---|
276 | for (k=0; k<N; k++) { |
---|
277 | hc = wc[k] - zc[k]; |
---|
278 | |
---|
279 | if (hc < minimum_allowed_height) { |
---|
280 | wc[k] = zc[k]; |
---|
281 | xmomc[k] = 0.0; |
---|
282 | ymomc[k] = 0.0; |
---|
283 | } |
---|
284 | |
---|
285 | } |
---|
286 | return 0; |
---|
287 | } |
---|
288 | |
---|
289 | |
---|
290 | |
---|
291 | int _assign_wind_field_values(int N, |
---|
292 | double* xmom_update, |
---|
293 | double* ymom_update, |
---|
294 | double* s_vec, |
---|
295 | double* phi_vec, |
---|
296 | double cw) { |
---|
297 | |
---|
298 | //Assign windfield values to momentum updates |
---|
299 | |
---|
300 | int k; |
---|
301 | double S, s, phi, u, v; |
---|
302 | |
---|
303 | for (k=0; k<N; k++) { |
---|
304 | |
---|
305 | s = s_vec[k]; |
---|
306 | phi = phi_vec[k]; |
---|
307 | |
---|
308 | //Convert to radians |
---|
309 | phi = phi*pi/180; |
---|
310 | |
---|
311 | //Compute velocity vector (u, v) |
---|
312 | u = s*cos(phi); |
---|
313 | v = s*sin(phi); |
---|
314 | |
---|
315 | //Compute wind stress |
---|
316 | S = cw * sqrt(u*u + v*v); |
---|
317 | xmom_update[k] += S*u; |
---|
318 | ymom_update[k] += S*v; |
---|
319 | } |
---|
320 | return 0; |
---|
321 | } |
---|
322 | |
---|
323 | |
---|
324 | |
---|
325 | /////////////////////////////////////////////////////////////////// |
---|
326 | // Gateways to Python |
---|
327 | |
---|
328 | PyObject *gravity(PyObject *self, PyObject *args) { |
---|
329 | // |
---|
330 | // gravity(g, h, v, x, xmom, ymom) |
---|
331 | // |
---|
332 | |
---|
333 | |
---|
334 | PyArrayObject *h, *v, *x, *xmom, *ymom; |
---|
335 | int k, i, N, k3, k6; |
---|
336 | double g, avg_h, zx, zy; |
---|
337 | double x0, y0, x1, y1, x2, y2, z0, z1, z2; |
---|
338 | |
---|
339 | if (!PyArg_ParseTuple(args, "dOOOOO", |
---|
340 | &g, &h, &v, &x, |
---|
341 | &xmom, &ymom)) |
---|
342 | return NULL; |
---|
343 | |
---|
344 | N = h -> dimensions[0]; |
---|
345 | for (k=0; k<N; k++) { |
---|
346 | k3 = 3*k; // base index |
---|
347 | k6 = 6*k; // base index |
---|
348 | |
---|
349 | avg_h = 0.0; |
---|
350 | for (i=0; i<3; i++) { |
---|
351 | avg_h += ((double *) h -> data)[k3+i]; |
---|
352 | } |
---|
353 | avg_h /= 3; |
---|
354 | |
---|
355 | |
---|
356 | //Compute bed slope |
---|
357 | x0 = ((double*) x -> data)[k6 + 0]; |
---|
358 | y0 = ((double*) x -> data)[k6 + 1]; |
---|
359 | x1 = ((double*) x -> data)[k6 + 2]; |
---|
360 | y1 = ((double*) x -> data)[k6 + 3]; |
---|
361 | x2 = ((double*) x -> data)[k6 + 4]; |
---|
362 | y2 = ((double*) x -> data)[k6 + 5]; |
---|
363 | |
---|
364 | |
---|
365 | z0 = ((double*) v -> data)[k3 + 0]; |
---|
366 | z1 = ((double*) v -> data)[k3 + 1]; |
---|
367 | z2 = ((double*) v -> data)[k3 + 2]; |
---|
368 | |
---|
369 | _gradient(x0, y0, x1, y1, x2, y2, z0, z1, z2, &zx, &zy); |
---|
370 | |
---|
371 | //Update momentum |
---|
372 | ((double*) xmom -> data)[k] += -g*zx*avg_h; |
---|
373 | ((double*) ymom -> data)[k] += -g*zy*avg_h; |
---|
374 | } |
---|
375 | |
---|
376 | return Py_BuildValue(""); |
---|
377 | } |
---|
378 | |
---|
379 | |
---|
380 | PyObject *manning_friction(PyObject *self, PyObject *args) { |
---|
381 | // |
---|
382 | // manning_friction(g, eps, h, uh, vh, eta, xmom_update, ymom_update) |
---|
383 | // |
---|
384 | |
---|
385 | |
---|
386 | PyArrayObject *w, *z, *uh, *vh, *eta, *xmom, *ymom; |
---|
387 | int N; |
---|
388 | double g, eps; |
---|
389 | |
---|
390 | if (!PyArg_ParseTuple(args, "ddOOOOOOO", |
---|
391 | &g, &eps, &w, &z, &uh, &vh, &eta, |
---|
392 | &xmom, &ymom)) |
---|
393 | return NULL; |
---|
394 | |
---|
395 | N = w -> dimensions[0]; |
---|
396 | _manning_friction(g, eps, N, |
---|
397 | (double*) w -> data, |
---|
398 | (double*) z -> data, |
---|
399 | (double*) uh -> data, |
---|
400 | (double*) vh -> data, |
---|
401 | (double*) eta -> data, |
---|
402 | (double*) xmom -> data, |
---|
403 | (double*) ymom -> data); |
---|
404 | |
---|
405 | return Py_BuildValue(""); |
---|
406 | } |
---|
407 | |
---|
408 | PyObject *rotate(PyObject *self, PyObject *args, PyObject *kwargs) { |
---|
409 | // |
---|
410 | // r = rotate(q, normal, direction=1) |
---|
411 | // |
---|
412 | // Where q is assumed to be a Float numeric array of length 3 and |
---|
413 | // normal a Float numeric array of length 2. |
---|
414 | |
---|
415 | |
---|
416 | PyObject *Q, *Normal; |
---|
417 | PyArrayObject *q, *r, *normal; |
---|
418 | |
---|
419 | static char *argnames[] = {"q", "normal", "direction", NULL}; |
---|
420 | int dimensions[1], i, direction=1; |
---|
421 | double n1, n2; |
---|
422 | |
---|
423 | // Convert Python arguments to C |
---|
424 | if (!PyArg_ParseTupleAndKeywords(args, kwargs, "OO|i", argnames, |
---|
425 | &Q, &Normal, &direction)) |
---|
426 | return NULL; |
---|
427 | |
---|
428 | //Input checks (convert sequences into numeric arrays) |
---|
429 | q = (PyArrayObject *) |
---|
430 | PyArray_ContiguousFromObject(Q, PyArray_DOUBLE, 0, 0); |
---|
431 | normal = (PyArrayObject *) |
---|
432 | PyArray_ContiguousFromObject(Normal, PyArray_DOUBLE, 0, 0); |
---|
433 | |
---|
434 | //Allocate space for return vector r (don't DECREF) |
---|
435 | dimensions[0] = 3; |
---|
436 | r = (PyArrayObject *) PyArray_FromDims(1, dimensions, PyArray_DOUBLE); |
---|
437 | |
---|
438 | //Copy |
---|
439 | for (i=0; i<3; i++) { |
---|
440 | ((double *) (r -> data))[i] = ((double *) (q -> data))[i]; |
---|
441 | } |
---|
442 | |
---|
443 | //Get normal and direction |
---|
444 | n1 = ((double *) normal -> data)[0]; |
---|
445 | n2 = ((double *) normal -> data)[1]; |
---|
446 | if (direction == -1) n2 = -n2; |
---|
447 | |
---|
448 | //Rotate |
---|
449 | _rotate((double *) r -> data, n1, n2); |
---|
450 | |
---|
451 | //Release numeric arrays |
---|
452 | Py_DECREF(q); |
---|
453 | Py_DECREF(normal); |
---|
454 | |
---|
455 | //return result using PyArray to avoid memory leak |
---|
456 | return PyArray_Return(r); |
---|
457 | } |
---|
458 | |
---|
459 | |
---|
460 | PyObject *compute_fluxes(PyObject *self, PyObject *args) { |
---|
461 | /*Compute all fluxes and the timestep suitable for all volumes |
---|
462 | in domain. |
---|
463 | |
---|
464 | Compute total flux for each conserved quantity using "flux_function" |
---|
465 | |
---|
466 | Fluxes across each edge are scaled by edgelengths and summed up |
---|
467 | Resulting flux is then scaled by area and stored in |
---|
468 | explicit_update for each of the three conserved quantities |
---|
469 | level, xmomentum and ymomentum |
---|
470 | |
---|
471 | The maximal allowable speed computed by the flux_function for each volume |
---|
472 | is converted to a timestep that must not be exceeded. The minimum of |
---|
473 | those is computed as the next overall timestep. |
---|
474 | |
---|
475 | Python call: |
---|
476 | domain.timestep = compute_fluxes(timestep, |
---|
477 | domain.epsilon, |
---|
478 | domain.g, |
---|
479 | domain.neighbours, |
---|
480 | domain.neighbour_edges, |
---|
481 | domain.normals, |
---|
482 | domain.edgelengths, |
---|
483 | domain.radii, |
---|
484 | domain.areas, |
---|
485 | Level.edge_values, |
---|
486 | Xmom.edge_values, |
---|
487 | Ymom.edge_values, |
---|
488 | Bed.edge_values, |
---|
489 | Level.boundary_values, |
---|
490 | Xmom.boundary_values, |
---|
491 | Ymom.boundary_values, |
---|
492 | Level.explicit_update, |
---|
493 | Xmom.explicit_update, |
---|
494 | Ymom.explicit_update) |
---|
495 | |
---|
496 | |
---|
497 | Post conditions: |
---|
498 | domain.explicit_update is reset to computed flux values |
---|
499 | domain.timestep is set to the largest step satisfying all volumes. |
---|
500 | |
---|
501 | |
---|
502 | */ |
---|
503 | |
---|
504 | |
---|
505 | PyArrayObject *neighbours, *neighbour_edges, |
---|
506 | *normals, *edgelengths, *radii, *areas, |
---|
507 | *level_edge_values, |
---|
508 | *xmom_edge_values, |
---|
509 | *ymom_edge_values, |
---|
510 | *bed_edge_values, |
---|
511 | *level_boundary_values, |
---|
512 | *xmom_boundary_values, |
---|
513 | *ymom_boundary_values, |
---|
514 | *level_explicit_update, |
---|
515 | *xmom_explicit_update, |
---|
516 | *ymom_explicit_update; |
---|
517 | |
---|
518 | |
---|
519 | //Local variables |
---|
520 | double timestep, max_speed, epsilon, g; |
---|
521 | double normal[2], ql[3], qr[3], zl, zr; |
---|
522 | double flux[3], edgeflux[3]; //Work arrays for summing up fluxes |
---|
523 | |
---|
524 | int number_of_elements, k, i, j, m, n; |
---|
525 | int ki, nm, ki2; //Index shorthands |
---|
526 | |
---|
527 | |
---|
528 | // Convert Python arguments to C |
---|
529 | if (!PyArg_ParseTuple(args, "dddOOOOOOOOOOOOOOOO", |
---|
530 | ×tep, |
---|
531 | &epsilon, |
---|
532 | &g, |
---|
533 | &neighbours, |
---|
534 | &neighbour_edges, |
---|
535 | &normals, |
---|
536 | &edgelengths, &radii, &areas, |
---|
537 | &level_edge_values, |
---|
538 | &xmom_edge_values, |
---|
539 | &ymom_edge_values, |
---|
540 | &bed_edge_values, |
---|
541 | &level_boundary_values, |
---|
542 | &xmom_boundary_values, |
---|
543 | &ymom_boundary_values, |
---|
544 | &level_explicit_update, |
---|
545 | &xmom_explicit_update, |
---|
546 | &ymom_explicit_update)) { |
---|
547 | PyErr_SetString(PyExc_RuntimeError, "Input arguments failed"); |
---|
548 | return NULL; |
---|
549 | } |
---|
550 | |
---|
551 | number_of_elements = level_edge_values -> dimensions[0]; |
---|
552 | |
---|
553 | |
---|
554 | for (k=0; k<number_of_elements; k++) { |
---|
555 | |
---|
556 | //Reset work array |
---|
557 | for (j=0; j<3; j++) flux[j] = 0.0; |
---|
558 | |
---|
559 | //Loop through neighbours and compute edge flux for each |
---|
560 | for (i=0; i<3; i++) { |
---|
561 | ki = k*3+i; |
---|
562 | ql[0] = ((double *) level_edge_values -> data)[ki]; |
---|
563 | ql[1] = ((double *) xmom_edge_values -> data)[ki]; |
---|
564 | ql[2] = ((double *) ymom_edge_values -> data)[ki]; |
---|
565 | zl = ((double *) bed_edge_values -> data)[ki]; |
---|
566 | |
---|
567 | //Quantities at neighbour on nearest face |
---|
568 | n = ((int *) neighbours -> data)[ki]; |
---|
569 | if (n < 0) { |
---|
570 | m = -n-1; //Convert negative flag to index |
---|
571 | qr[0] = ((double *) level_boundary_values -> data)[m]; |
---|
572 | qr[1] = ((double *) xmom_boundary_values -> data)[m]; |
---|
573 | qr[2] = ((double *) ymom_boundary_values -> data)[m]; |
---|
574 | zr = zl; //Extend bed elevation to boundary |
---|
575 | } else { |
---|
576 | m = ((int *) neighbour_edges -> data)[ki]; |
---|
577 | |
---|
578 | nm = n*3+m; |
---|
579 | qr[0] = ((double *) level_edge_values -> data)[nm]; |
---|
580 | qr[1] = ((double *) xmom_edge_values -> data)[nm]; |
---|
581 | qr[2] = ((double *) ymom_edge_values -> data)[nm]; |
---|
582 | zr = ((double *) bed_edge_values -> data)[nm]; |
---|
583 | } |
---|
584 | |
---|
585 | // Outward pointing normal vector |
---|
586 | // normal = domain.normals[k, 2*i:2*i+2] |
---|
587 | ki2 = 2*ki; //k*6 + i*2 |
---|
588 | normal[0] = ((double *) normals -> data)[ki2]; |
---|
589 | normal[1] = ((double *) normals -> data)[ki2+1]; |
---|
590 | |
---|
591 | //Edge flux computation |
---|
592 | flux_function(ql, qr, zl, zr, |
---|
593 | normal[0], normal[1], |
---|
594 | epsilon, g, |
---|
595 | edgeflux, &max_speed); |
---|
596 | |
---|
597 | |
---|
598 | //flux -= edgeflux * edgelengths[k,i] |
---|
599 | for (j=0; j<3; j++) { |
---|
600 | flux[j] -= edgeflux[j]*((double *) edgelengths -> data)[ki]; |
---|
601 | } |
---|
602 | |
---|
603 | //Update timestep |
---|
604 | //timestep = min(timestep, domain.radii[k]/max_speed) |
---|
605 | if (max_speed > epsilon) { |
---|
606 | timestep = min(timestep, ((double *) radii -> data)[k]/max_speed); |
---|
607 | } |
---|
608 | } // end for i |
---|
609 | |
---|
610 | //Normalise by area and store for when all conserved |
---|
611 | //quantities get updated |
---|
612 | // flux /= areas[k] |
---|
613 | for (j=0; j<3; j++) { |
---|
614 | flux[j] /= ((double *) areas -> data)[k]; |
---|
615 | } |
---|
616 | |
---|
617 | ((double *) level_explicit_update -> data)[k] = flux[0]; |
---|
618 | ((double *) xmom_explicit_update -> data)[k] = flux[1]; |
---|
619 | ((double *) ymom_explicit_update -> data)[k] = flux[2]; |
---|
620 | |
---|
621 | } //end for k |
---|
622 | |
---|
623 | return Py_BuildValue("d", timestep); |
---|
624 | } |
---|
625 | |
---|
626 | |
---|
627 | |
---|
628 | PyObject *protect(PyObject *self, PyObject *args) { |
---|
629 | // |
---|
630 | // protect(minimum_allowed_height, wc, zc, xmomc, ymomc) |
---|
631 | |
---|
632 | |
---|
633 | PyArrayObject |
---|
634 | *wc, //Level at centroids |
---|
635 | *zc, //Elevation at centroids |
---|
636 | *xmomc, //Momentums at centroids |
---|
637 | *ymomc; |
---|
638 | |
---|
639 | |
---|
640 | int N; |
---|
641 | double minimum_allowed_height; |
---|
642 | |
---|
643 | // Convert Python arguments to C |
---|
644 | if (!PyArg_ParseTuple(args, "dOOOO", |
---|
645 | &minimum_allowed_height, |
---|
646 | &wc, &zc, &xmomc, &ymomc)) |
---|
647 | return NULL; |
---|
648 | |
---|
649 | N = wc -> dimensions[0]; |
---|
650 | |
---|
651 | _protect(N, |
---|
652 | minimum_allowed_height, |
---|
653 | (double*) wc -> data, |
---|
654 | (double*) zc -> data, |
---|
655 | (double*) xmomc -> data, |
---|
656 | (double*) ymomc -> data); |
---|
657 | |
---|
658 | return Py_BuildValue(""); |
---|
659 | } |
---|
660 | |
---|
661 | |
---|
662 | |
---|
663 | PyObject *balance_deep_and_shallow(PyObject *self, PyObject *args) { |
---|
664 | // |
---|
665 | // balance_deep_and_shallow(wc, zc, hc, wv, zv, hv, |
---|
666 | // xmomc, ymomc, xmomv, ymomv) |
---|
667 | |
---|
668 | |
---|
669 | PyArrayObject |
---|
670 | *wc, //Level at centroids |
---|
671 | *zc, //Elevation at centroids |
---|
672 | *hc, //Height at centroids |
---|
673 | *wv, //Level at vertices |
---|
674 | *zv, //Elevation at vertices |
---|
675 | *hv, //Heights at vertices |
---|
676 | *xmomc, //Momentums at centroids and vertices |
---|
677 | *ymomc, |
---|
678 | *xmomv, |
---|
679 | *ymomv; |
---|
680 | |
---|
681 | int N; //, err; |
---|
682 | |
---|
683 | // Convert Python arguments to C |
---|
684 | if (!PyArg_ParseTuple(args, "OOOOOOOOOO", |
---|
685 | &wc, &zc, &hc, |
---|
686 | &wv, &zv, &hv, |
---|
687 | &xmomc, &ymomc, &xmomv, &ymomv)) |
---|
688 | return NULL; |
---|
689 | |
---|
690 | N = wc -> dimensions[0]; |
---|
691 | |
---|
692 | _balance_deep_and_shallow(N, |
---|
693 | (double*) wc -> data, |
---|
694 | (double*) zc -> data, |
---|
695 | (double*) hc -> data, |
---|
696 | (double*) wv -> data, |
---|
697 | (double*) zv -> data, |
---|
698 | (double*) hv -> data, |
---|
699 | (double*) xmomc -> data, |
---|
700 | (double*) ymomc -> data, |
---|
701 | (double*) xmomv -> data, |
---|
702 | (double*) ymomv -> data); |
---|
703 | |
---|
704 | |
---|
705 | return Py_BuildValue(""); |
---|
706 | } |
---|
707 | |
---|
708 | |
---|
709 | |
---|
710 | |
---|
711 | |
---|
712 | |
---|
713 | |
---|
714 | PyObject *assign_windfield_values(PyObject *self, PyObject *args) { |
---|
715 | // |
---|
716 | // assign_windfield_values(xmom_update, ymom_update, |
---|
717 | // s_vec, phi_vec, self.const) |
---|
718 | |
---|
719 | |
---|
720 | |
---|
721 | PyArrayObject //(one element per triangle) |
---|
722 | *s_vec, //Speeds |
---|
723 | *phi_vec, //Bearings |
---|
724 | *xmom_update, //Momentum updates |
---|
725 | *ymom_update; |
---|
726 | |
---|
727 | |
---|
728 | int N; |
---|
729 | double cw; |
---|
730 | |
---|
731 | // Convert Python arguments to C |
---|
732 | if (!PyArg_ParseTuple(args, "OOOOd", |
---|
733 | &xmom_update, |
---|
734 | &ymom_update, |
---|
735 | &s_vec, &phi_vec, |
---|
736 | &cw)) |
---|
737 | return NULL; |
---|
738 | |
---|
739 | N = xmom_update -> dimensions[0]; |
---|
740 | |
---|
741 | _assign_wind_field_values(N, |
---|
742 | (double*) xmom_update -> data, |
---|
743 | (double*) ymom_update -> data, |
---|
744 | (double*) s_vec -> data, |
---|
745 | (double*) phi_vec -> data, |
---|
746 | cw); |
---|
747 | |
---|
748 | return Py_BuildValue(""); |
---|
749 | } |
---|
750 | |
---|
751 | |
---|
752 | |
---|
753 | |
---|
754 | ////////////////////////////////////////// |
---|
755 | // Method table for python module |
---|
756 | static struct PyMethodDef MethodTable[] = { |
---|
757 | /* The cast of the function is necessary since PyCFunction values |
---|
758 | * only take two PyObject* parameters, and rotate() takes |
---|
759 | * three. |
---|
760 | */ |
---|
761 | |
---|
762 | {"rotate", (PyCFunction)rotate, METH_VARARGS | METH_KEYWORDS, "Print out"}, |
---|
763 | {"compute_fluxes", compute_fluxes, METH_VARARGS, "Print out"}, |
---|
764 | {"gravity", gravity, METH_VARARGS, "Print out"}, |
---|
765 | {"manning_friction", manning_friction, METH_VARARGS, "Print out"}, |
---|
766 | {"balance_deep_and_shallow", balance_deep_and_shallow, |
---|
767 | METH_VARARGS, "Print out"}, |
---|
768 | {"protect", protect, METH_VARARGS | METH_KEYWORDS, "Print out"}, |
---|
769 | {"assign_windfield_values", assign_windfield_values, |
---|
770 | METH_VARARGS | METH_KEYWORDS, "Print out"}, |
---|
771 | //{"distribute_to_vertices_and_edges", |
---|
772 | // distribute_to_vertices_and_edges, METH_VARARGS}, |
---|
773 | //{"update_conserved_quantities", |
---|
774 | // update_conserved_quantities, METH_VARARGS}, |
---|
775 | //{"set_initialcondition", |
---|
776 | // set_initialcondition, METH_VARARGS}, |
---|
777 | {NULL, NULL} |
---|
778 | }; |
---|
779 | |
---|
780 | // Module initialisation |
---|
781 | void initshallow_water_ext(void){ |
---|
782 | Py_InitModule("shallow_water_ext", MethodTable); |
---|
783 | |
---|
784 | import_array(); //Necessary for handling of NumPY structures |
---|
785 | } |
---|