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$ |
---|
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 | |
---|
36 | module remap_vars |
---|
37 | |
---|
38 | use kinds_mod |
---|
39 | use constants |
---|
40 | use grids |
---|
41 | |
---|
42 | implicit none |
---|
43 | |
---|
44 | !----------------------------------------------------------------------- |
---|
45 | ! |
---|
46 | ! module variables |
---|
47 | ! |
---|
48 | !----------------------------------------------------------------------- |
---|
49 | |
---|
50 | integer (kind=int_kind), parameter :: & |
---|
51 | norm_opt_none = 1 & |
---|
52 | , norm_opt_dstarea = 2 & |
---|
53 | , norm_opt_frcarea = 3 |
---|
54 | |
---|
55 | integer (kind=int_kind), parameter :: & |
---|
56 | map_type_conserv = 1 & |
---|
57 | , map_type_bilinear = 2 & |
---|
58 | , map_type_bicubic = 3 & |
---|
59 | , map_type_distwgt = 4 |
---|
60 | |
---|
61 | integer (kind=int_kind), save :: & |
---|
62 | max_links_map1 & ! current size of link arrays |
---|
63 | , num_links_map1 & ! actual number of links for remapping |
---|
64 | , max_links_map2 & ! current size of link arrays |
---|
65 | , num_links_map2 & ! actual number of links for remapping |
---|
66 | , num_maps & ! num of remappings for this grid pair |
---|
67 | , num_wts & ! num of weights used in remapping |
---|
68 | , map_type & ! identifier for remapping method |
---|
69 | , norm_opt & ! option for normalization (conserv only) |
---|
70 | , resize_increment ! default amount to increase array size |
---|
71 | |
---|
72 | integer (kind=int_kind), dimension(:), allocatable, save :: & |
---|
73 | grid1_add_map1, & ! grid1 address for each link in mapping 1 |
---|
74 | grid2_add_map1, & ! grid2 address for each link in mapping 1 |
---|
75 | grid1_add_map2, & ! grid1 address for each link in mapping 2 |
---|
76 | grid2_add_map2 ! grid2 address for each link in mapping 2 |
---|
77 | |
---|
78 | real (kind=dbl_kind), dimension(:,:), allocatable, save :: & |
---|
79 | wts_map1, & ! map weights for each link (num_wts,max_links) |
---|
80 | wts_map2 ! map weights for each link (num_wts,max_links) |
---|
81 | |
---|
82 | !*********************************************************************** |
---|
83 | |
---|
84 | contains |
---|
85 | |
---|
86 | !*********************************************************************** |
---|
87 | |
---|
88 | subroutine init_remap_vars |
---|
89 | |
---|
90 | !----------------------------------------------------------------------- |
---|
91 | ! |
---|
92 | ! this routine initializes some variables and provides an initial |
---|
93 | ! allocation of arrays (fairly large so frequent resizing |
---|
94 | ! unnecessary). |
---|
95 | ! |
---|
96 | !----------------------------------------------------------------------- |
---|
97 | |
---|
98 | !----------------------------------------------------------------------- |
---|
99 | ! |
---|
100 | ! determine the number of weights |
---|
101 | ! |
---|
102 | !----------------------------------------------------------------------- |
---|
103 | |
---|
104 | select case (map_type) |
---|
105 | case(map_type_conserv) |
---|
106 | num_wts = 3 |
---|
107 | case(map_type_bilinear) |
---|
108 | num_wts = 1 |
---|
109 | case(map_type_bicubic) |
---|
110 | num_wts = 4 |
---|
111 | case(map_type_distwgt) |
---|
112 | num_wts = 1 |
---|
113 | end select |
---|
114 | |
---|
115 | !----------------------------------------------------------------------- |
---|
116 | ! |
---|
117 | ! initialize num_links and set max_links to four times the largest |
---|
118 | ! of the destination grid sizes initially (can be changed later). |
---|
119 | ! set a default resize increment to increase the size of link |
---|
120 | ! arrays if the number of links exceeds the initial size |
---|
121 | ! |
---|
122 | !----------------------------------------------------------------------- |
---|
123 | |
---|
124 | num_links_map1 = 0 |
---|
125 | max_links_map1 = 4*grid2_size |
---|
126 | if (num_maps > 1) then |
---|
127 | num_links_map2 = 0 |
---|
128 | max_links_map1 = max(4*grid1_size,4*grid2_size) |
---|
129 | max_links_map2 = max_links_map1 |
---|
130 | endif |
---|
131 | |
---|
132 | resize_increment = 0.1*max(grid1_size,grid2_size) |
---|
133 | |
---|
134 | !----------------------------------------------------------------------- |
---|
135 | ! |
---|
136 | ! allocate address and weight arrays for mapping 1 |
---|
137 | ! |
---|
138 | !----------------------------------------------------------------------- |
---|
139 | |
---|
140 | allocate (grid1_add_map1(max_links_map1), & |
---|
141 | grid2_add_map1(max_links_map1), & |
---|
142 | wts_map1(num_wts, max_links_map1)) |
---|
143 | |
---|
144 | !----------------------------------------------------------------------- |
---|
145 | ! |
---|
146 | ! allocate address and weight arrays for mapping 2 if necessary |
---|
147 | ! |
---|
148 | !----------------------------------------------------------------------- |
---|
149 | |
---|
150 | if (num_maps > 1) then |
---|
151 | allocate (grid1_add_map2(max_links_map2), & |
---|
152 | grid2_add_map2(max_links_map2), & |
---|
153 | wts_map2(num_wts, max_links_map2)) |
---|
154 | endif |
---|
155 | |
---|
156 | !----------------------------------------------------------------------- |
---|
157 | |
---|
158 | end subroutine init_remap_vars |
---|
159 | |
---|
160 | !*********************************************************************** |
---|
161 | |
---|
162 | subroutine resize_remap_vars(nmap, increment) |
---|
163 | |
---|
164 | !----------------------------------------------------------------------- |
---|
165 | ! |
---|
166 | ! this routine resizes remapping arrays by increasing(decreasing) |
---|
167 | ! the max_links by increment |
---|
168 | ! |
---|
169 | !----------------------------------------------------------------------- |
---|
170 | |
---|
171 | !----------------------------------------------------------------------- |
---|
172 | ! |
---|
173 | ! input variables |
---|
174 | ! |
---|
175 | !----------------------------------------------------------------------- |
---|
176 | |
---|
177 | integer (kind=int_kind), intent(in) :: & |
---|
178 | nmap, & ! identifies which mapping array to resize |
---|
179 | increment ! the number of links to add(subtract) to arrays |
---|
180 | |
---|
181 | !----------------------------------------------------------------------- |
---|
182 | ! |
---|
183 | ! local variables |
---|
184 | ! |
---|
185 | !----------------------------------------------------------------------- |
---|
186 | |
---|
187 | integer (kind=int_kind) :: & |
---|
188 | ierr, & ! error flag |
---|
189 | mxlinks ! size of link arrays |
---|
190 | |
---|
191 | integer (kind=int_kind), dimension(:), allocatable :: & |
---|
192 | add1_tmp, & ! temp array for resizing address arrays |
---|
193 | add2_tmp ! temp array for resizing address arrays |
---|
194 | |
---|
195 | real (kind=dbl_kind), dimension(:,:), allocatable :: & |
---|
196 | wts_tmp ! temp array for resizing weight arrays |
---|
197 | |
---|
198 | !----------------------------------------------------------------------- |
---|
199 | ! |
---|
200 | ! resize map 1 arrays if required. |
---|
201 | ! |
---|
202 | !----------------------------------------------------------------------- |
---|
203 | |
---|
204 | select case (nmap) |
---|
205 | case(1) |
---|
206 | |
---|
207 | !*** |
---|
208 | !*** allocate temporaries to hold original values |
---|
209 | !*** |
---|
210 | |
---|
211 | mxlinks = size(grid1_add_map1) |
---|
212 | allocate (add1_tmp(mxlinks), add2_tmp(mxlinks), & |
---|
213 | wts_tmp(num_wts,mxlinks)) |
---|
214 | |
---|
215 | add1_tmp = grid1_add_map1 |
---|
216 | add2_tmp = grid2_add_map1 |
---|
217 | wts_tmp = wts_map1 |
---|
218 | |
---|
219 | !*** |
---|
220 | !*** deallocate originals and increment max_links then |
---|
221 | !*** reallocate arrays at new size |
---|
222 | !*** |
---|
223 | |
---|
224 | deallocate (grid1_add_map1, grid2_add_map1, wts_map1) |
---|
225 | max_links_map1 = mxlinks + increment |
---|
226 | allocate (grid1_add_map1(max_links_map1), & |
---|
227 | grid2_add_map1(max_links_map1), & |
---|
228 | wts_map1(num_wts,max_links_map1)) |
---|
229 | |
---|
230 | !*** |
---|
231 | !*** restore original values from temp arrays and |
---|
232 | !*** deallocate temps |
---|
233 | !*** |
---|
234 | |
---|
235 | mxlinks = min(mxlinks, max_links_map1) |
---|
236 | grid1_add_map1(1:mxlinks) = add1_tmp (1:mxlinks) |
---|
237 | grid2_add_map1(1:mxlinks) = add2_tmp (1:mxlinks) |
---|
238 | wts_map1 (:,1:mxlinks) = wts_tmp(:,1:mxlinks) |
---|
239 | deallocate(add1_tmp, add2_tmp, wts_tmp) |
---|
240 | |
---|
241 | !----------------------------------------------------------------------- |
---|
242 | ! |
---|
243 | ! resize map 2 arrays if required. |
---|
244 | ! |
---|
245 | !----------------------------------------------------------------------- |
---|
246 | |
---|
247 | case(2) |
---|
248 | |
---|
249 | !*** |
---|
250 | !*** allocate temporaries to hold original values |
---|
251 | !*** |
---|
252 | |
---|
253 | mxlinks = size(grid1_add_map2) |
---|
254 | allocate (add1_tmp(mxlinks), add2_tmp(mxlinks), & |
---|
255 | wts_tmp(num_wts,mxlinks),stat=ierr) |
---|
256 | if (ierr .ne. 0) then |
---|
257 | print *,'error allocating temps in resize: ',ierr |
---|
258 | stop |
---|
259 | endif |
---|
260 | |
---|
261 | add1_tmp = grid1_add_map2 |
---|
262 | add2_tmp = grid2_add_map2 |
---|
263 | wts_tmp = wts_map2 |
---|
264 | |
---|
265 | !*** |
---|
266 | !*** deallocate originals and increment max_links then |
---|
267 | !*** reallocate arrays at new size |
---|
268 | !*** |
---|
269 | |
---|
270 | deallocate (grid1_add_map2, grid2_add_map2, wts_map2) |
---|
271 | max_links_map2 = mxlinks + increment |
---|
272 | allocate (grid1_add_map2(max_links_map2), & |
---|
273 | grid2_add_map2(max_links_map2), & |
---|
274 | wts_map2(num_wts,max_links_map2),stat=ierr) |
---|
275 | if (ierr .ne. 0) then |
---|
276 | print *,'error allocating new arrays in resize: ',ierr |
---|
277 | stop |
---|
278 | endif |
---|
279 | |
---|
280 | |
---|
281 | !*** |
---|
282 | !*** restore original values from temp arrays and |
---|
283 | !*** deallocate temps |
---|
284 | !*** |
---|
285 | |
---|
286 | mxlinks = min(mxlinks, max_links_map2) |
---|
287 | grid1_add_map2(1:mxlinks) = add1_tmp (1:mxlinks) |
---|
288 | grid2_add_map2(1:mxlinks) = add2_tmp (1:mxlinks) |
---|
289 | wts_map2 (:,1:mxlinks) = wts_tmp(:,1:mxlinks) |
---|
290 | deallocate(add1_tmp, add2_tmp, wts_tmp) |
---|
291 | |
---|
292 | end select |
---|
293 | |
---|
294 | !----------------------------------------------------------------------- |
---|
295 | |
---|
296 | end subroutine resize_remap_vars |
---|
297 | |
---|
298 | !*********************************************************************** |
---|
299 | |
---|
300 | end module remap_vars |
---|
301 | |
---|
302 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|