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