scrip.F90 8.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! This routine is the driver for computing the addresses and weights
  4. ! for interpolating between two grids on a sphere.
  5. !
  6. ! Modified slightly to get name of namelist file from command line - sga 2/12/05
  7. !
  8. !-----------------------------------------------------------------------
  9. !
  10. ! CVS:$Id: scrip.f,v 1.6 2001/08/21 21:06:44 pwjones Exp $
  11. !
  12. ! Copyright (c) 1997, 1998 the Regents of the University of
  13. ! California.
  14. !
  15. ! This software and ancillary information (herein called software)
  16. ! called SCRIP is made available under the terms described here.
  17. ! The software has been approved for release with associated
  18. ! LA-CC Number 98-45.
  19. !
  20. ! Unless otherwise indicated, this software has been authored
  21. ! by an employee or employees of the University of California,
  22. ! operator of the Los Alamos National Laboratory under Contract
  23. ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
  24. ! Government has rights to use, reproduce, and distribute this
  25. ! software. The public may copy and use this software without
  26. ! charge, provided that this Notice and any statement of authorship
  27. ! are reproduced on all copies. Neither the Government nor the
  28. ! University makes any warranty, express or implied, or assumes
  29. ! any liability or responsibility for the use of this software.
  30. !
  31. ! If software is modified to produce derivative works, such modified
  32. ! software should be clearly marked, so as not to confuse it with
  33. ! the version available from Los Alamos National Laboratory.
  34. !
  35. !***********************************************************************
  36. program scrip
  37. !-----------------------------------------------------------------------
  38. use kinds_mod ! module defining data types
  39. use constants ! module for common constants
  40. use iounits ! I/O unit manager
  41. use timers ! CPU timers
  42. use grids ! module with grid information
  43. use remap_vars ! common remapping variables
  44. use remap_conservative ! routines for conservative remap
  45. use remap_distance_weight ! routines for dist-weight remap
  46. use remap_bilinear ! routines for bilinear interp
  47. use remap_bicubic ! routines for bicubic interp
  48. use remap_write ! routines for remap output
  49. implicit none
  50. !-----------------------------------------------------------------------
  51. !
  52. ! input namelist variables
  53. !
  54. !-----------------------------------------------------------------------
  55. character (char_len) :: &
  56. grid1_file, & ! filename of grid file containing grid1
  57. grid2_file, & ! filename of grid file containing grid2
  58. interp_file1, & ! filename for output remap data (map1)
  59. interp_file2, & ! filename for output remap data (map2)
  60. map1_name, & ! name for mapping from grid1 to grid2
  61. map2_name, & ! name for mapping from grid2 to grid1
  62. map_method, & ! choice for mapping method
  63. normalize_opt, & ! option for normalizing weights
  64. output_opt ! option for output conventions
  65. integer (kind=int_kind) :: &
  66. nmap ! number of mappings to compute (1 or 2)
  67. namelist /remap_inputs/ grid1_file, grid2_file, &
  68. interp_file1, interp_file2, &
  69. map1_name, map2_name, num_maps, &
  70. luse_grid1_area, luse_grid2_area, &
  71. map_method, normalize_opt, output_opt, &
  72. restrict_type, num_srch_bins
  73. !-----------------------------------------------------------------------
  74. !
  75. ! local variables
  76. !
  77. !-----------------------------------------------------------------------
  78. integer (kind=int_kind) :: n, & ! dummy counter
  79. iunit ! unit number for namelist file
  80. character (char_len) :: nm_in
  81. #if defined ARGC
  82. integer :: iargc
  83. external iargc
  84. if (iargc() == 1) then
  85. call getarg(1, nm_in)
  86. else
  87. write(6,*) 'need name of namelist file'
  88. stop
  89. endif
  90. #else
  91. write(6,*) 'enter name for namelist file'
  92. read(5,*) nm_in
  93. #endif
  94. !-----------------------------------------------------------------------
  95. !
  96. ! initialize timers
  97. !
  98. !-----------------------------------------------------------------------
  99. call timers_init
  100. do n=1,max_timers
  101. call timer_clear(n)
  102. end do
  103. !-----------------------------------------------------------------------
  104. !
  105. ! read input namelist
  106. !
  107. !-----------------------------------------------------------------------
  108. grid1_file = 'unknown'
  109. grid2_file = 'unknown'
  110. interp_file1 = 'unknown'
  111. interp_file2 = 'unknown'
  112. map1_name = 'unknown'
  113. map2_name = 'unknown'
  114. luse_grid1_area = .false.
  115. luse_grid2_area = .false.
  116. num_maps = 2
  117. map_type = 1
  118. normalize_opt = 'fracarea'
  119. output_opt = 'scrip'
  120. restrict_type = 'latitude'
  121. num_srch_bins = 900
  122. call get_unit(iunit)
  123. open(iunit, file=nm_in, status='old', form='formatted')
  124. read(iunit, nml=remap_inputs)
  125. call release_unit(iunit)
  126. select case(map_method)
  127. case ('conservative')
  128. map_type = map_type_conserv
  129. luse_grid_centers = .false.
  130. case ('bilinear')
  131. map_type = map_type_bilinear
  132. luse_grid_centers = .true.
  133. case ('bicubic')
  134. map_type = map_type_bicubic
  135. luse_grid_centers = .true.
  136. case ('distwgt')
  137. map_type = map_type_distwgt
  138. luse_grid_centers = .true.
  139. case default
  140. stop 'unknown mapping method'
  141. end select
  142. select case(normalize_opt(1:4))
  143. case ('none')
  144. norm_opt = norm_opt_none
  145. case ('frac')
  146. norm_opt = norm_opt_frcarea
  147. case ('dest')
  148. norm_opt = norm_opt_dstarea
  149. case default
  150. stop 'unknown normalization option'
  151. end select
  152. !-----------------------------------------------------------------------
  153. !
  154. ! initialize grid information for both grids
  155. !
  156. !-----------------------------------------------------------------------
  157. call grid_init(grid1_file, grid2_file)
  158. write(stdout, *) ' Computing remappings between: ',grid1_name
  159. write(stdout, *) ' and ',grid2_name
  160. !-----------------------------------------------------------------------
  161. !
  162. ! initialize some remapping variables.
  163. !
  164. !-----------------------------------------------------------------------
  165. call init_remap_vars
  166. !-----------------------------------------------------------------------
  167. !
  168. ! call appropriate interpolation setup routine based on type of
  169. ! remapping requested.
  170. !
  171. !-----------------------------------------------------------------------
  172. select case(map_type)
  173. case(map_type_conserv)
  174. call remap_conserv
  175. case(map_type_bilinear)
  176. call remap_bilin
  177. case(map_type_distwgt)
  178. call remap_distwgt
  179. case(map_type_bicubic)
  180. call remap_bicub
  181. case default
  182. stop 'Invalid Map Type'
  183. end select
  184. !-----------------------------------------------------------------------
  185. !
  186. ! reduce size of remapping arrays and then write remapping info
  187. ! to a file.
  188. !
  189. !-----------------------------------------------------------------------
  190. if (num_links_map1 /= max_links_map1) then
  191. call resize_remap_vars(1, num_links_map1-max_links_map1)
  192. endif
  193. if ((num_maps > 1) .and. (num_links_map2 /= max_links_map2)) then
  194. call resize_remap_vars(2, num_links_map2-max_links_map2)
  195. endif
  196. call write_remap(map1_name, map2_name, &
  197. interp_file1, interp_file2, output_opt)
  198. !-----------------------------------------------------------------------
  199. end program scrip
  200. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!