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

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

devel/Python : Implemented Davies relaxation in Baroclinic_3D_ullrich.py

File size: 9.6 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
14from numpy import pi, log, exp, sin, cos
15
16# Baroclinic instability test based on Ullrich et al. 2015, QJRMS
17
18parser = argparse.ArgumentParser()
19
20parser.add_argument("--mpi_ni", type=int, default=1,
21                    help="number of x processors")
22parser.add_argument("--mpi_nj", type=int, default=1,
23                    help="number of y processors")
24parser.add_argument("--T", type=float, default=5.,
25                    help="Length of time slice in seconds")
26parser.add_argument("--Davies_N1", type=int, default=5)
27parser.add_argument("--Davies_N2", type=int, default=5)
28
29
30args = parser.parse_args()
31
32
33def baroclinic_3D(Lx,nx,Ly,ny,llm,ztop=25000.):
34    Rd    = 287.0      # Gas constant for dryy air (j kg^-1 K^-1)
35    T0    = 288.0      # Reference temperature (K)
36    lap   = 0.005      # Lapse rate (K m^-1)
37    b     = 2.         # Non dimensional vertical width parameter
38    u0    = 35.        # Reference zonal wind speed (m s^-1)
39    a     = 6.371229e6 # Radius of the Earth (m)
40    ptop  = 2000.
41    y0    = Ly*0.5
42    Cpd   = 1004.5
43    p0    = 1e5
44
45    omega  = 7.292e-5  # Angular velocity of the Earth (s^-1)   
46    phi0   = 90.*np.pi/180.0 # Reference latitude North pi/4 (deg)
47    f0     = 2*omega*np.sin(phi0) 
48    beta0  = 2*omega*np.cos(phi0)/a
49    fb     = 2*omega*np.sin(phi0) - y0*2*omega*np.cos(phi0)/a
50
51    def Phi_xy(y):
52        fc = y*y - (Ly*y/pi)*sin(2*pi*y/Ly)
53        fd = Ly*Ly/(2*pi*pi)*cos(2*pi*y/Ly)
54        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)) )
55
56    def Phi_xyeta(y,eta): return T0*g/lap*(1-eta**(Rd*lap/g)) + Phi_xy(y)*log(eta)*exp(-((log(eta)/b)**2))
57    def ulon(x,y,eta): return -u0*(sin(pi*y/Ly)**2)*log(eta)*(eta**(-log(eta)/b/b))
58    def tmean(eta) : return T0*eta**(Rd*lap/g)
59    def T(y,eta) : return tmean(eta)+(Phi_xy(y)/Rd)*(((2/(b*b))*(log(eta))**2)-1)*exp(-((0.5*log(eta))**2)) 
60    def p(eta): return p0*eta     # eta  = p/p_s
61
62    def eta(alpha) : return (1-(lap*ztop*alpha/(T0)))**(g/(Rd*lap)) # roughly equispaced levels
63
64    filename = 'cart_%03d_%03d.nc'%(nx,ny)
65    print 'Reading Cartesian mesh ...'
66    def coriolis(x,y): return f0+beta0*(y+.5*Ly)
67    meshfile = meshes.DYNAMICO_Format(filename)
68    pmesh = meshes.Unstructured_PMesh(comm,meshfile)
69    pmesh.partition_curvilinear(args.mpi_ni,args.mpi_nj)
70    nqdyn, radius = 1, None
71    mesh  = meshes.Local_Mesh(pmesh, llm, nqdyn, radius, coriolis)
72
73
74    alpha_k = (np.arange(llm)  +.5)/llm
75    alpha_l = (np.arange(llm+1)+ 0.)/llm
76    x_ik, alpha_ik = np.meshgrid(mesh.lon_i, alpha_k, indexing='ij')
77    y_ik, alpha_ik = np.meshgrid(mesh.lat_i, alpha_k, indexing='ij')
78    x_il, alpha_il = np.meshgrid(mesh.lon_i, alpha_l, indexing='ij')
79    y_il, alpha_il = np.meshgrid(mesh.lat_i, alpha_l, indexing='ij')
80    x_ek, alpha_ek = np.meshgrid(mesh.lon_e, alpha_k, indexing='ij')
81    y_ek, alpha_ek = np.meshgrid(mesh.lat_e, alpha_k, indexing='ij')
82
83    print('----------------')
84    print 'ztop(ptop) according to Eq. 7:', T0/lap*(1.-(ptop/p0)**(Rd*lap/g)) 
85    print(np.shape(alpha_k),np.shape(alpha_l))
86    thermo = dyn.Ideal_perfect(Cpd, Rd, p0, T0)
87
88    eta_il = eta(alpha_il)
89    eta_ik = eta(alpha_ik)
90    eta_ek = eta(alpha_ek)
91    print('min max eta_il', np.min(eta_il),np.max(eta_il))
92
93    Phi_il = Phi_xyeta(y_il, eta_il)
94    Phi_ik = Phi_xyeta(y_ik, eta_ik)
95    p_ik = p(eta_ik)
96    T_ik = T(y_ik, eta_ik)  #ik full level(40), il(41)
97
98    gas = thermo.set_pT(p_ik,T_ik)
99    mass_ik = mesh.field_mass()
100    for l in range(llm): 
101        mass_ik[:,l]=(Phi_il[:,l+1]-Phi_il[:,l])/(g*gas.v[:,l])
102    Sik, ujk, Wil  = gas.s*mass_ik, mesh.field_u(), mesh.field_w()
103   
104    u_ek = mesh.ucov3D(ulon(x_ek, y_ek, eta_ek), 0.*eta_ek)
105
106    print 'ztop (m) = ', Phi_il[0,-1]/g, ztop
107    ptop = p(eta(1.))
108    print 'ptop (Pa) = ', gas.p[0,-1], ptop
109
110    params=dyn.Struct()
111    params.ptop=ptop
112    params.dx=dx
113    params.dx_g0=dx/g
114    params.g = g
115
116    # define parameters for lower BC
117    pbot = p(eta_il[:,0])
118    print 'min p, T :', pbot.min(), tmean(pbot/p0)
119    gas_bot = thermo.set_pT(pbot, tmean(pbot/p0))
120    params.pbot = gas_bot.p
121    params.rho_bot = 1e6/gas_bot.v
122   
123    return thermo, mesh, params, prec.asnum([mass_ik,Sik,ujk,Phi_il,Wil]), gas
124
125def diagnose(Phi,S,m,W):
126    s=S/m ; s=.5*(s+abs(s))
127    for l in range(llm):
128        v[:,l]=(Phi[:,l+1]-Phi[:,l])/(g*m[:,l])
129        w[:,l]=.5*params.g*(W[:,l+1]+W[:,l])/m[:,l]
130        z[:,l]=.5*(Phi[:,l+1]+Phi[:,l])/params.g
131    gas = thermo.set_vs(v,s)
132    return gas, w, z
133
134class Davies:
135    def __init__(self,N1,N2,x_i,y_i,x_e,y_e):
136        self.N1, self.N2 = N1, N2
137        self.beta_i = self.mask(x_i,y_i)
138        self.beta_e = self.mask(x_e,y_e)
139    def mask0(self,c,c0,delta): # 1D building block for Davies relaxation
140        N1, N2 = self.N1, self.N2
141        N3=N1+N2
142        m = np.zeros(c.size)
143
144        for i in range(c.size):
145            ci=c[i]
146            m[i] = (1.+np.cos((ci-c0+N3*delta)*np.pi/(N2*delta)))/2.0
147            if ci < c0-N3*delta : m[i]=1.
148            if ci > c0-N1*delta : m[i]=0.
149        return m
150    def relax(self, llm, step, flow):
151        beta_i, beta_e = self.beta_i, self.beta_e
152        m,S,u,Phi,W=flow
153        for l in range(llm):
154            step.mass[:,l]        = beta_i*step.mass[:,l]        + (1.-beta_i)*m[:,l]
155            step.theta_rhodz[:,l] = beta_i*step.theta_rhodz[:,l] + (1.-beta_i)*S[:,l]
156            step.u[:,l]           = beta_e*step.u[:,l]           + (1.-beta_e)*u[:,l]
157        for l in range(llm+1):
158            step.geopot[:,l]      = beta_i*step.geopot[:,l]      + (1.-beta_i)*Phi[:,l]
159            step.W[:,l]           = beta_i*step.W[:,l]           + (1.-beta_i)*W[:,l]
160
161
162class myDavies(Davies):
163    def mask(self,x,y):
164#        return self.mask0(y,Ly,dy)
165        return self.mask0(y,-.5*Ly,dy)*self.mask0(-y,-.5*Ly,dy)
166       
167with xios.Client() as client: # setup XIOS which creates the DYNAMICO communicator
168    comm = client.comm
169    mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size()
170    print '%d/%d starting'%(mpi_rank,mpi_size)
171
172    g, Lx, Ly = 9.81, 4e7, 6e6
173    nx, ny, llm = 200, 30, 25
174    dx,dy=Lx/nx,Ly/ny
175
176    unst.setvar('g',g)
177
178    thermo, mesh, params, flow0, gas0 = baroclinic_3D(Lx,nx,Ly,ny,llm)
179
180    mass_bl,mass_dak,mass_dbk = meshes.compute_hybrid_coefs(flow0[0])
181    print 'Type of mass_bl, mass_dak, mass_dbk : ', [x.dtype for x in mass_bl, mass_dak, mass_dbk]
182    unst.ker.dynamico_init_hybrid(mass_bl,mass_dak,mass_dbk)
183
184    T=3600.
185    dt=360.
186    dz = flow0[3].max()/(params.g*llm)
187    nt = int(math.ceil(T/dt))
188    dt = T/nt
189    print 'Time step : %d x %g = %g s' % (nt,dt,nt*dt)
190   
191    caldyn_thermo, caldyn_eta = unst.thermo_entropy, unst.eta_mass
192
193    if False: # time stepping in Python
194        caldyn = unst.Caldyn_NH(caldyn_thermo,caldyn_eta, mesh,thermo,params,params.g)
195        scheme = time_step.ARK2(caldyn.bwd_fast_slow, dt)
196        def next_flow(m,S,u,Phi,W):
197            return scheme.advance((m,S,u,Phi,W),nt)
198
199    else: # time stepping in Fortran
200        scheme = time_step.ARK2(None, dt)
201        caldyn_step = unst.caldyn_step_NH(mesh,scheme,1, caldyn_thermo,caldyn_eta,thermo,params,params.g)
202        davies = myDavies(args.Davies_N1, args.Davies_N2, mesh.lon_i, mesh.lat_i, mesh.lon_e,mesh.lat_e)
203        def next_flow(m,S,u,Phi,W):
204            # junk,fast,slow = caldyn.bwd_fast_slow(flow, 0.)
205            caldyn_step.mass[:,:], caldyn_step.theta_rhodz[:,:], caldyn_step.u[:,:] = m,S,u
206            caldyn_step.geopot[:,:], caldyn_step.W[:,:] = Phi,W
207            for i in range(nt):
208                caldyn_step.next()
209                davies.relax(llm, caldyn_step, flow0)
210            return (caldyn_step.mass.copy(), caldyn_step.theta_rhodz.copy(), caldyn_step.u.copy(),
211                    caldyn_step.geopot.copy(), caldyn_step.W.copy())
212
213    m,S,u,Phi,W=flow0
214    if caldyn_thermo == unst.thermo_theta:
215        s=S/m
216        theta = thermo.T0*np.exp(s/thermo.Cpd)
217        S=m*theta
218        title_format = 'Potential temperature at t=%g s (K)'
219    else:
220        title_format = 'Specific entropy at t=%g s (J/K/kg)'
221
222    w=mesh.field_mass()
223    z=mesh.field_mass()
224
225
226    #T, nslice, dt = 3600., 1, 3600.
227    Nslice=24
228
229    with xios.Context_Curvilinear(mesh,1, 24*3600) as context:
230        # now XIOS knows about the mesh and we can write to disk
231        v = mesh.field_mass() # specific volume (diagnosed)
232       
233        for i in range(Nslice):
234            context.update_calendar(i)
235
236            # Diagnose quantities of interest from prognostic variables m,S,u,Phi,W
237            gas, w, z = diagnose(Phi,S,m,W)
238
239            # write to disk
240            context.send_field_primal('temp', gas0.T)
241            context.send_field_primal('p', gas.p)
242            context.send_field_primal('theta', gas.s)
243            context.send_field_primal('uz', w)
244
245            print 'ptop, model top (m) :', unst.getvar('ptop'), Phi.max()/unst.getvar('g')
246            #if args.mpi_ni*args.mpi_nj==1: plot()
247
248            time1, elapsed1 =time.time(), unst.getvar('elapsed')
249            m,S,u,Phi,W = next_flow(m,S,u,Phi,W)
250            time2, elapsed2 =time.time(), unst.getvar('elapsed')
251            factor = 1000./nt
252            print 'ms per full time step : ', factor*(time2-time1), factor*(elapsed2-elapsed1)
253            factor = 1e9/(4*nt*nx*ny*llm)
254            print 'nanosec per gridpoint per full time step : ', factor*(time2-time1), factor*(elapsed2-elapsed1)
255
256        context.update_calendar(Nslice+1) # make sure XIOS writes last iteration
257print('************DONE************')
Note: See TracBrowser for help on using the repository browser.