remap_read.f 36 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! This routine reads remapping information from files written
  4. ! by remap_setup. If remapping in both directions are required,
  5. ! two input files must be specified.
  6. !
  7. !-----------------------------------------------------------------------
  8. !
  9. ! CVS:$Id: remap_read.f,v 1.6 2000/04/19 21:56:26 pwjones Exp $
  10. !
  11. ! Copyright (c) 1997, 1998 the Regents of the University of
  12. ! California.
  13. !
  14. ! This software and ancillary information (herein called software)
  15. ! called SCRIP is made available under the terms described here.
  16. ! The software has been approved for release with associated
  17. ! LA-CC Number 98-45.
  18. !
  19. ! Unless otherwise indicated, this software has been authored
  20. ! by an employee or employees of the University of California,
  21. ! operator of the Los Alamos National Laboratory under Contract
  22. ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
  23. ! Government has rights to use, reproduce, and distribute this
  24. ! software. The public may copy and use this software without
  25. ! charge, provided that this Notice and any statement of authorship
  26. ! are reproduced on all copies. Neither the Government nor the
  27. ! University makes any warranty, express or implied, or assumes
  28. ! any liability or responsibility for the use of this software.
  29. !
  30. ! If software is modified to produce derivative works, such modified
  31. ! software should be clearly marked, so as not to confuse it with
  32. ! the version available from Los Alamos National Laboratory.
  33. !
  34. !***********************************************************************
  35. module remap_read
  36. !-----------------------------------------------------------------------
  37. !
  38. ! contains routines for reading a remap file
  39. !
  40. !-----------------------------------------------------------------------
  41. use kinds_mod ! defines common data types
  42. use constants ! defines useful constants
  43. use grids ! includes all grid information
  44. use netcdf_mod ! module with netcdf vars and utilities
  45. use remap_vars ! module for all required remapping variables
  46. implicit none
  47. !-----------------------------------------------------------------------
  48. !
  49. ! module variables
  50. !
  51. !-----------------------------------------------------------------------
  52. !-----------------------------------------------------------------------
  53. !
  54. ! various netCDF ids for files variables
  55. !
  56. !-----------------------------------------------------------------------
  57. integer (kind=int_kind), private :: ! netCDF ids
  58. & ncstat, nc_file_id,
  59. & nc_srcgrdsize_id, nc_dstgrdsize_id,
  60. & nc_srcgrdcorn_id, nc_dstgrdcorn_id,
  61. & nc_srcgrdrank_id, nc_dstgrdrank_id,
  62. & nc_srcgrddims_id, nc_dstgrddims_id,
  63. & nc_numlinks_id, nc_numwgts_id,
  64. & nc_srcgrdimask_id, nc_dstgrdimask_id,
  65. & nc_srcgrdcntrlat_id, nc_srcgrdcntrlon_id,
  66. & nc_srcgrdcrnrlat_id, nc_srcgrdcrnrlon_id,
  67. & nc_srcgrdarea_id, nc_srcgrdfrac_id,
  68. & nc_dstgrdcntrlat_id, nc_dstgrdcntrlon_id,
  69. & nc_dstgrdcrnrlat_id, nc_dstgrdcrnrlon_id,
  70. & nc_dstgrdarea_id, nc_dstgrdfrac_id,
  71. & nc_srcgrdadd_id, nc_dstgrdadd_id, nc_rmpmatrix_id
  72. !***********************************************************************
  73. contains
  74. !***********************************************************************
  75. subroutine read_remap(map_name, interp_file)
  76. !-----------------------------------------------------------------------
  77. !
  78. ! this driver routine reads some global attributes and then
  79. ! calls a specific read routine based on file conventions
  80. !
  81. !-----------------------------------------------------------------------
  82. !-----------------------------------------------------------------------
  83. !
  84. ! input variables
  85. !
  86. !-----------------------------------------------------------------------
  87. character(char_len), intent(in) ::
  88. & interp_file ! filename for remap data
  89. !-----------------------------------------------------------------------
  90. !
  91. ! output variables
  92. !
  93. !-----------------------------------------------------------------------
  94. character(char_len), intent(out) ::
  95. & map_name ! name for mapping
  96. !-----------------------------------------------------------------------
  97. !
  98. ! local variables
  99. !
  100. !-----------------------------------------------------------------------
  101. character(char_len) ::
  102. & map_method ! character string for map_type
  103. &, normalize_opt ! character string for normalization option
  104. &, convention ! character string for output convention
  105. !-----------------------------------------------------------------------
  106. !
  107. ! open file and read some global information
  108. !
  109. !-----------------------------------------------------------------------
  110. ncstat = nf_open(interp_file, NF_NOWRITE, nc_file_id)
  111. call netcdf_error_handler(ncstat)
  112. !***
  113. !*** map name
  114. !***
  115. map_name = ' '
  116. ncstat = nf_get_att_text(nc_file_id, NF_GLOBAL, 'title',
  117. & map_name)
  118. call netcdf_error_handler(ncstat)
  119. print *,'Reading remapping:',trim(map_name)
  120. print *,'From file:',trim(interp_file)
  121. !***
  122. !*** normalization option
  123. !***
  124. normalize_opt = ' '
  125. ncstat = nf_get_att_text(nc_file_id, NF_GLOBAL, 'normalization',
  126. & normalize_opt)
  127. call netcdf_error_handler(ncstat)
  128. select case(normalize_opt)
  129. case ('none')
  130. norm_opt = norm_opt_none
  131. case ('fracarea')
  132. norm_opt = norm_opt_frcarea
  133. case ('destarea')
  134. norm_opt = norm_opt_dstarea
  135. case default
  136. print *,'normalize_opt = ',normalize_opt
  137. stop 'Invalid normalization option'
  138. end select
  139. !***
  140. !*** map method
  141. !***
  142. map_method = ' '
  143. ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'map_method',
  144. & map_method)
  145. call netcdf_error_handler(ncstat)
  146. select case(map_method)
  147. case('Conservative remapping')
  148. map_type = map_type_conserv
  149. case('Bilinear remapping')
  150. map_type = map_type_bilinear
  151. case('Distance weighted avg of nearest neighbors')
  152. map_type = map_type_distwgt
  153. case('Bicubic remapping')
  154. map_type = map_type_bicubic
  155. case default
  156. print *,'map_type = ',map_method
  157. stop 'Invalid Map Type'
  158. end select
  159. !***
  160. !*** file convention
  161. !***
  162. convention = ' '
  163. ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'conventions',
  164. & convention)
  165. call netcdf_error_handler(ncstat)
  166. !-----------------------------------------------------------------------
  167. !
  168. ! call appropriate read routine based on output convention
  169. !
  170. !-----------------------------------------------------------------------
  171. select case(convention)
  172. case ('SCRIP')
  173. call read_remap_scrip
  174. case ('NCAR-CSM')
  175. call read_remap_csm
  176. case default
  177. print *,'convention = ',convention
  178. stop 'unknown output file convention'
  179. end select
  180. !-----------------------------------------------------------------------
  181. end subroutine read_remap
  182. !***********************************************************************
  183. subroutine read_remap_scrip
  184. !-----------------------------------------------------------------------
  185. !
  186. ! the routine reads a netCDF file to extract remapping info
  187. ! in SCRIP format
  188. !
  189. !-----------------------------------------------------------------------
  190. !-----------------------------------------------------------------------
  191. !
  192. ! local variables
  193. !
  194. !-----------------------------------------------------------------------
  195. character (char_len) ::
  196. & grid1_name ! grid name for source grid
  197. &, grid2_name ! grid name for dest grid
  198. integer (kind=int_kind) ::
  199. & n ! dummy index
  200. integer (kind=int_kind), dimension(:), allocatable ::
  201. & grid1_mask_int, ! integer masks to determine
  202. & grid2_mask_int ! cells that participate in map
  203. !-----------------------------------------------------------------------
  204. !
  205. ! read some additional global attributes
  206. !
  207. !-----------------------------------------------------------------------
  208. !***
  209. !*** source and destination grid names
  210. !***
  211. grid1_name = ' '
  212. grid2_name = ' '
  213. ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'source_grid',
  214. & grid1_name)
  215. call netcdf_error_handler(ncstat)
  216. ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'dest_grid',
  217. & grid2_name)
  218. call netcdf_error_handler(ncstat)
  219. print *,' '
  220. print *,'Remapping between:',trim(grid1_name)
  221. print *,'and ',trim(grid2_name)
  222. print *,' '
  223. !-----------------------------------------------------------------------
  224. !
  225. ! read dimension information
  226. !
  227. !-----------------------------------------------------------------------
  228. ncstat = nf_inq_dimid(nc_file_id, 'src_grid_size',
  229. & nc_srcgrdsize_id)
  230. call netcdf_error_handler(ncstat)
  231. ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdsize_id, grid1_size)
  232. call netcdf_error_handler(ncstat)
  233. ncstat = nf_inq_dimid(nc_file_id, 'dst_grid_size',
  234. & nc_dstgrdsize_id)
  235. call netcdf_error_handler(ncstat)
  236. ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdsize_id, grid2_size)
  237. call netcdf_error_handler(ncstat)
  238. ncstat = nf_inq_dimid(nc_file_id, 'src_grid_corners',
  239. & nc_srcgrdcorn_id)
  240. call netcdf_error_handler(ncstat)
  241. ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdcorn_id,
  242. & grid1_corners)
  243. call netcdf_error_handler(ncstat)
  244. ncstat = nf_inq_dimid(nc_file_id, 'dst_grid_corners',
  245. & nc_dstgrdcorn_id)
  246. call netcdf_error_handler(ncstat)
  247. ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdcorn_id,
  248. & grid2_corners)
  249. call netcdf_error_handler(ncstat)
  250. ncstat = nf_inq_dimid(nc_file_id, 'src_grid_rank',
  251. & nc_srcgrdrank_id)
  252. call netcdf_error_handler(ncstat)
  253. ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdrank_id,
  254. & grid1_rank)
  255. call netcdf_error_handler(ncstat)
  256. ncstat = nf_inq_dimid(nc_file_id, 'dst_grid_rank',
  257. & nc_dstgrdrank_id)
  258. call netcdf_error_handler(ncstat)
  259. ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdrank_id,
  260. & grid2_rank)
  261. call netcdf_error_handler(ncstat)
  262. ncstat = nf_inq_dimid(nc_file_id, 'num_links',
  263. & nc_numlinks_id)
  264. call netcdf_error_handler(ncstat)
  265. ncstat = nf_inq_dimlen(nc_file_id, nc_numlinks_id,
  266. & num_links_map1)
  267. call netcdf_error_handler(ncstat)
  268. ncstat = nf_inq_dimid(nc_file_id, 'num_wgts',
  269. & nc_numwgts_id)
  270. call netcdf_error_handler(ncstat)
  271. ncstat = nf_inq_dimlen(nc_file_id, nc_numwgts_id, num_wts)
  272. call netcdf_error_handler(ncstat)
  273. !-----------------------------------------------------------------------
  274. !
  275. ! allocate arrays
  276. !
  277. !-----------------------------------------------------------------------
  278. allocate( grid1_dims (grid1_rank),
  279. & grid1_center_lat(grid1_size),
  280. & grid1_center_lon(grid1_size),
  281. & grid1_area (grid1_size),
  282. & grid1_frac (grid1_size),
  283. & grid1_mask (grid1_size),
  284. & grid1_mask_int (grid1_size),
  285. & grid1_corner_lat(grid1_corners, grid1_size),
  286. & grid1_corner_lon(grid1_corners, grid1_size) )
  287. allocate( grid2_dims (grid2_rank),
  288. & grid2_center_lat(grid2_size),
  289. & grid2_center_lon(grid2_size),
  290. & grid2_area (grid2_size),
  291. & grid2_frac (grid2_size),
  292. & grid2_mask (grid2_size),
  293. & grid2_mask_int (grid2_size),
  294. & grid2_corner_lat(grid2_corners, grid2_size),
  295. & grid2_corner_lon(grid2_corners, grid2_size) )
  296. allocate( grid1_add_map1(num_links_map1),
  297. & grid2_add_map1(num_links_map1),
  298. & wts_map1(num_wts,num_links_map1) )
  299. !-----------------------------------------------------------------------
  300. !
  301. ! get variable ids
  302. !
  303. !-----------------------------------------------------------------------
  304. ncstat = nf_inq_varid(nc_file_id, 'src_grid_dims',
  305. & nc_srcgrddims_id)
  306. call netcdf_error_handler(ncstat)
  307. ncstat = nf_inq_varid(nc_file_id, 'src_grid_imask',
  308. & nc_srcgrdimask_id)
  309. call netcdf_error_handler(ncstat)
  310. ncstat = nf_inq_varid(nc_file_id, 'src_grid_center_lat',
  311. & nc_srcgrdcntrlat_id)
  312. call netcdf_error_handler(ncstat)
  313. ncstat = nf_inq_varid(nc_file_id, 'src_grid_center_lon',
  314. & nc_srcgrdcntrlon_id)
  315. call netcdf_error_handler(ncstat)
  316. ncstat = nf_inq_varid(nc_file_id, 'src_grid_corner_lat',
  317. & nc_srcgrdcrnrlat_id)
  318. call netcdf_error_handler(ncstat)
  319. ncstat = nf_inq_varid(nc_file_id, 'src_grid_corner_lon',
  320. & nc_srcgrdcrnrlon_id)
  321. call netcdf_error_handler(ncstat)
  322. ncstat = nf_inq_varid(nc_file_id, 'src_grid_area',
  323. & nc_srcgrdarea_id)
  324. call netcdf_error_handler(ncstat)
  325. ncstat = nf_inq_varid(nc_file_id, 'src_grid_frac',
  326. & nc_srcgrdfrac_id)
  327. call netcdf_error_handler(ncstat)
  328. ncstat = nf_inq_varid(nc_file_id, 'dst_grid_dims',
  329. & nc_dstgrddims_id)
  330. call netcdf_error_handler(ncstat)
  331. ncstat = nf_inq_varid(nc_file_id, 'dst_grid_imask',
  332. & nc_dstgrdimask_id)
  333. call netcdf_error_handler(ncstat)
  334. ncstat = nf_inq_varid(nc_file_id, 'dst_grid_center_lat',
  335. & nc_dstgrdcntrlat_id)
  336. call netcdf_error_handler(ncstat)
  337. ncstat = nf_inq_varid(nc_file_id, 'dst_grid_center_lon',
  338. & nc_dstgrdcntrlon_id)
  339. call netcdf_error_handler(ncstat)
  340. ncstat = nf_inq_varid(nc_file_id, 'dst_grid_corner_lat',
  341. & nc_dstgrdcrnrlat_id)
  342. call netcdf_error_handler(ncstat)
  343. ncstat = nf_inq_varid(nc_file_id, 'dst_grid_corner_lon',
  344. & nc_dstgrdcrnrlon_id)
  345. call netcdf_error_handler(ncstat)
  346. ncstat = nf_inq_varid(nc_file_id, 'dst_grid_area',
  347. & nc_dstgrdarea_id)
  348. call netcdf_error_handler(ncstat)
  349. ncstat = nf_inq_varid(nc_file_id, 'dst_grid_frac',
  350. & nc_dstgrdfrac_id)
  351. call netcdf_error_handler(ncstat)
  352. ncstat = nf_inq_varid(nc_file_id, 'src_address',
  353. & nc_srcgrdadd_id)
  354. call netcdf_error_handler(ncstat)
  355. ncstat = nf_inq_varid(nc_file_id, 'dst_address',
  356. & nc_dstgrdadd_id)
  357. call netcdf_error_handler(ncstat)
  358. ncstat = nf_inq_varid(nc_file_id, 'remap_matrix',
  359. & nc_rmpmatrix_id)
  360. call netcdf_error_handler(ncstat)
  361. !-----------------------------------------------------------------------
  362. !
  363. ! read all variables
  364. !
  365. !-----------------------------------------------------------------------
  366. ncstat = nf_get_var_int(nc_file_id, nc_srcgrddims_id,
  367. & grid1_dims)
  368. call netcdf_error_handler(ncstat)
  369. ncstat = nf_get_var_int(nc_file_id, nc_srcgrdimask_id,
  370. & grid1_mask_int)
  371. call netcdf_error_handler(ncstat)
  372. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcntrlat_id,
  373. & grid1_center_lat)
  374. call netcdf_error_handler(ncstat)
  375. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcntrlon_id,
  376. & grid1_center_lon)
  377. call netcdf_error_handler(ncstat)
  378. grid1_units = ' '
  379. ncstat = nf_get_att_text(nc_file_id, nc_srcgrdcntrlat_id, 'units',
  380. & grid1_units)
  381. call netcdf_error_handler(ncstat)
  382. select case (grid1_units(1:7))
  383. case ('degrees')
  384. grid1_center_lat = grid1_center_lat*deg2rad
  385. grid1_center_lon = grid1_center_lon*deg2rad
  386. case ('radians')
  387. !*** no conversion necessary
  388. case default
  389. print *,'unknown units supplied for grid1 center lat/lon: '
  390. print *,'proceeding assuming radians'
  391. end select
  392. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcrnrlat_id,
  393. & grid1_corner_lat)
  394. call netcdf_error_handler(ncstat)
  395. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcrnrlon_id,
  396. & grid1_corner_lon)
  397. call netcdf_error_handler(ncstat)
  398. grid1_units = ' '
  399. ncstat = nf_get_att_text(nc_file_id, nc_srcgrdcrnrlat_id, 'units',
  400. & grid1_units)
  401. call netcdf_error_handler(ncstat)
  402. select case (grid1_units(1:7))
  403. case ('degrees')
  404. grid1_corner_lat = grid1_corner_lat*deg2rad
  405. grid1_corner_lon = grid1_corner_lon*deg2rad
  406. case ('radians')
  407. !*** no conversion necessary
  408. case default
  409. print *,'unknown units supplied for grid1 corner lat/lon: '
  410. print *,'proceeding assuming radians'
  411. end select
  412. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdarea_id,
  413. & grid1_area)
  414. call netcdf_error_handler(ncstat)
  415. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdfrac_id,
  416. & grid1_frac)
  417. call netcdf_error_handler(ncstat)
  418. ncstat = nf_get_var_int(nc_file_id, nc_dstgrddims_id,
  419. & grid2_dims)
  420. call netcdf_error_handler(ncstat)
  421. ncstat = nf_get_var_int(nc_file_id, nc_dstgrdimask_id,
  422. & grid2_mask_int)
  423. call netcdf_error_handler(ncstat)
  424. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcntrlat_id,
  425. & grid2_center_lat)
  426. call netcdf_error_handler(ncstat)
  427. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcntrlon_id,
  428. & grid2_center_lon)
  429. call netcdf_error_handler(ncstat)
  430. grid2_units = ' '
  431. ncstat = nf_get_att_text(nc_file_id, nc_dstgrdcntrlat_id, 'units',
  432. & grid2_units)
  433. call netcdf_error_handler(ncstat)
  434. select case (grid2_units(1:7))
  435. case ('degrees')
  436. grid2_center_lat = grid2_center_lat*deg2rad
  437. grid2_center_lon = grid2_center_lon*deg2rad
  438. case ('radians')
  439. !*** no conversion necessary
  440. case default
  441. print *,'unknown units supplied for grid2 center lat/lon: '
  442. print *,'proceeding assuming radians'
  443. end select
  444. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcrnrlat_id,
  445. & grid2_corner_lat)
  446. call netcdf_error_handler(ncstat)
  447. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcrnrlon_id,
  448. & grid2_corner_lon)
  449. call netcdf_error_handler(ncstat)
  450. grid2_units = ' '
  451. ncstat = nf_get_att_text(nc_file_id, nc_dstgrdcrnrlat_id, 'units',
  452. & grid2_units)
  453. call netcdf_error_handler(ncstat)
  454. select case (grid2_units(1:7))
  455. case ('degrees')
  456. grid2_corner_lat = grid2_corner_lat*deg2rad
  457. grid2_corner_lon = grid2_corner_lon*deg2rad
  458. case ('radians')
  459. !*** no conversion necessary
  460. case default
  461. print *,'unknown units supplied for grid2 corner lat/lon: '
  462. print *,'proceeding assuming radians'
  463. end select
  464. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdarea_id,
  465. & grid2_area)
  466. call netcdf_error_handler(ncstat)
  467. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdfrac_id,
  468. & grid2_frac)
  469. call netcdf_error_handler(ncstat)
  470. ncstat = nf_get_var_int(nc_file_id, nc_srcgrdadd_id,
  471. & grid1_add_map1)
  472. call netcdf_error_handler(ncstat)
  473. ncstat = nf_get_var_int(nc_file_id, nc_dstgrdadd_id,
  474. & grid2_add_map1)
  475. call netcdf_error_handler(ncstat)
  476. ncstat = nf_get_var_double(nc_file_id, nc_rmpmatrix_id,
  477. & wts_map1)
  478. call netcdf_error_handler(ncstat)
  479. !-----------------------------------------------------------------------
  480. !
  481. ! initialize logical mask
  482. !
  483. !-----------------------------------------------------------------------
  484. where (grid1_mask_int == 1)
  485. grid1_mask = .true.
  486. elsewhere
  487. grid1_mask = .false.
  488. endwhere
  489. where (grid2_mask_int == 1)
  490. grid2_mask = .true.
  491. elsewhere
  492. grid2_mask = .false.
  493. endwhere
  494. deallocate(grid1_mask_int, grid2_mask_int)
  495. !-----------------------------------------------------------------------
  496. !
  497. ! close input file
  498. !
  499. !-----------------------------------------------------------------------
  500. ncstat = nf_close(nc_file_id)
  501. call netcdf_error_handler(ncstat)
  502. !-----------------------------------------------------------------------
  503. end subroutine read_remap_scrip
  504. !***********************************************************************
  505. subroutine read_remap_csm
  506. !-----------------------------------------------------------------------
  507. !
  508. ! the routine reads a netCDF file to extract remapping info
  509. ! in NCAR-CSM format
  510. !
  511. !-----------------------------------------------------------------------
  512. !-----------------------------------------------------------------------
  513. !
  514. ! local variables
  515. !
  516. !-----------------------------------------------------------------------
  517. character (char_len) ::
  518. & grid1_name ! grid name for source grid
  519. &, grid2_name ! grid name for dest grid
  520. integer (kind=int_kind) ::
  521. & nc_numwgts1_id ! extra netCDF id for num_wgts > 1
  522. &, nc_rmpmatrix2_id ! extra netCDF id for high-order remap matrix
  523. real (kind=dbl_kind), dimension(:),allocatable ::
  524. & wts1 ! CSM wants single array for 1st-order wts
  525. real (kind=dbl_kind), dimension(:,:),allocatable ::
  526. & wts2 ! write remaining weights in different array
  527. integer (kind=int_kind) ::
  528. & n ! dummy index
  529. integer (kind=int_kind), dimension(:), allocatable ::
  530. & grid1_mask_int, ! integer masks to determine
  531. & grid2_mask_int ! cells that participate in map
  532. !-----------------------------------------------------------------------
  533. !
  534. ! read some additional global attributes
  535. !
  536. !-----------------------------------------------------------------------
  537. !***
  538. !*** source and destination grid names
  539. !***
  540. grid1_name = ' '
  541. grid2_name = ' '
  542. ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'domain_a',
  543. & grid1_name)
  544. call netcdf_error_handler(ncstat)
  545. ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'domain_b',
  546. & grid2_name)
  547. call netcdf_error_handler(ncstat)
  548. print *,' '
  549. print *,'Remapping between:',trim(grid1_name)
  550. print *,'and ',trim(grid2_name)
  551. print *,' '
  552. !-----------------------------------------------------------------------
  553. !
  554. ! read dimension information
  555. !
  556. !-----------------------------------------------------------------------
  557. ncstat = nf_inq_dimid(nc_file_id, 'n_a', nc_srcgrdsize_id)
  558. call netcdf_error_handler(ncstat)
  559. ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdsize_id, grid1_size)
  560. call netcdf_error_handler(ncstat)
  561. ncstat = nf_inq_dimid(nc_file_id, 'n_b', nc_dstgrdsize_id)
  562. call netcdf_error_handler(ncstat)
  563. ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdsize_id, grid2_size)
  564. call netcdf_error_handler(ncstat)
  565. ncstat = nf_inq_dimid(nc_file_id, 'nv_a', nc_srcgrdcorn_id)
  566. call netcdf_error_handler(ncstat)
  567. ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdcorn_id,
  568. & grid1_corners)
  569. call netcdf_error_handler(ncstat)
  570. ncstat = nf_inq_dimid(nc_file_id, 'nv_b', nc_dstgrdcorn_id)
  571. call netcdf_error_handler(ncstat)
  572. ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdcorn_id,
  573. & grid2_corners)
  574. call netcdf_error_handler(ncstat)
  575. ncstat = nf_inq_dimid(nc_file_id, 'src_grid_rank',
  576. & nc_srcgrdrank_id)
  577. call netcdf_error_handler(ncstat)
  578. ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdrank_id,
  579. & grid1_rank)
  580. call netcdf_error_handler(ncstat)
  581. ncstat = nf_inq_dimid(nc_file_id, 'dst_grid_rank',
  582. & nc_dstgrdrank_id)
  583. call netcdf_error_handler(ncstat)
  584. ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdrank_id,
  585. & grid2_rank)
  586. call netcdf_error_handler(ncstat)
  587. ncstat = nf_inq_dimid(nc_file_id, 'n_s',
  588. & nc_numlinks_id)
  589. call netcdf_error_handler(ncstat)
  590. ncstat = nf_inq_dimlen(nc_file_id, nc_numlinks_id,
  591. & num_links_map1)
  592. call netcdf_error_handler(ncstat)
  593. ncstat = nf_inq_dimid(nc_file_id, 'num_wgts',
  594. & nc_numwgts_id)
  595. call netcdf_error_handler(ncstat)
  596. ncstat = nf_inq_dimlen(nc_file_id, nc_numwgts_id, num_wts)
  597. call netcdf_error_handler(ncstat)
  598. if (num_wts > 1) then
  599. ncstat = nf_inq_dimid(nc_file_id, 'num_wgts1',
  600. & nc_numwgts1_id)
  601. call netcdf_error_handler(ncstat)
  602. endif
  603. !-----------------------------------------------------------------------
  604. !
  605. ! allocate arrays
  606. !
  607. !-----------------------------------------------------------------------
  608. allocate( grid1_dims (grid1_rank),
  609. & grid1_center_lat(grid1_size),
  610. & grid1_center_lon(grid1_size),
  611. & grid1_area (grid1_size),
  612. & grid1_frac (grid1_size),
  613. & grid1_mask (grid1_size),
  614. & grid1_mask_int (grid1_size),
  615. & grid1_corner_lat(grid1_corners, grid1_size),
  616. & grid1_corner_lon(grid1_corners, grid1_size) )
  617. allocate( grid2_dims (grid2_rank),
  618. & grid2_center_lat(grid2_size),
  619. & grid2_center_lon(grid2_size),
  620. & grid2_area (grid2_size),
  621. & grid2_frac (grid2_size),
  622. & grid2_mask (grid2_size),
  623. & grid2_mask_int (grid2_size),
  624. & grid2_corner_lat(grid2_corners, grid2_size),
  625. & grid2_corner_lon(grid2_corners, grid2_size) )
  626. allocate( grid1_add_map1(num_links_map1),
  627. & grid2_add_map1(num_links_map1),
  628. & wts_map1(num_wts,num_links_map1),
  629. & wts1(num_links_map1),
  630. & wts2(num_wts-1,num_links_map1) )
  631. !-----------------------------------------------------------------------
  632. !
  633. ! get variable ids
  634. !
  635. !-----------------------------------------------------------------------
  636. ncstat = nf_inq_varid(nc_file_id, 'src_grid_dims',
  637. & nc_srcgrddims_id)
  638. call netcdf_error_handler(ncstat)
  639. ncstat = nf_inq_varid(nc_file_id, 'mask_a',
  640. & nc_srcgrdimask_id)
  641. call netcdf_error_handler(ncstat)
  642. ncstat = nf_inq_varid(nc_file_id, 'yc_a', nc_srcgrdcntrlat_id)
  643. call netcdf_error_handler(ncstat)
  644. ncstat = nf_inq_varid(nc_file_id, 'xc_a', nc_srcgrdcntrlon_id)
  645. call netcdf_error_handler(ncstat)
  646. ncstat = nf_inq_varid(nc_file_id, 'yv_a', nc_srcgrdcrnrlat_id)
  647. call netcdf_error_handler(ncstat)
  648. ncstat = nf_inq_varid(nc_file_id, 'xv_a', nc_srcgrdcrnrlon_id)
  649. call netcdf_error_handler(ncstat)
  650. ncstat = nf_inq_varid(nc_file_id, 'area_a', nc_srcgrdarea_id)
  651. call netcdf_error_handler(ncstat)
  652. ncstat = nf_inq_varid(nc_file_id, 'frac_a', nc_srcgrdfrac_id)
  653. call netcdf_error_handler(ncstat)
  654. ncstat = nf_inq_varid(nc_file_id, 'dst_grid_dims',
  655. & nc_dstgrddims_id)
  656. call netcdf_error_handler(ncstat)
  657. ncstat = nf_inq_varid(nc_file_id, 'mask_b',
  658. & nc_dstgrdimask_id)
  659. call netcdf_error_handler(ncstat)
  660. ncstat = nf_inq_varid(nc_file_id, 'yc_b', nc_dstgrdcntrlat_id)
  661. call netcdf_error_handler(ncstat)
  662. ncstat = nf_inq_varid(nc_file_id, 'xc_b', nc_dstgrdcntrlon_id)
  663. call netcdf_error_handler(ncstat)
  664. ncstat = nf_inq_varid(nc_file_id, 'yv_b', nc_dstgrdcrnrlat_id)
  665. call netcdf_error_handler(ncstat)
  666. ncstat = nf_inq_varid(nc_file_id, 'xv_b', nc_dstgrdcrnrlon_id)
  667. call netcdf_error_handler(ncstat)
  668. ncstat = nf_inq_varid(nc_file_id, 'area_b', nc_dstgrdarea_id)
  669. call netcdf_error_handler(ncstat)
  670. ncstat = nf_inq_varid(nc_file_id, 'frac_b', nc_dstgrdfrac_id)
  671. call netcdf_error_handler(ncstat)
  672. ncstat = nf_inq_varid(nc_file_id, 'col', nc_srcgrdadd_id)
  673. call netcdf_error_handler(ncstat)
  674. ncstat = nf_inq_varid(nc_file_id, 'row', nc_dstgrdadd_id)
  675. call netcdf_error_handler(ncstat)
  676. ncstat = nf_inq_varid(nc_file_id, 'S', nc_rmpmatrix_id)
  677. call netcdf_error_handler(ncstat)
  678. if (num_wts > 1) then
  679. ncstat = nf_inq_varid(nc_file_id, 'S2', nc_rmpmatrix2_id)
  680. call netcdf_error_handler(ncstat)
  681. endif
  682. !-----------------------------------------------------------------------
  683. !
  684. ! read all variables
  685. !
  686. !-----------------------------------------------------------------------
  687. ncstat = nf_get_var_int(nc_file_id, nc_srcgrddims_id,
  688. & grid1_dims)
  689. call netcdf_error_handler(ncstat)
  690. ncstat = nf_get_var_int(nc_file_id, nc_srcgrdimask_id,
  691. & grid1_mask_int)
  692. call netcdf_error_handler(ncstat)
  693. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcntrlat_id,
  694. & grid1_center_lat)
  695. call netcdf_error_handler(ncstat)
  696. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcntrlon_id,
  697. & grid1_center_lon)
  698. call netcdf_error_handler(ncstat)
  699. ncstat = nf_get_att_text(nc_file_id, nc_srcgrdcntrlat_id, 'units',
  700. & grid1_units)
  701. call netcdf_error_handler(ncstat)
  702. select case (grid1_units(1:7))
  703. case ('degrees')
  704. grid1_center_lat = grid1_center_lat*deg2rad
  705. grid1_center_lon = grid1_center_lon*deg2rad
  706. case ('radians')
  707. !*** no conversion necessary
  708. case default
  709. print *,'unknown units supplied for grid1 center lat/lon: '
  710. print *,'proceeding assuming radians'
  711. end select
  712. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcrnrlat_id,
  713. & grid1_corner_lat)
  714. call netcdf_error_handler(ncstat)
  715. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcrnrlon_id,
  716. & grid1_corner_lon)
  717. call netcdf_error_handler(ncstat)
  718. ncstat = nf_get_att_text(nc_file_id, nc_srcgrdcrnrlat_id, 'units',
  719. & grid1_units)
  720. call netcdf_error_handler(ncstat)
  721. select case (grid1_units(1:7))
  722. case ('degrees')
  723. grid1_corner_lat = grid1_corner_lat*deg2rad
  724. grid1_corner_lon = grid1_corner_lon*deg2rad
  725. case ('radians')
  726. !*** no conversion necessary
  727. case default
  728. print *,'unknown units supplied for grid1 corner lat/lon: '
  729. print *,'proceeding assuming radians'
  730. end select
  731. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdarea_id,
  732. & grid1_area)
  733. call netcdf_error_handler(ncstat)
  734. ncstat = nf_get_var_double(nc_file_id, nc_srcgrdfrac_id,
  735. & grid1_frac)
  736. call netcdf_error_handler(ncstat)
  737. ncstat = nf_get_var_int(nc_file_id, nc_dstgrddims_id,
  738. & grid2_dims)
  739. call netcdf_error_handler(ncstat)
  740. ncstat = nf_get_var_int(nc_file_id, nc_dstgrdimask_id,
  741. & grid2_mask_int)
  742. call netcdf_error_handler(ncstat)
  743. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcntrlat_id,
  744. & grid2_center_lat)
  745. call netcdf_error_handler(ncstat)
  746. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcntrlon_id,
  747. & grid2_center_lon)
  748. call netcdf_error_handler(ncstat)
  749. ncstat = nf_get_att_text(nc_file_id, nc_dstgrdcntrlat_id, 'units',
  750. & grid2_units)
  751. call netcdf_error_handler(ncstat)
  752. select case (grid2_units(1:7))
  753. case ('degrees')
  754. grid2_center_lat = grid2_center_lat*deg2rad
  755. grid2_center_lon = grid2_center_lon*deg2rad
  756. case ('radians')
  757. !*** no conversion necessary
  758. case default
  759. print *,'unknown units supplied for grid2 center lat/lon: '
  760. print *,'proceeding assuming radians'
  761. end select
  762. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcrnrlat_id,
  763. & grid2_corner_lat)
  764. call netcdf_error_handler(ncstat)
  765. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcrnrlon_id,
  766. & grid2_corner_lon)
  767. call netcdf_error_handler(ncstat)
  768. ncstat = nf_get_att_text(nc_file_id, nc_dstgrdcrnrlat_id, 'units',
  769. & grid2_units)
  770. call netcdf_error_handler(ncstat)
  771. select case (grid2_units(1:7))
  772. case ('degrees')
  773. grid2_corner_lat = grid2_corner_lat*deg2rad
  774. grid2_corner_lon = grid2_corner_lon*deg2rad
  775. case ('radians')
  776. !*** no conversion necessary
  777. case default
  778. print *,'unknown units supplied for grid2 corner lat/lon: '
  779. print *,'proceeding assuming radians'
  780. end select
  781. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdarea_id,
  782. & grid2_area)
  783. call netcdf_error_handler(ncstat)
  784. ncstat = nf_get_var_double(nc_file_id, nc_dstgrdfrac_id,
  785. & grid2_frac)
  786. call netcdf_error_handler(ncstat)
  787. ncstat = nf_get_var_int(nc_file_id, nc_srcgrdadd_id,
  788. & grid1_add_map1)
  789. call netcdf_error_handler(ncstat)
  790. ncstat = nf_get_var_int(nc_file_id, nc_dstgrdadd_id,
  791. & grid2_add_map1)
  792. call netcdf_error_handler(ncstat)
  793. ncstat = nf_get_var_double(nc_file_id, nc_rmpmatrix_id,
  794. & wts1)
  795. wts_map1(1,:) = wts1
  796. deallocate(wts1)
  797. if (num_wts > 1) then
  798. ncstat = nf_get_var_double(nc_file_id, nc_rmpmatrix2_id,
  799. & wts2)
  800. wts_map1(2:,:) = wts2
  801. deallocate(wts2)
  802. endif
  803. call netcdf_error_handler(ncstat)
  804. !-----------------------------------------------------------------------
  805. !
  806. ! initialize logical mask
  807. !
  808. !-----------------------------------------------------------------------
  809. where (grid1_mask_int == 1)
  810. grid1_mask = .true.
  811. elsewhere
  812. grid1_mask = .false.
  813. endwhere
  814. where (grid2_mask_int == 1)
  815. grid2_mask = .true.
  816. elsewhere
  817. grid2_mask = .false.
  818. endwhere
  819. deallocate(grid1_mask_int, grid2_mask_int)
  820. !-----------------------------------------------------------------------
  821. !
  822. ! close input file
  823. !
  824. !-----------------------------------------------------------------------
  825. ncstat = nf_close(nc_file_id)
  826. call netcdf_error_handler(ncstat)
  827. !-----------------------------------------------------------------------
  828. end subroutine read_remap_csm
  829. !***********************************************************************
  830. end module remap_read
  831. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!