-
Notifications
You must be signed in to change notification settings - Fork 0
/
roms_run_t.py
319 lines (263 loc) · 11.4 KB
/
roms_run_t.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
# -*- coding: utf-8 -*-
"""
Script for creating the input files for a time-varying seamount problem.
Last Modified 29 July 17
@author: bperfect
"""
import hgrid
import vgrid
from create_ini import create_ini
from grid import *
import numpy as np
from nc_create_grid import *
import netCDF4 as netCDF
from scipy import integrate
dirname = '/home/bperfect/seamounts/ROMS/seamount/'
restart = False
el = 120000 #90000 #y-domain width
xl = 180000 #180000 #x-domain width
u_inf = 0.1
u_offset = 3000.0
u_scale = 450.0
grid_x = 3*180 #210
grid_y = 3*120 #140
grid_z = 80 #60
rst_file = dirname + 'ocean_rst.nc'
ini_file = dirname + 'ocean_ini_f3e5n1e3.nc'
grd_file = dirname + 'ocean_grd_f3e5n1e3.nc'
bry_file = dirname + 'ocean_bry_f3e5n1e3.nc'
tstart = -1
tend = 401
dt = 10 # must reduce this if ramping up
time = np.arange(tstart,tend,dt)
time_fac = 1.0+0.0*time
#time_fac = 0.5*np.tanh(time/4.0-2.5)+0.5
inifac = time_fac[-tstart]
ntimes = time.size
# f and N, in radians per second
f = 3.0e-5
N = 1.0e-3
#Constants
rho_to_t=5.73
rho0=1025.0
g=9.81
depth=5000
sigma=10000
height = 3500
Fr = u_inf/N/height
print('Froude number is ' + str(Fr))
Ro = u_inf/f/2.355/sigma
print(('Rossby number is ' + str(Ro)))
def calcH(horGrid):
height=3500
sigma=10000
return depth-height*np.exp(-(horGrid.x_rho-xl*1/4)**2/sigma**2-(horGrid.y_rho-el/2)**2/sigma**2)
def calcU(z):
#return z*0+0.1
return u_inf+z*0#/2.0+u_inf/2.0*np.tanh((z+u_offset)/u_scale)
def calcDUDZ(z):
return 0*z#u_inf/2.0/u_scale*(1-(np.tanh((z+u_offset)/u_scale))**2)
def calcBackground(z):
''' We would like the background stratification to be determined by N, but we need to ensure the thermal wind
doesn't create an unstable stratification. We need to create a function that transforms a vector z into a
density. z does not necessarily start from the depth (because of the seamount), so it becomes necessary to
calculate the density at z[0] and then do a cumulative integration of the minimum allowable stratification
to obtain the density at the desired z levels
'''
z_min = z[0]
terrain_elevation = np.arange(-depth,z_min,10)
# N_sq_min_terrain = f_overall_max*el*u_inf/2.0/u_scale**2*np.tanh((terrain_elevation+u_offset)/u_scale)*(1.0-(np.tanh((terrain_elevation+u_offset)/u_scale))**2)
return 17.0+integrate.trapz(rho_to_t*rho0/g*np.asarray(N**2+terrain_elevation*0),terrain_elevation)+integrate.cumtrapz(rho_to_t*rho0/g*np.asarray(N**2+z*0),z,initial=0)
# z_min = z[0]
# terrain_elevation = np.arange(-depth,z_min,10)
# f_overall_max = 0.00013
# # obtain minimum allowable N^2 based on thermal wind restriction
# N_sq_min = f_overall_max*el*u_inf/2.0/u_scale**2*np.tanh((z+u_offset)/u_scale)*(1.0-(np.tanh((z+u_offset)/u_scale))**2)
# N_sq = [float(max(N**2,np.abs(x)*1.05)) for x in N_sq_min]
# #Compute starting density at the lowest Z (because terrain is raised off of the sea floor)
# N_sq_min_terrain = f_overall_max*el*u_inf/2.0/u_scale**2*np.tanh((terrain_elevation+u_offset)/u_scale)*(1.0-(np.tanh((terrain_elevation+u_offset)/u_scale))**2)
# N_sq_terrain = [float(max(N**2,np.abs(x)*1.05)) for x in N_sq_min_terrain]
# return 17.0+integrate.trapz(rho_to_t*rho0/g*np.asarray(N_sq_terrain),terrain_elevation)+integrate.cumtrapz(rho_to_t*rho0/g*np.asarray(N_sq),z,initial=0)
def calcThermalWind(y,z):
return -rho_to_t*rho0/g*f*(y-el/2)*calcDUDZ(z)
def calcFS(y,u):
return -f/9.81*y*u
def calcV(arg):
return arg*0.0
theta_b = 3
theta_s = 0.65
#Follow the same code, for the refined grid
#Create the original array
grid1 = np.mgrid[-1:grid_y+1:(grid_y+3)*1j,-1:grid_x+1:(grid_x+3)*1j]
xgrid=grid1[0]
ygrid=grid1[1]
#generate horizontal grid object
horGrid1 = CGrid(np.around(grid1[1])*xl/grid_x,np.around(grid1[0])*el/grid_y)
h=calcH(horGrid1)
#generate vertical grid object
s4 = s_coordinate_4(h, theta_b,theta_s, 100, grid_z) #theta_b, theta_s, Tcline, N,
#combine the horizontal and vertical grid to make a ROMS_Grid object
grd = ROMS_Grid("Seamount grid", horGrid1, s4)
#write grid object to a netCDF file
nc_create_grid(grd_file, grd, f, True)
#write initial condition netCDF file
create_ini(ini_file, grd)
ini = netCDF.Dataset(ini_file,'a',format='NETCDF3_64BIT')
# initial file variables
create_ini(bry_file,grd)
bry = netCDF.Dataset(bry_file,'a',format='NETCDF3_64BIT')
# Boundary-specific variables
bry.createDimension('v2d_time', np.size(time))
bry.createDimension('v3d_time', np.size(time))
bry.createDimension('zeta_time', np.size(time))
bry.createDimension('temp_time', np.size(time))
bry.createVariable('v2d_time', 'f8', ('v2d_time'))
bry.variables['v2d_time'][:]=time
bry.variables['v2d_time'].units = 'days'
bry.createVariable('v3d_time', 'f8', ('v3d_time'))
bry.variables['v3d_time'][:]=time
bry.variables['v3d_time'].units = 'days'
bry.createVariable('zeta_time', 'f8', ('zeta_time'))
bry.variables['zeta_time'][:]=time
bry.variables['zeta_time'].units = 'mdays'
bry.createVariable('temp_time', 'f8', ('temp_time'))
bry.variables['temp_time'][:]=time
bry.variables['temp_time'].units = 'days'
bry.variables['ocean_time'][:]=time
bry.createVariable('ubar_west', 'f8', ('v2d_time', 'eta_u'))
bry.createVariable('ubar_east', 'f8', ('v2d_time', 'eta_u'))
bry.createVariable('ubar_north', 'f8', ('v2d_time', 'xi_u'))
bry.createVariable('ubar_south', 'f8', ('v2d_time', 'xi_u'))
bry.createVariable('vbar_west', 'f8', ('v2d_time', 'eta_v'))
bry.createVariable('vbar_east', 'f8', ('v2d_time', 'eta_v'))
bry.createVariable('vbar_north', 'f8', ('v2d_time', 'xi_v'))
bry.createVariable('vbar_south', 'f8', ('v2d_time', 'xi_v'))
bry.createVariable('u_west', 'f8', ('v3d_time', 's_rho', 'eta_u'))
bry.createVariable('u_east', 'f8', ('v3d_time', 's_rho', 'eta_u'))
bry.createVariable('u_north', 'f8', ('v3d_time', 's_rho', 'xi_u'))
bry.createVariable('u_south', 'f8', ('v3d_time', 's_rho', 'xi_u'))
bry.createVariable('v_west', 'f8', ('v3d_time', 's_rho', 'eta_v'))
bry.createVariable('v_east', 'f8', ('v3d_time', 's_rho', 'eta_v'))
bry.createVariable('v_south', 'f8', ('v3d_time', 's_rho', 'xi_v'))
bry.createVariable('v_north', 'f8', ('v3d_time', 's_rho', 'xi_v'))
bry.createVariable('temp_west', 'f8', ('temp_time', 's_rho', 'eta_rho'))
bry.createVariable('temp_east', 'f8', ('temp_time', 's_rho', 'eta_rho'))
bry.createVariable('temp_north', 'f8', ('temp_time', 's_rho', 'xi_rho'))
bry.createVariable('temp_south', 'f8', ('temp_time', 's_rho', 'xi_rho'))
bry.createVariable('zeta_west', 'f8', ('zeta_time', 'eta_rho'))
bry.createVariable('zeta_east', 'f8', ('zeta_time', 'eta_rho'))
bry.createVariable('zeta_north', 'f8', ('zeta_time', 'xi_rho'))
bry.createVariable('zeta_south', 'f8', ('zeta_time', 'xi_rho'))
Cs_r = bry.variables['Cs_r'][:]
h = grd.vgrid.h
y_rho = grd.hgrid.y_rho
# h might have to be switched to the correct coordinates
u = ini.variables['u'][:]
ubar = ini.variables['ubar'][:]
for i in range(u.shape[1]):
for j in range(u.shape[2]):
ubar[i,j] = np.trapz(calcU(Cs_r*h[i,j]),Cs_r)
for k in range(u.shape[0]):
u[k,i,j] = calcU(Cs_r[k]*h[i,j])
ini.variables['u'][:] = u*inifac
ini.variables['ubar'][:] = ubar*inifac
#v = ini.variables['v'][:]
#vbar = ini.variables['vbar'][:]
#for i in range(v.shape[1]):
# for j in range(v.shape[2]):
# vbar[i,j] = np.trapz(calcV(Cs_r*h[i,j]),Cs_r)
# for k in range(v.shape[0]):
# v[k,i,j] = calcV(Cs_r[k]*h[i,j])
#ini.variables['v'][:] = v
#ini.variables['vbar'][:] = vbar
''' M2 Boundaries'''
#West
temp = ubar[:,1]
fac = np.tile(time_fac[:,np.newaxis],[1,temp.size])
value = np.transpose(np.repeat(temp[:,np.newaxis],np.size(time),axis=1))
bry.variables['ubar_west'][:] = np.multiply(fac,value)
bry.variables['vbar_west'][:] = 0.0*bry.variables['vbar_west'][:]
#East
temp = ubar[:,-1]
value = np.transpose(np.repeat(temp[:,np.newaxis],np.size(time),axis=1))
bry.variables['ubar_east'][:] = np.multiply(fac,value)
bry.variables['vbar_east'][:] = 0.0*bry.variables['vbar_east'][:]
#South
temp = ubar[1,:]
fac = np.tile(time_fac[:,np.newaxis],[1,temp.size])
value = np.transpose(np.repeat(temp[:,np.newaxis],np.size(time),axis=1))
bry.variables['ubar_south'][:] = np.multiply(fac,value)
bry.variables['vbar_south'][:] = 0.0*bry.variables['vbar_south'][:]
temp = ubar[-1,:]
value = np.transpose(np.repeat(temp[:,np.newaxis],np.size(time),axis=1))
bry.variables['ubar_north'][:] = np.multiply(fac,value)
bry.variables['vbar_north'][:] = 0.0*bry.variables['vbar_north'][:]
''' M3 Boundaries '''
temp = u[:,:,1]
fac = np.tile(time_fac[:,np.newaxis,np.newaxis],[1,temp.shape[0],temp.shape[1]])
value = np.repeat(temp[np.newaxis,:,:],np.size(time),axis=0)
bry.variables['u_west'][:] = np.multiply(fac,value)
bry.variables['v_west'][:] = 0.0*bry.variables['v_west'][:]
temp = u[:,:,-1]
value = np.repeat(temp[np.newaxis,:,:],np.size(time),axis=0)
bry.variables['u_east'][:] = np.multiply(fac,value)
bry.variables['v_east'][:] = 0.0*bry.variables['v_east'][:]
temp = u[:,1,:]
fac = np.tile(time_fac[:,np.newaxis,np.newaxis],[1,temp.shape[0],temp.shape[1]])
value = np.repeat(temp[np.newaxis,:,:],np.size(time),axis=0)
bry.variables['u_south'][:] = np.multiply(fac,value)
bry.variables['v_south'][:] = 0.0*bry.variables['v_south'][:]
temp = u[:,-1,:]
value = np.repeat(temp[np.newaxis,:,:],np.size(time),axis=0)
bry.variables['u_north'][:] = np.multiply(fac,value)
bry.variables['v_north'][:] = 0.0*bry.variables['v_north'][:]
''' Temp Boundaries '''
background = ini.variables['temp'][:]
thermalWind = deepcopy(background)
for i in range(background.shape[1]):
for j in range(background.shape[2]):
z=Cs_r*h[i,j]
y=y_rho[i,j]
background[:,i,j] = calcBackground(z)
thermalWind[:,i,j] = calcThermalWind(y,z)
temperature = background + thermalWind*inifac
ini.variables['temp'][:] = temperature
temp2 = background[:,:,1]
temp = thermalWind[:,:,1]
fac = np.tile(time_fac[:,np.newaxis,np.newaxis],[1,temp.shape[0],temp.shape[1]])
value = np.repeat(temp[np.newaxis,:,:],np.size(time),axis=0)
bry.variables['temp_west'][:] = np.multiply(fac,value)+temp2
temp = thermalWind[:,:,-1]
temp2 = background[:,:,-1]
value = np.repeat(temp[np.newaxis,:,:],np.size(time),axis=0)
bry.variables['temp_east'][:] = np.multiply(fac,value)+temp2
temp = thermalWind[:,1,:]
temp2 = background[:,1,:]
fac = np.tile(time_fac[:,np.newaxis,np.newaxis],[1,temp.shape[0],temp.shape[1]])
value = np.repeat(temp[np.newaxis,:,:],np.size(time),axis=0)
bry.variables['temp_south'][:] = np.multiply(fac,value)+temp2
temp = thermalWind[:,-1,:]
temp2 = background[:,-1,:]
value = np.repeat(temp[np.newaxis,:,:],np.size(time),axis=0)
bry.variables['temp_north'][:] = np.multiply(fac,value)+temp2
''' Zeta Boundaries'''
surfU = u[-1,1,1]
zeta = calcFS(y_rho,surfU)
ini.variables['zeta'][:] = zeta[:]*inifac
ini.variables['zeta'][:] = zeta[:]*inifac
temp = zeta[:,1]
fac = np.tile(time_fac[:,np.newaxis],[1,temp.shape[0]])
value = np.repeat(temp[np.newaxis,:],np.size(time),axis=0)
bry.variables['zeta_west'][:] = np.multiply(fac,value)
temp = zeta[:,-1]
value = np.repeat(temp[np.newaxis,:],np.size(time),axis=0)
bry.variables['zeta_east'][:] = np.multiply(fac,value)
temp = zeta[1,:]
fac = np.tile(time_fac[:,np.newaxis],[1,temp.shape[0]])
value = np.repeat(temp[np.newaxis,:],np.size(time),axis=0)
bry.variables['zeta_south'][:] = np.multiply(fac,value)
temp = zeta[-1,:]
value = np.repeat(temp[np.newaxis,:],np.size(time),axis=0)
bry.variables['zeta_north'][:] = np.multiply(fac,value)
if restart == True:
nc_ini_from_small_rst()