convertPOPT.f 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! This file converts a POP grid.dat file to a remapping grid file
  4. ! in netCDF format.
  5. !
  6. !-----------------------------------------------------------------------
  7. !
  8. ! CVS:$Id: convertPOPT.f,v 1.4 2001/08/21 21:22:56 pwjones Exp $
  9. !
  10. ! Copyright (c) 1997, 1998 the Regents of the University of
  11. ! California.
  12. !
  13. ! Unless otherwise indicated, this software has been authored
  14. ! by an employee or employees of the University of California,
  15. ! operator of the Los Alamos National Laboratory under Contract
  16. ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
  17. ! Government has rights to use, reproduce, and distribute this
  18. ! software. The public may copy and use this software without
  19. ! charge, provided that this Notice and any statement of authorship
  20. ! are reproduced on all copies. Neither the Government nor the
  21. ! University makes any warranty, express or implied, or assumes
  22. ! any liability or responsibility for the use of this software.
  23. !
  24. !***********************************************************************
  25. program convertPOPT
  26. !-----------------------------------------------------------------------
  27. !
  28. ! This file converts a POP grid.dat file to a remapping grid file.
  29. !
  30. !-----------------------------------------------------------------------
  31. use kinds_mod
  32. use constants
  33. use iounits
  34. use netcdf_mod
  35. implicit none
  36. !-----------------------------------------------------------------------
  37. !
  38. ! variables that describe the grid
  39. ! 4/3 nx = 192, ny = 128
  40. ! 2/3 (mod) nx = 384, ny = 288
  41. ! x3p Greenland DP nx = 100, ny = 116
  42. ! x2p Greenland DP nx = 160, ny = 192
  43. ! x1p Greenland DP nx = 320, ny = 384
  44. !
  45. !-----------------------------------------------------------------------
  46. integer (kind=int_kind), parameter ::
  47. & nx = 320, ny = 384,
  48. & grid_size = nx*ny,
  49. & grid_rank = 2,
  50. & grid_corners = 4
  51. integer (kind=int_kind), dimension(2) ::
  52. & grid_dims ! size of each dimension
  53. character(char_len), parameter ::
  54. & grid_name = 'Greenland DP x1p',
  55. & grid_file_in = '/scratch/pwjones/grid.320x384.da',
  56. & grid_topo_in = '/scratch/pwjones/kmt.320x384.da',
  57. & grid_file_out = '/scratch/pwjones/Greenland_DP_x1p.nc'
  58. real (kind=dbl_kind), parameter ::
  59. & radius = 6370.0e5_dbl_kind ! radius of Earth (cm)
  60. &, area_norm = one/(radius*radius)
  61. !-----------------------------------------------------------------------
  62. !
  63. ! grid coordinates and masks
  64. !
  65. !-----------------------------------------------------------------------
  66. integer (kind=int_kind), dimension(grid_size) ::
  67. & grid_imask
  68. real (kind=dbl_kind), dimension(grid_size) ::
  69. & grid_area , ! area as computed in POP
  70. & grid_center_lat, ! lat/lon coordinates for
  71. & grid_center_lon ! each grid center in radians
  72. real (kind=dbl_kind), dimension(grid_corners,grid_size) ::
  73. & grid_corner_lat, ! lat/lon coordinates for
  74. & grid_corner_lon ! each grid corner in radians
  75. real (kind=dbl_kind), dimension(nx,ny) ::
  76. & HTN, HTE ! T-cell grid lengths
  77. !-----------------------------------------------------------------------
  78. !
  79. ! other local variables
  80. !
  81. !-----------------------------------------------------------------------
  82. integer (kind=int_kind) :: i, j, n, iunit, ocn_add, im1, jm1
  83. integer (kind=int_kind) ::
  84. & ncstat, ! general netCDF status variable
  85. & nc_grid_id, ! netCDF grid dataset id
  86. & nc_gridsize_id, ! netCDF grid size dim id
  87. & nc_gridcorn_id, ! netCDF grid corner dim id
  88. & nc_gridrank_id, ! netCDF grid rank dim id
  89. & nc_griddims_id, ! netCDF grid dimensions id
  90. & nc_grdcntrlat_id, ! netCDF grid center lat id
  91. & nc_grdcntrlon_id, ! netCDF grid center lon id
  92. & nc_grdimask_id, ! netCDF grid mask id
  93. & nc_gridarea_id, ! netCDF grid area id
  94. & nc_grdcrnrlat_id, ! netCDF grid corner lat id
  95. & nc_grdcrnrlon_id ! netCDF grid corner lon id
  96. integer (kind=int_kind), dimension(2) ::
  97. & nc_dims2_id ! netCDF dim id array for 2-d arrays
  98. real (kind=dbl_kind) :: tmplon, dxt, dyt
  99. !-----------------------------------------------------------------------
  100. !
  101. ! read in grid info
  102. ! lat/lon info is on velocity points which correspond
  103. ! to the NE corner (in logical space) of the grid cell.
  104. !
  105. !-----------------------------------------------------------------------
  106. call get_unit(iunit)
  107. open(unit=iunit, file=grid_topo_in, status='old',
  108. & form='unformatted', access='direct', recl=grid_size*4)
  109. read (unit=iunit,rec=1) grid_imask
  110. call release_unit(iunit)
  111. call get_unit(iunit)
  112. open(unit=iunit, file=grid_file_in, status='old',
  113. & form='unformatted', access='direct', recl=grid_size*8)
  114. read (unit=iunit, rec=1) grid_corner_lat(3,:)
  115. read (unit=iunit, rec=2) grid_corner_lon(3,:)
  116. read (unit=iunit, rec=3) HTN
  117. read (unit=iunit, rec=4) HTE
  118. call release_unit(iunit)
  119. grid_dims(1) = nx
  120. grid_dims(2) = ny
  121. !-----------------------------------------------------------------------
  122. !
  123. ! convert KMT field to integer grid mask
  124. !
  125. !-----------------------------------------------------------------------
  126. grid_imask = min(grid_imask, 1)
  127. !-----------------------------------------------------------------------
  128. !
  129. ! compute remaining corners
  130. !
  131. !-----------------------------------------------------------------------
  132. do j=1,ny
  133. do i=1,nx
  134. ocn_add = (j-1)*nx + i
  135. if (i .ne. 1) then
  136. im1 = ocn_add - 1
  137. else
  138. im1 = ocn_add + nx - 1
  139. endif
  140. grid_corner_lat(4,ocn_add) = grid_corner_lat(3,im1)
  141. grid_corner_lon(4,ocn_add) = grid_corner_lon(3,im1)
  142. end do
  143. end do
  144. do j=2,ny
  145. do i=1,nx
  146. ocn_add = (j-1)*nx + i
  147. jm1 = (j-2)*nx + i
  148. grid_corner_lat(2,ocn_add) = grid_corner_lat(3,jm1)
  149. grid_corner_lat(1,ocn_add) = grid_corner_lat(4,jm1)
  150. grid_corner_lon(2,ocn_add) = grid_corner_lon(3,jm1)
  151. grid_corner_lon(1,ocn_add) = grid_corner_lon(4,jm1)
  152. end do
  153. end do
  154. !-----------------------------------------------------------------------
  155. !
  156. ! mock up the lower row boundaries
  157. !
  158. !-----------------------------------------------------------------------
  159. do i=1,nx
  160. grid_corner_lat(1,i) = -pih + tiny
  161. grid_corner_lat(2,i) = -pih + tiny
  162. grid_corner_lon(1,i) = grid_corner_lon(4,i)
  163. grid_corner_lon(2,i) = grid_corner_lon(3,i)
  164. end do
  165. !-----------------------------------------------------------------------
  166. !
  167. ! correct for 0,2pi longitude crossings
  168. !
  169. !-----------------------------------------------------------------------
  170. do ocn_add=1,grid_size
  171. if (grid_corner_lon(1,ocn_add) > pi2)
  172. & grid_corner_lon(1,ocn_add) =
  173. & grid_corner_lon(1,ocn_add) - pi2
  174. if (grid_corner_lon(1,ocn_add) < 0.0)
  175. & grid_corner_lon(1,ocn_add) =
  176. & grid_corner_lon(1,ocn_add) + pi2
  177. do n=2,grid_corners
  178. tmplon = grid_corner_lon(n ,ocn_add) -
  179. & grid_corner_lon(n-1,ocn_add)
  180. if (tmplon < -three*pih) grid_corner_lon(n,ocn_add) =
  181. & grid_corner_lon(n,ocn_add) + pi2
  182. if (tmplon > three*pih) grid_corner_lon(n,ocn_add) =
  183. & grid_corner_lon(n,ocn_add) - pi2
  184. end do
  185. end do
  186. !-----------------------------------------------------------------------
  187. !
  188. ! compute ocean cell centers by averaging corner values
  189. !
  190. !-----------------------------------------------------------------------
  191. do ocn_add=1,grid_size
  192. grid_center_lat(ocn_add) = grid_corner_lat(1,ocn_add)
  193. grid_center_lon(ocn_add) = grid_corner_lon(1,ocn_add)
  194. do n=2,grid_corners
  195. grid_center_lat(ocn_add) = grid_center_lat(ocn_add) +
  196. & grid_corner_lat(n,ocn_add)
  197. grid_center_lon(ocn_add) = grid_center_lon(ocn_add) +
  198. & grid_corner_lon(n,ocn_add)
  199. end do
  200. grid_center_lat(ocn_add) = grid_center_lat(ocn_add)/
  201. & float(grid_corners)
  202. grid_center_lon(ocn_add) = grid_center_lon(ocn_add)/
  203. & float(grid_corners)
  204. if (grid_center_lon(ocn_add) > pi2)
  205. & grid_center_lon(ocn_add) = grid_center_lon(ocn_add) - pi2
  206. if (grid_center_lon(ocn_add) < 0.0)
  207. & grid_center_lon(ocn_add) = grid_center_lon(ocn_add) + pi2
  208. end do
  209. !-----------------------------------------------------------------------
  210. !
  211. ! compute cell areas in same way as POP
  212. !
  213. !-----------------------------------------------------------------------
  214. n = 0
  215. do j=1,ny
  216. if (j > 1) then
  217. jm1 = j-1
  218. else
  219. jm1 = 1
  220. endif
  221. do i=1,nx
  222. if (i > 1) then
  223. im1 = i-1
  224. else
  225. im1 = nx
  226. endif
  227. n = n+1
  228. dxt = half*(HTN(i,j) + HTN(i,jm1))
  229. dyt = half*(HTE(i,j) + HTE(im1,j))
  230. if (dxt == zero) dxt=one
  231. if (dyt == zero) dyt=one
  232. grid_area(n) = dxt*dyt*area_norm
  233. end do
  234. end do
  235. !-----------------------------------------------------------------------
  236. !
  237. ! set up attributes for netCDF file
  238. !
  239. !-----------------------------------------------------------------------
  240. !***
  241. !*** create netCDF dataset for this grid
  242. !***
  243. ncstat = nf_create (grid_file_out, NF_CLOBBER,
  244. & nc_grid_id)
  245. call netcdf_error_handler(ncstat)
  246. ncstat = nf_put_att_text (nc_grid_id, NF_GLOBAL, 'title',
  247. & len_trim(grid_name), grid_name)
  248. call netcdf_error_handler(ncstat)
  249. !***
  250. !*** define grid size dimension
  251. !***
  252. ncstat = nf_def_dim (nc_grid_id, 'grid_size', grid_size,
  253. & nc_gridsize_id)
  254. call netcdf_error_handler(ncstat)
  255. !***
  256. !*** define grid rank dimension
  257. !***
  258. ncstat = nf_def_dim (nc_grid_id, 'grid_rank', grid_rank,
  259. & nc_gridrank_id)
  260. call netcdf_error_handler(ncstat)
  261. !***
  262. !*** define grid corner dimension
  263. !***
  264. ncstat = nf_def_dim (nc_grid_id, 'grid_corners', grid_corners,
  265. & nc_gridcorn_id)
  266. call netcdf_error_handler(ncstat)
  267. !***
  268. !*** define grid dim size array
  269. !***
  270. ncstat = nf_def_var (nc_grid_id, 'grid_dims', NF_INT,
  271. & 1, nc_gridrank_id, nc_griddims_id)
  272. call netcdf_error_handler(ncstat)
  273. !***
  274. !*** define grid center latitude array
  275. !***
  276. ncstat = nf_def_var (nc_grid_id, 'grid_center_lat', NF_DOUBLE,
  277. & 1, nc_gridsize_id, nc_grdcntrlat_id)
  278. call netcdf_error_handler(ncstat)
  279. ncstat = nf_put_att_text (nc_grid_id, nc_grdcntrlat_id, 'units',
  280. & 7, 'radians')
  281. call netcdf_error_handler(ncstat)
  282. !***
  283. !*** define grid center longitude array
  284. !***
  285. ncstat = nf_def_var (nc_grid_id, 'grid_center_lon', NF_DOUBLE,
  286. & 1, nc_gridsize_id, nc_grdcntrlon_id)
  287. call netcdf_error_handler(ncstat)
  288. ncstat = nf_put_att_text (nc_grid_id, nc_grdcntrlon_id, 'units',
  289. & 7, 'radians')
  290. call netcdf_error_handler(ncstat)
  291. !***
  292. !*** define grid area array
  293. !***
  294. ncstat = nf_def_var (nc_grid_id, 'grid_area', NF_DOUBLE,
  295. & 1, nc_gridsize_id, nc_gridarea_id)
  296. call netcdf_error_handler(ncstat)
  297. !***
  298. !*** define grid mask
  299. !***
  300. ncstat = nf_def_var (nc_grid_id, 'grid_imask', NF_INT,
  301. & 1, nc_gridsize_id, nc_grdimask_id)
  302. call netcdf_error_handler(ncstat)
  303. ncstat = nf_put_att_text (nc_grid_id, nc_grdimask_id, 'units',
  304. & 8, 'unitless')
  305. call netcdf_error_handler(ncstat)
  306. !***
  307. !*** define grid corner latitude array
  308. !***
  309. nc_dims2_id(1) = nc_gridcorn_id
  310. nc_dims2_id(2) = nc_gridsize_id
  311. ncstat = nf_def_var (nc_grid_id, 'grid_corner_lat', NF_DOUBLE,
  312. & 2, nc_dims2_id, nc_grdcrnrlat_id)
  313. call netcdf_error_handler(ncstat)
  314. ncstat = nf_put_att_text (nc_grid_id, nc_grdcrnrlat_id, 'units',
  315. & 7, 'radians')
  316. call netcdf_error_handler(ncstat)
  317. !***
  318. !*** define grid corner longitude array
  319. !***
  320. ncstat = nf_def_var (nc_grid_id, 'grid_corner_lon', NF_DOUBLE,
  321. & 2, nc_dims2_id, nc_grdcrnrlon_id)
  322. call netcdf_error_handler(ncstat)
  323. ncstat = nf_put_att_text (nc_grid_id, nc_grdcrnrlon_id, 'units',
  324. & 7, 'radians')
  325. call netcdf_error_handler(ncstat)
  326. !***
  327. !*** end definition stage
  328. !***
  329. ncstat = nf_enddef(nc_grid_id)
  330. call netcdf_error_handler(ncstat)
  331. !-----------------------------------------------------------------------
  332. !
  333. ! write grid data
  334. !
  335. !-----------------------------------------------------------------------
  336. ncstat = nf_put_var_int(nc_grid_id, nc_griddims_id, grid_dims)
  337. call netcdf_error_handler(ncstat)
  338. ncstat = nf_put_var_double(nc_grid_id, nc_grdcntrlat_id,
  339. & grid_center_lat)
  340. ncstat = nf_put_var_int(nc_grid_id, nc_grdimask_id, grid_imask)
  341. call netcdf_error_handler(ncstat)
  342. ncstat = nf_put_var_double(nc_grid_id, nc_gridarea_id,
  343. & grid_area)
  344. call netcdf_error_handler(ncstat)
  345. ncstat = nf_put_var_double(nc_grid_id, nc_grdcntrlat_id,
  346. & grid_center_lat)
  347. call netcdf_error_handler(ncstat)
  348. ncstat = nf_put_var_double(nc_grid_id, nc_grdcntrlon_id,
  349. & grid_center_lon)
  350. call netcdf_error_handler(ncstat)
  351. ncstat = nf_put_var_double(nc_grid_id, nc_grdcrnrlat_id,
  352. & grid_corner_lat)
  353. call netcdf_error_handler(ncstat)
  354. ncstat = nf_put_var_double(nc_grid_id, nc_grdcrnrlon_id,
  355. & grid_corner_lon)
  356. call netcdf_error_handler(ncstat)
  357. ncstat = nf_close(nc_grid_id)
  358. !***********************************************************************
  359. end program convertPOPT
  360. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!