scrip_test_repeat.f 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709
  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. !-----------------------------------------------------------------------
  8. !
  9. ! CVS:$Id: scrip_test_repeat.f,v 1.3 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. program remap_test_repeat
  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. implicit none
  44. !-----------------------------------------------------------------------
  45. !
  46. ! interface for remap routine
  47. !
  48. !-----------------------------------------------------------------------
  49. interface
  50. subroutine remap(dst_array, map_wts, dst_add, src_add,
  51. & src_array, src_grad1, src_grad2, src_grad3)
  52. use kinds_mod
  53. use constants
  54. implicit none
  55. integer (kind=int_kind), dimension(:), intent(in) ::
  56. & dst_add, src_add
  57. real (kind=dbl_kind), dimension(:,:), intent(in) :: map_wts
  58. real (kind=dbl_kind), dimension(:), intent(in) :: src_array
  59. real (kind=dbl_kind), dimension(:), intent(in), optional ::
  60. & src_grad1, src_grad2, src_grad3
  61. real (kind=dbl_kind), dimension(:), intent(inout) :: dst_array
  62. end subroutine remap
  63. end interface
  64. !-----------------------------------------------------------------------
  65. !
  66. ! input namelist variables
  67. !
  68. !-----------------------------------------------------------------------
  69. integer (kind=int_kind) ::
  70. & field_choice, ! choice of field to be interpolated
  71. & num_repeats ! number of times to repeat remappings
  72. character (char_len) ::
  73. & interp_file1, ! filename containing remap data (map1)
  74. & interp_file2, ! filename containing remap data (map2)
  75. & output_file ! filename containing output test data
  76. namelist /remap_inputs/ field_choice, num_repeats,
  77. & interp_file1, interp_file2, output_file
  78. !-----------------------------------------------------------------------
  79. !
  80. ! local variables
  81. !
  82. !-----------------------------------------------------------------------
  83. character (char_len) ::
  84. & map_name1, ! name for mapping from grid1 to grid2
  85. & map_name2 ! name for mapping from grid2 to grid1
  86. integer (kind=int_kind) :: ! netCDF ids for files and arrays
  87. & n, ncstat, nc_outfile_id,
  88. & nc_srcgrdcntrlat_id, nc_srcgrdcntrlon_id,
  89. & nc_dstgrdcntrlat_id, nc_dstgrdcntrlon_id,
  90. & nc_srcgrdrank_id, nc_dstgrdrank_id,
  91. & nc_srcgrdimask_id, nc_dstgrdimask_id,
  92. & nc_srcgrdarea_id, nc_dstgrdarea_id,
  93. & nc_srcgrdfrac_id, nc_dstgrdfrac_id,
  94. & nc_srcarray_id, nc_srcgradlat_id, nc_srcgradlon_id,
  95. & nc_dstarray1_id, nc_dstarray2_id,
  96. & nc_dsterror1_id, nc_dsterror2_id
  97. integer (kind=int_kind), dimension(:), allocatable ::
  98. & nc_grid1size_id, nc_grid2size_id
  99. !-----------------------------------------------------------------------
  100. character (char_len) ::
  101. & dim_name ! netCDF dimension name
  102. integer (kind=int_kind) :: i,j,n,
  103. & iunit ! unit number for namelist file
  104. integer (kind=int_kind), dimension(:), allocatable ::
  105. & grid1_imask, grid2_imask
  106. real (kind=dbl_kind) ::
  107. & length ! length scale for cosine hill test field
  108. real (kind=dbl_kind), dimension(:), allocatable ::
  109. & grid1_array,
  110. & grid1_tmp,
  111. & grad1_lat,
  112. & grad1_lon,
  113. & grid2_array,
  114. & grid2_err,
  115. & grid2_tmp,
  116. & grad2_lat,
  117. & grad2_lon
  118. !-----------------------------------------------------------------------
  119. !
  120. ! read namelist for file and mapping info
  121. !
  122. !-----------------------------------------------------------------------
  123. call get_unit(iunit)
  124. open(iunit, file='repeat_test_in', status='old', form='formatted')
  125. read(iunit, nml=remap_inputs)
  126. call release_unit(iunit)
  127. write(*,nml=remap_inputs)
  128. !-----------------------------------------------------------------------
  129. !
  130. ! read remapping data
  131. !
  132. !-----------------------------------------------------------------------
  133. num_maps = 2
  134. call read_remap(map_name1, map_name2, interp_file1, interp_file2)
  135. !-----------------------------------------------------------------------
  136. !
  137. ! allocate arrays
  138. !
  139. !-----------------------------------------------------------------------
  140. allocate (grid1_array (grid1_size),
  141. & grid1_tmp (grid1_size),
  142. & grad1_lat (grid1_size),
  143. & grad1_lon (grid1_size),
  144. & grid1_imask (grid1_size),
  145. & grid2_array (grid2_size),
  146. & grid2_err (grid2_size),
  147. & grid2_tmp (grid2_size),
  148. & grad2_lat (grid2_size),
  149. & grad2_lon (grid2_size),
  150. & grid2_imask (grid2_size))
  151. where (grid1_mask)
  152. grid1_imask = 1
  153. elsewhere
  154. grid1_imask = 0
  155. endwhere
  156. where (grid2_mask)
  157. grid2_imask = 1
  158. elsewhere
  159. grid2_imask = 0
  160. endwhere
  161. !-----------------------------------------------------------------------
  162. !
  163. ! setup a NetCDF file for output
  164. !
  165. !-----------------------------------------------------------------------
  166. !***
  167. !*** create netCDF dataset
  168. !***
  169. ncstat = nf_create (output_file, NF_CLOBBER, nc_outfile_id)
  170. call netcdf_error_handler(ncstat)
  171. ncstat = nf_put_att_text (nc_outfile_id, NF_GLOBAL, 'title',
  172. & len_trim(map_name1), map_name1)
  173. call netcdf_error_handler(ncstat)
  174. !***
  175. !*** define grid size dimensions
  176. !***
  177. allocate( nc_grid1size_id(grid1_rank),
  178. & nc_grid2size_id(grid2_rank))
  179. do n=1,grid1_rank
  180. write(dim_name,1000) 'grid1_dim',n
  181. ncstat = nf_def_dim (nc_outfile_id, dim_name,
  182. & grid1_dims(n), nc_grid1size_id(n))
  183. call netcdf_error_handler(ncstat)
  184. end do
  185. do n=1,grid2_rank
  186. write(dim_name,1000) 'grid2_dim',n
  187. ncstat = nf_def_dim (nc_outfile_id, dim_name,
  188. & grid2_dims(n), nc_grid2size_id(n))
  189. call netcdf_error_handler(ncstat)
  190. end do
  191. 1000 format(a9,i1)
  192. !***
  193. !*** define grid center latitude array
  194. !***
  195. ncstat = nf_def_var (nc_outfile_id, 'src_grid_center_lat',
  196. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  197. & nc_srcgrdcntrlat_id)
  198. call netcdf_error_handler(ncstat)
  199. ncstat = nf_put_att_text (nc_outfile_id, nc_srcgrdcntrlat_id,
  200. & 'units', 7, 'radians')
  201. call netcdf_error_handler(ncstat)
  202. ncstat = nf_def_var (nc_outfile_id, 'dst_grid_center_lat',
  203. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  204. & nc_dstgrdcntrlat_id)
  205. call netcdf_error_handler(ncstat)
  206. ncstat = nf_put_att_text (nc_outfile_id, nc_dstgrdcntrlat_id,
  207. & 'units', 7, 'radians')
  208. call netcdf_error_handler(ncstat)
  209. !***
  210. !*** define grid center longitude array
  211. !***
  212. ncstat = nf_def_var (nc_outfile_id, 'src_grid_center_lon',
  213. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  214. & nc_srcgrdcntrlon_id)
  215. call netcdf_error_handler(ncstat)
  216. ncstat = nf_put_att_text (nc_outfile_id, nc_srcgrdcntrlon_id,
  217. & 'units', 7, 'radians')
  218. call netcdf_error_handler(ncstat)
  219. ncstat = nf_def_var (nc_outfile_id, 'dst_grid_center_lon',
  220. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  221. & nc_dstgrdcntrlon_id)
  222. call netcdf_error_handler(ncstat)
  223. ncstat = nf_put_att_text (nc_outfile_id, nc_dstgrdcntrlon_id,
  224. & 'units', 7, 'radians')
  225. call netcdf_error_handler(ncstat)
  226. !***
  227. !*** define grid mask
  228. !***
  229. ncstat = nf_def_var (nc_outfile_id, 'src_grid_imask', NF_INT,
  230. & grid1_rank, nc_grid1size_id, nc_srcgrdimask_id)
  231. call netcdf_error_handler(ncstat)
  232. ncstat = nf_put_att_text (nc_outfile_id, nc_srcgrdimask_id,
  233. & 'units', 8, 'unitless')
  234. call netcdf_error_handler(ncstat)
  235. ncstat = nf_def_var (nc_outfile_id, 'dst_grid_imask', NF_INT,
  236. & grid2_rank, nc_grid2size_id, nc_dstgrdimask_id)
  237. call netcdf_error_handler(ncstat)
  238. ncstat = nf_put_att_text (nc_outfile_id, nc_dstgrdimask_id,
  239. & 'units', 8, 'unitless')
  240. call netcdf_error_handler(ncstat)
  241. !***
  242. !*** define grid area arrays
  243. !***
  244. ncstat = nf_def_var (nc_outfile_id, 'src_grid_area',
  245. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  246. & nc_srcgrdarea_id)
  247. call netcdf_error_handler(ncstat)
  248. ncstat = nf_def_var (nc_outfile_id, 'dst_grid_area',
  249. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  250. & nc_dstgrdarea_id)
  251. call netcdf_error_handler(ncstat)
  252. !***
  253. !*** define grid fraction arrays
  254. !***
  255. ncstat = nf_def_var (nc_outfile_id, 'src_grid_frac',
  256. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  257. & nc_srcgrdfrac_id)
  258. call netcdf_error_handler(ncstat)
  259. ncstat = nf_def_var (nc_outfile_id, 'dst_grid_frac',
  260. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  261. & nc_dstgrdfrac_id)
  262. call netcdf_error_handler(ncstat)
  263. !***
  264. !*** define source array
  265. !***
  266. ncstat = nf_def_var (nc_outfile_id, 'src_array',
  267. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  268. & nc_srcarray_id)
  269. call netcdf_error_handler(ncstat)
  270. !***
  271. !*** define gradient arrays
  272. !***
  273. ncstat = nf_def_var (nc_outfile_id, 'src_grad_lat',
  274. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  275. & nc_srcgradlat_id)
  276. call netcdf_error_handler(ncstat)
  277. ncstat = nf_def_var (nc_outfile_id, 'src_grad_lon',
  278. & NF_DOUBLE, grid1_rank, nc_grid1size_id,
  279. & nc_srcgradlon_id)
  280. call netcdf_error_handler(ncstat)
  281. !***
  282. !*** define destination arrays
  283. !***
  284. ncstat = nf_def_var (nc_outfile_id, 'dst_array1',
  285. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  286. & nc_dstarray1_id)
  287. call netcdf_error_handler(ncstat)
  288. ncstat = nf_def_var (nc_outfile_id, 'dst_array2',
  289. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  290. & nc_dstarray2_id)
  291. call netcdf_error_handler(ncstat)
  292. !***
  293. !*** define error arrays
  294. !***
  295. ncstat = nf_def_var (nc_outfile_id, 'dst_error1',
  296. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  297. & nc_dsterror1_id)
  298. call netcdf_error_handler(ncstat)
  299. ncstat = nf_def_var (nc_outfile_id, 'dst_error2',
  300. & NF_DOUBLE, grid2_rank, nc_grid2size_id,
  301. & nc_dsterror2_id)
  302. call netcdf_error_handler(ncstat)
  303. !***
  304. !*** end definition stage
  305. !***
  306. ncstat = nf_enddef(nc_outfile_id)
  307. call netcdf_error_handler(ncstat)
  308. !-----------------------------------------------------------------------
  309. !
  310. ! write some grid info
  311. !
  312. !-----------------------------------------------------------------------
  313. !***
  314. !*** write grid center latitude array
  315. !***
  316. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgrdcntrlat_id,
  317. & grid1_center_lat)
  318. call netcdf_error_handler(ncstat)
  319. ncstat = nf_put_var_double(nc_outfile_id, nc_dstgrdcntrlat_id,
  320. & grid2_center_lat)
  321. call netcdf_error_handler(ncstat)
  322. !***
  323. !*** write grid center longitude array
  324. !***
  325. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgrdcntrlon_id,
  326. & grid1_center_lon)
  327. call netcdf_error_handler(ncstat)
  328. ncstat = nf_put_var_double(nc_outfile_id, nc_dstgrdcntrlon_id,
  329. & grid2_center_lon)
  330. call netcdf_error_handler(ncstat)
  331. !***
  332. !*** write grid mask
  333. !***
  334. ncstat = nf_put_var_int(nc_outfile_id, nc_srcgrdimask_id,
  335. & grid1_imask)
  336. call netcdf_error_handler(ncstat)
  337. ncstat = nf_put_var_int(nc_outfile_id, nc_dstgrdimask_id,
  338. & grid2_imask)
  339. call netcdf_error_handler(ncstat)
  340. !***
  341. !*** define grid area arrays
  342. !***
  343. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgrdarea_id,
  344. & grid1_area)
  345. call netcdf_error_handler(ncstat)
  346. ncstat = nf_put_var_double(nc_outfile_id, nc_dstgrdarea_id,
  347. & grid2_area)
  348. call netcdf_error_handler(ncstat)
  349. !***
  350. !*** define grid fraction arrays
  351. !***
  352. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgrdfrac_id,
  353. & grid1_frac)
  354. call netcdf_error_handler(ncstat)
  355. ncstat = nf_put_var_double(nc_outfile_id, nc_dstgrdfrac_id,
  356. & grid2_frac)
  357. call netcdf_error_handler(ncstat)
  358. !-----------------------------------------------------------------------
  359. !
  360. ! set up fields for test cases based on user choice
  361. !
  362. !-----------------------------------------------------------------------
  363. select case (field_choice)
  364. case(1) !*** cosine hill at lon=pi and lat=0
  365. length = 0.1*pi2
  366. grid1_array = cos(grid1_center_lat)*cos(grid1_center_lon)
  367. grid2_array = cos(grid2_center_lat)*cos(grid2_center_lon)
  368. grid1_tmp = acos(-grid1_array)/length
  369. grid2_tmp = acos(-grid2_array)/length
  370. where (grid1_tmp <= one)
  371. grad1_lat = (pi/length)*sin(pi*grid1_tmp)*
  372. & sin(grid1_center_lat)*cos(grid1_center_lon)/
  373. & sqrt(one-grid1_array**2)
  374. grad1_lon = (pi/length)*sin(pi*grid1_tmp)*
  375. & sin(grid1_center_lon)/
  376. & sqrt(one-grid1_array**2)
  377. grid1_array = two + cos(pi*grid1_tmp)
  378. elsewhere
  379. grid1_array = one
  380. grad1_lat = zero
  381. grad1_lon = zero
  382. endwhere
  383. where (grid2_tmp <= one)
  384. grad2_lat = (pi/length)*sin(pi*grid2_tmp)*
  385. & sin(grid2_center_lat)*cos(grid2_center_lon)/
  386. & sqrt(one-grid2_array**2)
  387. grad2_lon = (pi/length)*sin(pi*grid2_tmp)*
  388. & sin(grid2_center_lon)/
  389. & sqrt(one-grid2_array**2)
  390. grid2_array = two + cos(pi*grid2_tmp)
  391. elsewhere
  392. grid2_array = one
  393. grad2_lat = zero
  394. grad2_lon = zero
  395. endwhere
  396. where (.not. grid1_mask)
  397. grid1_array = zero
  398. grad1_lat = zero
  399. grad1_lon = zero
  400. end where
  401. where (.not. grid2_mask)
  402. grid2_array = zero
  403. grad2_lat = zero
  404. grad2_lon = zero
  405. end where
  406. case(2) !*** pseudo-spherical harmonic l=2,m=2
  407. where (grid1_mask)
  408. grid1_array = two + cos(grid1_center_lat)**2*
  409. & cos(two*grid1_center_lon)
  410. grad1_lat = -sin(two*grid1_center_lat)*
  411. & cos(two*grid1_center_lon)
  412. grad1_lon = -two*cos(grid1_center_lat)*
  413. & sin(two*grid1_center_lon)
  414. elsewhere
  415. grid1_array = zero
  416. grad1_lat = zero
  417. grad1_lon = zero
  418. end where
  419. where (grid2_mask)
  420. grid2_array = two + cos(grid2_center_lat)**2*
  421. & cos(two*grid2_center_lon)
  422. grad2_lat = -sin(two*grid2_center_lat)*
  423. & cos(two*grid2_center_lon)
  424. grad2_lon = -two*cos(grid2_center_lat)*
  425. & sin(two*grid2_center_lon)
  426. elsewhere
  427. grid2_array = zero
  428. grad2_lat = zero
  429. grad2_lon = zero
  430. end where
  431. case(3) !*** pseudo-spherical harmonic l=32, m=16
  432. where (grid1_mask)
  433. grid1_array = two + sin(two*grid1_center_lat)**16*
  434. & cos(16.*grid1_center_lon)
  435. grad1_lat = 32.*sin(two*grid1_center_lat)**15*
  436. & cos(two*grid1_center_lat)*
  437. & cos(16.*grid1_center_lon)
  438. grad1_lon = -32.*sin(two*grid1_center_lat)**15*
  439. & sin(grid1_center_lat)*
  440. & sin(16.*grid1_center_lon)
  441. elsewhere
  442. grid1_array = zero
  443. grad1_lat = zero
  444. grad1_lon = zero
  445. end where
  446. where (grid2_mask)
  447. grid2_array = two + sin(two*grid2_center_lat)**16*
  448. & cos(16.*grid2_center_lon)
  449. grad2_lat = 32.*sin(two*grid2_center_lat)**15*
  450. & cos(two*grid2_center_lat)*
  451. & cos(16.*grid2_center_lon)
  452. grad2_lon = -32.*sin(two*grid2_center_lat)**15*
  453. & sin(grid2_center_lat)*
  454. & sin(16.*grid2_center_lon)
  455. elsewhere
  456. grid2_array = zero
  457. grad2_lat = zero
  458. grad2_lon = zero
  459. end where
  460. case default
  461. stop 'Bad choice for field to interpolate'
  462. end select
  463. !-----------------------------------------------------------------------
  464. !
  465. ! test repeated first-order maps between grid1 and grid2
  466. !
  467. !-----------------------------------------------------------------------
  468. call remap(grid2_tmp, wts_map1, grid2_add_map1, grid1_add_map1,
  469. & grid1_array)
  470. do n=1,num_repeats
  471. call remap(grid1_tmp, wts_map2, grid1_add_map2, grid2_add_map2,
  472. & grid2_tmp)
  473. call remap(grid2_tmp, wts_map1, grid2_add_map1, grid1_add_map1,
  474. & grid1_tmp)
  475. end do
  476. where (grid2_frac > .999)
  477. grid2_err = (grid2_tmp - grid2_array)/grid2_array
  478. elsewhere
  479. grid2_err = zero
  480. end where
  481. print *,'First order mapping from grid1 to grid2:'
  482. print *,'----------------------------------------'
  483. print *,'Grid1 min,max: ',minval(grid1_array),maxval(grid1_array)
  484. print *,'Grid2 min,max: ',minval(grid2_tmp ),maxval(grid2_tmp )
  485. print *,' Err2 min,max: ',minval(grid2_err),maxval(grid2_err)
  486. print *,' Err2 mean: ',sum(abs(grid2_err))/
  487. & count(grid2_frac > .001)
  488. !***
  489. !*** Conservation Test
  490. !***
  491. print *,'Conservation:'
  492. print *,'Grid1 Integral = ',sum(grid1_array*grid1_area*grid1_frac)
  493. print *,'Grid2 Integral = ',sum(grid2_tmp *grid2_area*grid2_frac)
  494. !-----------------------------------------------------------------------
  495. !
  496. ! write results to NetCDF file
  497. !
  498. !-----------------------------------------------------------------------
  499. ncstat = nf_put_var_double(nc_outfile_id, nc_srcarray_id,
  500. & grid1_array)
  501. call netcdf_error_handler(ncstat)
  502. ncstat = nf_put_var_double(nc_outfile_id, nc_dstarray1_id,
  503. & grid2_tmp )
  504. call netcdf_error_handler(ncstat)
  505. ncstat = nf_put_var_double(nc_outfile_id, nc_dsterror1_id,
  506. & grid2_err)
  507. call netcdf_error_handler(ncstat)
  508. !-----------------------------------------------------------------------
  509. !
  510. ! test repeated second-order mapppings between grid1 and grid2
  511. !
  512. !-----------------------------------------------------------------------
  513. if (num_wts > 1) then
  514. call remap(grid2_tmp , wts_map1, grid2_add_map1, grid1_add_map1,
  515. & grid1_array, src_grad1=grad1_lat,
  516. & src_grad2=grad1_lon)
  517. do n=1,num_repeats
  518. call remap(grid1_tmp, wts_map2, grid1_add_map2, grid2_add_map2,
  519. & grid2_tmp, src_grad1=grad2_lat,
  520. & src_grad2=grad2_lon)
  521. call remap(grid2_tmp, wts_map1, grid2_add_map1, grid1_add_map1,
  522. & grid1_tmp, src_grad1=grad1_lat,
  523. & src_grad2=grad1_lon)
  524. end do
  525. where (grid2_frac > .999)
  526. grid2_err = (grid2_tmp - grid2_array)/grid2_array
  527. elsewhere
  528. grid2_err = zero
  529. end where
  530. print *,'Second order mapping from grid1 to grid2:'
  531. print *,'-----------------------------------------'
  532. print *,'Grid1 min,max: ',minval(grid1_array),maxval(grid1_array)
  533. print *,'Grid2 min,max: ',minval(grid2_tmp ),maxval(grid2_tmp )
  534. print *,' Err2 min,max: ',minval(grid2_err),maxval(grid2_err)
  535. print *,' Err2 mean: ',sum(abs(grid2_err))/
  536. & count(grid2_frac > .001)
  537. !***
  538. !*** Conservation Test
  539. !***
  540. print *,'Conservation:'
  541. print *,'Grid1 Integral = ',sum(grid1_array*grid1_area*grid1_frac)
  542. print *,'Grid2 Integral = ',sum(grid2_tmp *grid2_area*grid2_frac)
  543. !-----------------------------------------------------------------------
  544. !
  545. ! write results to NetCDF file
  546. !
  547. !-----------------------------------------------------------------------
  548. ncstat = nf_put_var_double(nc_outfile_id, nc_srcgradlon_id,
  549. & grad1_lon)
  550. call netcdf_error_handler(ncstat)
  551. ncstat = nf_put_var_double(nc_outfile_id, nc_dstarray2_id,
  552. & grid2_tmp )
  553. call netcdf_error_handler(ncstat)
  554. ncstat = nf_put_var_double(nc_outfile_id, nc_dsterror2_id,
  555. & grid2_err)
  556. call netcdf_error_handler(ncstat)
  557. endif
  558. !-----------------------------------------------------------------------
  559. !
  560. ! close netCDF file
  561. !
  562. !-----------------------------------------------------------------------
  563. ncstat = nf_close(nc_outfile_id)
  564. call netcdf_error_handler(ncstat)
  565. !-----------------------------------------------------------------------
  566. end program remap_test_repeat
  567. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!