source: codes/icosagcm/devel/Python/test/py/Baroclinic_3D_ullrich.py @ 785

Last change on this file since 785 was 785, checked in by jisesh, 6 years ago

devel/Python : added bilaplacian dissipation for theta + diagnostic of curl

File size: 12.7 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 xios
7from dynamico import precision as prec
8from dynamico.meshes import Cartesian_mesh as Mesh
9
10import math as math
11import numpy as np
12import time
13import argparse
14import logging
15from numpy import pi, log, exp, sin, cos
16
17# Baroclinic instability test based on Ullrich et al. 2015, QJRMS
18
19parser = argparse.ArgumentParser()
20
21parser.add_argument("--mpi_ni", type=int, default=1,
22                    help="number of x processors")
23parser.add_argument("--mpi_nj", type=int, default=1,
24                    help="number of y processors")
25parser.add_argument("--T", type=float, default=3600.,
26                    help="Length of time slice in seconds")
27parser.add_argument("--N", type=int, default=48,
28                    help="Number of time slices")
29parser.add_argument("--log",type=str,default='DEBUG',choices=('DEBUG', 'INFO', 'WARNING', 'ERROR', 'CRITICAL'),
30                    help="Set verbosity level of logging")
31parser.add_argument("--Davies_N1", type=int, default=3)
32parser.add_argument("--Davies_N2", type=int, default=3)
33
34
35args = parser.parse_args()
36
37
38def baroclinic_3D(Lx,nx,Ly,ny,llm,ztop=25000.):
39    Rd    = 287.0      # Gas constant for dryy air (j kg^-1 K^-1)
40    T0    = 288.0      # Reference temperature (K)
41    lap   = 0.005      # Lapse rate (K m^-1)
42    b     = 2.         # Non dimensional vertical width parameter
43    u0    = 35.        # Reference zonal wind speed (m s^-1)
44    a     = 6.371229e6 # Radius of the Earth (m)
45    ptop  = 2000.
46    y0    = Ly*0.5
47    Cpd   = 1004.5
48    p0    = 1e5
49    up    = 1. # amplitude (m/s)
50    xc,yc,lp = 0.,Ly/2.,600000.
51
52    omega  = 7.292e-5  # Angular velocity of the Earth (s^-1)   
53    phi0   = 90.*np.pi/180.0 # Reference latitude North pi/4 (deg)
54    f0     = 2*omega*np.sin(phi0) 
55    beta0  = 2*omega*np.cos(phi0)/a
56    fb     = 2*omega*np.sin(phi0) - y0*2*omega*np.cos(phi0)/a
57
58    def Phi_xy(y):
59        fc = y*y - (Ly*y/pi)*sin(2*pi*y/Ly)
60        fd = Ly*Ly/(2*pi*pi)*cos(2*pi*y/Ly)
61        return .5*u0*( fb*(y-y0-Ly/(2*pi)*sin(2*pi*y/Ly)) + .5*beta0*(fc-fd-(Ly*Ly/3.)- Ly*Ly/(2*pi*pi)) )
62
63    def Phi_xyeta(y,eta): return T0*g/lap*(1-eta**(Rd*lap/g)) + Phi_xy(y)*log(eta)*exp(-((log(eta)/b)**2))
64    def ulon(x,y,eta):
65        u = -u0*(sin(pi*y/Ly)**2)*log(eta)*(eta**(-log(eta)/(b*b))) 
66        u = u + up*exp(-(((x-xc)**2+(y-yc)**2)/lp**2))
67        return  u
68
69    def tmean(eta) : return T0*eta**(Rd*lap/g)
70    def T(y,eta) : return tmean(eta)+(Phi_xy(y)/Rd)*(((2/(b*b))*(log(eta))**2)-1)*exp(-((0.5*log(eta))**2)) 
71    def p(eta): return p0*eta     # eta  = p/p_s
72
73    def eta(alpha) : return (1-(lap*ztop*alpha/(T0)))**(g/(Rd*lap)) # roughly equispaced levels
74    def coriolis(x,y): return f0+beta0*y
75
76    filename = 'cart_%03d_%03d.nc'%(nx,ny)
77    logging.info('Reading Cartesian mesh ...')
78    meshfile = meshes.DYNAMICO_Format(filename)
79    pmesh = meshes.Unstructured_PMesh(comm,meshfile)
80    pmesh.partition_curvilinear(args.mpi_ni,args.mpi_nj)
81    nqdyn, radius = 1, None
82    mesh  = meshes.Local_Mesh(pmesh, llm, nqdyn, radius, coriolis)
83
84
85    alpha_k = (np.arange(llm)  +.5)/llm
86    alpha_l = (np.arange(llm+1)+ 0.)/llm
87    x_ik, alpha_ik = np.meshgrid(mesh.lon_i,       alpha_k, indexing='ij')
88    y_ik, alpha_ik = np.meshgrid(mesh.lat_i+.5*Ly, alpha_k, indexing='ij')
89    x_il, alpha_il = np.meshgrid(mesh.lon_i,       alpha_l, indexing='ij')
90    y_il, alpha_il = np.meshgrid(mesh.lat_i+.5*Ly, alpha_l, indexing='ij')
91    x_ek, alpha_ek = np.meshgrid(mesh.lon_e,       alpha_k, indexing='ij')
92    y_ek, alpha_ek = np.meshgrid(mesh.lat_e+.5*Ly, alpha_k, indexing='ij')
93
94    print('ztop(ptop) according to Eq. 7:', T0/lap*(1.-(ptop/p0)**(Rd*lap/g))) 
95    print(np.shape(alpha_k),np.shape(alpha_l))
96    thermo = dyn.Ideal_perfect(Cpd, Rd, p0, T0)
97
98    eta_il = eta(alpha_il)
99    eta_ik = eta(alpha_ik)
100    eta_ek = eta(alpha_ek)
101    print('min max eta_il', np.min(eta_il),np.max(eta_il))
102
103    Phi_il = Phi_xyeta(y_il, eta_il)
104    Phi_ik = Phi_xyeta(y_ik, eta_ik)
105    p_ik = p(eta_ik)
106    T_ik = T(y_ik, eta_ik)  #ik full level(40), il(41)
107
108    gas = thermo.set_pT(p_ik,T_ik)
109    mass_ik = mesh.field_mass()
110    for l in range(llm): 
111        mass_ik[:,l]=(Phi_il[:,l+1]-Phi_il[:,l])/(g*gas.v[:,l])
112    Sik, Wil  = gas.s*mass_ik, mesh.field_w()
113   
114    u_ek = mesh.ucov3D(ulon(x_ek, y_ek, eta_ek), 0.*eta_ek)
115
116    print('ztop (m) = ', Phi_il[0,-1]/g, ztop)
117    ptop = p(eta(1.))
118    print( 'ptop (Pa) = ', gas.p[0,-1], ptop)
119
120    params=dyn.Struct()
121    params.ptop=ptop
122    params.dx=dx
123    params.dx_g0=dx/g
124    params.g = g
125
126    # define parameters for lower BC
127    pbot = p(eta_il[:,0])
128    print( 'min p, T :', pbot.min(), tmean(pbot/p0))
129    gas_bot = thermo.set_pT(pbot, tmean(pbot/p0))
130    params.pbot = gas_bot.p
131    params.rho_bot = 1e6/gas_bot.v
132   
133    return thermo, mesh, params, prec.asnum([mass_ik,Sik,u_ek,Phi_il,Wil]), gas
134
135def diagnose(Phi,S,m,W):
136    s=S/m
137    for l in range(llm):
138        v[:,l]=(Phi[:,l+1]-Phi[:,l])/(g*m[:,l])
139        w[:,l]=.5*params.g*(W[:,l+1]+W[:,l])/m[:,l]
140        z[:,l]=.5*(Phi[:,l+1]+Phi[:,l])/params.g
141    gas = thermo.set_vs(v,s)
142    return gas, w, z
143
144def grad(mesh, s_ik):
145    llm=s_ik.shape[1]
146    edge_num, left, right = mesh.edge_num, mesh.left, mesh.right
147    grad_ek = prec.zeros((edge_num, llm))
148    for e in range(edge_num):
149        left_e = left[e]-1 # Fortran => Python indexing
150        right_e = right[e]-1 # Fortran => Python indexing
151        grad_ek[e,:] = s_ik[right_e,:]-s_ik[left_e,:]
152    return grad_ek
153
154def div(mesh, u_ek):
155    llm=u_ek.shape[1]
156    primal_num, primal_deg, primal_edge = mesh.primal_num, mesh.primal_deg, mesh.primal_edge
157    primal_ne, le_de, Ai = mesh.primal_ne, mesh.le_de, mesh.Ai
158    div_ik = prec.zeros((primal_num, llm))
159    for i in range(primal_num):
160        div_i = 0.
161        for iedge in range(primal_deg[i]):
162            e = primal_edge[i,iedge]-1 # Fortran => Python indexing
163            div_i = div_i + u_ek[e,:]*primal_ne[i,iedge]*le_de[e]
164        div_ik[i,:] = div_i / Ai[i]
165    return div_ik
166
167def curl(mesh, u_ek):
168    llm=u_ek.shape[1]
169    dual_num, dual_deg, dual_edge, dual_ne, Av = mesh.dual_num, mesh.dual_deg, mesh.dual_edge, mesh.dual_ne, mesh.Av
170    curl_vk = prec.zeros((dual_num, llm))
171    for v in range(dual_num):
172        curl_v = 0.
173        for iedge in range(dual_deg[v]):
174            e = dual_edge[v,iedge]-1 # Fortran => Python indexing
175            curl_v = curl_v + u_ek[e,:]*dual_ne[v,iedge]
176        curl_vk[v,:] = curl_v / Av[v]
177    return curl_vk
178
179class Davies:
180    def __init__(self,N1,N2,x_i,y_i,x_e,y_e):
181        self.N1, self.N2 = N1, N2
182        self.beta_i = self.mask(x_i,y_i)
183        self.beta_e = self.mask(x_e,y_e)
184        self.x_e,self.y_e = x_e,y_e
185    def mask0(self,c,c0,delta): # 1D building block for Davies relaxation
186        N1, N2 = self.N1, self.N2
187        m = np.zeros(c.size)
188        c1,c2 = c0-N1*delta, c0-(N1+N2)*delta
189        for i in range(c.size):
190            ci=c[i]
191            m[i] = (1.+np.cos((ci-c2)*np.pi/(N2*delta)))/2.0
192            if ci < c2 : m[i]=1.
193            if ci > c1 : m[i]=0.
194        return m
195    def relax(self, llm, step, flow):
196        beta_i, beta_e = self.beta_i, self.beta_e
197        m,S,u,Phi,W=flow
198        for l in range(llm):
199            step.mass[:,l]        = beta_i*step.mass[:,l]        + (1.-beta_i)*m[:,l]
200            step.theta_rhodz[:,l] = beta_i*step.theta_rhodz[:,l] + (1.-beta_i)*S[:,l]
201            step.u[:,l]           = beta_e*step.u[:,l]           + (1.-beta_e)*u[:,l]
202        for l in range(llm+1):
203            step.geopot[:,l]      = beta_i*step.geopot[:,l]      + (1.-beta_i)*Phi[:,l]
204            step.W[:,l]           = beta_i*step.W[:,l]           + (1.-beta_i)*W[:,l]
205
206class myDavies(Davies):
207    def mask(self,x,y):
208        logging.debug('davies dy = %f' % dy)
209        return self.mask0(y,.5*Ly,dy)*self.mask0(-y,.5*Ly,dy)
210
211
212       
213with xios.Client() as client: # setup XIOS which creates the DYNAMICO communicator
214    comm = client.comm
215    mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size()
216
217    #-------------Logging verbosity configuration---------------------------
218    myloglevel = args.log
219    myloglevel = getattr(logging, myloglevel.upper())
220    logging.basicConfig(filename='out.log',filemode='w',level=myloglevel, 
221        format='%(processName)s:%(name)s:%(filename)s:%(module)s:%(funcName)s:%(lineno)d:%(levelname)s:%(message)s' )
222
223    logging.debug('%d/%d starting'%(mpi_rank,mpi_size))
224
225    g, Lx, Ly = 9.81, 4e7, 6e6
226    nx, ny, llm = 200, 30, 30
227    dx,dy = Lx/nx,Ly/ny
228
229    unst.setvar('g',g)
230
231    thermo, mesh, params, flow0, gas0 = baroclinic_3D(Lx,nx,Ly,ny,llm)
232
233    mass_bl,mass_dak,mass_dbk = meshes.compute_hybrid_coefs(flow0[0])
234    print( 'Type of mass_bl, mass_dak, mass_dbk : ', [x.dtype for x in mass_bl, mass_dak, mass_dbk])
235    unst.ker.dynamico_init_hybrid(mass_bl,mass_dak,mass_dbk)
236
237    T  = args.T
238    dt = 100.#180.#360.
239    dz = flow0[3].max()/(params.g*llm)
240    nt = int(math.ceil(T/dt))
241    dt = T/nt
242    logging.info( 'Time step : %d x %g = %g s' % (nt,dt,nt*dt))
243   
244    caldyn_thermo, caldyn_eta = unst.thermo_entropy, unst.eta_mass
245
246    if False: # time stepping in Python
247        caldyn = unst.Caldyn_NH(caldyn_thermo,caldyn_eta, mesh,thermo,params,params.g)
248        scheme = time_step.ARK2(caldyn.bwd_fast_slow, dt)
249        def next_flow(m,S,u,Phi,W):
250            return scheme.advance((m,S,u,Phi,W),nt)
251
252    else: # time stepping in Fortran
253        scheme = time_step.ARK2(None, dt)
254        caldyn_step = unst.caldyn_step_NH(mesh,scheme,1, caldyn_thermo,caldyn_eta,thermo,params,params.g)
255        davies = myDavies(args.Davies_N1, args.Davies_N2, mesh.lon_i, mesh.lat_i, mesh.lon_e,mesh.lat_e)
256        print('mask min/max :', davies.beta_i.min(), davies.beta_i.max() )
257        logging.debug('mask min/max :')
258        logging.debug(davies.beta_i.min())
259        logging.debug(davies.beta_i.max())
260        def next_flow(m,S,u,Phi,W):
261            # junk,fast,slow = caldyn.bwd_fast_slow(flow, 0.)
262            caldyn_step.mass[:,:], caldyn_step.theta_rhodz[:,:], caldyn_step.u[:,:] = m,S,u
263            caldyn_step.geopot[:,:], caldyn_step.W[:,:] = Phi,W
264            kappa = dx**4/43200.
265            for i in range(nt):
266                caldyn_step.next()
267                davies.relax(llm, caldyn_step, flow0)
268                m,S = caldyn_step.mass, caldyn_step.theta_rhodz
269                s = S/m
270                laps = mesh.field_mass()
271                bilaps = mesh.field_mass()
272                unst.ker.dynamico_scalar_laplacian(s,laps)
273                unst.ker.dynamico_scalar_laplacian(laps,bilaps)
274#                grads = grad(mesh,s)
275#                laps = div(mesh,grads)
276#                gradlaps = grad(mesh,laps) # bilaplacian
277#                bilaps = div(mesh,gradlaps)               
278                caldyn_step.theta_rhodz[:] = S - dt*kappa*bilaps*m # Euler step
279            return (caldyn_step.mass.copy(), caldyn_step.theta_rhodz.copy(), caldyn_step.u.copy(),
280                    caldyn_step.geopot.copy(), caldyn_step.W.copy())
281
282    m,S,u,Phi,W = flow0
283    if caldyn_thermo == unst.thermo_theta:
284        s=S/m
285        theta = thermo.T0*np.exp(s/thermo.Cpd)
286        S=m*theta
287        title_format = 'Potential temperature at t=%g s (K)'
288    else:
289        title_format = 'Specific entropy at t=%g s (J/K/kg)'
290
291    w=mesh.field_mass()
292    z=mesh.field_mass()
293
294
295    #T, nslice, dt = 3600., 1, 3600.
296    Nslice = args.N
297
298    with xios.Context_Curvilinear(mesh,1, 24*3600) as context:
299        # now XIOS knows about the mesh and we can write to disk
300        v = mesh.field_mass() # specific volume (diagnosed)
301       
302        for i in range(Nslice):
303            context.update_calendar(i+1)
304
305            # Diagnose quantities of interest from prognostic variables m,S,u,Phi,W
306            gas, w, z = diagnose(Phi,S,m,W)
307            curl_vk = curl(mesh, u)
308            # write to disk
309            context.send_field_primal('ps', davies.beta_i)
310            context.send_field_primal('temp', gas.T)
311            context.send_field_primal('p', gas.p)
312            context.send_field_primal('theta', gas.s)
313            context.send_field_primal('uz', w)
314            context.send_field_primal('z', z)
315            context.send_field_primal('curl', curl_vk)
316
317            print( 'ptop, model top (m) :', unst.getvar('ptop'), Phi.max()/unst.getvar('g'))
318
319            time1, elapsed1 =time.time(), unst.getvar('elapsed')
320            m,S,u,Phi,W = next_flow(m,S,u,Phi,W)
321            time2, elapsed2 =time.time(), unst.getvar('elapsed')
322            factor = 1000./nt
323            print( 'ms per full time step : ', factor*(time2-time1), factor*(elapsed2-elapsed1))
324            factor = 1e9/(4*nt*nx*ny*llm)
325            print( 'nanosec per gridpoint per full time step : ', factor*(time2-time1), factor*(elapsed2-elapsed1))
326
327logging.info('************DONE************')
Note: See TracBrowser for help on using the repository browser.