from dynamico.meshes import zeros import numpy as np import netCDF4 as cdf import argparse #----------------------- Cartesian mesh ----------------------- def concat(x,y): z = np.asarray([x,y]) return z.transpose() # arrays is a list of arrays # vals is a list of tuples # each tuple is stored in each array def put(ij, deg, arrays, vals): k=0 for vv in vals: # vv is a tuple of values to be stored in arrays for array,v in zip(arrays,vv): array[ij,k]=v k=k+1 deg[ij]=k class Cartesian_mesh: def __init__(self,nx,ny,llm,nqdyn,Lx,Ly): dx,dy = Lx/nx, Ly/ny self.dx, self.dy = dx,dy self.nx, self.ny, self.llm, self.nqdyn = nx,ny,llm,nqdyn self.field_z = self.field_mass # 1D coordinate arrays x=(np.arange(nx)-nx/2.)*Lx/nx y=(np.arange(ny)-ny/2.)*Ly/ny lev=np.arange(llm) levp1=np.arange(llm+1) self.x, self.y, self.lev, self.levp1 = x,y,lev,levp1 # 3D coordinate arrays self.xx,self.yy,self.ll = np.meshgrid(x,y,lev, indexing='ij') self.xxp1,self.yyp1,self.llp1 = np.meshgrid(x,y,levp1, indexing='ij') # beware conventions for indexing # Fortran order : llm,nx*ny,nqdyn / indices start at 1 # Python order : nqdyn,ny,nx,llm / indices start at 0 # indices below follow Fortran while x,y follow Python/C index=lambda x,y : ((x+(nx*(y+2*ny)))%(nx*ny))+1 indexu=lambda x,y : 2*index(x,y)-1 indexv=lambda x,y : 2*index(x,y) indices = lambda shape : np.zeros(shape,np.int32) primal_deg = indices(nx*ny) primal_edge = indices((nx*ny,4)) primal_ne = indices((nx*ny,4)) dual_deg = indices(nx*ny) dual_edge = indices((nx*ny,4)) dual_ne = indices((nx*ny,4)) primal_vertex = indices((nx*ny,4)) dual_vertex = indices((nx*ny,4)) Riv2 = np.zeros((nx*ny,4)) left = indices(2*nx*ny) right = indices(2*nx*ny) up = indices(2*nx*ny) down = indices(2*nx*ny) le_de = np.zeros(2*nx*ny) le = np.zeros(2*nx*ny) de = np.zeros(2*nx*ny) angle_e = np.zeros(2*nx*ny) trisk_deg = indices(2*nx*ny) trisk = indices((2*nx*ny,4)) # 4 TRiSK coefs per edge wee = np.zeros((2*nx*ny,4)) lon_i = indices(nx*ny) lon_v = indices(nx*ny) lon_e = indices(2*nx*ny) lat_i = indices(nx*ny) lat_v = indices(nx*ny) lat_e = indices(2*nx*ny) for x in range(nx): for y in range(ny): # NB : Fortran indices start at 1 # primal cells put(index(x,y)-1, primal_deg,(primal_edge,primal_vertex,primal_ne), ((indexu(x,y),index(x,y),1), (indexv(x,y),index(x-1,y),1), (indexu(x-1,y),index(x-1,y-1),-1), (indexv(x,y-1),index(x,y-1),-1) )) # dual cells put(index(x,y)-1,dual_deg,(dual_edge,dual_vertex,dual_ne,Riv2), ((indexv(x+1,y),index(x,y),1,.25), (indexu(x,y+1),index(x+1,y),-1,.25), (indexv(x,y),index(x+1,y+1),-1,.25), (indexu(x,y),index(x,y+1),1,.25) )) # edges : # left and right are adjacent primal cells # flux is positive when going from left to right # up and down are adjacent dual cells (vertices) # circulation is positive when going from down to up # u-points ij =indexu(x,y)-1 left[ij]=index(x,y) left[ij]=index(x,y) right[ij]=index(x+1,y) down[ij]=index(x,y-1) up[ij]=index(x,y) #le_de[ij]=dy/dx le[ij]=dy de[ij]=dx le_de[ij]=le[ij]/de[ij] angle_e[ij]=0. put(ij,trisk_deg,(trisk,wee),( (indexv(x,y),.25), (indexv(x+1,y),.25), (indexv(x,y-1),.25), (indexv(x+1,y-1),.25))) # v-points ij = indexv(x,y)-1 left[ij]=index(x,y) right[ij]=index(x,y+1) down[ij]=index(x,y) up[ij]=index(x-1,y) #le_de[ij]=dx/dy le[ij]=dx de[ij]=dy le_de[ij]=le[ij]/de[ij] angle_e[ij]=.5*np.pi put(ij,trisk_deg,(trisk,wee),( (indexu(x,y),-.25), (indexu(x-1,y),-.25), (indexu(x,y+1),-.25), (indexu(x-1,y+1),-.25))) ij = index(x,y)-1 lon_i[ij], lat_i[ij] = x, y ij = index(x,y)-1 lon_v[ij], lat_v[ij] = x+.5, y+.5 ij = indexu(x,y)-1 lon_e[ij], lat_e[ij] = x+.5, y ij = indexv(x,y)-1 lon_e[ij], lat_e[ij] = x, y+.5 Aiv=np.zeros(nx*ny)+dx*dy Ai=Av=np.zeros(nx*ny)+dx*dy self.llm=llm self.nqdyn=nqdyn self.nx=nx self.ny=ny self.primal_deg=primal_deg self.primal_edge=primal_edge self.primal_ne=primal_ne self.dual_deg=dual_deg self.dual_edge=dual_edge self.dual_ne=dual_ne self.dual_vertex=dual_vertex self.primal_vertex=primal_vertex self.left=left self.right=right self.down=down self.up=up self.trisk_deg=trisk_deg self.trisk=trisk self.Aiv=Aiv self.Ai=Ai self.Av=Av self.le=le self.de=de self.angle_e=angle_e self.Riv2=Riv2 self.wee=wee self.lon_i = lon_i self.lon_v = lon_v self.lon_e = lon_e #self.lon_ev = indices(2*nx*ny) self.lat_i = lat_i self.lat_v = lat_v self.lat_e = lat_e #self.lat_ev = indices(2*nx*ny) def field_theta(self,n=1): return zeros((n,self.nqdyn,self.ny,self.nx,self.llm)) def field_mass(self,n=1): return zeros((n,self.ny,self.nx,self.llm)) def field_z(self,n=1): return zeros((n,self.ny,self.nx,self.llm)) def field_w(self,n=1): return zeros((n,self.ny,self.nx,self.llm+1)) def field_u(self,n=1): if n==1: return np.zeros((self.ny,2*self.nx,self.llm)) else: return np.zeros((n,self.ny,2*self.nx,self.llm)) def field_ps(self,n=1): return zeros((n,self.ny,self.nx)) def ucomp(self,u): return u[:,range(0,2*self.nx,2),:] def set_ucomp(self,uv,u): uv[:,range(0,2*self.nx,2),:]=u def vcomp(self,u): return u[:,range(1,2*self.nx,2),:] #----------------------------WRITING MESH---------------------- def ncwrite(self, name): """The method writes Cartesian mesh on the disc. Args: Mesh parameters""" #----------------OPENING NETCDF OUTPUT FILE------------ try: f = cdf.Dataset(name, 'w', format='NETCDF4') except: print("CartesianMesh.ncwrite : Error occurred while opening new netCDF mesh file") raise #----------------DEFINING DIMENSIONS-------------------- for dimname, dimsize in [("primal_cell", self.primal_deg.size), ("dual_cell", self.dual_deg.size), ("edge", self.left.size), ("primal_edge_or_vertex", self.primal_edge.shape[1]), ("dual_edge_or_vertex", self.dual_edge.shape[1]), ("TWO", 2), ("trisk_edge", 4)]: f.createDimension(dimname,dimsize) #----------------DEFINING VARIABLES--------------------- f.description = "Cartesian_mesh" f.nx, f.ny = self.nx, self.ny def create_vars(dimname, info): for vname, vtype, vdata in info: # print('create_vars', dimname, vname, vdata.shape) var = f.createVariable(vname,vtype,dimname) var[:] = vdata create_vars("primal_cell", [("primal_deg","i4",self.primal_deg), ("Ai","f8",self.Aiv), ("lon_i","f8",self.lon_i), ("lat_i","f8",self.lat_i)] ) create_vars("dual_cell", [("dual_deg","i4",self.dual_deg), ("Av","f8",self.Aiv), ("lon_v","f8",self.lon_v), ("lat_v","f8",self.lat_v), ] ) create_vars( ("primal_cell","primal_edge_or_vertex"), [("primal_edge", "i4", self.primal_edge), ("primal_ne", "i4", self.primal_ne), ("primal_vertex", "i4", self.primal_vertex)] ) create_vars( ("dual_cell","dual_edge_or_vertex"), [("dual_edge", "i4", self.dual_edge), ("dual_vertex","i4",self.dual_vertex), ("dual_ne","i4",self.dual_ne), ("Riv2","f8",self.Riv2)] ) create_vars("edge", [("trisk_deg","i4",self.trisk_deg), ("le","f8",self.le), ("de","f8",self.de), ("lon_e","f8", self.lon_e), ("lat_e","f8", self.lat_e), ("angle_e","f8", self.angle_e)] ) create_vars( ("edge","TWO"), [("edge_left_right","i4", concat(self.left,self.right)), ("edge_down_up","i4", concat(self.down,self.up))] ) create_vars( ("edge","trisk_edge"), [("trisk","i4",self.trisk), ("wee","f8",self.wee)] ) f.close() #--------------------------------- Main program ------------------------------------- parser = argparse.ArgumentParser() parser.add_argument("-nx", type=int, default=128, choices=None, help="number of x points") parser.add_argument("-ny", type=int, default=128, choices=None, help="number of y points") parser.add_argument("-Lx", type=float, default=8., choices=None, help="Lx") parser.add_argument("-Ly", type=float, default=8., choices=None, help="Ly") parser.add_argument("-llm", type=int, default=1, choices=[1], help="number of vertical levels") args = parser.parse_args() nx, ny, Lx, Ly, llm, nqdyn = args.nx, args.ny,args.Lx, args.Ly, args.llm, 1 dx,dy=Lx/nx,Ly/ny mesh = Cartesian_mesh(nx,ny,llm,nqdyn,Lx,Ly) print('Successfully initialized Cartesian mesh') mesh.ncwrite('cart_'+str(nx)+'.nc') print('Successfully written Cartesian mesh to NetCDF File')