[681] | 1 | print 'Starting' |
---|
| 2 | |
---|
| 3 | from mpi4py import MPI |
---|
| 4 | comm = MPI.COMM_WORLD |
---|
| 5 | mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size() |
---|
| 6 | print '%d/%d starting'%(mpi_rank,mpi_size) |
---|
| 7 | prefix='fig_RSW2_MPAS_W02/%02d'%mpi_rank |
---|
| 8 | |
---|
[666] | 9 | print 'Loading DYNAMICO modules ...' |
---|
| 10 | from dynamico import unstructured as unst |
---|
[681] | 11 | from dynamico.meshes import MPAS_Format, Unstructured_PMesh as PMesh, Local_Mesh as Mesh |
---|
[666] | 12 | from dynamico import time_step |
---|
| 13 | print '...Done' |
---|
| 14 | |
---|
[641] | 15 | print 'Loading modules ...' |
---|
| 16 | import math as math |
---|
| 17 | import matplotlib.pyplot as plt |
---|
| 18 | import numpy as np |
---|
| 19 | print '...Done' |
---|
| 20 | |
---|
[681] | 21 | grid, llm, nqdyn = 2562, 1,1 # 2562, 10242, 40962 |
---|
[641] | 22 | Omega, radius, g, gh0 = 2.*np.pi/86400., 6.4e6, 1., 2.94e4 |
---|
[681] | 23 | N, T, courant = 40, 10400., 1.2 # simulation length = N*T |
---|
[641] | 24 | |
---|
| 25 | print 'Omega, planetary PV', Omega, 2*Omega/gh0 |
---|
| 26 | |
---|
| 27 | def f(lon,lat): return 2*Omega*np.sin(lat) # Coriolis parameter |
---|
| 28 | print 'Reading MPAS mesh ...' |
---|
[680] | 29 | meshfile = MPAS_Format('grids/x1.%d.grid.nc'%grid) |
---|
[681] | 30 | pmesh = PMesh(comm,meshfile) |
---|
[760] | 31 | pmesh.partition_metis() |
---|
[681] | 32 | mesh = Mesh(pmesh, llm, nqdyn, radius, f) |
---|
[641] | 33 | print '...Done' |
---|
| 34 | lon, lat = mesh.lon_i, mesh.lat_i |
---|
| 35 | x,y,z = np.cos(lat)*np.cos(lon), np.cos(lat)*np.sin(lon), np.sin(lat) |
---|
| 36 | |
---|
| 37 | unst.setvar('g',g) |
---|
| 38 | |
---|
| 39 | c0 = math.sqrt(gh0) # phase speed of barotropic mode |
---|
| 40 | dx = mesh.de.min() |
---|
| 41 | dt = courant*dx/c0 |
---|
[681] | 42 | print dx, dt |
---|
[641] | 43 | nt = int(math.ceil(T/dt)) |
---|
| 44 | dt = T/nt |
---|
| 45 | print dx, dt, dt*c0/dx |
---|
| 46 | print T, nt |
---|
| 47 | |
---|
[749] | 48 | scheme = time_step.RK4(None, dt) |
---|
[641] | 49 | print dt, scheme.csjl, scheme.cfjl |
---|
[642] | 50 | step = unst.caldyn_step_TRSW(mesh,scheme,nt) |
---|
[641] | 51 | |
---|
| 52 | u0 = Omega*radius/12. # cf Williamson (1991), p.13 |
---|
| 53 | gh1 = radius*Omega*u0+.5*u0**2 |
---|
| 54 | print 'Williamson (1991) test 2, u0=', u0 |
---|
| 55 | ulon = u0*np.cos(mesh.lat_e) |
---|
| 56 | Phi0 = gh0 - gh1*(np.sin(mesh.lat_i)**2) |
---|
| 57 | zeta0 = (2*u0/radius+2*Omega)*np.sin(mesh.lat_v) |
---|
| 58 | Phi0v = gh0 - (radius*Omega*u0+.5*u0**2)*(np.sin(mesh.lat_v)**2) |
---|
| 59 | q0 = zeta0/Phi0v |
---|
| 60 | fu_perp = mesh.ucov2D(0.*ulon,2*Omega*np.sin(mesh.lat_e)*ulon) |
---|
| 61 | |
---|
| 62 | # initial condition |
---|
| 63 | step.mass[:]=Phi0 |
---|
| 64 | step.theta_rhodz[:]=Phi0 |
---|
| 65 | step.u[:]=mesh.ucov2D(ulon,0.*ulon) |
---|
| 66 | |
---|
| 67 | gh,u = step.mass.copy(), step.u.copy() |
---|
| 68 | |
---|
| 69 | print gh.shape, gh.min(), gh.max(), u.min(), u.max() |
---|
| 70 | flow=gh,u |
---|
| 71 | |
---|
[681] | 72 | unst.ker.dynamico_update_halo(mesh.com_edges.index, 1, u.shape[0], u) |
---|
| 73 | |
---|
[641] | 74 | for i in range(20): |
---|
| 75 | if True: |
---|
[642] | 76 | step.next() |
---|
[641] | 77 | gh,u = step.mass, step.u |
---|
| 78 | print i, gh.shape, gh.min(), gh.max(), u.min(), u.max() |
---|
| 79 | mesh.plot_i(gh-Phi0) ; plt.title('err(gh)'); |
---|
[681] | 80 | plt.savefig('%s_err_gh_%02d.png'%(prefix,i)) |
---|
[641] | 81 | plt.close() |
---|
| 82 | else: |
---|
| 83 | # advance by one time step using dynamico_ARK |
---|
| 84 | gh,u = step.mass.copy(), step.u.copy() |
---|
| 85 | print i, gh.shape, gh.min(), gh.max(), u.min(), u.max() |
---|
| 86 | step.next() |
---|
| 87 | gh_now,u_now = step.mass.copy(), step.u.copy() |
---|
| 88 | dmass,dhs,du_slow,du_fast = step.drhodz[:], step.dtheta_rhodz[:], step.du_slow[:], step.du_fast[:] |
---|
| 89 | du=du_slow+du_fast |
---|
| 90 | |
---|
| 91 | # do the same using caldyn, obtain fast/slow tendencies first |
---|
| 92 | gh_ref, uref=flow |
---|
| 93 | junk, fast, slow = caldyn.bwd_fast_slow(flow,0.) |
---|
| 94 | dmass_ref, duslow_ref = slow |
---|
| 95 | junk, dufast_ref = fast |
---|
| 96 | du_ref = duslow_ref+dufast_ref |
---|
| 97 | flow=scheme.next(flow) |
---|
| 98 | gh_nowref, u_nowref=flow |
---|
| 99 | # mesh.plot_i(gh) ; plt.title('gh'); |
---|
| 100 | # mesh.plot_i(Phi0) ; plt.title('gh'); |
---|
| 101 | fig = plt.figure(figsize=(6, 8)) |
---|
| 102 | f, ((ax1, ax2), (ax3, ax4), (ax5,ax6), (ax7,ax8)) = plt.subplots(4,2) |
---|
| 103 | # mesh.plot_i(dmass-dhs) ; plt.title('dt*d(gh)/dt'); |
---|
| 104 | ax1.scatter(duslow_ref,du_slow-duslow_ref) ; ax1.set_title('du_slow'); |
---|
| 105 | ax2.scatter(dufast_ref,du_fast-dufast_ref) ; ax2.set_title('du_fast'); |
---|
| 106 | ax3.scatter(du_ref,du-du_ref) ; ax3.set_title('du'); |
---|
| 107 | ax4.scatter(gh_ref,gh-gh_ref) ; ax4.set_title('gh'); |
---|
| 108 | ax5.scatter(uref, u-uref) ; ax5.set_title('u'); |
---|
| 109 | ax6.scatter(dmass_ref,dmass-dmass_ref) ; ax6.set_title('dmass'); |
---|
| 110 | ax7.scatter(gh_nowref,gh_now-gh_nowref) ; ax7.set_title('gh_now'); |
---|
| 111 | ax8.scatter(u_nowref,u_now-u_nowref) ; ax8.set_title('u_now'); |
---|
| 112 | f.subplots_adjust(hspace=1.) |
---|
[681] | 113 | plt.savefig('%s_scatter_%02d.png'%(prefix,i)) |
---|
[641] | 114 | plt.close() |
---|
| 115 | |
---|
| 116 | print 'Time spent in DYNAMICO (s) : ', unst.getvar('elapsed') |
---|