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

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

devel/unstructured : implement reference vs physical mesh for spherical meshes

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