source: codes/icosagcm/devel/Python/test/py/NH_3D_bubble.py @ 798

Last change on this file since 798 was 753, checked in by dubos, 6 years ago

devel/unstructured : fixed indexing bug in Cartesian mesh

File size: 5.5 KB
Line 
1from dynamico import unstructured as unst
2from dynamico import dyn
3from dynamico import time_step
4from dynamico import DCMIP
5from dynamico import meshes
6from dynamico import precision as prec
7from dynamico.meshes import Cartesian_mesh as Mesh
8import math as math
9import matplotlib.pyplot as plt
10import numpy as np
11import time
12
13def thermal_bubble_3D(Lx,nx,Ly,ny,llm,ztop=1000., zc=350.,
14                      rc=250, thetac=0.5, x0=0., y0=0.):
15    Cpd, Rd, g, p0,theta0, T0 = 1004.5, 287.,9.81, 1e5, 300., 300.
16    nqdyn = 1
17   
18    Phi = lambda eta : g*ztop*eta
19    p=lambda Phi : p0*np.exp(-Phi/(Rd*T0))
20    zz = lambda p: -(Rd*T0*np.log(p/p0))/g
21    rr = lambda x,y,p: np.sqrt((x-x0)**2 + (y-y0)**2 + (zz(p)-zc)**2)
22    sa = lambda x,y,p: rr(x,y,p) < rc
23    deform = lambda x,y,p: (0.5*thetac*(1+np.cos(np.pi*rr(x,y,p)/rc)))*sa(x,y,p)
24    temp =  lambda p: theta0*(p/p0)**(Rd/Cpd)
25    T = lambda x,y,p: deform(x,y,p) + temp(p)
26
27    mesh = Mesh(nx,ny,llm,nqdyn,Lx,Ly,0.)
28    thermo = dyn.Ideal_perfect(Cpd, Rd, p0, T0)
29
30    Phi_il = Phi(mesh.llp1/float(llm))
31    Phi_ik = Phi((mesh.ll+.5)/llm)
32    p_ik = p(Phi_ik)
33    T_ik = T(mesh.xx, mesh.yy, p_ik)
34
35    gas = thermo.set_pT(p_ik,T_ik)
36    mass_ik = mesh.field_mass()
37    for l in range(llm): 
38        mass_ik[:,:,l]=(Phi_il[:,:,l+1]-Phi_il[:,:,l])/(g*gas.v[:,:,l])
39    Sik, ujk, Wil  = gas.s*mass_ik, mesh.field_u(), mesh.field_w()
40   
41    print 'ztop (m) = ', Phi_il[0,0,-1]/g, ztop
42    ptop = p(g*ztop)
43    print 'ptop (Pa) = ', gas.p[0,0,-1], ptop
44    params=dyn.Struct()
45    params.ptop=ptop
46    params.dx=dx
47    params.dx_g0=dx/g
48    params.g = g
49    pbot = p(Phi_il[:,:,0])
50    gas_bot = thermo.set_pT(pbot, temp(pbot))
51    params.pbot = gas_bot.p
52    params.rho_bot = 1e6/gas_bot.v
53   
54    return thermo, mesh, params, prec.asnum([mass_ik,Sik,ujk,Phi_il,Wil]), gas
55
56#Lx, nx, llm, thetac, T, Nslice, courant = 2000., 100, 50, 30., 5., 10, 2.8
57Lx, nx, llm, thetac, T, Nslice, courant = 2000., 20, 79, 30, 5., 10, 2.8
58#Lx, nx, llm, thetac, T, Nslice, courant = 3000., 75, 25, -30, 5., 10, 2.8
59
60nqdyn, dx = 1, Lx/nx
61Ly,ny,dy = Lx,nx,dx
62thermo, mesh, params, flow0, gas0 = thermal_bubble_3D(Lx,nx,Ly,ny,llm, thetac=thetac)
63
64# compute hybrid coefs from initial distribution of mass
65mass_bl,mass_dak,mass_dbk = meshes.compute_hybrid_coefs(flow0[0])
66print 'Type of mass_bl, mass_dak, mass_dbk : ', [x.dtype for x in mass_bl, mass_dak, mass_dbk]
67unst.ker.dynamico_init_hybrid(mass_bl,mass_dak,mass_dbk)
68
69dz = flow0[3].max()/(params.g*llm)
70#dt = courant*.5/np.sqrt(gas0.c2.max()*(dx**-2+dy**-2+dz**-2))
71dt = courant*.5/np.sqrt(gas0.c2.max()*(dx**-2+dy**-2))
72nt = int(math.ceil(T/dt))
73dt = T/nt
74print 'Time step : %d x %g s' % (nt,dt)
75
76#caldyn_thermo, caldyn_eta = unst.thermo_theta, unst.eta_mass
77caldyn_thermo, caldyn_eta = unst.thermo_entropy, unst.eta_mass
78#caldyn_thermo, caldyn_eta = unst.thermo_entropy, unst.eta_lag
79
80if False: # time stepping in Python
81    caldyn = unst.Caldyn_NH(caldyn_thermo,caldyn_eta, mesh,thermo,params,params.g)
82    scheme = time_step.ARK2(caldyn.bwd_fast_slow, dt)
83    def next_flow(m,S,u,Phi,W):
84        # junk,fast,slow = caldyn.bwd_fast_slow(flow, 0.)
85        return scheme.advance((m,S,u,Phi,W),nt)
86
87else: # time stepping in Fortran
88    scheme = time_step.ARK2(None, dt)
89    caldyn_step = unst.caldyn_step_NH(mesh,scheme,nt, caldyn_thermo,caldyn_eta,
90                                      thermo,params,params.g)
91    def next_flow(m,S,u,Phi,W):
92        # junk,fast,slow = caldyn.bwd_fast_slow(flow, 0.)
93        caldyn_step.mass[:,:,:], caldyn_step.theta_rhodz[:,:,:], caldyn_step.u[:,:,:] = m,S,u
94        caldyn_step.geopot[:,:,:], caldyn_step.W[:,:,:] = Phi,W
95        caldyn_step.next()
96        return (caldyn_step.mass.copy(), caldyn_step.theta_rhodz.copy(), caldyn_step.u.copy(),
97                caldyn_step.geopot.copy(), caldyn_step.W.copy())
98   
99m,S,u,Phi,W=flow0
100if caldyn_thermo == unst.thermo_theta:
101    s=S/m
102    theta = thermo.T0*np.exp(s/thermo.Cpd)
103    S=m*theta
104    title_format = 'Potential temperature at t=%g s (K)'
105else:
106    title_format = 'Specific entropy at t=%g s (J/K/kg)'
107
108w=mesh.field_mass()
109z=mesh.field_mass()
110
111for it in range(Nslice):
112    s=S/m ; s=.5*(s+abs(s))
113    for l in range(llm):
114        w[:,:,l]=.5*params.g*(W[:,:,l+1]+W[:,:,l])/m[:,:,l]
115        z[:,:,l]=.5*(Phi[:,:,l+1]+Phi[:,:,l])/params.g
116       
117    print 'ptop, model top (m) :', unst.getvar('ptop'), Phi.max()/unst.getvar('g')
118    jj=ny/2
119    xx,zz,ss,ww = mesh.xx[jj,:,:]/1000., z[jj,:,:]/1000., s[jj,:,:], w[jj,:,:]
120
121    f, (ax1, ax2) = plt.subplots(1, 2, sharey=True, figsize=(12,4))
122   
123    c=ax1.contourf(xx,zz,ss,20)
124    ax1.set_xlim((-.5,.5)), ax1.set_xlabel('x (km)')
125    ax1.set_ylim((0.,1.)), ax1.set_ylabel('z (km)')       
126    plt.colorbar(c,ax=ax1)
127    ax1.set_title(title_format % (it*T,))
128
129    #    plt.show()
130   
131    #    plt.figure(figsize=(12,5))
132    c=ax2.contourf(xx,zz,ww,20)
133    ax2.set_xlim((-.5,.5)), ax2.set_xlabel('x (km)')
134    ax2.set_ylim((0.,1.)), ax2.set_ylabel('z (km)')       
135    plt.colorbar(c,ax=ax2)
136    ax2.set_title('Vertical velocity at t=%g s (m/s)' % (it*T,))
137    #    plt.tight_layout()
138    #    plt.show()
139    plt.savefig('fig_NH_3D_bubble/%02d.png'%it)
140   
141    time1, elapsed1 =time.time(), unst.getvar('elapsed')
142    m,S,u,Phi,W = next_flow(m,S,u,Phi,W)
143    time2, elapsed2 =time.time(), unst.getvar('elapsed')
144    factor = 1000./nt
145    print 'ms per full time step : ', factor*(time2-time1), factor*(elapsed2-elapsed1)
146    factor = 1e9/(4*nt*nx*ny*llm)
147    print 'nanosec per gridpoint per full time step : ', factor*(time2-time1), factor*(elapsed2-elapsed1)
148       
Note: See TracBrowser for help on using the repository browser.