1 | import time |
---|
2 | import math |
---|
3 | import numpy as np |
---|
4 | import netCDF4 as cdf |
---|
5 | import matplotlib.tri as tri |
---|
6 | import matplotlib.pyplot as plt |
---|
7 | from matplotlib.patches import Polygon |
---|
8 | from matplotlib.collections import PatchCollection |
---|
9 | |
---|
10 | import unstructured as unst |
---|
11 | import parallel |
---|
12 | from util import list_stencil, Base_class, inverse_list |
---|
13 | from precision import zeros |
---|
14 | |
---|
15 | radian=180/math.pi |
---|
16 | |
---|
17 | #------------------- Hybrid mass-based coordinate ------------- |
---|
18 | |
---|
19 | # compute hybrid coefs from distribution of mass |
---|
20 | |
---|
21 | def 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) |
---|
48 | def 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 |
---|
75 | def 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 | |
---|
83 | class Cartesian_mesh: |
---|
84 | def __init__(self,nx,ny,llm,nqdyn,Lx,Ly,f): |
---|
85 | dx,dy = Lx/nx, Ly/ny |
---|
86 | self.dx, self.dy, self.f = dx,dy,f |
---|
87 | self.nx, self.ny, self.llm, self.nqdyn = nx,ny,llm,nqdyn |
---|
88 | self.field_z = self.field_mass |
---|
89 | # 1D coordinate arrays |
---|
90 | x=(np.arange(nx)-nx/2.)*Lx/nx |
---|
91 | y=(np.arange(ny)-ny/2.)*Ly/ny |
---|
92 | lev=np.arange(llm) |
---|
93 | levp1=np.arange(llm+1) |
---|
94 | self.x, self.y, self.lev, self.levp1 = x,y,lev,levp1 |
---|
95 | # 3D coordinate arrays |
---|
96 | self.xx,self.yy,self.ll = np.meshgrid(x,y,lev, indexing='ij') |
---|
97 | self.xxp1,self.yyp1,self.llp1 = np.meshgrid(x,y,levp1, indexing='ij') |
---|
98 | # beware conventions for indexing |
---|
99 | # Fortran order : llm,nx*ny,nqdyn / indices start at 1 |
---|
100 | # Python order : nqdyn,ny,nx,llm / indices start at 0 |
---|
101 | # indices below follow Fortran while x,y follow Python/C |
---|
102 | index=lambda x,y : ((x+(nx*(y+2*ny)))%(nx*ny))+1 |
---|
103 | indexu=lambda x,y : 2*index(x,y)-1 |
---|
104 | indexv=lambda x,y : 2*index(x,y) |
---|
105 | indices = lambda shape : np.zeros(shape,np.int32) |
---|
106 | |
---|
107 | primal_nb = indices(nx*ny) |
---|
108 | primal_edge = indices((nx*ny,4)) |
---|
109 | primal_ne = indices((nx*ny,4)) |
---|
110 | |
---|
111 | dual_nb = indices(nx*ny) |
---|
112 | dual_edge = indices((nx*ny,4)) |
---|
113 | dual_ne = indices((nx*ny,4)) |
---|
114 | dual_vertex = indices((nx*ny,4)) |
---|
115 | Riv2 = np.zeros((nx*ny,4)) |
---|
116 | |
---|
117 | left = indices(2*nx*ny) |
---|
118 | right = indices(2*nx*ny) |
---|
119 | up = indices(2*nx*ny) |
---|
120 | down = indices(2*nx*ny) |
---|
121 | le_de = np.zeros(2*nx*ny) |
---|
122 | |
---|
123 | trisk_deg = indices(2*nx*ny) |
---|
124 | trisk = indices((2*nx*ny,4)) # 4 TRiSK coefs per edge |
---|
125 | wee = np.zeros((2*nx*ny,4)) |
---|
126 | |
---|
127 | for x in range(nx): |
---|
128 | for y in range(ny): |
---|
129 | # NB : Fortran indices start at 1 |
---|
130 | # primal cells |
---|
131 | put(index(x,y)-1,primal_nb,(primal_edge,primal_ne), |
---|
132 | ((indexu(x,y),1), |
---|
133 | (indexv(x,y),1), |
---|
134 | (indexu(x-1,y),-1), |
---|
135 | (indexv(x,y-1),-1) )) |
---|
136 | # dual cells |
---|
137 | put(index(x,y)-1,dual_nb,(dual_edge,dual_vertex,dual_ne,Riv2), |
---|
138 | ((indexv(x+1,y),index(x,y),1,.25), |
---|
139 | (indexu(x,y+1),index(x+1,y),-1,.25), |
---|
140 | (indexv(x,y),index(x+1,y+1),-1,.25), |
---|
141 | (indexu(x,y),index(x,y+1),1,.25) )) |
---|
142 | # edges : |
---|
143 | # left and right are adjacent primal cells |
---|
144 | # flux is positive when going from left to right |
---|
145 | # up and down are adjacent dual cells (vertices) |
---|
146 | # circulation is positive when going from down to up |
---|
147 | # u-points |
---|
148 | ij =indexu(x,y)-1 |
---|
149 | left[ij]=index(x,y) |
---|
150 | right[ij]=index(x+1,y) |
---|
151 | down[ij]=index(x,y-1) |
---|
152 | up[ij]=index(x,y) |
---|
153 | le_de[ij]=dy/dx |
---|
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 | left[ij]=index(x,y) |
---|
162 | right[ij]=index(x,y+1) |
---|
163 | down[ij]=index(x,y) |
---|
164 | up[ij]=index(x-1,y) |
---|
165 | le_de[ij]=dx/dy |
---|
166 | put(ij,trisk_deg,(trisk,wee),( |
---|
167 | (indexu(x,y),.25), |
---|
168 | (indexu(x-1,y),.25), |
---|
169 | (indexu(x,y+1),.25), |
---|
170 | (indexu(x-1,y+1),.25))) |
---|
171 | Aiv=np.zeros(nx*ny) |
---|
172 | Aiv[:]=dx*dy # Ai=Av=dx*dy |
---|
173 | unst.init_mesh(llm,nqdyn,2*nx*ny,nx*ny,nx*ny,4,4,4, |
---|
174 | primal_nb,primal_edge,primal_ne, |
---|
175 | dual_nb,dual_edge,dual_ne,dual_vertex, |
---|
176 | left,right,down,up,trisk_deg,trisk, |
---|
177 | Aiv,Aiv,f+0.*Aiv,le_de,Riv2,-wee) |
---|
178 | |
---|
179 | def field_theta(self,n=1): return zeros((n,self.nqdyn,self.ny,self.nx,self.llm)) |
---|
180 | def field_mass(self,n=1): return zeros((n,self.ny,self.nx,self.llm)) |
---|
181 | def field_z(self,n=1): return zeros((n,self.ny,self.nx,self.llm)) |
---|
182 | def field_w(self,n=1): return zeros((n,self.ny,self.nx,self.llm+1)) |
---|
183 | def field_u(self,n=1): return zeros((n,self.ny,2*self.nx,self.llm)) |
---|
184 | def field_ps(self,n=1): return zeros((n,self.ny,self.nx)) |
---|
185 | def ucomp(self,u): |
---|
186 | return u[range(0,2*self.nx,2),:] if self.ny==1 else u[:,range(0,2*self.nx,2),:] |
---|
187 | def set_ucomp(self,uv,u): |
---|
188 | if self.ny==1 : uv[range(0,2*self.nx,2),:]=u |
---|
189 | else : uv[:,range(0,2*self.nx,2),:]=u |
---|
190 | def vcomp(self,u): return u[:,range(1,2*self.nx,2),:] |
---|
191 | |
---|
192 | #---------------------- MPAS fully unstructured mesh ------------------------ |
---|
193 | |
---|
194 | class Mesh_Format(Base_class): |
---|
195 | def getdims(self,*names): return [len(self.nc.dimensions[name]) for name in names] |
---|
196 | def get_pdims(self,comm,*names): return [self.get_pdim(comm,name) for name in names] |
---|
197 | def get_pvars(self,pdim,*names): return [self.get_pvar(pdim,name) for name in names] |
---|
198 | def getvars(self,*names): |
---|
199 | for name in names : print "getvar %s ..."%name |
---|
200 | time1=time.time() |
---|
201 | ret=[self.getvar(name) for name in names] |
---|
202 | print "... Done in %f seconds"%(time.time()-time1) |
---|
203 | return ret |
---|
204 | |
---|
205 | class MPAS_Format(Mesh_Format): |
---|
206 | start_index=1 # indexing starts at 1 as in Fortran |
---|
207 | translate= { |
---|
208 | 'primal_num':'nCells', |
---|
209 | 'edge_num':'nEdges', |
---|
210 | 'dual_num':'nVertices', |
---|
211 | 'primal_deg':'nEdgesOnCell', 'trisk_deg':'nEdgesOnEdge', |
---|
212 | 'primal_edge':'edgesOnCell', 'left_right':'cellsOnEdge', |
---|
213 | 'dual_edge':'edgesOnVertex', 'down_up':'verticesOnEdge', |
---|
214 | 'primal_vertex':'verticesOnCell', 'dual_vertex':'cellsOnVertex', |
---|
215 | 'trisk':'edgesOnEdge', 'down_up':'verticesOnEdge', 'left_right':'cellsOnEdge', |
---|
216 | 'le':'dvEdge', 'de':'dcEdge', 'Ai':'areaCell', 'Av':'areaTriangle', |
---|
217 | 'lat_i':'latCell','lon_i':'lonCell','lat_v':'latVertex','lon_v':'lonVertex', |
---|
218 | 'lat_e':'latEdge','lon_e':'lonEdge','angle_e':'angleEdge', |
---|
219 | 'wee':'weightsOnEdge','Riv2':'kiteAreasOnVertex'} |
---|
220 | def __init__(self,gridfilename): |
---|
221 | try: |
---|
222 | self.nc = cdf.Dataset(gridfilename, "r") |
---|
223 | except RuntimeError: |
---|
224 | print """ |
---|
225 | Unable to open grid file %s, maybe you forgot to download it ? |
---|
226 | To do so, go to the 'Python/' dir and execute './get_MPAS_grids.sh'. |
---|
227 | """ % gridfile |
---|
228 | raise |
---|
229 | def get_pdim(self,comm,name): |
---|
230 | name = self.translate[name] |
---|
231 | return parallel.PDim(self.nc.dimensions[name], comm) |
---|
232 | def getvar(self,name): |
---|
233 | # first deal with special cases |
---|
234 | if name == 'dual_deg': |
---|
235 | dual_deg = [3 for i in range(self.dual_num) ] |
---|
236 | dual_deg = np.ascontiguousarray(dual_deg,dtype=np.int32) |
---|
237 | # NB : Fortran code expects 32-bit ints |
---|
238 | return dual_deg |
---|
239 | # general case |
---|
240 | name=self.translate[name] |
---|
241 | return self.nc.variables[name][:] |
---|
242 | def get_pvar(self,pdim,name): |
---|
243 | # provides parallel access to a NetCDF variable, with name translation |
---|
244 | # first deal with special cases |
---|
245 | if name == 'dual_deg': return parallel.CstPArray1D(pdim, np.int32, 3) |
---|
246 | if name == 'edge_deg': return parallel.CstPArray1D(pdim, np.int32, 2) |
---|
247 | # general case |
---|
248 | name=self.translate[name] |
---|
249 | return parallel.PArray(pdim, self.nc.variables[name]) |
---|
250 | def normalize(self, mesh, radius): |
---|
251 | (trisk_deg, trisk, Ai, |
---|
252 | de, le, wee, Av, Riv2) = (mesh.trisk_deg, mesh.trisk, mesh.Ai, |
---|
253 | mesh.de, mesh.le, mesh.wee, mesh.Av, mesh.Riv2) |
---|
254 | edge_num, dual_num, r2 = de.size, Av.size, radius**2 |
---|
255 | # fix normalization of wee and Riv2 weights |
---|
256 | for edge1 in range(edge_num): |
---|
257 | for i in range(trisk_deg[edge1]): |
---|
258 | edge2=trisk[edge1,i]-1 # NB Fortran vs C/Python indexing |
---|
259 | wee[edge1,i] = de[edge1]*wee[edge1,i]/le[edge2] |
---|
260 | for ivertex in range(dual_num): |
---|
261 | for j in range(3): |
---|
262 | Riv2[ivertex,j]=Riv2[ivertex,j]/Av[ivertex] |
---|
263 | r=Riv2[ivertex,:] |
---|
264 | r=sum(r) |
---|
265 | if abs(r-1.)>1e-6 : print 'error with Riv2 at vertex ', ivertex |
---|
266 | # scale lengths and areas |
---|
267 | mesh.le, mesh.de, mesh.Av, mesh.Ai = radius*le, radius*de, r2*Av, r2*Ai |
---|
268 | |
---|
269 | def compute_ne(num,deg,edges,left,right): |
---|
270 | ne = np.zeros_like(edges) |
---|
271 | for cell in range(num): |
---|
272 | for iedge in range(deg[cell]): |
---|
273 | edge = edges[cell,iedge]-1 |
---|
274 | if left[edge]==cell+1: ne[cell,iedge]+=1 |
---|
275 | if right[edge]==cell+1: ne[cell,iedge]-=1 |
---|
276 | # if ne[cell,iedge]==0 : print 'compute_ne error at cell,iedge', cell, iedge |
---|
277 | return ne |
---|
278 | |
---|
279 | class Abstract_Mesh(Base_class): |
---|
280 | def field_theta(self,n=1): return zeros((n,self.nqdyn,self.primal_num,self.llm)) |
---|
281 | def field_mass(self,n=1): return zeros((n,self.primal_num,self.llm)) |
---|
282 | def field_z(self,n=1): return zeros((n,self.dual_num,self.llm)) |
---|
283 | def field_w(self,n=1): return zeros((n,self.primal_num,self.llm+1)) |
---|
284 | def field_u(self,n=1): return zeros((n,self.edge_num,self.llm)) |
---|
285 | def field_ps(self,n=1): return zeros((n,self.primal_num,)) |
---|
286 | def ucov2D(self, ulon, ulat): |
---|
287 | return self.de*(ulon*np.cos(self.angle_e)+ulat*np.sin(self.angle_e)) |
---|
288 | def ucov3D(self, ulon, ulat): |
---|
289 | ucov = zeros((self.edge_num,ulon.shape[1])) |
---|
290 | for edge in range(self.edge_num): |
---|
291 | angle=self.angle_e[edge] |
---|
292 | ucov[edge,:] = self.de[edge]*(ulon[edge,:]*math.cos(angle)+ulat[edge,:]*math.sin(angle)) |
---|
293 | return ucov |
---|
294 | def plot(self, tri,data, **kwargs): |
---|
295 | plt.figure(figsize=(10,4)) |
---|
296 | plt.gca().set_aspect('equal') |
---|
297 | plt.tricontourf(tri, data, 20, **kwargs) |
---|
298 | plt.colorbar() |
---|
299 | plt.ylim((-90,90)) |
---|
300 | plt.xlim((-180,180)) |
---|
301 | def plot_i(self,data, **kwargs): |
---|
302 | self.plot(self.primal,data, **kwargs) |
---|
303 | def plot_v(self,data, **kwargs): |
---|
304 | self.plot(self.dual,data, **kwargs) |
---|
305 | def plot_e(self,data, **kwargs): |
---|
306 | self.plot(self.triedge,data, **kwargs) |
---|
307 | |
---|
308 | class Unstructured_Mesh(Abstract_Mesh): |
---|
309 | def __init__(self, gridfile, llm, nqdyn, radius, f): |
---|
310 | getvars = gridfile.getvars |
---|
311 | # get positions, lengths, surfaces and weights |
---|
312 | le,de,Ai,Av,wee,Riv2 = getvars('le','de','Ai','Av','wee','Riv2') |
---|
313 | lat_i,lon_i,lat_v,lon_v = getvars('lat_i','lon_i','lat_v','lon_v') |
---|
314 | lat_e,lon_e,angle_e = getvars('lat_e','lon_e','angle_e') |
---|
315 | fv = f(lon_v,lat_v) # Coriolis parameter |
---|
316 | primal_num, dual_num, edge_num = [ x.size for x in Ai, Av, le ] |
---|
317 | gridfile.primal_num, gridfile.dual_num = primal_num, dual_num |
---|
318 | primal_deg, dual_deg, trisk_deg = getvars('primal_deg','dual_deg','trisk_deg') |
---|
319 | # get indices for stencils |
---|
320 | # primal <-> edges |
---|
321 | primal_edge, left_right = getvars('primal_edge','left_right') |
---|
322 | left,right = left_right[:,0], left_right[:,1] |
---|
323 | primal_ne = compute_ne(primal_num,primal_deg,primal_edge,left,right) |
---|
324 | # dual <-> edges |
---|
325 | dual_edge, down_up = getvars('dual_edge','down_up') |
---|
326 | down,up = down_up[:,0], down_up[:,1] |
---|
327 | dual_ne = compute_ne(dual_num,dual_deg,dual_edge,up,down) |
---|
328 | # primal <-> dual, edges <-> edges (TRiSK) |
---|
329 | primal_vertex, dual_vertex, trisk = getvars('primal_vertex','dual_vertex','trisk') |
---|
330 | # CRITICAL : force arrays left, etc. to be contiguous in memory: |
---|
331 | left,right,up,down = [np.ascontiguousarray(x) for x in (left,right,up,down)] |
---|
332 | trisk,wee,primal_edge = [np.ascontiguousarray(x) for x in (trisk,wee,primal_edge)] |
---|
333 | self.set_members(locals(), 'gridfile','llm','nqdyn','radius','f', |
---|
334 | 'primal_num', 'dual_num', 'edge_num', 'primal_deg', 'primal_vertex', |
---|
335 | 'lat_i', 'lon_i', 'lat_v', 'lon_v', 'lat_e', 'lon_e', 'angle_e', |
---|
336 | 'Ai', 'Av', 'fv', 'le', 'de', |
---|
337 | 'trisk_deg', 'trisk', 'wee', 'Riv2') |
---|
338 | gridfile.normalize(self, radius) |
---|
339 | self.le_de = self.le/self.de |
---|
340 | max_primal_deg, max_dual_deg, max_trisk_deg = [ x.shape[1] for x in primal_edge, dual_edge, trisk] |
---|
341 | unst.init_mesh(llm, nqdyn, edge_num, primal_num,dual_num, |
---|
342 | max_trisk_deg, max_primal_deg, max_dual_deg, |
---|
343 | primal_deg,primal_edge,primal_ne, |
---|
344 | dual_deg,dual_edge,dual_ne,dual_vertex, |
---|
345 | left,right,down,up,trisk_deg,trisk, |
---|
346 | self.Ai,self.Av,self.fv,le/de,Riv2,wee) |
---|
347 | self.primal = tri.Triangulation(lon_i*radian, lat_i*radian) |
---|
348 | self.dual = tri.Triangulation(lon_v*radian, lat_v*radian) |
---|
349 | self.triedge = tri.Triangulation(lon_e*radian, lat_e*radian) |
---|
350 | |
---|
351 | class Unstructured_PMesh(Abstract_Mesh): # Mesh data distributed across MPI processes |
---|
352 | def __init__(self,comm, gridfile): |
---|
353 | self.comm, self.gridfile, get_pvars = comm, gridfile, gridfile.get_pvars |
---|
354 | dim_primal, dim_edge, dim_dual = gridfile.get_pdims(comm, |
---|
355 | 'primal_num','edge_num','dual_num') |
---|
356 | primal_deg, primal_vertex, primal_edge, lon_i, lat_i, Ai = get_pvars( |
---|
357 | dim_primal, 'primal_deg', 'primal_vertex', 'primal_edge', 'lon_i', 'lat_i', 'Ai' ) |
---|
358 | edge_deg, trisk_deg, left_right, down_up, trisk, le, de, wee, lon_e, lat_e, angle_e = get_pvars( |
---|
359 | dim_edge, 'edge_deg', 'trisk_deg', 'left_right', 'down_up', |
---|
360 | 'trisk', 'le', 'de', 'wee', 'lon_e', 'lat_e', 'angle_e') |
---|
361 | dual_deg, dual_vertex, dual_edge, Riv2, lon_v, lat_v, Av = get_pvars( |
---|
362 | dim_dual, 'dual_deg', 'dual_vertex', 'dual_edge', 'Riv2', 'lon_v', 'lat_v', 'Av') |
---|
363 | # Let indices start at 0 |
---|
364 | for x in (primal_vertex, primal_edge, dual_vertex,dual_edge, |
---|
365 | left_right, down_up, trisk ) : x.data = x.data-gridfile.start_index |
---|
366 | self.edge2cell = Stencil_glob(edge_deg, left_right) |
---|
367 | self.cell2edge = Stencil_glob(primal_deg, primal_edge) |
---|
368 | self.cell2vertex = Stencil_glob(primal_deg, primal_vertex) |
---|
369 | self.edge2vertex = Stencil_glob(edge_deg, down_up) |
---|
370 | self.vertex2edge = Stencil_glob(dual_deg, dual_edge) |
---|
371 | self.vertex2cell = Stencil_glob(dual_deg, dual_vertex, Riv2) |
---|
372 | self.edge2edge = Stencil_glob(trisk_deg, trisk, wee) |
---|
373 | print 'Partitioning ...' |
---|
374 | edge_owner = unst.partition_mesh(trisk_deg, trisk, comm.Get_size()) |
---|
375 | edge_owner = parallel.LocPArray1D(dim_edge, edge_owner) |
---|
376 | primal_owner = partition_from_stencil(edge_owner, primal_deg, primal_edge) |
---|
377 | primal_owner = parallel.LocPArray1D(dim_primal, primal_owner) |
---|
378 | self.set_members(locals(), 'dim_primal', 'dim_dual', 'dim_edge', |
---|
379 | 'edge_owner', 'primal_owner', 'primal_deg', 'primal_vertex', |
---|
380 | 'lon_v', 'lat_v', 'Av', 'lon_i', 'lat_i', 'Ai', |
---|
381 | 'le', 'de', 'lon_e', 'lat_e', 'angle_e') |
---|
382 | def get(self, getter, *names): return [ getter(self.__dict__[x]) for x in names ] |
---|
383 | |
---|
384 | def chain(start, links): |
---|
385 | for link in links: |
---|
386 | start = link(start).neigh_set |
---|
387 | yield start |
---|
388 | |
---|
389 | def patches(degree, bounds, lon, lat): |
---|
390 | for i in range(degree.size): |
---|
391 | nb_edge=degree[i] |
---|
392 | bounds_cell = bounds[i,0:nb_edge] |
---|
393 | lat_cell = lat[bounds_cell] |
---|
394 | lon_cell = lon[bounds_cell] |
---|
395 | orig=lon_cell[0] |
---|
396 | lon_cell = lon_cell-orig+180. |
---|
397 | lon_cell = np.mod(lon_cell,360.) |
---|
398 | lon_cell = lon_cell+orig-180. |
---|
399 | # if np.abs(lon_cell-orig).max()>10. : |
---|
400 | # print '%d patches :'%mpi_rank, lon_cell |
---|
401 | lonlat_cell = np.zeros((nb_edge,2)) |
---|
402 | lonlat_cell[:,0],lonlat_cell[:,1] = lon_cell,lat_cell |
---|
403 | polygon = Polygon(lonlat_cell, True) |
---|
404 | yield polygon |
---|
405 | |
---|
406 | class Local_Mesh(Abstract_Mesh): |
---|
407 | def __init__(self, pmesh, llm, nqdyn, radius, f): |
---|
408 | comm, dim_primal, primal_owner, dim_edge, edge_owner, dim_dual = pmesh.members( |
---|
409 | 'comm', 'dim_primal', 'primal_owner', 'dim_edge', 'edge_owner', 'dim_dual') |
---|
410 | print 'Constructing halos ...' |
---|
411 | edges_E0 = find_my_cells(edge_owner) |
---|
412 | cells_C0, edges_E1, vertices_V1, edges_E2, cells_C1 = chain( |
---|
413 | edges_E0, ( pmesh.edge2cell, pmesh.cell2edge, pmesh.edge2vertex, pmesh.vertex2edge, pmesh.edge2cell) ) |
---|
414 | edges_E0, edges_E1, edges_E2 = progressive_list(edges_E0, edges_E1, edges_E2) |
---|
415 | cells_C0, cells_C1 = progressive_list(cells_C0, cells_C1) |
---|
416 | print 'E2,E1,E0 ; C1,C0 : ', map(len, (edges_E2, edges_E1, edges_E0, cells_C1, cells_C0)) |
---|
417 | get_own_cells, get_all_cells = dim_primal.getter(cells_C0), dim_primal.getter(cells_C1) |
---|
418 | get_all_duals, get_all_edges = dim_dual.getter(vertices_V1), dim_edge.getter(edges_E2) |
---|
419 | |
---|
420 | self.lon_i, self.lat_i, self.Ai = pmesh.get(get_all_cells, 'lon_i', 'lat_i', 'Ai') |
---|
421 | self.lon_v, self.lat_v, self.Av = pmesh.get(get_all_duals, 'lon_v', 'lat_v', 'Av') |
---|
422 | self.fv = f(self.lon_v,self.lat_v) # Coriolis parameter |
---|
423 | self.le, self.de, self.lon_e, self.lat_e, self.angle_e = pmesh.get( |
---|
424 | get_all_edges, 'le', 'de', 'lon_e', 'lat_e', 'angle_e') |
---|
425 | primal_num, dual_num, edge_num = map(len, (cells_C1, vertices_V1, edges_E2)) |
---|
426 | |
---|
427 | # construct local stencils from global stencils |
---|
428 | dict_E2, dict_C1, dict_V1 = map(inverse_list, (edges_E2, cells_C1, vertices_V1)) |
---|
429 | down_up = pmesh.edge2vertex( edges_E2, dict_V1) |
---|
430 | left_right = pmesh.edge2cell( edges_E2, dict_C1) |
---|
431 | primal_edge = pmesh.cell2edge( cells_C1, dict_E2) |
---|
432 | dual_edge = pmesh.vertex2edge( vertices_V1, dict_E2) |
---|
433 | trisk = pmesh.edge2edge( edges_E2, dict_E2) |
---|
434 | Riv2 = pmesh.vertex2cell( vertices_V1, dict_C1) |
---|
435 | print 'Edge-vertex degree min/max %d %d'%(down_up.degree.min(), down_up.degree.max() ) |
---|
436 | print 'Edge-cell degree min/max %d %d'%(left_right.degree.min(), left_right.degree.max() ) |
---|
437 | print 'Primal degree min/max : %d %d'%(primal_edge.degree.min(), primal_edge.degree.max() ) |
---|
438 | print 'Dual degree min/max : %d %d'%(dual_edge.degree.min(), dual_edge.degree.max() ) |
---|
439 | print 'TRiSK degree min/max : %d %d'%(trisk.degree.min(), trisk.degree.max() ) |
---|
440 | print 'Riv2 degree min/max : %d %d'%(Riv2.degree.min(), Riv2.degree.max() ) |
---|
441 | primal_deg, primal_edge = primal_edge.degree, primal_edge.neigh_loc |
---|
442 | dual_deg, dual_edge = dual_edge.degree, dual_edge.neigh_loc |
---|
443 | trisk_deg, trisk, wee = trisk.degree, trisk.neigh_loc, trisk.weight |
---|
444 | dual_deg, dual_vertex, Riv2 = Riv2.degree, Riv2.neigh_loc, Riv2.weight |
---|
445 | edge_deg, left, right = left_right.degree, left_right.neigh_loc[:,0], left_right.neigh_loc[:,1] |
---|
446 | edge_deg, down, up = down_up.degree, down_up.neigh_loc[:,0], down_up.neigh_loc[:,1] |
---|
447 | for edge in range(edge_num): |
---|
448 | if edge_deg[edge]<2: up[edge]=down[edge] |
---|
449 | max_primal_deg, max_dual_deg, max_trisk_deg = [x.shape[1] for x in primal_edge, dual_edge, trisk] |
---|
450 | |
---|
451 | # construct own primal mesh for XIOS output |
---|
452 | primal_own_glo = find_my_cells(primal_owner) |
---|
453 | primal_vertex = pmesh.cell2vertex(primal_own_glo, dict_V1) |
---|
454 | primal_own_deg, primal_vertex = primal_vertex.degree, primal_vertex.neigh_loc |
---|
455 | primal_own_glo = np.asarray(primal_own_glo, dtype=np.int32) |
---|
456 | primal_own_loc = [dict_C1[i] for i in primal_own_glo] |
---|
457 | primal_own_num = len(primal_own_glo) |
---|
458 | primal_own_all = comm.allgather(primal_own_num) |
---|
459 | displ = sum(primal_own_all[0:comm.Get_rank()]) # required by XIOS |
---|
460 | print 'local primal mesh :', primal_vertex.min(), primal_vertex.max(), primal_own_all, displ |
---|
461 | |
---|
462 | # compute_ne(...) and normalize(...) expect indices starting at 1 |
---|
463 | trisk, dual_vertex, primal_edge, dual_edge = trisk+1, dual_vertex+1, primal_edge+1, dual_edge+1 |
---|
464 | left, right, down, up = left+1, right+1, down+1, up+1 |
---|
465 | primal_ne = compute_ne(primal_num,primal_deg,primal_edge,left,right) |
---|
466 | dual_ne = compute_ne(dual_num,dual_deg,dual_edge,up,down) |
---|
467 | |
---|
468 | # store what we have computed |
---|
469 | self.set_members(locals(), 'nqdyn', 'llm', 'primal_num', 'dual_num', 'edge_num', |
---|
470 | 'primal_own_glo', 'primal_own_loc', 'primal_own_deg', |
---|
471 | 'primal_deg', 'primal_vertex', 'displ', |
---|
472 | 'trisk_deg', 'trisk', 'wee', 'dual_deg', 'dual_vertex', 'Riv2') |
---|
473 | |
---|
474 | # normalize and propagate info to Fortran side |
---|
475 | pmesh.gridfile.normalize(self,radius) |
---|
476 | unst.init_mesh(llm, nqdyn, self.edge_num, self.primal_num, self.dual_num, |
---|
477 | max_trisk_deg, max_primal_deg, max_dual_deg, |
---|
478 | primal_deg,primal_edge,primal_ne, |
---|
479 | dual_deg,dual_edge,dual_ne,dual_vertex, |
---|
480 | left,right,down,up,trisk_deg,trisk, |
---|
481 | self.Ai,self.Av,self.fv,self.le/self.de,Riv2,wee) |
---|
482 | |
---|
483 | # setup halo transfers - NB : llm must be set before we call set_dynamico_transfer |
---|
484 | self.com_primal = parallel.Halo_Xchange( |
---|
485 | 42, dim_primal, cells_C1, get_all_cells(primal_owner)) |
---|
486 | self.com_edges = parallel.Halo_Xchange( |
---|
487 | 73, dim_edge, edges_E2, get_all_edges(edge_owner)) |
---|
488 | self.com_primal.set_dynamico_transfer('primal') |
---|
489 | self.com_edges.set_dynamico_transfer('edge') |
---|
490 | |
---|
491 | self.primal = tri.Triangulation(self.lon_i*radian, self.lat_i*radian) |
---|
492 | self.dual = tri.Triangulation(self.lon_v*radian, self.lat_v*radian) |
---|
493 | self.triedge = tri.Triangulation(self.lon_e*radian, self.lat_e*radian) |
---|
494 | |
---|
495 | def plot_patches(self, ax, clim, degree, bounds, lon, lat, data): |
---|
496 | nb_vertex = lon.size # global |
---|
497 | p = list(patches(degree, bounds, lon, lat)) |
---|
498 | p = PatchCollection(p, linewidth=0.01) |
---|
499 | p.set_array(data) # set values at each polygon (cell) |
---|
500 | p.set_clim(clim) |
---|
501 | ax.add_collection(p) |
---|
502 | |
---|
503 | #-------------------------------------- Mesh reordering ------------------------------------------# |
---|
504 | |
---|
505 | # NB : Indices start at 0 |
---|
506 | # permutations map an old index to a new index : perm[old_index]=new_index |
---|
507 | # enumerate(perm) is the list [ (0,perm[0]) , (1,perm[1]), ...] = [(old, new) ... ] |
---|
508 | |
---|
509 | def reorder_values(perm, *datas): |
---|
510 | # data in datas : ndarray containing indices ; perm : permutation to apply to values |
---|
511 | for data in datas: |
---|
512 | for x in np.nditer(data, op_flags=['readwrite']): |
---|
513 | x[...]=perm[x] |
---|
514 | |
---|
515 | def reorder_indices(perm, *datas): |
---|
516 | # data in datas : ndarray containing values ; perm : permutation to apply to indices |
---|
517 | for data in datas: |
---|
518 | data_out, ndim = data.copy(), data.ndim |
---|
519 | if ndim == 1: |
---|
520 | for old,new in enumerate(perm): |
---|
521 | data_out[new]=data[old] |
---|
522 | if ndim == 2: |
---|
523 | for old,new in enumerate(perm): |
---|
524 | data_out[new,:]=data[old,:] |
---|
525 | data[:]=data_out[:] |
---|
526 | |
---|
527 | def key_to_perm(keys): |
---|
528 | # keys : maps an old index to a key |
---|
529 | # returns : a list mapping old indices to new indices obtained by sorting keys in ascending order |
---|
530 | ordered = [old for key,old in sorted(zip(keys,range(len(keys))))] |
---|
531 | # ordered maps new indices to old indices |
---|
532 | # but we want to map old indices to new indices |
---|
533 | perm = list(range(len(keys))) |
---|
534 | for new,old in enumerate(ordered): |
---|
535 | perm[old]=new |
---|
536 | return perm # maps old indices to new indices |
---|
537 | |
---|
538 | def toint(x): return np.int32((x+1.)*65536) |
---|
539 | |
---|
540 | def morton_index(x,y,z): |
---|
541 | nn, xx, yy, zz = x.size, toint(x), toint(y), toint(z) |
---|
542 | idx=np.zeros(nn, dtype=np.int64) |
---|
543 | unst.ker.dynamico_morton_encode(nn, xx,yy,zz, idx) |
---|
544 | return idx |
---|
545 | |
---|
546 | #-------------------------------------- Mesh partitioning ------------------------------------------# |
---|
547 | |
---|
548 | # NB : Indices start at 0 |
---|
549 | |
---|
550 | def partition_from_stencil(owner2, degree, stencil): |
---|
551 | # 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 |
---|
552 | dim1, dim2= degree.dim, owner2.dim |
---|
553 | degree, stencil, n = degree.data, stencil.data, dim1.end-dim1.start |
---|
554 | cells2 = list_stencil(degree, stencil) |
---|
555 | cells2 = sorted(list(set(list(cells2)))) # list of cells for which we need to know owner2 |
---|
556 | get2 = dim2.getter(cells2) |
---|
557 | owner1, owner2, glob2loc = np.zeros(n, dtype=np.int32), get2(owner2), get2.dict |
---|
558 | for i in range(n): |
---|
559 | owners = [ owner2[glob2loc[stencil[i,j]]] for j in range(degree[i]) ] |
---|
560 | if i%2 == 0 : owner1[i] = min(owners) |
---|
561 | else : owner1[i] = max(owners) |
---|
562 | return owner1 |
---|
563 | |
---|
564 | def find_my_cells(owner): # owner : a PArray1D containing the data returned by partition_mesh |
---|
565 | dim, comm, owner = owner.dim, owner.dim.comm, owner.data |
---|
566 | mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size() |
---|
567 | cells=[set() for i in range(mpi_size)] |
---|
568 | for i in range(len(owner)): |
---|
569 | cells[owner[i]].add(dim.start+i) |
---|
570 | cells = [sorted(list(cells[i])) for i in range(mpi_size)] |
---|
571 | mycells = comm.alltoall(cells) |
---|
572 | mycells = sorted(sum(mycells, [])) # concatenate into a single list |
---|
573 | return mycells |
---|
574 | |
---|
575 | #---------------------------------- Stencil management ---------------------------------------# |
---|
576 | |
---|
577 | # NB : Indices start at 0 |
---|
578 | |
---|
579 | # Class Stencil represents an adjacency relationship (e.g. cell=>edges) |
---|
580 | # using adjacency information read from PArrays |
---|
581 | # It computes a list of "edges" adjacent to a given list of "cells" |
---|
582 | # This is used to form the sets E0 -> C0 -> E1 -> V1 -> E2 -> C1 |
---|
583 | # which are then used to form lists of global indices for V1,C1,E2 such that |
---|
584 | # C0, E0, E1 form contiguous subsets of C1, E2 starting from 0 |
---|
585 | |
---|
586 | def reindex(vertex_dict, degree, bounds): |
---|
587 | for i in range(degree.size): |
---|
588 | for j in range(degree[i]): |
---|
589 | bounds[i,j] = vertex_dict[bounds[i,j]] |
---|
590 | return bounds |
---|
591 | |
---|
592 | class Stencil_glob: |
---|
593 | def __init__(self, degree, neigh, weight=None): |
---|
594 | self.degree, self.neigh, self.weight = degree, neigh, weight |
---|
595 | def __call__(self, cells, neigh_dict=None): |
---|
596 | return Stencil(cells, self.degree, self.neigh, neigh_dict, self.weight) |
---|
597 | |
---|
598 | class Stencil: |
---|
599 | def __init__(self, cells, degree, neigh, neigh_dict, weight=None): |
---|
600 | get = degree.dim.getter(cells) |
---|
601 | mydegree, myneigh = [get(x) for x in (degree, neigh)] |
---|
602 | if not weight is None : myweight = get(weight) |
---|
603 | if neigh_dict is None : |
---|
604 | keep = lambda n : True |
---|
605 | else : # if neigh_dict is present, only neighbours in dict are retained |
---|
606 | keep = lambda n : n in neigh_dict |
---|
607 | neigh_set = list_stencil(mydegree, myneigh, keep) |
---|
608 | self.neigh_set = list(set(list(neigh_set) )) |
---|
609 | rej=0 |
---|
610 | for i in range(len(mydegree)): # keep only elements in neigh_dict, in-place |
---|
611 | k=0 |
---|
612 | for j in range(mydegree[i]): |
---|
613 | n=myneigh[i,j] |
---|
614 | if keep(n): |
---|
615 | myneigh[i,k]=n |
---|
616 | if not weight is None : myweight[i,k]=myweight[i,j] |
---|
617 | k=k+1 |
---|
618 | else: |
---|
619 | rej=rej+1 |
---|
620 | mydegree[i]=k |
---|
621 | if neigh_dict is None: |
---|
622 | neigh_dict = {j:i for i,j in enumerate(self.neigh_set)} |
---|
623 | myneigh_loc = reindex(neigh_dict, mydegree, myneigh) |
---|
624 | self.degree, self.neigh_glob, self.neigh_loc = mydegree, myneigh, myneigh_loc |
---|
625 | if weight is not None: self.weight = myweight |
---|
626 | |
---|
627 | def progressive_iter(mylist, cell_lists): |
---|
628 | for thelist in cell_lists: |
---|
629 | mylist = mylist + list(set(thelist)-set(mylist)) |
---|
630 | yield mylist |
---|
631 | def progressive_list(*cell_lists): |
---|
632 | # cell_lists : a tuple of lists of global indices, with each list a subset of the next |
---|
633 | # returns : a list 'mylist' such that for each list 'thelist' in cell_lists, thelist = mylist[0:len(thelist)] |
---|
634 | # example : edge_list = progressive_list(E0,E1,E2) with E0,E1,E2 increasing lists of cell edges |
---|
635 | return list(progressive_iter([], cell_lists)) |
---|