source: codes/icosagcm/devel/Python/dynamico/meshes.py @ 984

Last change on this file since 984 was 984, checked in by dubos, 5 years ago

devel/Python : bugfix dynamico/meshes.py

File size: 41.8 KB
Line 
1from __future__ import absolute_import     
2
3import time
4import math
5import numpy as np
6import netCDF4 as cdf
7import matplotlib.tri as tri
8import matplotlib.pyplot as plt
9from matplotlib.patches import Polygon
10from matplotlib.collections import PatchCollection
11
12from dynamico import parallel
13from dynamico.util import list_stencil, BaseClass, inverse_list
14from dynamico import getargs
15from dynamico.partition import partition_mesh
16
17log_master, log_world = getargs.getLogger(__name__)
18INFO, DEBUG, ERROR = log_master.info, log_master.debug, log_world.error
19INFO_ALL, DEBUG_ALL = log_world.info, log_world.debug
20
21radian=180/math.pi # convert from radians to degrees
22
23#------------------- Hybrid mass-based coordinate -------------
24
25# compute hybrid coefs from distribution of mass
26
27def compute_hybrid_coefs(mass):
28    if mass.ndim==2 :
29        nx,llm=mass.shape
30        mass_dak = np.zeros((nx,llm))
31        mass_dbk = np.zeros((nx,llm))
32        mass_bl = np.zeros((nx,llm+1))+1.
33        for i in range(nx):
34            m_i = mass[i,:]
35            dbk_i = m_i/sum(m_i)
36            mass_dbk[i,:] = dbk_i
37            mass_bl[i,1:]= 1.-dbk_i.cumsum()
38    if mass.ndim==3 :
39        nx,ny,llm=mass.shape
40        mass_dak = np.zeros((nx,ny,llm))
41        mass_dbk = np.zeros((nx,ny,llm))
42        mass_bl = np.zeros((nx,ny,llm+1))+1.
43        for i in range(nx):
44            for j in range(ny):
45                m_ij = mass[i,j,:]
46                dbk_ij = m_ij/sum(m_ij)
47                mass_dbk[i,j,:] = dbk_ij
48                mass_bl[i,j,1:]= 1.-dbk_ij.cumsum()
49 
50    return mass_bl, mass_dak, mass_dbk
51
52# from theta.k90       : rhodz(l,ij) = mass_dak(l,ij) + mass_col(ij)*mass_dbk(l,ij)
53# from caldyn_wflux.k90: wflux(l,ij) = mass_bl(l,ij) * dmass_col(ij) - convm(l,ij)
54def mass_coefs_from_pressure_coefs(g, ap_il, bp_il):
55    # p = Mg + ptop = a + b.ps
56    # ps = Mcol.g+ptop
57    # M = mass_a + mass_b.Mcol
58    # M = (a+b.ps-ptop)/g
59    #   = (a+b.ptop-ptop)/g + b.Mcol
60    # mass_a  = (a+b.ptop)/g
61    # mass_da = (da+db.ptop)/g
62    # mass_b  = b
63    # mass_db = db
64    n, llm, ptop = ap_il.shape[0], ap_il.shape[1]-1, ap_il[0,-1]
65    DEBUG('hybrid ptop : %f' % ptop) 
66    mass_al = (ap_il+ptop*bp_il)/g
67    mass_dak, mass_dbk = np.zeros((n,llm)), np.zeros((n,llm))
68    for k in range(llm):
69        mass_dak[:,k] = mass_al[:,k]-mass_al[:,k+1]
70        mass_dbk[:,k] = bp_il[:,k]-bp_il[:,k+1]
71    DEBUG('%s'%mass_al[0,:])
72    DEBUG('%s'%mass_dak[0,:])
73    DEBUG('%s'%mass_dbk[0,:])
74    return [ np.asarray(x) for x in bp_il, mass_dak, mass_dbk ]
75
76#----------------------- Base class for "user" meshes ---------
77
78class BaseMesh(BaseClass):
79    def zeros(self,dims): return np.zeros([n for n in dims if n>1]) # overriden by meshes.dev.DevMesh
80
81#----------------------- Cartesian mesh -----------------------
82
83# arrays is a list of arrays
84# vals is a list of tuples
85# each tuple is stored in each array
86def put(ij, deg, arrays, vals):
87    k=0
88    for vv in vals: # vv is a tuple of values to be stored in arrays
89        for array,v in zip(arrays,vv):
90            array[ij,k]=v
91        k=k+1       
92    deg[ij]=k
93
94def concat(x,y):
95    z = np.asarray([x,y])
96    return z.transpose()
97
98class Cartesian_Mesh(BaseMesh):
99    def __init__(self,nx,ny,llm,nqdyn,Lx,Ly,f):
100        dx,dy = Lx/nx, Ly/ny
101        self.dx, self.dy, self.f = dx,dy,f
102        self.nx, self.ny, self.llm, self.nqdyn = nx,ny,llm,nqdyn
103        self.field_z = self.field_mass
104        # 1D coordinate arrays
105        x=(np.arange(nx)-nx/2.)*Lx/nx
106        y=(np.arange(ny)-ny/2.)*Ly/ny
107        lev=np.arange(llm)
108        levp1=np.arange(llm+1)
109        self.x, self.y, self.lev, self.levp1 = x,y,lev,levp1
110        # 3D coordinate arrays
111        self.yy,self.xx,self.ll = np.meshgrid(y,x,lev, indexing='ij')       
112        self.yyp1,self.xxp1,self.llp1 = np.meshgrid(y,x,levp1, indexing='ij')
113        # beware conventions for indexing
114        # Fortran order : llm,nx*ny,nqdyn  / indices start at 1
115        # Python order : nqdyn,ny,nx,llm   / indices start at 0
116        # indices below follow Fortran while x,y follow Python/C
117        def periodize(z,nz): return (z+2*nz)%nz
118        def index(x,y):      return 1+periodize(x,nx)+nx*periodize(y,ny)
119        def indexu(x,y):     return 2*index(x,y)-1
120        def indexv(x,y):     return 2*index(x,y)
121        def zeros(shape,n=1): return [np.zeros(shape) for i in range(n)] if n>1 else np.zeros(shape)
122        def indices(shape,n=1): return [np.zeros(shape,np.int32) for i in range(n)] if n>1 else np.zeros(shape,np.int32)
123
124        primal_num, dual_num, edge_num = nx*ny, nx*ny, 2*nx*ny
125
126        Aiv, lon_i,lon_v,lat_i,lat_v   = zeros( nx*ny,   5)
127        lon_e,lat_e,le,de,angle_e      = zeros( 2*nx*ny, 5)
128        Riv2                           = zeros((nx*ny,4)  )
129        wee                            = zeros((2*nx*ny,4))
130        primal_bounds_lon, primal_bounds_lat, dual_bounds_lon, dual_bounds_lat = zeros( (nx*ny,4), 4)
131
132        primal_deg,dual_deg                 = indices( nx*ny,   2)
133        primal_edge,primal_ne,primal_vertex = indices((nx*ny,4),3)
134        dual_edge,dual_ne,dual_vertex       = indices((nx*ny,4),3)
135        trisk_deg,left,right,up,down        = indices( 2*nx*ny, 5)
136        trisk                               = indices((2*nx*ny,4)) # 4 TRiSK coefs per edge
137   
138        def putval(ij, arrays, vals): 
139            for array,v in zip(arrays,vals): array[ij]=v
140        for x in range(nx):
141            for y in range(ny):
142                # NB : Fortran indices start at 1
143                # primal cells
144                ij=index(x,y)-1
145                lon_i[ij], lat_i[ij] = x, y
146                put(ij, primal_deg,(primal_edge,primal_vertex,primal_ne, primal_bounds_lon, primal_bounds_lat),
147                    ((indexu(x,y),   index(x,y),      1, x+.5, y+.5), 
148                     (indexv(x,y),   index(x-1,y),    1, x-.5, y+.5),
149                     (indexu(x-1,y), index(x-1,y-1), -1, x-.5, y-.5),
150                     (indexv(x,y-1), index(x,y-1),   -1, x+.5, y-.5) ))
151                # dual cells
152                ij=index(x,y)-1
153                lon_v[ij], lat_v[ij] = x+.5, y+.5
154                put(ij,dual_deg,(dual_edge,dual_vertex,dual_ne,Riv2, dual_bounds_lon, dual_bounds_lat),
155                    ((indexv(x+1,y),index(x,y),      1, .25, x,    y),
156                     (indexu(x,y+1),index(x+1,y),   -1, .25, x+1., y),
157                     (indexv(x,y),  index(x+1,y+1), -1, .25, x+1., y+1.),
158                     (indexu(x,y),  index(x,y+1),    1, .25, x,    y+1.) ))
159                # edges :
160                # left and right are adjacent primal cells
161                # flux is positive when going from left to right
162                # up and down are adjacent dual cells (vertices)
163                # circulation is positive when going from down to up
164                # u-points
165                ij =indexu(x,y)-1 
166                putval(ij, (left, right, down, up, le, de, angle_e, lon_e, lat_e),
167                       (index(x,y), index(x+1,y), index(x,y-1), index(x,y), dy, dx, 0., x+.5, y ) )
168                put(ij,trisk_deg,(trisk,wee),(
169                    (indexv(x,y),.25),
170                    (indexv(x+1,y),.25),
171                    (indexv(x,y-1),.25),
172                    (indexv(x+1,y-1),.25)))
173                # v-points
174                ij = indexv(x,y)-1
175                putval(ij, (left, right, down, up, le, de, angle_e, lon_e, lat_e),
176                       (index(x,y), index(x,y+1), index(x,y), index(x-1,y), dx, dy, .5*math.pi, x, y+.5 ) )
177                put(ij,trisk_deg,(trisk,wee),(
178                    (indexu(x,y),-.25),
179                    (indexu(x-1,y),-.25),
180                    (indexu(x,y+1),-.25),
181                    (indexu(x-1,y+1),-.25)))
182
183        primal_i, primal_j = [ x.astype(np.int32) for x in lon_i, lat_i]
184        edge_i, edge_j = [ x.astype(np.int32) for x in lon_e, lat_e]
185        lon_i, lon_v, lon_e = [(x+.5)*dx-.5*Lx for x in lon_i, lon_v, lon_e]
186        lat_i, lat_v, lat_e = [(y+.5)*dy-.5*Ly for y in lat_i, lat_v, lat_e]
187
188        Aiv[:]=dx*dy # Ai=Av=dx*dy
189        le_de, fv = le/de, f+0.*Aiv
190       
191        self.set_members(locals(), 'primal_num',  'dual_num', 'edge_num',
192                         'primal_deg', 'primal_edge', 'primal_ne', 'primal_vertex',
193                         'primal_bounds_lon', 'primal_bounds_lat',
194                         'dual_deg', 'dual_edge', 'dual_ne', 'dual_vertex',
195                         'dual_bounds_lon', 'dual_bounds_lat',
196                         'left','right','down','up',
197                         'trisk_deg','trisk','Riv2','wee',
198                         'lon_i','lat_i','lon_v','lat_v', 'fv',
199                         'le_de','le','de','lon_e','lat_e','angle_e',
200                         'primal_i','primal_j','edge_i','edge_j')
201        self.Ai, self.Av = Aiv, Aiv
202        self.to_dynamico()
203    def to_dynamico(self): pass
204    def ncwrite(self, name):
205        """The method writes Cartesian mesh on the disc.
206            Args: Mesh parameters"""
207
208        #----------------OPENING NETCDF OUTPUT FILE------------
209
210        try:
211            f = cdf.Dataset(name, 'w', format='NETCDF4')
212        except:
213            ERROR("Cartesian_mesh.ncwrite : Error occurred while opening new netCDF mesh file")
214            raise
215
216        #----------------DEFINING DIMENSIONS--------------------
217
218        for dimname, dimsize in [("primal_cell", self.primal_deg.size),
219                                 ("dual_cell", self.dual_deg.size),
220                                 ("edge", self.left.size),
221                                 ("primal_edge_or_vertex", self.primal_edge.shape[1]),
222                                 ("dual_edge_or_vertex", self.dual_edge.shape[1]),
223                                 ("TWO", 2),
224                                 ("trisk_edge", 4)]:
225            f.createDimension(dimname,dimsize)
226
227        #----------------DEFINING VARIABLES---------------------
228       
229        f.description = "Cartesian_mesh"
230        f.meshtype, f.nx, f.ny = 'curvilinear', self.nx, self.ny
231
232        def create_vars(dimname, info):
233            for vname, vtype, vdata in info:
234                var = f.createVariable(vname,vtype,dimname)
235                var[:] = vdata
236               
237        create_vars("primal_cell", 
238                    [("primal_deg","i4",self.primal_deg), 
239                     ("primal_i","i4",self.primal_i),
240                     ("primal_j","i4",self.primal_j),
241                     ("Ai","f8",self.Ai),
242                     ("lon_i","f8",self.lon_i),
243                     ("lat_i","f8",self.lat_i)] )
244
245        create_vars("dual_cell", 
246                    [("dual_deg","i4",self.dual_deg), 
247                     ("Av","f8",self.Av),
248                     ("lon_v","f8",self.lon_v),
249                     ("lat_v","f8",self.lat_v),
250                     ] )
251
252        create_vars( ("primal_cell","primal_edge_or_vertex"), 
253                     [("primal_edge", "i4", self.primal_edge),
254                      ("primal_ne", "i4", self.primal_ne),
255                      ("primal_vertex", "i4", self.primal_vertex),
256                      ("primal_bounds_lon", "f8", self.primal_bounds_lon),
257                      ("primal_bounds_lat", "f8", self.primal_bounds_lat)] )
258
259        create_vars( ("dual_cell","dual_edge_or_vertex"), 
260                     [("dual_edge", "i4", self.dual_edge),
261                      ("dual_vertex","i4",self.dual_vertex),
262                      ("dual_bounds_lon", "f8", self.dual_bounds_lon),
263                      ("dual_bounds_lat", "f8", self.dual_bounds_lat),
264                      ("dual_ne","i4",self.dual_ne),
265                      ("Riv2","f8",self.Riv2)] )
266
267        create_vars("edge",
268                    [("trisk_deg","i4",self.trisk_deg),
269                     ("le","f8",self.le),
270                     ("de","f8",self.de),
271                     ("lon_e","f8", self.lon_e),
272                     ("lat_e","f8", self.lat_e),
273                     ("angle_e","f8", self.angle_e),
274                     ("edge_i","i4",self.edge_i),
275                     ("edge_j","i4",self.edge_j)] )
276
277        create_vars( ("edge","TWO"),
278                     [("edge_left_right","i4", concat(self.left,self.right)),
279                      ("edge_down_up","i4", concat(self.down,self.up))] )
280                     
281
282        create_vars( ("edge","trisk_edge"),
283                    [("trisk","i4",self.trisk),
284                     ("wee","f8",self.wee)] )
285
286        f.close()
287
288    def field_theta(self,n=1): return self.zeros((n,self.nqdyn,self.ny,self.nx,self.llm))
289    def field_mass(self,n=1): return self.zeros((n,self.ny,self.nx,self.llm))
290    def field_z(self,n=1): return self.zeros((n,self.ny,self.nx,self.llm))
291    def field_w(self,n=1): return self.zeros((n,self.ny,self.nx,self.llm+1))
292    def field_u(self,n=1): return self.zeros((n,self.ny,2*self.nx,self.llm))
293    def field_ps(self,n=1): return self.zeros((n,self.ny,self.nx))
294    def field_ucomp(self,n=1): return self.zeros((n,self.ny,self.nx,self.llm))
295    def ucomp(self,u): 
296        return u[range(0,2*self.nx,2),:] if self.ny==1 else u[:,range(0,2*self.nx,2),:]
297    def set_ucomp(self,uv,u):
298        if self.ny==1 : uv[range(0,2*self.nx,2),:]=u
299        else : uv[:,range(0,2*self.nx,2),:]=u
300    def vcomp(self,u):
301        return u[range(1,2*self.nx,2),:] if self.ny==1 else u[:,range(1,2*self.nx,2),:]
302    def set_vcomp(self,uv,v):
303        if self.ny==1 : uv[range(1,2*self.nx,2),:]=v
304        else : uv[:,range(1,2*self.nx,2),:]=v
305
306#---------------------- DYNAMICO format for fully unstructured mesh and curvilinear meshes ----------------------
307
308class Mesh_Format(BaseMesh):
309    def getdims(self,*names): return [len(self.nc.dimensions[name]) for name in names]
310    def get_pdims(self,comm,*names): return [self.get_pdim(comm,name) for name in names]
311    def get_pvars(self,pdim,*names): return [self.get_pvar(pdim,name) for name in names]
312    def getvars(self,*names): 
313            for name in names : DEBUG("getvar %s ..."%name)
314            time1=time.time()
315            ret=[self.getvar(name) for name in names]
316            DEBUG("... Done in %f seconds"%(time.time()-time1))
317            return ret
318   
319class DYNAMICO_Format(Mesh_Format):
320    """ Helps class Unstructured_Mesh to read a DYNAMICO mesh file."""
321    start_index=1 # indexing starts at 1 as in Fortran
322    def __init__(self,gridfilename):
323        nc = cdf.Dataset(gridfilename, "r")
324        self.nc, self.meshtype = nc, nc.meshtype
325        if self.meshtype == 'curvilinear':
326            self.nx, self.ny = nc.nx, nc.ny
327    def get_pdim(self,comm,name): return parallel.PDim(self.nc.dimensions[name], comm)
328    def getvar(self,name): return self.nc.variables[name][:]
329    def get_pvar(self,pdim,name):
330        # provides parallel access to a NetCDF variable, with name translation
331        # first deal with special cases
332        if name == 'edge_deg': return parallel.CstPArray1D(pdim, np.int32, 2)
333        # general case
334        return parallel.PArray(pdim, self.nc.variables[name])
335    def normalize(self, mesh): pass
336
337#---------------------- MPAS fully unstructured mesh ------------------------
338
339class MPAS_Format(Mesh_Format):
340    """ Helps class Unstructured_Mesh to read a MPAS mesh file."""
341    start_index=1 # indexing starts at 1 as in Fortran
342    meshtype='unstructured'
343    translate= {
344        'primal_cell':'nCells', 'edge':'nEdges', 'dual_cell':'nVertices',
345        'primal_deg':'nEdgesOnCell', 'trisk_deg':'nEdgesOnEdge',
346        'primal_edge':'edgesOnCell', 'dual_edge':'edgesOnVertex', 
347        'primal_vertex':'verticesOnCell', 'dual_vertex':'cellsOnVertex',
348        'trisk':'edgesOnEdge', 'edge_down_up':'verticesOnEdge', 'edge_left_right':'cellsOnEdge',
349        'le':'dvEdge', 'de':'dcEdge', 'Ai':'areaCell', 'Av':'areaTriangle',
350        'lat_i':'latCell','lon_i':'lonCell','lat_v':'latVertex','lon_v':'lonVertex',
351        'lat_e':'latEdge','lon_e':'lonEdge','angle_e':'angleEdge',
352        'wee':'weightsOnEdge','Riv2':'kiteAreasOnVertex',
353        'primal_ne':'signOnCell', 'dual_ne':'signOnVertex',
354        'primal_bounds_lon':None, 'primal_bounds_lat':None,
355        'dual_bounds_lon':None, 'dual_bounds_lat':None }
356    def __init__(self,gridfilename):
357        try:
358            self.nc = cdf.Dataset(gridfilename, "r")
359        except RuntimeError:
360            ERROR("""
361Unable to open grid file %s, maybe you forgot to download it ?
362To do so, go to the 'Python/' dir and execute './get_MPAS_grids.sh'.
363            """ % gridfile)
364            raise
365    def get_pdim(self,comm,name): 
366        name = self.translate[name]
367        return parallel.PDim(self.nc.dimensions[name], comm)
368    def getvar(self,name):
369        # first deal with special cases
370        if name == 'dual_deg':           
371            dual_deg = [3 for i in range(self.dual_num) ]
372            dual_deg = np.ascontiguousarray(dual_deg,dtype=np.int32) 
373            # NB : Fortran code expects 32-bit ints
374            return dual_deg
375        # general case
376        name=self.translate[name]
377        return self.nc.variables[name][:]
378    def get_pvar(self,pdim,name): 
379        # provides parallel access to a NetCDF variable, with name translation
380        # first deal with special cases
381        if name == 'dual_deg': return parallel.CstPArray1D(pdim, np.int32, 3)
382        if name == 'edge_deg': return parallel.CstPArray1D(pdim, np.int32, 2)
383        # general case
384        name=self.translate[name]
385        return None if name is None else parallel.PArray(pdim, self.nc.variables[name])
386    def normalize(self, mesh):
387        (trisk_deg, trisk, Ai,
388         de, le, wee, Av, Riv2) = (mesh.trisk_deg, mesh.trisk, mesh.Ai, 
389                                       mesh.de, mesh.le, mesh.wee, mesh.Av, mesh.Riv2)
390        edge_num, dual_num = de.size, Av.size
391        # fix normalization of wee and Riv2 weights
392        for edge1 in range(edge_num):
393            for i in range(trisk_deg[edge1]):
394                edge2=trisk[edge1,i]-1  # NB Fortran vs C/Python indexing
395                wee[edge1,i] = de[edge1]*wee[edge1,i]/le[edge2]
396        for ivertex in range(dual_num):
397            for j in range(3):
398                Riv2[ivertex,j]=Riv2[ivertex,j]/Av[ivertex]
399            r=Riv2[ivertex,:]
400            r=sum(r)
401            if abs(r-1.)>1e-6 : WARNING('error with Riv2 at vertex %d'%ivertex)
402        # scale lengths and areas
403        # this is now done by apply_map
404        # mesh.le, mesh.de, mesh.Av, mesh.Ai = radius*le, radius*de, r2*Av, r2*Ai
405    def get_bounds(self, mesh):
406        def bounds(deg, vertex, lon, lat):
407            num, maxdeg = deg.size, deg.max()
408            vertex_lon, vertex_lat = np.zeros((num, maxdeg)), np.zeros((num, maxdeg))
409            for cell in range(num):
410                cdeg=deg[cell]
411                for ivertex in range(maxdeg):
412                    if ivertex<cdeg:
413                        vtx = vertex[cell,ivertex] - 1
414                        vlon, vlat = lon[vtx], lat[vtx]
415                    vertex_lon[cell, ivertex] = vlon
416                    vertex_lat[cell, ivertex] = vlat
417            return vertex_lon, vertex_lat
418        mesh.primal_bounds_lon, mesh.primal_bounds_lat = bounds(mesh.primal_own_deg, mesh.primal_vertex, mesh.lon_v, mesh.lat_v)
419        mesh.dual_bounds_lon, mesh.dual_bounds_lat = bounds(mesh.dual_deg, mesh.dual_vertex, mesh.lon_i, mesh.lat_i)
420
421def compute_ne(num,deg,edges,left,right):
422    ne = np.zeros_like(edges)
423    for cell in range(num):
424        for iedge in range(deg[cell]):
425            edge = edges[cell,iedge]-1
426            if left[edge]==cell+1: ne[cell,iedge]+=1
427            if right[edge]==cell+1: ne[cell,iedge]-=1
428#            if ne[cell,iedge]==0 : print 'compute_ne error at cell,iedge', cell, iedge
429    return ne
430
431class Abstract_Mesh(BaseMesh):
432    def to_dynamico(self): pass
433    def field_theta(self,n=1): return self.zeros((n,self.nqdyn,self.primal_num,self.llm))
434    def field_mass(self,n=1):  return self.zeros((n,self.primal_num,self.llm))
435    def field_z(self,n=1):     return self.zeros((n,self.dual_num,self.llm))
436    def field_w(self,n=1):     return self.zeros((n,self.primal_num,self.llm+1))
437    def field_u(self,n=1):     return self.zeros((n,self.edge_num,self.llm))
438    def field_ps(self,n=1):    return self.zeros((n,self.primal_num,))
439    def ucov2D(self, ulon, ulat): 
440        return self.de*(ulon*np.cos(self.angle_e)+ulat*np.sin(self.angle_e))
441    def ucov3D(self, ulon, ulat):
442        ucov = self.zeros((self.edge_num,ulon.shape[1]))
443        for edge in range(self.edge_num):
444            angle=self.angle_e[edge]
445            ucov[edge,:] = self.de[edge]*(ulon[edge,:]*math.cos(angle)+ulat[edge,:]*math.sin(angle))
446        return ucov
447    def plot(self, tri,data, **kwargs):
448        plt.figure(figsize=(10,4))
449        plt.gca().set_aspect('equal')
450        plt.tricontourf(tri, data, 20, **kwargs)
451        plt.colorbar()
452        plt.ylim((-90,90))
453        plt.xlim((-180,180))
454    def plot_i(self,data, **kwargs):
455        self.plot(self.primal,data, **kwargs)
456    def plot_v(self,data, **kwargs):
457        self.plot(self.dual,data, **kwargs)
458    def plot_e(self,data, **kwargs):
459        self.plot(self.triedge,data, **kwargs)
460    def apply_map(self, mapping):
461        """Before calling apply_map, lat and lon are coordinates in the reference domain.
462        After calling apply_map, lat and lon are coordinates in the physical domain."""
463        lat_i, lat_e, lat_v = self.lat_i, self.lat_e, self.lat_v
464        lon_i, lon_e, lon_v = self.lon_i, self.lon_e, self.lon_v
465        factor_e = mapping.map_factor(lon_e, lat_e)
466        factor2_i = mapping.map_factor(lon_i, lat_i)**2
467        factor2_v = mapping.map_factor(lon_v, lat_v)**2
468        self.le, self.de, self.Av, self.Ai = self.le*factor_e, self.de*factor_e, self.Av*factor2_v, self.Ai*factor2_i
469        self.angle_e += mapping.map_angle(lon_e, lat_e)
470        self.lon_i, self.lat_i = mapping.map(lon_i, lat_i)
471        self.lon_v, self.lat_v = mapping.map(lon_v, lat_v)
472        self.lon_e, self.lat_e = mapping.map(lon_e, lat_e)
473        self.primal_bounds_lon, self.primal_bounds_lat = mapping.map( self.primal_bounds_lon, self.primal_bounds_lat)
474        self.dual_bounds_lon, self.dual_bounds_lat = mapping.map( self.dual_bounds_lon, self.dual_bounds_lat)
475        self.ref_lon_i, self.ref_lat_i = lon_i, lat_i
476        self.ref_lon_e, self.ref_lat_e = lon_e, lat_e
477        self.ref_lon_v, self.ref_lat_v = lon_v, lat_v
478
479class Unstructured_Mesh(Abstract_Mesh):
480    def __init__(self, gridfile, llm, nqdyn, radius, f):
481        getvars = gridfile.getvars
482        # get positions, lengths, surfaces and weights
483        le,de,Ai,Av,wee,Riv2 = getvars('le','de','Ai','Av','wee','Riv2')
484        le_de = le/de
485        lat_i,lon_i,lat_v,lon_v = getvars('lat_i','lon_i','lat_v','lon_v')
486        lat_e,lon_e,angle_e = getvars('lat_e','lon_e','angle_e')
487        fv = f(lon_v,lat_v) # Coriolis parameter
488        primal_num, dual_num, edge_num = [ x.size for x in Ai, Av, le ]
489        gridfile.primal_num, gridfile.dual_num = primal_num, dual_num
490        primal_deg, dual_deg, trisk_deg = getvars('primal_deg','dual_deg','trisk_deg')       
491        # get indices for stencils
492        # primal <-> edges
493        primal_edge, left_right = getvars('primal_edge','edge_left_right')
494        left,right = left_right[:,0], left_right[:,1]
495        primal_ne = compute_ne(primal_num,primal_deg,primal_edge,left,right)
496        # dual <-> edges
497        dual_edge, down_up = getvars('dual_edge','edge_down_up')
498        down,up = down_up[:,0], down_up[:,1]
499        dual_ne = compute_ne(dual_num,dual_deg,dual_edge,up,down)
500        # primal <-> dual, edges <-> edges (TRiSK)
501        primal_vertex, dual_vertex, trisk = getvars('primal_vertex','dual_vertex','trisk')
502        # CRITICAL : force arrays left, etc. to be contiguous in memory:
503        left,right,up,down = [np.ascontiguousarray(x) for x in (left,right,up,down)]
504        trisk,wee,primal_edge = [np.ascontiguousarray(x) for x in (trisk,wee,primal_edge)]
505        self.set_members(locals(), 'gridfile','llm','nqdyn','radius','f',
506                         'primal_num', 'dual_num', 'edge_num', 'primal_deg', 'primal_vertex', 'primal_ne',
507                         'lat_i', 'lon_i', 'lat_v', 'lon_v', 'lat_e', 'lon_e', 'angle_e',
508                         'Ai', 'Av', 'fv', 'le', 'de', 'down', 'up', 'le_de',
509                         'trisk_deg', 'trisk', 'wee', 'Riv2', 'dual_ne',
510                         'left', 'right',
511                         'primal_edge')
512        gridfile.normalize(self)
513        self.to_dynamico()
514        self.primal  = tri.Triangulation(lon_i*radian, lat_i*radian)
515        self.dual    = tri.Triangulation(lon_v*radian, lat_v*radian)
516        self.triedge = tri.Triangulation(lon_e*radian, lat_e*radian)       
517       
518class Unstructured_PMesh(Abstract_Mesh): # Mesh data distributed across MPI processes
519    def __init__(self,comm, gridfile):
520        self.comm, self.gridfile, get_pvars = comm, gridfile, gridfile.get_pvars
521        dim_primal, dim_edge, dim_dual = gridfile.get_pdims(comm,
522            'primal_cell','edge','dual_cell')
523        primal_deg, primal_vertex, primal_edge, primal_ne, lon_i, lat_i, Ai = get_pvars(
524            dim_primal, 'primal_deg', 'primal_vertex', 'primal_edge', 'primal_ne', 'lon_i', 'lat_i', 'Ai' )
525        primal_bounds_lon, primal_bounds_lat = get_pvars(dim_primal, 'primal_bounds_lon', 'primal_bounds_lat')
526        edge_deg, trisk_deg, left_right, down_up, trisk, le, de, wee, lon_e, lat_e, angle_e = get_pvars(
527             dim_edge, 'edge_deg', 'trisk_deg', 'edge_left_right', 'edge_down_up', 
528             'trisk', 'le', 'de', 'wee', 'lon_e', 'lat_e', 'angle_e')
529        dual_deg, dual_vertex, dual_edge, dual_ne, Riv2, lon_v, lat_v, Av = get_pvars(
530            dim_dual, 'dual_deg', 'dual_vertex', 'dual_edge', 'dual_ne', 'Riv2', 'lon_v', 'lat_v', 'Av')
531        dual_bounds_lon, dual_bounds_lat = get_pvars(dim_dual, 'dual_bounds_lon', 'dual_bounds_lat')
532        if gridfile.meshtype == 'curvilinear':
533            self.primal_i, self.primal_j = get_pvars(dim_primal, 'primal_i','primal_j')
534            self.dual_i, self.dual_j = get_pvars(dim_dual, 'primal_i','primal_j')
535            self.edge_i, self.edge_j = get_pvars(dim_edge, 'edge_i','edge_j')
536
537        # Let indices start at 0
538        for x in (primal_vertex, primal_edge, dual_vertex,dual_edge,
539            left_right, down_up, trisk ) : x.data = x.data-gridfile.start_index
540        self.edge2cell   = Stencil_glob(edge_deg, left_right)
541        self.cell2edge   = Stencil_glob(primal_deg, primal_edge, primal_ne)
542        self.cell2vertex = Stencil_glob(primal_deg, primal_vertex)
543        self.edge2vertex = Stencil_glob(edge_deg, down_up)
544        self.vertex2edge = Stencil_glob(dual_deg, dual_edge, dual_ne)
545        self.vertex2cell = Stencil_glob(dual_deg, dual_vertex, Riv2)
546        self.edge2edge   = Stencil_glob(trisk_deg, trisk, wee)
547        self.set_members(locals(), 'dim_primal', 'dim_dual', 'dim_edge',
548                         'primal_deg', 'primal_vertex', 'primal_edge',
549                         'primal_bounds_lon', 'primal_bounds_lat',
550                         'trisk_deg', 'trisk', 'dual_deg', 'dual_edge',
551                         'dual_bounds_lon', 'dual_bounds_lat',
552                         'lon_v', 'lat_v', 'Av', 'lon_i', 'lat_i', 'Ai', 
553                         'le', 'de', 'lon_e', 'lat_e', 'angle_e')
554    def partition_metis(self):
555        INFO('Partitioning with ParMetis...')
556# TODO        edge_owner = partition_mesh(self.comm, self.trisk_deg, self.trisk)
557        edge_owner = partition_mesh(self.trisk_deg, self.trisk, self.comm.Get_size())
558        self.edge_owner = parallel.LocPArray1D(self.dim_edge, edge_owner)
559        primal_owner = partition_from_stencil(self.edge_owner, self.primal_deg, self.primal_edge)
560        self.primal_owner = parallel.LocPArray1D(self.dim_primal, primal_owner)
561        dual_owner = partition_from_stencil(self.edge_owner, self.dual_deg, self.dual_edge)
562        self.dual_owner = parallel.LocPArray1D(self.dim_dual, dual_owner)
563
564    def partition_curvilinear(self, ni,nj):
565        def owner(dim,i,j):
566            owner_i, owner_j =  (ni*i.data)/nx, (nj*j.data)/ny
567            return parallel.LocPArray1D(dim, owner_i + ni*owner_j)
568        nx, ny = self.gridfile.nx, self.gridfile.ny
569        INFO('Partitioning %d x %d cells in %d x %d blocks ...'%(nx,ny,ni,nj))
570        n = self.comm.Get_size()
571        assert n == ni*nj, 'Mismatch between ni x nj = %d and MPI_SIZE = %d.'%(ni*nj, n)
572        self.edge_owner   = owner(self.dim_edge,   self.edge_i,   self.edge_j)
573        self.primal_owner = owner(self.dim_primal, self.primal_i, self.primal_j)
574        self.dual_owner = self.primal_owner
575        DEBUG('partition_curvilinear %d : %s'%(self.comm.Get_rank(), self.primal_owner.data) )
576
577    def get(self, getter, *names): return [ getter(self.__dict__[x]) for x in names ]
578
579def chain(start, links):
580    for link in links:
581        start = link(start).neigh_set
582        yield start
583
584def patches(degree, bounds, lon, lat):
585    for i in range(degree.size):
586        nb_edge=degree[i]
587        bounds_cell = bounds[i,0:nb_edge]
588        lat_cell    = lat[bounds_cell]
589        lon_cell    = lon[bounds_cell]
590        orig=lon_cell[0]
591        lon_cell    = lon_cell-orig+180.
592        lon_cell    = np.mod(lon_cell,360.)
593        lon_cell    = lon_cell+orig-180.
594#        if np.abs(lon_cell-orig).max()>10. :
595#            print '%d patches :'%mpi_rank, lon_cell
596        lonlat_cell = np.zeros((nb_edge,2))
597        lonlat_cell[:,0],lonlat_cell[:,1] = lon_cell,lat_cell
598        polygon = Polygon(lonlat_cell, True)
599        yield polygon
600
601class Local_Mesh(Abstract_Mesh):
602    Halo_Xchange = parallel.Halo_Xchange # overriden in dev.meshes
603    def __init__(self, pmesh, llm, nqdyn, mapping):
604        comm, dim_primal, primal_owner, dim_edge, edge_owner, dim_dual, dual_owner = pmesh.members(
605            'comm', 'dim_primal', 'primal_owner', 'dim_edge', 'edge_owner', 'dim_dual', 'dual_owner' )
606        INFO('Constructing halos ...')
607        edges_E0 = find_my_cells(edge_owner)
608        cells_C0, edges_E1, vertices_V1, edges_E2, cells_C1 = chain(
609            edges_E0, ( pmesh.edge2cell, pmesh.cell2edge, pmesh.edge2vertex, pmesh.vertex2edge, pmesh.edge2cell) )
610        edges_E0, edges_E1, edges_E2 = progressive_list(edges_E0, edges_E1, edges_E2)
611        cells_C0, cells_C1 = progressive_list(cells_C0, cells_C1)
612        DEBUG('E2,E1,E0 ; C1,C0 : %s' % map(len, (edges_E2, edges_E1, edges_E0, cells_C1, cells_C0)) )
613        get_own_cells, get_all_cells = dim_primal.getter(cells_C0), dim_primal.getter(cells_C1)
614        get_all_duals, get_all_edges = dim_dual.getter(vertices_V1), dim_edge.getter(edges_E2)
615
616        primal_num, dual_num, edge_num = map(len, (cells_C1, vertices_V1, edges_E2))
617
618        self.meshtype = pmesh.gridfile.meshtype
619        if self.meshtype == 'curvilinear' :
620            self.nx, self.ny = pmesh.gridfile.nx, pmesh.gridfile.ny
621            self.primal_i, self.primal_j = pmesh.get(get_all_cells, 'primal_i', 'primal_j')
622            self.dual_i, self.dual_j = pmesh.get(get_all_duals, 'dual_i', 'dual_j')
623
624        # construct local stencils from global stencils
625        dict_E2, dict_C1, dict_V1 = map(inverse_list, (edges_E2, cells_C1, vertices_V1))
626        down_up       = pmesh.edge2vertex( edges_E2,    dict_V1)
627        left_right    = pmesh.edge2cell(   edges_E2,    dict_C1)
628        primal_edge   = pmesh.cell2edge(   cells_C1,    dict_E2)
629        dual_edge     = pmesh.vertex2edge( vertices_V1, dict_E2)
630        trisk         = pmesh.edge2edge(   edges_E2,    dict_E2)
631        Riv2          = pmesh.vertex2cell( vertices_V1, dict_C1)
632        DEBUG('Edge-vertex degree min/max %d %d'%(down_up.degree.min(), down_up.degree.max()) )
633        DEBUG('Edge-cell degree min/max %d %d'%(left_right.degree.min(), left_right.degree.max()) )
634        DEBUG('Primal degree min/max : %d %d'%(primal_edge.degree.min(), primal_edge.degree.max()) )
635        DEBUG('Dual degree min/max : %d %d'%(dual_edge.degree.min(), dual_edge.degree.max()) )
636        DEBUG('TRiSK degree min/max : %d %d'%(trisk.degree.min(), trisk.degree.max()) )
637        DEBUG('Riv2 degree min/max : %d %d'%(Riv2.degree.min(), Riv2.degree.max()) )
638        primal_deg, primal_edge, primal_ne = primal_edge.degree,   primal_edge.neigh_loc, primal_edge.weight
639        dual_deg,   dual_edge, dual_ne  = dual_edge.degree, dual_edge.neigh_loc, dual_edge.weight
640        trisk_deg,  trisk, wee      = trisk.degree,         trisk.neigh_loc, trisk.weight
641        dual_deg, dual_vertex, Riv2 = Riv2.degree,          Riv2.neigh_loc,  Riv2.weight
642        edge_deg, left, right = left_right.degree, left_right.neigh_loc[:,0], left_right.neigh_loc[:,1]
643        edge_deg, down, up    = down_up.degree,    down_up.neigh_loc[:,0],    down_up.neigh_loc[:,1]
644        for edge in range(edge_num): 
645            if edge_deg[edge]<2: up[edge]=down[edge]
646
647        # construct own primal mesh for XIOS output
648        primal_own_glo = find_my_cells(primal_owner)
649        primal_vertex = pmesh.cell2vertex(primal_own_glo, dict_V1)
650        primal_own_deg, primal_vertex = primal_vertex.degree, primal_vertex.neigh_loc
651        primal_own_glo = np.asarray(primal_own_glo, dtype=np.int32)
652        primal_own_loc = [dict_C1[i] for i in primal_own_glo]
653        primal_own_num = len(primal_own_glo)
654        primal_own_all = comm.allgather(primal_own_num)
655        displ = sum(primal_own_all[0:comm.Get_rank()]) # required by XIOS
656#        print 'local primal mesh :', primal_vertex.min(), primal_vertex.max(), primal_own_all, displ
657
658        # construct own dual mesh for XIOS output
659        dual_own_glo = find_my_cells(dual_owner)
660        dual_own_glo = np.asarray(dual_own_glo, dtype=np.int32)
661        dual_own_loc = [dict_V1[i] for i in dual_own_glo]
662       
663        # compute_ne(...) and normalize(...) expect indices starting at 1
664        for ind in (trisk, primal_vertex, dual_vertex, primal_edge, dual_edge) : ind+=1
665        left, right, down, up = left+1, right+1, down+1, up+1 # these ones are subarrays, we need deep copies
666
667        # read from mesh file : coordinates, lengths and areas in reference domain
668        lon_i, lat_i, Ai = pmesh.get(get_all_cells, 'lon_i', 'lat_i', 'Ai')
669
670
671        lon_v, lat_v, Av = pmesh.get(get_all_duals, 'lon_v', 'lat_v', 'Av')
672        le, de, lon_e, lat_e, angle_e = pmesh.get(get_all_edges, 'le', 'de', 'lon_e', 'lat_e', 'angle_e')
673        le_de = le/de
674
675        # store what we have read/computed
676        self.set_members(locals(), 'nqdyn', 'llm', 'primal_num', 'dual_num', 'edge_num',
677                         'primal_own_glo', 'primal_own_loc', 'primal_own_deg',
678                         'primal_deg', 'primal_vertex', 'displ', 'dual_own_loc', 
679                         'trisk_deg', 'trisk', 'wee',
680                         'dual_deg', 'dual_edge', 'dual_ne', 'dual_vertex', 'Riv2',
681                         'left', 'right', 'down', 'up', 
682                         'primal_deg', 'primal_edge', 'primal_ne',
683                         'vertices_V1', 'edges_E2',
684                         'lon_i', 'lat_i', 'Ai',
685                         'lon_v', 'lat_v', 'Av',
686                         'le', 'de', 'le_de', 'lon_e', 'lat_e', 'angle_e')
687
688        # obtain bounds of primal and dual cells (for XIOS)
689        if pmesh.primal_bounds_lon is not None:
690            INFO('Reading cell bounds from mesh file.')
691            self.primal_bounds_lon, self.primal_bounds_lat = pmesh.get(get_all_cells, 'primal_bounds_lon', 'primal_bounds_lat')
692            self.dual_bounds_lon, self.dual_bounds_lat = pmesh.get(get_all_duals, 'dual_bounds_lon', 'dual_bounds_lat')
693        else:
694            INFO('Mesh file does not contain information abound cell bounds.')
695            pmesh.gridfile.get_bounds(self)
696
697        # normalize data read from file depending on file format
698        pmesh.gridfile.normalize(self)
699        # map from reference domain to physical domain
700        self.apply_map(mapping)
701        # now latitudes and longitudes refer to the physical domain
702
703        self.fv = mapping.coriolis(self.lon_v,self.lat_v) # Coriolis parameter
704        # propagate info to Fortran side
705        self.to_dynamico()
706        # setup halo transfers - NB : llm must be set before we call set_dynamico_transfer
707        self.com_primal = self.Halo_Xchange(
708            42, dim_primal, cells_C1, get_all_cells(primal_owner))
709        self.com_edges = self.Halo_Xchange(
710            73, dim_edge, edges_E2, get_all_edges(edge_owner))
711        self.com_primal.set_dynamico_transfer('primal')
712        self.com_edges.set_dynamico_transfer('edge')
713
714        self.primal  = tri.Triangulation(self.lon_i*radian, self.lat_i*radian)
715        self.dual    = tri.Triangulation(self.lon_v*radian, self.lat_v*radian)
716        self.triedge = tri.Triangulation(self.lon_e*radian, self.lat_e*radian)       
717
718    def plot_patches(self, ax, clim, degree, bounds, lon, lat, data):
719        nb_vertex = lon.size # global
720        p = list(patches(degree, bounds, lon, lat))
721        p = PatchCollection(p, linewidth=0.01)
722        p.set_array(data) # set values at each polygon (cell)
723        p.set_clim(clim)
724        ax.add_collection(p)
725
726#--------------------------------------  Mesh reordering  ------------------------------------------#
727
728# NB : Indices start at 0
729# permutations map an old index to a new index : perm[old_index]=new_index
730# enumerate(perm) is the list [ (0,perm[0]) , (1,perm[1]), ...] = [(old, new) ... ]
731
732def reorder_values(perm, *datas): 
733    # data in datas : ndarray containing indices ; perm : permutation to apply to values
734    for data in datas:
735        for x in np.nditer(data, op_flags=['readwrite']):
736            x[...]=perm[x]
737
738def reorder_indices(perm, *datas): 
739    # data in datas : ndarray containing values ; perm : permutation to apply to indices
740    for data in datas:
741        data_out, ndim = data.copy(), data.ndim
742        if ndim == 1:
743            for old,new in enumerate(perm):
744                data_out[new]=data[old]
745        if ndim == 2:
746            for old,new in enumerate(perm):
747                data_out[new,:]=data[old,:]
748        data[:]=data_out[:]
749
750def key_to_perm(keys): 
751    # keys : maps an old index to a key
752    # returns : a list mapping old indices to new indices obtained by sorting keys in ascending order
753    ordered = [old for key,old in sorted(zip(keys,range(len(keys))))] 
754    # ordered maps new indices to old indices
755    # but we want to map old indices to new indices
756    perm = list(range(len(keys)))
757    for new,old in enumerate(ordered):
758        perm[old]=new
759    return perm  # maps old indices to new indices
760
761#-------------------------------------- Mesh partitioning ------------------------------------------#
762
763# NB : Indices start at 0
764
765def partition_from_stencil(owner2, degree, stencil):
766    # given a stencil dim1->dim2 and owner2 on dim2, define owner[i] on dim1 as min(stencil[i,:] if i is even, max if odd
767    dim1, dim2= degree.dim, owner2.dim
768    degree, stencil, n = degree.data, stencil.data, dim1.end-dim1.start
769    cells2 = list_stencil(degree, stencil)
770    cells2 = sorted(list(set(list(cells2)))) # list of cells for which we need to know owner2
771    get2 = dim2.getter(cells2)
772    owner1, owner2, glob2loc = np.zeros(n, dtype=np.int32), get2(owner2), get2.dict
773    for i in range(n):
774        owners = [ owner2[glob2loc[stencil[i,j]]] for j in range(degree[i]) ]
775        if i%2 == 0 : owner1[i] = min(owners)
776        else : owner1[i] = max(owners)
777    return owner1
778           
779def find_my_cells(owner): # owner : a PArray1D containing the data returned by partition_mesh
780    dim, comm, owner = owner.dim, owner.dim.comm, owner.data
781    mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size()
782    cells=[set() for i in range(mpi_size)]
783    for i in range(len(owner)):
784        cells[owner[i]].add(dim.start+i)
785    cells = [sorted(list(cells[i])) for i in range(mpi_size)]
786    mycells = comm.alltoall(cells)
787    mycells = sorted(sum(mycells, [])) # concatenate into a single list
788    return mycells
789
790#---------------------------------- Stencil management ---------------------------------------#
791
792# NB : Indices start at 0
793
794# Class Stencil represents an adjacency relationship (e.g. cell=>edges)
795# using adjacency information read from PArrays
796# It computes a list of "edges" adjacent to a given list of "cells"
797# This is used to form the sets E0 -> C0 -> E1 -> V1 -> E2 -> C1
798# which are then used to form lists of global indices for V1,C1,E2 such that
799# C0, E0, E1 form contiguous subsets of C1, E2 starting from 0
800
801def reindex(vertex_dict, degree, bounds):
802    for i in range(degree.size):
803        for j in range(degree[i]):
804            bounds[i,j] = vertex_dict[bounds[i,j]]
805    return bounds
806
807class Stencil_glob(object):
808    def __init__(self, degree, neigh, weight=None):
809        self.degree, self.neigh, self.weight = degree, neigh, weight
810    def __call__(self, cells, neigh_dict=None):
811        return Stencil(cells, self.degree, self.neigh, neigh_dict, self.weight)
812           
813class Stencil(object):
814    def __init__(self, cells, degree, neigh, neigh_dict, weight=None):
815        get = degree.dim.getter(cells)
816        mydegree, myneigh = [get(x) for x in (degree, neigh)]
817        if not weight is None : myweight = get(weight)
818        if neigh_dict is None : 
819            keep = lambda n : True
820        else : # if neigh_dict is present, only neighbours in dict are retained
821            keep = lambda n : n in neigh_dict
822        neigh_set = list_stencil(mydegree, myneigh, keep)
823        self.neigh_set = list(set(list(neigh_set)     ))                         
824        rej=0
825        for i in range(len(mydegree)): # keep only elements in neigh_dict, in-place
826            k=0
827            for j in range(mydegree[i]):
828                n=myneigh[i,j]
829                if keep(n):
830                    myneigh[i,k]=n
831                    if not weight is None : myweight[i,k]=myweight[i,j]
832                    k=k+1
833                else:
834                    rej=rej+1
835            mydegree[i]=k
836        if neigh_dict is None:
837            neigh_dict = {j:i for i,j in enumerate(self.neigh_set)}
838        myneigh_loc = reindex(neigh_dict, mydegree, myneigh)
839        self.degree, self.neigh_glob, self.neigh_loc = mydegree, myneigh, myneigh_loc
840        if weight is not None: self.weight = myweight
841           
842def progressive_iter(mylist, cell_lists):
843    for thelist in cell_lists: 
844        mylist = mylist + list(set(thelist)-set(mylist))
845        yield mylist
846def progressive_list(*cell_lists):
847    # cell_lists : a tuple of lists of global indices, with each list a subset of the next
848    # returns : a list 'mylist' such that for each list 'thelist' in cell_lists, thelist = mylist[0:len(thelist)]
849    # example : edge_list = progressive_list(E0,E1,E2) with E0,E1,E2 increasing lists of cell edges
850    return list(progressive_iter([], cell_lists))
Note: See TracBrowser for help on using the repository browser.