[1677] | 1 | C**** |
---|
| 2 | C ***************************** |
---|
| 3 | C * OASIS MODULE - LEVEL ? * |
---|
| 4 | C * ------------- ------- * |
---|
| 5 | C ***************************** |
---|
| 6 | C |
---|
| 7 | C**** remap_bilinear - calculate bilinear remapping |
---|
| 8 | C |
---|
| 9 | C Purpose: |
---|
| 10 | C ------- |
---|
| 11 | C Adaptation of SCRIP 1.4 remap_bilinear module to calculate |
---|
| 12 | C bilinear remapping. |
---|
| 13 | C |
---|
| 14 | C** Interface: |
---|
| 15 | C --------- |
---|
| 16 | C *CALL* *remap_bilin** |
---|
| 17 | C |
---|
| 18 | C Input: |
---|
| 19 | C ----- |
---|
| 20 | C |
---|
| 21 | C Output: |
---|
| 22 | C ------ |
---|
| 23 | C |
---|
| 24 | C History: |
---|
| 25 | C ------- |
---|
| 26 | C Version Programmer Date Description |
---|
| 27 | C ------- ---------- ---- ----------- |
---|
| 28 | C 2.5 S. Valcke 2002/09 Treatment for masked point |
---|
| 29 | C |
---|
| 30 | C %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 31 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 32 | ! |
---|
| 33 | ! this module contains necessary routines for performing an |
---|
| 34 | ! bilinear interpolation. |
---|
| 35 | ! |
---|
| 36 | !----------------------------------------------------------------------- |
---|
| 37 | ! |
---|
| 38 | ! CVS:$Id: remap_bilinear.f,v 1.1.1.1 2005/03/23 16:01:11 adm Exp $ |
---|
| 39 | ! |
---|
| 40 | ! Copyright (c) 1997, 1998 the Regents of the University of |
---|
| 41 | ! California. |
---|
| 42 | ! |
---|
| 43 | ! This software and ancillary information (herein called software) |
---|
| 44 | ! called SCRIP is made available under the terms described here. |
---|
| 45 | ! The software has been approved for release with associated |
---|
| 46 | ! LA-CC Number 98-45. |
---|
| 47 | ! |
---|
| 48 | ! Unless otherwise indicated, this software has been authored |
---|
| 49 | ! by an employee or employees of the University of California, |
---|
| 50 | ! operator of the Los Alamos National Laboratory under Contract |
---|
| 51 | ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S. |
---|
| 52 | ! Government has rights to use, reproduce, and distribute this |
---|
| 53 | ! software. The public may copy and use this software without |
---|
| 54 | ! charge, provided that this Notice and any statement of authorship |
---|
| 55 | ! are reproduced on all copies. Neither the Government nor the |
---|
| 56 | ! University makes any warranty, express or implied, or assumes |
---|
| 57 | ! any liability or responsibility for the use of this software. |
---|
| 58 | ! |
---|
| 59 | ! If software is modified to produce derivative works, such modified |
---|
| 60 | ! software should be clearly marked, so as not to confuse it with |
---|
| 61 | ! the version available from Los Alamos National Laboratory. |
---|
| 62 | ! |
---|
| 63 | !*********************************************************************** |
---|
| 64 | |
---|
| 65 | module remap_bilinear |
---|
| 66 | |
---|
| 67 | !----------------------------------------------------------------------- |
---|
| 68 | |
---|
| 69 | use kinds_mod ! defines common data types |
---|
| 70 | use constants ! defines common constants |
---|
| 71 | use grids ! module containing grid info |
---|
| 72 | use remap_vars ! module containing remap info |
---|
| 73 | USE mod_unit |
---|
| 74 | USE mod_printing |
---|
| 75 | implicit NONE |
---|
| 76 | |
---|
| 77 | !----------------------------------------------------------------------- |
---|
| 78 | |
---|
| 79 | integer (kind=int_kind), parameter :: |
---|
| 80 | & max_iter = 100 ! max iteration count for i,j iteration |
---|
| 81 | |
---|
| 82 | real (kind=dbl_kind), parameter :: |
---|
| 83 | & converge = 1.e-10_dbl_kind ! convergence criterion |
---|
| 84 | |
---|
| 85 | !*********************************************************************** |
---|
| 86 | |
---|
| 87 | contains |
---|
| 88 | |
---|
| 89 | !*********************************************************************** |
---|
| 90 | |
---|
| 91 | subroutine remap_bilin (lextrapdone) |
---|
| 92 | |
---|
| 93 | !----------------------------------------------------------------------- |
---|
| 94 | ! |
---|
| 95 | ! this routine computes the weights for a bilinear interpolation. |
---|
| 96 | ! |
---|
| 97 | !----------------------------------------------------------------------- |
---|
| 98 | |
---|
| 99 | LOGICAL :: |
---|
| 100 | & lextrapdone ! logical, true if EXTRAP done on field |
---|
| 101 | |
---|
| 102 | !----------------------------------------------------------------------- |
---|
| 103 | ! |
---|
| 104 | ! local variables |
---|
| 105 | ! |
---|
| 106 | !----------------------------------------------------------------------- |
---|
| 107 | |
---|
| 108 | integer (kind=int_kind) :: n,icount, i, |
---|
| 109 | & dst_add, ! destination address |
---|
| 110 | & iter, ! iteration counter |
---|
| 111 | & nmap ! index of current map being computed |
---|
| 112 | |
---|
| 113 | integer (kind=int_kind), dimension(4) :: |
---|
| 114 | & src_add ! address for the four source points |
---|
| 115 | |
---|
| 116 | real (kind=dbl_kind), dimension(4) :: |
---|
| 117 | & src_lats, ! latitudes of four bilinear corners |
---|
| 118 | & src_lons, ! longitudes of four bilinear corners |
---|
| 119 | & wgts ! bilinear weights for four corners |
---|
| 120 | |
---|
| 121 | real (kind=dbl_kind) :: |
---|
| 122 | & plat, plon, ! lat/lon coords of destination point |
---|
| 123 | & iguess, jguess, ! current guess for bilinear coordinate |
---|
| 124 | & deli, delj, ! corrections to i,j |
---|
| 125 | & dth1, dth2, dth3, ! some latitude differences |
---|
| 126 | & dph1, dph2, dph3, ! some longitude differences |
---|
| 127 | & dthp, dphp, ! difference between point and sw corner |
---|
| 128 | & mat1, mat2, mat3, mat4, ! matrix elements |
---|
| 129 | & determinant, ! matrix determinant |
---|
| 130 | & sum_wgts ! sum of weights for normalization |
---|
| 131 | |
---|
| 132 | real (kind=dbl_kind) :: |
---|
| 133 | & coslat_dst, sinlat_dst, coslon_dst, sinlon_dst, |
---|
| 134 | & dist_min, distance, ! for computing dist-weighted avg |
---|
| 135 | & src_latsnn |
---|
| 136 | |
---|
| 137 | integer (kind=int_kind) :: min_add, max_add, srch_add, src_addnn |
---|
| 138 | |
---|
| 139 | !----------------------------------------------------------------------- |
---|
| 140 | ! |
---|
| 141 | ! compute mappings from grid1 to grid2 |
---|
| 142 | ! |
---|
| 143 | !----------------------------------------------------------------------- |
---|
| 144 | |
---|
| 145 | nmap = 1 |
---|
| 146 | if (grid1_rank /= 2) then |
---|
| 147 | stop 'Can not do bilinear interpolation when grid_rank /= 2' |
---|
| 148 | endif |
---|
| 149 | |
---|
| 150 | !*** |
---|
| 151 | !*** loop over destination grid |
---|
| 152 | !*** |
---|
| 153 | |
---|
| 154 | grid_loop1: do dst_add = 1, grid2_size |
---|
| 155 | |
---|
| 156 | if (.not. grid2_mask(dst_add)) cycle grid_loop1 |
---|
| 157 | |
---|
| 158 | plat = grid2_center_lat(dst_add) |
---|
| 159 | plon = grid2_center_lon(dst_add) |
---|
| 160 | |
---|
| 161 | !*** |
---|
| 162 | !*** find nearest square of grid points on source grid |
---|
| 163 | !*** |
---|
| 164 | |
---|
| 165 | call grid_search_bilin(src_add, src_lats, src_lons, |
---|
| 166 | & min_add, max_add, |
---|
| 167 | & plat, plon, grid1_dims, |
---|
| 168 | & grid1_center_lat, grid1_center_lon, |
---|
| 169 | & grid1_bound_box, bin_addr1, lextrapdone) |
---|
| 170 | |
---|
| 171 | |
---|
| 172 | if (src_add(1) > 0) THEN |
---|
| 173 | |
---|
| 174 | !*** |
---|
| 175 | !*** if the 4 surrounding points have been found and are |
---|
| 176 | !*** non-masked, do bilinear interpolation |
---|
| 177 | !*** |
---|
| 178 | |
---|
| 179 | grid2_frac(dst_add) = one |
---|
| 180 | |
---|
| 181 | |
---|
| 182 | dth1 = src_lats(2) - src_lats(1) |
---|
| 183 | dth2 = src_lats(4) - src_lats(1) |
---|
| 184 | dth3 = src_lats(3) - src_lats(2) - dth2 |
---|
| 185 | |
---|
| 186 | dph1 = src_lons(2) - src_lons(1) |
---|
| 187 | dph2 = src_lons(4) - src_lons(1) |
---|
| 188 | dph3 = src_lons(3) - src_lons(2) |
---|
| 189 | |
---|
| 190 | if (dph1 > three*pih) dph1 = dph1 - pi2 |
---|
| 191 | if (dph2 > three*pih) dph2 = dph2 - pi2 |
---|
| 192 | if (dph3 > three*pih) dph3 = dph3 - pi2 |
---|
| 193 | if (dph1 < -three*pih) dph1 = dph1 + pi2 |
---|
| 194 | if (dph2 < -three*pih) dph2 = dph2 + pi2 |
---|
| 195 | if (dph3 < -three*pih) dph3 = dph3 + pi2 |
---|
| 196 | |
---|
| 197 | dph3 = dph3 - dph2 |
---|
| 198 | |
---|
| 199 | iguess = half |
---|
| 200 | jguess = half |
---|
| 201 | |
---|
| 202 | iter_loop1: do iter=1,max_iter |
---|
| 203 | |
---|
| 204 | dthp = plat - src_lats(1) - dth1*iguess - |
---|
| 205 | & dth2*jguess - dth3*iguess*jguess |
---|
| 206 | dphp = plon - src_lons(1) |
---|
| 207 | |
---|
| 208 | if (dphp > three*pih) dphp = dphp - pi2 |
---|
| 209 | if (dphp < -three*pih) dphp = dphp + pi2 |
---|
| 210 | |
---|
| 211 | dphp = dphp - dph1*iguess - dph2*jguess - |
---|
| 212 | & dph3*iguess*jguess |
---|
| 213 | |
---|
| 214 | mat1 = dth1 + dth3*jguess |
---|
| 215 | mat2 = dth2 + dth3*iguess |
---|
| 216 | mat3 = dph1 + dph3*jguess |
---|
| 217 | mat4 = dph2 + dph3*iguess |
---|
| 218 | |
---|
| 219 | determinant = mat1*mat4 - mat2*mat3 |
---|
| 220 | |
---|
| 221 | deli = (dthp*mat4 - mat2*dphp)/determinant |
---|
| 222 | delj = (mat1*dphp - dthp*mat3)/determinant |
---|
| 223 | |
---|
| 224 | if (abs(deli) < converge .and. |
---|
| 225 | & abs(delj) < converge) exit iter_loop1 |
---|
| 226 | |
---|
| 227 | iguess = iguess + deli |
---|
| 228 | jguess = jguess + delj |
---|
| 229 | |
---|
| 230 | end do iter_loop1 |
---|
| 231 | |
---|
| 232 | if (iter <= max_iter) then |
---|
| 233 | |
---|
| 234 | !*** |
---|
| 235 | !*** successfully found i,j - compute weights |
---|
| 236 | !*** |
---|
| 237 | |
---|
| 238 | wgts(1) = (one-iguess)*(one-jguess) |
---|
| 239 | wgts(2) = iguess*(one-jguess) |
---|
| 240 | wgts(3) = iguess*jguess |
---|
| 241 | wgts(4) = (one-iguess)*jguess |
---|
| 242 | |
---|
| 243 | call store_link_bilin(dst_add, src_add, wgts, nmap) |
---|
| 244 | |
---|
| 245 | else |
---|
| 246 | print *,'Point coords: ',plat,plon |
---|
| 247 | print *,'Dest grid lats: ',src_lats |
---|
| 248 | print *,'Dest grid lons: ',src_lons |
---|
| 249 | print *,'Dest grid addresses: ',src_add |
---|
| 250 | print *,'Current i,j : ',iguess, jguess |
---|
| 251 | stop 'Iteration for i,j exceed max iteration count' |
---|
| 252 | endif |
---|
| 253 | |
---|
| 254 | else if (src_add(1) < 0) THEN |
---|
| 255 | |
---|
| 256 | !*** |
---|
| 257 | !*** Search for bilinear failed or at least one of the 4 |
---|
| 258 | !*** neighbours was masked. Do distance-weighted average using |
---|
| 259 | !*** the non-masked points among the 4 closest ones. |
---|
| 260 | !*** |
---|
| 261 | |
---|
| 262 | IF (nlogprt .eq. 2) then |
---|
| 263 | WRITE(nulou,*) ' ' |
---|
| 264 | WRITE(nulou,*) |
---|
| 265 | & 'WARNING: Cannot make bilinear interpolation for target point' |
---|
| 266 | & ,dst_add |
---|
| 267 | WRITE(nulou,*) |
---|
| 268 | & 'Using non-masked points among 4 nearest neighbors.' |
---|
| 269 | WRITE(nulou,*) ' ' |
---|
| 270 | ENDIF |
---|
| 271 | |
---|
| 272 | !*** |
---|
| 273 | !*** Find the 4 closest points |
---|
| 274 | !*** |
---|
| 275 | |
---|
| 276 | coslat_dst = cos(plat) |
---|
| 277 | sinlat_dst = sin(plat) |
---|
| 278 | coslon_dst = cos(plon) |
---|
| 279 | sinlon_dst = sin(plon) |
---|
| 280 | src_add = 0 |
---|
| 281 | dist_min = bignum |
---|
| 282 | src_lats = bignum |
---|
| 283 | do srch_add = min_add,max_add |
---|
| 284 | distance=acos(coslat_dst*cos(grid1_center_lat(srch_add))* |
---|
| 285 | & (coslon_dst*cos(grid1_center_lon(srch_add)) + |
---|
| 286 | & sinlon_dst*sin(grid1_center_lon(srch_add)))+ |
---|
| 287 | & sinlat_dst*sin(grid1_center_lat(srch_add))) |
---|
| 288 | |
---|
| 289 | if (distance < dist_min) then |
---|
| 290 | sort_loop: do n=1,4 |
---|
| 291 | if (distance < src_lats(n)) then |
---|
| 292 | do i=4,n+1,-1 |
---|
| 293 | src_add (i) = src_add (i-1) |
---|
| 294 | src_lats(i) = src_lats(i-1) |
---|
| 295 | end do |
---|
| 296 | src_add (n) = srch_add |
---|
| 297 | src_lats(n) = distance |
---|
| 298 | dist_min = src_lats(4) |
---|
| 299 | exit sort_loop |
---|
| 300 | endif |
---|
| 301 | end do sort_loop |
---|
| 302 | endif |
---|
| 303 | end do |
---|
| 304 | |
---|
| 305 | src_lons = one/(src_lats + tiny) |
---|
| 306 | distance = sum(src_lons) |
---|
| 307 | src_lats = src_lons/distance |
---|
| 308 | |
---|
| 309 | !*** |
---|
| 310 | !*** Among 4 closest points, keep only the non-masked ones |
---|
| 311 | !*** |
---|
| 312 | |
---|
| 313 | icount = 0 |
---|
| 314 | do n=1,4 |
---|
| 315 | if (grid1_mask(src_add(n)) .or. |
---|
| 316 | & (.not. grid1_mask(src_add(n)) .and. lextrapdone)) then |
---|
| 317 | icount = icount + 1 |
---|
| 318 | else |
---|
| 319 | src_lats(n) = zero |
---|
| 320 | endif |
---|
| 321 | end do |
---|
| 322 | |
---|
| 323 | if (icount > 0) then |
---|
| 324 | !*** renormalize weights |
---|
| 325 | sum_wgts = sum(src_lats) |
---|
| 326 | wgts(1) = src_lats(1)/sum_wgts |
---|
| 327 | wgts(2) = src_lats(2)/sum_wgts |
---|
| 328 | wgts(3) = src_lats(3)/sum_wgts |
---|
| 329 | wgts(4) = src_lats(4)/sum_wgts |
---|
| 330 | |
---|
| 331 | grid2_frac(dst_add) = one |
---|
| 332 | call store_link_bilin(dst_add, src_add, wgts, nmap) |
---|
| 333 | ELSE |
---|
| 334 | WRITE(nulou,*) ' ' |
---|
| 335 | WRITE(nulou,*) |
---|
| 336 | & 'All 4 surrounding source grid points are masked' |
---|
| 337 | WRITE(nulou,*) 'for target point ',dst_add |
---|
| 338 | WRITE(nulou,*) 'with longitude and latitude', plon, plat |
---|
| 339 | WRITE(nulou,*) 'Using the nearest non-masked neighbour.' |
---|
| 340 | WRITE(nulou,*) ' ' |
---|
| 341 | |
---|
| 342 | src_latsnn = bignum |
---|
| 343 | do srch_add = min_add,max_add |
---|
| 344 | if (grid1_mask(srch_add) .or. |
---|
| 345 | & (.not. grid1_mask(srch_add) .and. lextrapdone)) then |
---|
| 346 | distance=acos(coslat_dst |
---|
| 347 | & *cos(grid1_center_lat(srch_add))* |
---|
| 348 | & (coslon_dst*cos(grid1_center_lon(srch_add)) + |
---|
| 349 | & sinlon_dst*sin(grid1_center_lon(srch_add)))+ |
---|
| 350 | & sinlat_dst*sin(grid1_center_lat(srch_add))) |
---|
| 351 | |
---|
| 352 | if (distance < src_latsnn) then |
---|
| 353 | src_addnn = srch_add |
---|
| 354 | src_latsnn = distance |
---|
| 355 | endif |
---|
| 356 | endif |
---|
| 357 | end DO |
---|
| 358 | |
---|
| 359 | WRITE(nulou,*) |
---|
| 360 | & 'Nearest non masked neighbour is source point ',src_addnn |
---|
| 361 | WRITE(nulou,*) 'with longitude and latitude', |
---|
| 362 | & grid1_center_lon(src_addnn), grid1_center_lat(src_addnn) |
---|
| 363 | WRITE(nulou,*) ' ' |
---|
| 364 | |
---|
| 365 | wgts(1) = 1. |
---|
| 366 | wgts(2) = 0. |
---|
| 367 | wgts(3) = 0. |
---|
| 368 | wgts(4) = 0. |
---|
| 369 | src_add(1) = src_addnn |
---|
| 370 | src_add(2) = 0 |
---|
| 371 | src_add(3) = 0 |
---|
| 372 | src_add(4) = 0 |
---|
| 373 | |
---|
| 374 | grid2_frac(dst_add) = one |
---|
| 375 | call store_link_bilin(dst_add, src_add, wgts, nmap) |
---|
| 376 | |
---|
| 377 | ENDIF |
---|
| 378 | ENDIF |
---|
| 379 | end do grid_loop1 |
---|
| 380 | |
---|
| 381 | end subroutine remap_bilin |
---|
| 382 | |
---|
| 383 | !*********************************************************************** |
---|
| 384 | |
---|
| 385 | subroutine grid_search_bilin(src_add, src_lats, src_lons, |
---|
| 386 | & min_add, max_add, |
---|
| 387 | & plat, plon, src_grid_dims, |
---|
| 388 | & src_center_lat, src_center_lon, |
---|
| 389 | & src_grid_bound_box, |
---|
| 390 | & src_bin_add, lextrapdone) |
---|
| 391 | |
---|
| 392 | !----------------------------------------------------------------------- |
---|
| 393 | ! |
---|
| 394 | ! this routine finds the location of the search point plat, plon |
---|
| 395 | ! in the source grid and returns the corners needed for a bilinear |
---|
| 396 | ! interpolation. |
---|
| 397 | ! |
---|
| 398 | !----------------------------------------------------------------------- |
---|
| 399 | |
---|
| 400 | !----------------------------------------------------------------------- |
---|
| 401 | ! |
---|
| 402 | ! output variables |
---|
| 403 | ! |
---|
| 404 | !----------------------------------------------------------------------- |
---|
| 405 | |
---|
| 406 | integer (kind=int_kind), dimension(4), intent(out) :: |
---|
| 407 | & src_add ! address of each corner point enclosing P |
---|
| 408 | |
---|
| 409 | real (kind=dbl_kind), dimension(4), intent(out) :: |
---|
| 410 | & src_lats, ! latitudes of the four corner points |
---|
| 411 | & src_lons ! longitudes of the four corner points |
---|
| 412 | |
---|
| 413 | integer (kind=int_kind) :: min_add, max_add |
---|
| 414 | |
---|
| 415 | !----------------------------------------------------------------------- |
---|
| 416 | ! |
---|
| 417 | ! input variables |
---|
| 418 | ! |
---|
| 419 | !----------------------------------------------------------------------- |
---|
| 420 | |
---|
| 421 | real (kind=dbl_kind), intent(in) :: |
---|
| 422 | & plat, ! latitude of the search point |
---|
| 423 | & plon ! longitude of the search point |
---|
| 424 | |
---|
| 425 | integer (kind=int_kind), dimension(2), intent(in) :: |
---|
| 426 | & src_grid_dims ! size of each src grid dimension |
---|
| 427 | |
---|
| 428 | real (kind=dbl_kind), dimension(:), intent(in) :: |
---|
| 429 | & src_center_lat, ! latitude of each src grid center |
---|
| 430 | & src_center_lon ! longitude of each src grid center |
---|
| 431 | |
---|
| 432 | real (kind=dbl_kind), dimension(:,:), intent(in) :: |
---|
| 433 | & src_grid_bound_box ! bound box for source grid |
---|
| 434 | |
---|
| 435 | integer (kind=int_kind), dimension(:,:), intent(in) :: |
---|
| 436 | & src_bin_add ! latitude bins for restricting |
---|
| 437 | |
---|
| 438 | LOGICAL :: |
---|
| 439 | & lextrapdone ! logical, true if EXTRAP done on field |
---|
| 440 | !----------------------------------------------------------------------- |
---|
| 441 | ! |
---|
| 442 | ! local variables |
---|
| 443 | ! |
---|
| 444 | !----------------------------------------------------------------------- |
---|
| 445 | |
---|
| 446 | integer (kind=int_kind) :: n, next_n, srch_add, ni, ! dummy indices |
---|
| 447 | & nx, ny, ntotmask, ! dimensions of src grid |
---|
| 448 | & i, j, jp1, ip1, n_add, e_add, ne_add ! addresses |
---|
| 449 | |
---|
| 450 | real (kind=dbl_kind) :: ! vectors for cross-product check |
---|
| 451 | & vec1_lat, vec1_lon, |
---|
| 452 | & vec2_lat, vec2_lon, cross_product, cross_product_last |
---|
| 453 | |
---|
| 454 | !----------------------------------------------------------------------- |
---|
| 455 | ! |
---|
| 456 | ! restrict search first using bins |
---|
| 457 | ! |
---|
| 458 | !----------------------------------------------------------------------- |
---|
| 459 | |
---|
| 460 | src_add = 0 |
---|
| 461 | |
---|
| 462 | min_add = size(src_center_lat) |
---|
| 463 | max_add = 1 |
---|
| 464 | do n=1,num_srch_bins |
---|
| 465 | if (plat >= bin_lats(1,n) .and. plat <= bin_lats(2,n) .and. |
---|
| 466 | & plon >= bin_lons(1,n) .and. plon <= bin_lons(2,n)) then |
---|
| 467 | min_add = min(min_add, src_bin_add(1,n)) |
---|
| 468 | max_add = max(max_add, src_bin_add(2,n)) |
---|
| 469 | endif |
---|
| 470 | end do |
---|
| 471 | |
---|
| 472 | !----------------------------------------------------------------------- |
---|
| 473 | ! |
---|
| 474 | ! now perform a more detailed search |
---|
| 475 | ! |
---|
| 476 | !----------------------------------------------------------------------- |
---|
| 477 | |
---|
| 478 | nx = src_grid_dims(1) |
---|
| 479 | ny = src_grid_dims(2) |
---|
| 480 | |
---|
| 481 | srch_loop: do srch_add = min_add,max_add |
---|
| 482 | |
---|
| 483 | !*** first check bounding box |
---|
| 484 | |
---|
| 485 | if (plat <= src_grid_bound_box(2,srch_add) .and. |
---|
| 486 | & plat >= src_grid_bound_box(1,srch_add) .and. |
---|
| 487 | & plon <= src_grid_bound_box(4,srch_add) .and. |
---|
| 488 | & plon >= src_grid_bound_box(3,srch_add)) then |
---|
| 489 | |
---|
| 490 | !*** |
---|
| 491 | !*** we are within bounding box so get really serious |
---|
| 492 | !*** |
---|
| 493 | |
---|
| 494 | !*** determine neighbor addresses |
---|
| 495 | |
---|
| 496 | j = (srch_add - 1)/nx +1 |
---|
| 497 | i = srch_add - (j-1)*nx |
---|
| 498 | |
---|
| 499 | !*** find ip1 |
---|
| 500 | !*** Note: I do not want to take into account the number |
---|
| 501 | !*** of overlapping grid points, as the underlying cell |
---|
| 502 | !*** will be found in all cases if the grid overlaps. |
---|
| 503 | |
---|
| 504 | if (i < nx) then |
---|
| 505 | ip1 = i + 1 |
---|
| 506 | else |
---|
| 507 | ip1 = 1 |
---|
| 508 | endif |
---|
| 509 | |
---|
| 510 | if (j < ny) then |
---|
| 511 | jp1 = j+1 |
---|
| 512 | else |
---|
| 513 | jp1 = 1 |
---|
| 514 | endif |
---|
| 515 | |
---|
| 516 | n_add = (jp1 - 1)*nx + i |
---|
| 517 | e_add = (j - 1)*nx + ip1 |
---|
| 518 | ne_add = (jp1 - 1)*nx + ip1 |
---|
| 519 | |
---|
| 520 | src_lats(1) = src_center_lat(srch_add) |
---|
| 521 | src_lats(2) = src_center_lat(e_add) |
---|
| 522 | src_lats(3) = src_center_lat(ne_add) |
---|
| 523 | src_lats(4) = src_center_lat(n_add) |
---|
| 524 | |
---|
| 525 | src_lons(1) = src_center_lon(srch_add) |
---|
| 526 | src_lons(2) = src_center_lon(e_add) |
---|
| 527 | src_lons(3) = src_center_lon(ne_add) |
---|
| 528 | src_lons(4) = src_center_lon(n_add) |
---|
| 529 | |
---|
| 530 | !*** |
---|
| 531 | !*** for consistency, we must make sure all lons are in |
---|
| 532 | !*** same 2pi interval |
---|
| 533 | !*** |
---|
| 534 | |
---|
| 535 | vec1_lon = src_lons(1) - plon |
---|
| 536 | if (vec1_lon > pi) then |
---|
| 537 | src_lons(1) = src_lons(1) - pi2 |
---|
| 538 | else if (vec1_lon < -pi) then |
---|
| 539 | src_lons(1) = src_lons(1) + pi2 |
---|
| 540 | endif |
---|
| 541 | do n=2,4 |
---|
| 542 | vec1_lon = src_lons(n) - src_lons(1) |
---|
| 543 | if (vec1_lon > pi) then |
---|
| 544 | src_lons(n) = src_lons(n) - pi2 |
---|
| 545 | else if (vec1_lon < -pi) then |
---|
| 546 | src_lons(n) = src_lons(n) + pi2 |
---|
| 547 | endif |
---|
| 548 | end do |
---|
| 549 | |
---|
| 550 | corner_loop: do n=1,4 |
---|
| 551 | next_n = MOD(n,4) + 1 |
---|
| 552 | |
---|
| 553 | !*** |
---|
| 554 | !*** here we take the cross product of the vector making |
---|
| 555 | !*** up each box side with the vector formed by the vertex |
---|
| 556 | !*** and search point. if all the cross products are |
---|
| 557 | !*** positive, the point is contained in the box. |
---|
| 558 | !*** |
---|
| 559 | |
---|
| 560 | vec1_lat = src_lats(next_n) - src_lats(n) |
---|
| 561 | vec1_lon = src_lons(next_n) - src_lons(n) |
---|
| 562 | vec2_lat = plat - src_lats(n) |
---|
| 563 | vec2_lon = plon - src_lons(n) |
---|
| 564 | |
---|
| 565 | !*** |
---|
| 566 | !*** check for 0,2pi crossings |
---|
| 567 | !*** |
---|
| 568 | |
---|
| 569 | if (vec1_lon > three*pih) then |
---|
| 570 | vec1_lon = vec1_lon - pi2 |
---|
| 571 | else if (vec1_lon < -three*pih) then |
---|
| 572 | vec1_lon = vec1_lon + pi2 |
---|
| 573 | endif |
---|
| 574 | if (vec2_lon > three*pih) then |
---|
| 575 | vec2_lon = vec2_lon - pi2 |
---|
| 576 | else if (vec2_lon < -three*pih) then |
---|
| 577 | vec2_lon = vec2_lon + pi2 |
---|
| 578 | endif |
---|
| 579 | |
---|
| 580 | cross_product = vec1_lon*vec2_lat - vec2_lon*vec1_lat |
---|
| 581 | |
---|
| 582 | !*** |
---|
| 583 | !*** if cross product is less than zero, this cell |
---|
| 584 | !*** doesn't work |
---|
| 585 | !*** |
---|
| 586 | |
---|
| 587 | if (n == 1) cross_product_last = cross_product |
---|
| 588 | if (cross_product*cross_product_last < zero) |
---|
| 589 | & exit corner_loop |
---|
| 590 | cross_product_last = cross_product |
---|
| 591 | |
---|
| 592 | end do corner_loop |
---|
| 593 | |
---|
| 594 | !*** |
---|
| 595 | !*** if cross products all same sign, we found the location |
---|
| 596 | !*** |
---|
| 597 | |
---|
| 598 | if (n > 4) then |
---|
| 599 | src_add(1) = srch_add |
---|
| 600 | src_add(2) = e_add |
---|
| 601 | src_add(3) = ne_add |
---|
| 602 | src_add(4) = n_add |
---|
| 603 | |
---|
| 604 | ! Check if one point is masked; IF so, nearest-neighbour |
---|
| 605 | ! interpolation will be used |
---|
| 606 | |
---|
| 607 | ntotmask = 0 |
---|
| 608 | do ni=1,4 |
---|
| 609 | if (.not. grid1_mask(src_add(ni)).and. |
---|
| 610 | & .not. lextrapdone) |
---|
| 611 | & ntotmask = ntotmask + 1 |
---|
| 612 | end DO |
---|
| 613 | IF (ntotmask .gt. 0) src_add(1) = -src_add(1) |
---|
| 614 | |
---|
| 615 | return |
---|
| 616 | endif |
---|
| 617 | |
---|
| 618 | !*** |
---|
| 619 | !*** otherwise move on to next cell |
---|
| 620 | !*** |
---|
| 621 | |
---|
| 622 | endif !bounding box check |
---|
| 623 | end do srch_loop |
---|
| 624 | |
---|
| 625 | !*** |
---|
| 626 | !*** if no cell found, point is likely either in a box that straddles |
---|
| 627 | !*** either pole or is outside the grid. Put src_add = -1 so that |
---|
| 628 | !*** distance-weighted average of the 4 non-masked closest points |
---|
| 629 | !*** is done in calling routine. |
---|
| 630 | |
---|
| 631 | src_add = -1 |
---|
| 632 | |
---|
| 633 | !----------------------------------------------------------------------- |
---|
| 634 | |
---|
| 635 | end subroutine grid_search_bilin |
---|
| 636 | |
---|
| 637 | !*********************************************************************** |
---|
| 638 | |
---|
| 639 | subroutine store_link_bilin(dst_add, src_add, weights, nmap) |
---|
| 640 | |
---|
| 641 | !----------------------------------------------------------------------- |
---|
| 642 | ! |
---|
| 643 | ! this routine stores the address and weight for four links |
---|
| 644 | ! associated with one destination point in the appropriate address |
---|
| 645 | ! and weight arrays and resizes those arrays if necessary. |
---|
| 646 | ! |
---|
| 647 | !----------------------------------------------------------------------- |
---|
| 648 | |
---|
| 649 | !----------------------------------------------------------------------- |
---|
| 650 | ! |
---|
| 651 | ! input variables |
---|
| 652 | ! |
---|
| 653 | !----------------------------------------------------------------------- |
---|
| 654 | |
---|
| 655 | integer (kind=int_kind), intent(in) :: |
---|
| 656 | & dst_add, ! address on destination grid |
---|
| 657 | & nmap ! identifies which direction for mapping |
---|
| 658 | |
---|
| 659 | integer (kind=int_kind), dimension(4), intent(in) :: |
---|
| 660 | & src_add ! addresses on source grid |
---|
| 661 | |
---|
| 662 | real (kind=dbl_kind), dimension(4), intent(in) :: |
---|
| 663 | & weights ! array of remapping weights for these links |
---|
| 664 | |
---|
| 665 | !----------------------------------------------------------------------- |
---|
| 666 | ! |
---|
| 667 | ! local variables |
---|
| 668 | ! |
---|
| 669 | !----------------------------------------------------------------------- |
---|
| 670 | |
---|
| 671 | integer (kind=int_kind) :: n, ! dummy index |
---|
| 672 | & num_links_old ! placeholder for old link number |
---|
| 673 | |
---|
| 674 | !----------------------------------------------------------------------- |
---|
| 675 | ! |
---|
| 676 | ! increment number of links and check to see if remap arrays need |
---|
| 677 | ! to be increased to accomodate the new link. then store the |
---|
| 678 | ! link. |
---|
| 679 | ! |
---|
| 680 | !----------------------------------------------------------------------- |
---|
| 681 | |
---|
| 682 | select case (nmap) |
---|
| 683 | case(1) |
---|
| 684 | |
---|
| 685 | num_links_old = num_links_map1 |
---|
| 686 | num_links_map1 = num_links_old + 4 |
---|
| 687 | |
---|
| 688 | if (num_links_map1 > max_links_map1) |
---|
| 689 | & call resize_remap_vars(1,resize_increment) |
---|
| 690 | |
---|
| 691 | do n=1,4 |
---|
| 692 | grid1_add_map1(num_links_old+n) = src_add(n) |
---|
| 693 | grid2_add_map1(num_links_old+n) = dst_add |
---|
| 694 | wts_map1 (1,num_links_old+n) = weights(n) |
---|
| 695 | end do |
---|
| 696 | |
---|
| 697 | case(2) |
---|
| 698 | |
---|
| 699 | num_links_old = num_links_map2 |
---|
| 700 | num_links_map2 = num_links_old + 4 |
---|
| 701 | |
---|
| 702 | if (num_links_map2 > max_links_map2) |
---|
| 703 | & call resize_remap_vars(2,resize_increment) |
---|
| 704 | |
---|
| 705 | do n=1,4 |
---|
| 706 | grid1_add_map2(num_links_old+n) = dst_add |
---|
| 707 | grid2_add_map2(num_links_old+n) = src_add(n) |
---|
| 708 | wts_map2 (1,num_links_old+n) = weights(n) |
---|
| 709 | end do |
---|
| 710 | |
---|
| 711 | end select |
---|
| 712 | |
---|
| 713 | !----------------------------------------------------------------------- |
---|
| 714 | |
---|
| 715 | end subroutine store_link_bilin |
---|
| 716 | |
---|
| 717 | !*********************************************************************** |
---|
| 718 | |
---|
| 719 | end module remap_bilinear |
---|
| 720 | |
---|
| 721 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|