1 | from __future__ import print_function |
---|
2 | from dynamico import getargs |
---|
3 | |
---|
4 | getargs.add("--T", type=float, default=3600., |
---|
5 | help="Length of time slice in seconds") |
---|
6 | getargs.add("--perturb", type=float, default=1., |
---|
7 | help="Amplitude of wind perturbation in m/s") |
---|
8 | getargs.add("--N", type=int, default=48, |
---|
9 | help="Number of time slices") |
---|
10 | getargs.add("--Davies_N1", type=int, default=3) |
---|
11 | getargs.add("--Davies_N2", type=int, default=3) |
---|
12 | getargs.add("--nx", type=int, default=200) |
---|
13 | getargs.add("--ny", type=int, default=30) |
---|
14 | getargs.add("--betaplane", action='store_true') |
---|
15 | getargs.defaults(dt=360., llm=50) |
---|
16 | |
---|
17 | logging = getargs.Logger('main') |
---|
18 | args = getargs.parse() |
---|
19 | |
---|
20 | from dynamico import unstructured as unst |
---|
21 | from dynamico import dyn |
---|
22 | from dynamico import time_step |
---|
23 | from dynamico import DCMIP |
---|
24 | from dynamico import meshes |
---|
25 | from dynamico import xios |
---|
26 | from dynamico import precision as prec |
---|
27 | from dynamico.meshes import Cartesian_mesh as Mesh |
---|
28 | from dynamico.kernels import grad, curl, div, KE |
---|
29 | from dynamico.LAM import Davies |
---|
30 | |
---|
31 | import math as math |
---|
32 | import numpy as np |
---|
33 | import time |
---|
34 | from numpy import pi, log, exp, sin, cos |
---|
35 | |
---|
36 | # Baroclinic instability test based on Ullrich et al. 2015, QJRMS |
---|
37 | |
---|
38 | |
---|
39 | def 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. |
---|
47 | y0 = .5*Ly |
---|
48 | Cpd = 1004.5 |
---|
49 | p0 = 1e5 |
---|
50 | up = args.perturb # amplitude (m/s) |
---|
51 | xc,yc,lp = 0.,Ly/2.,600000. |
---|
52 | |
---|
53 | omega = 7.292e-5 # Angular velocity of the Earth (s^-1) |
---|
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. |
---|
57 | fb = f0 - beta0*y0 |
---|
58 | |
---|
59 | def Phi_xy(y): |
---|
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 | |
---|
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)) |
---|
65 | def ulon(x,y,eta): |
---|
66 | u = -u0*(sin(pi*y/Ly)**2)*log(eta)*(eta**(-log(eta)/(b*b))) |
---|
67 | u = u + up*exp(-(((x-xc)**2+(y-yc)**2)/lp**2)) |
---|
68 | return u |
---|
69 | |
---|
70 | def tmean(eta) : return T0*eta**(Rd*lap/g) |
---|
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)) |
---|
72 | def p(eta): return p0*eta # eta = p/p_s |
---|
73 | |
---|
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 |
---|
76 | |
---|
77 | filename = 'cart_%03d_%03d.nc'%(nx,ny) |
---|
78 | logging.info('Reading Cartesian mesh ...') |
---|
79 | meshfile = meshes.DYNAMICO_Format(filename) |
---|
80 | pmesh = meshes.Unstructured_PMesh(comm,meshfile) |
---|
81 | pmesh.partition_curvilinear(args.mpi_ni,args.mpi_nj) |
---|
82 | nqdyn, radius = 1, None |
---|
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 |
---|
87 | x_ik, alpha_ik = np.meshgrid(mesh.lon_i, alpha_k, indexing='ij') |
---|
88 | y_ik, alpha_ik = np.meshgrid(mesh.lat_i+y0, 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+y0, 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+y0, 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_il = p(eta_ik), p(eta_il) |
---|
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 | # mass_ik[:,l]=(p_il[:,l]-p_il[:,l+1])/g |
---|
113 | Sik, Wil = gas.s*mass_ik, mesh.field_w() |
---|
114 | |
---|
115 | u_ek = mesh.ucov3D(ulon(x_ek, y_ek, eta_ek), 0.*eta_ek) |
---|
116 | |
---|
117 | print('ztop (m) = ', Phi_il[0,-1]/g, ztop) |
---|
118 | ptop = p(eta(1.)) |
---|
119 | print( 'ptop (Pa) = ', gas.p[0,-1], ptop) |
---|
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]) |
---|
129 | print( 'min p, T :', pbot.min(), tmean(pbot/p0)) |
---|
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 | |
---|
134 | return thermo, mesh, params, prec.asnum([mass_ik,Sik,u_ek,Phi_il,Wil]), gas |
---|
135 | |
---|
136 | def diagnose(Phi,S,m,W): |
---|
137 | s=S/m |
---|
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 | |
---|
145 | class myDavies(Davies): |
---|
146 | def mask(self,x,y): |
---|
147 | logging.debug('davies dy = %f'% dy) |
---|
148 | return self.mask0(y,.5*Ly,dy)*self.mask0(-y,.5*Ly,dy) |
---|
149 | |
---|
150 | with 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() |
---|
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' ) |
---|
155 | |
---|
156 | logging.info('%d/%d starting'%(mpi_rank,mpi_size)) |
---|
157 | |
---|
158 | g, Lx, Ly = 9.81, 4e7, 6e6 |
---|
159 | nx, ny, llm = args.nx, args.ny, args.llm |
---|
160 | dx,dy = Lx/nx,Ly/ny |
---|
161 | |
---|
162 | unst.setvar('g',g) |
---|
163 | |
---|
164 | thermo, mesh, params, flow0, gas0 = baroclinic_3D(Lx,nx,Ly,ny,llm) |
---|
165 | |
---|
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) |
---|
168 | |
---|
169 | T = args.T |
---|
170 | dt = args.dt |
---|
171 | dz = flow0[3].max()/(params.g*llm) |
---|
172 | nt = int(math.ceil(T/dt)) |
---|
173 | dt = T/nt |
---|
174 | logging.info( 'Time step : %d x %g = %g s' % (nt,dt,nt*dt)) |
---|
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) |
---|
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) |
---|
190 | davies = myDavies(args.Davies_N1, args.Davies_N2, mesh.lon_i, mesh.lat_i, mesh.lon_e,mesh.lat_e) |
---|
191 | # print('mask min/max :', davies.beta_i.min(), davies.beta_i.max() ) |
---|
192 | logging.debug('mask min/max :') |
---|
193 | logging.debug('%f'% davies.beta_i.min()) |
---|
194 | logging.debug('%f'% davies.beta_i.max()) |
---|
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 |
---|
199 | kappa = dx**4/432000. |
---|
200 | for i in range(nt): |
---|
201 | caldyn_step.next() |
---|
202 | davies.relax(llm, caldyn_step, flow0) |
---|
203 | m,S,u = caldyn_step.mass, caldyn_step.theta_rhodz, caldyn_step.u |
---|
204 | s = S/m |
---|
205 | laps, bilaps = mesh.field_mass(), mesh.field_mass() |
---|
206 | lapu, bilapu = mesh.field_u(), mesh.field_u() |
---|
207 | unst.ker.dynamico_scalar_laplacian(s,laps) |
---|
208 | unst.ker.dynamico_scalar_laplacian(laps,bilaps) |
---|
209 | unst.ker.dynamico_curl_laplacian(u,lapu) |
---|
210 | unst.ker.dynamico_curl_laplacian(lapu,bilapu) |
---|
211 | caldyn_step.theta_rhodz[:] = S - dt*kappa*bilaps*m # Euler step |
---|
212 | caldyn_step.u[:] = u - dt*kappa*bilapu # Euler step |
---|
213 | |
---|
214 | return (caldyn_step.mass.copy(), caldyn_step.theta_rhodz.copy(), caldyn_step.u.copy(), |
---|
215 | caldyn_step.geopot.copy(), caldyn_step.W.copy()) |
---|
216 | |
---|
217 | m,S,u,Phi,W = flow0 |
---|
218 | if caldyn_thermo == unst.thermo_theta: |
---|
219 | s=S/m |
---|
220 | theta = thermo.T0*exp(s/thermo.Cpd) |
---|
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. |
---|
230 | Nslice = args.N |
---|
231 | |
---|
232 | with xios.Context_Curvilinear(mesh,1, 24*3600) as context: |
---|
233 | # now XIOS knows about the mesh and we can write to disk |
---|
234 | v = mesh.field_mass() # specific volume (diagnosed) |
---|
235 | |
---|
236 | for i in range(Nslice+1): |
---|
237 | context.update_calendar(i+1) |
---|
238 | |
---|
239 | # Diagnose quantities of interest from prognostic variables m,S,u,Phi,W |
---|
240 | gas, w, z = diagnose(Phi,S,m,W) |
---|
241 | ps = caldyn_step.ps |
---|
242 | curl_vk = curl(mesh, u) |
---|
243 | KE_ik = KE(mesh,u) |
---|
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) |
---|
246 | drhodz, hflux = caldyn_step.drhodz[0,:,:], caldyn_step.hflux[:,:] |
---|
247 | |
---|
248 | # write to disk |
---|
249 | context.send_field_primal('ps', ps) |
---|
250 | context.send_field_primal('temp', gas.T) |
---|
251 | |
---|
252 | context.send_field_primal('p', gas.p) |
---|
253 | # context.send_field_primal('p', KE_ik) |
---|
254 | # context.send_field_primal('p', drhodz) |
---|
255 | |
---|
256 | context.send_field_primal('theta', gas.s) |
---|
257 | context.send_field_primal('uz', w) |
---|
258 | context.send_field_primal('z', z) |
---|
259 | context.send_field_primal('div_fast', div_fast) |
---|
260 | |
---|
261 | context.send_field_primal('div_slow', div_slow) |
---|
262 | # context.send_field_primal('div_slow', div(mesh,hflux) ) |
---|
263 | |
---|
264 | context.send_field_dual('curl', curl_vk) |
---|
265 | |
---|
266 | print( 'ptop, model top (m) :', unst.getvar('ptop'), Phi.max()/unst.getvar('g')) |
---|
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 |
---|
272 | print( 'ms per full time step : ', factor*(time2-time1), factor*(elapsed2-elapsed1)) |
---|
273 | factor = 1e9/(4*nt*nx*ny*llm) |
---|
274 | print( 'nanosec per gridpoint per full time step : ', factor*(time2-time1), factor*(elapsed2-elapsed1)) |
---|
275 | |
---|
276 | logging.info('************DONE************') |
---|