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