source: inundation/ga/storm_surge/pyvolution/shallow_water.py @ 620

Last change on this file since 620 was 620, checked in by ole, 20 years ago

Added functionality for storing momentum (x and y) as well as water levels in sww files.

File size: 42.7 KB
Line 
1"""Class Domain -
22D triangular domains for finite-volume computations of
3the shallow water wave equation.
4
5This module contains a specialisation of class Domain from module domain.py
6consisting of methods specific to the Shallow Water Wave Equation
7
8FIXME: Write equations here!
9
10
11Conserved quantities are w (water level or stage), uh (x momentum)
12and vh (y momentum).
13
14
15Ole Nielsen, Stephen Roberts, Duncan Gray, Christopher Zoppou
16Geoscience Australia, 2004   
17"""
18
19from domain import *
20Generic_domain = Domain #Rename
21
22class Domain(Generic_domain):
23
24    def __init__(self, coordinates, vertices, boundary = None,
25                 tagged_elements = None):
26
27        conserved_quantities = ['level', 'xmomentum', 'ymomentum']
28        other_quantities = ['elevation', 'friction'] 
29       
30        Generic_domain.__init__(self, coordinates, vertices, boundary,
31                                conserved_quantities, other_quantities,
32                                tagged_elements)
33
34        from config import minimum_allowed_height, g
35        self.minimum_allowed_height = minimum_allowed_height
36        self.g = g
37
38        self.forcing_terms.append(gravity)
39        self.forcing_terms.append(manning_friction)
40
41        #Realtime visualisation
42        self.visualise = False
43
44        #Stored output
45        self.store = False
46        self.format = 'sww'   
47        self.smooth = True
48
49        #Reduction operation for get_vertex_values             
50        #from util import mean
51        #self.reduction = mean
52        self.reduction = min  #Looks better near steep slopes
53       
54        self.quantities_to_be_stored = ['level']
55
56       
57        #Establish shortcuts to relevant quantities (for efficiency)
58        #self.w = self.quantities['level']
59        #self.uh = self.quantities['xmomentum']                 
60        #self.vh = self.quantities['ymomentum']                         
61        #self.z = self.quantities['elevation']         
62        #self.eta = self.quantities['friction']                 
63
64    def check_integrity(self):
65        Generic_domain.check_integrity(self)
66
67        #Check that we are solving the shallow water wave equation
68
69        msg = 'First conserved quantity must be "level"'
70        assert self.conserved_quantities[0] == 'level', msg
71        msg = 'Second conserved quantity must be "xmomentum"'
72        assert self.conserved_quantities[1] == 'xmomentum', msg
73        msg = 'Third conserved quantity must be "ymomentum"'
74        assert self.conserved_quantities[2] == 'ymomentum', msg
75       
76
77        #Check that levels are >= bed elevation
78        from Numeric import alltrue, greater_equal
79       
80        level = self.quantities['level']
81        bed = self.quantities['elevation']       
82
83        msg = 'All water levels must be greater than the bed elevation'
84        assert alltrue( greater_equal(
85            level.vertex_values, bed.vertex_values )), msg
86
87        assert alltrue( greater_equal(
88            level.edge_values, bed.edge_values )), msg
89
90        assert alltrue( greater_equal(
91            level.centroid_values, bed.centroid_values )), msg       
92
93
94    def compute_fluxes(self):
95        #Call correct module function
96        #(either from this module or C-extension)
97        compute_fluxes(self)
98
99    def distribute_to_vertices_and_edges(self):
100        #Call correct module function
101        #(either from this module or C-extension)       
102        distribute_to_vertices_and_edges(self)
103
104
105    #FIXME: Under construction   
106#     def set_defaults(self):
107#         """Set default values for uninitialised quantities.
108#         This is specific to the shallow water wave equation
109#         Defaults for 'elevation', 'friction', 'xmomentum' and 'ymomentum'
110#         are 0.0. Default for 'level' is whatever the value of 'elevation'.
111#         """
112
113#         for name in self.other_quantities + self.conserved_quantities:
114#             print name
115#             print self.quantities.keys()
116#             if not self.quantities.has_key(name):
117#                 if name == 'level':
118                   
119#                     if self.quantities.has_key('elevation'):
120#                         z = self.quantities['elevation'].vertex_values
121#                         self.set_quantity(name, z)
122#                     else:
123#                         self.set_quantity(name, 0.0)                       
124#                 else:   
125#                     self.set_quantity(name, 0.0)
126                   
127
128
129#         #Lift negative heights up
130#         #z = self.quantities['elevation'].vertex_values
131#         #w = self.quantities['level'].vertex_values
132
133#         #h = w-z
134
135#         #for k in range(h.shape[0]):
136#         #    for i in range(3):
137#         #        if h[k, i] < 0.0:
138#         #            w[k, i] = z[k, i]
139
140                   
141#         #self.quantities['level'].interpolate()
142               
143
144
145    def evolve(self, yieldstep = None, finaltime = None):
146        """Specialisation of basic evolve method from parent class
147        """
148
149        #Call check integrity here rather than from user scripts
150        #self.check_integrity()
151       
152        #Initialise real time viz if requested
153        if self.visualise is True and self.time == 0.0:
154            import realtime_visualisation as visualise
155            visualise.create_surface(self)
156
157        #Store model data, e.g. for visualisation   
158        if self.store is True and self.time == 0.0:
159            self.initialise_storage()
160            #print 'Storing results in ' + self.writer.filename
161        else:
162            #print 'Results will not be stored.'
163            #print 'To store results set domain.store = True'
164            pass
165            #FIXME: Diagnostic output should be controlled by
166            # a 'verbose' flag living in domain (or in a parent class)
167
168        #Call basic machinery from parent class
169        for t in Generic_domain.evolve(self, yieldstep, finaltime):
170            #Real time viz
171            if self.visualise is True:
172                visualise.update(self)
173
174            #Store model data, e.g. for subsequent visualisation   
175            if self.store is True:
176                #self.store_timestep(['level', 'xmomentum', 'ymomentum'])
177                self.store_timestep(self.quantities_to_be_stored)
178               
179                #FIXME: Could maybe be taken from specified list
180                #of 'store every step' quantities       
181
182            #Pass control on to outer loop for more specific actions   
183            yield(t)
184       
185
186    def initialise_storage(self):       
187        """Create and initialise self.writer object for storing data.
188        Also, save x,y and bed elevation
189        """
190
191        import data_manager
192       
193        #Initialise writer
194        self.writer = data_manager.get_dataobject(self, mode = 'w')
195
196        #Store vertices and connectivity
197        self.writer.store_connectivity()                   
198
199    def store_timestep(self, name):
200        """Store named quantity and time.
201
202        Precondition:
203           self.write has been initialised
204        """
205        self.writer.store_timestep(name)
206       
207#Rotation of momentum vector
208def rotate(q, normal, direction = 1):
209    """Rotate the momentum component q (q[1], q[2])
210    from x,y coordinates to coordinates based on normal vector.
211
212    If direction is negative the rotation is inverted.
213   
214    Input vector is preserved
215
216    This function is specific to the shallow water wave equation
217    """
218
219    #FIXME: Needs to be tested
220
221    from Numeric import zeros, Float
222   
223    assert len(q) == 3,\
224           'Vector of conserved quantities must have length 3'\
225           'for 2D shallow water equation'
226
227    try:
228        l = len(normal)
229    except:
230        raise 'Normal vector must be an Numeric array'
231
232    #FIXME: Put this test into C-extension as well
233    assert l == 2, 'Normal vector must have 2 components'
234
235 
236    n1 = normal[0]
237    n2 = normal[1]   
238   
239    r = zeros(len(q), Float) #Rotated quantities
240    r[0] = q[0]              #First quantity, height, is not rotated
241
242    if direction == -1:
243        n2 = -n2
244
245
246    r[1] =  n1*q[1] + n2*q[2]
247    r[2] = -n2*q[1] + n1*q[2]
248   
249    return r
250
251
252
253####################################
254# Flux computation       
255def flux_function(normal, ql, qr, zl, zr): 
256    """Compute fluxes between volumes for the shallow water wave equation
257    cast in terms of w = h+z using the 'central scheme' as described in
258
259    Kurganov, Noelle, Petrova. 'Semidiscrete Central-Upwind Schemes For
260    Hyperbolic Conservation Laws and Hamilton-Jacobi Equations'.
261    Siam J. Sci. Comput. Vol. 23, No. 3, pp. 707-740.
262
263    The implemented formula is given in equation (3.15) on page 714
264
265    Conserved quantities w, uh, vh are stored as elements 0, 1 and 2
266    in the numerical vectors ql an qr.
267
268    Bed elevations zl and zr.
269    """
270
271    from config import g, epsilon
272    from math import sqrt
273    from Numeric import array
274
275    #Align momentums with x-axis
276    q_left  = rotate(ql, normal, direction = 1)
277    q_right = rotate(qr, normal, direction = 1)
278
279    z = (zl+zr)/2 #Take average of field values
280
281    w_left  = q_left[0]   #w=h+z
282    h_left  = w_left-z
283    uh_left = q_left[1]
284
285    if h_left < epsilon:
286        u_left = 0.0  #Could have been negative
287        h_left = 0.0
288    else:   
289        u_left  = uh_left/h_left
290
291
292    w_right  = q_right[0]  #w=h+z
293    h_right  = w_right-z
294    uh_right = q_right[1]
295
296
297    if h_right < epsilon:
298        u_right = 0.0  #Could have been negative
299        h_right = 0.0
300    else:   
301        u_right  = uh_right/h_right
302
303    vh_left  = q_left[2]
304    vh_right = q_right[2]       
305
306    soundspeed_left  = sqrt(g*h_left) 
307    soundspeed_right = sqrt(g*h_right)
308
309    #Maximal wave speed
310    s_max = max(u_left+soundspeed_left, u_right+soundspeed_right, 0)
311   
312    #Minimal wave speed       
313    s_min = min(u_left-soundspeed_left, u_right-soundspeed_right, 0)
314
315    #Flux computation
316    flux_left  = array([u_left*h_left,
317                        u_left*uh_left + 0.5*g*h_left**2,
318                        u_left*vh_left])
319    flux_right = array([u_right*h_right,
320                        u_right*uh_right + 0.5*g*h_right**2,
321                        u_right*vh_right])   
322
323    denom = s_max-s_min
324    if denom == 0.0:
325        edgeflux = array([0.0, 0.0, 0.0])
326        max_speed = 0.0
327    else:   
328        edgeflux = (s_max*flux_left - s_min*flux_right)/denom
329        edgeflux += s_max*s_min*(q_right-q_left)/denom
330
331        edgeflux = rotate(edgeflux, normal, direction=-1)
332        max_speed = max(abs(s_max), abs(s_min))
333
334    return edgeflux, max_speed       
335
336
337def compute_fluxes(domain):
338    """Compute all fluxes and the timestep suitable for all volumes
339    in domain.
340
341    Compute total flux for each conserved quantity using "flux_function"
342
343    Fluxes across each edge are scaled by edgelengths and summed up
344    Resulting flux is then scaled by area and stored in
345    explicit_update for each of the three conserved quantities
346    level, xmomentum and ymomentum   
347
348    The maximal allowable speed computed by the flux_function for each volume
349    is converted to a timestep that must not be exceeded. The minimum of
350    those is computed as the next overall timestep.
351
352    Post conditions:
353      domain.explicit_update is reset to computed flux values
354      domain.timestep is set to the largest step satisfying all volumes.
355    """
356
357    import sys
358    from Numeric import zeros, Float
359
360    N = domain.number_of_elements
361   
362    #Shortcuts
363    Level = domain.quantities['level']
364    Xmom = domain.quantities['xmomentum']
365    Ymom = domain.quantities['ymomentum']
366    Bed = domain.quantities['elevation']       
367
368    #Arrays
369    level = Level.edge_values
370    xmom =  Xmom.edge_values
371    ymom =  Ymom.edge_values
372    bed =   Bed.edge_values   
373
374    level_bdry = Level.boundary_values
375    xmom_bdry =  Xmom.boundary_values
376    ymom_bdry =  Ymom.boundary_values
377   
378    flux = zeros(3, Float) #Work array for summing up fluxes
379
380    #Loop
381    timestep = float(sys.maxint)   
382    for k in range(N):
383
384        flux[:] = 0.  #Reset work array
385        for i in range(3):
386            #Quantities inside volume facing neighbour i
387            ql = [level[k, i], xmom[k, i], ymom[k, i]]
388            zl = bed[k, i]
389
390            #Quantities at neighbour on nearest face
391            n = domain.neighbours[k,i] 
392            if n < 0:
393                m = -n-1 #Convert negative flag to index
394                qr = [level_bdry[m], xmom_bdry[m], ymom_bdry[m]]
395                zr = zl #Extend bed elevation to boundary
396            else:   
397                m = domain.neighbour_edges[k,i]
398                qr = [level[n, m], xmom[n, m], ymom[n, m]]
399                zr = bed[n, m]               
400
401               
402            #Outward pointing normal vector   
403            normal = domain.normals[k, 2*i:2*i+2]
404
405            #Flux computation using provided function
406            edgeflux, max_speed = flux_function(normal, ql, qr, zl, zr)
407            flux -= edgeflux * domain.edgelengths[k,i]
408
409            #Update optimal_timestep
410            try:
411                timestep = min(timestep, domain.radii[k]/max_speed)
412            except ZeroDivisionError:
413                pass
414
415        #Normalise by area and store for when all conserved
416        #quantities get updated
417        flux /= domain.areas[k]
418        Level.explicit_update[k] = flux[0]
419        Xmom.explicit_update[k] = flux[1]
420        Ymom.explicit_update[k] = flux[2]
421
422    #print 'FLUX l', Level.explicit_update
423    #print 'FLUX x', Xmom.explicit_update
424    #print 'FLUX y', Ymom.explicit_update   
425   
426    domain.timestep = timestep   
427
428
429def compute_fluxes_c(domain):
430    """Wrapper calling C version of compute fluxes
431    """
432
433    import sys
434    from Numeric import zeros, Float
435
436    N = domain.number_of_elements
437   
438    #Shortcuts
439    Level = domain.quantities['level']
440    Xmom = domain.quantities['xmomentum']
441    Ymom = domain.quantities['ymomentum']
442    Bed = domain.quantities['elevation']       
443
444    timestep = float(sys.maxint)   
445    from shallow_water_ext import compute_fluxes
446    domain.timestep = compute_fluxes(timestep, domain.epsilon, domain.g,
447                                     domain.neighbours,
448                                     domain.neighbour_edges,
449                                     domain.normals,
450                                     domain.edgelengths,                       
451                                     domain.radii,
452                                     domain.areas,
453                                     Level.edge_values, 
454                                     Xmom.edge_values, 
455                                     Ymom.edge_values, 
456                                     Bed.edge_values,   
457                                     Level.boundary_values,
458                                     Xmom.boundary_values,
459                                     Ymom.boundary_values,
460                                     Level.explicit_update,
461                                     Xmom.explicit_update,
462                                     Ymom.explicit_update)
463                                     
464
465####################################
466# Module functions for gradient limiting (distribute_to_vertices_and_edges)
467
468def distribute_to_vertices_and_edges(domain):
469    """Distribution from centroids to vertices specific to the
470    shallow water wave
471    equation.
472
473    It will ensure that h (w-z) is always non-negative even in the
474    presence of steep bed-slopes by taking a weighted average between shallow
475    and deep cases.
476   
477    In addition, all conserved quantities get distributed as per either a
478    constant (order==1) or a piecewise linear function (order==2).
479   
480    FIXME: more explanation about removal of artificial variability etc
481
482    Precondition:
483      All quantities defined at centroids and bed elevation defined at
484      vertices.
485
486    Postcondition
487      Conserved quantities defined at vertices
488   
489    """
490
491    #Remove very thin layers of water
492    protect_against_infinitesimal_and_negative_heights(domain)
493
494    #Extrapolate all conserved quantities
495    for name in domain.conserved_quantities:
496        Q = domain.quantities[name]
497        if domain.order == 1:
498            Q.extrapolate_first_order()
499        elif domain.order == 2:
500            Q.extrapolate_second_order()
501            Q.limit()                               
502        else:
503            raise 'Unknown order'
504
505    #Take bed elevation into account when water heights are small   
506    balance_deep_and_shallow(domain)
507
508    #Compute edge values by interpolation
509    for name in domain.conserved_quantities:
510        Q = domain.quantities[name]
511        Q.interpolate_from_vertices_to_edges()
512                   
513
514
515def dry(domain):
516    """Protect against infinitesimal heights and associated high velocities
517    at vertices
518    """
519
520    #FIXME: Experimental
521   
522    #Shortcuts
523    wv = domain.quantities['level'].vertex_values
524    zv = domain.quantities['elevation'].vertex_values
525    xmomv = domain.quantities['xmomentum'].vertex_values
526    ymomv = domain.quantities['ymomentum'].vertex_values           
527    hv = wv - zv  #Water depths at vertices
528
529    #Update
530    for k in range(domain.number_of_elements):
531        hmax = max(hv[k, :])
532       
533        if hmax < domain.minimum_allowed_height:
534            #Control level
535            wv[k, :] = zv[k, :]
536
537            #Control momentum
538            xmomv[k,:] = ymomv[k,:] = 0.0
539   
540
541def protect_against_infinitesimal_and_negative_heights(domain):
542    """Protect against infinitesimal heights and associated high velocities
543    """
544   
545    #Shortcuts
546    wc = domain.quantities['level'].centroid_values
547    zc = domain.quantities['elevation'].centroid_values
548    xmomc = domain.quantities['xmomentum'].centroid_values
549    ymomc = domain.quantities['ymomentum'].centroid_values           
550    hc = wc - zc  #Water depths at centroids   
551
552    #Update
553    for k in range(domain.number_of_elements):
554
555        if hc[k] < domain.minimum_allowed_height: 
556            #Control level
557            wc[k] = zc[k]
558
559            #Control momentum
560            xmomc[k] = ymomc[k] = 0.0
561       
562        #From 'newstyle
563        #if hc[k] < domain.minimum_allowed_height:       
564        #     if hc[k] < 0.0:
565        #            #Control level and height
566        #            wc[k] = zc[k]
567        #
568        #        #Control momentum
569        #        xmomc[k] = ymomc[k] = 0.0
570        #else:
571           
572
573
574def protect_against_infinitesimal_and_negative_heights_c(domain):
575    """Protect against infinitesimal heights and associated high velocities
576    """
577   
578    #Shortcuts
579    wc = domain.quantities['level'].centroid_values
580    zc = domain.quantities['elevation'].centroid_values
581    xmomc = domain.quantities['xmomentum'].centroid_values
582    ymomc = domain.quantities['ymomentum'].centroid_values           
583
584    from shallow_water_ext import protect
585
586    protect(domain.minimum_allowed_height, wc, zc, xmomc, ymomc)
587   
588   
589def balance_deep_and_shallow(domain):
590    """Compute linear combination between stage as computed by
591    gradient-limiters and stage computed as constant height above bed.
592    The former takes precedence when heights are large compared to the
593    bed slope while the latter takes precedence when heights are
594    relatively small.  Anything in between is computed as a balanced
595    linear combination in order to avoid numerical disturbances which
596    would otherwise appear as a result of hard switching between
597    modes.
598    """
599   
600    #Shortcuts
601    wc = domain.quantities['level'].centroid_values
602    zc = domain.quantities['elevation'].centroid_values
603    hc = wc - zc
604   
605    wv = domain.quantities['level'].vertex_values
606    zv = domain.quantities['elevation'].vertex_values
607    hv = wv-zv
608
609
610    #Computed linear combination between constant levels and and
611    #levels parallel to the bed elevation.     
612    for k in range(domain.number_of_elements):
613        #Compute maximal variation in bed elevation
614        #  This quantitiy is
615        #    dz = max_i abs(z_i - z_c)
616        #  and it is independent of dimension
617        #  In the 1d case zc = (z0+z1)/2
618        #  In the 2d case zc = (z0+z1+z2)/3
619
620        dz = max(abs(zv[k,0]-zc[k]),
621                 abs(zv[k,1]-zc[k]),
622                 abs(zv[k,2]-zc[k]))
623
624       
625        hmin = min( hv[k,:] )
626
627       
628        #Create alpha in [0,1], where alpha==0 means using shallow
629        #first order scheme and alpha==1 means using the stage w as
630        #computed by the gradient limiter (1st or 2nd order)
631        #
632        #If hmin > dz/2 then alpha = 1 and the bed will have no effect
633        #If hmin < 0 then alpha = 0 reverting to constant height above bed.
634     
635        if dz > 0.0:
636            alpha = max( min( 2*hmin/dz, 1.0), 0.0 )
637        else:
638            #Flat bed
639            alpha = 1.0
640
641
642        #Weighted balance between stage parallel to bed elevation
643        #(wvi = zvi + hc) and stage as computed by 1st or 2nd
644        #order gradient limiter
645        #(wvi = zvi + hvi) where i=0,1,2 denotes the vertex ids
646        #
647        #It follows that the updated wvi is
648        #  wvi := (1-alpha)*(zvi+hc) + alpha*(zvi+hvi) =
649        #  zvi + hc + alpha*(hvi - hc)
650        #
651        #Note that hvi = zc+hc-zvi in the first order case (constant).
652
653        if alpha < 1:         
654            for i in range(3):
655                wv[k,i] = zv[k,i] + hc[k] + alpha*(hv[k,i]-hc[k])
656           
657
658            #Momentums at centroids
659            xmomc = domain.quantities['xmomentum'].centroid_values
660            ymomc = domain.quantities['ymomentum'].centroid_values       
661
662            #Momentums at vertices
663            xmomv = domain.quantities['xmomentum'].vertex_values
664            ymomv = domain.quantities['ymomentum'].vertex_values         
665
666            # Update momentum as a linear combination of
667            # xmomc and ymomc (shallow) and momentum
668            # from extrapolator xmomv and ymomv (deep).
669            xmomv[k,:] = (1-alpha)*xmomc[k] + alpha*xmomv[k,:]           
670            ymomv[k,:] = (1-alpha)*ymomc[k] + alpha*ymomv[k,:]
671                                   
672           
673
674def balance_deep_and_shallow_c(domain):
675    """Wrapper for C implementation
676    """
677   
678    #Shortcuts
679    wc = domain.quantities['level'].centroid_values
680    zc = domain.quantities['elevation'].centroid_values
681    hc = wc - zc
682   
683    wv = domain.quantities['level'].vertex_values
684    zv = domain.quantities['elevation'].vertex_values
685    hv = wv-zv
686
687    #Momentums at centroids
688    xmomc = domain.quantities['xmomentum'].centroid_values
689    ymomc = domain.quantities['ymomentum'].centroid_values       
690
691    #Momentums at vertices
692    xmomv = domain.quantities['xmomentum'].vertex_values
693    ymomv = domain.quantities['ymomentum'].vertex_values         
694
695   
696
697    from shallow_water_ext import balance_deep_and_shallow
698    balance_deep_and_shallow(wc, zc, hc, wv, zv, hv,
699                             xmomc, ymomc, xmomv, ymomv)
700
701   
702
703
704###############################################
705#Boundary - specific to the shallow water wave equation
706class Reflective_boundary(Boundary):
707    """Reflective boundary returns same conserved quantities as
708    those present in its neighbour volume but reflected.
709
710    This class is specific to the shallow water equation as it
711    works with the momentum quantities assumed to be the second
712    and third conserved quantities.
713    """
714   
715    def __init__(self, domain = None):       
716        Boundary.__init__(self)
717
718        if domain is None:
719            msg = 'Domain must be specified for reflective boundary'
720            raise msg
721       
722        #Handy shorthands
723        self.level = domain.quantities['level'].edge_values
724        self.xmom = domain.quantities['xmomentum'].edge_values
725        self.ymom = domain.quantities['ymomentum'].edge_values         
726        self.normals = domain.normals
727       
728        from Numeric import zeros, Float       
729        self.conserved_quantities = zeros(3, Float) 
730       
731    def __repr__(self):
732        return 'Reflective_boundary'
733
734
735    def evaluate(self, vol_id, edge_id):
736        """Reflective boundaries reverses the outward momentum
737        of the volume they serve.
738        """
739
740        q = self.conserved_quantities
741        q[0] = self.level[vol_id, edge_id]
742        q[1] = self.xmom[vol_id, edge_id] 
743        q[2] = self.ymom[vol_id, edge_id] 
744       
745        normal = self.normals[vol_id, 2*edge_id:2*edge_id+2]     
746
747       
748        r = rotate(q, normal, direction = 1)
749        r[1] = -r[1]
750        q = rotate(r, normal, direction = -1)
751
752        return q
753
754
755#########################
756#Standard forcing terms:
757#
758def gravity(domain):
759    """Apply gravitational pull in the presence of bed slope
760    """
761
762    from util import gradient
763    from Numeric import zeros, Float, array, sum
764
765    xmom = domain.quantities['xmomentum'].explicit_update
766    ymom = domain.quantities['ymomentum'].explicit_update
767
768    Level = domain.quantities['level']
769    Elevation = domain.quantities['elevation']   
770    h = Level.edge_values - Elevation.edge_values
771    v = Elevation.vertex_values
772   
773    x = domain.get_vertex_coordinates()
774    g = domain.g
775   
776    for k in range(domain.number_of_elements):   
777        avg_h = sum( h[k,:] )/3
778       
779        #Compute bed slope
780        x0, y0, x1, y1, x2, y2 = x[k,:]   
781        z0, z1, z2 = v[k,:]
782
783        zx, zy = gradient(x0, y0, x1, y1, x2, y2, z0, z1, z2)
784
785        #Update momentum
786        xmom[k] += -g*zx*avg_h
787        ymom[k] += -g*zy*avg_h       
788
789
790def gravity_c(domain):
791    """Wrapper calling C version
792    """   
793
794    xmom = domain.quantities['xmomentum'].explicit_update
795    ymom = domain.quantities['ymomentum'].explicit_update
796
797    Level = domain.quantities['level']   
798    Elevation = domain.quantities['elevation']   
799    h = Level.edge_values - Elevation.edge_values
800    v = Elevation.vertex_values
801   
802    x = domain.get_vertex_coordinates()
803    g = domain.g
804   
805   
806    from shallow_water_ext import gravity
807    gravity(g, h, v, x, xmom, ymom)
808   
809       
810def manning_friction(domain):
811    """Apply (Manning) friction to water momentum
812    """
813
814    from math import sqrt
815
816    w = domain.quantities['level'].centroid_values
817    z = domain.quantities['elevation'].centroid_values
818    h = w-z
819   
820    uh = domain.quantities['xmomentum'].centroid_values
821    vh = domain.quantities['ymomentum'].centroid_values
822    eta = domain.quantities['friction'].centroid_values   
823   
824    xmom_update = domain.quantities['xmomentum'].semi_implicit_update
825    ymom_update = domain.quantities['ymomentum'].semi_implicit_update   
826   
827    N = domain.number_of_elements
828    eps = domain.minimum_allowed_height
829    g = domain.g
830   
831    for k in range(N):
832        if eta[k] >= eps:
833            if h[k] >= eps:
834                S = -g * eta[k]**2 * sqrt((uh[k]**2 + vh[k]**2))
835                S /= h[k]**(7.0/3)
836
837                #Update momentum
838                xmom_update[k] += S*uh[k]
839                ymom_update[k] += S*vh[k]
840           
841       
842def manning_friction_c(domain):
843    """Wrapper for c version
844    """
845
846
847    xmom = domain.quantities['xmomentum']
848    ymom = domain.quantities['ymomentum']   
849       
850    w = domain.quantities['level'].centroid_values
851    z = domain.quantities['elevation'].centroid_values
852   
853    uh = xmom.centroid_values
854    vh = ymom.centroid_values
855    eta = domain.quantities['friction'].centroid_values   
856   
857    xmom_update = xmom.semi_implicit_update
858    ymom_update = ymom.semi_implicit_update   
859   
860    N = domain.number_of_elements
861    eps = domain.minimum_allowed_height
862    g = domain.g
863
864    from shallow_water_ext import manning_friction
865    manning_friction(g, eps, w, z, uh, vh, eta, xmom_update, ymom_update)
866
867
868def linear_friction(domain):
869    """Apply linear friction to water momentum
870
871    Assumes quantity: 'linear_friction' to be present
872    """
873
874    from math import sqrt
875
876    w = domain.quantities['level'].centroid_values
877    z = domain.quantities['elevation'].centroid_values
878    h = w-z
879
880    uh = domain.quantities['xmomentum'].centroid_values
881    vh = domain.quantities['ymomentum'].centroid_values
882    tau = domain.quantities['linear_friction'].centroid_values   
883   
884    xmom_update = domain.quantities['xmomentum'].semi_implicit_update
885    ymom_update = domain.quantities['ymomentum'].semi_implicit_update   
886   
887    N = domain.number_of_elements
888    eps = domain.minimum_allowed_height
889   
890    for k in range(N):
891        if tau[k] >= eps:
892            if h[k] >= eps:
893                S = -tau[k]/h[k]
894
895                #Update momentum
896                xmom_update[k] += S*uh[k]
897                ymom_update[k] += S*vh[k]
898           
899       
900       
901def check_forcefield(f):
902    """Check that f is either
903    1: a callable object of x,y,t, where x and y are vectors
904       and that it returns an array or a list of same length
905       as x and y
906    2: a scalar   
907    """
908
909    from Numeric import ones, Float, array
910
911   
912    if callable(f):
913        N = 3
914        x = ones(3, Float)
915        y = ones(3, Float)           
916        try:
917            q = f(x, y, 1.0)
918        except Exception, e:
919            msg = 'Function %s could not be executed:\n%s' %(f, e)
920            raise msg
921                   
922        try:
923            q = array(q).astype(Float)
924        except:
925            msg = 'Return value from vector function %s could ' %f
926            msg += 'not be converted into a Numeric array of floats.\n'
927            msg += 'Specified function should return either list or array.' 
928            raise msg
929
930        msg = 'Return vector from function %s' %f
931        msg += 'must have same lenght as input vectors'
932        assert len(q) == N, msg
933       
934    else:       
935        try:
936            f = float(f)
937        except:
938            msg = 'Force field %s must be either a scalar' %f
939            msg += ' or a vector function'
940            raise msg
941    return f   
942
943
944class Wind_stress:
945    """Apply wind stress to water momentum
946    """
947
948    def __init__(self, s, phi):
949        """Initialise windfield from wind speed s [m/s]
950        and wind direction phi [degrees]
951       
952        Inputs v and phi can be either scalars or Python functions, e.g.
953
954        W = Wind_stress(10, 178)
955
956        #FIXME - 'normal' degrees are assumed for now, i.e. the
957        vector (1,0) has zero degrees.
958        We may need to convert from 'compass' degrees later on and also
959        map from True north to grid north.
960
961        Arguments can also be Python functions of x,y,t as in
962
963        def speed(x,y,t):
964            ...
965            return s
966       
967        def angle(x,y,t):
968            ...
969            return phi
970
971        where x and y are vectors.
972
973        and then pass the functions in
974
975        W = Wind_stress(speed, angle)
976
977        The instantiated object W can be appended to the list of
978        forcing_terms as in
979
980        domain.forcing_terms.append(W)
981       
982        """
983
984        from config import rho_a, rho_w, eta_w
985        from Numeric import array, Float               
986
987        self.speed = check_forcefield(s)
988        self.phi = check_forcefield(phi)
989
990        self.const = eta_w*rho_a/rho_w
991       
992
993    def __call__(self, domain):
994        """Evaluate windfield based on values found in domain
995        """
996
997        from math import pi, cos, sin, sqrt
998        from Numeric import Float, ones, ArrayType
999       
1000        xmom_update = domain.quantities['xmomentum'].explicit_update
1001        ymom_update = domain.quantities['ymomentum'].explicit_update   
1002       
1003        N = domain.number_of_elements
1004        t = domain.time       
1005       
1006        if callable(self.speed):
1007            xc = domain.get_centroid_coordinates()           
1008            s_vec = self.speed(xc[:,0], xc[:,1], t)
1009        else:
1010            #Assume s is a scalar
1011
1012            try:
1013                s_vec = self.speed * ones(N, Float)
1014            except:
1015                msg = 'Speed must be either callable or a scalar: %s' %self.s
1016                raise msg
1017
1018
1019        if callable(self.phi):
1020            xc = domain.get_centroid_coordinates()           
1021            phi_vec = self.phi(xc[:,0], xc[:,1], t)
1022        else:
1023            #Assume phi is a scalar
1024
1025            try:
1026                phi_vec = self.phi * ones(N, Float)           
1027            except:
1028                msg = 'Angle must be either callable or a scalar: %s' %self.phi
1029                raise msg
1030
1031        assign_windfield_values(xmom_update, ymom_update,
1032                                s_vec, phi_vec, self.const)
1033
1034
1035def assign_windfield_values(xmom_update, ymom_update,
1036                            s_vec, phi_vec, const):
1037    """Python version of assigning wind field to update vectors.
1038    A c version also exists (for speed)
1039    """
1040    from math import pi, cos, sin, sqrt
1041   
1042    N = len(s_vec)
1043    for k in range(N):
1044        s = s_vec[k]
1045        phi = phi_vec[k]
1046
1047        #Convert to radians
1048        phi = phi*pi/180
1049
1050        #Compute velocity vector (u, v)
1051        u = s*cos(phi)
1052        v = s*sin(phi)
1053       
1054        #Compute wind stress
1055        S = const * sqrt(u**2 + v**2)
1056        xmom_update[k] += S*u
1057        ymom_update[k] += S*v       
1058           
1059
1060class Wind_stress_from_file:
1061    """Apply wind stress read from a file to water momentum
1062
1063    At this stage the wind field is assumed to depend on time only, i.e
1064    no spatial dependency.
1065   
1066    FIXME: This class may be incorporated in the generic Wind_stress class
1067    Superseded
1068    """
1069
1070   
1071    def __init__(self, filename):
1072        """Initialise windfield from a file with the following format
1073
1074        time [DD/MM/YY hh:mm:ss], speed [m/s] direction [degrees]
1075       
1076        e.g.
1077       
1078        03/09/04 19:15:00, 9.53 39
1079       
1080        domain.forcing_terms.append(W)
1081       
1082        """
1083
1084        import time
1085        from Numeric import array
1086        from config import time_format
1087
1088
1089        from config import rho_a, rho_w, eta_w
1090        self.const = eta_w*rho_a/rho_w
1091       
1092
1093        try:
1094            fid = open(filename)
1095        except Exception, e:
1096            msg = 'Wind stress file %s could not be opened: %s\n'\
1097                  %(filename, e)
1098            raise msg
1099
1100
1101        line = fid.readline()
1102        fid.close()
1103        fields = line.split(',')
1104        msg = 'File %s must have the format date, values'
1105        assert len(fields) == 2, msg
1106
1107        try:
1108            starttime = time.mktime(time.strptime(fields[0], time_format))
1109        except ValueError:   
1110            msg = 'First field in file %s must be' %filename
1111            msg += ' date-time with format %s.\n' %time_format
1112            msg += 'I got %s instead.' %fields[0] 
1113            raise msg
1114
1115        #Split values
1116        values = []
1117        for value in fields[1].split():
1118            values.append(float(value))
1119
1120        q = array(values)
1121       
1122        msg = 'ERROR: File function must return a 1d list or array '
1123        assert len(q.shape) == 1, msg
1124           
1125
1126        msg = 'Values specified in file must convert to an array of length 2'
1127        assert len(q) == 2, msg
1128
1129        self.filename = filename
1130        self.domain = domain
1131        self.starttime = starttime
1132
1133
1134        if domain.starttime is None:
1135            domain.starttime = starttime
1136        else:
1137            msg = 'Start time as specified in domain (%s) is earlier '
1138            'than the starttime of file %s: %s'\
1139                  %(domain.starttime, self.filename, self.starttime) 
1140            if self.starttime > domain.starttime:
1141                raise msg
1142
1143
1144        self.read_times() #Now read all times in.
1145
1146       
1147    def read_times(self):
1148        from Numeric import zeros, Float, alltrue
1149        from config import time_format
1150        import time
1151       
1152        fid = open(self.filename)       
1153        lines = fid.readlines()
1154        fid.close()
1155       
1156        N = len(lines)
1157
1158        T = zeros(N, Float)       #Time
1159        Q = zeros((N, 2), Float)  #Wind field (s, phi)
1160       
1161        for i, line in enumerate(lines):
1162            fields = line.split(',')
1163            real_time = time.mktime(time.strptime(fields[0], time_format))
1164
1165            T[i] = real_time - self.start_time
1166
1167            for j, value in enumerate(fields[1].split()):
1168                Q[i, j] = float(value)
1169
1170        msg = 'File %s must list time as a monotonuosly ' %self.filename
1171        msg += 'increasing sequence'
1172        assert alltrue( T[1:] - T[:-1] > 0 ), msg
1173
1174        self.T = T     #Time
1175        self.Q = Q     #Windfied
1176        self.index = 0 #Initial index
1177
1178       
1179    def __repr__(self):
1180        return 'Wind field from file'
1181
1182       
1183
1184    def __call__(self, domain):
1185        """Evaluate windfield based on values found in domain
1186        """
1187
1188        from math import pi, cos, sin, sqrt
1189       
1190        xmom_update = domain.quantities['xmomentum'].explicit_update
1191        ymom_update = domain.quantities['ymomentum'].explicit_update   
1192       
1193        N = domain.number_of_elements
1194        t = self.domain.time
1195       
1196        msg = 'Time given in File %s does not match model time'\
1197              %self.filename
1198       
1199        if t < self.T[0]: raise msg
1200        if t > self.T[-1]: raise msg       
1201
1202        oldindex = self.index
1203
1204        #Find slot
1205        while t > self.T[self.index]: self.index += 1
1206        while t < self.T[self.index]: self.index -= 1
1207
1208        #t is now between index and index+1
1209        ratio = (t - self.T[self.index])/\
1210                (self.T[self.index+1] - self.T[self.index])
1211
1212        #Compute interpolated values for s and phi
1213        q = self.Q[self.index,:] +\
1214            ratio*(self.Q[self.index+1,:] - self.Q[self.index,:]) 
1215
1216        s = q[0]
1217       
1218        #Convert to radians
1219        phi = q[1]*pi/180
1220       
1221        #Compute velocity vector (u, v)
1222        u = s*cos(phi)
1223        v = s*sin(phi)
1224
1225        #Compute wind stress for this time step
1226        S = self.const * sqrt(u**2 + v**2)
1227        xmom_update += S*u
1228        ymom_update += S*v       
1229       
1230
1231
1232
1233
1234###########################
1235###########################
1236#Geometries
1237
1238
1239#FIXME: Rethink this way of creating values.
1240
1241
1242class Weir:
1243    """Set a bathymetry for weir with a hole and a downstream gutter
1244    x,y are assumed to be in the unit square
1245    """
1246   
1247    def __init__(self, stage):
1248        self.inflow_stage = stage
1249
1250    def __call__(self, x, y):   
1251        from Numeric import zeros, Float
1252        from math import sqrt
1253   
1254        N = len(x)
1255        assert N == len(y)   
1256
1257        z = zeros(N, Float)
1258        for i in range(N):
1259            z[i] = -x[i]/2  #General slope
1260
1261            #Flattish bit to the left
1262            if x[i] < 0.3:
1263                z[i] = -x[i]/10 
1264           
1265            #Weir
1266            if x[i] >= 0.3 and x[i] < 0.4:
1267                z[i] = -x[i]+0.9
1268
1269            #Dip
1270            x0 = 0.6
1271            #depth = -1.3
1272            depth = -1.0
1273            #plateaux = -0.9
1274            plateaux = -0.6           
1275            if y[i] < 0.7:
1276                if x[i] > x0 and x[i] < 0.9:
1277                    z[i] = depth
1278
1279                #RHS plateaux
1280                if x[i] >= 0.9:
1281                    z[i] = plateaux
1282                   
1283                   
1284            elif y[i] >= 0.7 and y[i] < 1.5:
1285                #Restrict and deepen
1286                if x[i] >= x0 and x[i] < 0.8:
1287                    z[i] = depth-(y[i]/3-0.3)
1288                    #z[i] = depth-y[i]/5
1289                    #z[i] = depth
1290                elif x[i] >= 0.8:   
1291                    #RHS plateaux
1292                    z[i] = plateaux
1293                   
1294            elif y[i] >= 1.5:
1295                if x[i] >= x0 and x[i] < 0.8 + (y[i]-1.5)/1.2:
1296                    #Widen up and stay at constant depth
1297                    z[i] = depth-1.5/5
1298                elif x[i] >= 0.8 + (y[i]-1.5)/1.2:                       
1299                    #RHS plateaux
1300                    z[i] = plateaux                                   
1301                   
1302
1303            #Hole in weir (slightly higher than inflow condition)
1304            if x[i] >= 0.3 and x[i] < 0.4 and y[i] > 0.2 and y[i] < 0.4:
1305                z[i] = -x[i]+self.inflow_stage + 0.02
1306
1307            #Channel behind weir
1308            x0 = 0.5
1309            if x[i] >= 0.4 and x[i] < x0 and y[i] > 0.2 and y[i] < 0.4:
1310                z[i] = -x[i]+self.inflow_stage + 0.02
1311               
1312            if x[i] >= x0 and x[i] < 0.6 and y[i] > 0.2 and y[i] < 0.4:
1313                #Flatten it out towards the end
1314                z[i] = -x0+self.inflow_stage + 0.02 + (x0-x[i])/5
1315
1316            #Hole to the east
1317            x0 = 1.1; y0 = 0.35
1318            #if x[i] < -0.2 and y < 0.5:
1319            if sqrt((2*(x[i]-x0))**2 + (2*(y[i]-y0))**2) < 0.2:
1320                z[i] = sqrt(((x[i]-x0))**2 + ((y[i]-y0))**2)-1.0
1321
1322            #Tiny channel draining hole
1323            if x[i] >= 1.14 and x[i] < 1.2 and y[i] >= 0.4 and y[i] < 0.6:
1324                z[i] = -0.9 #North south
1325               
1326            if x[i] >= 0.9 and x[i] < 1.18 and y[i] >= 0.58 and y[i] < 0.65:
1327                z[i] = -1.0 + (x[i]-0.9)/3 #East west
1328           
1329           
1330
1331            #Stuff not in use
1332           
1333            #Upward slope at inlet to the north west
1334            #if x[i] < 0.0: # and y[i] > 0.5:
1335            #    #z[i] = -y[i]+0.5  #-x[i]/2
1336            #    z[i] = x[i]/4 - y[i]**2 + 0.5
1337
1338            #Hole to the west
1339            #x0 = -0.4; y0 = 0.35 # center
1340            #if sqrt((2*(x[i]-x0))**2 + (2*(y[i]-y0))**2) < 0.2:
1341            #    z[i] = sqrt(((x[i]-x0))**2 + ((y[i]-y0))**2)-0.2
1342
1343
1344
1345
1346
1347        return z/2
1348
1349class Weir_simple:
1350    """Set a bathymetry for weir with a hole and a downstream gutter
1351    x,y are assumed to be in the unit square
1352    """
1353   
1354    def __init__(self, stage):
1355        self.inflow_stage = stage
1356
1357    def __call__(self, x, y):   
1358        from Numeric import zeros, Float
1359   
1360        N = len(x)
1361        assert N == len(y)   
1362
1363        z = zeros(N, Float)
1364        for i in range(N):
1365            z[i] = -x[i]  #General slope
1366
1367            #Flat bit to the left
1368            if x[i] < 0.3:
1369                z[i] = -x[i]/10  #General slope
1370           
1371            #Weir
1372            if x[i] > 0.3 and x[i] < 0.4:
1373                z[i] = -x[i]+0.9
1374
1375            #Dip
1376            if x[i] > 0.6 and x[i] < 0.9:
1377                z[i] = -x[i]-0.5  #-y[i]/5
1378
1379            #Hole in weir (slightly higher than inflow condition)
1380            if x[i] > 0.3 and x[i] < 0.4 and y[i] > 0.2 and y[i] < 0.4:
1381                z[i] = -x[i]+self.inflow_stage + 0.05       
1382           
1383
1384        return z/2
1385       
1386
1387           
1388class Constant_height:
1389    """Set an initial condition with constant water height, e.g
1390    stage s = z+h
1391    """
1392    def __init__(self, W, h):
1393        self.W = W
1394        self.h = h
1395
1396    def __call__(self, x, y):
1397        if self.W is None:
1398            from Numeric import ones, Float
1399            return self.h*ones(len(x), Float)
1400        else:
1401            return self.W(x,y) + self.h
1402
1403
1404
1405##############################################
1406#Initialise module
1407
1408
1409import compile
1410if compile.can_use_C_extension('shallow_water_ext.c'):
1411    #Replace python version with c implementations
1412   
1413    from shallow_water_ext import rotate, assign_windfield_values
1414    compute_fluxes = compute_fluxes_c
1415    gravity = gravity_c
1416    manning_friction = manning_friction_c
1417    balance_deep_and_shallow = balance_deep_and_shallow_c
1418    protect_against_infinitesimal_and_negative_heights = protect_against_infinitesimal_and_negative_heights_c
1419
1420   
1421    #distribute_to_vertices_and_edges = distribute_to_vertices_and_edges_c
1422
1423
1424#Optimisation with psyco
1425from config import use_psyco
1426if use_psyco:
1427    try:
1428        import psyco
1429    except:
1430        msg = 'WARNING: psyco (speedup) could not import'+\
1431              ', you may want to consider installing it'
1432        print msg
1433    else:
1434        psyco.bind(Domain.distribute_to_vertices_and_edges)
1435        psyco.bind(Domain.compute_fluxes)
1436       
1437if __name__ == "__main__":
1438    pass
Note: See TracBrowser for help on using the repository browser.