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

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

devel/Python : temporary fix to meshes.py

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