[229] | 1 | """Class Domain - |
---|
| 2 | 2D triangular domains for finite-volume computations of |
---|
| 3 | the shallow water wave equation. |
---|
| 4 | |
---|
| 5 | This module contains a specialisation of class Domain from module domain.py |
---|
| 6 | consisting of methods specific to the Shallow Water Wave Equation |
---|
| 7 | |
---|
| 8 | FIXME: Write equations here! |
---|
| 9 | |
---|
| 10 | |
---|
| 11 | Conserved quantities are w (water level or stage), uh (x momentum) |
---|
| 12 | and vh (y momentum). |
---|
| 13 | |
---|
| 14 | |
---|
| 15 | Ole Nielsen, Stephen Roberts, Duncan Gray, Christopher Zoppou |
---|
| 16 | Geoscience Australia, 2004 |
---|
| 17 | """ |
---|
| 18 | |
---|
| 19 | from domain import * |
---|
| 20 | Generic_domain = Domain #Rename |
---|
| 21 | |
---|
| 22 | class Domain(Generic_domain): |
---|
| 23 | |
---|
[415] | 24 | def __init__(self, coordinates, vertices, boundary = None, |
---|
| 25 | tagged_elements = None): |
---|
[229] | 26 | |
---|
| 27 | conserved_quantities = ['level', 'xmomentum', 'ymomentum'] |
---|
| 28 | other_quantities = ['elevation', 'friction'] |
---|
| 29 | |
---|
| 30 | Generic_domain.__init__(self, coordinates, vertices, boundary, |
---|
[415] | 31 | conserved_quantities, other_quantities, |
---|
| 32 | tagged_elements) |
---|
[229] | 33 | |
---|
[240] | 34 | from config import minimum_allowed_height, g |
---|
[229] | 35 | self.minimum_allowed_height = minimum_allowed_height |
---|
[240] | 36 | self.g = g |
---|
[229] | 37 | |
---|
| 38 | self.forcing_terms.append(gravity) |
---|
| 39 | self.forcing_terms.append(manning_friction) |
---|
| 40 | |
---|
[269] | 41 | #Realtime visualisation |
---|
| 42 | self.visualise = False |
---|
| 43 | |
---|
| 44 | #Stored output |
---|
| 45 | self.store=False |
---|
[281] | 46 | self.format = 'sww' |
---|
[269] | 47 | self.smooth = True |
---|
| 48 | |
---|
| 49 | #Reduction operation for get_vertex_values |
---|
| 50 | #from pytools.stats import mean |
---|
| 51 | #self.reduction = mean |
---|
| 52 | self.reduction = min #Looks better near steep slopes |
---|
| 53 | |
---|
| 54 | |
---|
[268] | 55 | #Establish shortcuts to relevant quantities (for efficiency) |
---|
[462] | 56 | #self.w = self.quantities['level'] |
---|
| 57 | #self.uh = self.quantities['xmomentum'] |
---|
| 58 | #self.vh = self.quantities['ymomentum'] |
---|
| 59 | #self.z = self.quantities['elevation'] |
---|
| 60 | #self.eta = self.quantities['friction'] |
---|
[272] | 61 | |
---|
[229] | 62 | def check_integrity(self): |
---|
| 63 | Generic_domain.check_integrity(self) |
---|
| 64 | |
---|
| 65 | #Check that we are solving the shallow water wave equation |
---|
| 66 | |
---|
| 67 | msg = 'First conserved quantity must be "level"' |
---|
| 68 | assert self.conserved_quantities[0] == 'level', msg |
---|
| 69 | msg = 'Second conserved quantity must be "xmomentum"' |
---|
| 70 | assert self.conserved_quantities[1] == 'xmomentum', msg |
---|
| 71 | msg = 'Third conserved quantity must be "ymomentum"' |
---|
| 72 | assert self.conserved_quantities[2] == 'ymomentum', msg |
---|
| 73 | |
---|
| 74 | |
---|
| 75 | #Check that levels are >= bed elevation |
---|
| 76 | from Numeric import alltrue, greater_equal |
---|
| 77 | |
---|
| 78 | level = self.quantities['level'] |
---|
| 79 | bed = self.quantities['elevation'] |
---|
| 80 | |
---|
| 81 | msg = 'All water levels must be greater than the bed elevation' |
---|
| 82 | assert alltrue( greater_equal( |
---|
| 83 | level.vertex_values, bed.vertex_values )), msg |
---|
| 84 | |
---|
| 85 | assert alltrue( greater_equal( |
---|
| 86 | level.edge_values, bed.edge_values )), msg |
---|
| 87 | |
---|
| 88 | assert alltrue( greater_equal( |
---|
| 89 | level.centroid_values, bed.centroid_values )), msg |
---|
| 90 | |
---|
| 91 | |
---|
| 92 | def compute_fluxes(self): |
---|
| 93 | #Call correct module function |
---|
| 94 | #(either from this module or C-extension) |
---|
| 95 | compute_fluxes(self) |
---|
| 96 | |
---|
| 97 | def distribute_to_vertices_and_edges(self): |
---|
| 98 | #Call correct module function |
---|
| 99 | #(either from this module or C-extension) |
---|
| 100 | distribute_to_vertices_and_edges(self) |
---|
[269] | 101 | |
---|
| 102 | |
---|
[462] | 103 | #FIXME: Under construction |
---|
| 104 | # def set_defaults(self): |
---|
| 105 | # """Set default values for uninitialised quantities. |
---|
| 106 | # This is specific to the shallow water wave equation |
---|
| 107 | # Defaults for 'elevation', 'friction', 'xmomentum' and 'ymomentum' |
---|
| 108 | # are 0.0. Default for 'level' is whatever the value of 'elevation'. |
---|
| 109 | # """ |
---|
| 110 | |
---|
| 111 | # for name in self.other_quantities + self.conserved_quantities: |
---|
| 112 | # print name |
---|
| 113 | # print self.quantities.keys() |
---|
| 114 | # if not self.quantities.has_key(name): |
---|
| 115 | # if name == 'level': |
---|
| 116 | |
---|
| 117 | # if self.quantities.has_key('elevation'): |
---|
| 118 | # z = self.quantities['elevation'].vertex_values |
---|
| 119 | # self.set_quantity(name, z) |
---|
| 120 | # else: |
---|
| 121 | # self.set_quantity(name, 0.0) |
---|
| 122 | # else: |
---|
| 123 | # self.set_quantity(name, 0.0) |
---|
| 124 | |
---|
| 125 | |
---|
| 126 | |
---|
| 127 | # #Lift negative heights up |
---|
| 128 | # #z = self.quantities['elevation'].vertex_values |
---|
| 129 | # #w = self.quantities['level'].vertex_values |
---|
| 130 | |
---|
| 131 | # #h = w-z |
---|
| 132 | |
---|
| 133 | # #for k in range(h.shape[0]): |
---|
| 134 | # # for i in range(3): |
---|
| 135 | # # if h[k, i] < 0.0: |
---|
| 136 | # # w[k, i] = z[k, i] |
---|
| 137 | |
---|
| 138 | |
---|
| 139 | # #self.quantities['level'].interpolate() |
---|
| 140 | |
---|
| 141 | |
---|
| 142 | |
---|
[269] | 143 | def evolve(self, yieldstep = None, finaltime = None): |
---|
[271] | 144 | """Specialisation of basic evolve method from parent class |
---|
| 145 | """ |
---|
| 146 | |
---|
| 147 | #Initialise real time viz if requested |
---|
[269] | 148 | if self.visualise is True and self.time == 0.0: |
---|
| 149 | import realtime_visualisation as visualise |
---|
| 150 | visualise.create_surface(self) |
---|
| 151 | |
---|
| 152 | #Store model data, e.g. for visualisation |
---|
| 153 | if self.store is True and self.time == 0.0: |
---|
[281] | 154 | self.initialise_storage() |
---|
| 155 | |
---|
[271] | 156 | |
---|
| 157 | #Call basic machinery from parent class |
---|
[269] | 158 | for t in Generic_domain.evolve(self, yieldstep, finaltime): |
---|
| 159 | #Real time viz |
---|
| 160 | if self.visualise is True: |
---|
| 161 | visualise.update(self) |
---|
| 162 | |
---|
| 163 | #Store model data, e.g. for subsequent visualisation |
---|
| 164 | if self.store is True: |
---|
[287] | 165 | self.store_timestep('level') |
---|
[281] | 166 | #FIXME: Could maybe be taken from specified list |
---|
| 167 | #of 'store every step' quantities |
---|
[269] | 168 | |
---|
[280] | 169 | |
---|
[269] | 170 | #Pass control on to outer loop for more specific actions |
---|
| 171 | yield(t) |
---|
[229] | 172 | |
---|
[240] | 173 | |
---|
[281] | 174 | def initialise_storage(self): |
---|
[287] | 175 | """Create and initialise self.writer object for storing data. |
---|
| 176 | Also, save x,y and bed elevation |
---|
[281] | 177 | """ |
---|
[240] | 178 | |
---|
[281] | 179 | import data_manager |
---|
| 180 | |
---|
| 181 | #Initialise writer |
---|
| 182 | self.writer = data_manager.get_dataobject(self, mode = 'w') |
---|
| 183 | |
---|
| 184 | #Store vertices and connectivity |
---|
| 185 | self.writer.store_connectivity() |
---|
| 186 | |
---|
[287] | 187 | def store_timestep(self, name): |
---|
| 188 | """Store named quantity and time. |
---|
| 189 | |
---|
| 190 | Precondition: |
---|
| 191 | self.write has been initialised |
---|
[280] | 192 | """ |
---|
[287] | 193 | self.writer.store_timestep(name) |
---|
[281] | 194 | |
---|
[240] | 195 | #Rotation of momentum vector |
---|
| 196 | def rotate(q, normal, direction = 1): |
---|
| 197 | """Rotate the momentum component q (q[1], q[2]) |
---|
| 198 | from x,y coordinates to coordinates based on normal vector. |
---|
| 199 | |
---|
| 200 | If direction is negative the rotation is inverted. |
---|
| 201 | |
---|
| 202 | Input vector is preserved |
---|
| 203 | |
---|
| 204 | This function is specific to the shallow water wave equation |
---|
| 205 | """ |
---|
| 206 | |
---|
| 207 | #FIXME: Needs to be tested |
---|
| 208 | |
---|
| 209 | from Numeric import zeros, Float |
---|
| 210 | |
---|
| 211 | assert len(q) == 3,\ |
---|
| 212 | 'Vector of conserved quantities must have length 3'\ |
---|
| 213 | 'for 2D shallow water equation' |
---|
| 214 | |
---|
| 215 | try: |
---|
| 216 | l = len(normal) |
---|
| 217 | except: |
---|
| 218 | raise 'Normal vector must be an Numeric array' |
---|
| 219 | |
---|
| 220 | #FIXME: Put this test into C-extension as well |
---|
| 221 | assert l == 2, 'Normal vector must have 2 components' |
---|
| 222 | |
---|
| 223 | |
---|
| 224 | n1 = normal[0] |
---|
| 225 | n2 = normal[1] |
---|
| 226 | |
---|
| 227 | r = zeros(len(q), Float) #Rotated quantities |
---|
| 228 | r[0] = q[0] #First quantity, height, is not rotated |
---|
| 229 | |
---|
| 230 | if direction == -1: |
---|
| 231 | n2 = -n2 |
---|
| 232 | |
---|
| 233 | |
---|
| 234 | r[1] = n1*q[1] + n2*q[2] |
---|
| 235 | r[2] = -n2*q[1] + n1*q[2] |
---|
| 236 | |
---|
| 237 | return r |
---|
| 238 | |
---|
| 239 | |
---|
| 240 | |
---|
[229] | 241 | #################################### |
---|
| 242 | # Flux computation |
---|
| 243 | def flux_function(normal, ql, qr, zl, zr): |
---|
| 244 | """Compute fluxes between volumes for the shallow water wave equation |
---|
| 245 | cast in terms of w = h+z using the 'central scheme' as described in |
---|
| 246 | |
---|
| 247 | Kurganov, Noelle, Petrova. 'Semidiscrete Central-Upwind Schemes For |
---|
| 248 | Hyperbolic Conservation Laws and Hamilton-Jacobi Equations'. |
---|
| 249 | Siam J. Sci. Comput. Vol. 23, No. 3, pp. 707-740. |
---|
| 250 | |
---|
| 251 | The implemented formula is given in equation (3.15) on page 714 |
---|
| 252 | |
---|
| 253 | Conserved quantities w, uh, vh are stored as elements 0, 1 and 2 |
---|
| 254 | in the numerical vectors ql an qr. |
---|
| 255 | |
---|
| 256 | Bed elevations zl and zr. |
---|
| 257 | """ |
---|
| 258 | |
---|
[232] | 259 | from config import g, epsilon |
---|
[229] | 260 | from math import sqrt |
---|
| 261 | from Numeric import array |
---|
| 262 | |
---|
| 263 | #Align momentums with x-axis |
---|
| 264 | q_left = rotate(ql, normal, direction = 1) |
---|
| 265 | q_right = rotate(qr, normal, direction = 1) |
---|
| 266 | |
---|
| 267 | z = (zl+zr)/2 #Take average of field values |
---|
| 268 | |
---|
| 269 | w_left = q_left[0] #w=h+z |
---|
| 270 | h_left = w_left-z |
---|
| 271 | uh_left = q_left[1] |
---|
| 272 | |
---|
[232] | 273 | if h_left < epsilon: |
---|
| 274 | u_left = 0.0 #Could have been negative |
---|
| 275 | h_left = 0.0 |
---|
| 276 | else: |
---|
[229] | 277 | u_left = uh_left/h_left |
---|
| 278 | |
---|
| 279 | |
---|
| 280 | w_right = q_right[0] #w=h+z |
---|
| 281 | h_right = w_right-z |
---|
| 282 | uh_right = q_right[1] |
---|
| 283 | |
---|
| 284 | |
---|
[232] | 285 | if h_right < epsilon: |
---|
| 286 | u_right = 0.0 #Could have been negative |
---|
| 287 | h_right = 0.0 |
---|
| 288 | else: |
---|
[229] | 289 | u_right = uh_right/h_right |
---|
| 290 | |
---|
| 291 | vh_left = q_left[2] |
---|
| 292 | vh_right = q_right[2] |
---|
| 293 | |
---|
| 294 | soundspeed_left = sqrt(g*h_left) |
---|
| 295 | soundspeed_right = sqrt(g*h_right) |
---|
| 296 | |
---|
| 297 | #Maximal wave speed |
---|
| 298 | s_max = max(u_left+soundspeed_left, u_right+soundspeed_right, 0) |
---|
| 299 | |
---|
| 300 | #Minimal wave speed |
---|
| 301 | s_min = min(u_left-soundspeed_left, u_right-soundspeed_right, 0) |
---|
| 302 | |
---|
| 303 | #Flux computation |
---|
[232] | 304 | flux_left = array([u_left*h_left, |
---|
| 305 | u_left*uh_left + 0.5*g*h_left**2, |
---|
| 306 | u_left*vh_left]) |
---|
| 307 | flux_right = array([u_right*h_right, |
---|
| 308 | u_right*uh_right + 0.5*g*h_right**2, |
---|
| 309 | u_right*vh_right]) |
---|
[229] | 310 | |
---|
| 311 | denom = s_max-s_min |
---|
| 312 | if denom == 0.0: |
---|
[232] | 313 | edgeflux = array([0.0, 0.0, 0.0]) |
---|
[229] | 314 | max_speed = 0.0 |
---|
| 315 | else: |
---|
[232] | 316 | edgeflux = (s_max*flux_left - s_min*flux_right)/denom |
---|
| 317 | edgeflux += s_max*s_min*(q_right-q_left)/denom |
---|
[229] | 318 | |
---|
[232] | 319 | edgeflux = rotate(edgeflux, normal, direction=-1) |
---|
[229] | 320 | max_speed = max(abs(s_max), abs(s_min)) |
---|
| 321 | |
---|
[232] | 322 | return edgeflux, max_speed |
---|
[229] | 323 | |
---|
| 324 | |
---|
| 325 | def compute_fluxes(domain): |
---|
| 326 | """Compute all fluxes and the timestep suitable for all volumes |
---|
| 327 | in domain. |
---|
| 328 | |
---|
| 329 | Compute total flux for each conserved quantity using "flux_function" |
---|
| 330 | |
---|
| 331 | Fluxes across each edge are scaled by edgelengths and summed up |
---|
| 332 | Resulting flux is then scaled by area and stored in |
---|
[240] | 333 | explicit_update for each of the three conserved quantities |
---|
| 334 | level, xmomentum and ymomentum |
---|
[229] | 335 | |
---|
| 336 | The maximal allowable speed computed by the flux_function for each volume |
---|
| 337 | is converted to a timestep that must not be exceeded. The minimum of |
---|
| 338 | those is computed as the next overall timestep. |
---|
| 339 | |
---|
| 340 | Post conditions: |
---|
| 341 | domain.explicit_update is reset to computed flux values |
---|
| 342 | domain.timestep is set to the largest step satisfying all volumes. |
---|
| 343 | """ |
---|
| 344 | |
---|
| 345 | import sys |
---|
| 346 | from Numeric import zeros, Float |
---|
| 347 | |
---|
| 348 | N = domain.number_of_elements |
---|
| 349 | |
---|
| 350 | #Shortcuts |
---|
| 351 | Level = domain.quantities['level'] |
---|
| 352 | Xmom = domain.quantities['xmomentum'] |
---|
| 353 | Ymom = domain.quantities['ymomentum'] |
---|
| 354 | Bed = domain.quantities['elevation'] |
---|
| 355 | |
---|
| 356 | #Arrays |
---|
| 357 | level = Level.edge_values |
---|
| 358 | xmom = Xmom.edge_values |
---|
| 359 | ymom = Ymom.edge_values |
---|
| 360 | bed = Bed.edge_values |
---|
| 361 | |
---|
| 362 | level_bdry = Level.boundary_values |
---|
| 363 | xmom_bdry = Xmom.boundary_values |
---|
| 364 | ymom_bdry = Ymom.boundary_values |
---|
| 365 | |
---|
| 366 | flux = zeros(3, Float) #Work array for summing up fluxes |
---|
| 367 | |
---|
| 368 | #Loop |
---|
[240] | 369 | timestep = float(sys.maxint) |
---|
[229] | 370 | for k in range(N): |
---|
| 371 | |
---|
| 372 | flux[:] = 0. #Reset work array |
---|
| 373 | for i in range(3): |
---|
| 374 | #Quantities inside volume facing neighbour i |
---|
| 375 | ql = [level[k, i], xmom[k, i], ymom[k, i]] |
---|
| 376 | zl = bed[k, i] |
---|
| 377 | |
---|
| 378 | #Quantities at neighbour on nearest face |
---|
[240] | 379 | n = domain.neighbours[k,i] |
---|
[229] | 380 | if n < 0: |
---|
[240] | 381 | m = -n-1 #Convert negative flag to index |
---|
[229] | 382 | qr = [level_bdry[m], xmom_bdry[m], ymom_bdry[m]] |
---|
| 383 | zr = zl #Extend bed elevation to boundary |
---|
| 384 | else: |
---|
[240] | 385 | m = domain.neighbour_edges[k,i] |
---|
[229] | 386 | qr = [level[n, m], xmom[n, m], ymom[n, m]] |
---|
| 387 | zr = bed[n, m] |
---|
| 388 | |
---|
| 389 | |
---|
| 390 | #Outward pointing normal vector |
---|
[240] | 391 | normal = domain.normals[k, 2*i:2*i+2] |
---|
[229] | 392 | |
---|
| 393 | #Flux computation using provided function |
---|
| 394 | edgeflux, max_speed = flux_function(normal, ql, qr, zl, zr) |
---|
[240] | 395 | flux -= edgeflux * domain.edgelengths[k,i] |
---|
| 396 | |
---|
[229] | 397 | #Update optimal_timestep |
---|
| 398 | try: |
---|
[240] | 399 | timestep = min(timestep, domain.radii[k]/max_speed) |
---|
[229] | 400 | except ZeroDivisionError: |
---|
| 401 | pass |
---|
| 402 | |
---|
| 403 | #Normalise by area and store for when all conserved |
---|
| 404 | #quantities get updated |
---|
[240] | 405 | flux /= domain.areas[k] |
---|
[229] | 406 | Level.explicit_update[k] = flux[0] |
---|
| 407 | Xmom.explicit_update[k] = flux[1] |
---|
| 408 | Ymom.explicit_update[k] = flux[2] |
---|
| 409 | |
---|
[458] | 410 | #print 'FLUX l', Level.explicit_update |
---|
| 411 | #print 'FLUX x', Xmom.explicit_update |
---|
| 412 | #print 'FLUX y', Ymom.explicit_update |
---|
| 413 | |
---|
[229] | 414 | domain.timestep = timestep |
---|
| 415 | |
---|
| 416 | |
---|
[240] | 417 | def compute_fluxes_c(domain): |
---|
[246] | 418 | """Wrapper calling C version of compute fluxes |
---|
[240] | 419 | """ |
---|
| 420 | |
---|
| 421 | import sys |
---|
| 422 | from Numeric import zeros, Float |
---|
| 423 | |
---|
| 424 | N = domain.number_of_elements |
---|
| 425 | |
---|
| 426 | #Shortcuts |
---|
| 427 | Level = domain.quantities['level'] |
---|
| 428 | Xmom = domain.quantities['xmomentum'] |
---|
| 429 | Ymom = domain.quantities['ymomentum'] |
---|
| 430 | Bed = domain.quantities['elevation'] |
---|
| 431 | |
---|
| 432 | timestep = float(sys.maxint) |
---|
| 433 | from shallow_water_ext import compute_fluxes |
---|
| 434 | domain.timestep = compute_fluxes(timestep, domain.epsilon, domain.g, |
---|
| 435 | domain.neighbours, |
---|
| 436 | domain.neighbour_edges, |
---|
| 437 | domain.normals, |
---|
| 438 | domain.edgelengths, |
---|
| 439 | domain.radii, |
---|
| 440 | domain.areas, |
---|
| 441 | Level.edge_values, |
---|
| 442 | Xmom.edge_values, |
---|
| 443 | Ymom.edge_values, |
---|
| 444 | Bed.edge_values, |
---|
| 445 | Level.boundary_values, |
---|
| 446 | Xmom.boundary_values, |
---|
| 447 | Ymom.boundary_values, |
---|
| 448 | Level.explicit_update, |
---|
| 449 | Xmom.explicit_update, |
---|
| 450 | Ymom.explicit_update) |
---|
[268] | 451 | |
---|
[240] | 452 | |
---|
[229] | 453 | #################################### |
---|
| 454 | # Module functions for gradient limiting (distribute_to_vertices_and_edges) |
---|
| 455 | |
---|
| 456 | def distribute_to_vertices_and_edges(domain): |
---|
[274] | 457 | """Distribution from centroids to vertices specific to the |
---|
| 458 | shallow water wave |
---|
| 459 | equation. |
---|
[229] | 460 | |
---|
[274] | 461 | It will ensure that h (w-z) is always non-negative even in the |
---|
| 462 | presence of steep bed-slopes by taking a weighted average between shallow |
---|
| 463 | and deep cases. |
---|
[234] | 464 | |
---|
[274] | 465 | In addition, all conserved quantities get distributed as per either a |
---|
| 466 | constant (order==1) or a piecewise linear function (order==2). |
---|
| 467 | |
---|
| 468 | FIXME: more explanation about removal of artificial variability etc |
---|
| 469 | |
---|
| 470 | Precondition: |
---|
| 471 | All quantities defined at centroids and bed elevation defined at |
---|
| 472 | vertices. |
---|
| 473 | |
---|
| 474 | Postcondition |
---|
| 475 | Conserved quantities defined at vertices |
---|
| 476 | |
---|
| 477 | """ |
---|
| 478 | |
---|
| 479 | #Remove very thin layers of water |
---|
[443] | 480 | protect_against_infinitesimal_and_negative_heights(domain) |
---|
[274] | 481 | |
---|
| 482 | #Extrapolate all conserved quantities |
---|
[229] | 483 | for name in domain.conserved_quantities: |
---|
| 484 | Q = domain.quantities[name] |
---|
[274] | 485 | if domain.order == 1: |
---|
| 486 | Q.extrapolate_first_order() |
---|
| 487 | elif domain.order == 2: |
---|
| 488 | Q.extrapolate_second_order() |
---|
| 489 | Q.limit() |
---|
| 490 | else: |
---|
| 491 | raise 'Unknown order' |
---|
| 492 | |
---|
[443] | 493 | #Take bed elevation into account when water heights are small |
---|
[274] | 494 | balance_deep_and_shallow(domain) |
---|
[443] | 495 | |
---|
[274] | 496 | #Compute edge values by interpolation |
---|
| 497 | for name in domain.conserved_quantities: |
---|
| 498 | Q = domain.quantities[name] |
---|
[268] | 499 | Q.interpolate_from_vertices_to_edges() |
---|
| 500 | |
---|
[443] | 501 | |
---|
| 502 | |
---|
| 503 | def dry(domain): |
---|
| 504 | """Protect against infinitesimal heights and associated high velocities |
---|
| 505 | at vertices |
---|
| 506 | """ |
---|
| 507 | |
---|
| 508 | #FIXME: Experimental |
---|
[268] | 509 | |
---|
[443] | 510 | #Shortcuts |
---|
| 511 | wv = domain.quantities['level'].vertex_values |
---|
| 512 | zv = domain.quantities['elevation'].vertex_values |
---|
| 513 | xmomv = domain.quantities['xmomentum'].vertex_values |
---|
| 514 | ymomv = domain.quantities['ymomentum'].vertex_values |
---|
| 515 | hv = wv - zv #Water depths at vertices |
---|
[229] | 516 | |
---|
[443] | 517 | #Update |
---|
| 518 | for k in range(domain.number_of_elements): |
---|
| 519 | hmax = max(hv[k, :]) |
---|
| 520 | |
---|
| 521 | if hmax < domain.minimum_allowed_height: |
---|
| 522 | #Control level |
---|
| 523 | wv[k, :] = zv[k, :] |
---|
[242] | 524 | |
---|
[443] | 525 | #Control momentum |
---|
| 526 | xmomv[k,:] = ymomv[k,:] = 0.0 |
---|
| 527 | |
---|
| 528 | |
---|
| 529 | def protect_against_infinitesimal_and_negative_heights(domain): |
---|
[234] | 530 | """Protect against infinitesimal heights and associated high velocities |
---|
[229] | 531 | """ |
---|
| 532 | |
---|
[234] | 533 | #Shortcuts |
---|
[229] | 534 | wc = domain.quantities['level'].centroid_values |
---|
| 535 | zc = domain.quantities['elevation'].centroid_values |
---|
| 536 | xmomc = domain.quantities['xmomentum'].centroid_values |
---|
[234] | 537 | ymomc = domain.quantities['ymomentum'].centroid_values |
---|
| 538 | hc = wc - zc #Water depths at centroids |
---|
[229] | 539 | |
---|
[234] | 540 | #Update |
---|
[229] | 541 | for k in range(domain.number_of_elements): |
---|
| 542 | |
---|
[443] | 543 | if hc[k] < domain.minimum_allowed_height: |
---|
| 544 | #Control level |
---|
| 545 | wc[k] = zc[k] |
---|
| 546 | |
---|
[272] | 547 | #Control momentum |
---|
| 548 | xmomc[k] = ymomc[k] = 0.0 |
---|
[443] | 549 | |
---|
| 550 | #From 'newstyle |
---|
| 551 | #if hc[k] < domain.minimum_allowed_height: |
---|
| 552 | # if hc[k] < 0.0: |
---|
| 553 | # #Control level and height |
---|
| 554 | # wc[k] = zc[k] |
---|
| 555 | # |
---|
| 556 | # #Control momentum |
---|
| 557 | # xmomc[k] = ymomc[k] = 0.0 |
---|
| 558 | #else: |
---|
[273] | 559 | |
---|
[234] | 560 | |
---|
[273] | 561 | |
---|
[443] | 562 | def protect_against_infinitesimal_and_negative_heights_c(domain): |
---|
[273] | 563 | """Protect against infinitesimal heights and associated high velocities |
---|
| 564 | """ |
---|
| 565 | |
---|
| 566 | #Shortcuts |
---|
| 567 | wc = domain.quantities['level'].centroid_values |
---|
| 568 | zc = domain.quantities['elevation'].centroid_values |
---|
| 569 | xmomc = domain.quantities['xmomentum'].centroid_values |
---|
| 570 | ymomc = domain.quantities['ymomentum'].centroid_values |
---|
| 571 | |
---|
| 572 | from shallow_water_ext import protect |
---|
| 573 | |
---|
| 574 | protect(domain.minimum_allowed_height, wc, zc, xmomc, ymomc) |
---|
| 575 | |
---|
[443] | 576 | |
---|
[234] | 577 | def balance_deep_and_shallow(domain): |
---|
[266] | 578 | """Compute linear combination between stage as computed by |
---|
| 579 | gradient-limiters and stage computed as constant height above bed. |
---|
| 580 | The former takes precedence when heights are large compared to the |
---|
| 581 | bed slope while the latter takes precedence when heights are |
---|
| 582 | relatively small. Anything in between is computed as a balanced |
---|
| 583 | linear combination in order to avoid numerical disturbances which |
---|
| 584 | would otherwise appear as a result of hard switching between |
---|
| 585 | modes. |
---|
[244] | 586 | """ |
---|
[229] | 587 | |
---|
[234] | 588 | #Shortcuts |
---|
[229] | 589 | wc = domain.quantities['level'].centroid_values |
---|
| 590 | zc = domain.quantities['elevation'].centroid_values |
---|
| 591 | hc = wc - zc |
---|
[234] | 592 | |
---|
[229] | 593 | wv = domain.quantities['level'].vertex_values |
---|
| 594 | zv = domain.quantities['elevation'].vertex_values |
---|
| 595 | hv = wv-zv |
---|
| 596 | |
---|
[244] | 597 | |
---|
[229] | 598 | #Computed linear combination between constant levels and and |
---|
| 599 | #levels parallel to the bed elevation. |
---|
[244] | 600 | for k in range(domain.number_of_elements): |
---|
[229] | 601 | #Compute maximal variation in bed elevation |
---|
| 602 | # This quantitiy is |
---|
| 603 | # dz = max_i abs(z_i - z_c) |
---|
| 604 | # and it is independent of dimension |
---|
| 605 | # In the 1d case zc = (z0+z1)/2 |
---|
| 606 | # In the 2d case zc = (z0+z1+z2)/3 |
---|
| 607 | |
---|
[244] | 608 | dz = max(abs(zv[k,0]-zc[k]), |
---|
| 609 | abs(zv[k,1]-zc[k]), |
---|
| 610 | abs(zv[k,2]-zc[k])) |
---|
[229] | 611 | |
---|
[244] | 612 | |
---|
| 613 | hmin = min( hv[k,:] ) |
---|
[229] | 614 | |
---|
[267] | 615 | |
---|
[229] | 616 | #Create alpha in [0,1], where alpha==0 means using shallow |
---|
[244] | 617 | #first order scheme and alpha==1 means using the stage w as |
---|
| 618 | #computed by the gradient limiter (1st or 2nd order) |
---|
| 619 | # |
---|
| 620 | #If hmin > dz/2 then alpha = 1 and the bed will have no effect |
---|
| 621 | #If hmin < 0 then alpha = 0 reverting to constant height above bed. |
---|
[229] | 622 | |
---|
[244] | 623 | if dz > 0.0: |
---|
| 624 | alpha = max( min( 2*hmin/dz, 1.0), 0.0 ) |
---|
| 625 | else: |
---|
| 626 | #Flat bed |
---|
[229] | 627 | alpha = 1.0 |
---|
[244] | 628 | |
---|
| 629 | |
---|
| 630 | #Weighted balance between stage parallel to bed elevation |
---|
[246] | 631 | #(wvi = zvi + hc) and stage as computed by 1st or 2nd |
---|
| 632 | #order gradient limiter |
---|
[244] | 633 | #(wvi = zvi + hvi) where i=0,1,2 denotes the vertex ids |
---|
| 634 | # |
---|
| 635 | #It follows that the updated wvi is |
---|
| 636 | # wvi := (1-alpha)*(zvi+hc) + alpha*(zvi+hvi) = |
---|
| 637 | # zvi + hc + alpha*(hvi - hc) |
---|
| 638 | # |
---|
| 639 | #Note that hvi = zc+hc-zvi in the first order case (constant). |
---|
| 640 | |
---|
[229] | 641 | if alpha < 1: |
---|
| 642 | for i in range(3): |
---|
| 643 | wv[k,i] = zv[k,i] + hc[k] + alpha*(hv[k,i]-hc[k]) |
---|
| 644 | |
---|
| 645 | |
---|
| 646 | #Momentums at centroids |
---|
| 647 | xmomc = domain.quantities['xmomentum'].centroid_values |
---|
| 648 | ymomc = domain.quantities['ymomentum'].centroid_values |
---|
| 649 | |
---|
| 650 | #Momentums at vertices |
---|
| 651 | xmomv = domain.quantities['xmomentum'].vertex_values |
---|
| 652 | ymomv = domain.quantities['ymomentum'].vertex_values |
---|
| 653 | |
---|
| 654 | # Update momentum as a linear combination of |
---|
| 655 | # xmomc and ymomc (shallow) and momentum |
---|
| 656 | # from extrapolator xmomv and ymomv (deep). |
---|
[267] | 657 | xmomv[k,:] = (1-alpha)*xmomc[k] + alpha*xmomv[k,:] |
---|
| 658 | ymomv[k,:] = (1-alpha)*ymomc[k] + alpha*ymomv[k,:] |
---|
| 659 | |
---|
[229] | 660 | |
---|
[266] | 661 | |
---|
| 662 | def balance_deep_and_shallow_c(domain): |
---|
| 663 | """Wrapper for C implementation |
---|
| 664 | """ |
---|
[229] | 665 | |
---|
[266] | 666 | #Shortcuts |
---|
| 667 | wc = domain.quantities['level'].centroid_values |
---|
| 668 | zc = domain.quantities['elevation'].centroid_values |
---|
| 669 | hc = wc - zc |
---|
| 670 | |
---|
| 671 | wv = domain.quantities['level'].vertex_values |
---|
| 672 | zv = domain.quantities['elevation'].vertex_values |
---|
| 673 | hv = wv-zv |
---|
[229] | 674 | |
---|
[266] | 675 | #Momentums at centroids |
---|
| 676 | xmomc = domain.quantities['xmomentum'].centroid_values |
---|
| 677 | ymomc = domain.quantities['ymomentum'].centroid_values |
---|
[229] | 678 | |
---|
[266] | 679 | #Momentums at vertices |
---|
| 680 | xmomv = domain.quantities['xmomentum'].vertex_values |
---|
| 681 | ymomv = domain.quantities['ymomentum'].vertex_values |
---|
| 682 | |
---|
| 683 | |
---|
| 684 | |
---|
| 685 | from shallow_water_ext import balance_deep_and_shallow |
---|
[267] | 686 | balance_deep_and_shallow(wc, zc, hc, wv, zv, hv, |
---|
[266] | 687 | xmomc, ymomc, xmomv, ymomv) |
---|
| 688 | |
---|
| 689 | |
---|
| 690 | |
---|
| 691 | |
---|
[229] | 692 | ############################################### |
---|
| 693 | #Boundary - specific to the shallow water wave equation |
---|
| 694 | class Reflective_boundary(Boundary): |
---|
| 695 | """Reflective boundary returns same conserved quantities as |
---|
| 696 | those present in its neighbour volume but reflected. |
---|
| 697 | |
---|
| 698 | This class is specific to the shallow water equation as it |
---|
| 699 | works with the momentum quantities assumed to be the second |
---|
| 700 | and third conserved quantities. |
---|
| 701 | """ |
---|
| 702 | |
---|
| 703 | def __init__(self, domain = None): |
---|
| 704 | Boundary.__init__(self) |
---|
| 705 | |
---|
| 706 | if domain is None: |
---|
| 707 | msg = 'Domain must be specified for reflective boundary' |
---|
| 708 | raise msg |
---|
| 709 | |
---|
[263] | 710 | #Handy shorthands |
---|
| 711 | self.level = domain.quantities['level'].edge_values |
---|
| 712 | self.xmom = domain.quantities['xmomentum'].edge_values |
---|
| 713 | self.ymom = domain.quantities['ymomentum'].edge_values |
---|
| 714 | self.normals = domain.normals |
---|
| 715 | |
---|
| 716 | from Numeric import zeros, Float |
---|
| 717 | self.conserved_quantities = zeros(3, Float) |
---|
[229] | 718 | |
---|
| 719 | def __repr__(self): |
---|
| 720 | return 'Reflective_boundary' |
---|
| 721 | |
---|
| 722 | |
---|
| 723 | def evaluate(self, vol_id, edge_id): |
---|
| 724 | """Reflective boundaries reverses the outward momentum |
---|
| 725 | of the volume they serve. |
---|
| 726 | """ |
---|
[263] | 727 | |
---|
| 728 | q = self.conserved_quantities |
---|
| 729 | q[0] = self.level[vol_id, edge_id] |
---|
| 730 | q[1] = self.xmom[vol_id, edge_id] |
---|
| 731 | q[2] = self.ymom[vol_id, edge_id] |
---|
| 732 | |
---|
| 733 | normal = self.normals[vol_id, 2*edge_id:2*edge_id+2] |
---|
| 734 | |
---|
[229] | 735 | |
---|
| 736 | r = rotate(q, normal, direction = 1) |
---|
| 737 | r[1] = -r[1] |
---|
| 738 | q = rotate(r, normal, direction = -1) |
---|
| 739 | |
---|
| 740 | return q |
---|
| 741 | |
---|
| 742 | |
---|
| 743 | ######################### |
---|
| 744 | #Standard forcing terms: |
---|
| 745 | # |
---|
| 746 | def gravity(domain): |
---|
| 747 | """Implement forcing function for bed slope working with |
---|
| 748 | consecutive data structures of class Volume |
---|
| 749 | """ |
---|
| 750 | |
---|
| 751 | from util import gradient |
---|
| 752 | from Numeric import zeros, Float, array, sum |
---|
| 753 | |
---|
[246] | 754 | xmom = domain.quantities['xmomentum'].explicit_update |
---|
| 755 | ymom = domain.quantities['ymomentum'].explicit_update |
---|
| 756 | |
---|
[229] | 757 | Level = domain.quantities['level'] |
---|
| 758 | Elevation = domain.quantities['elevation'] |
---|
| 759 | h = Level.edge_values - Elevation.edge_values |
---|
[246] | 760 | v = Elevation.vertex_values |
---|
[229] | 761 | |
---|
[246] | 762 | x = domain.get_vertex_coordinates() |
---|
| 763 | g = domain.g |
---|
| 764 | |
---|
[229] | 765 | for k in range(domain.number_of_elements): |
---|
| 766 | avg_h = sum( h[k,:] )/3 |
---|
| 767 | |
---|
| 768 | #Compute bed slope |
---|
[246] | 769 | x0, y0, x1, y1, x2, y2 = x[k,:] |
---|
| 770 | z0, z1, z2 = v[k,:] |
---|
[229] | 771 | |
---|
| 772 | zx, zy = gradient(x0, y0, x1, y1, x2, y2, z0, z1, z2) |
---|
| 773 | |
---|
| 774 | #Update momentum |
---|
[246] | 775 | xmom[k] += -g*zx*avg_h |
---|
| 776 | ymom[k] += -g*zy*avg_h |
---|
[229] | 777 | |
---|
| 778 | |
---|
[246] | 779 | def gravity_c(domain): |
---|
| 780 | """Wrapper calling C version |
---|
| 781 | """ |
---|
| 782 | |
---|
| 783 | xmom = domain.quantities['xmomentum'].explicit_update |
---|
| 784 | ymom = domain.quantities['ymomentum'].explicit_update |
---|
| 785 | |
---|
| 786 | Level = domain.quantities['level'] |
---|
| 787 | Elevation = domain.quantities['elevation'] |
---|
| 788 | h = Level.edge_values - Elevation.edge_values |
---|
| 789 | v = Elevation.vertex_values |
---|
| 790 | |
---|
| 791 | x = domain.get_vertex_coordinates() |
---|
| 792 | g = domain.g |
---|
| 793 | |
---|
| 794 | |
---|
| 795 | from shallow_water_ext import gravity |
---|
| 796 | gravity(g, h, v, x, xmom, ymom) |
---|
| 797 | |
---|
| 798 | |
---|
[229] | 799 | def manning_friction(domain): |
---|
| 800 | """Apply (Manning) friction to water momentum |
---|
| 801 | """ |
---|
| 802 | |
---|
| 803 | from math import sqrt |
---|
| 804 | |
---|
[246] | 805 | w = domain.quantities['level'].centroid_values |
---|
| 806 | uh = domain.quantities['xmomentum'].centroid_values |
---|
| 807 | vh = domain.quantities['ymomentum'].centroid_values |
---|
| 808 | eta = domain.quantities['friction'].centroid_values |
---|
[229] | 809 | |
---|
[246] | 810 | xmom_update = domain.quantities['xmomentum'].semi_implicit_update |
---|
| 811 | ymom_update = domain.quantities['ymomentum'].semi_implicit_update |
---|
| 812 | |
---|
| 813 | N = domain.number_of_elements |
---|
| 814 | eps = domain.minimum_allowed_height |
---|
| 815 | g = domain.g |
---|
| 816 | |
---|
| 817 | for k in range(N): |
---|
| 818 | if w[k] >= eps: |
---|
| 819 | S = -g * eta[k]**2 * sqrt((uh[k]**2 + vh[k]**2)) |
---|
| 820 | S /= w[k]**(7.0/3) |
---|
[229] | 821 | |
---|
| 822 | #Update momentum |
---|
[246] | 823 | xmom_update[k] += S*uh[k] |
---|
| 824 | ymom_update[k] += S*vh[k] |
---|
[229] | 825 | |
---|
[246] | 826 | |
---|
| 827 | def manning_friction_c(domain): |
---|
| 828 | """Wrapper for c version |
---|
| 829 | """ |
---|
[229] | 830 | |
---|
[246] | 831 | |
---|
[263] | 832 | xmom = domain.quantities['xmomentum'] |
---|
| 833 | ymom = domain.quantities['ymomentum'] |
---|
| 834 | |
---|
[246] | 835 | w = domain.quantities['level'].centroid_values |
---|
[263] | 836 | uh = xmom.centroid_values |
---|
| 837 | vh = ymom.centroid_values |
---|
[246] | 838 | eta = domain.quantities['friction'].centroid_values |
---|
| 839 | |
---|
[263] | 840 | xmom_update = xmom.semi_implicit_update |
---|
| 841 | ymom_update = ymom.semi_implicit_update |
---|
[246] | 842 | |
---|
| 843 | N = domain.number_of_elements |
---|
| 844 | eps = domain.minimum_allowed_height |
---|
| 845 | g = domain.g |
---|
| 846 | |
---|
| 847 | from shallow_water_ext import manning_friction |
---|
| 848 | manning_friction(g, eps, w, uh, vh, eta, xmom_update, ymom_update) |
---|
| 849 | |
---|
| 850 | |
---|
[229] | 851 | ########################### |
---|
| 852 | ########################### |
---|
| 853 | #Geometries |
---|
| 854 | |
---|
| 855 | |
---|
| 856 | #FIXME: Rethink this way of creating values. |
---|
| 857 | |
---|
| 858 | |
---|
| 859 | class Weir: |
---|
| 860 | """Set a bathymetry for weir with a hole and a downstream gutter |
---|
| 861 | x,y are assumed to be in the unit square |
---|
| 862 | """ |
---|
| 863 | |
---|
| 864 | def __init__(self, stage): |
---|
| 865 | self.inflow_stage = stage |
---|
| 866 | |
---|
| 867 | def __call__(self, x, y): |
---|
| 868 | from Numeric import zeros, Float |
---|
| 869 | from math import sqrt |
---|
| 870 | |
---|
| 871 | N = len(x) |
---|
| 872 | assert N == len(y) |
---|
| 873 | |
---|
| 874 | z = zeros(N, Float) |
---|
| 875 | for i in range(N): |
---|
| 876 | z[i] = -x[i]/2 #General slope |
---|
| 877 | |
---|
| 878 | #Flattish bit to the left |
---|
| 879 | if x[i] < 0.3: |
---|
| 880 | z[i] = -x[i]/10 |
---|
| 881 | |
---|
| 882 | #Weir |
---|
| 883 | if x[i] >= 0.3 and x[i] < 0.4: |
---|
| 884 | z[i] = -x[i]+0.9 |
---|
| 885 | |
---|
| 886 | #Dip |
---|
| 887 | x0 = 0.6 |
---|
| 888 | #depth = -1.3 |
---|
| 889 | depth = -1.0 |
---|
| 890 | #plateaux = -0.9 |
---|
| 891 | plateaux = -0.6 |
---|
| 892 | if y[i] < 0.7: |
---|
| 893 | if x[i] > x0 and x[i] < 0.9: |
---|
| 894 | z[i] = depth |
---|
| 895 | |
---|
| 896 | #RHS plateaux |
---|
| 897 | if x[i] >= 0.9: |
---|
| 898 | z[i] = plateaux |
---|
| 899 | |
---|
| 900 | |
---|
| 901 | elif y[i] >= 0.7 and y[i] < 1.5: |
---|
| 902 | #Restrict and deepen |
---|
| 903 | if x[i] >= x0 and x[i] < 0.8: |
---|
| 904 | z[i] = depth-(y[i]/3-0.3) |
---|
| 905 | #z[i] = depth-y[i]/5 |
---|
| 906 | #z[i] = depth |
---|
| 907 | elif x[i] >= 0.8: |
---|
| 908 | #RHS plateaux |
---|
| 909 | z[i] = plateaux |
---|
| 910 | |
---|
| 911 | elif y[i] >= 1.5: |
---|
| 912 | if x[i] >= x0 and x[i] < 0.8 + (y[i]-1.5)/1.2: |
---|
| 913 | #Widen up and stay at constant depth |
---|
| 914 | z[i] = depth-1.5/5 |
---|
| 915 | elif x[i] >= 0.8 + (y[i]-1.5)/1.2: |
---|
| 916 | #RHS plateaux |
---|
| 917 | z[i] = plateaux |
---|
| 918 | |
---|
| 919 | |
---|
| 920 | #Hole in weir (slightly higher than inflow condition) |
---|
| 921 | if x[i] >= 0.3 and x[i] < 0.4 and y[i] > 0.2 and y[i] < 0.4: |
---|
| 922 | z[i] = -x[i]+self.inflow_stage + 0.02 |
---|
| 923 | |
---|
| 924 | #Channel behind weir |
---|
| 925 | x0 = 0.5 |
---|
| 926 | if x[i] >= 0.4 and x[i] < x0 and y[i] > 0.2 and y[i] < 0.4: |
---|
| 927 | z[i] = -x[i]+self.inflow_stage + 0.02 |
---|
| 928 | |
---|
| 929 | if x[i] >= x0 and x[i] < 0.6 and y[i] > 0.2 and y[i] < 0.4: |
---|
| 930 | #Flatten it out towards the end |
---|
| 931 | z[i] = -x0+self.inflow_stage + 0.02 + (x0-x[i])/5 |
---|
| 932 | |
---|
| 933 | #Hole to the east |
---|
| 934 | x0 = 1.1; y0 = 0.35 |
---|
| 935 | #if x[i] < -0.2 and y < 0.5: |
---|
| 936 | if sqrt((2*(x[i]-x0))**2 + (2*(y[i]-y0))**2) < 0.2: |
---|
| 937 | z[i] = sqrt(((x[i]-x0))**2 + ((y[i]-y0))**2)-1.0 |
---|
| 938 | |
---|
| 939 | #Tiny channel draining hole |
---|
| 940 | if x[i] >= 1.14 and x[i] < 1.2 and y[i] >= 0.4 and y[i] < 0.6: |
---|
| 941 | z[i] = -0.9 #North south |
---|
| 942 | |
---|
| 943 | if x[i] >= 0.9 and x[i] < 1.18 and y[i] >= 0.58 and y[i] < 0.65: |
---|
| 944 | z[i] = -1.0 + (x[i]-0.9)/3 #East west |
---|
| 945 | |
---|
| 946 | |
---|
| 947 | |
---|
| 948 | #Stuff not in use |
---|
| 949 | |
---|
| 950 | #Upward slope at inlet to the north west |
---|
| 951 | #if x[i] < 0.0: # and y[i] > 0.5: |
---|
| 952 | # #z[i] = -y[i]+0.5 #-x[i]/2 |
---|
| 953 | # z[i] = x[i]/4 - y[i]**2 + 0.5 |
---|
| 954 | |
---|
| 955 | #Hole to the west |
---|
| 956 | #x0 = -0.4; y0 = 0.35 # center |
---|
| 957 | #if sqrt((2*(x[i]-x0))**2 + (2*(y[i]-y0))**2) < 0.2: |
---|
| 958 | # z[i] = sqrt(((x[i]-x0))**2 + ((y[i]-y0))**2)-0.2 |
---|
| 959 | |
---|
| 960 | |
---|
| 961 | |
---|
| 962 | |
---|
| 963 | |
---|
| 964 | return z/2 |
---|
| 965 | |
---|
| 966 | class Weir_simple: |
---|
| 967 | """Set a bathymetry for weir with a hole and a downstream gutter |
---|
| 968 | x,y are assumed to be in the unit square |
---|
| 969 | """ |
---|
| 970 | |
---|
| 971 | def __init__(self, stage): |
---|
| 972 | self.inflow_stage = stage |
---|
| 973 | |
---|
| 974 | def __call__(self, x, y): |
---|
| 975 | from Numeric import zeros, Float |
---|
| 976 | |
---|
| 977 | N = len(x) |
---|
| 978 | assert N == len(y) |
---|
| 979 | |
---|
| 980 | z = zeros(N, Float) |
---|
| 981 | for i in range(N): |
---|
| 982 | z[i] = -x[i] #General slope |
---|
| 983 | |
---|
| 984 | #Flat bit to the left |
---|
| 985 | if x[i] < 0.3: |
---|
| 986 | z[i] = -x[i]/10 #General slope |
---|
| 987 | |
---|
| 988 | #Weir |
---|
| 989 | if x[i] > 0.3 and x[i] < 0.4: |
---|
| 990 | z[i] = -x[i]+0.9 |
---|
| 991 | |
---|
| 992 | #Dip |
---|
| 993 | if x[i] > 0.6 and x[i] < 0.9: |
---|
| 994 | z[i] = -x[i]-0.5 #-y[i]/5 |
---|
| 995 | |
---|
| 996 | #Hole in weir (slightly higher than inflow condition) |
---|
| 997 | if x[i] > 0.3 and x[i] < 0.4 and y[i] > 0.2 and y[i] < 0.4: |
---|
| 998 | z[i] = -x[i]+self.inflow_stage + 0.05 |
---|
| 999 | |
---|
| 1000 | |
---|
| 1001 | return z/2 |
---|
| 1002 | |
---|
| 1003 | |
---|
| 1004 | |
---|
| 1005 | class Constant_height: |
---|
| 1006 | """Set an initial condition with constant water height, e.g |
---|
| 1007 | stage s = z+h |
---|
| 1008 | """ |
---|
| 1009 | def __init__(self, W, h): |
---|
| 1010 | self.W = W |
---|
| 1011 | self.h = h |
---|
| 1012 | |
---|
| 1013 | def __call__(self, x, y): |
---|
| 1014 | if self.W is None: |
---|
| 1015 | from Numeric import ones, Float |
---|
| 1016 | return self.h*ones(len(x), Float) |
---|
| 1017 | else: |
---|
| 1018 | return self.W(x,y) + self.h |
---|
| 1019 | |
---|
| 1020 | |
---|
[443] | 1021 | |
---|
[229] | 1022 | ############################################## |
---|
| 1023 | #Initialise module |
---|
| 1024 | |
---|
[240] | 1025 | |
---|
| 1026 | import compile |
---|
| 1027 | if compile.can_use_C_extension('shallow_water_ext.c'): |
---|
| 1028 | #Replace python version with c implementations |
---|
[259] | 1029 | |
---|
[240] | 1030 | from shallow_water_ext import rotate |
---|
| 1031 | compute_fluxes = compute_fluxes_c |
---|
[246] | 1032 | gravity = gravity_c |
---|
| 1033 | manning_friction = manning_friction_c |
---|
[273] | 1034 | balance_deep_and_shallow = balance_deep_and_shallow_c |
---|
[443] | 1035 | protect_against_infinitesimal_and_negative_heights = protect_against_infinitesimal_and_negative_heights_c |
---|
| 1036 | |
---|
[259] | 1037 | |
---|
[240] | 1038 | #distribute_to_vertices_and_edges = distribute_to_vertices_and_edges_c |
---|
| 1039 | |
---|
| 1040 | |
---|
| 1041 | #Optimisation with psyco |
---|
| 1042 | from config import use_psyco |
---|
| 1043 | if use_psyco: |
---|
| 1044 | try: |
---|
| 1045 | import psyco |
---|
| 1046 | except: |
---|
| 1047 | msg = 'WARNING: psyco (speedup) could not import'+\ |
---|
| 1048 | ', you may want to consider installing it' |
---|
| 1049 | print msg |
---|
| 1050 | else: |
---|
[278] | 1051 | psyco.bind(Domain.distribute_to_vertices_and_edges) |
---|
| 1052 | psyco.bind(Domain.compute_fluxes) |
---|
| 1053 | |
---|
[240] | 1054 | if __name__ == "__main__": |
---|
| 1055 | pass |
---|