[632] | 1 | from dynamico import unstructured as unst |
---|
| 2 | from dynamico import dyn |
---|
| 3 | from dynamico import time_step |
---|
| 4 | from dynamico import DCMIP |
---|
| 5 | from dynamico import meshes |
---|
| 6 | from dynamico.meshes import MPAS_Mesh as Mesh |
---|
| 7 | import math as math |
---|
| 8 | import matplotlib.pyplot as plt |
---|
| 9 | import numpy as np |
---|
| 10 | import time |
---|
| 11 | |
---|
| 12 | #------------------------ initial condition ------------------------- |
---|
| 13 | |
---|
| 14 | Cpd, Rd, g = 1004.5, 287., 9.81 |
---|
| 15 | N, Teq, peq = 0.01, 300., 1e5 # background |
---|
| 16 | lonc, d2 = 2.*math.pi/3., 25e6 # perturbation |
---|
| 17 | N2, g2 = N*N, g*g |
---|
| 18 | N2_g2, g2_N2 = N2/g2, g2/N2 |
---|
| 19 | G, kappa = g2/(N2*Cpd), Rd/Cpd |
---|
| 20 | |
---|
| 21 | def psTs(lat) : |
---|
| 22 | cos2latm1 = np.cos(2*lat)-1. |
---|
| 23 | Ts = G+(Teq-G)*np.exp(-.25*N2_g2*u0*u0*cos2latm1) |
---|
| 24 | ps = peq*np.exp(u0/(4.*G*Rd)*u0*cos2latm1)*((Ts/Teq)**(1./kappa)) |
---|
| 25 | return ps, Ts |
---|
| 26 | |
---|
| 27 | def Phi_top(lat) : |
---|
| 28 | ps, Ts = psTs(lat) |
---|
| 29 | return -g2_N2*np.log( 1+(Ts/G)*( ((ptop/ps)**kappa) - 1 ) ) |
---|
| 30 | |
---|
| 31 | def T0_p0(lon,lat,Phi): |
---|
| 32 | ps, Ts = psTs(lat) |
---|
| 33 | pb = ps*( (G/Ts)*(np.exp(-N2_g2*Phi)-1.)+1. )**(1./kappa) |
---|
| 34 | Tb = G + (Ts-G)*np.exp(N2_g2*Phi) |
---|
| 35 | r = radius*np.arccos(np.cos(lat)*np.cos(lon-lonc)) |
---|
| 36 | Theta = dT*np.sin(math.pi*Phi/(g*ztop))*d2/(d2+r*r) |
---|
| 37 | T = Tb + Theta*(pb/peq)**kappa |
---|
| 38 | return T, pb |
---|
| 39 | |
---|
| 40 | def DCMIP31_3D(grid): |
---|
| 41 | def Phi(eta, lat): return eta*Phi_top(lat) |
---|
| 42 | def ulon(lat): return u0*np.cos(lat) |
---|
| 43 | def ulat(lat): return 0.*lat |
---|
| 44 | def f(lon,lat): return 0.*np.sin(lat) # Coriolis parameter |
---|
| 45 | |
---|
| 46 | |
---|
| 47 | nqdyn, preff, Treff = 1, 1e5, 300. |
---|
| 48 | thermo = dyn.Ideal_perfect(Cpd, Rd, preff, Treff) |
---|
| 49 | mesh = Mesh('grids/x1.%d.grid.nc'%grid, llm, nqdyn, radius, f) |
---|
| 50 | lev,levp1 = np.arange(llm),np.arange(llm+1) |
---|
| 51 | lat_ik,k_i = np.meshgrid(mesh.lat_i,lev, indexing='ij') |
---|
| 52 | lon_ik,k_i = np.meshgrid(mesh.lon_i,lev, indexing='ij') |
---|
| 53 | lat_il,l_i = np.meshgrid(mesh.lat_i,levp1, indexing='ij') |
---|
| 54 | lon_il,l_i = np.meshgrid(mesh.lon_i,levp1, indexing='ij') |
---|
| 55 | lat_ek,k_e = np.meshgrid(mesh.lat_e,lev, indexing='ij') |
---|
| 56 | |
---|
| 57 | Phi_ik, Phi_il = Phi(k_i/float(llm),lat_ik), Phi(l_i/float(llm),lat_il) |
---|
| 58 | Tb_ik, pb_ik = T0_p0(lon_ik, lat_ik, Phi_ik) |
---|
| 59 | Tb_il, pb_il = T0_p0(lon_il, lat_il, Phi_il) |
---|
| 60 | gas = thermo.set_pT(pb_ik,Tb_ik) |
---|
| 61 | mass_ik = mesh.field_mass() |
---|
| 62 | for l in range(llm): |
---|
| 63 | mass_ik[:,l]=(pb_il[:,l]-pb_il[:,l+1])/g |
---|
| 64 | Sik = gas.s*mass_ik |
---|
| 65 | |
---|
| 66 | # start at hydrostatic geopotential |
---|
| 67 | pb_ik = .5*(pb_il[:,1:]+pb_il[:,0:-1]) # pressure at full layers |
---|
| 68 | gas = thermo.set_ps(pb_ik,gas.s) |
---|
| 69 | for l in range(llm): |
---|
| 70 | Phi_il[:,l+1] = Phi_il[:,l] + g*mass_ik[:,l]*gas.v[:,l] |
---|
| 71 | |
---|
| 72 | ujk, Wil = mesh.ucov3D(ulon(lat_ek),ulat(lat_ek)), mesh.field_w() |
---|
| 73 | |
---|
| 74 | print 'ztop (m) = ', Phi_il[0,-1]/g, ztop |
---|
| 75 | print 'ptop (Pa) = ', gas.p[0,-1], ptop |
---|
| 76 | dx=mesh.de.min() |
---|
| 77 | params=dyn.Struct() |
---|
| 78 | params.g, params.ztop, params.ptop = g, ztop, ptop |
---|
| 79 | params.dx, params.dx_g0 = dx, dx/g |
---|
| 80 | params.pbot, params.rho_bot = 1e5+0.*mesh.lat_i, 1e6+0.*mesh.lat_i |
---|
| 81 | return thermo, mesh, params, (mass_ik,Sik,ujk,Phi_il,Wil), gas |
---|
| 82 | |
---|
| 83 | #------------------------ main program ------------------------- |
---|
| 84 | |
---|
| 85 | ztop, u0, dT, Xfactor = 1e4, 20., 1., 125. |
---|
| 86 | radius, grid, llm = 6.37122e6/Xfactor, 10242, 20 |
---|
| 87 | T, Nslice, courant = 360., 10, 3.0 |
---|
| 88 | |
---|
| 89 | junk, ptop = T0_p0(0.,0.,g*ztop) |
---|
| 90 | |
---|
| 91 | thermo, mesh, params, flow0, gas0 = DCMIP31_3D(grid) |
---|
| 92 | llm, dx = mesh.llm, params.dx |
---|
| 93 | dz = params.ztop/llm |
---|
| 94 | print 'llm, nb_hex, dx, dz =', llm, mesh.Ai.size, dx, dz |
---|
| 95 | |
---|
| 96 | caldyn_thermo, caldyn_eta = unst.thermo_entropy, unst.eta_lag |
---|
| 97 | caldyn = unst.Caldyn_NH(caldyn_thermo,caldyn_eta, mesh,params,thermo,params,params.g) |
---|
| 98 | |
---|
| 99 | mesh.plot_e(mesh.le/mesh.de) ; plt.title('le/de') |
---|
| 100 | plt.savefig('fig_NH_3D_DCMIP31/le_de.png') ;plt.close() |
---|
| 101 | |
---|
| 102 | mesh.plot_i(mesh.Ai) ; plt.title('Ai') |
---|
| 103 | plt.savefig('fig_NH_3D_DCMIP31/Ai.png'); plt.close() |
---|
| 104 | |
---|
| 105 | #dt = courant*.5/np.sqrt(gas0.c2.max()*(dx**-2+dz**-2)) |
---|
| 106 | dt = courant*.5*dx/np.sqrt(gas0.c2.max()) |
---|
| 107 | |
---|
| 108 | nt = int(math.ceil(T/dt)) |
---|
| 109 | dt = T/nt |
---|
| 110 | print 'Time step : %g s' % dt |
---|
| 111 | |
---|
| 112 | #scheme = time_step.RK4(caldyn.bwd_fast_slow, dt) |
---|
| 113 | scheme = time_step.ARK2(caldyn.bwd_fast_slow, dt, a32=0.7) |
---|
| 114 | |
---|
| 115 | flow=flow0 |
---|
| 116 | w=mesh.field_mass() |
---|
| 117 | z=mesh.field_mass() |
---|
| 118 | |
---|
| 119 | def plots(it): |
---|
| 120 | junk,fast,slow = caldyn.bwd_fast_slow(flow, 0.) |
---|
| 121 | m,S,u,Phi,W = flow |
---|
| 122 | s=S/m ; s=.5*(s+abs(s)) |
---|
| 123 | for l in range(llm): |
---|
| 124 | w[:,l]=.5*params.g*(W[:,l+1]+W[:,l])/m[:,l] |
---|
| 125 | z[:,l]=.5*(Phi[:,l+1]+Phi[:,l])/params.g |
---|
| 126 | |
---|
| 127 | print( 'ptop, model top (m), max(w), max(W) (m/s) :', |
---|
| 128 | unst.getvar('ptop'), Phi.max()/unst.getvar('g'), w.max(), W.max() ) |
---|
| 129 | |
---|
| 130 | mesh.plot_i(w[:,llm/2]) |
---|
| 131 | plt.title('vertical velocity at t=%d'%(it*T)) |
---|
| 132 | plt.savefig('fig_NH_3D_DCMIP31/w%02d.png'%it) |
---|
| 133 | plt.close() |
---|
| 134 | |
---|
| 135 | plots(0) |
---|
| 136 | |
---|
| 137 | for it in range(Nslice): |
---|
| 138 | unst.setvar('debug_hevi_solver',True) |
---|
| 139 | time1, elapsed1 =time.time(), unst.getvar('elapsed') |
---|
| 140 | flow = scheme.advance(flow,nt) |
---|
| 141 | time2, elapsed2 =time.time(), unst.getvar('elapsed') |
---|
| 142 | factor = 1000./(3*nt) |
---|
| 143 | print 'ms per call to caldyn_hevi : ', factor*(time2-time1), factor*(elapsed2-elapsed1) |
---|
| 144 | factor = 1e9/(4*nt*w.size) |
---|
| 145 | print 'nanosec per gridpoint per call to caldyn_hevi : ', factor*(time2-time1), factor*(elapsed2-elapsed1) |
---|
| 146 | plots(it+1) |
---|
| 147 | |
---|