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

Last change on this file was 7164, checked in by ole, 16 years ago

Ported optimisation of shallow_water_ext.c
from changeset:7143 into numpy branch.

Changes to test_shallow_water.py and test_data_manager.py
implemented in this changeset still
need to be ported.

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