remap_vars.f 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! this module contains necessary variables for remapping between
  4. ! two grids. also routines for resizing and initializing these
  5. ! variables.
  6. !
  7. !-----------------------------------------------------------------------
  8. !
  9. ! CVS:$Id: remap_vars.f,v 1.5 2000/04/19 21:56:26 pwjones Exp $
  10. !
  11. ! Copyright (c) 1997, 1998 the Regents of the University of
  12. ! California.
  13. !
  14. ! This software and ancillary information (herein called software)
  15. ! called SCRIP is made available under the terms described here.
  16. ! The software has been approved for release with associated
  17. ! LA-CC Number 98-45.
  18. !
  19. ! Unless otherwise indicated, this software has been authored
  20. ! by an employee or employees of the University of California,
  21. ! operator of the Los Alamos National Laboratory under Contract
  22. ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
  23. ! Government has rights to use, reproduce, and distribute this
  24. ! software. The public may copy and use this software without
  25. ! charge, provided that this Notice and any statement of authorship
  26. ! are reproduced on all copies. Neither the Government nor the
  27. ! University makes any warranty, express or implied, or assumes
  28. ! any liability or responsibility for the use of this software.
  29. !
  30. ! If software is modified to produce derivative works, such modified
  31. ! software should be clearly marked, so as not to confuse it with
  32. ! the version available from Los Alamos National Laboratory.
  33. !
  34. !***********************************************************************
  35. module remap_vars
  36. use kinds_mod
  37. use constants
  38. use grids
  39. implicit none
  40. !-----------------------------------------------------------------------
  41. !
  42. ! module variables
  43. !
  44. !-----------------------------------------------------------------------
  45. integer (kind=int_kind), parameter ::
  46. & norm_opt_none = 1
  47. &, norm_opt_dstarea = 2
  48. &, norm_opt_frcarea = 3
  49. integer (kind=int_kind), parameter ::
  50. & map_type_conserv = 1
  51. &, map_type_bilinear = 2
  52. &, map_type_bicubic = 3
  53. &, map_type_distwgt = 4
  54. integer (kind=int_kind), save ::
  55. & max_links_map1 ! current size of link arrays
  56. &, num_links_map1 ! actual number of links for remapping
  57. &, max_links_map2 ! current size of link arrays
  58. &, num_links_map2 ! actual number of links for remapping
  59. &, num_maps ! num of remappings for this grid pair
  60. &, num_wts ! num of weights used in remapping
  61. &, map_type ! identifier for remapping method
  62. &, norm_opt ! option for normalization (conserv only)
  63. &, resize_increment ! default amount to increase array size
  64. integer (kind=int_kind), dimension(:), allocatable, save ::
  65. & grid1_add_map1, ! grid1 address for each link in mapping 1
  66. & grid2_add_map1, ! grid2 address for each link in mapping 1
  67. & grid1_add_map2, ! grid1 address for each link in mapping 2
  68. & grid2_add_map2 ! grid2 address for each link in mapping 2
  69. real (kind=dbl_kind), dimension(:,:), allocatable, save ::
  70. & wts_map1, ! map weights for each link (num_wts,max_links)
  71. & wts_map2 ! map weights for each link (num_wts,max_links)
  72. !***********************************************************************
  73. contains
  74. !***********************************************************************
  75. subroutine init_remap_vars
  76. !-----------------------------------------------------------------------
  77. !
  78. ! this routine initializes some variables and provides an initial
  79. ! allocation of arrays (fairly large so frequent resizing
  80. ! unnecessary).
  81. !
  82. !-----------------------------------------------------------------------
  83. !-----------------------------------------------------------------------
  84. !
  85. ! determine the number of weights
  86. !
  87. !-----------------------------------------------------------------------
  88. select case (map_type)
  89. case(map_type_conserv)
  90. num_wts = 3
  91. case(map_type_bilinear)
  92. num_wts = 1
  93. case(map_type_bicubic)
  94. num_wts = 4
  95. case(map_type_distwgt)
  96. num_wts = 1
  97. end select
  98. !-----------------------------------------------------------------------
  99. !
  100. ! initialize num_links and set max_links to four times the largest
  101. ! of the destination grid sizes initially (can be changed later).
  102. ! set a default resize increment to increase the size of link
  103. ! arrays if the number of links exceeds the initial size
  104. !
  105. !-----------------------------------------------------------------------
  106. num_links_map1 = 0
  107. max_links_map1 = 4*grid2_size
  108. if (num_maps > 1) then
  109. num_links_map2 = 0
  110. max_links_map1 = max(4*grid1_size,4*grid2_size)
  111. max_links_map2 = max_links_map1
  112. endif
  113. resize_increment = 0.1*max(grid1_size,grid2_size)
  114. !-----------------------------------------------------------------------
  115. !
  116. ! allocate address and weight arrays for mapping 1
  117. !
  118. !-----------------------------------------------------------------------
  119. allocate (grid1_add_map1(max_links_map1),
  120. & grid2_add_map1(max_links_map1),
  121. & wts_map1(num_wts, max_links_map1))
  122. !-----------------------------------------------------------------------
  123. !
  124. ! allocate address and weight arrays for mapping 2 if necessary
  125. !
  126. !-----------------------------------------------------------------------
  127. if (num_maps > 1) then
  128. allocate (grid1_add_map2(max_links_map2),
  129. & grid2_add_map2(max_links_map2),
  130. & wts_map2(num_wts, max_links_map2))
  131. endif
  132. !-----------------------------------------------------------------------
  133. end subroutine init_remap_vars
  134. !***********************************************************************
  135. subroutine resize_remap_vars(nmap, increment)
  136. !-----------------------------------------------------------------------
  137. !
  138. ! this routine resizes remapping arrays by increasing(decreasing)
  139. ! the max_links by increment
  140. !
  141. !-----------------------------------------------------------------------
  142. !-----------------------------------------------------------------------
  143. !
  144. ! input variables
  145. !
  146. !-----------------------------------------------------------------------
  147. integer (kind=int_kind), intent(in) ::
  148. & nmap, ! identifies which mapping array to resize
  149. & increment ! the number of links to add(subtract) to arrays
  150. !-----------------------------------------------------------------------
  151. !
  152. ! local variables
  153. !
  154. !-----------------------------------------------------------------------
  155. integer (kind=int_kind) ::
  156. & ierr, ! error flag
  157. & mxlinks ! size of link arrays
  158. integer (kind=int_kind), dimension(:), allocatable ::
  159. & add1_tmp, ! temp array for resizing address arrays
  160. & add2_tmp ! temp array for resizing address arrays
  161. real (kind=dbl_kind), dimension(:,:), allocatable ::
  162. & wts_tmp ! temp array for resizing weight arrays
  163. !-----------------------------------------------------------------------
  164. !
  165. ! resize map 1 arrays if required.
  166. !
  167. !-----------------------------------------------------------------------
  168. select case (nmap)
  169. case(1)
  170. !***
  171. !*** allocate temporaries to hold original values
  172. !***
  173. mxlinks = size(grid1_add_map1)
  174. allocate (add1_tmp(mxlinks), add2_tmp(mxlinks),
  175. & wts_tmp(num_wts,mxlinks))
  176. add1_tmp = grid1_add_map1
  177. add2_tmp = grid2_add_map1
  178. wts_tmp = wts_map1
  179. !***
  180. !*** deallocate originals and increment max_links then
  181. !*** reallocate arrays at new size
  182. !***
  183. deallocate (grid1_add_map1, grid2_add_map1, wts_map1)
  184. max_links_map1 = mxlinks + increment
  185. allocate (grid1_add_map1(max_links_map1),
  186. & grid2_add_map1(max_links_map1),
  187. & wts_map1(num_wts,max_links_map1))
  188. !***
  189. !*** restore original values from temp arrays and
  190. !*** deallocate temps
  191. !***
  192. mxlinks = min(mxlinks, max_links_map1)
  193. grid1_add_map1(1:mxlinks) = add1_tmp (1:mxlinks)
  194. grid2_add_map1(1:mxlinks) = add2_tmp (1:mxlinks)
  195. wts_map1 (:,1:mxlinks) = wts_tmp(:,1:mxlinks)
  196. deallocate(add1_tmp, add2_tmp, wts_tmp)
  197. !-----------------------------------------------------------------------
  198. !
  199. ! resize map 2 arrays if required.
  200. !
  201. !-----------------------------------------------------------------------
  202. case(2)
  203. !***
  204. !*** allocate temporaries to hold original values
  205. !***
  206. mxlinks = size(grid1_add_map2)
  207. allocate (add1_tmp(mxlinks), add2_tmp(mxlinks),
  208. & wts_tmp(num_wts,mxlinks),stat=ierr)
  209. if (ierr .ne. 0) then
  210. print *,'error allocating temps in resize: ',ierr
  211. stop
  212. endif
  213. add1_tmp = grid1_add_map2
  214. add2_tmp = grid2_add_map2
  215. wts_tmp = wts_map2
  216. !***
  217. !*** deallocate originals and increment max_links then
  218. !*** reallocate arrays at new size
  219. !***
  220. deallocate (grid1_add_map2, grid2_add_map2, wts_map2)
  221. max_links_map2 = mxlinks + increment
  222. allocate (grid1_add_map2(max_links_map2),
  223. & grid2_add_map2(max_links_map2),
  224. & wts_map2(num_wts,max_links_map2),stat=ierr)
  225. if (ierr .ne. 0) then
  226. print *,'error allocating new arrays in resize: ',ierr
  227. stop
  228. endif
  229. !***
  230. !*** restore original values from temp arrays and
  231. !*** deallocate temps
  232. !***
  233. mxlinks = min(mxlinks, max_links_map2)
  234. grid1_add_map2(1:mxlinks) = add1_tmp (1:mxlinks)
  235. grid2_add_map2(1:mxlinks) = add2_tmp (1:mxlinks)
  236. wts_map2 (:,1:mxlinks) = wts_tmp(:,1:mxlinks)
  237. deallocate(add1_tmp, add2_tmp, wts_tmp)
  238. end select
  239. !-----------------------------------------------------------------------
  240. end subroutine resize_remap_vars
  241. !***********************************************************************
  242. end module remap_vars
  243. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!