scripshape.F90 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419
  1. PROGRAM scripshape
  2. !
  3. ! program to take output from the SCRIP weights generator
  4. ! and rearrange the data into a series of 2D fields suitable
  5. ! for reading with iom_get in NEMO configurations using the
  6. ! interpolation on the fly option
  7. !
  8. USE netcdf
  9. IMPLICIT none
  10. INTEGER :: ncId, VarId, status
  11. INTEGER :: start(4), count(4)
  12. CHARACTER(LEN=1) :: y
  13. INTEGER :: nd, ns, nl, nw, sx, sy, dx, dy
  14. INTEGER :: i, j, k, m, n, smax
  15. !
  16. ! ifort -O2 -o scripshape scripshape.f90 \
  17. ! -I/nerc/packages/netcdfifort/v3.6.0-pl1/include \
  18. ! -L/nerc/packages/netcdfifort/v3.6.0-pl1/lib -lnetcdf
  19. !
  20. !
  21. INTEGER(KIND=4), ALLOCATABLE :: src(:)
  22. INTEGER(KIND=4), ALLOCATABLE :: dst(:)
  23. REAL(KIND=8), ALLOCATABLE :: wgt(:,:)
  24. REAL(KIND=8), ALLOCATABLE :: src1(:,:),dst1(:,:),wgt1(:,:)
  25. LOGICAL :: around, verbose
  26. #if defined ARGC
  27. INTEGER(KIND=4) :: iargc
  28. EXTERNAL :: iargc
  29. #endif
  30. CHARACTER(LEN=256) :: interp_file, output_file, name_file
  31. INTEGER :: ew_wrap
  32. NAMELIST /shape_inputs/ interp_file, output_file, ew_wrap
  33. ! scripshape requires 1 arguments; the name of the file containing
  34. ! the input namelist.
  35. ! This namelist contains:
  36. ! the name of the input file containing the weights ! produced by SCRIP in its format;
  37. ! the name of the new output file which ! is to contain the reorganized fields ready for input to NEMO.
  38. ! the east-west wrapping of the input grid (-1, 0, 1 and 2 are accepted values)
  39. !
  40. ! E.g.
  41. ! interp_file = 'data_nemo_bilin.nc'
  42. ! output_file = 'weights_bilin.nc'
  43. ! ew_wrap = 2
  44. !
  45. #if defined ARGC
  46. IF (iargc() == 1) THEN
  47. CALL getarg(1, name_file)
  48. ELSE
  49. WRITE(*,*) 'Usage: scripshape namelist_file'
  50. STOP
  51. ENDIF
  52. #else
  53. WRITE(6,*) 'enter name of namelist file'
  54. READ(5,*) name_file
  55. #endif
  56. interp_file = 'none'
  57. output_file = 'none'
  58. ew_wrap = 0
  59. OPEN(12, FILE=name_file, STATUS='OLD', FORM='FORMATTED')
  60. READ(12, NML=shape_inputs)
  61. CLOSE(12)
  62. !
  63. INQUIRE(FILE = TRIM(interp_file), EXIST=around)
  64. IF (.not.around) THEN
  65. WRITE(*,*) 'Input file: '//TRIM(interp_file)//' not found'
  66. STOP
  67. ENDIF
  68. !
  69. INQUIRE(FILE = TRIM(output_file), EXIST=around)
  70. IF (around) THEN
  71. WRITE(*,*) 'Output file: '//TRIM(output_file)//' exists'
  72. WRITE(*,*) 'Ok to overwrite (y/n)?'
  73. READ(5,'(a)') y
  74. IF ( y .ne. 'y' .AND. y .ne. 'Y' ) STOP
  75. ENDIF
  76. !
  77. verbose = .true.
  78. !
  79. ! Obtain grid size information from interp_file
  80. !
  81. CALL ncgetsize
  82. !
  83. ! Allocate array spaces
  84. !
  85. ALLOCATE(src(nl), STAT=status)
  86. IF(status /= 0 ) CALL alloc_err('src')
  87. ALLOCATE(dst(nl), STAT=status)
  88. IF(status /= 0 ) CALL alloc_err('dst')
  89. ALLOCATE(wgt(nw,nl), STAT=status)
  90. IF(status /= 0 ) CALL alloc_err('wgt')
  91. ALLOCATE(src1(dx,dy), STAT=status)
  92. IF(status /= 0 ) CALL alloc_err('src1')
  93. ALLOCATE(dst1(dx,dy), STAT=status)
  94. IF(status /= 0 ) CALL alloc_err('dst1')
  95. ALLOCATE(wgt1(dx,dy), STAT=status)
  96. IF(status /= 0 ) CALL alloc_err('wgt1')
  97. !
  98. ! Read all required data from interp_file
  99. !
  100. CALL ncgetfields
  101. !
  102. ! Check that dst is monotonically increasing
  103. !
  104. DO k = 1,nl-1
  105. IF(dst(k+1).lt.dst(k)) THEN
  106. WRITE(*,*) 'non-monotonic at ',k
  107. WRITE(*,*) dst(k-4:k+16)
  108. STOP
  109. ENDIF
  110. END DO
  111. !
  112. ! Remove references to the top row of src
  113. !
  114. IF(verbose) WRITE(*,*) &
  115. 'Removing references to the top row of the source grid'
  116. smax = (sy-1)*sx
  117. n = 0
  118. DO k = 1,nl
  119. IF(src(k).gt.smax-1) THEN
  120. src(k) = src(k)-sx
  121. n = n + 1
  122. ENDIF
  123. END DO
  124. IF(verbose) WRITE(*,*) n,' values changed (',100.*n/nl,'%)'
  125. !
  126. ! Loop through weights for each corner in turn and
  127. ! rearrange weight fields into separate 2D fields for
  128. ! reading with iom_get in NEMO
  129. !
  130. DO k = 1,nw
  131. DO n = 1,4
  132. i = 0
  133. j = 1
  134. DO m = n,nl,4
  135. i = i+1
  136. IF(i.gt.dx) THEN
  137. i = 1
  138. j = j + 1
  139. ENDIF
  140. src1(i,j) = src(m)
  141. dst1(i,j) = dst(m)
  142. wgt1(i,j) = wgt(k,m)
  143. END DO
  144. !
  145. ! Write out this set which will be labelled with
  146. ! a 2 digit number equal to n+4*(k-1)
  147. !
  148. CALL wrwgts
  149. !
  150. END DO
  151. END DO
  152. STOP
  153. CONTAINS
  154. !
  155. !----------------------------------------------------------------------*
  156. SUBROUTINE ncgetsize
  157. !
  158. ! Access grid size information in interp_file and set the
  159. ! following integers:
  160. !
  161. ! nd = dst_grid_size
  162. ! ns = src_grid_size
  163. ! nl = num_links
  164. ! nw = num_wgts
  165. ! sx,sy = src_grid_dims
  166. ! dx,dy = dst_grid_dims
  167. !
  168. INTEGER idims(2)
  169. !
  170. status = nf90_open(interp_file, nf90_NoWrite, ncid)
  171. IF(status /= nf90_NoErr) CALL handle_err(status)
  172. !
  173. status = nf90_inq_dimid(ncid, 'dst_grid_size', VarId)
  174. IF(status /= nf90_NoErr) CALL handle_err(status)
  175. status = nf90_inquire_dimension(ncid, VarId, LEN = nd)
  176. IF(status /= nf90_NoErr) CALL handle_err(status)
  177. !
  178. status = nf90_inq_dimid(ncid, 'src_grid_size', VarId)
  179. IF(status /= nf90_NoErr) CALL handle_err(status)
  180. status = nf90_inquire_dimension(ncid, VarId, LEN = ns)
  181. IF(status /= nf90_NoErr) CALL handle_err(status)
  182. !
  183. status = nf90_inq_dimid(ncid, 'num_links', VarId)
  184. IF(status /= nf90_NoErr) CALL handle_err(status)
  185. status = nf90_inquire_dimension(ncid, VarId, LEN = nl)
  186. IF(status /= nf90_NoErr) CALL handle_err(status)
  187. !
  188. status = nf90_inq_dimid(ncid, 'num_wgts', VarId)
  189. IF(status /= nf90_NoErr) CALL handle_err(status)
  190. status = nf90_inquire_dimension(ncid, VarId, LEN = nw)
  191. IF(status /= nf90_NoErr) CALL handle_err(status)
  192. !
  193. start = 1
  194. count = 2
  195. status = nf90_inq_varid(ncid, 'src_grid_dims', VarId)
  196. IF(status /= nf90_NoErr) CALL handle_err(status)
  197. status = nf90_get_var(ncid, VarId, idims, start, count)
  198. IF(status /= nf90_NoErr) CALL handle_err(status)
  199. sx = idims(1) ; sy = idims(2)
  200. !
  201. status = nf90_inq_varid(ncid, 'dst_grid_dims', VarId)
  202. IF(status /= nf90_NoErr) CALL handle_err(status)
  203. status = nf90_get_var(ncid, VarId, idims, start, count)
  204. IF(status /= nf90_NoErr) CALL handle_err(status)
  205. dx = idims(1) ; dy = idims(2)
  206. !
  207. status = nf90_close(ncid)
  208. IF (status /= nf90_noerr) CALL handle_err(status)
  209. !
  210. IF(verbose) THEN
  211. WRITE(*,*) 'Detected sizes: '
  212. WRITE(*,*) 'dst_grid_size: ', nd
  213. WRITE(*,*) 'src_grid_size: ', ns
  214. WRITE(*,*) 'num_links : ', nl
  215. WRITE(*,*) 'num_wgts : ', nw
  216. WRITE(*,*) 'src_grid_dims: ', sx, ' x ', sy
  217. WRITE(*,*) 'dst_grid_dims: ', dx, ' x ', dy
  218. ENDIF
  219. !
  220. END SUBROUTINE ncgetsize
  221. !----------------------------------------------------------------------*
  222. SUBROUTINE ncgetfields
  223. !
  224. ! Read all required data from interp_file. The data read are:
  225. !
  226. ! netcdf variable size internal array
  227. !-----------------+-------+--------------
  228. ! src_address nl src
  229. ! dst_address nl dst
  230. ! remap_matrix (nw,nl) wgt
  231. !
  232. status = nf90_open(interp_file, nf90_NoWrite, ncid)
  233. IF(status /= nf90_NoErr) CALL handle_err(status)
  234. !
  235. status = nf90_inq_varid(ncid, 'src_address', VarId)
  236. IF(status /= nf90_NoErr) CALL handle_err(status)
  237. !
  238. ! Read the values for src
  239. status = nf90_get_var(ncid, VarId, src, &
  240. start = (/ 1 /), &
  241. count = (/ nl /))
  242. IF(status /= nf90_NoErr) CALL handle_err(status)
  243. !
  244. status = nf90_inq_varid(ncid, 'dst_address', VarId)
  245. IF(status /= nf90_NoErr) CALL handle_err(status)
  246. !
  247. ! Read the values for dst
  248. status = nf90_get_var(ncid, VarId, dst, &
  249. start = (/ 1 /), &
  250. count = (/ nl /))
  251. IF(status /= nf90_NoErr) CALL handle_err(status)
  252. !
  253. status = nf90_inq_varid(ncid, 'remap_matrix', VarId)
  254. IF(status /= nf90_NoErr) CALL handle_err(status)
  255. !
  256. ! Read the values for wgt
  257. status = nf90_get_var(ncid, VarId, wgt, &
  258. start = (/ 1, 1 /), &
  259. count = (/ nw, nl /))
  260. IF(status /= nf90_NoErr) CALL handle_err(status)
  261. !
  262. status = nf90_close(ncid)
  263. IF (status /= nf90_noerr) CALL handle_err(status)
  264. !
  265. END SUBROUTINE ncgetfields
  266. !----------------------------------------------------------------------*
  267. SUBROUTINE handle_err(status)
  268. !
  269. ! Simple netcdf error checking routine
  270. !
  271. INTEGER, intent ( in) :: status
  272. !
  273. IF(status /= nf90_noerr) THEN
  274. IF(trim(nf90_strerror(status)) .eq. 'Attribute not found') THEN
  275. ! ignore
  276. ELSE
  277. WRITE(*,*) trim(nf90_strerror(status))
  278. STOP "Stopped"
  279. END IF
  280. END IF
  281. END SUBROUTINE handle_err
  282. !----------------------------------------------------------------------*
  283. SUBROUTINE alloc_err(arname)
  284. !
  285. ! Simple allocation error checking routine
  286. !
  287. CHARACTER(LEN=*) :: arname
  288. !
  289. WRITE(*,*) 'Allocation error attempting to ALLOCATE '//arname
  290. STOP "Stopped"
  291. END SUBROUTINE alloc_err
  292. !
  293. !----------------------------------------------------------------------*
  294. SUBROUTINE wrwgts
  295. !
  296. ! Write out each set of 2D fields to output_file.
  297. ! Each call will write out a set of srcXX, dstXX and wgtXX fields
  298. ! where XX is a two digit number equal to n + 4*(k-1). The first
  299. ! and last calls to this routine initialise and close the output
  300. ! file respectively. The first call is detected when k*n=1 and the
  301. ! last call is detected when k*n=4*nw. The outfile file remains
  302. ! open between the first and last calls.
  303. !
  304. INTEGER :: status, ncid, ncin
  305. INTEGER :: Lontdid, Lattdid
  306. INTEGER :: tvid, tvid2, tvid3
  307. INTEGER :: ioldfill
  308. CHARACTER(LEN=2) :: cs
  309. SAVE ncid, Lontdid, Lattdid
  310. !
  311. IF(k*n.eq.1) THEN
  312. !
  313. ! Create output_file and set the dimensions
  314. !
  315. status = nf90_create(output_file, nf90_Clobber, ncid)
  316. IF(status /= nf90_NoErr) CALL handle_err(status)
  317. status = nf90_set_fill(ncid, nf90_NoFill, ioldfill)
  318. IF(status /= nf90_NoErr) CALL handle_err(status)
  319. !
  320. status = nf90_def_dim(ncid, "lon", dx, Lontdid)
  321. IF(status /= nf90_NoErr) CALL handle_err(status)
  322. status = nf90_def_dim(ncid, "lat", dy, Lattdid)
  323. IF(status /= nf90_NoErr) CALL handle_err(status)
  324. !
  325. status = nf90_put_att(ncid, nf90_global, 'ew_wrap', ew_wrap)
  326. IF(status /= nf90_NoErr) CALL handle_err(status)
  327. ELSE
  328. !
  329. ! Reenter define mode
  330. !
  331. status = nf90_redef(ncid)
  332. IF(status /= nf90_NoErr) CALL handle_err(status)
  333. ENDIF
  334. !
  335. WRITE(cs,'(i2.2)') n + 4*(k-1)
  336. !
  337. ! Define new variables
  338. !
  339. status = nf90_def_var(ncid, "src"//cs, nf90_double, &
  340. (/ Lontdid, Lattdid /), tvid)
  341. IF(status /= nf90_NoErr) CALL handle_err(status)
  342. !
  343. status = nf90_def_var(ncid, "dst"//cs, nf90_double, &
  344. (/ Lontdid, Lattdid /), tvid2)
  345. IF(status /= nf90_NoErr) CALL handle_err(status)
  346. !
  347. status = nf90_def_var(ncid, "wgt"//cs, nf90_double, &
  348. (/ Lontdid, Lattdid /), tvid3)
  349. IF(status /= nf90_NoErr) CALL handle_err(status)
  350. !
  351. ! Leave define mode
  352. !
  353. status = nf90_enddef(ncid)
  354. IF(status /= nf90_NoErr) CALL handle_err(status)
  355. !
  356. ! Write the data
  357. !
  358. status = nf90_put_var(ncid, tvid, src1, &
  359. start = (/ 1, 1 /), &
  360. count = (/ dx, dy /) )
  361. IF(status /= nf90_NoErr) CALL handle_err(status)
  362. !
  363. status = nf90_put_var(ncid, tvid2, dst1, &
  364. start = (/ 1, 1 /), &
  365. count = (/ dx, dy /) )
  366. IF(status /= nf90_NoErr) CALL handle_err(status)
  367. !
  368. status = nf90_put_var(ncid, tvid3, wgt1, &
  369. start = (/ 1, 1 /), &
  370. count = (/ dx, dy /) )
  371. IF(status /= nf90_NoErr) CALL handle_err(status)
  372. !
  373. IF(k*n.eq.4*nw) THEN
  374. !
  375. ! -- Reenter define mode
  376. !
  377. status = nf90_redef(ncid)
  378. IF(status /= nf90_NoErr) CALL handle_err(status)
  379. !
  380. ! -- Reopen interp_file and transfer some global attributes
  381. !
  382. status = nf90_open(interp_file, nf90_NoWrite, ncin)
  383. IF(status /= nf90_NoErr) CALL handle_err(status)
  384. !
  385. status = nf90_copy_att(ncin,NF90_GLOBAL,'title', ncid,NF90_GLOBAL)
  386. IF(status /= nf90_NoErr) CALL handle_err(status)
  387. !
  388. status = nf90_copy_att(ncin,NF90_GLOBAL,'normalization',ncid,NF90_GLOBAL)
  389. IF(status /= nf90_NoErr) CALL handle_err(status)
  390. !
  391. status = nf90_copy_att(ncin,NF90_GLOBAL,'map_method', ncid,NF90_GLOBAL)
  392. IF(status /= nf90_NoErr) CALL handle_err(status)
  393. !
  394. status = nf90_copy_att(ncin,NF90_GLOBAL,'conventions', ncid,NF90_GLOBAL)
  395. IF(status /= nf90_NoErr) CALL handle_err(status)
  396. !
  397. status = nf90_copy_att(ncin,NF90_GLOBAL,'history', ncid,NF90_GLOBAL)
  398. IF(status /= nf90_NoErr) CALL handle_err(status)
  399. !
  400. ! -- Close interp_file
  401. !
  402. status = nf90_close(ncin)
  403. IF(status /= nf90_NoErr) CALL handle_err(status)
  404. !
  405. ! -- Close output_file
  406. !
  407. status = nf90_close(ncid)
  408. IF(status /= nf90_NoErr) CALL handle_err(status)
  409. ENDIF
  410. END SUBROUTINE wrwgts
  411. END PROGRAM scripshape