from dynamico import getargs getargs.add("--LAM", action='store_true') # Args for global mesh getargs.add("--grid", type=int, help='Number of hexagons', default=2562) # Args for LAM getargs.add("--nx", type=int, help='Zonal dimension of LAM mesh', default=100) getargs.add("--ny", type=int, help='Meridional dimension of LAM mesh', default=100) getargs.add("--dx", type=float, help='Resolution at center of LAM domain', default=1e5) getargs.add("--center_lat", type=float, help='Latitude in degrees of LAM center', default=0.) getargs.add("--Davies_N1", type=int, default=3) getargs.add("--Davies_N2", type=int, default=3) args = getargs.parse() log_master, log_world = getargs.getLogger() INFO, DEBUG, ERROR = log_master.info, log_master.debug, log_world.error INFO('Starting') from mpi4py import MPI comm = MPI.COMM_WORLD mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size() INFO('%d/%d starting'%(mpi_rank,mpi_size)) prefix='fig_RSW2_MPAS_W02/%02d'%mpi_rank INFO('Loading DYNAMICO modules ...') from dynamico.dev import unstructured as unst from dynamico.dev.meshes import MPAS_Format, Unstructured_PMesh as PMesh, Local_Mesh as Mesh from dynamico.dev import meshes from dynamico import time_step from dynamico import maps from dynamico.LAM import Davies print '...Done' print 'Loading modules ...' import math as math import matplotlib.pyplot as plt import numpy as np print '...Done' #--------------------------- functions and classes ----------------------------- class myDavies(Davies): def mask(self,X,Y): # X and Y are coordinates in the reference domain (cell = unit square) # numerical domain extends from -nx/2 ... nx/2 and -ny/2 ... ny/2 # useful domain extends from -X0 ... X0 and -Y0 ... Y0 N3 = args.Davies_N1+args.Davies_N2 X0 = args.nx/2. - N3 Y0 = args.ny/2. - N3 mask = self.mask0( X,X0,1.) # Western boundary mask *= self.mask0(-X,X0,1.) # Eastern boundary mask *= self.mask0( Y,Y0,1.) # Northern boundary mask *= self.mask0(-Y,Y0,1.) # Southern boundary return mask #----------------------------- main program -------------------------------- llm, nqdyn = 1,1 # 2562, 10242, 40962 Omega, radius, g, gh0 = 2.*np.pi/86400., 6.4e6, 1., 2.94e4 N, T, courant = 40, 10400., 1.2 # simulation length = N*T print 'Omega, planetary PV', Omega, 2*Omega/gh0 if args.LAM: nx, ny = args.nx, args.ny filename = 'cart_%03d_%03d.nc'%(nx,ny) INFO('Reading Cartesian mesh ...') meshfile = meshes.DYNAMICO_Format(filename) pmesh = meshes.Unstructured_PMesh(comm,meshfile) pmesh.partition_curvilinear(args.mpi_ni,args.mpi_nj) planet = maps.PolarStereoMap(radius,Omega, args.dx, args.center_lat*np.pi/180.) else: planet = maps.SphereMap(radius, Omega) print 'Reading MPAS mesh ...' meshfile = MPAS_Format('grids/x1.%d.grid.nc'%args.grid) pmesh = PMesh(comm,meshfile) pmesh.partition_metis() mesh = Mesh(pmesh, llm, nqdyn, planet) print '...Done' lon, lat = mesh.lon_i, mesh.lat_i x,y,z = np.cos(lat)*np.cos(lon), np.cos(lat)*np.sin(lon), np.sin(lat) unst.setvar('g',g) c0 = math.sqrt(gh0) # phase speed of barotropic mode dx = mesh.de.min() dt = courant*dx/c0 print dx, dt nt = int(math.ceil(T/dt)) dt = T/nt print dx, dt, dt*c0/dx print T, nt scheme = time_step.RK4(None, dt) print dt, scheme.csjl, scheme.cfjl step = unst.caldyn_step_TRSW(mesh,scheme,1) u0 = Omega*radius/12. # cf Williamson (1991), p.13 gh1 = radius*Omega*u0+.5*u0**2 print 'Williamson (1991) test 2, u0=', u0 ulon = u0*np.cos(mesh.lat_e) gh_ini = gh0 - gh1*(np.sin(mesh.lat_i)**2) u_ini = mesh.ucov2D(ulon,0.*ulon) # unst.ker.dynamico_update_halo(mesh.com_edges.index, 1, u.shape[0], u) if args.LAM: davies = myDavies(args.Davies_N1, args.Davies_N2, mesh.ref_lon_i, mesh.ref_lat_i, mesh.ref_lon_e,mesh.ref_lat_e) def relax(): davies.relax_RSW(llm, step, (gh_ini, u_ini) ) else: def relax(): pass def next_flow(Phi, u): step.mass[:], step.theta_rhodz[:], step.u[:] = Phi, Phi, u for i in range(nt): step.next() relax() return step.mass.copy(), step.u.copy() gh, u = gh_ini, u_ini print gh.shape, gh.min(), gh.max(), u.min(), u.max() for i in range(20): if True: gh, u = next_flow(gh,u) # step.next() # gh,u = step.mass, step.u print i, gh.shape, gh.min(), gh.max(), u.min(), u.max() mesh.plot_i(gh-gh_ini) ; plt.title('err(gh)'); plt.savefig('%s_err_gh_%02d.png'%(prefix,i)) plt.close() else: # advance by one time step using dynamico_ARK gh,u = step.mass.copy(), step.u.copy() print i, gh.shape, gh.min(), gh.max(), u.min(), u.max() step.next() gh_now,u_now = step.mass.copy(), step.u.copy() dmass,dhs,du_slow,du_fast = step.drhodz[:], step.dtheta_rhodz[:], step.du_slow[:], step.du_fast[:] du=du_slow+du_fast # do the same using caldyn, obtain fast/slow tendencies first gh_ref, uref=flow junk, fast, slow = caldyn.bwd_fast_slow(flow,0.) dmass_ref, duslow_ref = slow junk, dufast_ref = fast du_ref = duslow_ref+dufast_ref flow=scheme.next(flow) gh_nowref, u_nowref=flow # mesh.plot_i(gh) ; plt.title('gh'); # mesh.plot_i(Phi0) ; plt.title('gh'); fig = plt.figure(figsize=(6, 8)) f, ((ax1, ax2), (ax3, ax4), (ax5,ax6), (ax7,ax8)) = plt.subplots(4,2) # mesh.plot_i(dmass-dhs) ; plt.title('dt*d(gh)/dt'); ax1.scatter(duslow_ref,du_slow-duslow_ref) ; ax1.set_title('du_slow'); ax2.scatter(dufast_ref,du_fast-dufast_ref) ; ax2.set_title('du_fast'); ax3.scatter(du_ref,du-du_ref) ; ax3.set_title('du'); ax4.scatter(gh_ref,gh-gh_ref) ; ax4.set_title('gh'); ax5.scatter(uref, u-uref) ; ax5.set_title('u'); ax6.scatter(dmass_ref,dmass-dmass_ref) ; ax6.set_title('dmass'); ax7.scatter(gh_nowref,gh_now-gh_nowref) ; ax7.set_title('gh_now'); ax8.scatter(u_nowref,u_now-u_nowref) ; ax8.set_title('u_now'); f.subplots_adjust(hspace=1.) plt.savefig('%s_scatter_%02d.png'%(prefix,i)) plt.close() print 'Time spent in DYNAMICO (s) : ', unst.getvar('elapsed')