prof_TS_z_box.py 5.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212
  1. #!/usr/bin/env python
  2. # L. Brodeau, June 2014
  3. import sys
  4. import os
  5. import numpy as nmp
  6. from netCDF4 import Dataset
  7. from os.path import basename
  8. import barakuda_tool as bt
  9. import barakuda_physics as bphy
  10. venv_needed = {'ORCA','EXP','DIAG_D','FILE_DEF_BOXES','CPREF','MM_FILE',
  11. 'NN_T','NN_S'}
  12. vdic = bt.check_env_var(sys.argv[0], venv_needed)
  13. CONFEXP = vdic['ORCA']+'-'+vdic['EXP']
  14. cname_script = basename(sys.argv[0])
  15. print '\n'+cname_script
  16. narg = len(sys.argv)
  17. if narg < 2: print 'Usage: '+cname_script+' <year>'; sys.exit(0)
  18. cy = sys.argv[1] ; jy=int(cy)
  19. # First will read name and coordinates of rectangular boxes to treat into file FILE_DEF_BOXES
  20. ##############################################################################################
  21. vboxes, vi1, vj1, vi2, vj2 = bt.read_coor(vdic['FILE_DEF_BOXES'])
  22. nbb = len(vboxes)
  23. print ''
  24. # Opening mesh-mask and T-grid file of NEMO:
  25. ############################################
  26. bt.chck4f(vdic['MM_FILE'])
  27. id_mm = Dataset(vdic['MM_FILE'])
  28. Xmask = id_mm.variables['tmask'][0,:,:,:]
  29. ze1t = id_mm.variables['e1t'] [0,:,:]
  30. ze2t = id_mm.variables['e2t'] [0,:,:]
  31. id_mm.close()
  32. cf_in = vdic['CPREF']+cy+'0101_'+cy+'1231_grid_T.nc' ; bt.chck4f(cf_in)
  33. id_in = Dataset(cf_in)
  34. Vtime = id_in.variables['time_counter'][:]
  35. Vdepth = id_in.variables['deptht'][:]
  36. Xlon = id_in.variables['nav_lon'][:,:]
  37. Xlat = id_in.variables['nav_lat'][:,:]
  38. Xtheta = id_in.variables[vdic['NN_T']][:,:,:,:]
  39. Xsali = id_in.variables[vdic['NN_S']][:,:,:,:]
  40. print '(has ',Xtheta.shape[0],' time snapshots)\n'
  41. id_in.close()
  42. [ Nt, nk, nj, ni ] = Xtheta.shape
  43. print 'Nt, nk, nj, ni =', Nt, nk, nj, ni
  44. if Nt == 12:
  45. for jt in range(Nt): Vtime[jt] = float(jy) + (float(jt) + 0.5)/float(Nt)
  46. # Loop along boxes:
  47. for jb in range(nbb):
  48. cbox = vboxes[jb]
  49. i1 = vi1[jb]
  50. j1 = vj1[jb]
  51. i2 = vi2[jb]+1
  52. j2 = vj2[jb]+1
  53. print '\n *** Treating box '+cbox+' => ', i1, j1, i2-1, j2-1
  54. # Filling box arrays:
  55. # ~~~~~~~~~~~~~~~~~~~
  56. nx_b = i2 - i1
  57. ny_b = j2 - j1
  58. shape_array = [ Nt, nk, ny_b, nx_b ]
  59. # Temporary arrays:
  60. Ztmp = nmp.zeros(ny_b*nx_b) ; Ztmp.shape = [ ny_b, nx_b ]
  61. Zar2 = nmp.zeros(ny_b*nx_b) ; Zar2.shape = [ ny_b, nx_b ]
  62. # Time series for each level of the mean sigma0 on the box:
  63. Zm1 = nmp.zeros(Nt*nk) ; Zm1.shape = [ Nt, nk ]
  64. Tm1 = nmp.zeros(Nt*nk) ; Tm1.shape = [ Nt, nk ]
  65. Sm1 = nmp.zeros(Nt*nk) ; Sm1.shape = [ Nt, nk ]
  66. Ztmp[:,:] = ze1t[j1:j2, i1:i2]*ze2t[j1:j2, i1:i2]
  67. for jk in range(nk):
  68. rarea = nmp.sum(Ztmp[:,:]*Xmask[jk, j1:j2, i1:i2])
  69. if rarea > 0.:
  70. Zar2[:,:] = Ztmp[:,:]*Xmask[jk, j1:j2, i1:i2]
  71. for jt in range(Nt):
  72. # Sigma0
  73. Zm1[jt,jk] = nmp.sum(bphy.sigma0(Xtheta[jt,jk, j1:j2, i1:i2],Xsali[jt,jk, j1:j2, i1:i2])*Zar2[:,:]) / rarea
  74. # Theta
  75. Tm1[jt,jk] = nmp.sum(Xtheta[jt,jk, j1:j2, i1:i2]*Zar2[:,:]) / rarea
  76. # Salinity
  77. Sm1[jt,jk] = nmp.sum( Xsali[jt,jk, j1:j2, i1:i2]*Zar2[:,:]) / rarea
  78. else:
  79. Zm1[:,jk] = nmp.nan ; Tm1[:,jk] = nmp.nan ; Sm1[:,jk] = nmp.nan
  80. # Writing in output file
  81. ########################
  82. cf_out = vdic['DIAG_D']+'/TS_z_box_'+cbox+'_'+CONFEXP+'.nc'
  83. l_nc_is_new = not os.path.exists(cf_out)
  84. # Creating/Opening output Netcdf file:
  85. if l_nc_is_new:
  86. f_out = Dataset(cf_out, 'w', format='NETCDF3_CLASSIC')
  87. else:
  88. f_out = Dataset(cf_out, 'a', format='NETCDF3_CLASSIC')
  89. if l_nc_is_new:
  90. jrec2write = 0
  91. # Creating Dimensions:
  92. f_out.createDimension('time', None)
  93. f_out.createDimension('deptht', nk)
  94. # Creating variables:
  95. id_t = f_out.createVariable('time','f4',('time',)) ; id_t.units = 'year'
  96. id_z = f_out.createVariable('deptht','f4',('deptht',)) ; id_z.units = 'm'
  97. #lili
  98. id_v01 = f_out.createVariable('sigma0', 'f4',('time','deptht',))
  99. id_v01.unit = ''; id_v01.long_name = 'sigma0 density averaged on box '+cbox
  100. id_v02 = f_out.createVariable('theta','f4',('time','deptht',))
  101. id_v02.unit = 'deg.C'; id_v02.long_name = 'potential temperature on box '+cbox
  102. id_v03 = f_out.createVariable('S','f4',('time','deptht',))
  103. id_v03.unit = 'PSU'; id_v03.long_name = 'salinity on box '+cbox
  104. id_z[:] = Vdepth[:]
  105. for jm in range(Nt):
  106. id_t[jrec2write+jm] = Vtime[jm]
  107. id_v01[jrec2write+jm,:] = Zm1[jm,:]
  108. id_v02[jrec2write+jm,:] = Tm1[jm,:]
  109. id_v03[jrec2write+jm,:] = Sm1[jm,:]
  110. f_out.box_coordinates = cbox+' => '+str(i1)+','+str(j1)+' -> '+str(i2-1)+','+str(j2-1)
  111. f_out.box_file = vdic['FILE_DEF_BOXES']
  112. f_out.Author = 'L. Brodeau ('+cname_script+' of Barakuda)'
  113. else:
  114. vt = f_out.variables['time']
  115. jrec2write = len(vt)
  116. v01 = f_out.variables['sigma0']
  117. v02 = f_out.variables['theta']
  118. v03 = f_out.variables['S']
  119. for jm in range(Nt):
  120. vt [jrec2write+jm] = Vtime[jm]
  121. v01[jrec2write+jm,:] = Zm1[jm,:]
  122. v02[jrec2write+jm,:] = Tm1[jm,:]
  123. v03[jrec2write+jm,:] = Sm1[jm,:]
  124. f_out.close()
  125. print cf_out+' written!\n'