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

Last change on this file since 800 was 800, checked in by jisesh, 5 years ago

devel/unstructured : introduced mapping from reference domain to physical domain

File size: 11.4 KB
RevLine 
[790]1from __future__ import print_function
[795]2from __future__ import division
3
[790]4from dynamico import getargs
5
6getargs.add("--T", type=float, default=3600.,
7                    help="Length of time slice in seconds")
8getargs.add("--perturb", type=float, default=1.,
9                    help="Amplitude of wind perturbation in m/s")
10getargs.add("--N", type=int, default=48,
11                    help="Number of time slices")
12getargs.add("--Davies_N1", type=int, default=3)
13getargs.add("--Davies_N2", type=int, default=3)
14getargs.add("--nx", type=int, default=200)
15getargs.add("--ny", type=int, default=30)
16getargs.add("--betaplane", action='store_true')
[795]17
18getargs.add("--kappa_divgrad",  type=float, default=3.0e15)
19getargs.add("--kappa_curlcurl", type=float, default=3.0e15)
[797]20getargs.add("--ztop", type=float, default=30000.)
[795]21
[790]22getargs.defaults(dt=360., llm=50)
23
[795]24getargs.config_log(filename='out.log',filemode='w') # must be done before calling Logger()
25#    getargs.config_log(filename='out.log',filemode='w',
26#                    format='%(processName)s:%(name)s:%(filename)s:%(module)s:%(funcName)s:%(lineno)d:%(levelname)s:%(message)s' )
27
[790]28logging = getargs.Logger('main')
29args = getargs.parse()
30
[761]31from dynamico import unstructured as unst
32from dynamico import dyn
33from dynamico import time_step
34from dynamico import DCMIP
35from dynamico import meshes
[800]36from dynamico import maps
[761]37from dynamico import xios
38from dynamico import precision as prec
39from dynamico.meshes import Cartesian_mesh as Mesh
[791]40from dynamico.kernels import grad, curl, div, KE
41from dynamico.LAM import Davies
[761]42
43import math as math
44import numpy as np
45import time
46from numpy import pi, log, exp, sin, cos
47
48# Baroclinic instability test based on Ullrich et al. 2015, QJRMS
49
[794]50def create_pmesh(nx,ny):
51    filename = 'cart_%03d_%03d.nc'%(nx,ny)
52    logging.info('Reading Cartesian mesh ...')
53    meshfile = meshes.DYNAMICO_Format(filename)
54    pmesh = meshes.Unstructured_PMesh(comm,meshfile)
55    pmesh.partition_curvilinear(args.mpi_ni,args.mpi_nj)
56    return pmesh
[764]57
[800]58def baroclinic_3D(pmesh, dx,Lx,Ly,llm,ztop):
[761]59    Rd    = 287.0      # Gas constant for dryy air (j kg^-1 K^-1)
60    T0    = 288.0      # Reference temperature (K)
61    lap   = 0.005      # Lapse rate (K m^-1)
62    b     = 2.         # Non dimensional vertical width parameter
63    u0    = 35.        # Reference zonal wind speed (m s^-1)
64    a     = 6.371229e6 # Radius of the Earth (m)
[787]65    y0    = .5*Ly
[761]66    Cpd   = 1004.5
67    p0    = 1e5
[790]68    up    = args.perturb # amplitude (m/s)
[774]69    xc,yc,lp = 0.,Ly/2.,600000.
[761]70
71    omega  = 7.292e-5  # Angular velocity of the Earth (s^-1)   
[790]72    phi0   = 45.*pi/180.0 # Reference latitude North pi/4 (deg)
73    f0     = 2*omega*sin(phi0) 
74    beta0  = 2*omega*cos(phi0)/a if args.betaplane else 0.
[787]75    fb     = f0 - beta0*y0
[761]76
[769]77    def Phi_xy(y):
[761]78        fc = y*y - (Ly*y/pi)*sin(2*pi*y/Ly)
79        fd = Ly*Ly/(2*pi*pi)*cos(2*pi*y/Ly)
80        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)) )
81
[769]82    def Phi_xyeta(y,eta): return T0*g/lap*(1-eta**(Rd*lap/g)) + Phi_xy(y)*log(eta)*exp(-((log(eta)/b)**2))
[774]83    def ulon(x,y,eta):
84        u = -u0*(sin(pi*y/Ly)**2)*log(eta)*(eta**(-log(eta)/(b*b))) 
[790]85        u = u + up*exp(-(((x-xc)**2+(y-yc)**2)/lp**2))
[774]86        return  u
87
[761]88    def tmean(eta) : return T0*eta**(Rd*lap/g)
[769]89    def T(y,eta) : return tmean(eta)+(Phi_xy(y)/Rd)*(((2/(b*b))*(log(eta))**2)-1)*exp(-((0.5*log(eta))**2)) 
[761]90    def p(eta): return p0*eta     # eta  = p/p_s
91
[787]92    def eta(alpha) : return (1.-(lap*ztop*alpha/T0))**(g/(Rd*lap)) # roughly equispaced levels
[800]93   
94    betaplane = maps.BetaPlaneMap(dx, .5*Lx, .5*Ly, f0, beta0, .5*Ly)
95    nqdyn = 1
96    mesh = meshes.Local_Mesh(pmesh, llm, nqdyn, betaplane)
[761]97
98    alpha_k = (np.arange(llm)  +.5)/llm
99    alpha_l = (np.arange(llm+1)+ 0.)/llm
[800]100    x_ik, alpha_ik = np.meshgrid(mesh.lon_i, alpha_k, indexing='ij')
101    y_ik, alpha_ik = np.meshgrid(mesh.lat_i, alpha_k, indexing='ij')
102    x_il, alpha_il = np.meshgrid(mesh.lon_i, alpha_l, indexing='ij')
103    y_il, alpha_il = np.meshgrid(mesh.lat_i, alpha_l, indexing='ij')
104    x_ek, alpha_ek = np.meshgrid(mesh.lon_e, alpha_k, indexing='ij')
105    y_ek, alpha_ek = np.meshgrid(mesh.lat_e, alpha_k, indexing='ij')
[761]106
107    print(np.shape(alpha_k),np.shape(alpha_l))
108    thermo = dyn.Ideal_perfect(Cpd, Rd, p0, T0)
109
110    eta_il = eta(alpha_il)
111    eta_ik = eta(alpha_ik)
112    eta_ek = eta(alpha_ek)
[790]113#    print('min max eta_il', np.min(eta_il),np.max(eta_il))
[761]114
[769]115    Phi_il = Phi_xyeta(y_il, eta_il)
116    Phi_ik = Phi_xyeta(y_ik, eta_ik)
[787]117    p_ik, p_il = p(eta_ik), p(eta_il)
[769]118    T_ik = T(y_ik, eta_ik)  #ik full level(40), il(41)
[761]119
120    gas = thermo.set_pT(p_ik,T_ik)
121    mass_ik = mesh.field_mass()
122    for l in range(llm): 
123        mass_ik[:,l]=(Phi_il[:,l+1]-Phi_il[:,l])/(g*gas.v[:,l])
[787]124#        mass_ik[:,l]=(p_il[:,l]-p_il[:,l+1])/g
[774]125    Sik, Wil  = gas.s*mass_ik, mesh.field_w()
[761]126   
127    u_ek = mesh.ucov3D(ulon(x_ek, y_ek, eta_ek), 0.*eta_ek)
128
[774]129    print('ztop (m) = ', Phi_il[0,-1]/g, ztop)
[761]130    ptop = p(eta(1.))
[774]131    print( 'ptop (Pa) = ', gas.p[0,-1], ptop)
[797]132    print('ztop(ptop) according to Eq. 7:', T0/lap*(1.-(ptop/p0)**(Rd*lap/g))) 
[761]133
134    params=dyn.Struct()
135    params.ptop=ptop
136    params.dx=dx
137    params.dx_g0=dx/g
138    params.g = g
139
140    # define parameters for lower BC
141    pbot = p(eta_il[:,0])
[774]142    print( 'min p, T :', pbot.min(), tmean(pbot/p0))
[761]143    gas_bot = thermo.set_pT(pbot, tmean(pbot/p0))
144    params.pbot = gas_bot.p
145    params.rho_bot = 1e6/gas_bot.v
146   
[774]147    return thermo, mesh, params, prec.asnum([mass_ik,Sik,u_ek,Phi_il,Wil]), gas
[761]148
[795]149def diagnose(Phi,S,m,W,u):
[780]150    s=S/m
[795]151    curl_vk = curl(mesh, u)
152    abs_vort_vk = mesh.field_z() # absolute vorticity
153    un = mesh.field_u()
154    v = mesh.field_mass() # specific volume
155    w = mesh.field_mass()
156    z = mesh.field_mass()
[769]157    for l in range(llm):
158        v[:,l]=(Phi[:,l+1]-Phi[:,l])/(g*m[:,l])
159        w[:,l]=.5*params.g*(W[:,l+1]+W[:,l])/m[:,l]
160        z[:,l]=.5*(Phi[:,l+1]+Phi[:,l])/params.g
[795]161        un[:,l]=u[:,l]/mesh.de
162        abs_vort_vk[:,l]=curl_vk[:,l] + mesh.fv
[769]163    gas = thermo.set_vs(v,s)
[795]164    ps = gas.p[:,0]+ .5*g*m[:,0]
165    return gas, w, z, ps, un, curl_vk, abs_vort_vk
[769]166
[771]167class myDavies(Davies):
168    def mask(self,x,y):
[800]169        logging.debug('davies dy = %f'% dx)
170        return self.mask0(y,Ly,dx)*self.mask0(-y,0.,dx)
[771]171       
[763]172with xios.Client() as client: # setup XIOS which creates the DYNAMICO communicator
173    comm = client.comm
174    mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size()
[761]175
[790]176    logging.info('%d/%d starting'%(mpi_rank,mpi_size))
[774]177
[795]178    g, Lx, Ly = 9.80616, 4e7, 6e6
[787]179    nx, ny, llm = args.nx, args.ny, args.llm
[794]180    dx = Lx/nx
181   
[763]182    unst.setvar('g',g)
[794]183   
184    pmesh = create_pmesh(nx,ny)
[800]185    thermo, mesh, params, flow0, gas0 = baroclinic_3D(pmesh, dx,Lx,Ly,llm, args.ztop)
[761]186
[769]187    mass_bl,mass_dak,mass_dbk = meshes.compute_hybrid_coefs(flow0[0])
188    unst.ker.dynamico_init_hybrid(mass_bl,mass_dak,mass_dbk)
[761]189
[774]190    T  = args.T
[787]191    dt = args.dt
[769]192    dz = flow0[3].max()/(params.g*llm)
193    nt = int(math.ceil(T/dt))
194    dt = T/nt
[774]195    logging.info( 'Time step : %d x %g = %g s' % (nt,dt,nt*dt))
[794]196
[769]197   
198    caldyn_thermo, caldyn_eta = unst.thermo_entropy, unst.eta_mass
199
200    if False: # time stepping in Python
201        caldyn = unst.Caldyn_NH(caldyn_thermo,caldyn_eta, mesh,thermo,params,params.g)
202        scheme = time_step.ARK2(caldyn.bwd_fast_slow, dt)
203        def next_flow(m,S,u,Phi,W):
204            return scheme.advance((m,S,u,Phi,W),nt)
205
206    else: # time stepping in Fortran
207        scheme = time_step.ARK2(None, dt)
[787]208        if args.hydrostatic:
[795]209            caldyn_step = unst.caldyn_step_HPE(mesh,scheme,1, caldyn_thermo,caldyn_eta,
210                                               thermo,params,params.g)
[787]211        else:
[795]212            caldyn_step = unst.caldyn_step_NH(mesh,scheme,1, caldyn_thermo,caldyn_eta,
213                                              thermo,params,params.g)
214        davies = myDavies(args.Davies_N1, args.Davies_N2, 
215                          mesh.lon_i, mesh.lat_i, mesh.lon_e,mesh.lat_e)
[790]216#        print('mask min/max :', davies.beta_i.min(), davies.beta_i.max() )
[774]217        logging.debug('mask min/max :')
[790]218        logging.debug('%f'% davies.beta_i.min())
219        logging.debug('%f'% davies.beta_i.max())
[795]220
[769]221        def next_flow(m,S,u,Phi,W):
222            # junk,fast,slow = caldyn.bwd_fast_slow(flow, 0.)
223            caldyn_step.mass[:,:], caldyn_step.theta_rhodz[:,:], caldyn_step.u[:,:] = m,S,u
224            caldyn_step.geopot[:,:], caldyn_step.W[:,:] = Phi,W
[771]225            for i in range(nt):
226                caldyn_step.next()
227                davies.relax(llm, caldyn_step, flow0)
[792]228                m,S,u = caldyn_step.mass, caldyn_step.theta_rhodz, caldyn_step.u
[785]229                s = S/m
[792]230                laps, bilaps = mesh.field_mass(), mesh.field_mass()
231                lapu, bilapu = mesh.field_u(), mesh.field_u()
[785]232                unst.ker.dynamico_scalar_laplacian(s,laps)
233                unst.ker.dynamico_scalar_laplacian(laps,bilaps)
[792]234                unst.ker.dynamico_curl_laplacian(u,lapu)
235                unst.ker.dynamico_curl_laplacian(lapu,bilapu)
[795]236                caldyn_step.theta_rhodz[:] = S - dt*args.kappa_divgrad*bilaps*m # Euler step
237                caldyn_step.u[:] = u - dt*args.kappa_curlcurl*bilapu # Euler step
[792]238
[769]239            return (caldyn_step.mass.copy(), caldyn_step.theta_rhodz.copy(), caldyn_step.u.copy(),
[771]240                    caldyn_step.geopot.copy(), caldyn_step.W.copy())
[769]241
[774]242    m,S,u,Phi,W = flow0
[769]243    if caldyn_thermo == unst.thermo_theta:
244        s=S/m
[790]245        theta = thermo.T0*exp(s/thermo.Cpd)
[769]246        S=m*theta
247        title_format = 'Potential temperature at t=%g s (K)'
248    else:
249        title_format = 'Specific entropy at t=%g s (J/K/kg)'
250
251    #T, nslice, dt = 3600., 1, 3600.
[774]252    Nslice = args.N
[769]253
[795]254    temp_v = mesh.field_z(),     
255
[769]256    with xios.Context_Curvilinear(mesh,1, 24*3600) as context:
[763]257        # now XIOS knows about the mesh and we can write to disk
[792]258        for i in range(Nslice+1):
[785]259            context.update_calendar(i+1)
[769]260
261            # Diagnose quantities of interest from prognostic variables m,S,u,Phi,W
[795]262            gas, w, z, ps, un, zeta_vk, zeta_abs_vk = diagnose(Phi,S,m,W,u)
[790]263            KE_ik = KE(mesh,u)
[787]264            du_fast, du_slow = caldyn_step.du_fast[0,:,:], caldyn_step.du_slow[0,:,:] 
265            div_fast, div_slow = div(mesh,du_fast), div(mesh,du_slow)
[790]266            drhodz, hflux = caldyn_step.drhodz[0,:,:], caldyn_step.hflux[:,:]
267
[769]268            # write to disk
[787]269            context.send_field_primal('ps', ps)
[785]270            context.send_field_primal('temp', gas.T)
[790]271
[769]272            context.send_field_primal('p', gas.p)
[790]273#            context.send_field_primal('p', KE_ik)
274#            context.send_field_primal('p', drhodz)
275
[800]276            context.send_field_primal('theta', thermo.T0*exp(gas.s/thermo.Cpd))
277
[769]278            context.send_field_primal('uz', w)
[780]279            context.send_field_primal('z', z)
[787]280            context.send_field_primal('div_fast', div_fast)
[790]281
[787]282            context.send_field_primal('div_slow', div_slow)
[790]283
[795]284            context.send_field_dual('curl', zeta_vk)
285            context.send_field_dual('curl_abs', zeta_abs_vk)
[769]286
[774]287            print( 'ptop, model top (m) :', unst.getvar('ptop'), Phi.max()/unst.getvar('g'))
[769]288
289            time1, elapsed1 =time.time(), unst.getvar('elapsed')
290            m,S,u,Phi,W = next_flow(m,S,u,Phi,W)
291            time2, elapsed2 =time.time(), unst.getvar('elapsed')
292            factor = 1000./nt
[774]293            print( 'ms per full time step : ', factor*(time2-time1), factor*(elapsed2-elapsed1))
[769]294            factor = 1e9/(4*nt*nx*ny*llm)
[774]295            print( 'nanosec per gridpoint per full time step : ', factor*(time2-time1), factor*(elapsed2-elapsed1))
[769]296
[774]297logging.info('************DONE************')
Note: See TracBrowser for help on using the repository browser.