source: codes/icosagcm/devel/Python/test/py/RSW2_MPAS_W02.py @ 758

Last change on this file since 758 was 749, checked in by dubos, 6 years ago

devel/unstructured : fix recently introduced bugs

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