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

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

devel/unstructured : added kernel for curl curl ; used by Baroclinic_3D_ullrich

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