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