scrip.F 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582
  1. C****
  2. C *****************************
  3. C * OASIS ROUTINE - LEVEL 1 *
  4. C * ------------- ------- *
  5. C *****************************
  6. C****
  7. C***********************************************************************
  8. C This routine belongs to the SCRIP library. It is modified to run
  9. C within OASIS.
  10. C Modifications:
  11. C - routine does not read namelist but gets parameters from the
  12. C calling routine scriprmp.f
  13. C - map-method and noralize-option are written in capital letters
  14. C - routine grid_init is not called from scrip any more but was
  15. C called earlier from scriprmp
  16. C - call of two extra routines: free_grids and free_remap_vars to
  17. C allow multiple calls of SCRIP
  18. C - added case for GAUSWGT
  19. C - added 'REDUCED' case for bilinear and bicubic.
  20. C - hard coded num_maps=1 for USE in OASIS
  21. C - added lextrapdone argument
  22. C
  23. C Modified by V. Gayler, M&D 20.09.2001
  24. C Modified by D. Declat, CERFACS 27.06.2002
  25. C Modified by S. Valcke, CERFACS 27.08.2002
  26. C***********************************************************************
  27. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  28. !
  29. ! This routine is the driver for computing the addresses and weights
  30. ! for interpolating between two grids on a sphere.
  31. !
  32. !-----------------------------------------------------------------------
  33. !
  34. ! CVS:$Id: scrip.f 1831 2009-01-09 17:19:08Z valcke $
  35. !
  36. ! Copyright (c) 1997, 1998 the Regents of the University of
  37. ! California.
  38. !
  39. ! This software and ancillary information (herein called software)
  40. ! called SCRIP is made available under the terms described here.
  41. ! The software has been approved for release with associated
  42. ! LA-CC Number 98-45.
  43. !
  44. ! Unless otherwise indicated, this software has been authored
  45. ! by an employee or employees of the University of California,
  46. ! operator of the Los Alamos National Laboratory under Contract
  47. ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
  48. ! Government has rights to use, reproduce, and distribute this
  49. ! software. The public may copy and use this software without
  50. ! charge, provided that this Notice and any statement of authorship
  51. ! are reproduced on all copies. Neither the Government nor the
  52. ! University makes any warranty, express or implied, or assumes
  53. ! any liability or responsibility for the use of this software.
  54. !
  55. ! If software is modified to produce derivative works, such modified
  56. ! software should be clearly marked, so as not to confuse it with
  57. ! the version available from Los Alamos National Laboratory.
  58. !
  59. !***********************************************************************
  60. subroutine scrip
  61. $ (interp_file1, map1_name, m_method, n_opt,
  62. $ lextrapdone, rl_varmul, id_scripvoi)
  63. !-----------------------------------------------------------------------
  64. use kinds_mod ! module defining data types
  65. use constants ! module for common constants
  66. use iounits ! I/O unit manager
  67. use timers ! CPU timers
  68. use grids ! module with grid information
  69. use remap_vars ! common remapping variables
  70. use remap_conservative ! routines for conservative remap
  71. use remap_distance_weight ! routines for dist-weight remap
  72. use remap_gaussian_weight ! routines for gaus-weight remap
  73. use remap_bilinear ! routines for bilinear interp
  74. use remap_bicubic ! routines for bicubic interp
  75. use remap_bilinear_reduced ! routines for bilinear interp
  76. use remap_bicubic_reduced ! routines for bicubic interp
  77. use remap_write ! routines for remap output
  78. USE mod_oasis_flush
  79. implicit none
  80. !-----------------------------------------------------------------------
  81. !
  82. ! input variables
  83. !
  84. !-----------------------------------------------------------------------
  85. character (char_len), intent(in) ::
  86. & interp_file1, ! filename for output remap data (map1)
  87. & map1_name ! name for mapping from grid1 to grid2
  88. character*8, intent(in) ::
  89. & m_method, ! choice for mapping method
  90. & n_opt ! option for normalizing weights
  91. LOGICAL ::
  92. & lextrapdone ! logical, true if EXTRAP done on field
  93. REAL (kind=real_kind) ::
  94. & rl_varmul ! Gaussian variance (for GAUSWGT)
  95. INTEGER (kind=int_kind) ::
  96. & id_scripvoi ! number of neighbours for DISTWGT and GAUSWGT
  97. !-----------------------------------------------------------------------
  98. !
  99. ! local variables
  100. !
  101. !-----------------------------------------------------------------------
  102. integer (kind=int_kind) ::
  103. & n ! dummy counter
  104. character (char_len) ::
  105. & interp_file2, ! filename for output remap data (map2)
  106. & map2_name, ! name for mapping from grid2 to grid1
  107. & output_opt, ! option for output conventions
  108. & map_method, ! choice for mapping method
  109. & normalize_opt ! option for normalizing weights
  110. !---for the fracnnei options
  111. REAL (kind=dbl_kind),allocatable ::
  112. $weights_temp(:,:) ! temporary remapping weights
  113. INTEGER (kind=int_kind),allocatable ::
  114. $src_addr_temp(:), ! temporary remapping source addresses
  115. $dst_addr_temp(:) ! temporary remapping target addresses
  116. INTEGER (kind=int_kind) :: num_neigh !total number of Vmm
  117. logical (kind=log_kind) ::
  118. $ lnnei(grid2_size) ! flag for tricky points
  119. !-----------------------------------------------------------------------
  120. !
  121. IF (nlogprt .GE. 2) THEN
  122. WRITE (UNIT = nulou,FMT = *)' '
  123. WRITE (UNIT = nulou,FMT = *)'Entering routine scrip'
  124. CALL OASIS_FLUSH_SCRIP(nulou)
  125. ENDIF
  126. !
  127. !-----------------------------------------------------------------------
  128. !
  129. ! initialize timers
  130. !
  131. !-----------------------------------------------------------------------
  132. call timers_init
  133. do n=1,max_timers
  134. call timer_clear(n)
  135. end do
  136. !-----------------------------------------------------------------------
  137. !
  138. ! initialize variables of former SCRIP namelist
  139. !
  140. !-----------------------------------------------------------------------
  141. interp_file2 = 'unknown'
  142. map2_name = 'unknown'
  143. luse_grid1_area = .false.
  144. luse_grid2_area = .false.
  145. num_maps = 1
  146. output_opt = 'scrip'
  147. map_method = m_method
  148. normalize_opt = n_opt
  149. select case(map_method)
  150. case ('CONSERV')
  151. map_type = map_type_conserv
  152. case ('BILINEAR')
  153. map_type = map_type_bilinear
  154. case ('BICUBIC')
  155. map_type = map_type_bicubic
  156. case ('DISTWGT')
  157. map_type = map_type_distwgt
  158. case ('GAUSWGT')
  159. map_type = map_type_gauswgt
  160. case default
  161. stop 'unknown mapping method'
  162. end select
  163. SELECT CASE (normalize_opt)
  164. CASE ('FRACNNEI')
  165. lfracnnei = .true.
  166. END SELECT
  167. select case(normalize_opt(1:4))
  168. case ('NONE')
  169. norm_opt = norm_opt_none
  170. case ('FRAC')
  171. norm_opt = norm_opt_frcarea
  172. case ('DEST')
  173. norm_opt = norm_opt_dstarea
  174. CASE ('NONO')
  175. norm_opt = norm_opt_nonorm
  176. case default
  177. stop 'unknown normalization option'
  178. end select
  179. !
  180. IF (nlogprt .GE. 2) THEN
  181. WRITE (UNIT = nulou,FMT = *)' Computing remappings between: '
  182. & ,grid1_name
  183. WRITE (UNIT = nulou,FMT = *)' and '
  184. & ,grid2_name
  185. CALL OASIS_FLUSH_SCRIP(nulou)
  186. ENDIF
  187. !
  188. !-----------------------------------------------------------------------
  189. !
  190. ! initialize some remapping variables.
  191. !
  192. !-----------------------------------------------------------------------
  193. call init_remap_vars (id_scripvoi)
  194. !-----------------------------------------------------------------------
  195. !
  196. ! call appropriate interpolation setup routine based on type of
  197. ! remapping requested.
  198. !
  199. !-----------------------------------------------------------------------
  200. select case(map_type)
  201. case(map_type_conserv)
  202. call remap_conserv
  203. case(map_type_bilinear)
  204. IF (restrict_TYPE == 'REDUCED') then
  205. call remap_bilin_reduced (lextrapdone)
  206. ELSE
  207. call remap_bilin (lextrapdone)
  208. endif
  209. case(map_type_distwgt)
  210. call remap_distwgt (lextrapdone, id_scripvoi)
  211. case(map_type_gauswgt)
  212. call remap_gauswgt (lextrapdone, id_scripvoi, rl_varmul)
  213. case(map_type_bicubic)
  214. IF (restrict_TYPE == 'REDUCED') then
  215. call remap_bicub_reduced (lextrapdone)
  216. ELSE
  217. call remap_bicub (lextrapdone)
  218. endif
  219. case default
  220. stop 'Invalid Map Type'
  221. end select
  222. CALL sort_add(grid2_add_map1, grid1_add_map1, wts_map1,
  223. $ num_links_map1, num_wts)
  224. IF (lfracnnei) THEN
  225. CALL fracnnei_vmm(grid2_size, grid2_mask, num_links_map1,
  226. $ grid2_add_map1, grid1_add_map1, num_neigh, lnnei)
  227. allocate(weights_temp(num_wts,num_links_map1))
  228. allocate(src_addr_temp(num_links_map1))
  229. allocate(dst_addr_temp(num_links_map1))
  230. C* -- Store the weights, src_addr, dst_addr in temporary array
  231. DO n=1,num_links_map1
  232. weights_temp(:,n)= wts_map1(:,n)
  233. src_addr_temp(n) = grid1_add_map1(n)
  234. dst_addr_temp(n) = grid2_add_map1(n)
  235. enddo
  236. C -- Deallocate and reallocate the weights, src_addr, dst_addr
  237. deallocate(wts_map1,grid1_add_map1,grid2_add_map1)
  238. max_links_map1 = num_links_map1 + num_neigh
  239. allocate(wts_map1(num_wts,max_links_map1))
  240. allocate(grid1_add_map1(max_links_map1))
  241. allocate(grid2_add_map1(max_links_map1))
  242. CALL fracnnei(grid1_size, grid2_size, grid1_mask, grid2_mask,
  243. $ grid1_center_lon, grid1_center_lat,
  244. $ grid2_center_lon, grid2_center_lat,
  245. $ num_links_map1, num_wts, num_neigh, lnnei,
  246. $ weights_temp,src_addr_temp,dst_addr_temp,
  247. $ wts_map1, grid1_add_map1, grid2_add_map1)
  248. num_links_map1 = num_links_map1 + num_neigh
  249. deallocate(weights_temp,src_addr_temp,dst_addr_temp)
  250. lfracnnei = .false.
  251. ENDIF
  252. !
  253. #ifdef TREAT_OVERLAY
  254. !
  255. ! Change address if overlap point were found
  256. IF (map_type == 1) THEN
  257. DO n = 1, num_links_map1
  258. IF (grid1_add_map1(n) .ne. 0) THEN
  259. grid1_add_map1(n) = grid1_add_repl1(grid1_add_map1(n))
  260. ENDIF
  261. END DO
  262. ENDIF
  263. !
  264. #endif
  265. !
  266. DO n = 1, num_links_map1
  267. IF (.not. grid2_mask(grid2_add_map1(n))) wts_map1(:,n)=0.
  268. enddo
  269. !
  270. !-----------------------------------------------------------------------
  271. !
  272. ! reduce size of remapping arrays and then write remapping info
  273. ! to a file.
  274. !
  275. !-----------------------------------------------------------------------
  276. if (num_links_map1 /= max_links_map1) then
  277. call resize_remap_vars(1, num_links_map1-max_links_map1)
  278. endif
  279. if ((num_maps > 1) .and. (num_links_map2 /= max_links_map2)) then
  280. call resize_remap_vars(2, num_links_map2-max_links_map2)
  281. endif
  282. call write_remap(map1_name, map2_name,
  283. & interp_file1, interp_file2, output_opt)
  284. !-----------------------------------------------------------------------
  285. !
  286. ! deallocate allocatable arrays
  287. !
  288. !-----------------------------------------------------------------------
  289. call free_grids
  290. call free_remap_vars
  291. !
  292. IF (nlogprt .GE. 2) THEN
  293. WRITE (UNIT = nulou,FMT = *)' '
  294. WRITE (UNIT = nulou,FMT = *)'Leaving routine scrip'
  295. CALL OASIS_FLUSH_SCRIP(nulou)
  296. ENDIF
  297. !-----------------------------------------------------------------------!
  298. end subroutine scrip
  299. !
  300. subroutine sort_add(add1, add2, weights, num_links, num_wts)
  301. !-----------------------------------------------------------------------
  302. !
  303. ! this routine sorts address and weight arrays based on the
  304. ! destination address with the source address as a secondary
  305. ! sorting criterion. the method is a standard heap sort.
  306. !
  307. !-----------------------------------------------------------------------
  308. use kinds_mod ! defines common data types
  309. use constants ! defines common scalar constants
  310. USE mod_oasis_flush
  311. implicit none
  312. !-----------------------------------------------------------------------
  313. !
  314. ! Input and Output arrays
  315. !
  316. !-----------------------------------------------------------------------
  317. integer (kind=int_kind), intent(in) :: num_links, num_wts
  318. integer (kind=int_kind), intent(inout), dimension(num_links) ::
  319. & add1, ! destination address array (num_links)
  320. & add2 ! source address array
  321. real (kind=dbl_kind), intent(inout),
  322. & dimension(num_wts, num_links) ::
  323. & weights ! remapping weights (num_wts, num_links)
  324. !-----------------------------------------------------------------------
  325. !
  326. ! local variables
  327. !
  328. !-----------------------------------------------------------------------
  329. integer (kind=int_kind) ::
  330. ! & num_links, ! num of links for this mapping
  331. ! & num_wts, ! num of weights for this mapping
  332. & add1_tmp, add2_tmp, ! temp for addresses during swap
  333. & nwgt,
  334. & lvl, final_lvl, ! level indexes for heap sort levels
  335. & chk_lvl1, chk_lvl2, max_lvl
  336. real (kind=dbl_kind), dimension(SIZE(weights,DIM=1)) ::
  337. & wgttmp ! temp for holding wts during swap
  338. !-----------------------------------------------------------------------
  339. !
  340. IF (nlogprt .GE. 2) THEN
  341. WRITE (UNIT = nulou,FMT = *)' '
  342. WRITE (UNIT = nulou,FMT = *)'Entering routine sort_add'
  343. CALL OASIS_FLUSH_SCRIP(nulou)
  344. ENDIF
  345. !
  346. !-----------------------------------------------------------------------
  347. !
  348. ! determine total number of links to sort and number of weights
  349. !
  350. !-----------------------------------------------------------------------
  351. ! num_links = SIZE(add1)
  352. ! num_wts = SIZE(weights, DIM=1)
  353. !-----------------------------------------------------------------------
  354. !
  355. ! start at the lowest level (N/2) of the tree and sift lower
  356. ! values to the bottom of the tree, promoting the larger numbers
  357. !
  358. !-----------------------------------------------------------------------
  359. do lvl=num_links/2,1,-1
  360. final_lvl = lvl
  361. add1_tmp = add1(lvl)
  362. add2_tmp = add2(lvl)
  363. wgttmp(:) = weights(:,lvl)
  364. !***
  365. !*** loop until proper level is found for this link, or reach
  366. !*** bottom
  367. !***
  368. sift_loop1: do
  369. !***
  370. !*** find the largest of the two daughters
  371. !***
  372. chk_lvl1 = 2*final_lvl
  373. chk_lvl2 = 2*final_lvl+1
  374. if (chk_lvl1 .EQ. num_links) chk_lvl2 = chk_lvl1
  375. if ((add1(chk_lvl1) > add1(chk_lvl2)) .OR.
  376. & ((add1(chk_lvl1) == add1(chk_lvl2)) .AND.
  377. & (add2(chk_lvl1) > add2(chk_lvl2)))) then
  378. max_lvl = chk_lvl1
  379. else
  380. max_lvl = chk_lvl2
  381. endif
  382. !***
  383. !*** if the parent is greater than both daughters,
  384. !*** the correct level has been found
  385. !***
  386. if ((add1_tmp .GT. add1(max_lvl)) .OR.
  387. & ((add1_tmp .EQ. add1(max_lvl)) .AND.
  388. & (add2_tmp .GT. add2(max_lvl)))) then
  389. add1(final_lvl) = add1_tmp
  390. add2(final_lvl) = add2_tmp
  391. weights(:,final_lvl) = wgttmp(:)
  392. exit sift_loop1
  393. !***
  394. !*** otherwise, promote the largest daughter and push
  395. !*** down one level in the tree. if haven't reached
  396. !*** the end of the tree, repeat the process. otherwise
  397. !*** store last values and exit the loop
  398. !***
  399. else
  400. add1(final_lvl) = add1(max_lvl)
  401. add2(final_lvl) = add2(max_lvl)
  402. weights(:,final_lvl) = weights(:,max_lvl)
  403. final_lvl = max_lvl
  404. if (2*final_lvl > num_links) then
  405. add1(final_lvl) = add1_tmp
  406. add2(final_lvl) = add2_tmp
  407. weights(:,final_lvl) = wgttmp(:)
  408. exit sift_loop1
  409. endif
  410. endif
  411. end do sift_loop1
  412. end do
  413. !-----------------------------------------------------------------------
  414. !
  415. ! now that the heap has been sorted, strip off the top (largest)
  416. ! value and promote the values below
  417. !
  418. !-----------------------------------------------------------------------
  419. do lvl=num_links,3,-1
  420. !***
  421. !*** move the top value and insert it into the correct place
  422. !***
  423. add1_tmp = add1(lvl)
  424. add1(lvl) = add1(1)
  425. add2_tmp = add2(lvl)
  426. add2(lvl) = add2(1)
  427. wgttmp(:) = weights(:,lvl)
  428. weights(:,lvl) = weights(:,1)
  429. !***
  430. !*** as above this loop sifts the tmp values down until proper
  431. !*** level is reached
  432. !***
  433. final_lvl = 1
  434. sift_loop2: do
  435. !***
  436. !*** find the largest of the two daughters
  437. !***
  438. chk_lvl1 = 2*final_lvl
  439. chk_lvl2 = 2*final_lvl+1
  440. if (chk_lvl2 >= lvl) chk_lvl2 = chk_lvl1
  441. if ((add1(chk_lvl1) > add1(chk_lvl2)) .OR.
  442. & ((add1(chk_lvl1) == add1(chk_lvl2)) .AND.
  443. & (add2(chk_lvl1) > add2(chk_lvl2)))) then
  444. max_lvl = chk_lvl1
  445. else
  446. max_lvl = chk_lvl2
  447. endif
  448. !***
  449. !*** if the parent is greater than both daughters,
  450. !*** the correct level has been found
  451. !***
  452. if ((add1_tmp > add1(max_lvl)) .OR.
  453. & ((add1_tmp == add1(max_lvl)) .AND.
  454. & (add2_tmp > add2(max_lvl)))) then
  455. add1(final_lvl) = add1_tmp
  456. add2(final_lvl) = add2_tmp
  457. weights(:,final_lvl) = wgttmp(:)
  458. exit sift_loop2
  459. !***
  460. !*** otherwise, promote the largest daughter and push
  461. !*** down one level in the tree. if haven't reached
  462. !*** the end of the tree, repeat the process. otherwise
  463. !*** store last values and exit the loop
  464. !***
  465. else
  466. add1(final_lvl) = add1(max_lvl)
  467. add2(final_lvl) = add2(max_lvl)
  468. weights(:,final_lvl) = weights(:,max_lvl)
  469. final_lvl = max_lvl
  470. if (2*final_lvl >= lvl) then
  471. add1(final_lvl) = add1_tmp
  472. add2(final_lvl) = add2_tmp
  473. weights(:,final_lvl) = wgttmp(:)
  474. exit sift_loop2
  475. endif
  476. endif
  477. end do sift_loop2
  478. end do
  479. !***
  480. !*** swap the last two entries
  481. !***
  482. add1_tmp = add1(2)
  483. add1(2) = add1(1)
  484. add1(1) = add1_tmp
  485. add2_tmp = add2(2)
  486. add2(2) = add2(1)
  487. add2(1) = add2_tmp
  488. wgttmp (:) = weights(:,2)
  489. weights(:,2) = weights(:,1)
  490. weights(:,1) = wgttmp (:)
  491. !
  492. IF (nlogprt .GE. 2) THEN
  493. WRITE (UNIT = nulou,FMT = *)' '
  494. WRITE (UNIT = nulou,FMT = *)'Leaving routine sort_add'
  495. CALL OASIS_FLUSH_SCRIP(nulou)
  496. ENDIF
  497. !
  498. !-----------------------------------------------------------------------
  499. end subroutine sort_add
  500. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!