scrip_test.f 31 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! this program is a short driver that tests the remappings using
  4. ! a simple analytic field. the results are written in netCDF
  5. ! format.
  6. !
  7. ! CVS: $Id: scrip_test.f,v 1.6 2000/04/19 21:45:09 pwjones Exp $
  8. !
  9. !-----------------------------------------------------------------------
  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. program remap_test
  36. !-----------------------------------------------------------------------
  37. use kinds_mod ! defines common data types
  38. use constants ! defines common constants
  39. use iounits ! I/O unit manager
  40. use netcdf_mod ! netcdf I/O stuff
  41. use grids ! module containing grid info
  42. use remap_vars ! module containing remapping info
  43. use remap_mod ! module containing remapping routines
  44. use remap_read ! routines for reading remap files
  45. implicit none
  46. !-----------------------------------------------------------------------
  47. !
  48. ! input namelist variables
  49. !
  50. !-----------------------------------------------------------------------
  51. integer (kind=int_kind) ::
  52. & field_choice ! choice of field to be interpolated
  53. character (char_len) ::
  54. & interp_file, ! filename containing remap data (map1)
  55. & output_file ! filename for test results
  56. namelist /remap_inputs/ field_choice, interp_file, output_file
  57. !-----------------------------------------------------------------------
  58. !
  59. ! local variables
  60. !
  61. !-----------------------------------------------------------------------
  62. character (char_len) ::
  63. & map_name ! name for mapping from grid1 to grid2
  64. integer (kind=int_kind) :: ! netCDF ids for files and arrays
  65. & ncstat, nc_outfile_id,
  66. & nc_srcgrdcntrlat_id, nc_srcgrdcntrlon_id,
  67. & nc_dstgrdcntrlat_id, nc_dstgrdcntrlon_id,
  68. & nc_srcgrdrank_id, nc_dstgrdrank_id,
  69. & nc_srcgrdimask_id, nc_dstgrdimask_id,
  70. & nc_srcgrdarea_id, nc_dstgrdarea_id,
  71. & nc_srcgrdfrac_id, nc_dstgrdfrac_id,
  72. & nc_srcarray_id, nc_srcgradlat_id, nc_srcgradlon_id,
  73. & nc_dstarray1_id, nc_dstarray1a_id, nc_dstarray2_id,
  74. & nc_dsterror1_id, nc_dsterror1a_id, nc_dsterror2_id
  75. integer (kind=int_kind), dimension(:), allocatable ::
  76. & nc_grid1size_id, nc_grid2size_id
  77. !-----------------------------------------------------------------------
  78. character (char_len) ::
  79. & dim_name ! netCDF dimension name
  80. integer (kind=int_kind) :: i,j,n,imin,imax,idiff,
  81. & ip1,im1,jp1,jm1,nx,ny, ! for computing bicub gradients
  82. & in,is,ie,iw,ine,inw,ise,isw,
  83. & iunit ! unit number for namelist file
  84. integer (kind=int_kind), dimension(:), allocatable ::
  85. & grid1_imask, grid2_imask, grid2_count
  86. real (kind=dbl_kind) ::
  87. & delew, delns, ! variables for computing bicub gradients
  88. & length ! length scale for cosine hill test field
  89. real (kind=dbl_kind), dimension(:), allocatable ::
  90. & grid1_array,
  91. & grid1_tmp,
  92. & grad1_lat,
  93. & grad1_lon,
  94. & grad1_latlon,
  95. & grad1_lat_zero,
  96. & grad1_lon_zero,
  97. & grid2_array,
  98. & grid2_err,
  99. & grid2_tmp
  100. !-----------------------------------------------------------------------
  101. !
  102. ! read namelist for file and mapping info
  103. !
  104. !-----------------------------------------------------------------------
  105. call get_unit(iunit)
  106. open(iunit, file='scrip_test_in', status='old', form='formatted')
  107. read(iunit, nml=remap_inputs)
  108. call release_unit(iunit)
  109. write(*,nml=remap_inputs)
  110. !-----------------------------------------------------------------------
  111. !
  112. ! read remapping data
  113. !
  114. !-----------------------------------------------------------------------
  115. call read_remap(map_name, interp_file)
  116. !-----------------------------------------------------------------------
  117. !
  118. ! allocate arrays
  119. !
  120. !-----------------------------------------------------------------------
  121. allocate (grid1_array (grid1_size),
  122. & grid1_tmp (grid1_size),
  123. & grad1_lat (grid1_size),
  124. & grad1_lon (grid1_size),
  125. & grad1_lat_zero (grid1_size),
  126. & grad1_lon_zero (grid1_size),
  127. & grid1_imask (grid1_size),
  128. & grid2_array (grid2_size),
  129. & grid2_err (grid2_size),
  130. & grid2_tmp (grid2_size),
  131. & grid2_imask (grid2_size),
  132. & grid2_count (grid2_size))
  133. where (grid1_mask)
  134. grid1_imask = 1
  135. elsewhere
  136. grid1_imask = 0
  137. endwhere
  138. where (grid2_mask)
  139. grid2_imask = 1
  140. elsewhere
  141. grid2_imask = 0
  142. endwhere
  143. !-----------------------------------------------------------------------
  144. !
  145. ! setup a NetCDF file for output
  146. !
  147. !-----------------------------------------------------------------------
  148. !***
  149. !*** create netCDF dataset
  150. !***
  151. ncstat = nf_create (output_file, NF_CLOBBER, nc_outfile_id)
  152. call netcdf_error_handler(ncstat)
  153. ncstat = nf_put_att_text (nc_outfile_id, NF_GLOBAL, 'title',
  154. & len_trim(map_name), map_name)
  155. call netcdf_error_handler(ncstat)
  156. !***
  157. !*** define grid size dimensions
  158. !***
  159. allocate( nc_grid1size_id(grid1_rank),
  160. & nc_grid2size_id(grid2_rank))
  161. do n=1,grid1_rank
  162. write(dim_name,1000) 'grid1_dim',n
  163. ncstat = nf_def_dim (nc_outfile_id, dim_name,
  164. & grid1_dims(n), nc_grid1size_id(n))
  165. call netcdf_error_handler(ncstat)
  166. end do
  167. do n=1,grid2_rank
  168. write(dim_name,1000) 'grid2_dim',n
  169. ncstat = nf_def_dim (nc_outfile_id, dim_name,
  170. & grid2_dims(n), nc_grid2size_id(n))
  171. call netcdf_error_handler(ncstat)
  172. end do
  173. 1000 format(a9,i1)
  174. !***
  175. !*** define grid center latitude array
  176. !***
  177. ncstat = nf_def_var (nc_outfile_id, 'src_grid_center_lat',
  178. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  179. & nc_srcgrdcntrlat_id)
  180. call netcdf_error_handler(ncstat)
  181. ncstat = nf_put_att_text (nc_outfile_id, nc_srcgrdcntrlat_id,
  182. & 'units', 7, 'radians')
  183. call netcdf_error_handler(ncstat)
  184. ncstat = nf_def_var (nc_outfile_id, 'dst_grid_center_lat',
  185. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  186. & nc_dstgrdcntrlat_id)
  187. call netcdf_error_handler(ncstat)
  188. ncstat = nf_put_att_text (nc_outfile_id, nc_dstgrdcntrlat_id,
  189. & 'units', 7, 'radians')
  190. call netcdf_error_handler(ncstat)
  191. !***
  192. !*** define grid center longitude array
  193. !***
  194. ncstat = nf_def_var (nc_outfile_id, 'src_grid_center_lon',
  195. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  196. & nc_srcgrdcntrlon_id)
  197. call netcdf_error_handler(ncstat)
  198. ncstat = nf_put_att_text (nc_outfile_id, nc_srcgrdcntrlon_id,
  199. & 'units', 7, 'radians')
  200. call netcdf_error_handler(ncstat)
  201. ncstat = nf_def_var (nc_outfile_id, 'dst_grid_center_lon',
  202. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  203. & nc_dstgrdcntrlon_id)
  204. call netcdf_error_handler(ncstat)
  205. ncstat = nf_put_att_text (nc_outfile_id, nc_dstgrdcntrlon_id,
  206. & 'units', 7, 'radians')
  207. call netcdf_error_handler(ncstat)
  208. !***
  209. !*** define grid mask
  210. !***
  211. ncstat = nf_def_var (nc_outfile_id, 'src_grid_imask', NF_INT,
  212. & grid1_rank, nc_grid1size_id, nc_srcgrdimask_id)
  213. call netcdf_error_handler(ncstat)
  214. ncstat = nf_put_att_text (nc_outfile_id, nc_srcgrdimask_id,
  215. & 'units', 8, 'unitless')
  216. call netcdf_error_handler(ncstat)
  217. ncstat = nf_def_var (nc_outfile_id, 'dst_grid_imask', NF_INT,
  218. & grid2_rank, nc_grid2size_id, nc_dstgrdimask_id)
  219. call netcdf_error_handler(ncstat)
  220. ncstat = nf_put_att_text (nc_outfile_id, nc_dstgrdimask_id,
  221. & 'units', 8, 'unitless')
  222. call netcdf_error_handler(ncstat)
  223. !***
  224. !*** define grid area arrays
  225. !***
  226. ncstat = nf_def_var (nc_outfile_id, 'src_grid_area',
  227. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  228. & nc_srcgrdarea_id)
  229. call netcdf_error_handler(ncstat)
  230. ncstat = nf_def_var (nc_outfile_id, 'dst_grid_area',
  231. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  232. & nc_dstgrdarea_id)
  233. call netcdf_error_handler(ncstat)
  234. !***
  235. !*** define grid fraction arrays
  236. !***
  237. ncstat = nf_def_var (nc_outfile_id, 'src_grid_frac',
  238. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  239. & nc_srcgrdfrac_id)
  240. call netcdf_error_handler(ncstat)
  241. ncstat = nf_def_var (nc_outfile_id, 'dst_grid_frac',
  242. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  243. & nc_dstgrdfrac_id)
  244. call netcdf_error_handler(ncstat)
  245. !***
  246. !*** define source array
  247. !***
  248. ncstat = nf_def_var (nc_outfile_id, 'src_array',
  249. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  250. & nc_srcarray_id)
  251. call netcdf_error_handler(ncstat)
  252. !***
  253. !*** define gradient arrays
  254. !***
  255. ncstat = nf_def_var (nc_outfile_id, 'src_grad_lat',
  256. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  257. & nc_srcgradlat_id)
  258. call netcdf_error_handler(ncstat)
  259. ncstat = nf_def_var (nc_outfile_id, 'src_grad_lon',
  260. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  261. & nc_srcgradlon_id)
  262. call netcdf_error_handler(ncstat)
  263. !***
  264. !*** define destination arrays
  265. !***
  266. ncstat = nf_def_var (nc_outfile_id, 'dst_array1',
  267. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  268. & nc_dstarray1_id)
  269. call netcdf_error_handler(ncstat)
  270. ncstat = nf_def_var (nc_outfile_id, 'dst_array1a',
  271. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  272. & nc_dstarray1a_id)
  273. call netcdf_error_handler(ncstat)
  274. ncstat = nf_def_var (nc_outfile_id, 'dst_array2',
  275. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  276. & nc_dstarray2_id)
  277. call netcdf_error_handler(ncstat)
  278. !***
  279. !*** define error arrays
  280. !***
  281. ncstat = nf_def_var (nc_outfile_id, 'dst_error1',
  282. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  283. & nc_dsterror1_id)
  284. call netcdf_error_handler(ncstat)
  285. ncstat = nf_def_var (nc_outfile_id, 'dst_error1a',
  286. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  287. & nc_dsterror1a_id)
  288. call netcdf_error_handler(ncstat)
  289. ncstat = nf_def_var (nc_outfile_id, 'dst_error2',
  290. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  291. & nc_dsterror2_id)
  292. call netcdf_error_handler(ncstat)
  293. !***
  294. !*** end definition stage
  295. !***
  296. ncstat = nf_enddef(nc_outfile_id)
  297. call netcdf_error_handler(ncstat)
  298. !-----------------------------------------------------------------------
  299. !
  300. ! write some grid info
  301. !
  302. !-----------------------------------------------------------------------
  303. !***
  304. !*** write grid center latitude array
  305. !***
  306. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgrdcntrlat_id,
  307. & grid1_center_lat)
  308. call netcdf_error_handler(ncstat)
  309. ncstat = nf_put_var_double(nc_outfile_id, nc_dstgrdcntrlat_id,
  310. & grid2_center_lat)
  311. call netcdf_error_handler(ncstat)
  312. !***
  313. !*** write grid center longitude array
  314. !***
  315. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgrdcntrlon_id,
  316. & grid1_center_lon)
  317. call netcdf_error_handler(ncstat)
  318. ncstat = nf_put_var_double(nc_outfile_id, nc_dstgrdcntrlon_id,
  319. & grid2_center_lon)
  320. call netcdf_error_handler(ncstat)
  321. !***
  322. !*** write grid mask
  323. !***
  324. ncstat = nf_put_var_int(nc_outfile_id, nc_srcgrdimask_id,
  325. & grid1_imask)
  326. call netcdf_error_handler(ncstat)
  327. ncstat = nf_put_var_int(nc_outfile_id, nc_dstgrdimask_id,
  328. & grid2_imask)
  329. call netcdf_error_handler(ncstat)
  330. !***
  331. !*** define grid area arrays
  332. !***
  333. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgrdarea_id,
  334. & grid1_area)
  335. call netcdf_error_handler(ncstat)
  336. ncstat = nf_put_var_double(nc_outfile_id, nc_dstgrdarea_id,
  337. & grid2_area)
  338. call netcdf_error_handler(ncstat)
  339. !***
  340. !*** define grid fraction arrays
  341. !***
  342. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgrdfrac_id,
  343. & grid1_frac)
  344. call netcdf_error_handler(ncstat)
  345. ncstat = nf_put_var_double(nc_outfile_id, nc_dstgrdfrac_id,
  346. & grid2_frac)
  347. call netcdf_error_handler(ncstat)
  348. !-----------------------------------------------------------------------
  349. !
  350. ! set up fields for test cases based on user choice
  351. !
  352. !-----------------------------------------------------------------------
  353. select case (field_choice)
  354. case(1) !*** cosine hill at lon=pi and lat=0
  355. length = 0.1*pi2
  356. grid1_array = cos(grid1_center_lat)*cos(grid1_center_lon)
  357. grid2_array = cos(grid2_center_lat)*cos(grid2_center_lon)
  358. grid1_tmp = acos(-grid1_array)/length
  359. grid2_tmp = acos(-grid2_array)/length
  360. where (grid1_tmp <= one)
  361. grad1_lat = (pi/length)*sin(pi*grid1_tmp)*
  362. & sin(grid1_center_lat)*cos(grid1_center_lon)/
  363. & sqrt(one-grid1_array**2)
  364. grad1_lon = (pi/length)*sin(pi*grid1_tmp)*
  365. & sin(grid1_center_lon)/
  366. & sqrt(one-grid1_array**2)
  367. grid1_array = two + cos(pi*grid1_tmp)
  368. elsewhere
  369. grid1_array = one
  370. grad1_lat = zero
  371. grad1_lon = zero
  372. endwhere
  373. where (grid2_tmp <= one)
  374. grid2_array = two + cos(pi*grid2_tmp)
  375. elsewhere
  376. grid2_array = one
  377. endwhere
  378. where (.not. grid1_mask)
  379. grid1_array = zero
  380. grad1_lat = zero
  381. grad1_lon = zero
  382. end where
  383. where (grid2_frac < .001) grid2_array = zero
  384. case(2) !*** pseudo-spherical harmonic l=2,m=2
  385. where (grid1_mask)
  386. grid1_array = two + cos(grid1_center_lat)**2*
  387. & cos(two*grid1_center_lon)
  388. grad1_lat = -sin(two*grid1_center_lat)*
  389. & cos(two*grid1_center_lon)
  390. grad1_lon = -two*cos(grid1_center_lat)*
  391. & sin(two*grid1_center_lon)
  392. elsewhere
  393. grid1_array = zero
  394. grad1_lat = zero
  395. grad1_lon = zero
  396. end where
  397. where (grid2_frac > .001)
  398. grid2_array = two + cos(grid2_center_lat)**2*
  399. & cos(two*grid2_center_lon)
  400. elsewhere
  401. grid2_array = zero
  402. end where
  403. case(3) !*** pseudo-spherical harmonic l=32, m=16
  404. where (grid1_mask)
  405. grid1_array = two + sin(two*grid1_center_lat)**16*
  406. & cos(16.*grid1_center_lon)
  407. grad1_lat = 32.*sin(two*grid1_center_lat)**15*
  408. & cos(two*grid1_center_lat)*
  409. & cos(16.*grid1_center_lon)
  410. grad1_lon = -32.*sin(two*grid1_center_lat)**15*
  411. & sin(grid1_center_lat)*
  412. & sin(16.*grid1_center_lon)
  413. elsewhere
  414. grid1_array = zero
  415. grad1_lat = zero
  416. grad1_lon = zero
  417. end where
  418. where (grid2_frac > .001)
  419. grid2_array = two + sin(two*grid2_center_lat)**16*
  420. & cos(16.*grid2_center_lon)
  421. elsewhere
  422. grid2_array = zero
  423. end where
  424. case default
  425. stop 'Bad choice for field to interpolate'
  426. end select
  427. !-----------------------------------------------------------------------
  428. !
  429. ! if bicubic, we need 3 gradients in logical space
  430. !
  431. !-----------------------------------------------------------------------
  432. if (map_type == map_type_bicubic) then
  433. allocate (grad1_latlon (grid1_size))
  434. nx = grid1_dims(1)
  435. ny = grid1_dims(2)
  436. do n=1,grid1_size
  437. grad1_lat(n) = zero
  438. grad1_lon(n) = zero
  439. grad1_latlon(n) = zero
  440. if (grid1_mask(n)) then
  441. delew = half
  442. delns = half
  443. j = (n-1)/nx + 1
  444. i = n - (j-1)*nx
  445. ip1 = i+1
  446. im1 = i-1
  447. jp1 = j+1
  448. jm1 = j-1
  449. if (ip1 > nx) ip1 = ip1 - nx
  450. if (im1 < 1 ) im1 = nx
  451. if (jp1 > ny) then
  452. jp1 = j
  453. delns = one
  454. endif
  455. if (jm1 < 1 ) then
  456. jm1 = j
  457. delns = one
  458. endif
  459. in = (jp1-1)*nx + i
  460. is = (jm1-1)*nx + i
  461. ie = (j -1)*nx + ip1
  462. iw = (j -1)*nx + im1
  463. ine = (jp1-1)*nx + ip1
  464. inw = (jp1-1)*nx + im1
  465. ise = (jm1-1)*nx + ip1
  466. isw = (jm1-1)*nx + im1
  467. !*** compute i-gradient
  468. if (.not. grid1_mask(ie)) then
  469. ie = n
  470. delew = one
  471. endif
  472. if (.not. grid1_mask(iw)) then
  473. iw = n
  474. delew = one
  475. endif
  476. grad1_lat(n) = delew*(grid1_array(ie) - grid1_array(iw))
  477. !*** compute j-gradient
  478. if (.not. grid1_mask(in)) then
  479. in = n
  480. delns = one
  481. endif
  482. if (.not. grid1_mask(is)) then
  483. is = n
  484. delns = one
  485. endif
  486. grad1_lon(n) = delns*(grid1_array(in) - grid1_array(is))
  487. !*** compute ij-gradient
  488. delew = half
  489. if (jp1 == j .or. jm1 == j) then
  490. delns = one
  491. else
  492. delns = half
  493. endif
  494. if (.not. grid1_mask(ine)) then
  495. if (in /= n) then
  496. ine = in
  497. delew = one
  498. else if (ie /= n) then
  499. ine = ie
  500. inw = iw
  501. if (inw == n) delew = one
  502. delns = one
  503. else
  504. ine = n
  505. inw = iw
  506. delew = one
  507. delns = one
  508. endif
  509. endif
  510. if (.not. grid1_mask(inw)) then
  511. if (in /= n) then
  512. inw = in
  513. delew = one
  514. else if (iw /= n) then
  515. inw = iw
  516. ine = ie
  517. if (ie == n) delew = one
  518. delns = one
  519. else
  520. inw = n
  521. ine = ie
  522. delew = one
  523. delns = one
  524. endif
  525. endif
  526. grad1_lat_zero(n) = delew*(grid1_array(ine) -
  527. & grid1_array(inw))
  528. if (.not. grid1_mask(ise)) then
  529. if (is /= n) then
  530. ise = is
  531. delew = one
  532. else if (ie /= n) then
  533. ise = ie
  534. isw = iw
  535. if (isw == n) delew = one
  536. delns = one
  537. else
  538. ise = n
  539. isw = iw
  540. delew = one
  541. delns = one
  542. endif
  543. endif
  544. if (.not. grid1_mask(isw)) then
  545. if (is /= n) then
  546. isw = is
  547. delew = one
  548. else if (iw /= n) then
  549. isw = iw
  550. ise = ie
  551. if (ie == n) delew = one
  552. delns = one
  553. else
  554. isw = n
  555. ise = ie
  556. delew = one
  557. delns = one
  558. endif
  559. endif
  560. grad1_lon_zero(n) = delew*(grid1_array(ise) -
  561. & grid1_array(isw))
  562. grad1_latlon(n) = delns*(grad1_lat_zero(n) -
  563. & grad1_lon_zero(n))
  564. endif
  565. enddo
  566. endif
  567. !-----------------------------------------------------------------------
  568. !
  569. ! test a first-order map from grid1 to grid2
  570. !
  571. !-----------------------------------------------------------------------
  572. grad1_lat_zero = zero
  573. grad1_lon_zero = zero
  574. if (map_type /= map_type_bicubic) then
  575. call remap(grid2_tmp, wts_map1, grid2_add_map1, grid1_add_map1,
  576. & grid1_array)
  577. else
  578. call remap(grid2_tmp, wts_map1, grid2_add_map1, grid1_add_map1,
  579. & grid1_array, src_grad1=grad1_lat,
  580. & src_grad2=grad1_lon,
  581. & src_grad3=grad1_latlon)
  582. endif
  583. if (map_type == map_type_conserv) then
  584. select case (norm_opt)
  585. case (norm_opt_none)
  586. grid2_err = grid2_frac*grid2_area
  587. where (grid2_err /= zero)
  588. grid2_tmp = grid2_tmp/grid2_err
  589. else where
  590. grid2_tmp = zero
  591. end where
  592. case (norm_opt_frcarea)
  593. case (norm_opt_dstarea)
  594. where (grid2_frac /= zero)
  595. grid2_tmp = grid2_tmp/grid2_frac
  596. else where
  597. grid2_tmp = zero
  598. end where
  599. end select
  600. end if
  601. where (grid2_frac > .999)
  602. grid2_err = (grid2_tmp - grid2_array)/grid2_array
  603. elsewhere
  604. grid2_err = zero
  605. end where
  606. print *,'First order mapping from grid1 to grid2:'
  607. print *,'----------------------------------------'
  608. print *,'Grid1 min,max: ',minval(grid1_array),maxval(grid1_array)
  609. print *,'Grid2 min,max: ',minval(grid2_tmp ),maxval(grid2_tmp )
  610. print *,' Err2 min,max: ',minval(grid2_err),maxval(grid2_err)
  611. print *,' Err2 mean: ',sum(abs(grid2_err))/
  612. & count(grid2_frac > .999)
  613. !***
  614. !*** Conservation Test
  615. !***
  616. print *,'Conservation:'
  617. print *,'Grid1 Integral = ',sum(grid1_array*grid1_area*grid1_frac)
  618. print *,'Grid2 Integral = ',sum(grid2_tmp *grid2_area*grid2_frac)
  619. !-----------------------------------------------------------------------
  620. !
  621. ! write results to NetCDF file
  622. !
  623. !-----------------------------------------------------------------------
  624. ncstat = nf_put_var_double(nc_outfile_id, nc_srcarray_id,
  625. & grid1_array)
  626. call netcdf_error_handler(ncstat)
  627. ncstat = nf_put_var_double(nc_outfile_id, nc_dstarray1_id,
  628. & grid2_tmp )
  629. call netcdf_error_handler(ncstat)
  630. ncstat = nf_put_var_double(nc_outfile_id, nc_dsterror1_id,
  631. & grid2_err)
  632. call netcdf_error_handler(ncstat)
  633. !-----------------------------------------------------------------------
  634. !
  635. ! for conservative mappings:
  636. ! test a second-order map from grid1 to grid2 with only lat grads
  637. !
  638. !-----------------------------------------------------------------------
  639. if (map_type == map_type_conserv) then
  640. call remap(grid2_tmp, wts_map1, grid2_add_map1, grid1_add_map1,
  641. & grid1_array, src_grad1=grad1_lat,
  642. & src_grad2=grad1_lon_zero)
  643. select case (norm_opt)
  644. case (norm_opt_none)
  645. grid2_err = grid2_frac*grid2_area
  646. where (grid2_err /= zero)
  647. grid2_tmp = grid2_tmp/grid2_err
  648. else where
  649. grid2_tmp = zero
  650. end where
  651. case (norm_opt_frcarea)
  652. case (norm_opt_dstarea)
  653. where (grid2_frac /= zero)
  654. grid2_tmp = grid2_tmp/grid2_frac
  655. else where
  656. grid2_tmp = zero
  657. end where
  658. end select
  659. where (grid2_frac > .999)
  660. grid2_err = (grid2_tmp - grid2_array)/grid2_array
  661. elsewhere
  662. grid2_err = zero
  663. end where
  664. print *,'Second order mapping from grid1 to grid2 (lat only):'
  665. print *,'----------------------------------------'
  666. print *,'Grid1 min,max: ',minval(grid1_array),
  667. & maxval(grid1_array)
  668. print *,'Grid2 min,max: ',minval(grid2_tmp ),
  669. & maxval(grid2_tmp )
  670. print *,' Err2 min,max: ',minval(grid2_err),maxval(grid2_err)
  671. print *,' Err2 mean: ',sum(abs(grid2_err))/
  672. & count(grid2_frac > .999)
  673. !***
  674. !*** Conservation Test
  675. !***
  676. print *,'Conservation:'
  677. print *,'Grid1 Integral = ',
  678. & sum(grid1_array*grid1_area*grid1_frac)
  679. print *,'Grid2 Integral = ',
  680. & sum(grid2_tmp *grid2_area*grid2_frac)
  681. !-----------------------------------------------------------------------
  682. !
  683. ! write results to NetCDF file
  684. !
  685. !-----------------------------------------------------------------------
  686. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgradlat_id,
  687. & grad1_lat)
  688. call netcdf_error_handler(ncstat)
  689. ncstat = nf_put_var_double(nc_outfile_id, nc_dstarray1a_id,
  690. & grid2_tmp )
  691. call netcdf_error_handler(ncstat)
  692. ncstat = nf_put_var_double(nc_outfile_id, nc_dsterror1a_id,
  693. & grid2_err)
  694. call netcdf_error_handler(ncstat)
  695. !-----------------------------------------------------------------------
  696. !
  697. ! for conservative mappings:
  698. ! test a second-order map from grid1 to grid2
  699. !
  700. !-----------------------------------------------------------------------
  701. call remap(grid2_tmp,wts_map1,grid2_add_map1,grid1_add_map1,
  702. & grid1_array, src_grad1=grad1_lat,
  703. & src_grad2=grad1_lon)
  704. select case (norm_opt)
  705. case (norm_opt_none)
  706. grid2_err = grid2_frac*grid2_area
  707. where (grid2_err /= zero)
  708. grid2_tmp = grid2_tmp/grid2_err
  709. else where
  710. grid2_tmp = zero
  711. end where
  712. case (norm_opt_frcarea)
  713. case (norm_opt_dstarea)
  714. where (grid2_frac /= zero)
  715. grid2_tmp = grid2_tmp/grid2_frac
  716. else where
  717. grid2_tmp = zero
  718. end where
  719. end select
  720. where (grid2_frac > .999)
  721. grid2_err = (grid2_tmp - grid2_array)/grid2_array
  722. elsewhere
  723. grid2_err = zero
  724. end where
  725. print *,'Second order mapping from grid1 to grid2:'
  726. print *,'-----------------------------------------'
  727. print *,'Grid1 min,max: ',minval(grid1_array),
  728. & maxval(grid1_array)
  729. print *,'Grid2 min,max: ',minval(grid2_tmp ),
  730. & maxval(grid2_tmp )
  731. print *,' Err2 min,max: ',minval(grid2_err),maxval(grid2_err)
  732. print *,' Err2 mean: ',sum(abs(grid2_err))/
  733. & count(grid2_frac > .999)
  734. !***
  735. !*** Conservation Test
  736. !***
  737. print *,'Conservation:'
  738. print *,'Grid1 Integral = ',
  739. & sum(grid1_array*grid1_area*grid1_frac)
  740. print *,'Grid2 Integral = ',
  741. & sum(grid2_tmp *grid2_area*grid2_frac)
  742. !-----------------------------------------------------------------------
  743. !
  744. ! write results to NetCDF file
  745. !
  746. !-----------------------------------------------------------------------
  747. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgradlon_id,
  748. & grad1_lon)
  749. call netcdf_error_handler(ncstat)
  750. ncstat = nf_put_var_double(nc_outfile_id, nc_dstarray2_id,
  751. & grid2_tmp )
  752. call netcdf_error_handler(ncstat)
  753. ncstat = nf_put_var_double(nc_outfile_id, nc_dsterror2_id,
  754. & grid2_err)
  755. call netcdf_error_handler(ncstat)
  756. endif
  757. !-----------------------------------------------------------------------
  758. !
  759. ! close netCDF file
  760. !
  761. !-----------------------------------------------------------------------
  762. ncstat = nf_close(nc_outfile_id)
  763. call netcdf_error_handler(ncstat)
  764. !-----------------------------------------------------------------------
  765. !
  766. ! calculate some statistics
  767. !
  768. !-----------------------------------------------------------------------
  769. grid2_count = zero
  770. grid2_tmp = zero
  771. grid2_err = zero
  772. print *,'number of sparse matrix entries ',num_links_map1
  773. do n=1,num_links_map1
  774. grid2_count(grid2_add_map1(n)) =
  775. & grid2_count(grid2_add_map1(n)) + 1
  776. if (wts_map1(1,n) > one .or. wts_map1(1,n) < zero) then
  777. grid2_tmp(grid2_add_map1(n)) =
  778. & grid2_tmp(grid2_add_map1(n)) + 1
  779. grid2_err(grid2_add_map1(n)) = max(abs(wts_map1(1,n)),
  780. & grid2_err(grid2_add_map1(n)) )
  781. endif
  782. end do
  783. do n=1,grid2_size
  784. if (grid2_tmp(n) > zero) print *,n,grid2_err(n)
  785. end do
  786. imin = minval(grid2_count, mask=(grid2_count > 0))
  787. imax = maxval(grid2_count)
  788. idiff = (imax - imin)/10 + 1
  789. print *,'total number of dest cells ',grid2_size
  790. print *,'number of cells participating in remap ',
  791. & count(grid2_count > zero)
  792. print *,'min no of entries/row = ',imin
  793. print *,'max no of entries/row = ',imax
  794. imax = imin + idiff
  795. do n=1,10
  796. print *,'num of rows with entries between ',imin,' - ',imax-1,
  797. & count(grid2_count >= imin .and. grid2_count < imax)
  798. imin = imin + idiff
  799. imax = imax + idiff
  800. end do
  801. !-----------------------------------------------------------------------
  802. end program remap_test
  803. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!