convertgauss.f 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! This program creates a remapping grid file for Gaussian lat/lon
  4. ! grids (for spectral transform codes).
  5. !
  6. !-----------------------------------------------------------------------
  7. !
  8. ! CVS:$Id: convertgauss.f,v 1.3 2000/04/19 22:05:57 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 convert_gauss
  26. !-----------------------------------------------------------------------
  27. !
  28. ! This file creates a remapping grid file for a Gaussian grid
  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. !
  40. ! T42: nx=128 ny=64
  41. ! T62: nx=192 ny=94
  42. !
  43. !-----------------------------------------------------------------------
  44. integer (kind=int_kind), parameter ::
  45. & nx = 192, ny = 94,
  46. & grid_size = nx*ny,
  47. & grid_rank = 2,
  48. & grid_corners = 4
  49. character(char_len), parameter ::
  50. & grid_name = 'T62 Gaussian Grid',
  51. & grid_file_out = 'remap_grid_T62.nc'
  52. integer (kind=int_kind), dimension(grid_rank) ::
  53. & grid_dims
  54. !-----------------------------------------------------------------------
  55. !
  56. ! grid coordinates and masks
  57. !
  58. !-----------------------------------------------------------------------
  59. integer (kind=int_kind), dimension(grid_size) ::
  60. & grid_imask
  61. real (kind=dbl_kind), dimension(grid_size) ::
  62. & grid_center_lat, ! lat/lon coordinates for
  63. & grid_center_lon ! each grid center in degrees
  64. real (kind=dbl_kind), dimension(grid_corners,grid_size) ::
  65. & grid_corner_lat, ! lat/lon coordinates for
  66. & grid_corner_lon ! each grid corner in degrees
  67. !-----------------------------------------------------------------------
  68. !
  69. ! other local variables
  70. !
  71. !-----------------------------------------------------------------------
  72. integer (kind=int_kind) :: i, j, iunit, atm_add
  73. integer (kind=int_kind) ::
  74. & ncstat, ! general netCDF status variable
  75. & nc_grid_id, ! netCDF grid dataset id
  76. & nc_gridsize_id, ! netCDF grid size dim id
  77. & nc_gridcorn_id, ! netCDF grid corner dim id
  78. & nc_gridrank_id, ! netCDF grid rank dim id
  79. & nc_griddims_id, ! netCDF grid dimension size id
  80. & nc_grdcntrlat_id, ! netCDF grid center lat id
  81. & nc_grdcntrlon_id, ! netCDF grid center lon id
  82. & nc_grdimask_id, ! netCDF grid mask id
  83. & nc_grdcrnrlat_id, ! netCDF grid corner lat id
  84. & nc_grdcrnrlon_id ! netCDF grid corner lon id
  85. integer (kind=int_kind), dimension(2) ::
  86. & nc_dims2_id ! netCDF dim id array for 2-d arrays
  87. real (kind=dbl_kind) :: dlon, minlon, maxlon, centerlon,
  88. & minlat, maxlat, centerlat
  89. real (kind=dbl_kind), dimension(ny) :: gauss_root, gauss_wgt
  90. !-----------------------------------------------------------------------
  91. !
  92. ! compute longitudes of cell centers and corners. set up alon
  93. ! array for search routine.
  94. !
  95. !-----------------------------------------------------------------------
  96. grid_dims(1) = nx
  97. grid_dims(2) = ny
  98. dlon = 360./nx
  99. do i=1,nx
  100. centerlon = (i-1)*dlon
  101. minlon = centerlon - half*dlon
  102. maxlon = centerlon + half*dlon
  103. do j=1,ny
  104. atm_add = (j-1)*nx + i
  105. grid_center_lon(atm_add ) = centerlon
  106. grid_corner_lon(1,atm_add) = minlon
  107. grid_corner_lon(2,atm_add) = maxlon
  108. grid_corner_lon(3,atm_add) = maxlon
  109. grid_corner_lon(4,atm_add) = minlon
  110. end do
  111. end do
  112. !-----------------------------------------------------------------------
  113. !
  114. ! compute Gaussian latitudes and store in gauss_wgt.
  115. !
  116. !-----------------------------------------------------------------------
  117. call gquad(ny, gauss_root, gauss_wgt)
  118. do j=1,ny
  119. gauss_wgt(j) = pih - gauss_root(ny+1-j)
  120. end do
  121. !-----------------------------------------------------------------------
  122. !
  123. ! compute latitudes at cell centers and corners. set up alat
  124. ! array for search routine.
  125. !
  126. !-----------------------------------------------------------------------
  127. do j=1,ny
  128. centerlat = gauss_wgt(j)
  129. if (j .eq. 1) then
  130. minlat = -pih
  131. else
  132. minlat = ATAN((COS(gauss_wgt(j-1)) -
  133. & COS(gauss_wgt(j )))/
  134. & (SIN(gauss_wgt(j )) -
  135. & SIN(gauss_wgt(j-1))))
  136. endif
  137. if (j .eq. ny) then
  138. maxlat = pih
  139. else
  140. maxlat = ATAN((COS(gauss_wgt(j )) -
  141. & COS(gauss_wgt(j+1)))/
  142. & (SIN(gauss_wgt(j+1)) -
  143. & SIN(gauss_wgt(j ))))
  144. endif
  145. do i=1,nx
  146. atm_add = (j-1)*nx + i
  147. grid_center_lat(atm_add ) = centerlat*360./pi2
  148. grid_corner_lat(1,atm_add) = minlat*360./pi2
  149. grid_corner_lat(2,atm_add) = minlat*360./pi2
  150. grid_corner_lat(3,atm_add) = maxlat*360./pi2
  151. grid_corner_lat(4,atm_add) = maxlat*360./pi2
  152. end do
  153. end do
  154. !-----------------------------------------------------------------------
  155. !
  156. ! define mask
  157. !
  158. !-----------------------------------------------------------------------
  159. grid_imask = 1
  160. !-----------------------------------------------------------------------
  161. !
  162. ! set up attributes for netCDF file
  163. !
  164. !-----------------------------------------------------------------------
  165. !***
  166. !*** create netCDF dataset for this grid
  167. !***
  168. ncstat = nf_create (grid_file_out, NF_CLOBBER,
  169. & nc_grid_id)
  170. call netcdf_error_handler(ncstat)
  171. ncstat = nf_put_att_text (nc_grid_id, NF_GLOBAL, 'title',
  172. & len_trim(grid_name), grid_name)
  173. call netcdf_error_handler(ncstat)
  174. !***
  175. !*** define grid size dimension
  176. !***
  177. ncstat = nf_def_dim (nc_grid_id, 'grid_size', grid_size,
  178. & nc_gridsize_id)
  179. call netcdf_error_handler(ncstat)
  180. !***
  181. !*** define grid corner dimension
  182. !***
  183. ncstat = nf_def_dim (nc_grid_id, 'grid_corners', grid_corners,
  184. & nc_gridcorn_id)
  185. call netcdf_error_handler(ncstat)
  186. !***
  187. !*** define grid rank dimension
  188. !***
  189. ncstat = nf_def_dim (nc_grid_id, 'grid_rank', grid_rank,
  190. & nc_gridrank_id)
  191. call netcdf_error_handler(ncstat)
  192. !***
  193. !*** define grid dimension size array
  194. !***
  195. ncstat = nf_def_var (nc_grid_id, 'grid_dims', NF_INT,
  196. & 1, nc_gridrank_id, nc_griddims_id)
  197. call netcdf_error_handler(ncstat)
  198. !***
  199. !*** define grid center latitude array
  200. !***
  201. ncstat = nf_def_var (nc_grid_id, 'grid_center_lat', NF_DOUBLE,
  202. & 1, nc_gridsize_id, nc_grdcntrlat_id)
  203. call netcdf_error_handler(ncstat)
  204. ncstat = nf_put_att_text (nc_grid_id, nc_grdcntrlat_id, 'units',
  205. & 7, 'degrees')
  206. call netcdf_error_handler(ncstat)
  207. !***
  208. !*** define grid center longitude array
  209. !***
  210. ncstat = nf_def_var (nc_grid_id, 'grid_center_lon', NF_DOUBLE,
  211. & 1, nc_gridsize_id, nc_grdcntrlon_id)
  212. call netcdf_error_handler(ncstat)
  213. ncstat = nf_put_att_text (nc_grid_id, nc_grdcntrlon_id, 'units',
  214. & 7, 'degrees')
  215. call netcdf_error_handler(ncstat)
  216. !***
  217. !*** define grid mask
  218. !***
  219. ncstat = nf_def_var (nc_grid_id, 'grid_imask', NF_INT,
  220. & 1, nc_gridsize_id, nc_grdimask_id)
  221. call netcdf_error_handler(ncstat)
  222. ncstat = nf_put_att_text (nc_grid_id, nc_grdimask_id, 'units',
  223. & 8, 'unitless')
  224. call netcdf_error_handler(ncstat)
  225. !***
  226. !*** define grid corner latitude array
  227. !***
  228. nc_dims2_id(1) = nc_gridcorn_id
  229. nc_dims2_id(2) = nc_gridsize_id
  230. ncstat = nf_def_var (nc_grid_id, 'grid_corner_lat', NF_DOUBLE,
  231. & 2, nc_dims2_id, nc_grdcrnrlat_id)
  232. call netcdf_error_handler(ncstat)
  233. ncstat = nf_put_att_text (nc_grid_id, nc_grdcrnrlat_id, 'units',
  234. & 7, 'degrees')
  235. call netcdf_error_handler(ncstat)
  236. !***
  237. !*** define grid corner longitude array
  238. !***
  239. ncstat = nf_def_var (nc_grid_id, 'grid_corner_lon', NF_DOUBLE,
  240. & 2, nc_dims2_id, nc_grdcrnrlon_id)
  241. call netcdf_error_handler(ncstat)
  242. ncstat = nf_put_att_text (nc_grid_id, nc_grdcrnrlon_id, 'units',
  243. & 7, 'degrees')
  244. call netcdf_error_handler(ncstat)
  245. !***
  246. !*** end definition stage
  247. !***
  248. ncstat = nf_enddef(nc_grid_id)
  249. call netcdf_error_handler(ncstat)
  250. !-----------------------------------------------------------------------
  251. !
  252. ! write grid data
  253. !
  254. !-----------------------------------------------------------------------
  255. ncstat = nf_put_var_int(nc_grid_id, nc_griddims_id, grid_dims)
  256. call netcdf_error_handler(ncstat)
  257. ncstat = nf_put_var_int(nc_grid_id, nc_grdimask_id, grid_imask)
  258. call netcdf_error_handler(ncstat)
  259. ncstat = nf_put_var_double(nc_grid_id, nc_grdcntrlat_id,
  260. & grid_center_lat)
  261. call netcdf_error_handler(ncstat)
  262. ncstat = nf_put_var_double(nc_grid_id, nc_grdcntrlon_id,
  263. & grid_center_lon)
  264. call netcdf_error_handler(ncstat)
  265. ncstat = nf_put_var_double(nc_grid_id, nc_grdcrnrlat_id,
  266. & grid_corner_lat)
  267. call netcdf_error_handler(ncstat)
  268. ncstat = nf_put_var_double(nc_grid_id, nc_grdcrnrlon_id,
  269. & grid_corner_lon)
  270. call netcdf_error_handler(ncstat)
  271. ncstat = nf_close(nc_grid_id)
  272. call netcdf_error_handler(ncstat)
  273. !-----------------------------------------------------------------------
  274. end program convert_gauss
  275. !***********************************************************************
  276. subroutine gquad(l,root,w)
  277. !-----------------------------------------------------------------------
  278. !
  279. ! This subroutine finds the l roots (in theta) and gaussian weights
  280. ! associated with the legendre polynomial of degree l > 1.
  281. !
  282. !-----------------------------------------------------------------------
  283. use kinds_mod
  284. use constants
  285. implicit none
  286. !-----------------------------------------------------------------------
  287. !
  288. ! intent(in)
  289. !
  290. !-----------------------------------------------------------------------
  291. integer (kind=int_kind), intent(in) :: l
  292. !-----------------------------------------------------------------------
  293. !
  294. ! intent(out)
  295. !
  296. !-----------------------------------------------------------------------
  297. real (kind=dbl_kind), dimension(l), intent(out) ::
  298. & root, w
  299. !-----------------------------------------------------------------------
  300. !
  301. ! local
  302. !
  303. !-----------------------------------------------------------------------
  304. integer (kind=int_kind) :: l1, l2, l22, l3, k, i, j
  305. real (kind=dbl_kind) ::
  306. & del,co,p1,p2,p3,t1,t2,slope,s,c,pp1,pp2,p00
  307. !-----------------------------------------------------------------------
  308. !
  309. ! Define useful constants.
  310. !
  311. !-----------------------------------------------------------------------
  312. del= pi/float(4*l)
  313. l1 = l+1
  314. co = float(2*l+3)/float(l1**2)
  315. p2 = 1.0
  316. t2 = -del
  317. l2 = l/2
  318. k = 1
  319. p00 = one/sqrt(two)
  320. !-----------------------------------------------------------------------
  321. !
  322. ! Start search for each root by looking for crossing point.
  323. !
  324. !-----------------------------------------------------------------------
  325. do i=1,l2
  326. 10 t1 = t2
  327. t2 = t1+del
  328. p1 = p2
  329. s = sin(t2)
  330. c = cos(t2)
  331. pp1 = 1.0
  332. p3 = p00
  333. do j=1,l1
  334. pp2 = pp1
  335. pp1 = p3
  336. p3 = 2.0*sqrt((float(j**2)-0.250)/float(j**2))*c*pp1-
  337. & sqrt(float((2*j+1)*(j-1)*(j-1))/
  338. & float((2*j-3)*j*j))*pp2
  339. end do
  340. p2 = pp1
  341. if ((k*p2).gt.0) goto 10
  342. !-----------------------------------------------------------------------
  343. !
  344. ! Now converge using Newton-Raphson.
  345. !
  346. !-----------------------------------------------------------------------
  347. k = -k
  348. 20 continue
  349. slope = (t2-t1)/(p2-p1)
  350. t1 = t2
  351. t2 = t2-slope*p2
  352. p1 = p2
  353. s = sin(t2)
  354. c = cos(t2)
  355. pp1 = 1.0
  356. p3 = p00
  357. do j=1,l1
  358. pp2 = pp1
  359. pp1 = p3
  360. p3 = 2.0*sqrt((float(j**2)-0.250)/float(j**2))*c*pp1-
  361. & sqrt(float((2*j+1)*(j-1)*(j-1))/
  362. & float((2*j-3)*j*j))*pp2
  363. end do
  364. p2 = pp1
  365. if (abs(p2).gt.1.e-10) goto 20
  366. root(i) = t2
  367. w(i) = co*(sin(t2)/p3)**2
  368. end do
  369. !-----------------------------------------------------------------------
  370. !
  371. ! If l is odd, take care of odd point.
  372. !
  373. !-----------------------------------------------------------------------
  374. l22 = 2*l2
  375. if (l22 .ne. l) then
  376. l2 = l2+1
  377. t2 = pi/2.0
  378. root(l2) = t2
  379. s = sin(t2)
  380. c = cos(t2)
  381. pp1 = 1.0
  382. p3 = p00
  383. do j=1,l1
  384. pp2 = pp1
  385. pp1 = p3
  386. p3 = 2.0*sqrt((float(j**2)-0.250)/float(j**2))*c*pp1-
  387. & sqrt(float((2*j+1)*(j-1)*(j-1))/
  388. & float((2*j-3)*j*j))*pp2
  389. end do
  390. p2 = pp1
  391. w(l2) = co/p3**2
  392. endif
  393. !-----------------------------------------------------------------------
  394. !
  395. ! Use symmetry to compute remaining roots and weights.
  396. !
  397. !-----------------------------------------------------------------------
  398. l3 = l2+1
  399. do i=l3,l
  400. root(i) = pi-root(l-i+1)
  401. w(i) = w(l-i+1)
  402. end do
  403. !-----------------------------------------------------------------------
  404. end subroutine gquad
  405. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!