[615] | 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 | |
---|
| 8 | import dynamico.wrap as wrap |
---|
| 9 | from ctypes import c_void_p, c_int, c_double, c_bool |
---|
| 10 | radian=180/math.pi |
---|
| 11 | |
---|
| 12 | #------------- direct Cython interface to DYNAMICO routines -------------# |
---|
| 13 | |
---|
| 14 | cdef extern from "functions.h": |
---|
[618] | 15 | cdef void dynamico_init_params() |
---|
| 16 | cpdef void dynamico_setup_xios() |
---|
| 17 | cpdef void dynamico_xios_set_timestep(double) |
---|
| 18 | cpdef void dynamico_xios_update_calendar(int) |
---|
[615] | 19 | |
---|
| 20 | #------------- import and wrap DYNAMICO routines -------------# |
---|
| 21 | |
---|
| 22 | ker=wrap.Struct() # store imported fun X as funs.X |
---|
| 23 | |
---|
| 24 | check_args = False # use True instead of False for debugging, probably with some overhead |
---|
| 25 | |
---|
| 26 | try: |
---|
[618] | 27 | kernels = wrap.SharedLib(vars(ker), 'libicosa.so', check_args=check_args) |
---|
[615] | 28 | setvar, setvars, getvar, getvars = kernels.setvar, kernels.setvars, kernels.getvar, kernels.getvars |
---|
| 29 | except OSError: |
---|
| 30 | print """ |
---|
| 31 | Unable to load shared library 'libkernels.so' ! |
---|
| 32 | """ |
---|
| 33 | raise |
---|
| 34 | |
---|
| 35 | # providing a full prototype enables type-checking when calling |
---|
| 36 | # if a number n is present in the prototype, the previous type is repeated n times |
---|
| 37 | kernels.import_funs([ |
---|
| 38 | # ['setup_xios',None], |
---|
| 39 | # ['call_xios_set_timestep',c_double], |
---|
| 40 | # ['call_xios_update_calendar',c_int] |
---|
| 41 | # ['init_params',c_double], |
---|
[618] | 42 | ['dynamico_init_mesh',c_void_p,13], |
---|
| 43 | ['dynamico_init_metric', c_void_p,6], |
---|
| 44 | ['dynamico_init_hybrid', c_void_p,3], |
---|
| 45 | ['dynamico_caldyn_unstructured', c_double, c_void_p,20], |
---|
| 46 | # ['dynamico_partition_graph', c_int,2, c_void_p,3, c_int, c_void_p], |
---|
[615] | 47 | ]) |
---|
| 48 | |
---|
| 49 | # set/get global variables |
---|
| 50 | eta_mass,eta_lag=(1,2) |
---|
| 51 | thermo_theta,thermo_entropy,thermo_moist,thermo_boussinesq=(1,2,3,4) |
---|
| 52 | |
---|
| 53 | kernels.addvars( |
---|
| 54 | c_bool,'hydrostatic','debug_hevi_solver','rigid', |
---|
| 55 | c_int,'llm','nqdyn','primal_num','max_primal_deg', |
---|
| 56 | 'dual_num','max_dual_deg','edge_num','max_trisk_deg', |
---|
| 57 | 'caldyn_thermo','caldyn_eta','nb_threads', |
---|
| 58 | c_double,'elapsed','g', 'ptop', 'cpp', 'cppv', |
---|
| 59 | 'Rd', 'Rv', 'preff', 'Treff', 'pbot', 'rho_bot', 'Phi_bot') |
---|
| 60 | |
---|
| 61 | elapsed=0. |
---|
| 62 | |
---|
| 63 | #----------------------------- Base class for dynamics ------------------------ |
---|
| 64 | |
---|
| 65 | class Caldyn: |
---|
| 66 | def __init__(self,mesh): |
---|
| 67 | self.mesh=mesh |
---|
| 68 | fps, ftheta, fmass = mesh.field_ps, mesh.field_theta, mesh.field_mass |
---|
| 69 | fw, fu, fz = mesh.field_w, mesh.field_u, mesh.field_z |
---|
| 70 | self.ps, self.ms, self.dms = fps(), fps(), fps() |
---|
| 71 | self.s, self.hs, self.dhs = ftheta(), ftheta(), ftheta() |
---|
| 72 | self.pk, self.berni, self.geopot, self.hflux = fmass(),fmass(),fw(),fu() |
---|
| 73 | self.qu, self.qv = fu(),fz() |
---|
| 74 | self.fmass, self.ftheta, self.fu, self.fw = fmass, ftheta, fu, fw |
---|
| 75 | def bwd_fast_slow(self, flow, tau): |
---|
| 76 | global elapsed |
---|
| 77 | time1=time.time() |
---|
| 78 | flow,fast,slow = self._bwd_fast_slow_(flow,tau) |
---|
| 79 | time2=time.time() |
---|
| 80 | elapsed=elapsed+time2-time1 |
---|
| 81 | return flow,fast,slow |
---|
| 82 | |
---|
[618] | 83 | # when calling caldyn_unstructured, arrays for tendencies must be re-created each time |
---|
[615] | 84 | # to avoid overwriting in the same memory space when time scheme is multi-stage |
---|
| 85 | |
---|
| 86 | #-------------------------- Shallow-water dynamics --------------------- |
---|
| 87 | |
---|
| 88 | class Caldyn_RSW(Caldyn): |
---|
| 89 | def __init__(self,mesh): |
---|
| 90 | Caldyn.__init__(self,mesh) |
---|
| 91 | setvars(('hydrostatic','caldyn_thermo','caldyn_eta'), |
---|
| 92 | (True,thermo_boussinesq,eta_lag)) |
---|
| 93 | self.dhs = self.fmass() |
---|
[618] | 94 | dynamico_init_params() |
---|
[615] | 95 | def _bwd_fast_slow_(self, flow, tau): |
---|
| 96 | h,u = flow |
---|
| 97 | # h*s = h => uniform buoyancy s=1 => shallow-water |
---|
| 98 | dh, du_slow, du_fast, hs, buf = self.fmass(), self.fu(), self.fu(), h.copy(), self.geopot |
---|
[618] | 99 | ker.dynamico_caldyn_unstructured(tau, self.ms, h, hs, u, self.geopot, buf, |
---|
[615] | 100 | self.s, self.ps, self.pk, self.hflux, self.qv, |
---|
| 101 | self.dms, dh, self.dhs, du_fast, du_slow, |
---|
| 102 | buf, buf, buf, buf) |
---|
| 103 | return (h,u), (0.,du_fast), (dh,du_slow) |
---|
| 104 | |
---|
| 105 | #----------------------------------- HPE ------------------------------------ |
---|
| 106 | |
---|
| 107 | class Caldyn_HPE(Caldyn): |
---|
| 108 | def __init__(self,caldyn_thermo,caldyn_eta, mesh,metric,thermo,BC,g): |
---|
| 109 | Caldyn.__init__(self,mesh) |
---|
| 110 | setvars(('hydrostatic','caldyn_thermo','caldyn_eta', |
---|
| 111 | 'g','ptop','Rd','cpp','preff','Treff'), |
---|
| 112 | (True,caldyn_thermo,caldyn_eta, |
---|
| 113 | g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0)) |
---|
[618] | 114 | dynamico_init_params() |
---|
[615] | 115 | def _bwd_fast_slow_(self, flow, tau): |
---|
| 116 | dm, dS, du_slow, du_fast, buf = self.fmass(), self.ftheta(), self.fu(), self.fu(), self.geopot |
---|
| 117 | m,S,u = flow |
---|
[618] | 118 | ker.dynamico_caldyn_unstructured(tau, self.ms, m, S, u, self.geopot, buf, |
---|
[615] | 119 | self.s, self.ps, self.pk, self.hflux, self.qv, |
---|
| 120 | self.dms, dm, dS, du_fast, du_slow, |
---|
| 121 | buf, buf, buf, buf) |
---|
| 122 | return (m,S,u), (0.,0.,du_fast), (dm,dS,du_slow) |
---|
| 123 | |
---|
| 124 | #----------------------------------- NH ------------------------------------ |
---|
| 125 | |
---|
| 126 | class Caldyn_NH(Caldyn): |
---|
| 127 | def __init__(self,caldyn_thermo,caldyn_eta, mesh,metric,thermo,BC,g): |
---|
| 128 | Caldyn.__init__(self,mesh) |
---|
| 129 | setvars(('hydrostatic','caldyn_thermo','caldyn_eta', |
---|
| 130 | 'g','ptop','Rd','cpp','preff','Treff', |
---|
| 131 | 'pbot','rho_bot'), |
---|
| 132 | (False,caldyn_thermo,caldyn_eta, |
---|
| 133 | g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0, |
---|
| 134 | BC.pbot.max(), BC.rho_bot.max())) |
---|
[618] | 135 | dynamico_init_params() |
---|
[615] | 136 | def bwd_fast_slow(self, flow, tau): |
---|
| 137 | ftheta, fmass, fu, fw = self.ftheta, self.fmass, self.fu, self.fw |
---|
| 138 | dm, dS, du_slow, du_fast = fmass(), ftheta(), fu(), fu() |
---|
| 139 | dPhi_slow, dPhi_fast, dW_slow, dW_fast = fw(), fw(), fw(), fw() |
---|
| 140 | m,S,u,Phi,W = flow |
---|
[618] | 141 | ker.dynamico_caldyn_unstructured(tau, self.ms, m, S, u, Phi, W, |
---|
[615] | 142 | self.s, self.ps, self.pk, self.hflux, self.qv, |
---|
| 143 | self.dms, dm, dS, du_fast, du_slow, |
---|
| 144 | dPhi_fast, dPhi_slow, dW_fast, dW_slow) |
---|
| 145 | return ((m,S,u,Phi,W), (0.,0.,du_fast,dPhi_fast,dW_fast), |
---|
| 146 | (dm,dS,du_slow,dPhi_slow,dW_slow)) |
---|
| 147 | |
---|
| 148 | #-------------------------------- Hybrid mass-based coordinate ------------- |
---|
| 149 | |
---|
| 150 | # compute hybrid coefs from distribution of mass |
---|
| 151 | def compute_hybrid_coefs(mass): |
---|
| 152 | nx,llm=mass.shape |
---|
| 153 | mass_dak = np.zeros((nx,llm)) |
---|
| 154 | mass_dbk = np.zeros((nx,llm)) |
---|
| 155 | mass_bl = np.zeros((nx,llm+1))+1. |
---|
| 156 | for i in range(nx): |
---|
| 157 | m_i = mass[i,:] |
---|
| 158 | dbk_i = m_i/sum(m_i) |
---|
| 159 | mass_dbk[i,:] = dbk_i |
---|
| 160 | mass_bl[i,1:]= 1.-dbk_i.cumsum() |
---|
| 161 | return mass_bl, mass_dak, mass_dbk |
---|
| 162 | |
---|
| 163 | #----------------------- Cartesian mesh ----------------------- |
---|
| 164 | |
---|
| 165 | def squeeze(dims): |
---|
| 166 | # return np.zeros(dims) |
---|
| 167 | return np.zeros([n for n in dims if n>1]) |
---|
| 168 | |
---|
| 169 | # arrays is a list of arrays |
---|
| 170 | # vals is a list of tuples |
---|
| 171 | # each tuple is stored in each array |
---|
| 172 | def put(ij, deg, arrays, vals): |
---|
| 173 | k=0 |
---|
| 174 | for vv in vals: # vv is a tuple of values to be stored in arrays |
---|
| 175 | for array,v in zip(arrays,vv): |
---|
| 176 | array[ij,k]=v |
---|
| 177 | k=k+1 |
---|
| 178 | deg[ij]=k |
---|
| 179 | |
---|
| 180 | class Cartesian_mesh: |
---|
| 181 | def __init__(self,nx,ny,llm,nqdyn,Lx,Ly,f): |
---|
| 182 | dx,dy = Lx/nx, Ly/ny |
---|
| 183 | self.dx, self.dy, self.f = dx,dy,f |
---|
| 184 | self.nx, self.ny, self.llm, self.nqdyn = nx,ny,llm,nqdyn |
---|
| 185 | self.field_z = self.field_mass |
---|
| 186 | # 1D coordinate arrays |
---|
| 187 | x=(np.arange(nx)-nx/2.)*Lx/nx |
---|
| 188 | y=(np.arange(ny)-ny/2.)*Ly/ny |
---|
| 189 | lev=np.arange(llm) |
---|
| 190 | levp1=np.arange(llm+1) |
---|
| 191 | self.x, self.y, self.lev, self.levp1 = x,y,lev,levp1 |
---|
| 192 | # 3D coordinate arrays |
---|
| 193 | self.xx,self.yy,self.ll = np.meshgrid(x,y,lev, indexing='ij') |
---|
| 194 | self.xxp1,self.yyp1,self.llp1 = np.meshgrid(x,y,levp1, indexing='ij') |
---|
| 195 | # beware conventions for indexing |
---|
| 196 | # Fortran order : llm,nx*ny,nqdyn / indices start at 1 |
---|
| 197 | # Python order : nqdyn,ny,nx,llm / indices start at 0 |
---|
| 198 | # indices below follow Fortran while x,y follow Python/C |
---|
| 199 | index=lambda x,y : ((x+(nx*(y+2*ny)))%(nx*ny))+1 |
---|
| 200 | indexu=lambda x,y : 2*index(x,y)-1 |
---|
| 201 | indexv=lambda x,y : 2*index(x,y) |
---|
| 202 | indices = lambda shape : np.zeros(shape,np.int32) |
---|
| 203 | |
---|
| 204 | primal_nb = indices(nx*ny) |
---|
| 205 | primal_edge = indices((nx*ny,4)) |
---|
| 206 | primal_ne = indices((nx*ny,4)) |
---|
| 207 | |
---|
| 208 | dual_nb = indices(nx*ny) |
---|
| 209 | dual_edge = indices((nx*ny,4)) |
---|
| 210 | dual_ne = indices((nx*ny,4)) |
---|
| 211 | dual_vertex = indices((nx*ny,4)) |
---|
| 212 | Riv2 = np.zeros((nx*ny,4)) |
---|
| 213 | |
---|
| 214 | left = indices(2*nx*ny) |
---|
| 215 | right = indices(2*nx*ny) |
---|
| 216 | up = indices(2*nx*ny) |
---|
| 217 | down = indices(2*nx*ny) |
---|
| 218 | le_de = np.zeros(2*nx*ny) |
---|
| 219 | |
---|
| 220 | trisk_deg = indices(2*nx*ny) |
---|
| 221 | trisk = indices((2*nx*ny,4)) # 4 TRiSK coefs per edge |
---|
| 222 | wee = np.zeros((2*nx*ny,4)) |
---|
| 223 | |
---|
| 224 | for x in range(nx): |
---|
| 225 | for y in range(ny): |
---|
| 226 | # NB : Fortran indices start at 1 |
---|
| 227 | # primal cells |
---|
| 228 | put(index(x,y)-1,primal_nb,(primal_edge,primal_ne), |
---|
| 229 | ((indexu(x,y),1), |
---|
| 230 | (indexv(x,y),1), |
---|
| 231 | (indexu(x-1,y),-1), |
---|
| 232 | (indexv(x,y-1),-1) )) |
---|
| 233 | # dual cells |
---|
| 234 | put(index(x,y)-1,dual_nb,(dual_edge,dual_vertex,dual_ne,Riv2), |
---|
| 235 | ((indexv(x+1,y),index(x,y),1,.25), |
---|
| 236 | (indexu(x,y+1),index(x+1,y),-1,.25), |
---|
| 237 | (indexv(x,y),index(x+1,y+1),-1,.25), |
---|
| 238 | (indexu(x,y),index(x,y+1),1,.25) )) |
---|
| 239 | # edges : |
---|
| 240 | # left and right are adjacent primal cells |
---|
| 241 | # flux is positive when going from left to right |
---|
| 242 | # up and down are adjacent dual cells (vertices) |
---|
| 243 | # circulation is positive when going from down to up |
---|
| 244 | # u-points |
---|
| 245 | ij =indexu(x,y)-1 |
---|
| 246 | left[ij]=index(x,y) |
---|
| 247 | right[ij]=index(x+1,y) |
---|
| 248 | down[ij]=index(x,y-1) |
---|
| 249 | up[ij]=index(x,y) |
---|
| 250 | le_de[ij]=dy/dx |
---|
| 251 | put(ij,trisk_deg,(trisk,wee),( |
---|
| 252 | (indexv(x,y),-.25), |
---|
| 253 | (indexv(x+1,y),-.25), |
---|
| 254 | (indexv(x,y-1),-.25), |
---|
| 255 | (indexv(x+1,y-1),-.25))) |
---|
| 256 | # v-points |
---|
| 257 | ij = indexv(x,y)-1 |
---|
| 258 | left[ij]=index(x,y) |
---|
| 259 | right[ij]=index(x,y+1) |
---|
| 260 | down[ij]=index(x,y) |
---|
| 261 | up[ij]=index(x-1,y) |
---|
| 262 | le_de[ij]=dx/dy |
---|
| 263 | put(ij,trisk_deg,(trisk,wee),( |
---|
| 264 | (indexu(x,y),.25), |
---|
| 265 | (indexu(x-1,y),.25), |
---|
| 266 | (indexu(x,y+1),.25), |
---|
| 267 | (indexu(x-1,y+1),.25))) |
---|
| 268 | setvars(('llm','nqdyn','edge_num','primal_num','dual_num', |
---|
| 269 | 'max_trisk_deg','max_primal_deg','max_dual_deg'), |
---|
| 270 | (llm,nqdyn,2*nx*ny,nx*ny,nx*ny,4,4,4) ) |
---|
[618] | 271 | print('init_mesh ...') |
---|
| 272 | ker.dynamico_init_mesh(primal_nb,primal_edge,primal_ne, |
---|
[615] | 273 | dual_nb,dual_edge,dual_ne,dual_vertex, |
---|
| 274 | left,right,down,up,trisk_deg,trisk) |
---|
[618] | 275 | print ('...done') |
---|
| 276 | |
---|
[615] | 277 | Aiv=np.zeros(nx*ny)+dx*dy # Ai=Av=dx*dy |
---|
[618] | 278 | |
---|
| 279 | print('init_metric ...') |
---|
| 280 | ker.dynamico_init_metric(Aiv,Aiv,f+0.*Aiv,le_de,Riv2,-wee) |
---|
| 281 | print ('...done') |
---|
[615] | 282 | def field_theta(self): return squeeze((self.nqdyn,self.ny,self.nx,self.llm)) |
---|
| 283 | def field_mass(self): return squeeze((self.ny,self.nx,self.llm)) |
---|
| 284 | def field_z(self): return squeeze((self.ny,self.nx,self.llm)) |
---|
| 285 | def field_w(self): return squeeze((self.ny,self.nx,self.llm+1)) |
---|
| 286 | def field_u(self): return np.zeros((self.ny,2*self.nx,self.llm)) |
---|
| 287 | def field_ps(self): return squeeze((self.ny,self.nx)) |
---|
| 288 | def ucomp(self,u): return u[:,range(0,2*self.nx,2),:] |
---|
| 289 | def set_ucomp(self,uv,u): uv[:,range(0,2*self.nx,2),:]=u |
---|
| 290 | def vcomp(self,u): return u[:,range(1,2*self.nx,2),:] |
---|
| 291 | |
---|
| 292 | #---------------------- MPAS fully unstructured mesh ------------------------ |
---|
| 293 | |
---|
| 294 | def compute_ne(num,deg,edges,left,right): |
---|
| 295 | ne = np.zeros_like(edges) |
---|
| 296 | for cell in range(num): |
---|
| 297 | for iedge in range(deg[cell]): |
---|
| 298 | edge = edges[cell,iedge]-1 |
---|
| 299 | if left[edge]==cell+1: ne[cell,iedge]+=1 |
---|
| 300 | if right[edge]==cell+1: ne[cell,iedge]-=1 |
---|
| 301 | if ne[cell,iedge]==0 : print 'error at cell,iedge', cell, iedge |
---|
| 302 | return ne |
---|
| 303 | |
---|
| 304 | def plot(tri,data): |
---|
| 305 | plt.figure(figsize=(12,4)) |
---|
| 306 | plt.gca().set_aspect('equal') |
---|
| 307 | plt.tricontourf(tri, data, 20) |
---|
| 308 | plt.colorbar() |
---|
| 309 | plt.ylim((-90,90)) |
---|
| 310 | plt.xlim((0,360)) |
---|
| 311 | |
---|
| 312 | class MPAS_Mesh: |
---|
| 313 | def __init__(self, gridfile, llm, nqdyn, radius, f): |
---|
| 314 | self.gridfile, self.llm, self.nqdyn = gridfile,llm,nqdyn |
---|
| 315 | self.radius, self.f = radius, f |
---|
| 316 | # open mesh file, get main dimensions |
---|
| 317 | try: |
---|
| 318 | nc = cdf.Dataset(gridfile, "r") |
---|
| 319 | except RuntimeError: |
---|
| 320 | print """ |
---|
| 321 | Unable to open grid file %s, maybe you forgot to download it ? |
---|
| 322 | To do so, go to the 'Python/' dir and execute './get_MPAS_grids.sh'. |
---|
| 323 | """ % gridfile |
---|
| 324 | raise |
---|
| 325 | |
---|
| 326 | def getdims(*names): return [len(nc.dimensions[name]) for name in names] |
---|
| 327 | def getvars(*names): |
---|
| 328 | for name in names : print "getvar %s ..."%name |
---|
| 329 | time1=time.time() |
---|
| 330 | ret=[nc.variables[name][:] for name in names] |
---|
| 331 | print "... Done in %f seconds"%(time.time()-time1) |
---|
| 332 | return ret |
---|
| 333 | primal_num, edge_num, dual_num = getdims('nCells','nEdges','nVertices') |
---|
| 334 | print 'Number of primal cells, dual cells and edges : %d, %d, %d '%(primal_num,dual_num,edge_num) |
---|
| 335 | primal_deg, trisk_deg = getvars('nEdgesOnCell','nEdgesOnEdge') |
---|
| 336 | dual_deg = [3 for i in range(dual_num) ] |
---|
| 337 | dual_deg = np.ascontiguousarray(dual_deg,dtype=np.int32) # NB : Fortran code expects 32-bit ints |
---|
| 338 | |
---|
| 339 | # get indices for stencils |
---|
| 340 | # primal -> vertices (unused) |
---|
| 341 | primal_vertex, dual_vertex = getvars('verticesOnCell','cellsOnVertex') |
---|
| 342 | # primal <-> edges |
---|
| 343 | primal_edge, left_right = getvars('edgesOnCell','cellsOnEdge') |
---|
| 344 | left,right = left_right[:,0], left_right[:,1] |
---|
| 345 | primal_ne = compute_ne(primal_num,primal_deg,primal_edge,left,right) |
---|
| 346 | # dual <-> edges |
---|
| 347 | dual_edge, down_up = getvars('edgesOnVertex','verticesOnEdge') |
---|
| 348 | down,up = down_up[:,0], down_up[:,1] |
---|
| 349 | dual_ne = compute_ne(dual_num,dual_deg,dual_edge,up,down) |
---|
| 350 | # primal <-> dual, edges <-> edges (TRiSK) |
---|
| 351 | dual_vertex, trisk = getvars('cellsOnVertex','edgesOnEdge') |
---|
| 352 | # get positions, lengths, surfaces and weights |
---|
| 353 | le,de,Ai,Av = getvars('dvEdge','dcEdge','areaCell','areaTriangle') |
---|
| 354 | lat_i,lon_i = getvars('latCell','lonCell') |
---|
| 355 | lat_v,lon_v = getvars('latVertex','lonVertex') |
---|
| 356 | lat_e,lon_e,angle_e = getvars('latEdge','lonEdge','angleEdge') |
---|
| 357 | wee,Riv2 = getvars('weightsOnEdge','kiteAreasOnVertex') |
---|
| 358 | |
---|
| 359 | # fix normalization of wee and Riv2 weights |
---|
| 360 | for edge1 in range(edge_num): |
---|
| 361 | for i in range(trisk_deg[edge1]): |
---|
| 362 | edge2=trisk[edge1,i]-1 # NB Fortran vs C/Python indexing |
---|
| 363 | wee[edge1,i] = de[edge1]*wee[edge1,i]/le[edge2] |
---|
| 364 | for ivertex in range(dual_num): |
---|
| 365 | for j in range(3): |
---|
| 366 | Riv2[ivertex,j]=Riv2[ivertex,j]/Av[ivertex] |
---|
| 367 | r=Riv2[ivertex,:] |
---|
| 368 | r=sum(r) |
---|
| 369 | if abs(r-1.)>1e-6 : print 'error with Riv2 at vertex ', ivertex |
---|
| 370 | |
---|
| 371 | max_primal_deg, max_dual_deg, max_trisk_deg = getdims('maxEdges','vertexDegree','maxEdges2') |
---|
| 372 | |
---|
| 373 | # CRITICAL : force arrays left, etc. to be contiguous in memory: |
---|
| 374 | left,right,up,down = [np.ascontiguousarray(x) for x in (left,right,up,down)] |
---|
| 375 | trisk,wee,primal_edge = [np.ascontiguousarray(x) for x in (trisk,wee,primal_edge)] |
---|
| 376 | |
---|
| 377 | print ('Max stencil sizes (div,curl,trisk) : %d, %d, %d' |
---|
| 378 | % (max_primal_deg, max_dual_deg, max_trisk_deg) ) |
---|
| 379 | |
---|
| 380 | r2 = radius**2 |
---|
| 381 | Av = r2*Av |
---|
| 382 | fv = f(lon_v,lat_v) # Coriolis parameter |
---|
| 383 | self.Ai, self.Av, self.fv = r2*Ai,Av,fv |
---|
| 384 | self.le, self.de, self.le_de = radius*le, radius*de, le/de |
---|
| 385 | self.trisk_deg, self.trisk, self.wee = trisk_deg, trisk, wee |
---|
| 386 | setvars(('llm','nqdyn','edge_num','primal_num','dual_num', |
---|
| 387 | 'max_trisk_deg','max_primal_deg','max_dual_deg'), |
---|
| 388 | (llm, nqdyn, edge_num, primal_num,dual_num, |
---|
| 389 | max_trisk_deg, max_primal_deg, max_dual_deg) ) |
---|
| 390 | |
---|
[618] | 391 | ker.dynamico_init_mesh(primal_deg,primal_edge,primal_ne, |
---|
[615] | 392 | dual_deg,dual_edge,dual_ne,dual_vertex, |
---|
| 393 | left,right,down,up,trisk_deg,trisk) |
---|
[618] | 394 | ker.dynamico_init_metric(self.Ai,self.Av,self.fv,le/de,Riv2,wee) |
---|
[615] | 395 | |
---|
| 396 | for edge in range(edge_num): |
---|
| 397 | iedge = trisk[edge,0:trisk_deg[edge]] |
---|
| 398 | if iedge.min()<1 : |
---|
| 399 | print 'error' |
---|
| 400 | |
---|
| 401 | self.primal_num, self.edge_num, self.dual_num = primal_num, edge_num, dual_num |
---|
| 402 | def period(x) : return (x+2*math.pi)%(2*math.pi) |
---|
| 403 | lon_i, lon_v, lon_e = map(period, (lon_i,lon_v,lon_e)) |
---|
| 404 | self.lon_i, self.lat_i = lon_i, lat_i |
---|
| 405 | self.lon_v, self.lat_v = lon_v, lat_v |
---|
| 406 | self.lon_e, self.lat_e, self.angle_e = lon_e, lat_e, angle_e |
---|
| 407 | self.primal_deg, self.primal_vertex = primal_deg, primal_vertex |
---|
| 408 | self.primal = tri.Triangulation(lon_i*180./math.pi, lat_i*180./math.pi) |
---|
| 409 | self.dual_deg, self.dual_vertex = dual_deg, dual_vertex |
---|
| 410 | self.dual = tri.Triangulation(lon_v*180./math.pi, lat_v*180./math.pi) |
---|
| 411 | self.triedge = tri.Triangulation(lon_e*180./math.pi, lat_e*180./math.pi) |
---|
| 412 | self.dx = de.min() |
---|
| 413 | self.lon3D_i, self.ll3D = np.meshgrid(lon_i, range(llm)) |
---|
| 414 | self.lat3D_i, self.ll3D = np.meshgrid(lat_i, range(llm)) |
---|
| 415 | def field_theta(self): return squeeze((self.nqdyn,self.primal_num,self.llm)) |
---|
| 416 | def field_mass(self): return squeeze((self.primal_num,self.llm)) |
---|
| 417 | def field_z(self): return squeeze((self.dual_num,self.llm)) |
---|
| 418 | def field_w(self): return squeeze((self.primal_num,self.llm+1)) |
---|
| 419 | def field_u(self): return squeeze((self.edge_num,self.llm)) |
---|
| 420 | def field_ps(self): return squeeze((self.primal_num,)) |
---|
| 421 | def ucov2D(self, ulon, ulat): |
---|
| 422 | return self.de*(ulon*np.cos(self.angle_e)+ulat*np.sin(self.angle_e)) |
---|
| 423 | def ucov3D(self, ulon, ulat): |
---|
| 424 | ucov = np.zeros((self.edge_num,ulon.shape[1])) |
---|
| 425 | for edge in range(self.edge_num): |
---|
| 426 | angle=self.angle_e[edge] |
---|
| 427 | ucov[edge,:] = self.de[edge]*(ulon[edge,:]*math.cos(angle)+ulat[edge,:]*math.sin(angle)) |
---|
| 428 | return ucov |
---|
| 429 | def plot_i(self,data): |
---|
| 430 | plot(self.primal,data) |
---|
| 431 | def plot_v(self,data): |
---|
| 432 | plot(self.dual,data) |
---|
| 433 | def plot_e(self,data): |
---|
| 434 | plot(self.triedge,data) |
---|
| 435 | |
---|
| 436 | #-------------------------------------- Mesh partitioning ------------------------------------------# |
---|
| 437 | |
---|
| 438 | # Helper functions and interface to ParMETIS |
---|
| 439 | # list_stencil converts an adjacency graph from array format index[num_cells, MAX_EDGES] to compressed format |
---|
| 440 | # loc_stencil returns the start/end indices (vtxdist) expected by ParMETIS |
---|
| 441 | # i.e. index[start:end] with start=vtxdist[cell], end=vtxdist[cell+1] lists the edges of cell 'cell' |
---|
| 442 | |
---|
| 443 | def list_stencil(degree, stencil, cond=lambda x:True): |
---|
| 444 | for i in range(degree.size): |
---|
| 445 | for j in range(degree[i]): |
---|
| 446 | s=stencil[i,j] |
---|
| 447 | if cond(s): yield stencil[i,j] |
---|
| 448 | |
---|
| 449 | def loc_stencil(degree, stencil): |
---|
| 450 | loc=0 |
---|
| 451 | for i in range(degree.size): |
---|
| 452 | yield loc |
---|
| 453 | loc=loc+degree[i] |
---|
| 454 | yield loc |
---|
| 455 | |
---|
| 456 | def partition_mesh(degree, stencil, nparts): |
---|
| 457 | # arguments : PArray1D and PArray2D describing mesh, number of desired partitions |
---|
| 458 | dim_cell, degree, stencil = degree.dim, degree.data, stencil.data |
---|
| 459 | comm, vtxdist, idx_start, idx_end = dim_cell.comm, dim_cell.vtxdist, dim_cell.start, dim_cell.end |
---|
| 460 | mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size() |
---|
| 461 | adjncy_loc, xadj_loc = list_stencil(degree, stencil), loc_stencil(degree, stencil) |
---|
| 462 | adjncy_loc, xadj_loc = [np.asarray(list(x), dtype=np.int32) for x in (adjncy_loc, xadj_loc)] |
---|
| 463 | owner = np.zeros(idx_end-idx_start, dtype=np.int32); |
---|
[618] | 464 | ker.dynamico_partition_graph(mpi_rank, mpi_size, vtxdist, xadj_loc, adjncy_loc, nparts, owner) |
---|
[615] | 465 | return owner |
---|
| 466 | |
---|
| 467 | def partition_from_stencil(owner2, degree, stencil): |
---|
| 468 | # 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 |
---|
| 469 | dim1, dim2= degree.dim, owner2.dim |
---|
| 470 | degree, stencil, n = degree.data, stencil.data, dim1.end-dim1.start |
---|
| 471 | cells2 = list_stencil(degree, stencil) |
---|
| 472 | cells2 = sorted(list(set(list(cells2)))) # list of cells for which we need to know owner2 |
---|
| 473 | get2 = dim2.getter(cells2) |
---|
| 474 | owner1, owner2, glob2loc = np.zeros(n, dtype=np.int32), get2(owner2), get2.dict |
---|
| 475 | for i in range(n): |
---|
| 476 | owners = [ owner2[glob2loc[stencil[i,j]]] for j in range(degree[i]) ] |
---|
| 477 | if i%2 == 0 : owner1[i] = min(owners) |
---|
| 478 | else : owner1[i] = max(owners) |
---|
| 479 | return owner1 |
---|
| 480 | |
---|
| 481 | def find_my_cells(owner): # a PArray1D containing the data returned by partition_mesh |
---|
| 482 | dim, comm, owner = owner.dim, owner.dim.comm, owner.data |
---|
| 483 | mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size() |
---|
| 484 | cells=[set() for i in range(mpi_size)] |
---|
| 485 | for i in range(len(owner)): |
---|
| 486 | cells[owner[i]].add(dim.start+i) |
---|
| 487 | cells = [sorted(list(cells[i])) for i in range(mpi_size)] |
---|
| 488 | mycells = comm.alltoall(cells) |
---|
| 489 | mycells = sorted(sum(mycells, [])) # concatenate into a single list |
---|
| 490 | return mycells |
---|
| 491 | |
---|
| 492 | #---------------------------------- Stencil management ---------------------------------------# |
---|
| 493 | |
---|
| 494 | # Class Stencil represents an adjacency relationship (e.g. cell=>edges) |
---|
| 495 | # using adjacency information read from PArrays |
---|
| 496 | # It computes a list of "edges" adjacent to a given list of "cells" |
---|
| 497 | # This is used to form the sets E0 -> C0 -> E1 -> V1 -> E2 -> C1 |
---|
| 498 | # which are then used to form lists of global indices for V1,C1,E2 such that |
---|
| 499 | # C0, E0, E1 form contiguous subsets of C1, E2 starting from 0 |
---|
| 500 | |
---|
| 501 | def reindex(vertex_dict, degree, bounds): |
---|
| 502 | for i in range(degree.size): |
---|
| 503 | for j in range(degree[i]): |
---|
| 504 | bounds[i,j] = vertex_dict[bounds[i,j]] |
---|
| 505 | |
---|
| 506 | class Stencil_glob: |
---|
| 507 | def __init__(self, degree, neigh, weight=None): |
---|
| 508 | self.degree, self.neigh, self.weight = degree, neigh, weight |
---|
| 509 | def __call__(self, cells, neigh_dict=None): |
---|
| 510 | return Stencil(cells, self.degree, self.neigh, neigh_dict, self.weight) |
---|
| 511 | |
---|
| 512 | class Stencil: |
---|
| 513 | def __init__(self, cells, degree, neigh, neigh_dict, weight=None): |
---|
| 514 | get = degree.dim.getter(cells) |
---|
| 515 | mydegree, myneigh = [get(x) for x in (degree, neigh)] |
---|
| 516 | if not weight is None : myweight = get(weight) |
---|
| 517 | if neigh_dict is None : |
---|
| 518 | keep = lambda n : True |
---|
| 519 | else : # if neigh_dict is present, only neighbours in dict are retained |
---|
| 520 | keep = lambda n : n in neigh_dict |
---|
| 521 | neigh_set = list_stencil(mydegree, myneigh, keep) |
---|
| 522 | self.neigh_set = list(set(list(neigh_set) )) |
---|
| 523 | rej=0 |
---|
| 524 | for i in range(len(mydegree)): # keep only elements in neigh_dict, in-place |
---|
| 525 | k=0 |
---|
| 526 | for j in range(mydegree[i]): |
---|
| 527 | n=myneigh[i,j] |
---|
| 528 | if keep(n): |
---|
| 529 | myneigh[i,k]=n |
---|
| 530 | if not weight is None : myweight[i,k]=myweight[i,j] |
---|
| 531 | k=k+1 |
---|
| 532 | else: |
---|
| 533 | rej=rej+1 |
---|
| 534 | mydegree[i]=k |
---|
| 535 | if neigh_dict is None: |
---|
| 536 | neigh_dict = {j:i for i,j in enumerate(self.neigh_set)} |
---|
| 537 | myneigh_loc = reindex(neigh_dict, mydegree, myneigh) |
---|
| 538 | self.degree, self.neigh_glob, self.neigh_loc = mydegree, myneigh, myneigh_loc |
---|
| 539 | |
---|
| 540 | def progressive_iter(mylist, cell_lists): |
---|
| 541 | for thelist in cell_lists: |
---|
| 542 | mylist = mylist + list(set(thelist)-set(mylist)) |
---|
| 543 | yield mylist |
---|
| 544 | def progressive_list(*cell_lists): |
---|
| 545 | # cell_lists : a tuple of lists of global indices, with each list a subset of the next |
---|
| 546 | # returns : a list 'mylist' such that for each list 'thelist' in cell_lists, thelist = mylist[0:len(thelist)] |
---|
| 547 | # example : edge_list = progressive_list(E0,E1,E2) with E0,E1,E2 increasing lists of cell edges |
---|
| 548 | return list(progressive_iter([], cell_lists)) |
---|
| 549 | |
---|