00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_neigh_nearest_3d_real (grid_id, &
00012 coords1, coords2, coords3, &
00013 x_coords, y_coords, z_coords, coords_shape, &
00014 mask_array, mask_shape, mask_available, &
00015 sin_values, cos_values, &
00016 grid_valid_shape, &
00017 srcloc, nloc, nprev, nsearch, &
00018 neighbors_3d, num_neigh, &
00019 extra_search, ierror)
00020
00021
00022
00023 use PRISM_constants
00024
00025 use PSMILe, dummy_interface => PSMILe_Neigh_nearest_3d_real
00026
00027 Implicit none
00028
00029
00030
00031 Integer, Intent (In) :: grid_id
00032
00033
00034
00035
00036 Integer, Intent (In) :: nloc
00037
00038
00039
00040 Integer, Intent (In) :: nprev
00041
00042
00043
00044 Integer, Intent (In) :: nsearch
00045
00046
00047
00048
00049 Integer, Intent (In) :: srcloc (ndim_3d, nloc)
00050
00051
00052
00053
00054
00055
00056 Real, Intent (In) :: coords1 (nloc)
00057 Real, Intent (In) :: coords2 (nloc)
00058 Real, Intent (In) :: coords3 (nloc)
00059
00060
00061
00062 Integer, Intent (In) :: coords_shape (2, ndim_3d)
00063
00064
00065
00066 Real, Intent (In) ::
00067 x_coords( coords_shape(1,1):coords_shape(2,1),
00068 coords_shape(1,2):coords_shape(2,2),
00069 coords_shape(1,3):coords_shape(2,3))
00070 Real, Intent (In) ::
00071 y_coords( coords_shape(1,1):coords_shape(2,1),
00072 coords_shape(1,2):coords_shape(2,2),
00073 coords_shape(1,3):coords_shape(2,3))
00074 Real, Intent (In) ::
00075 z_coords( coords_shape(1,1):coords_shape(2,1),
00076 coords_shape(1,2):coords_shape(2,2),
00077 coords_shape(1,3):coords_shape(2,3))
00078
00079
00080
00081 Integer, Intent (In) :: mask_shape (2, ndim_3d)
00082
00083
00084
00085 Logical, Intent (In) ::
00086 mask_array ( mask_shape (1,1):mask_shape (2,1),
00087 mask_shape (1,2):mask_shape (2,2),
00088 mask_shape (1,3):mask_shape (2,3))
00089
00090
00091 Logical, Intent (In) :: mask_available
00092
00093
00094
00095
00096
00097 Integer, Intent (In) :: grid_valid_shape (2, ndim_3d)
00098
00099
00100
00101 Real, Intent (In) ::
00102 sin_values ( grid_valid_shape(1,1):grid_valid_shape(2,1),
00103 grid_valid_shape(1,2):grid_valid_shape(2,2),
00104 grid_valid_shape(1,3):grid_valid_shape(2,3), 2)
00105
00106
00107
00108 Real, Intent (In) ::
00109 cos_values ( grid_valid_shape(1,1):grid_valid_shape(2,1),
00110 grid_valid_shape(1,2):grid_valid_shape(2,2),
00111 grid_valid_shape(1,3):grid_valid_shape(2,3), 2)
00112
00113
00114
00115 Integer, Intent (In) :: num_neigh
00116
00117
00118
00119
00120
00121
00122 Type (Extra_search_info) :: extra_search
00123
00124
00125
00126
00127
00128
00129 Integer, Intent (Out) :: neighbors_3d (ndim_3d, nloc, num_neigh)
00130
00131
00132
00133 integer, Intent (Out) :: ierror
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148 Integer, Parameter :: nrefd = 4 * 4 * 4
00149 Integer, Parameter :: lon = 1
00150 Integer, Parameter :: lat = 2
00151
00152 Real, Parameter :: real_earth = 6400.0
00153
00154 Real, Parameter :: tol = 1d-6
00155 Real, Parameter :: eps = 1d-6
00156 Real, Parameter :: acosp1 = 1.0
00157 Real, Parameter :: acosm1 = -1.0
00158
00159
00160
00161
00162
00163 Integer :: i, ibeg, iend, n
00164
00165
00166
00167 Real, Pointer :: sin_search (:, :)
00168 Real, Pointer :: cos_search (:, :)
00169
00170 Integer :: ijkctl (ndim_3d, nrefd)
00171 Integer :: ijkdst (ndim_3d, nrefd)
00172 Integer :: irange (2, ndim_3d)
00173 Integer :: itemp (ndim_3d)
00174
00175 Logical :: mask_dist (nrefd)
00176 Real :: distance (nrefd), dist, fac, val
00177 Real :: dist_control
00178
00179 Integer :: ijk (ndim_3d), ijk0 (ndim_3d)
00180 Integer :: leni, lenij
00181 Integer :: ii, jj, kk, nref
00182
00183 Integer :: nfound
00184 Integer :: nmin (1)
00185
00186
00187
00188 Real :: real_huge
00189
00190
00191
00192 Integer, Parameter :: nerrp = 2
00193 Integer :: ierrp (nerrp)
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215 Character(len=len_cvs_string), save :: mycvs =
00216 '$Id: psmile_neigh_nearest_3d_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00217
00218
00219
00220
00221
00222 #ifdef OLD
00223 Integer, Parameter :: nrefd = 3 * 3 * 3
00224
00225 data ((ijkctl (i, n), i=1,ndim_3d), n = 1,nrefd) &
00226 / 0, 0, 0, &
00227 1, 0, 0, 0, 1, 0, 0, 0, 1, &
00228 -1, 0, 0, 0,-1, 0, 0, 0,-1, &
00229 -1,-1, 0, 1,-1, 0, &
00230 -1, 1, 0, 1, 1, 0, &
00231 -1, 0,-1, 1, 0,-1, &
00232 -1, 0, 1, 1, 0, 1, &
00233 0,-1,-1, 0, 1,-1, &
00234 0,-1, 1, 0, 1, 1, &
00235 -1,-1,-1, 1,-1,-1, &
00236 -1, 1,-1, 1, 1,-1, &
00237 -1,-1, 1, 1,-1, 1, &
00238 -1, 1, 1, 1, 1, 1/
00239 #endif
00240
00241
00242
00243 data ((ijkctl (i, n), i=1,ndim_3d), n = 1,nrefd) &
00244 / 0, 0, 0, 1, 0, 0, 2, 0, 0, 3, 0, 0, &
00245 0, 1, 0, 1, 1, 0, 2, 1, 0, 3, 1, 0, &
00246 0, 2, 0, 1, 2, 0, 2, 2, 0, 3, 2, 0, &
00247 0, 3, 0, 1, 3, 0, 2, 3, 0, 3, 3, 0, &
00248 0, 0, 1, 1, 0, 1, 2, 0, 1, 3, 0, 1, &
00249 0, 1, 1, 1, 1, 1, 2, 1, 1, 3, 1, 1, &
00250 0, 2, 1, 1, 2, 1, 2, 2, 1, 3, 2, 1, &
00251 0, 3, 1, 1, 3, 1, 2, 3, 1, 3, 3, 1, &
00252 0, 0, 2, 1, 0, 2, 2, 0, 2, 3, 0, 2, &
00253 0, 1, 2, 1, 1, 2, 2, 1, 2, 3, 1, 2, &
00254 0, 2, 2, 1, 2, 2, 2, 2, 2, 3, 2, 2, &
00255 0, 3, 2, 1, 3, 2, 2, 3, 2, 3, 3, 2, &
00256 0, 0, 3, 1, 0, 3, 2, 0, 3, 3, 0, 3, &
00257 0, 1, 3, 1, 1, 3, 2, 1, 3, 3, 1, 3, &
00258 0, 2, 3, 1, 2, 3, 2, 2, 3, 3, 2, 3, &
00259 0, 3, 3, 1, 3, 3, 2, 3, 3, 3, 3, 3/
00260
00261
00262
00263
00264
00265 #ifdef VERBOSE
00266 print 9990, trim(ch_id)
00267
00268 call psmile_flushstd
00269 #endif /* VERBOSE */
00270
00271 ibeg = nprev + 1
00272 iend = nprev + nsearch
00273
00274 #ifdef PRISM_ASSERTION
00275
00276 do i = ibeg, iend
00277
00278 if (srcloc(1,i) < coords_shape(1,1) .or. &
00279 srcloc(1,i) > coords_shape(2,1) .or. &
00280 srcloc(2,i) < coords_shape(1,2) .or. &
00281 srcloc(2,i) > coords_shape(2,2) .or. &
00282 srcloc(3,i) < coords_shape(1,3) .or. &
00283 srcloc(3,i) > coords_shape(2,3)) exit
00284 end do
00285
00286 if (i <= iend) then
00287 print *, i, srcloc(:, i), coords_shape
00288 call psmile_assert (__FILE__, __LINE__, &
00289 'wrong index')
00290 endif
00291 #endif
00292
00293 real_huge = huge (distance(1))
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306 Allocate (sin_search(ibeg:iend, lat), STAT = ierror)
00307
00308 if ( ierror > 0 ) then
00309 ierrp (1) = ierror
00310 ierrp (2) = nsearch * lat
00311
00312 ierror = PRISM_Error_Alloc
00313 call psmile_error ( ierror, 'sin_search', &
00314 ierrp, 2, __FILE__, __LINE__ )
00315 return
00316 endif
00317
00318 Allocate (cos_search(ibeg:iend, lat), STAT = ierror)
00319
00320 if ( ierror > 0 ) then
00321 ierrp (1) = ierror
00322 ierrp (2) = nsearch * lat
00323
00324 ierror = PRISM_Error_Alloc
00325 call psmile_error ( ierror, 'cos_search', &
00326 ierrp, 2, __FILE__, __LINE__ )
00327 return
00328 endif
00329
00330 sin_search (ibeg:iend, lon) = coords1 (ibeg:iend) * real_deg2rad
00331 sin_search (ibeg:iend, lat) = coords2 (ibeg:iend) * real_deg2rad
00332
00333 cos_search = cos (sin_search)
00334 sin_search = sin (sin_search)
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357 ijk0 = Grids(grid_id)%ijk0
00358
00359 do i = ibeg, iend
00360
00361
00362
00363 ijk = ijk0 + ((srcloc(:,i) - ijk0) / 4) * 4
00364
00365
00366
00367 irange (1, :) = max (ijk, grid_valid_shape (1, :))
00368 irange (2, :) = min (ijk+3, grid_valid_shape (2, :))
00369
00370 nref = (irange (2,1) - irange (1,1) + 1) * &
00371 (irange (2,2) - irange (1,2) + 1) * &
00372 (irange (2,3) - irange (1,3) + 1)
00373
00374 if (nref == nrefd) then
00375 do n = 1, nrefd
00376 ijkdst (:, n) = ijk + ijkctl (:, n)
00377 end do
00378 else
00379
00380 leni = irange (2, 1) - irange (1, 1) + 1
00381 lenij = (irange (2, 2) - irange (1, 2) + 1) * leni
00382
00383 do kk = 1, irange (2,3)-irange(1,3) + 1
00384 n = (kk-1)*lenij
00385 do jj = 1, irange (2,2)-irange(1,2) + 1
00386
00387 do ii = 1, leni
00388 ijkdst (1, n+ii) = irange(1,1) + ii - 1
00389 end do
00390
00391 ijkdst (2, n+1:n+leni) = irange (1,2) + jj - 1
00392 n = n + leni
00393 end do
00394
00395 ijkdst (3, (kk-1)*lenij+1:kk*lenij) = irange (1,3) + kk - 1
00396 end do
00397 endif
00398
00399
00400
00401 if (mask_available) then
00402
00403 do n = 1, nref
00404 mask_dist (n) = mask_array (ijkdst (1,n), ijkdst (2,n), &
00405 ijkdst (3,n))
00406 end do
00407
00408 n = count (mask_dist(1:nref))
00409
00410 if (n == 0) then
00411
00412
00413 nref = 0
00414 go to 10
00415 endif
00416
00417 if (n < nref) then
00418 n = 0
00419 do ii = 1, nref
00420 if (mask_dist (ii)) then
00421 n = n + 1
00422 if (n /= ii) ijkdst (:, n) = ijkdst (:, ii)
00423 endif
00424 end do
00425 nref = n
00426 endif
00427 endif
00428
00429
00430
00431 fac = real_earth + coords3 (i)
00432
00433 do n = 1, nref
00434 val = ( sin_values(ijkdst (1,n), &
00435 ijkdst (2,n), &
00436 ijkdst (3,n), lat) * sin_search(i, lat) &
00437 + cos_values(ijkdst (1,n), &
00438 ijkdst (2,n), &
00439 ijkdst (3,n), lat) * cos_search(i, lat) &
00440 * (cos_values(ijkdst (1,n), &
00441 ijkdst (2,n), &
00442 ijkdst (3,n), lon) * cos_search(i, lon) &
00443 +sin_values(ijkdst (1,n), &
00444 ijkdst (2,n), &
00445 ijkdst (3,n), lon) * sin_search(i, lon)) )
00446 #ifdef PRISM_ASSERTION
00447 if (val < acosm1 - eps .or. val > acosp1 + eps) then
00448 print *, i, val
00449 call psmile_assert (__FILE__, __LINE__, &
00450 'invalid value for acos')
00451 endif
00452 #endif
00453 val = min (acosp1, val)
00454 val = max (acosm1, val)
00455
00456 dist = fac * acos (val)
00457
00458 distance (n) = dist*dist + &
00459 ( coords3 (i)-z_coords(ijkdst (1,n), &
00460 ijkdst (2,n), &
00461 ijkdst (3,n)))**2
00462 #ifdef DEBUG2
00463 if (i == 2) then
00464 print *, 'dist, diff', dist, &
00465 coords3 (i)-z_coords(ijkdst (1,n), ijkdst (2,n), ijkdst (3,n))
00466 endif
00467 #endif
00468 end do
00469
00470 nmin = minloc (distance(1:nref))
00471 #ifdef MINLOCFIX
00472 if (nmin(1).eq.0) nmin=1
00473 #endif /* MINLOCFIX */
00474
00475 #ifdef DEBUG2
00476 if (i == 2) then
00477 print *, 'psmile_neigh_nearest_3d_real: ictl =', i
00478 print *, 'ijkdst and distance:'
00479 do n = 1, nref
00480 print *, ijkdst (:, n), distance (n)
00481 end do
00482 endif
00483 #endif
00484
00485
00486
00487 if (distance(nmin(1)) <= tol) then
00488 neighbors_3d (:, i, 1) = ijkdst (:, nmin(1))
00489
00490 do n = 2, num_neigh
00491 neighbors_3d (:, i, n) = grid_valid_shape (1, :) - 1
00492 end do
00493
00494 cycle
00495 endif
00496
00497
00498
00499
00500 10 continue
00501
00502 nfound = min (num_neigh, nref)
00503 if (nfound < num_neigh) then
00504 distance(nfound+1:num_neigh) = real_huge
00505
00506 do n = nfound+1, num_neigh
00507 ijkdst (:, n) = grid_valid_shape (1, :) - 1
00508 end do
00509 endif
00510
00511 do n = 1, nfound
00512
00513 nmin = minloc (distance(n:nref))
00514 #ifdef MINLOCFIX
00515 if (nmin(1).eq.0) nmin=1
00516 #endif /* MINLOCFIX */
00517
00518 if (nmin(1) /= 1) then
00519 jj = n + nmin (1) - 1
00520
00521 itemp (:) = ijkdst (:, n)
00522 ijkdst (:, n) = ijkdst (:, jj)
00523 ijkdst (:, jj) = itemp (:)
00524
00525 dist = distance (n)
00526 distance (n) = distance (jj)
00527 distance (jj) = dist
00528 endif
00529 end do
00530
00531 #ifdef PRISM_ASSERTION
00532 do n = 1, nfound-1
00533 if (distance(n) > distance (n+1)) exit
00534 end do
00535
00536 if (n <= nfound-1) then
00537 print *, 'n =' , n
00538 print *, 'distance (n) ', distance (n)
00539 print *, 'distance (n+1)', distance (n+1)
00540 print *, 'distance', distance (1:nfound)
00541 call psmile_assert (__FILE__, __LINE__, &
00542 'incorrect sort')
00543 endif
00544 #endif
00545
00546
00547
00548 dist_control = distance (num_neigh)
00549
00550 #ifdef TODO
00551 call psmile_mg_search_nneigh_3d_real (distance, ijkdst, num_neigh, &
00552 grid_id, ierror)
00553 #endif /* TODO */
00554
00555
00556
00557 do n = 1, num_neigh
00558 neighbors_3d (:, i, n) = ijkdst (:, n)
00559 end do
00560
00561 end do
00562
00563
00564
00565 Deallocate (cos_search)
00566 Deallocate (sin_search)
00567
00568
00569
00570 #ifdef VERBOSE
00571 print 9980, trim(ch_id), ierror
00572
00573 call psmile_flushstd
00574 #endif /* VERBOSE */
00575
00576
00577
00578 9990 format (1x, a, ': psmile_neigh_nearest_3d_real')
00579 9980 format (1x, a, ': psmile_neigh_nearest_3d_real: eof, ierror =', i3)
00580
00581 end subroutine PSMILe_Neigh_nearest_3d_real