source: branches/numpy/anuga/shallow_water/shallow_water_ext.c @ 6780

Last change on this file since 6780 was 6780, checked in by rwilson, 16 years ago

Changes to make legacy Numeric API work in numpy.

File size: 74.1 KB
Line 
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_domain.py for more documentation on
10// how to use this module
11//
12//
13// Ole Nielsen, GA 2004
14
15
16#include "Python.h"
17#include "numpy/arrayobject.h"
18#include "math.h"
19#include <stdio.h>
20
21#include "numpy_shim.h"
22
23// Shared code snippets
24#include "util_ext.h"
25
26
27const double pi = 3.14159265358979;
28
29// Computational function for rotation
30int _rotate(double *q, double n1, double n2) {
31  /*Rotate the momentum component q (q[1], q[2])
32    from x,y coordinates to coordinates based on normal vector (n1, n2).
33
34    Result is returned in array 3x1 r
35    To rotate in opposite direction, call rotate with (q, n1, -n2)
36
37    Contents of q are changed by this function */
38
39
40  double q1, q2;
41
42  // Shorthands
43  q1 = q[1];  // uh momentum
44  q2 = q[2];  // vh momentum
45
46  // Rotate
47  q[1] =  n1*q1 + n2*q2;
48  q[2] = -n2*q1 + n1*q2;
49
50  return 0;
51}
52
53int find_qmin_and_qmax(double dq0, double dq1, double dq2, 
54                       double *qmin, double *qmax){
55  // Considering the centroid of an FV triangle and the vertices of its
56  // auxiliary triangle, find
57  // qmin=min(q)-qc and qmax=max(q)-qc,
58  // where min(q) and max(q) are respectively min and max over the
59  // four values (at the centroid of the FV triangle and the auxiliary
60  // triangle vertices),
61  // and qc is the centroid
62  // dq0=q(vertex0)-q(centroid of FV triangle)
63  // dq1=q(vertex1)-q(vertex0)
64  // dq2=q(vertex2)-q(vertex0)
65 
66  if (dq0>=0.0){
67    if (dq1>=dq2){
68      if (dq1>=0.0)
69        *qmax=dq0+dq1;
70      else
71        *qmax=dq0;
72       
73      *qmin=dq0+dq2;
74      if (*qmin>=0.0) *qmin = 0.0;
75    }
76    else{// dq1<dq2
77      if (dq2>0)
78        *qmax=dq0+dq2;
79      else
80        *qmax=dq0;
81       
82      *qmin=dq0+dq1;   
83      if (*qmin>=0.0) *qmin=0.0;
84    }
85  }
86  else{//dq0<0
87    if (dq1<=dq2){
88      if (dq1<0.0)
89        *qmin=dq0+dq1;
90      else
91        *qmin=dq0;
92       
93      *qmax=dq0+dq2;   
94      if (*qmax<=0.0) *qmax=0.0;
95    }
96    else{// dq1>dq2
97      if (dq2<0.0)
98        *qmin=dq0+dq2;
99      else
100        *qmin=dq0;
101       
102      *qmax=dq0+dq1;
103      if (*qmax<=0.0) *qmax=0.0;
104    }
105  }
106  return 0;
107}
108
109int limit_gradient(double *dqv, double qmin, double qmax, double beta_w){
110  // Given provisional jumps dqv from the FV triangle centroid to its
111  // vertices and jumps qmin (qmax) between the centroid of the FV
112  // triangle and the minimum (maximum) of the values at the centroid of
113  // the FV triangle and the auxiliary triangle vertices,
114  // calculate a multiplicative factor phi by which the provisional
115  // vertex jumps are to be limited
116 
117  int i;
118  double r=1000.0, r0=1.0, phi=1.0;
119  static double TINY = 1.0e-100; // to avoid machine accuracy problems.
120  // FIXME: Perhaps use the epsilon used elsewhere.
121 
122  // Any provisional jump with magnitude < TINY does not contribute to
123  // the limiting process.
124  for (i=0;i<3;i++){
125    if (dqv[i]<-TINY)
126      r0=qmin/dqv[i];
127     
128    if (dqv[i]>TINY)
129      r0=qmax/dqv[i];
130     
131    r=min(r0,r);
132  }
133 
134  phi=min(r*beta_w,1.0);
135  for (i=0;i<3;i++)
136    dqv[i]=dqv[i]*phi;
137   
138  return 0;
139}
140
141
142double compute_froude_number(double uh,
143                             double h, 
144                             double g,
145                             double epsilon) {
146                         
147  // Compute Froude number; v/sqrt(gh)
148  // FIXME (Ole): Not currently in use
149   
150  double froude_number;
151 
152  //Compute Froude number (stability diagnostics)
153  if (h > epsilon) { 
154    froude_number = uh/sqrt(g*h)/h;
155  } else {
156    froude_number = 0.0;
157    // FIXME (Ole): What should it be when dry??
158  }
159 
160  return froude_number;
161}
162
163
164
165// Function to obtain speed from momentum and depth.
166// This is used by flux functions
167// Input parameters uh and h may be modified by this function.
168double _compute_speed(double *uh, 
169                      double *h, 
170                      double epsilon, 
171                      double h0) {
172 
173  double u;
174
175 
176  if (*h < epsilon) {
177    *h = 0.0;  //Could have been negative
178    u = 0.0;
179  } else {
180    u = *uh/(*h + h0/ *h);   
181  }
182 
183
184  // Adjust momentum to be consistent with speed
185  *uh = u * *h;
186 
187  return u;
188}
189
190// Innermost flux function (using stage w=z+h)
191int _flux_function_central(double *q_left, double *q_right,
192                           double z_left, double z_right,
193                           double n1, double n2,
194                           double epsilon, double H0, double g,
195                           double *edgeflux, double *max_speed) {
196
197  /*Compute fluxes between volumes for the shallow water wave equation
198    cast in terms of the 'stage', w = h+z using
199    the 'central scheme' as described in
200
201    Kurganov, Noelle, Petrova. 'Semidiscrete Central-Upwind Schemes For
202    Hyperbolic Conservation Laws and Hamilton-Jacobi Equations'.
203    Siam J. Sci. Comput. Vol. 23, No. 3, pp. 707-740.
204
205    The implemented formula is given in equation (3.15) on page 714
206  */
207
208  int i;
209
210  double w_left, h_left, uh_left, vh_left, u_left;
211  double w_right, h_right, uh_right, vh_right, u_right;
212  double v_left, v_right; 
213  double s_min, s_max, soundspeed_left, soundspeed_right;
214  double denom, z;
215
216  // Workspace (allocate once, use many)
217  static double q_left_rotated[3], q_right_rotated[3], flux_right[3], flux_left[3];
218
219  double h0 = H0*H0; // This ensures a good balance when h approaches H0.
220                     // But evidence suggests that h0 can be as little as
221                     // epsilon!
222 
223  // Copy conserved quantities to protect from modification
224  for (i=0; i<3; i++) {
225    q_left_rotated[i] = q_left[i];
226    q_right_rotated[i] = q_right[i];
227  }
228
229  // Align x- and y-momentum with x-axis
230  _rotate(q_left_rotated, n1, n2);
231  _rotate(q_right_rotated, n1, n2);
232
233  z = (z_left+z_right)/2; // Average elevation values.
234                          // Even though this will nominally allow for discontinuities
235                          // in the elevation data, there is currently no numerical
236                          // support for this so results may be strange near jumps in the bed.
237
238  // Compute speeds in x-direction
239  w_left = q_left_rotated[0];         
240  h_left = w_left-z;
241  uh_left = q_left_rotated[1];
242  u_left = _compute_speed(&uh_left, &h_left, epsilon, h0);
243
244  w_right = q_right_rotated[0];
245  h_right = w_right-z;
246  uh_right = q_right_rotated[1];
247  u_right = _compute_speed(&uh_right, &h_right, epsilon, h0); 
248
249  // Momentum in y-direction
250  vh_left  = q_left_rotated[2];
251  vh_right = q_right_rotated[2];
252
253  // Limit y-momentum if necessary 
254  v_left = _compute_speed(&vh_left, &h_left, epsilon, h0);
255  v_right = _compute_speed(&vh_right, &h_right, epsilon, h0);
256
257  // Maximal and minimal wave speeds
258  soundspeed_left  = sqrt(g*h_left);
259  soundspeed_right = sqrt(g*h_right);
260
261  s_max = max(u_left+soundspeed_left, u_right+soundspeed_right);
262  if (s_max < 0.0) s_max = 0.0;
263
264  s_min = min(u_left-soundspeed_left, u_right-soundspeed_right);
265  if (s_min > 0.0) s_min = 0.0;
266
267 
268  // Flux formulas
269  flux_left[0] = u_left*h_left;
270  flux_left[1] = u_left*uh_left + 0.5*g*h_left*h_left;
271  flux_left[2] = u_left*vh_left;
272
273  flux_right[0] = u_right*h_right;
274  flux_right[1] = u_right*uh_right + 0.5*g*h_right*h_right;
275  flux_right[2] = u_right*vh_right;
276
277 
278  // Flux computation
279  denom = s_max-s_min;
280  if (denom < epsilon) { // FIXME (Ole): Try using H0 here
281    for (i=0; i<3; i++) edgeflux[i] = 0.0;
282    *max_speed = 0.0;
283  } else {
284    for (i=0; i<3; i++) {
285      edgeflux[i] = s_max*flux_left[i] - s_min*flux_right[i];
286      edgeflux[i] += s_max*s_min*(q_right_rotated[i]-q_left_rotated[i]);
287      edgeflux[i] /= denom;
288    }
289
290    // Maximal wavespeed
291    *max_speed = max(fabs(s_max), fabs(s_min));
292
293    // Rotate back
294    _rotate(edgeflux, n1, -n2);
295  }
296 
297  return 0;
298}
299
300double erfcc(double x){
301    double z,t,result;
302
303    z=fabs(x);
304    t=1.0/(1.0+0.5*z);
305    result=t*exp(-z*z-1.26551223+t*(1.00002368+t*(.37409196+
306         t*(.09678418+t*(-.18628806+t*(.27886807+t*(-1.13520398+
307         t*(1.48851587+t*(-.82215223+t*.17087277)))))))));
308    if (x < 0.0) result = 2.0-result;
309
310    return result;
311    }
312
313
314
315// Computational function for flux computation (using stage w=z+h)
316int flux_function_kinetic(double *q_left, double *q_right,
317                  double z_left, double z_right,
318                  double n1, double n2,
319                  double epsilon, double H0, double g,
320                  double *edgeflux, double *max_speed) {
321
322  /*Compute fluxes between volumes for the shallow water wave equation
323    cast in terms of the 'stage', w = h+z using
324    the 'central scheme' as described in
325
326    Zhang et. al., Advances in Water Resources, 26(6), 2003, 635-647.
327  */
328
329  int i;
330
331  double w_left, h_left, uh_left, vh_left, u_left, F_left;
332  double w_right, h_right, uh_right, vh_right, u_right, F_right;
333  double s_min, s_max, soundspeed_left, soundspeed_right;
334  double z;
335  double q_left_rotated[3], q_right_rotated[3];
336
337  double h0 = H0*H0; //This ensures a good balance when h approaches H0.
338
339  //Copy conserved quantities to protect from modification
340  for (i=0; i<3; i++) {
341    q_left_rotated[i] = q_left[i];
342    q_right_rotated[i] = q_right[i];
343  }
344
345  //Align x- and y-momentum with x-axis
346  _rotate(q_left_rotated, n1, n2);
347  _rotate(q_right_rotated, n1, n2);
348
349  z = (z_left+z_right)/2; //Take average of field values
350
351  //Compute speeds in x-direction
352  w_left = q_left_rotated[0];         
353  h_left = w_left-z;
354  uh_left = q_left_rotated[1];
355  u_left =_compute_speed(&uh_left, &h_left, epsilon, h0);
356
357  w_right = q_right_rotated[0];
358  h_right = w_right-z;
359  uh_right = q_right_rotated[1];
360  u_right =_compute_speed(&uh_right, &h_right, epsilon, h0); 
361
362
363  //Momentum in y-direction
364  vh_left  = q_left_rotated[2];
365  vh_right = q_right_rotated[2];
366
367
368  //Maximal and minimal wave speeds
369  soundspeed_left  = sqrt(g*h_left);
370  soundspeed_right = sqrt(g*h_right);
371
372  s_max = max(u_left+soundspeed_left, u_right+soundspeed_right);
373  if (s_max < 0.0) s_max = 0.0;
374
375  s_min = min(u_left-soundspeed_left, u_right-soundspeed_right);
376  if (s_min > 0.0) s_min = 0.0;
377
378
379  F_left  = 0.0;
380  F_right = 0.0;
381  if (h_left > 0.0) F_left = u_left/sqrt(g*h_left);
382  if (h_right > 0.0) F_right = u_right/sqrt(g*h_right);
383
384  for (i=0; i<3; i++) edgeflux[i] = 0.0;
385  *max_speed = 0.0;
386
387  edgeflux[0] = h_left*u_left/2.0*erfcc(-F_left) +  \
388          h_left*sqrt(g*h_left)/2.0/sqrt(pi)*exp(-(F_left*F_left)) + \
389          h_right*u_right/2.0*erfcc(F_right) -  \
390          h_right*sqrt(g*h_right)/2.0/sqrt(pi)*exp(-(F_right*F_right));
391
392  edgeflux[1] = (h_left*u_left*u_left + g/2.0*h_left*h_left)/2.0*erfcc(-F_left) + \
393          u_left*h_left*sqrt(g*h_left)/2.0/sqrt(pi)*exp(-(F_left*F_left)) + \
394          (h_right*u_right*u_right + g/2.0*h_right*h_right)/2.0*erfcc(F_right) -  \
395          u_right*h_right*sqrt(g*h_right)/2.0/sqrt(pi)*exp(-(F_right*F_right));
396
397  edgeflux[2] = vh_left*u_left/2.0*erfcc(-F_left) + \
398          vh_left*sqrt(g*h_left)/2.0/sqrt(pi)*exp(-(F_left*F_left)) + \
399          vh_right*u_right/2.0*erfcc(F_right) - \
400          vh_right*sqrt(g*h_right)/2.0/sqrt(pi)*exp(-(F_right*F_right));
401
402  //Maximal wavespeed
403  *max_speed = max(fabs(s_max), fabs(s_min));
404
405  //Rotate back
406  _rotate(edgeflux, n1, -n2);
407
408  return 0;
409}
410
411
412
413
414void _manning_friction(double g, double eps, int N,
415                       double* w, double* z,
416                       double* uh, double* vh,
417                       double* eta, double* xmom, double* ymom) {
418
419  int k;
420  double S, h;
421
422  for (k=0; k<N; k++) {
423    if (eta[k] > eps) {
424      h = w[k]-z[k];
425      if (h >= eps) {
426        S = -g * eta[k]*eta[k] * sqrt((uh[k]*uh[k] + vh[k]*vh[k]));
427        S /= pow(h, 7.0/3);      //Expensive (on Ole's home computer)
428        //S /= exp(7.0/3.0*log(h));      //seems to save about 15% over manning_friction
429        //S /= h*h*(1 + h/3.0 - h*h/9.0); //FIXME: Could use a Taylor expansion
430
431
432        //Update momentum
433        xmom[k] += S*uh[k];
434        ymom[k] += S*vh[k];
435      }
436    }
437  }
438}
439
440
441/*
442void _manning_friction_explicit(double g, double eps, int N,
443                       double* w, double* z,
444                       double* uh, double* vh,
445                       double* eta, double* xmom, double* ymom) {
446
447  int k;
448  double S, h;
449
450  for (k=0; k<N; k++) {
451    if (eta[k] > eps) {
452      h = w[k]-z[k];
453      if (h >= eps) {
454        S = -g * eta[k]*eta[k] * sqrt((uh[k]*uh[k] + vh[k]*vh[k]));
455        S /= pow(h, 7.0/3);      //Expensive (on Ole's home computer)
456        //S /= exp(7.0/3.0*log(h));      //seems to save about 15% over manning_friction
457        //S /= h*h*(1 + h/3.0 - h*h/9.0); //FIXME: Could use a Taylor expansion
458
459
460        //Update momentum
461        xmom[k] += S*uh[k];
462        ymom[k] += S*vh[k];
463      }
464    }
465  }
466}
467*/
468
469
470
471
472double velocity_balance(double uh_i, 
473                        double uh,
474                        double h_i, 
475                        double h, 
476                        double alpha,
477                        double epsilon) {
478  // Find alpha such that speed at the vertex is within one
479  // order of magnitude of the centroid speed   
480
481  // FIXME(Ole): Work in progress
482 
483  double a, b, estimate;
484  double m=10; // One order of magnitude - allow for velocity deviations at vertices
485
486 
487  printf("alpha = %f, uh_i=%f, uh=%f, h_i=%f, h=%f\n",
488         alpha, uh_i, uh, h_i, h);
489     
490 
491 
492   
493  // Shorthands and determine inequality
494  if (fabs(uh) < epsilon) {
495    a = 1.0e10; // Limit
496  } else {
497    a = fabs(uh_i - uh)/fabs(uh);
498  }
499     
500  if (h < epsilon) {
501    b = 1.0e10; // Limit
502  } else {       
503    b = m*fabs(h_i - h)/h;
504  }
505
506  printf("a %f, b %f\n", a, b); 
507
508  if (a > b) {
509    estimate = (m-1)/(a-b);             
510   
511    printf("Alpha %f, estimate %f\n", 
512           alpha, estimate);   
513           
514    if (alpha < estimate) {
515      printf("Adjusting alpha from %f to %f\n", 
516             alpha, estimate);
517      alpha = estimate;
518    }
519  } else {
520 
521    if (h < h_i) {
522      estimate = (m-1)/(a-b);                 
523   
524      printf("Alpha %f, estimate %f\n", 
525             alpha, estimate);   
526           
527      if (alpha < estimate) {
528        printf("Adjusting alpha from %f to %f\n", 
529               alpha, estimate);
530        alpha = estimate;
531      }   
532    }
533    // Always fulfilled as alpha and m-1 are non negative
534  }
535 
536 
537  return alpha;
538}
539
540
541int _balance_deep_and_shallow(int N,
542                              double* wc,
543                              double* zc,
544                              double* wv,
545                              double* zv,
546                              double* hvbar, // Retire this
547                              double* xmomc,
548                              double* ymomc,
549                              double* xmomv,
550                              double* ymomv,
551                              double H0,
552                              int tight_slope_limiters,
553                              int use_centroid_velocities,
554                              double alpha_balance) {
555
556  int k, k3, i; 
557
558  double dz, hmin, alpha, h_diff, hc_k;
559  double epsilon = 1.0e-6; // FIXME: Temporary measure
560  double hv[3]; // Depths at vertices
561  double uc, vc; // Centroid speeds
562
563  // Compute linear combination between w-limited stages and
564  // h-limited stages close to the bed elevation.
565
566  for (k=0; k<N; k++) {
567    // Compute maximal variation in bed elevation
568    // This quantitiy is
569    //     dz = max_i abs(z_i - z_c)
570    // and it is independent of dimension
571    // In the 1d case zc = (z0+z1)/2
572    // In the 2d case zc = (z0+z1+z2)/3
573
574    k3 = 3*k;
575    hc_k = wc[k] - zc[k]; // Centroid value at triangle k
576   
577    dz = 0.0;
578    if (tight_slope_limiters == 0) {     
579      // FIXME: Try with this one precomputed
580      for (i=0; i<3; i++) {
581        dz = max(dz, fabs(zv[k3+i]-zc[k]));
582      }
583    }
584
585    // Calculate depth at vertices (possibly negative here!)
586    hv[0] = wv[k3] -   zv[k3];
587    hv[1] = wv[k3+1] - zv[k3+1];
588    hv[2] = wv[k3+2] - zv[k3+2];       
589   
590    // Calculate minimal depth across all three vertices
591    hmin = min(hv[0], min(hv[1], hv[2]));
592
593    //if (hmin < 0.0 ) {
594    //  printf("hmin = %f\n",hmin);
595    //}
596   
597
598    // Create alpha in [0,1], where alpha==0 means using the h-limited
599    // stage and alpha==1 means using the w-limited stage as
600    // computed by the gradient limiter (both 1st or 2nd order)
601    if (tight_slope_limiters == 0) {     
602      // If hmin > dz/alpha_balance then alpha = 1 and the bed will have no
603      // effect
604      // If hmin < 0 then alpha = 0 reverting to constant height above bed.
605      // The parameter alpha_balance==2 by default
606
607     
608      if (dz > 0.0) {
609        alpha = max( min( alpha_balance*hmin/dz, 1.0), 0.0 );     
610      } else {
611        alpha = 1.0;  // Flat bed
612      }
613      //printf("Using old style limiter\n");
614     
615    } else {
616
617      // Tight Slope Limiter (2007)
618   
619      // Make alpha as large as possible but still ensure that
620      // final depth is positive and that velocities at vertices
621      // are controlled
622   
623      if (hmin < H0) {
624        alpha = 1.0;
625        for (i=0; i<3; i++) {
626         
627          h_diff = hc_k - hv[i];         
628          if (h_diff <= 0) {
629            // Deep water triangle is further away from bed than
630            // shallow water (hbar < h). Any alpha will do
631         
632          } else { 
633            // Denominator is positive which means that we need some of the
634            // h-limited stage.
635           
636            alpha = min(alpha, (hc_k - H0)/h_diff);         
637          }
638        }
639
640        // Ensure alpha in [0,1]
641        if (alpha>1.0) alpha=1.0;
642        if (alpha<0.0) alpha=0.0;
643       
644      } else {
645        // Use w-limited stage exclusively in deeper water.
646        alpha = 1.0;       
647      }
648    }
649
650
651    //  Let
652    //
653    //    wvi be the w-limited stage (wvi = zvi + hvi)
654    //    wvi- be the h-limited state (wvi- = zvi + hvi-)
655    //
656    //
657    //  where i=0,1,2 denotes the vertex ids
658    //
659    //  Weighted balance between w-limited and h-limited stage is
660    //
661    //    wvi := (1-alpha)*(zvi+hvi-) + alpha*(zvi+hvi)
662    //
663    //  It follows that the updated wvi is
664    //    wvi := zvi + (1-alpha)*hvi- + alpha*hvi
665    //
666    //   Momentum is balanced between constant and limited
667
668
669    if (alpha < 1) {     
670      for (i=0; i<3; i++) {
671     
672        wv[k3+i] = zv[k3+i] + (1-alpha)*hc_k + alpha*hv[i];     
673
674        // Update momentum at vertices
675        if (use_centroid_velocities == 1) {
676          // This is a simple, efficient and robust option
677          // It uses first order approximation of velocities, but retains
678          // the order used by stage.
679       
680          // Speeds at centroids
681          if (hc_k > epsilon) {
682            uc = xmomc[k]/hc_k;
683            vc = ymomc[k]/hc_k;
684          } else {
685            uc = 0.0;
686            vc = 0.0;
687          }
688         
689          // Vertex momenta guaranteed to be consistent with depth guaranteeing
690          // controlled speed
691          hv[i] = wv[k3+i] - zv[k3+i]; // Recompute (balanced) vertex depth
692          xmomv[k3+i] = uc*hv[i];
693          ymomv[k3+i] = vc*hv[i];
694         
695        } else {
696          // Update momentum as a linear combination of
697          // xmomc and ymomc (shallow) and momentum
698          // from extrapolator xmomv and ymomv (deep).
699          // This assumes that values from xmomv and ymomv have
700          // been established e.g. by the gradient limiter.
701
702          // FIXME (Ole): I think this should be used with vertex momenta
703          // computed above using centroid_velocities instead of xmomc
704          // and ymomc as they'll be more representative first order
705          // values.
706         
707          xmomv[k3+i] = (1-alpha)*xmomc[k] + alpha*xmomv[k3+i];
708          ymomv[k3+i] = (1-alpha)*ymomc[k] + alpha*ymomv[k3+i];
709       
710        }
711      }
712    }
713  }
714  return 0;
715}
716 
717
718
719
720int _protect(int N,
721             double minimum_allowed_height,
722             double maximum_allowed_speed,
723             double epsilon,
724             double* wc,
725             double* zc,
726             double* xmomc,
727             double* ymomc) {
728
729  int k;
730  double hc;
731  double u, v, reduced_speed;
732
733
734  // Protect against initesimal and negative heights 
735  if (maximum_allowed_speed < epsilon) {
736    for (k=0; k<N; k++) {
737      hc = wc[k] - zc[k];
738
739      if (hc < minimum_allowed_height) {
740       
741        // Set momentum to zero and ensure h is non negative
742        xmomc[k] = 0.0;
743        ymomc[k] = 0.0;
744        if (hc <= 0.0) wc[k] = zc[k];
745      }
746    }
747  } else {
748   
749    // Protect against initesimal and negative heights
750    for (k=0; k<N; k++) {
751      hc = wc[k] - zc[k];
752
753      if (hc < minimum_allowed_height) {
754
755        //New code: Adjust momentum to guarantee speeds are physical
756        //          ensure h is non negative
757             
758        if (hc <= 0.0) {
759                wc[k] = zc[k];
760        xmomc[k] = 0.0;
761        ymomc[k] = 0.0;
762        } else {
763          //Reduce excessive speeds derived from division by small hc
764        //FIXME (Ole): This may be unnecessary with new slope limiters
765        //in effect.
766         
767          u = xmomc[k]/hc;
768          if (fabs(u) > maximum_allowed_speed) {
769            reduced_speed = maximum_allowed_speed * u/fabs(u);
770            //printf("Speed (u) has been reduced from %.3f to %.3f\n",
771            //   u, reduced_speed);
772            xmomc[k] = reduced_speed * hc;
773          }
774
775          v = ymomc[k]/hc;
776          if (fabs(v) > maximum_allowed_speed) {
777            reduced_speed = maximum_allowed_speed * v/fabs(v);
778            //printf("Speed (v) has been reduced from %.3f to %.3f\n",
779            //   v, reduced_speed);
780            ymomc[k] = reduced_speed * hc;
781          }
782        }
783      }
784    }
785  }
786  return 0;
787}
788
789
790
791int _assign_wind_field_values(int N,
792                              double* xmom_update,
793                              double* ymom_update,
794                              double* s_vec,
795                              double* phi_vec,
796                              double cw) {
797
798  // Assign windfield values to momentum updates
799
800  int k;
801  double S, s, phi, u, v;
802
803  for (k=0; k<N; k++) {
804
805    s = s_vec[k];
806    phi = phi_vec[k];
807
808    //Convert to radians
809    phi = phi*pi/180;
810
811    //Compute velocity vector (u, v)
812    u = s*cos(phi);
813    v = s*sin(phi);
814
815    //Compute wind stress
816    S = cw * sqrt(u*u + v*v);
817    xmom_update[k] += S*u;
818    ymom_update[k] += S*v;
819  }
820  return 0;
821}
822
823
824
825///////////////////////////////////////////////////////////////////
826// Gateways to Python
827
828
829
830PyObject *flux_function_central(PyObject *self, PyObject *args) {
831  //
832  // Gateway to innermost flux function.
833  // This is only used by the unit tests as the c implementation is
834  // normally called by compute_fluxes in this module.
835
836
837  PyArrayObject *normal, *ql, *qr,  *edgeflux;
838  double g, epsilon, max_speed, H0, zl, zr;
839
840  if (!PyArg_ParseTuple(args, "OOOddOddd",
841                        &normal, &ql, &qr, &zl, &zr, &edgeflux,
842                        &epsilon, &g, &H0)) {
843    PyErr_SetString(PyExc_RuntimeError, "shallow_water_ext.c: flux_function_central could not parse input arguments");
844    return NULL;
845  }
846
847 
848  _flux_function_central((double*) ql -> data, 
849                         (double*) qr -> data, 
850                         zl, 
851                         zr,                                             
852                         ((double*) normal -> data)[0],
853                         ((double*) normal -> data)[1],                 
854                         epsilon, H0, g,
855                         (double*) edgeflux -> data, 
856                         &max_speed);
857 
858  return Py_BuildValue("d", max_speed); 
859}
860
861
862
863PyObject *gravity(PyObject *self, PyObject *args) {
864  //
865  //  gravity(g, h, v, x, xmom, ymom)
866  //
867
868
869  PyArrayObject *h, *v, *x, *xmom, *ymom;
870  int k, N, k3, k6;
871  double g, avg_h, zx, zy;
872  double x0, y0, x1, y1, x2, y2, z0, z1, z2;
873  //double epsilon;
874
875  if (!PyArg_ParseTuple(args, "dOOOOO",
876                        &g, &h, &v, &x,
877                        &xmom, &ymom)) {
878    //&epsilon)) {
879    PyErr_SetString(PyExc_RuntimeError, "shallow_water_ext.c: gravity could not parse input arguments");
880    return NULL;
881  }
882
883  // check that numpy array objects arrays are C contiguous memory
884  CHECK_C_CONTIG(h);
885  CHECK_C_CONTIG(v);
886  CHECK_C_CONTIG(x);
887  CHECK_C_CONTIG(xmom);
888  CHECK_C_CONTIG(ymom);
889
890  N = h -> dimensions[0];
891  for (k=0; k<N; k++) {
892    k3 = 3*k;  // base index
893
894    // Get bathymetry
895    z0 = ((double*) v -> data)[k3 + 0];
896    z1 = ((double*) v -> data)[k3 + 1];
897    z2 = ((double*) v -> data)[k3 + 2];
898
899    // Optimise for flat bed
900    // Note (Ole): This didn't produce measurable speed up.
901    // Revisit later
902    //if (fabs(z0-z1)<epsilon && fabs(z1-z2)<epsilon) {
903    //  continue;
904    //}
905
906    // Get average depth from centroid values
907    avg_h = ((double *) h -> data)[k];
908
909    // Compute bed slope
910    k6 = 6*k;  // base index
911 
912    x0 = ((double*) x -> data)[k6 + 0];
913    y0 = ((double*) x -> data)[k6 + 1];
914    x1 = ((double*) x -> data)[k6 + 2];
915    y1 = ((double*) x -> data)[k6 + 3];
916    x2 = ((double*) x -> data)[k6 + 4];
917    y2 = ((double*) x -> data)[k6 + 5];
918
919
920    _gradient(x0, y0, x1, y1, x2, y2, z0, z1, z2, &zx, &zy);
921
922    // Update momentum
923    ((double*) xmom -> data)[k] += -g*zx*avg_h;
924    ((double*) ymom -> data)[k] += -g*zy*avg_h;
925  }
926
927  return Py_BuildValue("");
928}
929
930
931PyObject *manning_friction(PyObject *self, PyObject *args) {
932  //
933  // manning_friction(g, eps, h, uh, vh, eta, xmom_update, ymom_update)
934  //
935
936
937  PyArrayObject *w, *z, *uh, *vh, *eta, *xmom, *ymom;
938  int N;
939  double g, eps;
940
941  if (!PyArg_ParseTuple(args, "ddOOOOOOO",
942                        &g, &eps, &w, &z, &uh, &vh, &eta,
943                        &xmom, &ymom)) {
944    PyErr_SetString(PyExc_RuntimeError, "shallow_water_ext.c: manning_friction could not parse input arguments");
945    return NULL;
946  }
947
948  // check that numpy array objects arrays are C contiguous memory
949  CHECK_C_CONTIG(w);
950  CHECK_C_CONTIG(z);
951  CHECK_C_CONTIG(uh);
952  CHECK_C_CONTIG(vh);
953  CHECK_C_CONTIG(eta);
954  CHECK_C_CONTIG(xmom);
955  CHECK_C_CONTIG(ymom);
956
957  N = w -> dimensions[0];
958  _manning_friction(g, eps, N,
959                    (double*) w -> data,
960                    (double*) z -> data,
961                    (double*) uh -> data,
962                    (double*) vh -> data,
963                    (double*) eta -> data,
964                    (double*) xmom -> data,
965                    (double*) ymom -> data);
966
967  return Py_BuildValue("");
968}
969
970
971/*
972PyObject *manning_friction_explicit(PyObject *self, PyObject *args) {
973  //
974  // manning_friction_explicit(g, eps, h, uh, vh, eta, xmom_update, ymom_update)
975  //
976
977
978  PyArrayObject *w, *z, *uh, *vh, *eta, *xmom, *ymom;
979  int N;
980  double g, eps;
981
982  if (!PyArg_ParseTuple(args, "ddOOOOOOO",
983                        &g, &eps, &w, &z, &uh, &vh, &eta,
984                        &xmom, &ymom))
985    return NULL;
986
987  // check that numpy array objects arrays are C contiguous memory
988  CHECK_C_CONTIG(w);
989  CHECK_C_CONTIG(z);
990  CHECK_C_CONTIG(uh);
991  CHECK_C_CONTIG(vh);
992  CHECK_C_CONTIG(eta);
993  CHECK_C_CONTIG(xmom);
994  CHECK_C_CONTIG(ymom);
995
996  N = w -> dimensions[0];
997  _manning_friction_explicit(g, eps, N,
998                    (double*) w -> data,
999                    (double*) z -> data,
1000                    (double*) uh -> data,
1001                    (double*) vh -> data,
1002                    (double*) eta -> data,
1003                    (double*) xmom -> data,
1004                    (double*) ymom -> data);
1005
1006  return Py_BuildValue("");
1007}
1008*/
1009
1010
1011
1012// Computational routine
1013int _extrapolate_second_order_sw(int number_of_elements,
1014                                  double epsilon,
1015                                  double minimum_allowed_height,
1016                                  double beta_w,
1017                                  double beta_w_dry,
1018                                  double beta_uh,
1019                                  double beta_uh_dry,
1020                                  double beta_vh,
1021                                  double beta_vh_dry,
1022                                  long* surrogate_neighbours,
1023                                  long* number_of_boundaries,
1024                                  double* centroid_coordinates,
1025                                  double* stage_centroid_values,
1026                                  double* xmom_centroid_values,
1027                                  double* ymom_centroid_values,
1028                                  double* elevation_centroid_values,
1029                                  double* vertex_coordinates,
1030                                  double* stage_vertex_values,
1031                                  double* xmom_vertex_values,
1032                                  double* ymom_vertex_values,
1033                                  double* elevation_vertex_values,
1034                                  int optimise_dry_cells) {
1035                                 
1036                                 
1037
1038  // Local variables
1039  double a, b; // Gradient vector used to calculate vertex values from centroids
1040  int k,k0,k1,k2,k3,k6,coord_index,i;
1041  double x,y,x0,y0,x1,y1,x2,y2,xv0,yv0,xv1,yv1,xv2,yv2; // Vertices of the auxiliary triangle
1042  double dx1,dx2,dy1,dy2,dxv0,dxv1,dxv2,dyv0,dyv1,dyv2,dq0,dq1,dq2,area2;
1043  double dqv[3], qmin, qmax, hmin, hmax;
1044  double hc, h0, h1, h2, beta_tmp, hfactor;
1045
1046       
1047  for (k=0; k<number_of_elements; k++) {
1048    k3=k*3;
1049    k6=k*6;
1050
1051   
1052    if (number_of_boundaries[k]==3){
1053      // No neighbours, set gradient on the triangle to zero
1054     
1055      stage_vertex_values[k3]   = stage_centroid_values[k];
1056      stage_vertex_values[k3+1] = stage_centroid_values[k];
1057      stage_vertex_values[k3+2] = stage_centroid_values[k];
1058      xmom_vertex_values[k3]    = xmom_centroid_values[k];
1059      xmom_vertex_values[k3+1]  = xmom_centroid_values[k];
1060      xmom_vertex_values[k3+2]  = xmom_centroid_values[k];
1061      ymom_vertex_values[k3]    = ymom_centroid_values[k];
1062      ymom_vertex_values[k3+1]  = ymom_centroid_values[k];
1063      ymom_vertex_values[k3+2]  = ymom_centroid_values[k];
1064     
1065      continue;
1066    }
1067    else {
1068      // Triangle k has one or more neighbours.
1069      // Get centroid and vertex coordinates of the triangle
1070     
1071      // Get the vertex coordinates
1072      xv0 = vertex_coordinates[k6];   yv0=vertex_coordinates[k6+1];
1073      xv1 = vertex_coordinates[k6+2]; yv1=vertex_coordinates[k6+3];
1074      xv2 = vertex_coordinates[k6+4]; yv2=vertex_coordinates[k6+5];
1075     
1076      // Get the centroid coordinates
1077      coord_index=2*k;
1078      x=centroid_coordinates[coord_index];
1079      y=centroid_coordinates[coord_index+1];
1080     
1081      // Store x- and y- differentials for the vertices of
1082      // triangle k relative to the centroid
1083      dxv0=xv0-x; dxv1=xv1-x; dxv2=xv2-x;
1084      dyv0=yv0-y; dyv1=yv1-y; dyv2=yv2-y;
1085    }
1086
1087
1088           
1089    if (number_of_boundaries[k]<=1){
1090   
1091      //==============================================
1092      // Number of boundaries <= 1
1093      //==============================================   
1094   
1095   
1096      // If no boundaries, auxiliary triangle is formed
1097      // from the centroids of the three neighbours
1098      // If one boundary, auxiliary triangle is formed
1099      // from this centroid and its two neighbours
1100     
1101      k0=surrogate_neighbours[k3];
1102      k1=surrogate_neighbours[k3+1];
1103      k2=surrogate_neighbours[k3+2];
1104     
1105      // Get the auxiliary triangle's vertex coordinates
1106      // (really the centroids of neighbouring triangles)
1107      coord_index=2*k0;
1108      x0=centroid_coordinates[coord_index];
1109      y0=centroid_coordinates[coord_index+1];
1110     
1111      coord_index=2*k1;
1112      x1=centroid_coordinates[coord_index];
1113      y1=centroid_coordinates[coord_index+1];
1114     
1115      coord_index=2*k2;
1116      x2=centroid_coordinates[coord_index];
1117      y2=centroid_coordinates[coord_index+1];
1118     
1119      // Store x- and y- differentials for the vertices
1120      // of the auxiliary triangle
1121      dx1=x1-x0; dx2=x2-x0;
1122      dy1=y1-y0; dy2=y2-y0;
1123     
1124      // Calculate 2*area of the auxiliary triangle
1125      // The triangle is guaranteed to be counter-clockwise     
1126      area2 = dy2*dx1 - dy1*dx2; 
1127     
1128      // If the mesh is 'weird' near the boundary,
1129      // the triangle might be flat or clockwise:
1130      if (area2<=0) {
1131        PyErr_SetString(PyExc_RuntimeError, 
1132                        "shallow_water_ext.c: negative triangle area encountered");
1133        return -1;
1134      } 
1135     
1136      // Calculate heights of neighbouring cells
1137      hc = stage_centroid_values[k]  - elevation_centroid_values[k];
1138      h0 = stage_centroid_values[k0] - elevation_centroid_values[k0];
1139      h1 = stage_centroid_values[k1] - elevation_centroid_values[k1];
1140      h2 = stage_centroid_values[k2] - elevation_centroid_values[k2];
1141      hmin = min(min(h0,min(h1,h2)),hc);
1142      //hfactor = hc/(hc + 1.0);
1143
1144      hfactor = 0.0;
1145      if (hmin > 0.001 ) {
1146          hfactor = (hmin-0.001)/(hmin+0.004);
1147      }
1148     
1149      if (optimise_dry_cells) {     
1150        // Check if linear reconstruction is necessary for triangle k
1151        // This check will exclude dry cells.
1152
1153        hmax = max(h0,max(h1,h2));     
1154        if (hmax < epsilon) {
1155          continue;
1156        }
1157      }
1158
1159           
1160      //-----------------------------------
1161      // stage
1162      //-----------------------------------     
1163     
1164      // Calculate the difference between vertex 0 of the auxiliary
1165      // triangle and the centroid of triangle k
1166      dq0=stage_centroid_values[k0]-stage_centroid_values[k];
1167     
1168      // Calculate differentials between the vertices
1169      // of the auxiliary triangle (centroids of neighbouring triangles)
1170      dq1=stage_centroid_values[k1]-stage_centroid_values[k0];
1171      dq2=stage_centroid_values[k2]-stage_centroid_values[k0];
1172     
1173      // Calculate the gradient of stage on the auxiliary triangle
1174      a = dy2*dq1 - dy1*dq2;
1175      a /= area2;
1176      b = dx1*dq2 - dx2*dq1;
1177      b /= area2;
1178     
1179      // Calculate provisional jumps in stage from the centroid
1180      // of triangle k to its vertices, to be limited
1181      dqv[0]=a*dxv0+b*dyv0;
1182      dqv[1]=a*dxv1+b*dyv1;
1183      dqv[2]=a*dxv2+b*dyv2;
1184     
1185      // Now we want to find min and max of the centroid and the
1186      // vertices of the auxiliary triangle and compute jumps
1187      // from the centroid to the min and max
1188      find_qmin_and_qmax(dq0,dq1,dq2,&qmin,&qmax);
1189     
1190      // Playing with dry wet interface
1191      //hmin = qmin;
1192      //beta_tmp = beta_w_dry;
1193      //if (hmin>minimum_allowed_height)
1194      beta_tmp = beta_w_dry + (beta_w - beta_w_dry) * hfactor;
1195       
1196      //printf("min_alled_height = %f\n",minimum_allowed_height);
1197      //printf("hmin = %f\n",hmin);
1198      //printf("beta_w = %f\n",beta_w);
1199      //printf("beta_tmp = %f\n",beta_tmp);
1200      // Limit the gradient
1201      limit_gradient(dqv,qmin,qmax,beta_tmp); 
1202     
1203      for (i=0;i<3;i++)
1204        stage_vertex_values[k3+i]=stage_centroid_values[k]+dqv[i];
1205     
1206     
1207      //-----------------------------------
1208      // xmomentum
1209      //-----------------------------------           
1210
1211      // Calculate the difference between vertex 0 of the auxiliary
1212      // triangle and the centroid of triangle k     
1213      dq0=xmom_centroid_values[k0]-xmom_centroid_values[k];
1214     
1215      // Calculate differentials between the vertices
1216      // of the auxiliary triangle
1217      dq1=xmom_centroid_values[k1]-xmom_centroid_values[k0];
1218      dq2=xmom_centroid_values[k2]-xmom_centroid_values[k0];
1219     
1220      // Calculate the gradient of xmom on the auxiliary triangle
1221      a = dy2*dq1 - dy1*dq2;
1222      a /= area2;
1223      b = dx1*dq2 - dx2*dq1;
1224      b /= area2;
1225     
1226      // Calculate provisional jumps in stage from the centroid
1227      // of triangle k to its vertices, to be limited     
1228      dqv[0]=a*dxv0+b*dyv0;
1229      dqv[1]=a*dxv1+b*dyv1;
1230      dqv[2]=a*dxv2+b*dyv2;
1231     
1232      // Now we want to find min and max of the centroid and the
1233      // vertices of the auxiliary triangle and compute jumps
1234      // from the centroid to the min and max
1235      find_qmin_and_qmax(dq0,dq1,dq2,&qmin,&qmax);
1236      //beta_tmp = beta_uh;
1237      //if (hmin<minimum_allowed_height)
1238      //beta_tmp = beta_uh_dry;
1239      beta_tmp = beta_uh_dry + (beta_uh - beta_uh_dry) * hfactor;
1240
1241      // Limit the gradient
1242      limit_gradient(dqv,qmin,qmax,beta_tmp);
1243
1244      for (i=0;i<3;i++)
1245        xmom_vertex_values[k3+i]=xmom_centroid_values[k]+dqv[i];
1246     
1247     
1248      //-----------------------------------
1249      // ymomentum
1250      //-----------------------------------                 
1251
1252      // Calculate the difference between vertex 0 of the auxiliary
1253      // triangle and the centroid of triangle k
1254      dq0=ymom_centroid_values[k0]-ymom_centroid_values[k];
1255     
1256      // Calculate differentials between the vertices
1257      // of the auxiliary triangle
1258      dq1=ymom_centroid_values[k1]-ymom_centroid_values[k0];
1259      dq2=ymom_centroid_values[k2]-ymom_centroid_values[k0];
1260     
1261      // Calculate the gradient of xmom on the auxiliary triangle
1262      a = dy2*dq1 - dy1*dq2;
1263      a /= area2;
1264      b = dx1*dq2 - dx2*dq1;
1265      b /= area2;
1266     
1267      // Calculate provisional jumps in stage from the centroid
1268      // of triangle k to its vertices, to be limited
1269      dqv[0]=a*dxv0+b*dyv0;
1270      dqv[1]=a*dxv1+b*dyv1;
1271      dqv[2]=a*dxv2+b*dyv2;
1272     
1273      // Now we want to find min and max of the centroid and the
1274      // vertices of the auxiliary triangle and compute jumps
1275      // from the centroid to the min and max
1276      find_qmin_and_qmax(dq0,dq1,dq2,&qmin,&qmax);
1277     
1278      //beta_tmp = beta_vh;
1279      //
1280      //if (hmin<minimum_allowed_height)
1281      //beta_tmp = beta_vh_dry;
1282      beta_tmp = beta_vh_dry + (beta_vh - beta_vh_dry) * hfactor;       
1283
1284      // Limit the gradient
1285      limit_gradient(dqv,qmin,qmax,beta_tmp);
1286     
1287      for (i=0;i<3;i++)
1288        ymom_vertex_values[k3+i]=ymom_centroid_values[k]+dqv[i];
1289       
1290    } // End number_of_boundaries <=1
1291    else{
1292
1293      //==============================================
1294      // Number of boundaries == 2
1295      //==============================================       
1296       
1297      // One internal neighbour and gradient is in direction of the neighbour's centroid
1298     
1299      // Find the only internal neighbour (k1?)
1300      for (k2=k3;k2<k3+3;k2++){
1301        // Find internal neighbour of triangle k     
1302        // k2 indexes the edges of triangle k   
1303     
1304        if (surrogate_neighbours[k2]!=k)
1305          break;
1306      }
1307     
1308      if ((k2==k3+3)) {
1309        // If we didn't find an internal neighbour
1310        PyErr_SetString(PyExc_RuntimeError, 
1311                        "shallow_water_ext.c: Internal neighbour not found");     
1312        return -1;
1313      }
1314     
1315      k1=surrogate_neighbours[k2];
1316     
1317      // The coordinates of the triangle are already (x,y).
1318      // Get centroid of the neighbour (x1,y1)
1319      coord_index=2*k1;
1320      x1=centroid_coordinates[coord_index];
1321      y1=centroid_coordinates[coord_index+1];
1322     
1323      // Compute x- and y- distances between the centroid of
1324      // triangle k and that of its neighbour
1325      dx1=x1-x; dy1=y1-y;
1326     
1327      // Set area2 as the square of the distance
1328      area2=dx1*dx1+dy1*dy1;
1329     
1330      // Set dx2=(x1-x0)/((x1-x0)^2+(y1-y0)^2)
1331      // and dy2=(y1-y0)/((x1-x0)^2+(y1-y0)^2) which
1332      // respectively correspond to the x- and y- gradients
1333      // of the conserved quantities
1334      dx2=1.0/area2;
1335      dy2=dx2*dy1;
1336      dx2*=dx1;
1337     
1338     
1339      //-----------------------------------
1340      // stage
1341      //-----------------------------------           
1342
1343      // Compute differentials
1344      dq1=stage_centroid_values[k1]-stage_centroid_values[k];
1345     
1346      // Calculate the gradient between the centroid of triangle k
1347      // and that of its neighbour
1348      a=dq1*dx2;
1349      b=dq1*dy2;
1350     
1351      // Calculate provisional vertex jumps, to be limited
1352      dqv[0]=a*dxv0+b*dyv0;
1353      dqv[1]=a*dxv1+b*dyv1;
1354      dqv[2]=a*dxv2+b*dyv2;
1355     
1356      // Now limit the jumps
1357      if (dq1>=0.0){
1358        qmin=0.0;
1359        qmax=dq1;
1360      }
1361      else{
1362        qmin=dq1;
1363        qmax=0.0;
1364      }
1365     
1366      // Limit the gradient
1367      limit_gradient(dqv,qmin,qmax,beta_w);
1368     
1369      for (i=0;i<3;i++)
1370        stage_vertex_values[k3+i]=stage_centroid_values[k]+dqv[i];
1371     
1372     
1373      //-----------------------------------
1374      // xmomentum
1375      //-----------------------------------                       
1376     
1377      // Compute differentials
1378      dq1=xmom_centroid_values[k1]-xmom_centroid_values[k];
1379     
1380      // Calculate the gradient between the centroid of triangle k
1381      // and that of its neighbour
1382      a=dq1*dx2;
1383      b=dq1*dy2;
1384     
1385      // Calculate provisional vertex jumps, to be limited
1386      dqv[0]=a*dxv0+b*dyv0;
1387      dqv[1]=a*dxv1+b*dyv1;
1388      dqv[2]=a*dxv2+b*dyv2;
1389     
1390      // Now limit the jumps
1391      if (dq1>=0.0){
1392        qmin=0.0;
1393        qmax=dq1;
1394      }
1395      else{
1396        qmin=dq1;
1397        qmax=0.0;
1398      }
1399     
1400      // Limit the gradient     
1401      limit_gradient(dqv,qmin,qmax,beta_w);
1402     
1403      for (i=0;i<3;i++)
1404        xmom_vertex_values[k3+i]=xmom_centroid_values[k]+dqv[i];
1405     
1406     
1407      //-----------------------------------
1408      // ymomentum
1409      //-----------------------------------                       
1410
1411      // Compute differentials
1412      dq1=ymom_centroid_values[k1]-ymom_centroid_values[k];
1413     
1414      // Calculate the gradient between the centroid of triangle k
1415      // and that of its neighbour
1416      a=dq1*dx2;
1417      b=dq1*dy2;
1418     
1419      // Calculate provisional vertex jumps, to be limited
1420      dqv[0]=a*dxv0+b*dyv0;
1421      dqv[1]=a*dxv1+b*dyv1;
1422      dqv[2]=a*dxv2+b*dyv2;
1423     
1424      // Now limit the jumps
1425      if (dq1>=0.0){
1426        qmin=0.0;
1427        qmax=dq1;
1428      }
1429      else{
1430        qmin=dq1;
1431        qmax=0.0;
1432      }
1433     
1434      // Limit the gradient
1435      limit_gradient(dqv,qmin,qmax,beta_w);
1436     
1437      for (i=0;i<3;i++)
1438        ymom_vertex_values[k3+i]=ymom_centroid_values[k]+dqv[i];
1439       
1440    } // else [number_of_boundaries==2]
1441  } // for k=0 to number_of_elements-1
1442 
1443  return 0;
1444}                       
1445
1446
1447PyObject *extrapolate_second_order_sw(PyObject *self, PyObject *args) {
1448  /*Compute the vertex values based on a linear reconstruction
1449    on each triangle
1450   
1451    These values are calculated as follows:
1452    1) For each triangle not adjacent to a boundary, we consider the
1453       auxiliary triangle formed by the centroids of its three
1454       neighbours.
1455    2) For each conserved quantity, we integrate around the auxiliary
1456       triangle's boundary the product of the quantity and the outward
1457       normal vector. Dividing by the triangle area gives (a,b), the
1458       average of the vector (q_x,q_y) on the auxiliary triangle.
1459       We suppose that the linear reconstruction on the original
1460       triangle has gradient (a,b).
1461    3) Provisional vertex jumps dqv[0,1,2] are computed and these are
1462       then limited by calling the functions find_qmin_and_qmax and
1463       limit_gradient
1464
1465    Python call:
1466    extrapolate_second_order_sw(domain.surrogate_neighbours,
1467                                domain.number_of_boundaries
1468                                domain.centroid_coordinates,
1469                                Stage.centroid_values
1470                                Xmom.centroid_values
1471                                Ymom.centroid_values
1472                                domain.vertex_coordinates,
1473                                Stage.vertex_values,
1474                                Xmom.vertex_values,
1475                                Ymom.vertex_values)
1476
1477    Post conditions:
1478            The vertices of each triangle have values from a
1479            limited linear reconstruction
1480            based on centroid values
1481
1482  */
1483  PyArrayObject *surrogate_neighbours,
1484    *number_of_boundaries,
1485    *centroid_coordinates,
1486    *stage_centroid_values,
1487    *xmom_centroid_values,
1488    *ymom_centroid_values,
1489    *elevation_centroid_values,
1490    *vertex_coordinates,
1491    *stage_vertex_values,
1492    *xmom_vertex_values,
1493    *ymom_vertex_values,
1494    *elevation_vertex_values;
1495 
1496  PyObject *domain;
1497
1498 
1499  double beta_w, beta_w_dry, beta_uh, beta_uh_dry, beta_vh, beta_vh_dry;   
1500  double minimum_allowed_height, epsilon;
1501  int optimise_dry_cells, number_of_elements, e;
1502 
1503  // Provisional jumps from centroids to v'tices and safety factor re limiting
1504  // by which these jumps are limited
1505  // Convert Python arguments to C
1506  if (!PyArg_ParseTuple(args, "OOOOOOOOOOOOOi",
1507                        &domain,
1508                        &surrogate_neighbours,
1509                        &number_of_boundaries,
1510                        &centroid_coordinates,
1511                        &stage_centroid_values,
1512                        &xmom_centroid_values,
1513                        &ymom_centroid_values,
1514                        &elevation_centroid_values,
1515                        &vertex_coordinates,
1516                        &stage_vertex_values,
1517                        &xmom_vertex_values,
1518                        &ymom_vertex_values,
1519                        &elevation_vertex_values,
1520                        &optimise_dry_cells)) {                 
1521                       
1522    PyErr_SetString(PyExc_RuntimeError, 
1523                    "Input arguments to extrapolate_second_order_sw failed");
1524    return NULL;
1525  }
1526
1527  // check that numpy array objects arrays are C contiguous memory
1528  CHECK_C_CONTIG(surrogate_neighbours);
1529  CHECK_C_CONTIG(number_of_boundaries);
1530  CHECK_C_CONTIG(centroid_coordinates);
1531  CHECK_C_CONTIG(stage_centroid_values);
1532  CHECK_C_CONTIG(xmom_centroid_values);
1533  CHECK_C_CONTIG(ymom_centroid_values);
1534  CHECK_C_CONTIG(elevation_centroid_values);
1535  CHECK_C_CONTIG(vertex_coordinates);
1536  CHECK_C_CONTIG(stage_vertex_values);
1537  CHECK_C_CONTIG(xmom_vertex_values);
1538  CHECK_C_CONTIG(ymom_vertex_values);
1539  CHECK_C_CONTIG(elevation_vertex_values);
1540 
1541  // Get the safety factor beta_w, set in the config.py file.
1542  // This is used in the limiting process
1543 
1544
1545  beta_w                 = get_python_double(domain,"beta_w");
1546  beta_w_dry             = get_python_double(domain,"beta_w_dry");
1547  beta_uh                = get_python_double(domain,"beta_uh");
1548  beta_uh_dry            = get_python_double(domain,"beta_uh_dry");
1549  beta_vh                = get_python_double(domain,"beta_vh");
1550  beta_vh_dry            = get_python_double(domain,"beta_vh_dry"); 
1551
1552  minimum_allowed_height = get_python_double(domain,"minimum_allowed_height");
1553  epsilon                = get_python_double(domain,"epsilon");
1554
1555  number_of_elements = stage_centroid_values -> dimensions[0]; 
1556
1557  // Call underlying computational routine
1558  e = _extrapolate_second_order_sw(number_of_elements,
1559                                   epsilon,
1560                                   minimum_allowed_height,
1561                                   beta_w,
1562                                   beta_w_dry,
1563                                   beta_uh,
1564                                   beta_uh_dry,
1565                                   beta_vh,
1566                                   beta_vh_dry,
1567                                   (long*) surrogate_neighbours -> data,
1568                                   (long*) number_of_boundaries -> data,
1569                                   (double*) centroid_coordinates -> data,
1570                                   (double*) stage_centroid_values -> data,
1571                                   (double*) xmom_centroid_values -> data,
1572                                   (double*) ymom_centroid_values -> data,
1573                                   (double*) elevation_centroid_values -> data,
1574                                   (double*) vertex_coordinates -> data,
1575                                   (double*) stage_vertex_values -> data,
1576                                   (double*) xmom_vertex_values -> data,
1577                                   (double*) ymom_vertex_values -> data,
1578                                   (double*) elevation_vertex_values -> data,
1579                                   optimise_dry_cells);
1580  if (e == -1) {
1581    // Use error string set inside computational routine
1582    return NULL;
1583  }                               
1584 
1585 
1586  return Py_BuildValue("");
1587 
1588}// extrapolate_second-order_sw
1589
1590
1591
1592
1593PyObject *rotate(PyObject *self, PyObject *args, PyObject *kwargs) {
1594  //
1595  // r = rotate(q, normal, direction=1)
1596  //
1597  // Where q is assumed to be a Float numeric array of length 3 and
1598  // normal a Float numeric array of length 2.
1599
1600  // FIXME(Ole): I don't think this is used anymore
1601
1602  PyObject *Q, *Normal;
1603  PyArrayObject *q, *r, *normal;
1604
1605  static char *argnames[] = {"q", "normal", "direction", NULL};
1606  int dimensions[1], i, direction=1;
1607  double n1, n2;
1608
1609  // Convert Python arguments to C
1610  if (!PyArg_ParseTupleAndKeywords(args, kwargs, "OO|i", argnames,
1611                                   &Q, &Normal, &direction)) {
1612    PyErr_SetString(PyExc_RuntimeError, "shallow_water_ext.c: rotate could not parse input arguments");
1613    return NULL;
1614  } 
1615
1616  // Input checks (convert sequences into numeric arrays)
1617  q = (PyArrayObject *)
1618    PyArray_ContiguousFromObject(Q, PyArray_DOUBLE, 0, 0);
1619  normal = (PyArrayObject *)
1620    PyArray_ContiguousFromObject(Normal, PyArray_DOUBLE, 0, 0);
1621
1622
1623  if (normal -> dimensions[0] != 2) {
1624    PyErr_SetString(PyExc_RuntimeError, "Normal vector must have 2 components");
1625    return NULL;
1626  }
1627
1628  // Allocate space for return vector r (don't DECREF)
1629  dimensions[0] = 3;
1630  r = (PyArrayObject *) anuga_FromDims(1, dimensions, PyArray_DOUBLE);
1631
1632  // Copy
1633  for (i=0; i<3; i++) {
1634    ((double *) (r -> data))[i] = ((double *) (q -> data))[i];
1635  }
1636
1637  // Get normal and direction
1638  n1 = ((double *) normal -> data)[0];
1639  n2 = ((double *) normal -> data)[1];
1640  if (direction == -1) n2 = -n2;
1641
1642  // Rotate
1643  _rotate((double *) r -> data, n1, n2);
1644
1645  // Release numeric arrays
1646  Py_DECREF(q);
1647  Py_DECREF(normal);
1648
1649  // Return result using PyArray to avoid memory leak
1650  return PyArray_Return(r);
1651}
1652
1653
1654// Computational function for flux computation
1655double _compute_fluxes_central(int number_of_elements,
1656                               double timestep,
1657                               double epsilon,
1658                               double H0,
1659                               double g,
1660                               long* neighbours,
1661                               long* neighbour_edges,
1662                               double* normals,
1663                               double* edgelengths, 
1664                               double* radii, 
1665                               double* areas,
1666                               long* tri_full_flag,
1667                               double* stage_edge_values,
1668                               double* xmom_edge_values,
1669                               double* ymom_edge_values,
1670                               double* bed_edge_values,
1671                               double* stage_boundary_values,
1672                               double* xmom_boundary_values,
1673                               double* ymom_boundary_values,
1674                               double* stage_explicit_update,
1675                               double* xmom_explicit_update,
1676                               double* ymom_explicit_update,
1677                               long* already_computed_flux,
1678                               double* max_speed_array,
1679                               int optimise_dry_cells) {
1680                               
1681  // Local variables
1682  double max_speed, length, area, zl, zr;
1683  int k, i, m, n;
1684  int ki, nm=0, ki2; // Index shorthands
1685 
1686  // Workspace (making them static actually made function slightly slower (Ole)) 
1687  double ql[3], qr[3], edgeflux[3]; // Work array for summing up fluxes
1688
1689  static long call=1; // Static local variable flagging already computed flux
1690                               
1691       
1692  // Start computation
1693  call++; // Flag 'id' of flux calculation for this timestep
1694 
1695  // Set explicit_update to zero for all conserved_quantities.
1696  // This assumes compute_fluxes called before forcing terms
1697  for (k=0; k<number_of_elements; k++) {
1698    stage_explicit_update[k]=0.0;
1699    xmom_explicit_update[k]=0.0;
1700    ymom_explicit_update[k]=0.0; 
1701  }
1702
1703  // For all triangles
1704  for (k=0; k<number_of_elements; k++) {
1705   
1706    // Loop through neighbours and compute edge flux for each 
1707    for (i=0; i<3; i++) {
1708      ki = k*3+i; // Linear index (triangle k, edge i)
1709     
1710      if (already_computed_flux[ki] == call)
1711        // We've already computed the flux across this edge
1712        continue;
1713       
1714       
1715      ql[0] = stage_edge_values[ki];
1716      ql[1] = xmom_edge_values[ki];
1717      ql[2] = ymom_edge_values[ki];
1718      zl = bed_edge_values[ki];
1719
1720      // Quantities at neighbour on nearest face
1721      n = neighbours[ki];
1722      if (n < 0) {
1723        // Neighbour is a boundary condition
1724        m = -n-1; // Convert negative flag to boundary index
1725       
1726        qr[0] = stage_boundary_values[m];
1727        qr[1] = xmom_boundary_values[m];
1728        qr[2] = ymom_boundary_values[m];
1729        zr = zl; // Extend bed elevation to boundary
1730      } else {
1731        // Neighbour is a real element
1732        m = neighbour_edges[ki];
1733        nm = n*3+m; // Linear index (triangle n, edge m)
1734       
1735        qr[0] = stage_edge_values[nm];
1736        qr[1] = xmom_edge_values[nm];
1737        qr[2] = ymom_edge_values[nm];
1738        zr = bed_edge_values[nm];
1739      }
1740     
1741     
1742      if (optimise_dry_cells) {     
1743        // Check if flux calculation is necessary across this edge
1744        // This check will exclude dry cells.
1745        // This will also optimise cases where zl != zr as
1746        // long as both are dry
1747
1748        if ( fabs(ql[0] - zl) < epsilon && 
1749             fabs(qr[0] - zr) < epsilon ) {
1750          // Cell boundary is dry
1751         
1752          already_computed_flux[ki] = call; // #k Done 
1753          if (n>=0)
1754            already_computed_flux[nm] = call; // #n Done
1755       
1756          max_speed = 0.0;
1757          continue;
1758        }
1759      }
1760   
1761     
1762      // Outward pointing normal vector (domain.normals[k, 2*i:2*i+2])
1763      ki2 = 2*ki; //k*6 + i*2
1764
1765      // Edge flux computation (triangle k, edge i)
1766      _flux_function_central(ql, qr, zl, zr,
1767                             normals[ki2], normals[ki2+1],
1768                             epsilon, H0, g,
1769                             edgeflux, &max_speed);
1770     
1771     
1772      // Multiply edgeflux by edgelength
1773      length = edgelengths[ki];
1774      edgeflux[0] *= length;           
1775      edgeflux[1] *= length;           
1776      edgeflux[2] *= length;                       
1777     
1778     
1779      // Update triangle k with flux from edge i
1780      stage_explicit_update[k] -= edgeflux[0];
1781      xmom_explicit_update[k] -= edgeflux[1];
1782      ymom_explicit_update[k] -= edgeflux[2];
1783     
1784      already_computed_flux[ki] = call; // #k Done
1785     
1786     
1787      // Update neighbour n with same flux but reversed sign
1788      if (n>=0) {
1789        stage_explicit_update[n] += edgeflux[0];
1790        xmom_explicit_update[n] += edgeflux[1];
1791        ymom_explicit_update[n] += edgeflux[2];
1792       
1793        already_computed_flux[nm] = call; // #n Done
1794      }
1795
1796
1797      // Update timestep based on edge i and possibly neighbour n
1798      if (tri_full_flag[k] == 1) {
1799        if (max_speed > epsilon) {
1800       
1801          // Apply CFL condition for triangles joining this edge (triangle k and triangle n)
1802         
1803          // CFL for triangle k
1804          timestep = min(timestep, radii[k]/max_speed);
1805         
1806          if (n>=0) 
1807            // Apply CFL condition for neigbour n (which is on the ith edge of triangle k)
1808            timestep = min(timestep, radii[n]/max_speed);
1809         
1810          // Ted Rigby's suggested less conservative version
1811          //if (n>=0) {             
1812          //  timestep = min(timestep, (radii[k]+radii[n])/max_speed);
1813          //} else {
1814          //  timestep = min(timestep, radii[k]/max_speed);
1815          // }
1816        }
1817      }
1818     
1819    } // End edge i (and neighbour n)
1820   
1821   
1822    // Normalise triangle k by area and store for when all conserved
1823    // quantities get updated
1824    area = areas[k];
1825    stage_explicit_update[k] /= area;
1826    xmom_explicit_update[k] /= area;
1827    ymom_explicit_update[k] /= area;
1828   
1829   
1830    // Keep track of maximal speeds
1831    max_speed_array[k] = max_speed;
1832   
1833  } // End triangle k
1834       
1835                               
1836                               
1837  return timestep;
1838}
1839
1840//=========================================================================
1841// Python Glue
1842//=========================================================================
1843
1844PyObject *compute_fluxes_ext_central_new(PyObject *self, PyObject *args) {
1845  /*Compute all fluxes and the timestep suitable for all volumes
1846    in domain.
1847
1848    Compute total flux for each conserved quantity using "flux_function_central"
1849
1850    Fluxes across each edge are scaled by edgelengths and summed up
1851    Resulting flux is then scaled by area and stored in
1852    explicit_update for each of the three conserved quantities
1853    stage, xmomentum and ymomentum
1854
1855    The maximal allowable speed computed by the flux_function for each volume
1856    is converted to a timestep that must not be exceeded. The minimum of
1857    those is computed as the next overall timestep.
1858
1859    Python call:
1860    timestep = compute_fluxes(timestep, domain, stage, xmom, ymom, bed)
1861
1862
1863    Post conditions:
1864      domain.explicit_update is reset to computed flux values
1865      returns timestep which is the largest step satisfying all volumes.
1866
1867
1868  */
1869
1870    PyObject
1871        *domain,
1872        *stage, 
1873        *xmom, 
1874        *ymom, 
1875        *bed;
1876
1877    PyArrayObject
1878        *neighbours, 
1879        *neighbour_edges,
1880        *normals, 
1881        *edgelengths, 
1882        *radii, 
1883        *areas,
1884        *tri_full_flag,
1885        *stage_edge_values,
1886        *xmom_edge_values,
1887        *ymom_edge_values,
1888        *bed_edge_values,
1889        *stage_boundary_values,
1890        *xmom_boundary_values,
1891        *ymom_boundary_values,
1892        *stage_explicit_update,
1893        *xmom_explicit_update,
1894        *ymom_explicit_update,
1895        *already_computed_flux, //Tracks whether the flux across an edge has already been computed
1896        *max_speed_array; //Keeps track of max speeds for each triangle
1897
1898   
1899    double timestep, epsilon, H0, g;
1900    int optimise_dry_cells;
1901   
1902    // Convert Python arguments to C
1903    if (!PyArg_ParseTuple(args, "dOOOO", &timestep, &domain, &stage, &xmom, &ymom, &bed )) {
1904        PyErr_SetString(PyExc_RuntimeError, "Input arguments failed");
1905        return NULL;
1906    }
1907
1908    epsilon           = get_python_double(domain,"epsilon");
1909    H0                = get_python_double(domain,"H0");
1910    g                 = get_python_double(domain,"g");
1911    optimise_dry_cells = get_python_integer(domain,"optimse_dry_cells");
1912   
1913    neighbours        = get_consecutive_array(domain, "neighbours");
1914    neighbour_edges   = get_consecutive_array(domain, "neighbour_edges"); 
1915    normals           = get_consecutive_array(domain, "normals");
1916    edgelengths       = get_consecutive_array(domain, "edge_lengths");   
1917    radii             = get_consecutive_array(domain, "radii");   
1918    areas             = get_consecutive_array(domain, "areas");   
1919    tri_full_flag     = get_consecutive_array(domain, "tri_full_flag");
1920    already_computed_flux  = get_consecutive_array(domain, "already_computed_flux");
1921    max_speed_array   = get_consecutive_array(domain, "max_speed");
1922   
1923    stage_edge_values = get_consecutive_array(stage, "edge_values");   
1924    xmom_edge_values  = get_consecutive_array(xmom, "edge_values");   
1925    ymom_edge_values  = get_consecutive_array(ymom, "edge_values"); 
1926    bed_edge_values   = get_consecutive_array(bed, "edge_values");   
1927
1928    stage_boundary_values = get_consecutive_array(stage, "boundary_values");   
1929    xmom_boundary_values  = get_consecutive_array(xmom, "boundary_values");   
1930    ymom_boundary_values  = get_consecutive_array(ymom, "boundary_values"); 
1931
1932    stage_explicit_update = get_consecutive_array(stage, "explicit_update");   
1933    xmom_explicit_update  = get_consecutive_array(xmom, "explicit_update");   
1934    ymom_explicit_update  = get_consecutive_array(ymom, "explicit_update"); 
1935
1936
1937  int number_of_elements = stage_edge_values -> dimensions[0];
1938
1939  // Call underlying flux computation routine and update
1940  // the explicit update arrays
1941  timestep = _compute_fluxes_central(number_of_elements,
1942                                     timestep,
1943                                     epsilon,
1944                                     H0,
1945                                     g,
1946                                     (long*) neighbours -> data,
1947                                     (long*) neighbour_edges -> data,
1948                                     (double*) normals -> data,
1949                                     (double*) edgelengths -> data, 
1950                                     (double*) radii -> data, 
1951                                     (double*) areas -> data,
1952                                     (long*) tri_full_flag -> data,
1953                                     (double*) stage_edge_values -> data,
1954                                     (double*) xmom_edge_values -> data,
1955                                     (double*) ymom_edge_values -> data,
1956                                     (double*) bed_edge_values -> data,
1957                                     (double*) stage_boundary_values -> data,
1958                                     (double*) xmom_boundary_values -> data,
1959                                     (double*) ymom_boundary_values -> data,
1960                                     (double*) stage_explicit_update -> data,
1961                                     (double*) xmom_explicit_update -> data,
1962                                     (double*) ymom_explicit_update -> data,
1963                                     (long*) already_computed_flux -> data,
1964                                     (double*) max_speed_array -> data,
1965                                     optimise_dry_cells);
1966
1967  Py_DECREF(neighbours);
1968  Py_DECREF(neighbour_edges);
1969  Py_DECREF(normals);
1970  Py_DECREF(edgelengths);
1971  Py_DECREF(radii);
1972  Py_DECREF(areas);
1973  Py_DECREF(tri_full_flag);
1974  Py_DECREF(already_computed_flux);
1975  Py_DECREF(max_speed_array);
1976  Py_DECREF(stage_edge_values);
1977  Py_DECREF(xmom_edge_values);
1978  Py_DECREF(ymom_edge_values);
1979  Py_DECREF(bed_edge_values);
1980  Py_DECREF(stage_boundary_values);
1981  Py_DECREF(xmom_boundary_values);
1982  Py_DECREF(ymom_boundary_values);
1983  Py_DECREF(stage_explicit_update);
1984  Py_DECREF(xmom_explicit_update);
1985  Py_DECREF(ymom_explicit_update);
1986
1987 
1988  // Return updated flux timestep
1989  return Py_BuildValue("d", timestep);
1990}
1991
1992
1993
1994
1995
1996
1997PyObject *compute_fluxes_ext_central(PyObject *self, PyObject *args) {
1998  /*Compute all fluxes and the timestep suitable for all volumes
1999    in domain.
2000
2001    Compute total flux for each conserved quantity using "flux_function_central"
2002
2003    Fluxes across each edge are scaled by edgelengths and summed up
2004    Resulting flux is then scaled by area and stored in
2005    explicit_update for each of the three conserved quantities
2006    stage, xmomentum and ymomentum
2007
2008    The maximal allowable speed computed by the flux_function for each volume
2009    is converted to a timestep that must not be exceeded. The minimum of
2010    those is computed as the next overall timestep.
2011
2012    Python call:
2013    domain.timestep = compute_fluxes(timestep,
2014                                     domain.epsilon,
2015                                     domain.H0,
2016                                     domain.g,
2017                                     domain.neighbours,
2018                                     domain.neighbour_edges,
2019                                     domain.normals,
2020                                     domain.edgelengths,
2021                                     domain.radii,
2022                                     domain.areas,
2023                                     tri_full_flag,
2024                                     Stage.edge_values,
2025                                     Xmom.edge_values,
2026                                     Ymom.edge_values,
2027                                     Bed.edge_values,
2028                                     Stage.boundary_values,
2029                                     Xmom.boundary_values,
2030                                     Ymom.boundary_values,
2031                                     Stage.explicit_update,
2032                                     Xmom.explicit_update,
2033                                     Ymom.explicit_update,
2034                                     already_computed_flux,
2035                                     optimise_dry_cells)                                     
2036
2037
2038    Post conditions:
2039      domain.explicit_update is reset to computed flux values
2040      domain.timestep is set to the largest step satisfying all volumes.
2041
2042
2043  */
2044
2045
2046  PyArrayObject *neighbours, *neighbour_edges,
2047    *normals, *edgelengths, *radii, *areas,
2048    *tri_full_flag,
2049    *stage_edge_values,
2050    *xmom_edge_values,
2051    *ymom_edge_values,
2052    *bed_edge_values,
2053    *stage_boundary_values,
2054    *xmom_boundary_values,
2055    *ymom_boundary_values,
2056    *stage_explicit_update,
2057    *xmom_explicit_update,
2058    *ymom_explicit_update,
2059    *already_computed_flux, //Tracks whether the flux across an edge has already been computed
2060    *max_speed_array; //Keeps track of max speeds for each triangle
2061
2062   
2063  double timestep, epsilon, H0, g;
2064  int optimise_dry_cells;
2065   
2066  // Convert Python arguments to C
2067  if (!PyArg_ParseTuple(args, "ddddOOOOOOOOOOOOOOOOOOOi",
2068                        &timestep,
2069                        &epsilon,
2070                        &H0,
2071                        &g,
2072                        &neighbours,
2073                        &neighbour_edges,
2074                        &normals,
2075                        &edgelengths, &radii, &areas,
2076                        &tri_full_flag,
2077                        &stage_edge_values,
2078                        &xmom_edge_values,
2079                        &ymom_edge_values,
2080                        &bed_edge_values,
2081                        &stage_boundary_values,
2082                        &xmom_boundary_values,
2083                        &ymom_boundary_values,
2084                        &stage_explicit_update,
2085                        &xmom_explicit_update,
2086                        &ymom_explicit_update,
2087                        &already_computed_flux,
2088                        &max_speed_array,
2089                        &optimise_dry_cells)) {
2090    PyErr_SetString(PyExc_RuntimeError, "Input arguments failed");
2091    return NULL;
2092  }
2093
2094  // check that numpy array objects arrays are C contiguous memory
2095  CHECK_C_CONTIG(neighbours);
2096  CHECK_C_CONTIG(neighbour_edges);
2097  CHECK_C_CONTIG(normals);
2098  CHECK_C_CONTIG(edgelengths);
2099  CHECK_C_CONTIG(radii);
2100  CHECK_C_CONTIG(areas);
2101  CHECK_C_CONTIG(tri_full_flag);
2102  CHECK_C_CONTIG(stage_edge_values);
2103  CHECK_C_CONTIG(xmom_edge_values);
2104  CHECK_C_CONTIG(ymom_edge_values);
2105  CHECK_C_CONTIG(bed_edge_values);
2106  CHECK_C_CONTIG(stage_boundary_values);
2107  CHECK_C_CONTIG(xmom_boundary_values);
2108  CHECK_C_CONTIG(ymom_boundary_values);
2109  CHECK_C_CONTIG(stage_explicit_update);
2110  CHECK_C_CONTIG(xmom_explicit_update);
2111  CHECK_C_CONTIG(ymom_explicit_update);
2112  CHECK_C_CONTIG(already_computed_flux);
2113  CHECK_C_CONTIG(max_speed_array);
2114
2115  int number_of_elements = stage_edge_values -> dimensions[0];
2116
2117  // Call underlying flux computation routine and update
2118  // the explicit update arrays
2119  timestep = _compute_fluxes_central(number_of_elements,
2120                                     timestep,
2121                                     epsilon,
2122                                     H0,
2123                                     g,
2124                                     (long*) neighbours -> data,
2125                                     (long*) neighbour_edges -> data,
2126                                     (double*) normals -> data,
2127                                     (double*) edgelengths -> data, 
2128                                     (double*) radii -> data, 
2129                                     (double*) areas -> data,
2130                                     (long*) tri_full_flag -> data,
2131                                     (double*) stage_edge_values -> data,
2132                                     (double*) xmom_edge_values -> data,
2133                                     (double*) ymom_edge_values -> data,
2134                                     (double*) bed_edge_values -> data,
2135                                     (double*) stage_boundary_values -> data,
2136                                     (double*) xmom_boundary_values -> data,
2137                                     (double*) ymom_boundary_values -> data,
2138                                     (double*) stage_explicit_update -> data,
2139                                     (double*) xmom_explicit_update -> data,
2140                                     (double*) ymom_explicit_update -> data,
2141                                     (long*) already_computed_flux -> data,
2142                                     (double*) max_speed_array -> data,
2143                                     optimise_dry_cells);
2144 
2145  // Return updated flux timestep
2146  return Py_BuildValue("d", timestep);
2147}
2148
2149
2150
2151
2152
2153PyObject *compute_fluxes_ext_kinetic(PyObject *self, PyObject *args) {
2154  /*Compute all fluxes and the timestep suitable for all volumes
2155    in domain.
2156
2157    THIS IS AN EXPERIMENTAL FUNCTION - NORMALLY flux_function_central IS USED.
2158   
2159    Compute total flux for each conserved quantity using "flux_function_kinetic"
2160
2161    Fluxes across each edge are scaled by edgelengths and summed up
2162    Resulting flux is then scaled by area and stored in
2163    explicit_update for each of the three conserved quantities
2164    stage, xmomentum and ymomentum
2165
2166    The maximal allowable speed computed by the flux_function for each volume
2167    is converted to a timestep that must not be exceeded. The minimum of
2168    those is computed as the next overall timestep.
2169
2170    Python call:
2171    domain.timestep = compute_fluxes(timestep,
2172                                     domain.epsilon,
2173                                     domain.H0,
2174                                     domain.g,
2175                                     domain.neighbours,
2176                                     domain.neighbour_edges,
2177                                     domain.normals,
2178                                     domain.edgelengths,
2179                                     domain.radii,
2180                                     domain.areas,
2181                                     Stage.edge_values,
2182                                     Xmom.edge_values,
2183                                     Ymom.edge_values,
2184                                     Bed.edge_values,
2185                                     Stage.boundary_values,
2186                                     Xmom.boundary_values,
2187                                     Ymom.boundary_values,
2188                                     Stage.explicit_update,
2189                                     Xmom.explicit_update,
2190                                     Ymom.explicit_update,
2191                                     already_computed_flux)
2192
2193
2194    Post conditions:
2195      domain.explicit_update is reset to computed flux values
2196      domain.timestep is set to the largest step satisfying all volumes.
2197
2198
2199  */
2200
2201
2202  PyArrayObject *neighbours, *neighbour_edges,
2203    *normals, *edgelengths, *radii, *areas,
2204    *tri_full_flag,
2205    *stage_edge_values,
2206    *xmom_edge_values,
2207    *ymom_edge_values,
2208    *bed_edge_values,
2209    *stage_boundary_values,
2210    *xmom_boundary_values,
2211    *ymom_boundary_values,
2212    *stage_explicit_update,
2213    *xmom_explicit_update,
2214    *ymom_explicit_update,
2215    *already_computed_flux; // Tracks whether the flux across an edge has already been computed
2216
2217
2218  // Local variables
2219  double timestep, max_speed, epsilon, g, H0;
2220  double normal[2], ql[3], qr[3], zl, zr;
2221  double edgeflux[3]; //Work arrays for summing up fluxes
2222
2223  int number_of_elements, k, i, m, n;
2224  int ki, nm=0, ki2; //Index shorthands
2225  static long call=1;
2226
2227
2228  // Convert Python arguments to C
2229  if (!PyArg_ParseTuple(args, "ddddOOOOOOOOOOOOOOOOOO",
2230                        &timestep,
2231                        &epsilon,
2232                        &H0,
2233                        &g,
2234                        &neighbours,
2235                        &neighbour_edges,
2236                        &normals,
2237                        &edgelengths, &radii, &areas,
2238                        &tri_full_flag,
2239                        &stage_edge_values,
2240                        &xmom_edge_values,
2241                        &ymom_edge_values,
2242                        &bed_edge_values,
2243                        &stage_boundary_values,
2244                        &xmom_boundary_values,
2245                        &ymom_boundary_values,
2246                        &stage_explicit_update,
2247                        &xmom_explicit_update,
2248                        &ymom_explicit_update,
2249                        &already_computed_flux)) {
2250    PyErr_SetString(PyExc_RuntimeError, "Input arguments failed");
2251    return NULL;
2252  }
2253  number_of_elements = stage_edge_values -> dimensions[0];
2254  call++;//a static local variable to which already_computed_flux is compared
2255  //set explicit_update to zero for all conserved_quantities.
2256  //This assumes compute_fluxes called before forcing terms
2257  for (k=0; k<number_of_elements; k++) {
2258    ((double *) stage_explicit_update -> data)[k]=0.0;
2259    ((double *) xmom_explicit_update -> data)[k]=0.0;
2260    ((double *) ymom_explicit_update -> data)[k]=0.0;
2261  }
2262  //Loop through neighbours and compute edge flux for each
2263  for (k=0; k<number_of_elements; k++) {
2264    for (i=0; i<3; i++) {
2265      ki = k*3+i;
2266      if (((long *) already_computed_flux->data)[ki]==call)//we've already computed the flux across this edge
2267        continue;
2268      ql[0] = ((double *) stage_edge_values -> data)[ki];
2269      ql[1] = ((double *) xmom_edge_values -> data)[ki];
2270      ql[2] = ((double *) ymom_edge_values -> data)[ki];
2271      zl =    ((double *) bed_edge_values -> data)[ki];
2272
2273      //Quantities at neighbour on nearest face
2274      n = ((long *) neighbours -> data)[ki];
2275      if (n < 0) {
2276        m = -n-1; //Convert negative flag to index
2277        qr[0] = ((double *) stage_boundary_values -> data)[m];
2278        qr[1] = ((double *) xmom_boundary_values -> data)[m];
2279        qr[2] = ((double *) ymom_boundary_values -> data)[m];
2280        zr = zl; //Extend bed elevation to boundary
2281      } else {
2282        m = ((long *) neighbour_edges -> data)[ki];
2283        nm = n*3+m;
2284        qr[0] = ((double *) stage_edge_values -> data)[nm];
2285        qr[1] = ((double *) xmom_edge_values -> data)[nm];
2286        qr[2] = ((double *) ymom_edge_values -> data)[nm];
2287        zr =    ((double *) bed_edge_values -> data)[nm];
2288      }
2289      // Outward pointing normal vector
2290      // normal = domain.normals[k, 2*i:2*i+2]
2291      ki2 = 2*ki; //k*6 + i*2
2292      normal[0] = ((double *) normals -> data)[ki2];
2293      normal[1] = ((double *) normals -> data)[ki2+1];
2294      //Edge flux computation
2295      flux_function_kinetic(ql, qr, zl, zr,
2296                    normal[0], normal[1],
2297                    epsilon, H0, g,
2298                    edgeflux, &max_speed);
2299      //update triangle k
2300      ((long *) already_computed_flux->data)[ki]=call;
2301      ((double *) stage_explicit_update -> data)[k] -= edgeflux[0]*((double *) edgelengths -> data)[ki];
2302      ((double *) xmom_explicit_update -> data)[k] -= edgeflux[1]*((double *) edgelengths -> data)[ki];
2303      ((double *) ymom_explicit_update -> data)[k] -= edgeflux[2]*((double *) edgelengths -> data)[ki];
2304      //update the neighbour n
2305      if (n>=0){
2306        ((long *) already_computed_flux->data)[nm]=call;
2307        ((double *) stage_explicit_update -> data)[n] += edgeflux[0]*((double *) edgelengths -> data)[nm];
2308        ((double *) xmom_explicit_update -> data)[n] += edgeflux[1]*((double *) edgelengths -> data)[nm];
2309        ((double *) ymom_explicit_update -> data)[n] += edgeflux[2]*((double *) edgelengths -> data)[nm];
2310      }
2311      ///for (j=0; j<3; j++) {
2312        ///flux[j] -= edgeflux[j]*((double *) edgelengths -> data)[ki];
2313        ///}
2314        //Update timestep
2315        //timestep = min(timestep, domain.radii[k]/max_speed)
2316        //FIXME: SR Add parameter for CFL condition
2317    if ( ((long *) tri_full_flag -> data)[k] == 1) {
2318            if (max_speed > epsilon) {
2319                timestep = min(timestep, ((double *) radii -> data)[k]/max_speed);
2320                //maxspeed in flux_function is calculated as max(|u+a|,|u-a|)
2321                if (n>=0)
2322                    timestep = min(timestep, ((double *) radii -> data)[n]/max_speed);
2323            }
2324    }
2325    } // end for i
2326    //Normalise by area and store for when all conserved
2327    //quantities get updated
2328    ((double *) stage_explicit_update -> data)[k] /= ((double *) areas -> data)[k];
2329    ((double *) xmom_explicit_update -> data)[k] /= ((double *) areas -> data)[k];
2330    ((double *) ymom_explicit_update -> data)[k] /= ((double *) areas -> data)[k];
2331  } //end for k
2332  return Py_BuildValue("d", timestep);
2333}
2334
2335PyObject *protect(PyObject *self, PyObject *args) {
2336  //
2337  //    protect(minimum_allowed_height, maximum_allowed_speed, wc, zc, xmomc, ymomc)
2338
2339
2340  PyArrayObject
2341  *wc,            //Stage at centroids
2342  *zc,            //Elevation at centroids
2343  *xmomc,         //Momentums at centroids
2344  *ymomc;
2345
2346
2347  int N;
2348  double minimum_allowed_height, maximum_allowed_speed, epsilon;
2349
2350  // Convert Python arguments to C
2351  if (!PyArg_ParseTuple(args, "dddOOOO",
2352                        &minimum_allowed_height,
2353                        &maximum_allowed_speed,
2354                        &epsilon,
2355                        &wc, &zc, &xmomc, &ymomc)) {
2356    PyErr_SetString(PyExc_RuntimeError, "shallow_water_ext.c: protect could not parse input arguments");
2357    return NULL;
2358  } 
2359
2360  N = wc -> dimensions[0];
2361
2362  _protect(N,
2363           minimum_allowed_height,
2364           maximum_allowed_speed,
2365           epsilon,
2366           (double*) wc -> data,
2367           (double*) zc -> data,
2368           (double*) xmomc -> data,
2369           (double*) ymomc -> data);
2370
2371  return Py_BuildValue("");
2372}
2373
2374
2375
2376PyObject *balance_deep_and_shallow(PyObject *self, PyObject *args) {
2377  // Compute linear combination between stage as computed by
2378  // gradient-limiters limiting using w, and stage computed by
2379  // gradient-limiters limiting using h (h-limiter).
2380  // The former takes precedence when heights are large compared to the
2381  // bed slope while the latter takes precedence when heights are
2382  // relatively small.  Anything in between is computed as a balanced
2383  // linear combination in order to avoid numerical disturbances which
2384  // would otherwise appear as a result of hard switching between
2385  // modes.
2386  //
2387  //    balance_deep_and_shallow(wc, zc, wv, zv,
2388  //                             xmomc, ymomc, xmomv, ymomv)
2389
2390
2391  PyArrayObject
2392    *wc,            //Stage at centroids
2393    *zc,            //Elevation at centroids
2394    *wv,            //Stage at vertices
2395    *zv,            //Elevation at vertices
2396    *hvbar,         //h-Limited depths at vertices
2397    *xmomc,         //Momentums at centroids and vertices
2398    *ymomc,
2399    *xmomv,
2400    *ymomv;
2401 
2402  PyObject *domain, *Tmp;
2403   
2404  double alpha_balance = 2.0;
2405  double H0;
2406
2407  int N, tight_slope_limiters, use_centroid_velocities; //, err;
2408
2409  // Convert Python arguments to C
2410  if (!PyArg_ParseTuple(args, "OOOOOOOOOO",
2411                        &domain,
2412                        &wc, &zc, 
2413                        &wv, &zv, &hvbar,
2414                        &xmomc, &ymomc, &xmomv, &ymomv)) {
2415    PyErr_SetString(PyExc_RuntimeError, 
2416                    "shallow_water_ext.c: balance_deep_and_shallow could not parse input arguments");
2417    return NULL;
2418  } 
2419         
2420         
2421  // FIXME (Ole): I tested this without GetAttrString and got time down
2422  // marginally from 4.0s to 3.8s. Consider passing everything in
2423  // through ParseTuple and profile.
2424 
2425  // Pull out parameters
2426  Tmp = PyObject_GetAttrString(domain, "alpha_balance");
2427  if (!Tmp) {
2428    PyErr_SetString(PyExc_RuntimeError, "shallow_water_ext.c: balance_deep_and_shallow could not obtain object alpha_balance from domain");
2429    return NULL;
2430  } 
2431  alpha_balance = PyFloat_AsDouble(Tmp);
2432  Py_DECREF(Tmp);
2433
2434 
2435  Tmp = PyObject_GetAttrString(domain, "H0");
2436  if (!Tmp) {
2437    PyErr_SetString(PyExc_RuntimeError, "shallow_water_ext.c: balance_deep_and_shallow could not obtain object H0 from domain");
2438    return NULL;
2439  } 
2440  H0 = PyFloat_AsDouble(Tmp);
2441  Py_DECREF(Tmp);
2442
2443 
2444  Tmp = PyObject_GetAttrString(domain, "tight_slope_limiters");
2445  if (!Tmp) {
2446    PyErr_SetString(PyExc_RuntimeError, "shallow_water_ext.c: balance_deep_and_shallow could not obtain object tight_slope_limiters from domain");
2447    return NULL;
2448  } 
2449  tight_slope_limiters = PyInt_AsLong(Tmp);
2450  Py_DECREF(Tmp);
2451 
2452 
2453  Tmp = PyObject_GetAttrString(domain, "use_centroid_velocities");
2454  if (!Tmp) {
2455    PyErr_SetString(PyExc_RuntimeError, "shallow_water_ext.c: balance_deep_and_shallow could not obtain object use_centroid_velocities from domain");
2456    return NULL;
2457  } 
2458  use_centroid_velocities = PyInt_AsLong(Tmp);
2459  Py_DECREF(Tmp);
2460 
2461 
2462 
2463  N = wc -> dimensions[0];
2464  _balance_deep_and_shallow(N,
2465                            (double*) wc -> data,
2466                            (double*) zc -> data,
2467                            (double*) wv -> data,
2468                            (double*) zv -> data,
2469                            (double*) hvbar -> data,
2470                            (double*) xmomc -> data,
2471                            (double*) ymomc -> data,
2472                            (double*) xmomv -> data,
2473                            (double*) ymomv -> data,
2474                            H0,
2475                            (int) tight_slope_limiters,
2476                            (int) use_centroid_velocities,                         
2477                            alpha_balance);
2478
2479
2480  return Py_BuildValue("");
2481}
2482
2483
2484
2485
2486PyObject *assign_windfield_values(PyObject *self, PyObject *args) {
2487  //
2488  //      assign_windfield_values(xmom_update, ymom_update,
2489  //                              s_vec, phi_vec, self.const)
2490
2491
2492
2493  PyArrayObject   //(one element per triangle)
2494  *s_vec,         //Speeds
2495  *phi_vec,       //Bearings
2496  *xmom_update,   //Momentum updates
2497  *ymom_update;
2498
2499
2500  int N;
2501  double cw;
2502
2503  // Convert Python arguments to C
2504  if (!PyArg_ParseTuple(args, "OOOOd",
2505                        &xmom_update,
2506                        &ymom_update,
2507                        &s_vec, &phi_vec,
2508                        &cw)) {
2509    PyErr_SetString(PyExc_RuntimeError, "shallow_water_ext.c: assign_windfield_values could not parse input arguments");
2510    return NULL;
2511  } 
2512                       
2513
2514  N = xmom_update -> dimensions[0];
2515
2516  _assign_wind_field_values(N,
2517           (double*) xmom_update -> data,
2518           (double*) ymom_update -> data,
2519           (double*) s_vec -> data,
2520           (double*) phi_vec -> data,
2521           cw);
2522
2523  return Py_BuildValue("");
2524}
2525
2526
2527
2528
2529//-------------------------------
2530// Method table for python module
2531//-------------------------------
2532static struct PyMethodDef MethodTable[] = {
2533  /* The cast of the function is necessary since PyCFunction values
2534   * only take two PyObject* parameters, and rotate() takes
2535   * three.
2536   */
2537
2538  {"rotate", (PyCFunction)rotate, METH_VARARGS | METH_KEYWORDS, "Print out"},
2539  {"extrapolate_second_order_sw", extrapolate_second_order_sw, METH_VARARGS, "Print out"},
2540  {"compute_fluxes_ext_central", compute_fluxes_ext_central, METH_VARARGS, "Print out"},
2541  {"compute_fluxes_ext_central_new", compute_fluxes_ext_central_new, METH_VARARGS, "Print out"},
2542  {"compute_fluxes_ext_kinetic", compute_fluxes_ext_kinetic, METH_VARARGS, "Print out"},
2543  {"gravity", gravity, METH_VARARGS, "Print out"},
2544  {"manning_friction", manning_friction, METH_VARARGS, "Print out"},
2545  {"flux_function_central", flux_function_central, METH_VARARGS, "Print out"}, 
2546  {"balance_deep_and_shallow", balance_deep_and_shallow,
2547   METH_VARARGS, "Print out"},
2548  {"protect", protect, METH_VARARGS | METH_KEYWORDS, "Print out"},
2549  {"assign_windfield_values", assign_windfield_values,
2550   METH_VARARGS | METH_KEYWORDS, "Print out"},
2551  {NULL, NULL}
2552};
2553
2554// Module initialisation
2555void initshallow_water_ext(void){
2556  Py_InitModule("shallow_water_ext", MethodTable);
2557
2558  import_array(); // Necessary for handling of NumPY structures
2559}
Note: See TracBrowser for help on using the repository browser.