00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_neigh_nearx_3d_irr2_dble ( &
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, srclocz, nlocs, &
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_nearx_3d_irr2_dble
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) :: nlocs (ndim_2d)
00049
00050
00051
00052 Integer, Intent (In) :: srcloc (ndim_2d, nlocs (1))
00053 Integer, Intent (In) :: srclocz (nlocs (2))
00054
00055
00056
00057 Double Precision, Intent (In) :: coords1 (nlocs(1),nlocs(2))
00058 Double Precision, Intent (In) :: coords2 (nlocs(1),nlocs(2))
00059 Double Precision, Intent (In) :: coords3 (nlocs(1),nlocs(2))
00060
00061
00062
00063
00064
00065 Integer, Intent (In) :: coords_shape (2, ndim_3d)
00066
00067
00068
00069 Double Precision, Intent (In) :: x_coords(coords_shape(1,1):
00070 coords_shape(2,1))
00071 Double Precision, Intent (In) :: y_coords(coords_shape(1,2):
00072 coords_shape(2,2))
00073 Double Precision, Intent (In) :: z_coords(coords_shape(1,3):
00074 coords_shape(2,3))
00075
00076
00077
00078 Integer, Intent (In) :: mask_shape (2, ndim_3d)
00079
00080
00081
00082 Logical, Intent (In) :: mask_array (mask_shape (1,1):
00083 mask_shape (2,1),
00084 mask_shape (1,2):
00085 mask_shape (2,2),
00086 mask_shape (1,3):
00087 mask_shape (2,3))
00088
00089
00090 Logical, Intent (In) :: mask_available
00091
00092
00093
00094
00095
00096 Integer, Intent (In) :: grid_valid_shape (2, ndim_3d)
00097
00098
00099
00100 Double Precision, Intent (In) :: sin_values_lon (grid_valid_shape(1,1):
00101 grid_valid_shape(2,1))
00102 Double Precision, Intent (In) :: sin_values_lat (grid_valid_shape(1,2):
00103 grid_valid_shape(2,2))
00104
00105
00106
00107 Double Precision, Intent (In) :: cos_values_lon (grid_valid_shape(1,1):
00108 grid_valid_shape(2,1))
00109 Double Precision, Intent (In) :: cos_values_lat (grid_valid_shape(1,2):
00110 grid_valid_shape(2,2))
00111
00112
00113
00114 Integer, Intent (In) :: search_mode
00115
00116
00117
00118
00119
00120
00121 Integer, Intent (In) :: num_neigh
00122
00123
00124
00125
00126
00127
00128 Type (Extra_search_info) :: extra_search
00129
00130
00131
00132
00133
00134
00135 Integer, Intent (Out) :: neighbors_3d (ndim_3d, nloc, num_neigh)
00136
00137
00138
00139 Integer, Intent (Out) :: ierror
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150 Integer, Parameter :: lon = 1
00151 Integer, Parameter :: lat = 2
00152
00153
00154
00155
00156
00157 Double Precision, Pointer :: sin_search (:, :)
00158 Double Precision, Pointer :: cos_search (:, :)
00159 Double Precision, Pointer :: z_search (:)
00160
00161
00162
00163 Integer :: j
00164 Integer :: ibeg, iend
00165 Integer :: jbeg, jend
00166 Integer :: lenij
00167
00168 Integer, Pointer :: indices (:)
00169 Integer :: nprev1
00170
00171
00172
00173 Integer :: ijk0 (ndim_3d)
00174 Integer, Pointer :: ijk (:, :)
00175 Integer :: width (ndim_3d)
00176
00177
00178
00179 Integer, Parameter :: nerrp = 2
00180 Integer :: ierrp (nerrp)
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217 Character(len=len_cvs_string), save :: mycvs =
00218 '$Id: psmile_neigh_nearx_3d_irr2_dble.F90 2966 2011-02-18 09:47:30Z hanke $'
00219
00220
00221
00222
00223
00224 #ifdef VERBOSE
00225 print 9990, trim(ch_id)
00226
00227 call psmile_flushstd
00228 #endif /* VERBOSE */
00229
00230 ierror = 0
00231
00232 nprev1 = nprev + 1
00233
00234 ibeg = nprev + 1
00235 iend = nprev + nlocs(1)*nlocs(2)
00236
00237 indices => extra_search%indices
00238
00239
00240
00241
00242
00243
00244 do jbeg = 1, extra_search%n_extra
00245 if (indices(jbeg) >= ibeg) exit
00246 end do
00247
00248 if (jbeg > extra_search%n_extra) then
00249 #ifdef VERBOSE
00250 print 9980, trim(ch_id), ierror
00251 call psmile_flushstd
00252 #endif /* VERBOSE */
00253 return
00254 endif
00255
00256 do jend = jbeg, extra_search%n_extra
00257 if (indices(jend) > iend) exit
00258 end do
00259 jend = jend - 1
00260
00261 if (jbeg > jend) then
00262 #ifdef VERBOSE
00263 print 9980, trim(ch_id), ierror
00264 call psmile_flushstd
00265 #endif /* VERBOSE */
00266 return
00267 endif
00268
00269 ibeg = indices(jbeg)
00270 iend = indices(jend)
00271
00272 #ifdef PRISM_ASSERTION
00273
00274
00275
00276 if (search_mode < 1 .or. search_mode > ndim_3d) then
00277 print *, trim(ch_id), " search_mode ", search_mode
00278 call psmile_assert (__FILE__, __LINE__, &
00279 'search_mode is out of range ')
00280 endif
00281
00282 if (nloc < nprev + nlocs (1) * nlocs (2)) then
00283 print *, trim(ch_id), " nloc, nprev, nlocs ", nloc, nprev, nlocs
00284 call psmile_assert (__FILE__, __LINE__, &
00285 'nloc < nprev + PRODUCT (nlocs) ')
00286 endif
00287
00288
00289
00290
00291 do j = 1, nlocs (1)
00292
00293 if (srcloc(1,j) < grid_valid_shape(1,1) - 1 .or. &
00294 srcloc(1,j) > grid_valid_shape(2,1) .or. &
00295 srcloc(2,j) < grid_valid_shape(1,2) - 1 .or. &
00296 srcloc(2,j) > grid_valid_shape(2,2)) exit
00297 end do
00298
00299 if (j <= nlocs(1)) then
00300 print *, "Incorrect index in srcloc, j", j, srcloc(:, j)
00301 call psmile_assert (__FILE__, __LINE__, &
00302 'wrong index in srcloc')
00303 endif
00304
00305
00306 do j = 1, nlocs (2)
00307 if (srclocz(j) < grid_valid_shape(1,3) - 1 .or. &
00308 srclocz(j) > grid_valid_shape(2,3)) exit
00309 end do
00310
00311 if (j <= nlocs(2)) then
00312 print *, "Incorrect index in srclocz(j", j, srclocz(j)
00313 call psmile_assert (__FILE__, __LINE__, &
00314 'wrong index in srclocz')
00315 endif
00316 #endif
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326
00327
00328
00329 Allocate (sin_search(jbeg:jend, lat), STAT = ierror)
00330
00331 if ( ierror > 0 ) then
00332 ierrp (1) = ierror
00333 ierrp (2) = (jend-jbeg+1) * lat
00334
00335 ierror = PRISM_Error_Alloc
00336 call psmile_error ( ierror, 'sin_search', &
00337 ierrp, 2, __FILE__, __LINE__ )
00338 return
00339 endif
00340
00341 Allocate (cos_search(jbeg:jend, lat), STAT = ierror)
00342
00343 if ( ierror > 0 ) then
00344 ierrp (1) = ierror
00345 ierrp (2) = (jend-jbeg+1) * lat
00346
00347 ierror = PRISM_Error_Alloc
00348 call psmile_error ( ierror, 'cos_search', &
00349 ierrp, 2, __FILE__, __LINE__ )
00350 return
00351 endif
00352
00353 Allocate (z_search(jbeg:jend), STAT = ierror)
00354
00355 if ( ierror > 0 ) then
00356 ierrp (1) = ierror
00357 ierrp (2) = (jend-jbeg+1)
00358
00359 ierror = PRISM_Error_Alloc
00360 call psmile_error ( ierror, 'z_search', &
00361 ierrp, 2, __FILE__, __LINE__ )
00362 return
00363 endif
00364
00365 Allocate (ijk(ndim_3d, jbeg:jend), STAT = ierror)
00366
00367 if ( ierror > 0 ) then
00368 ierrp (1) = ierror
00369 ierrp (2) = (jend-jbeg+1) * ndim_3d
00370
00371 ierror = PRISM_Error_Alloc
00372 call psmile_error ( ierror, 'ijk', &
00373 ierrp, 2, __FILE__, __LINE__ )
00374 return
00375 endif
00376
00377
00378
00379 lenij = nlocs(1)
00380
00381 ijk (3, jbeg:jend) = (indices(jbeg:jend)-nprev1)/lenij + 1
00382 ijk (1, jbeg:jend) = mod(indices(jbeg:jend)-nprev1, nlocs(1)) + 1
00383
00384 #ifdef PRISM_ASSERTION
00385
00386 do j = jbeg, jend
00387 if (ijk(1,j) < 1 .or. ijk(1,j) > nlocs(1) .or. &
00388 ijk(3,j) < 1 .or. ijk(3,j) > nlocs(2)) exit
00389 end do
00390
00391 if (j <= jend) then
00392 print *, "Incorrect index in ijk", j, jbeg, jend, ijk(:, j)
00393 call psmile_assert (__FILE__, __LINE__, &
00394 'wrong index in ijk')
00395 endif
00396 #endif
00397
00398
00399
00400
00401
00402
00403 do j = jbeg, jend
00404 sin_search (j, lon) = coords1 (ijk(1,j),1) * dble_deg2rad
00405 end do
00406
00407
00408 do j = jbeg, jend
00409 sin_search (j, lat) = coords2 (ijk(1,j),1) * dble_deg2rad
00410 end do
00411
00412 cos_search = cos (sin_search)
00413 sin_search = sin (sin_search)
00414
00415
00416 do j = jbeg, jend
00417 z_search (j) = coords3 (1,ijk(3,j))
00418 end do
00419
00420
00421
00422
00423
00424 ijk0 = Grids(grid_id)%ijk0
00425
00426
00427
00428
00429 width (1:search_mode) = 3
00430 width (search_mode+1:ndim_3d) = 0
00431
00432 ijk (1, jbeg:jend) = ijk0(1) + ((srcloc(1, ijk(1, jbeg:jend)) - &
00433 ijk0(1)) / 4) * 4
00434
00435 if (search_mode >= 2) then
00436 ijk (2, jbeg:jend) = ijk0(2) + ((srcloc(2, ijk(1, jbeg:jend)) - &
00437 ijk0(2)) / 4) * 4
00438 else
00439 ijk (2, jbeg:jend) = srcloc(2, ijk(1, jbeg:jend))
00440 endif
00441
00442 if (search_mode >= 3) then
00443 ijk (3, jbeg:jend) = ijk0(3) + ((srclocz(ijk(3, jbeg:jend)) - &
00444 ijk0(3)) / 4) * 4
00445 else
00446 ijk (3, jbeg:jend) = srclocz(ijk(3, jbeg:jend))
00447 endif
00448
00449
00450
00451
00452 call psmile_neigh_nearx_sub_reg_dble ( &
00453 grid_id, &
00454 x_coords, y_coords, z_coords, &
00455 coords_shape, &
00456 mask_array, mask_shape, mask_available, &
00457 sin_values_lon, cos_values_lon, &
00458 sin_values_lat, cos_values_lat, &
00459 grid_valid_shape, search_mode, &
00460 neighbors_3d, num_neigh, nloc, &
00461 extra_search, &
00462 ijk, sin_search, cos_search, &
00463 z_search, jbeg, jend, width, ierror)
00464
00465
00466
00467 Deallocate (ijk)
00468 Deallocate (z_search)
00469 Deallocate (cos_search)
00470 Deallocate (sin_search)
00471
00472
00473
00474 #ifdef VERBOSE
00475 print 9980, trim(ch_id), ierror
00476
00477 call psmile_flushstd
00478 #endif /* VERBOSE */
00479
00480
00481
00482 9990 format (1x, a, ': psmile_neigh_nearx_3d_irr2_dble')
00483 9980 format (1x, a, ': psmile_neigh_nearx_3d_irr2_dble: eof, ierror =', i3)
00484
00485 end subroutine PSMILe_Neigh_nearx_3d_irr2_dble