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