scrip_use.f 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612
  1. !-----------------------------------------------------------------------
  2. !
  3. ! This script interpolates a (time, latitude, longitude) netcdf file
  4. ! using the interpolation weights computed by SCRIP :
  5. ! Spherical Coordinate Remapping and Interpolation Package
  6. !
  7. ! The arguments are passed through a namelist named scrip_use_in:
  8. !&remap_inputs
  9. ! remap_wgt = 'Weights from SCRIP'
  10. ! infile = 'input netcdf file'
  11. ! invertlat = TRUE/FALSE : should the latitudes be reverted ?
  12. ! var = 'netcdf variable name'
  13. ! fromregular = TRUE/FALSE : is the input grid regular ?
  14. ! outfile = 'output netcdf file'
  15. !/
  16. !
  17. ! History : Virginie Guemas - Initial version - 2011
  18. !-----------------------------------------------------------------------
  19. program scrip_use
  20. !-----------------------------------------------------------------------
  21. use kinds_mod ! defines common data types
  22. use constants ! defines common constants
  23. use iounits ! I/O unit manager
  24. use netcdf_mod ! netcdf I/O stuff
  25. use grids ! module containing grid info
  26. use remap_vars ! module containing remapping info
  27. use remap_mod ! module containing remapping routines
  28. use remap_read ! routines for reading remap files
  29. use read_input_file ! routines to read file to interpolate
  30. implicit none
  31. !-----------------------------------------------------------------------
  32. !
  33. ! input namelist variables
  34. !
  35. !-----------------------------------------------------------------------
  36. character (char_len) ::
  37. & remap_wgt, ! filename containing remap data (map1)
  38. & infile, ! filename containing input field
  39. & var, ! var name in infile
  40. & outfile ! filename to output interpolated field
  41. logical :: fromregular
  42. namelist /remap_inputs/ remap_wgt, infile, invertlat, var,
  43. & fromregular, outfile
  44. !-----------------------------------------------------------------------
  45. !
  46. ! local variables
  47. !
  48. !-----------------------------------------------------------------------
  49. character (char_len) ::
  50. & map_name ! name for mapping from grid1 to grid2
  51. integer (kind=int_kind) :: ! netCDF ids for files and arrays
  52. & ncstat, nc_outfile_id,
  53. & nc_dstgrdcntrlat_id, nc_dstgrdcntrlon_id,
  54. & nc_dstarray1_id, nc_dstarray1a_id, nc_dstarray2_id,
  55. & nc_vartime_id, nc_srcarray_id
  56. integer (kind=int_kind), dimension(:), allocatable ::
  57. & nc_grid2size_id, nc_grid1size_id
  58. !-----------------------------------------------------------------------
  59. character (char_len) ::
  60. & dim_name, attname ! netCDF dimension name
  61. integer (kind=int_kind) :: i,j,n,imin,imax,idiff,
  62. & ip1,im1,jp1,jm1,nx,ny, ! for computing bicub gradients
  63. & in,is,ie,iw,ine,inw,ise,isw,jatt,
  64. & iunit,jtime ! unit number for namelist file
  65. integer (kind=int_kind), dimension(:), allocatable ::
  66. & grid1_imask, grid2_imask, grid2_count
  67. real (kind=dbl_kind) ::
  68. & delew, delns, ! variables for computing bicub gradients
  69. & length ! length scale for cosine hill test field
  70. real (kind=dbl_kind), dimension(:), allocatable ::
  71. & grid1_tmp,
  72. & grad1_lat,
  73. & grad1_lon,
  74. & grad1_latlon,
  75. & grad1_lat_zero,
  76. & grad1_lon_zero,
  77. & grid2_array,
  78. & grid2_err,
  79. & grid2_tmp,
  80. & grid1_mask_grid2
  81. !-----------------------------------------------------------------------
  82. !
  83. ! read namelist for file and mapping info
  84. !
  85. !-----------------------------------------------------------------------
  86. call get_unit(iunit)
  87. open(iunit, file='scrip_use_in', status='old', form='formatted')
  88. read(iunit, nml=remap_inputs)
  89. call release_unit(iunit)
  90. write(*,nml=remap_inputs)
  91. !-----------------------------------------------------------------------
  92. !
  93. ! read remapping data
  94. !
  95. !-----------------------------------------------------------------------
  96. call read_remap(map_name, remap_wgt)
  97. !-----------------------------------------------------------------------
  98. !
  99. ! read input file
  100. !
  101. !-----------------------------------------------------------------------
  102. call read_input(infile, var)
  103. !-----------------------------------------------------------------------
  104. !
  105. ! allocate arrays
  106. !
  107. !-----------------------------------------------------------------------
  108. allocate (grid1_tmp (grid1_size),
  109. & grad1_lat (grid1_size),
  110. & grad1_lon (grid1_size),
  111. & grad1_latlon (grid1_size),
  112. & grad1_lat_zero (grid1_size),
  113. & grad1_lon_zero (grid1_size),
  114. & grid1_imask (grid1_size),
  115. & grid2_array (grid2_size),
  116. & grid2_err (grid2_size),
  117. & grid2_tmp (grid2_size),
  118. & grid2_imask (grid2_size),
  119. & grid2_count (grid2_size),
  120. & grid1_mask_grid2(grid2_size))
  121. where (grid1_mask)
  122. grid1_imask = 1
  123. elsewhere
  124. grid1_imask = 0
  125. endwhere
  126. where (grid2_mask)
  127. grid2_imask = 1
  128. elsewhere
  129. grid2_imask = 0
  130. endwhere
  131. !-----------------------------------------------------------------------
  132. !
  133. ! setup a NetCDF file for output
  134. !
  135. !-----------------------------------------------------------------------
  136. !***
  137. !*** create netCDF dataset
  138. !***
  139. ncstat = nf_create (outfile, NF_CLOBBER, nc_outfile_id)
  140. call netcdf_error_handler(ncstat)
  141. ncstat = nf_put_att_text (nc_outfile_id, NF_GLOBAL, 'title',
  142. & len_trim(map_name), map_name)
  143. call netcdf_error_handler(ncstat)
  144. !***
  145. !*** define grid size dimensions
  146. !***
  147. allocate( nc_grid2size_id(grid2_rank+1))
  148. allocate( nc_grid1size_id(grid1_rank+1))
  149. ncstat = nf_def_dim (nc_outfile_id, 'x',
  150. & grid2_dims(1), nc_grid2size_id(1))
  151. call netcdf_error_handler(ncstat)
  152. ncstat = nf_def_dim (nc_outfile_id, 'y',
  153. & grid2_dims(2), nc_grid2size_id(2))
  154. call netcdf_error_handler(ncstat)
  155. ncstat = nf_def_dim (nc_outfile_id, 'inx',
  156. & grid1_dims(1), nc_grid1size_id(1))
  157. call netcdf_error_handler(ncstat)
  158. if ( grid1_dims(2) > 0 ) then
  159. ncstat = nf_def_dim (nc_outfile_id, 'iny',
  160. & grid1_dims(2), nc_grid1size_id(2))
  161. call netcdf_error_handler(ncstat)
  162. endif
  163. !***
  164. !*** Create time dimension
  165. !***
  166. ncstat = nf_def_dim (nc_outfile_id, 'time',
  167. & NF_UNLIMITED, nc_grid2size_id(3))
  168. call netcdf_error_handler(ncstat)
  169. nc_grid1size_id(grid1_rank+1)=nc_grid2size_id(3)
  170. ncstat = nf_def_var (nc_outfile_id, 'time',
  171. & NF_DOUBLE, 1, nc_grid2size_id(3), nc_vartime_id)
  172. call netcdf_error_handler(ncstat)
  173. if ( nc_time_id > -1 ) then
  174. if ( natts >= 1 ) then
  175. do jatt = 1,natts
  176. ncstat = nf_inq_attname(nc_infile_id, nc_time_id, jatt,
  177. & attname)
  178. call netcdf_error_handler(ncstat)
  179. ncstat = nf_copy_att(nc_infile_id, nc_time_id, attname,
  180. & nc_outfile_id, nc_vartime_id)
  181. enddo
  182. endif
  183. endif
  184. !***
  185. !*** define grid center latitude array
  186. !***
  187. ncstat = nf_def_var (nc_outfile_id, 'latitude',
  188. & NF_DOUBLE, grid2_rank, nc_grid2size_id
  189. & (1:grid2_rank), nc_dstgrdcntrlat_id)
  190. call netcdf_error_handler(ncstat)
  191. ncstat = nf_put_att_text (nc_outfile_id, nc_dstgrdcntrlat_id,
  192. & 'units', 7, 'degrees')
  193. call netcdf_error_handler(ncstat)
  194. !***
  195. !*** define grid center longitude array
  196. !***
  197. ncstat = nf_def_var (nc_outfile_id, 'longitude',
  198. & NF_DOUBLE, grid2_rank, nc_grid2size_id
  199. & (1:grid2_rank), nc_dstgrdcntrlon_id)
  200. call netcdf_error_handler(ncstat)
  201. ncstat = nf_put_att_text (nc_outfile_id, nc_dstgrdcntrlon_id,
  202. & 'units', 7, 'degrees')
  203. call netcdf_error_handler(ncstat)
  204. !***
  205. !*** define source array
  206. !***
  207. ncstat = nf_def_var (nc_outfile_id, 'input',
  208. & NF_DOUBLE, (grid1_rank+1), nc_grid1size_id,
  209. & nc_srcarray_id)
  210. call netcdf_error_handler(ncstat)
  211. ncstat = nf_put_att_double (nc_outfile_id, nc_srcarray_id,
  212. & 'missing_value', NF_DOUBLE,1,dble(1e20))
  213. call netcdf_error_handler(ncstat)
  214. !***
  215. !*** define destination arrays
  216. !***
  217. ncstat = nf_def_var (nc_outfile_id, var,
  218. & NF_DOUBLE, ( grid2_rank + 1 ),
  219. & nc_grid2size_id, nc_dstarray1_id )
  220. call netcdf_error_handler(ncstat)
  221. ncstat = nf_put_att_double (nc_outfile_id, nc_dstarray1_id,
  222. & 'missing_value', NF_DOUBLE,1,dble(1e20))
  223. call netcdf_error_handler(ncstat)
  224. if ( nvaratts >= 1 ) then
  225. do jatt = 1,nvaratts
  226. ncstat = nf_inq_attname(nc_infile_id, nc_invar_id, jatt,
  227. & attname)
  228. call netcdf_error_handler(ncstat)
  229. if ((attname .ne. '_FillValue') .and. (attname .ne.
  230. & 'missing_value') ) then
  231. ncstat = nf_copy_att(nc_infile_id, nc_invar_id, attname,
  232. & nc_outfile_id, nc_dstarray1_id)
  233. call netcdf_error_handler(ncstat)
  234. endif
  235. enddo
  236. endif
  237. if ( nglobatts >= 1 ) then
  238. do jatt = 1,nglobatts
  239. ncstat = nf_inq_attname(nc_infile_id, NF_GLOBAL, jatt,
  240. & attname)
  241. call netcdf_error_handler(ncstat)
  242. ncstat = nf_copy_att(nc_infile_id, NF_GLOBAL, attname,
  243. & nc_outfile_id, NF_GLOBAL)
  244. call netcdf_error_handler(ncstat)
  245. enddo
  246. endif
  247. ncstat = nf_close(nc_infile_id)
  248. call netcdf_error_handler(ncstat)
  249. !***
  250. !*** end definition stage
  251. !***
  252. ncstat = nf_enddef(nc_outfile_id)
  253. call netcdf_error_handler(ncstat)
  254. !-----------------------------------------------------------------------
  255. !
  256. ! write some grid info
  257. !
  258. !-----------------------------------------------------------------------
  259. !***
  260. !*** write grid center latitude array
  261. !***
  262. ncstat = nf_put_var_double(nc_outfile_id, nc_dstgrdcntrlat_id,
  263. & grid2_center_lat*180./pi)
  264. call netcdf_error_handler(ncstat)
  265. !***
  266. !*** write grid center longitude array
  267. !***
  268. ncstat = nf_put_var_double(nc_outfile_id, nc_dstgrdcntrlon_id,
  269. & grid2_center_lon*180./pi)
  270. call netcdf_error_handler(ncstat)
  271. !-----------------------------------------------------------------------
  272. !
  273. ! Interpolate the input mask
  274. !
  275. !-----------------------------------------------------------------------
  276. call remap(grid1_mask_grid2, wts_map1, grid2_add_map1,
  277. & grid1_add_map1, real(grid1_imask,kind=dbl_kind))
  278. !-----------------------------------------------------------------------
  279. !
  280. ! Write time dimension
  281. !
  282. !-----------------------------------------------------------------------
  283. do jtime = 1,ntime
  284. ncstat = nf_put_vara_double(nc_outfile_id, nc_vartime_id,
  285. & jtime, 1, time(jtime))
  286. call netcdf_error_handler(ncstat)
  287. !-----------------------------------------------------------------------
  288. !
  289. ! if bicubic or 2nd-order conservative, 3 gradients needed in space
  290. !
  291. !-----------------------------------------------------------------------
  292. if ( fromregular .and. (map_type == map_type_bicubic .or.
  293. & map_type == map_type_conserv) ) then
  294. nx = grid1_dims(1)
  295. ny = grid1_dims(2)
  296. do n=1,grid1_size
  297. grad1_lat(n) = zero
  298. grad1_lon(n) = zero
  299. grad1_latlon(n) = zero
  300. if (grid1_mask(n)) then
  301. delew = half
  302. delns = half
  303. j = (n-1)/nx + 1
  304. i = n - (j-1)*nx
  305. ip1 = i+1
  306. im1 = i-1
  307. jp1 = j+1
  308. jm1 = j-1
  309. if (ip1 > nx) ip1 = ip1 - nx
  310. if (im1 < 1 ) im1 = nx
  311. if (jp1 > ny) then
  312. jp1 = j
  313. delns = one
  314. endif
  315. if (jm1 < 1 ) then
  316. jm1 = j
  317. delns = one
  318. endif
  319. in = (jp1-1)*nx + i
  320. is = (jm1-1)*nx + i
  321. ie = (j -1)*nx + ip1
  322. iw = (j -1)*nx + im1
  323. ine = (jp1-1)*nx + ip1
  324. inw = (jp1-1)*nx + im1
  325. ise = (jm1-1)*nx + ip1
  326. isw = (jm1-1)*nx + im1
  327. !*** compute i-gradient
  328. if (.not. grid1_mask(ie)) then
  329. ie = n
  330. delew = one
  331. endif
  332. if (.not. grid1_mask(iw)) then
  333. iw = n
  334. delew = one
  335. endif
  336. grad1_lat(n) = delew*(grid1_array(ie,jtime) -
  337. & grid1_array(iw,jtime))
  338. !*** compute j-gradient
  339. if (.not. grid1_mask(in)) then
  340. in = n
  341. delns = one
  342. endif
  343. if (.not. grid1_mask(is)) then
  344. is = n
  345. delns = one
  346. endif
  347. grad1_lon(n) = delns*(grid1_array(in,jtime) -
  348. & grid1_array(is,jtime))
  349. !*** compute ij-gradient
  350. delew = half
  351. if (jp1 == j .or. jm1 == j) then
  352. delns = one
  353. else
  354. delns = half
  355. endif
  356. if (.not. grid1_mask(ine)) then
  357. if (in /= n) then
  358. ine = in
  359. delew = one
  360. else if (ie /= n) then
  361. ine = ie
  362. inw = iw
  363. if (inw == n) delew = one
  364. delns = one
  365. else
  366. ine = n
  367. inw = iw
  368. delew = one
  369. delns = one
  370. endif
  371. endif
  372. if (.not. grid1_mask(inw)) then
  373. if (in /= n) then
  374. inw = in
  375. delew = one
  376. else if (iw /= n) then
  377. inw = iw
  378. ine = ie
  379. if (ie == n) delew = one
  380. delns = one
  381. else
  382. inw = n
  383. ine = ie
  384. delew = one
  385. delns = one
  386. endif
  387. endif
  388. grad1_lat_zero(n) = delew*(grid1_array(ine,jtime) -
  389. & grid1_array(inw,jtime))
  390. if (.not. grid1_mask(ise)) then
  391. if (is /= n) then
  392. ise = is
  393. delew = one
  394. else if (ie /= n) then
  395. ise = ie
  396. isw = iw
  397. if (isw == n) delew = one
  398. delns = one
  399. else
  400. ise = n
  401. isw = iw
  402. delew = one
  403. delns = one
  404. endif
  405. endif
  406. if (.not. grid1_mask(isw)) then
  407. if (is /= n) then
  408. isw = is
  409. delew = one
  410. else if (iw /= n) then
  411. isw = iw
  412. ise = ie
  413. if (ie == n) delew = one
  414. delns = one
  415. else
  416. isw = n
  417. ise = ie
  418. delew = one
  419. delns = one
  420. endif
  421. endif
  422. grad1_lon_zero(n) = delew*(grid1_array(ise,jtime) -
  423. & grid1_array(isw,jtime))
  424. grad1_latlon(n) = delns*(grad1_lat_zero(n) -
  425. & grad1_lon_zero(n))
  426. endif
  427. enddo
  428. endif
  429. !-----------------------------------------------------------------------
  430. !
  431. ! remapping from grid1 to grid2
  432. !
  433. !-----------------------------------------------------------------------
  434. grad1_lat_zero = zero
  435. grad1_lon_zero = zero
  436. if (map_type == map_type_bicubic) then
  437. if (fromregular) then
  438. call remap(grid2_tmp, wts_map1, grid2_add_map1,
  439. & grid1_add_map1, grid1_array(:,jtime),
  440. & src_grad1=grad1_lat,
  441. & src_grad2=grad1_lon,
  442. & src_grad3=grad1_latlon)
  443. else
  444. print*,"Input grid is not regular, bicubic interpolation "
  445. stop" is not possible : We stop "
  446. endif
  447. elseif (map_type == map_type_conserv .and. fromregular ) then
  448. call remap(grid2_tmp,wts_map1,grid2_add_map1,grid1_add_map1,
  449. & grid1_array(:,jtime), src_grad1=grad1_lat,
  450. & src_grad2=grad1_lon)
  451. else
  452. call remap(grid2_tmp, wts_map1, grid2_add_map1, grid1_add_map1,
  453. & grid1_array(:,jtime))
  454. endif
  455. if (map_type == map_type_conserv) then
  456. select case (norm_opt)
  457. case (norm_opt_none)
  458. grid2_err = grid2_frac*grid2_area
  459. where (grid2_err /= zero)
  460. grid2_tmp = grid2_tmp/grid2_err
  461. else where
  462. grid2_tmp = zero
  463. end where
  464. case (norm_opt_frcarea)
  465. case (norm_opt_dstarea)
  466. where (grid2_frac /= zero)
  467. grid2_tmp = grid2_tmp/grid2_frac
  468. else where
  469. grid2_tmp = zero
  470. end where
  471. end select
  472. end if
  473. !-----------------------------------------------------------------------
  474. !
  475. ! write results to NetCDF file
  476. !
  477. !-----------------------------------------------------------------------
  478. where (grid2_imask<0.5 .or. grid1_mask_grid2 == 0. )
  479. grid2_tmp=1e20
  480. end where
  481. ncstat = nf_put_vara_double(nc_outfile_id, nc_dstarray1_id,
  482. & (/1, 1, jtime/), (/grid2_dims, 1/), grid2_tmp )
  483. call netcdf_error_handler(ncstat)
  484. where (grid1_imask<0.5)
  485. grid1_array(:,jtime)=1e20
  486. end where
  487. ncstat = nf_put_vara_double(nc_outfile_id, nc_srcarray_id,
  488. & (/1, 1, jtime/), (/grid1_dims, 1/), grid1_array(:,jtime))
  489. call netcdf_error_handler(ncstat)
  490. enddo
  491. !-----------------------------------------------------------------------
  492. !
  493. ! close netCDF file
  494. !
  495. !-----------------------------------------------------------------------
  496. ncstat = nf_close(nc_outfile_id)
  497. call netcdf_error_handler(ncstat)
  498. !-----------------------------------------------------------------------
  499. end program scrip_use
  500. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!