remap_vars.F 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439
  1. C****
  2. C ************************
  3. C * OASIS MODULE *
  4. C * ------------ *
  5. C ************************
  6. C****
  7. C***********************************************************************
  8. C This module belongs to the SCRIP library. It is modified to run
  9. C within OASIS.
  10. C Modifications:
  11. C - introduction of logical flags to allow multiple calls of SCRIP
  12. C - deallocation of arrays not needed any more
  13. C - added CASE for GAUSWGT
  14. C
  15. C Modified by V. Gayler, M&D 20.09.2001
  16. C Modified by D. Declat, CERFACS 27.06.2002
  17. C***********************************************************************
  18. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  19. !
  20. ! this module contains necessary variables for remapping between
  21. ! two grids. also routines for resizing and initializing these
  22. ! variables.
  23. !
  24. !-----------------------------------------------------------------------
  25. !
  26. ! CVS:$Id: remap_vars.F 2826 2010-12-10 11:14:21Z valcke $
  27. !
  28. ! Copyright (c) 1997, 1998 the Regents of the University of
  29. ! California.
  30. !
  31. ! This software and ancillary information (herein called software)
  32. ! called SCRIP is made available under the terms described here.
  33. ! The software has been approved for release with associated
  34. ! LA-CC Number 98-45.
  35. !
  36. ! Unless otherwise indicated, this software has been authored
  37. ! by an employee or employees of the University of California,
  38. ! operator of the Los Alamos National Laboratory under Contract
  39. ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
  40. ! Government has rights to use, reproduce, and distribute this
  41. ! software. The public may copy and use this software without
  42. ! charge, provided that this Notice and any statement of authorship
  43. ! are reproduced on all copies. Neither the Government nor the
  44. ! University makes any warranty, express or implied, or assumes
  45. ! any liability or responsibility for the use of this software.
  46. !
  47. ! If software is modified to produce derivative works, such modified
  48. ! software should be clearly marked, so as not to confuse it with
  49. ! the version available from Los Alamos National Laboratory.
  50. !
  51. !***********************************************************************
  52. module remap_vars
  53. use kinds_mod
  54. use constants
  55. use grids
  56. USE mod_oasis_flush
  57. implicit none
  58. !-----------------------------------------------------------------------
  59. !
  60. ! module variables
  61. !
  62. !-----------------------------------------------------------------------
  63. integer (kind=int_kind), parameter ::
  64. & norm_opt_none = 1
  65. &, norm_opt_dstarea = 2
  66. &, norm_opt_frcarea = 3
  67. &, norm_opt_nonorm = 4
  68. integer (kind=int_kind), parameter ::
  69. & map_type_conserv = 1
  70. &, map_type_bilinear = 2
  71. &, map_type_bicubic = 3
  72. &, map_type_distwgt = 4
  73. &, map_type_gauswgt = 5
  74. integer (kind=int_kind), save ::
  75. & max_links_map1 ! current size of link arrays
  76. &, num_links_map1 ! actual number of links for remapping
  77. &, max_links_map2 ! current size of link arrays
  78. &, num_links_map2 ! actual number of links for remapping
  79. &, num_maps ! num of remappings for this grid pair
  80. &, num_wts ! num of weights used in remapping
  81. &, map_type ! identifier for remapping method
  82. &, norm_opt ! option for normalization (conserv only)
  83. &, resize_increment ! default amount to increase array size
  84. integer (kind=int_kind), dimension(:), allocatable, save ::
  85. & grid1_add_map1, ! grid1 address for each link in mapping 1
  86. & grid2_add_map1, ! grid2 address for each link in mapping 1
  87. & grid1_add_map2, ! grid1 address for each link in mapping 2
  88. & grid2_add_map2 ! grid2 address for each link in mapping 2
  89. #ifdef TREAT_OVERLAY
  90. INTEGER (kind=int_kind), dimension(:), allocatable, save ::
  91. & grid1_add_repl1 ! grid1 address to use after overlap calculation
  92. #endif TREAT_OVERLAY
  93. real (kind=dbl_kind), dimension(:,:), allocatable, save ::
  94. & wts_map1, ! map weights for each link (num_wts,max_links)
  95. & wts_map2 ! map weights for each link (num_wts,max_links)
  96. logical (kind=log_kind), save :: lfracnnei = .false.
  97. logical (kind=log_kind), save :: first_conserv = .true. ! flag to
  98. ! indicate, whether scrip is called from
  99. ! oasis for the first time
  100. logical (kind=log_kind), save :: first_call = .true. ! flag used in
  101. ! remap_conserve (store_link_cnsrv)
  102. !***********************************************************************
  103. contains
  104. !***********************************************************************
  105. subroutine init_remap_vars (id_scripvoi)
  106. !-----------------------------------------------------------------------
  107. !
  108. ! this routine initializes some variables and provides an initial
  109. ! allocation of arrays (fairly large so frequent resizing
  110. ! unnecessary).
  111. !
  112. !-----------------------------------------------------------------------
  113. !
  114. ! input variables
  115. !
  116. !-----------------------------------------------------------------------
  117. INTEGER (kind=int_kind)::
  118. & id_scripvoi ! number of neighbours for DISTWGT and GAUSWGT
  119. !-----------------------------------------------------------------------
  120. !
  121. IF (nlogprt .GE. 2) THEN
  122. WRITE (UNIT = nulou,FMT = *)' '
  123. WRITE (UNIT = nulou,FMT = *)'Entering routine init_remap_vars'
  124. CALL OASIS_FLUSH_SCRIP(nulou)
  125. ENDIF
  126. !
  127. !-----------------------------------------------------------------------
  128. !
  129. ! determine the number of weights
  130. !
  131. !-----------------------------------------------------------------------
  132. select case (map_type)
  133. case(map_type_conserv)
  134. num_wts = 3
  135. case(map_type_bilinear)
  136. num_wts = 1
  137. case(map_type_bicubic)
  138. IF (restrict_type == 'REDUCED') THEN
  139. num_wts = 1
  140. ELSE
  141. num_wts = 4
  142. ENDIF
  143. case(map_type_distwgt)
  144. num_wts = 1
  145. case(map_type_gauswgt)
  146. num_wts = 1
  147. end select
  148. !-----------------------------------------------------------------------
  149. !
  150. ! initialize num_links and set max_links to four times the largest
  151. ! of the destination grid sizes initially (can be changed later).
  152. ! set a default resize increment to increase the size of link
  153. ! arrays if the number of links exceeds the initial size
  154. !
  155. !-----------------------------------------------------------------------
  156. num_links_map1 = 0
  157. select case (map_type)
  158. case(map_type_conserv)
  159. max_links_map1 = 4*grid2_size
  160. case(map_type_bilinear)
  161. max_links_map1 = 4*grid2_size
  162. case(map_type_bicubic)
  163. IF (restrict_type == 'REDUCED') THEN
  164. max_links_map1 = 16*grid2_size
  165. ELSE
  166. max_links_map1 = 4*grid2_size
  167. ENDIF
  168. case(map_type_distwgt)
  169. max_links_map1 = id_scripvoi*grid2_size
  170. case(map_type_gauswgt)
  171. max_links_map1 = id_scripvoi*grid2_size
  172. END select
  173. if (num_maps > 1) then
  174. num_links_map2 = 0
  175. max_links_map1 = max(4*grid1_size,4*grid2_size)
  176. max_links_map2 = max_links_map1
  177. endif
  178. resize_increment = 0.1*max(grid1_size,grid2_size)
  179. !-----------------------------------------------------------------------
  180. !
  181. ! allocate address and weight arrays for mapping 1
  182. !
  183. !-----------------------------------------------------------------------
  184. allocate (grid1_add_map1(max_links_map1),
  185. & grid2_add_map1(max_links_map1),
  186. & wts_map1(num_wts, max_links_map1))
  187. #ifdef TREAT_OVERLAY
  188. allocate (grid1_add_repl1(grid1_size))
  189. #endif TREAT_OVERLAY
  190. !-----------------------------------------------------------------------
  191. !
  192. ! allocate address and weight arrays for mapping 2 if necessary
  193. !
  194. !-----------------------------------------------------------------------
  195. if (num_maps > 1) then
  196. allocate (grid1_add_map2(max_links_map2),
  197. & grid2_add_map2(max_links_map2),
  198. & wts_map2(num_wts, max_links_map2))
  199. endif
  200. !-----------------------------------------------------------------------
  201. !
  202. ! initialize flag for routine store_link_cnsrv
  203. !
  204. !-----------------------------------------------------------------------
  205. first_call = .true.
  206. !-----------------------------------------------------------------------
  207. !
  208. IF (nlogprt .GE. 2) THEN
  209. WRITE (UNIT = nulou,FMT = *)' '
  210. WRITE (UNIT = nulou,FMT = *)'Leaving routine init_remap_vars'
  211. CALL OASIS_FLUSH_SCRIP(nulou)
  212. ENDIF
  213. !
  214. end subroutine init_remap_vars
  215. !***********************************************************************
  216. subroutine resize_remap_vars(nmap, increment)
  217. !-----------------------------------------------------------------------
  218. !
  219. ! this routine resizes remapping arrays by increasing(decreasing)
  220. ! the max_links by increment
  221. !
  222. !-----------------------------------------------------------------------
  223. !-----------------------------------------------------------------------
  224. !
  225. ! input variables
  226. !
  227. !-----------------------------------------------------------------------
  228. integer (kind=int_kind), intent(in) ::
  229. & nmap, ! identifies which mapping array to resize
  230. & increment ! the number of links to add(subtract) to arrays
  231. !-----------------------------------------------------------------------
  232. !
  233. ! local variables
  234. !
  235. !-----------------------------------------------------------------------
  236. integer (kind=int_kind) ::
  237. & ierr, ! error flag
  238. & mxlinks ! size of link arrays
  239. integer (kind=int_kind), dimension(:), allocatable ::
  240. & add1_tmp, ! temp array for resizing address arrays
  241. & add2_tmp ! temp array for resizing address arrays
  242. real (kind=dbl_kind), dimension(:,:), allocatable ::
  243. & wts_tmp ! temp array for resizing weight arrays
  244. !
  245. IF (nlogprt .GE. 2) THEN
  246. WRITE (UNIT = nulou,FMT = *)' '
  247. WRITE (UNIT = nulou,FMT = *)
  248. & 'Entering routine resize_remap_vars'
  249. CALL OASIS_FLUSH_SCRIP(nulou)
  250. ENDIF
  251. !
  252. !-----------------------------------------------------------------------
  253. !
  254. ! resize map 1 arrays if required.
  255. !
  256. !-----------------------------------------------------------------------
  257. select case (nmap)
  258. case(1)
  259. !***
  260. !*** allocate temporaries to hold original values
  261. !***
  262. mxlinks = size(grid1_add_map1)
  263. allocate (add1_tmp(mxlinks), add2_tmp(mxlinks),
  264. & wts_tmp(num_wts,mxlinks))
  265. add1_tmp = grid1_add_map1
  266. add2_tmp = grid2_add_map1
  267. wts_tmp = wts_map1
  268. !***
  269. !*** deallocate originals and increment max_links then
  270. !*** reallocate arrays at new size
  271. !***
  272. deallocate (grid1_add_map1, grid2_add_map1, wts_map1)
  273. max_links_map1 = mxlinks + increment
  274. allocate (grid1_add_map1(max_links_map1),
  275. & grid2_add_map1(max_links_map1),
  276. & wts_map1(num_wts,max_links_map1))
  277. !***
  278. !*** restore original values from temp arrays and
  279. !*** deallocate temps
  280. !***
  281. mxlinks = min(mxlinks, max_links_map1)
  282. grid1_add_map1(1:mxlinks) = add1_tmp (1:mxlinks)
  283. grid2_add_map1(1:mxlinks) = add2_tmp (1:mxlinks)
  284. wts_map1 (:,1:mxlinks) = wts_tmp(:,1:mxlinks)
  285. deallocate(add1_tmp, add2_tmp, wts_tmp)
  286. !-----------------------------------------------------------------------
  287. !
  288. ! resize map 2 arrays if required.
  289. !
  290. !-----------------------------------------------------------------------
  291. case(2)
  292. !***
  293. !*** allocate temporaries to hold original values
  294. !***
  295. mxlinks = size(grid1_add_map2)
  296. allocate (add1_tmp(mxlinks), add2_tmp(mxlinks),
  297. & wts_tmp(num_wts,mxlinks),stat=ierr)
  298. if (ierr .ne. 0) THEN
  299. WRITE(nulou,*) ' '
  300. WRITE(nulou,*)'error allocating temps in resize: ',ierr
  301. CALL OASIS_FLUSH_SCRIP(nulou)
  302. stop
  303. endif
  304. add1_tmp = grid1_add_map2
  305. add2_tmp = grid2_add_map2
  306. wts_tmp = wts_map2
  307. !***
  308. !*** deallocate originals and increment max_links then
  309. !*** reallocate arrays at new size
  310. !***
  311. deallocate (grid1_add_map2, grid2_add_map2, wts_map2)
  312. max_links_map2 = mxlinks + increment
  313. allocate (grid1_add_map2(max_links_map2),
  314. & grid2_add_map2(max_links_map2),
  315. & wts_map2(num_wts,max_links_map2),stat=ierr)
  316. if (ierr .ne. 0) then
  317. WRITE(nulou,*) ' '
  318. WRITE(nulou,*)'error allocating new arrays in resize: ',ierr
  319. CALL OASIS_FLUSH_SCRIP(nulou)
  320. stop
  321. endif
  322. !***
  323. !*** restore original values from temp arrays and
  324. !*** deallocate temps
  325. !***
  326. mxlinks = min(mxlinks, max_links_map2)
  327. grid1_add_map2(1:mxlinks) = add1_tmp (1:mxlinks)
  328. grid2_add_map2(1:mxlinks) = add2_tmp (1:mxlinks)
  329. wts_map2 (:,1:mxlinks) = wts_tmp(:,1:mxlinks)
  330. deallocate(add1_tmp, add2_tmp, wts_tmp)
  331. end select
  332. !
  333. IF (nlogprt .GE. 2) THEN
  334. WRITE (UNIT = nulou,FMT = *)' '
  335. WRITE (UNIT = nulou,FMT = *)
  336. & 'Leaving routine resize_remap_vars'
  337. CALL OASIS_FLUSH_SCRIP(nulou)
  338. ENDIF
  339. !
  340. !-----------------------------------------------------------------------
  341. end subroutine resize_remap_vars
  342. !-----------------------------------------------------------------------
  343. !-----------------------------------------------------------------------
  344. subroutine free_remap_vars
  345. !-----------------------------------------------------------------------
  346. !
  347. ! subroutine to deallocate allocated arrays
  348. !
  349. !-----------------------------------------------------------------------
  350. !
  351. IF (nlogprt .GE. 2) THEN
  352. WRITE (UNIT = nulou,FMT = *)' '
  353. WRITE (UNIT = nulou,FMT = *)'Entering routine free_remap_vars'
  354. CALL OASIS_FLUSH_SCRIP(nulou)
  355. ENDIF
  356. !
  357. deallocate (grid1_add_map1, grid2_add_map1, wts_map1)
  358. #ifdef TREAT_OVERLAY
  359. deallocate (grid1_add_repl1)
  360. #endif TREAT_OVERLAY
  361. if (num_maps > 1) then
  362. deallocate (grid1_add_map2, grid2_add_map2, wts_map2)
  363. endif
  364. if (map_type == map_type_conserv) then
  365. first_call = .true.
  366. first_conserv = .false.
  367. endif
  368. !
  369. IF (nlogprt .GE. 2) THEN
  370. WRITE (UNIT = nulou,FMT = *)' '
  371. WRITE (UNIT = nulou,FMT = *)'Leaving routine free_remap_vars'
  372. CALL OASIS_FLUSH_SCRIP(nulou)
  373. ENDIF
  374. !
  375. !-----------------------------------------------------------------------
  376. end subroutine free_remap_vars
  377. !***********************************************************************
  378. end module remap_vars
  379. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!