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