00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_neigh_nearx_irreg2_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, 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_nearx_irreg2_dble
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 Integer, Intent (In) :: nprev
00048
00049
00050
00051 Integer, Intent (In) :: srcloc (ndim_2d, nlocs (1))
00052 Integer, Intent (In) :: srclocz (nlocs (2))
00053
00054
00055
00056 Double Precision, Intent (In) :: coords1 (nlocs(1))
00057 Double Precision, Intent (In) :: coords2 (nlocs(1))
00058 Double Precision, Intent (In) :: coords3 (nlocs(2))
00059
00060
00061
00062 Integer, Intent (In) :: coords_shape (2, ndim_3d)
00063
00064
00065
00066 Double Precision, Intent (In) :: x_coords(coords_shape(1,1):
00067 coords_shape(2,1),
00068 coords_shape(1,2):
00069 coords_shape(2,2))
00070 Double Precision, Intent (In) :: y_coords(coords_shape(1,1):
00071 coords_shape(2,1),
00072 coords_shape(1,2):
00073 coords_shape(2,2))
00074 Double Precision, Intent (In) :: z_coords(coords_shape(1,3):
00075 coords_shape(2,3))
00076
00077
00078
00079 Integer, Intent (In) :: mask_shape (2, ndim_3d)
00080
00081
00082
00083 Logical, Intent (In) :: mask_array (mask_shape (1,1):
00084 mask_shape (2,1),
00085 mask_shape (1,2):
00086 mask_shape (2,2),
00087 mask_shape (1,3):
00088 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 Double Precision, Intent (In) :: sin_values (grid_valid_shape(1,1):
00102 grid_valid_shape(2,1),
00103 grid_valid_shape(1,2):
00104 grid_valid_shape(2,2),
00105 2)
00106
00107
00108
00109 Double Precision, Intent (In) :: cos_values (grid_valid_shape(1,1):
00110 grid_valid_shape(2,1),
00111 grid_valid_shape(1,2):
00112 grid_valid_shape(2,2),
00113 2)
00114
00115
00116
00117 Integer, Intent (In) :: search_mode
00118
00119
00120
00121
00122
00123 Integer, Intent (In) :: num_neigh
00124
00125
00126
00127
00128
00129
00130 Type (Extra_search_info) :: extra_search
00131
00132
00133
00134
00135
00136
00137 Integer, Intent (Out) :: neighbors_3d (ndim_3d, nloc, num_neigh)
00138
00139
00140
00141 Integer, Intent (Out) :: ierror
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152 Integer, Parameter :: lon = 1
00153 Integer, Parameter :: lat = 2
00154
00155
00156
00157
00158
00159 Double Precision, Pointer :: sin_search (:, :)
00160 Double Precision, Pointer :: cos_search (:, :)
00161 Double Precision, Pointer :: z_search (:)
00162
00163
00164
00165 Integer :: i, ibeg, iend
00166 Integer :: j, jbeg, jend
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 Character(len=len_cvs_string), save :: mycvs =
00211 '$Id: psmile_neigh_nearx_irreg2_dble.F90 2325 2010-04-21 15:00:07Z valcke $'
00212
00213
00214
00215
00216
00217 #ifdef VERBOSE
00218 print 9990, trim(ch_id)
00219
00220 call psmile_flushstd
00221 #endif /* VERBOSE */
00222
00223 ierror = 0
00224
00225 nprev1 = nprev + 1
00226
00227 ibeg = nprev + 1
00228 iend = nprev + nlocs(1)*nlocs(2)
00229
00230 indices => extra_search%indices
00231
00232
00233
00234
00235
00236
00237 do jbeg = 1, extra_search%n_extra
00238 if (indices(jbeg) >= ibeg) exit
00239 end do
00240
00241 if (jbeg > extra_search%n_extra) return
00242
00243 do jend = jbeg, extra_search%n_extra
00244 if (indices(jend) > iend) exit
00245 end do
00246 jend = jend - 1
00247
00248 if (jbeg > jend) return
00249
00250 ibeg = indices(jbeg)
00251 iend = indices(jend)
00252
00253 #ifdef PRISM_ASSERTION
00254
00255
00256
00257 if (search_mode < 2 .or. search_mode > ndim_3d) then
00258 print *, trim(ch_id), " search_mode ", search_mode
00259 call psmile_assert (__FILE__, __LINE__, &
00260 'search_mode is out of range ')
00261 endif
00262
00263 if (nloc < nprev + nlocs (1) * nlocs (2) ) then
00264 print *, trim(ch_id), " nloc, nprev, nlocs ", nloc, nprev, nlocs
00265 call psmile_assert (__FILE__, __LINE__, &
00266 'nloc < nprev + PRODUCT (nlocs) ')
00267 endif
00268
00269
00270
00271
00272
00273
00274 do i = 1, nlocs (1)
00275
00276 if (srcloc(1,i) < grid_valid_shape(1,1)-1 .or. &
00277 srcloc(1,i) > grid_valid_shape(2,1) .or. &
00278 srcloc(2,i) < grid_valid_shape(1,2)-1 .or. &
00279 srcloc(2,i) > grid_valid_shape(2,2)) exit
00280 end do
00281
00282 if (i <= nlocs(1)) then
00283 print *, "Incorrect index in srcloc, i", i, srcloc(:, i)
00284 print *, "grid_valid_shape: ", grid_valid_shape
00285 call psmile_assert (__FILE__, __LINE__, &
00286 'wrong index in srcloc')
00287 endif
00288
00289
00290 do i = 1, nlocs (2)
00291 if (srclocz(i) < grid_valid_shape(1,3)-1 .or. &
00292 srclocz(i) > grid_valid_shape(2,3)) exit
00293 end do
00294
00295 if (i <= nlocs(2)) then
00296 print *, "Incorrect index in srclocz(2), i", i, srclocz(i)
00297 call psmile_assert (__FILE__, __LINE__, &
00298 'wrong index in srclocz')
00299 endif
00300 #endif
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314 Allocate (sin_search(jbeg:jend, lat), &
00315 cos_search(jbeg:jend, lat), &
00316 z_search(jbeg:jend), STAT = ierror)
00317
00318 if ( ierror > 0 ) then
00319 ierrp (1) = ierror
00320 ierrp (2) = (jend-jbeg+1) * (lat * 2 + 1)
00321
00322 ierror = PRISM_Error_Alloc
00323 call psmile_error ( ierror, 'sin_search, ..., z_search', &
00324 ierrp, 2, __FILE__, __LINE__ )
00325 return
00326 endif
00327
00328 Allocate (ijk(ndim_3d, jbeg:jend), STAT = ierror)
00329
00330 if ( ierror > 0 ) then
00331 ierrp (1) = ierror
00332 ierrp (2) = (jend-jbeg+1) * ndim_3d
00333
00334 ierror = PRISM_Error_Alloc
00335 call psmile_error ( ierror, 'ijk', &
00336 ierrp, 2, __FILE__, __LINE__ )
00337 return
00338 endif
00339
00340
00341
00342 ijk (3, jbeg:jend) = (indices(jbeg:jend)-nprev1)/nlocs(1) + 1
00343 ijk (2, jbeg:jend) = mod(indices(jbeg:jend)-nprev1, nlocs(1)) + 1
00344
00345 #ifdef PRISM_ASSERTION
00346
00347 do j = jbeg, jend
00348 if (ijk(2,j) < 1 .or. ijk(2,j) > nlocs(1) .or. &
00349 ijk(3,j) < 1 .or. ijk(3,j) > nlocs(2)) exit
00350 end do
00351
00352 if (j <= jend) then
00353 print *, "Incorrect index in ijk", j, jbeg, jend, ijk(2:3, j)
00354 call psmile_assert (__FILE__, __LINE__, &
00355 'wrong index in ijk')
00356 endif
00357 #endif
00358
00359
00360
00361
00362 do j = jbeg, jend
00363 sin_search (j, lon) = coords1 (ijk(2,j)) * dble_deg2rad
00364 end do
00365
00366
00367 do j = jbeg, jend
00368 sin_search (j, lat) = coords2 (ijk(2,j)) * dble_deg2rad
00369 end do
00370
00371 cos_search = cos (sin_search)
00372 sin_search = sin (sin_search)
00373
00374
00375 do j = jbeg, jend
00376 z_search (j) = coords3 (ijk(3,j))
00377 end do
00378
00379
00380
00381
00382
00383 ijk0 = Grids(grid_id)%ijk0
00384
00385
00386
00387
00388
00389
00390 width (1:search_mode) = 3
00391
00392 if (search_mode == 2) then
00393 ijk (ndim_3d, jbeg:jend) = max (srclocz (ijk(3,jbeg:jend)), &
00394 grid_valid_shape(1,ndim_3d))
00395
00396
00397 width (ndim_3d) = 0
00398
00399 else
00400
00401 ijk (ndim_3d, jbeg:jend) = ijk0(ndim_3d) + &
00402 ((max(srclocz(ijk(3,jbeg:jend)), grid_valid_shape(1,ndim_3d)) &
00403 - ijk0(ndim_3d)) / 4) * 4
00404
00405 endif
00406
00407 do i = 1, 2
00408 ijk (i, jbeg:jend) = ijk0(i) + &
00409 ((max(srcloc(i, ijk(2, jbeg:jend)), grid_valid_shape(1,i)) &
00410 - ijk0(i)) / 4) * 4
00411 end do
00412
00413
00414
00415
00416
00417
00418
00419 call psmile_neigh_nearx_sub_irr_dble (grid_id, &
00420 x_coords, y_coords, z_coords, coords_shape, &
00421 mask_array, mask_shape, mask_available, &
00422 sin_values, cos_values, grid_valid_shape, search_mode, &
00423 neighbors_3d, num_neigh, nloc, extra_search, &
00424 ijk, sin_search, cos_search, z_search, jbeg, jend, &
00425 width, ierror)
00426
00427
00428
00429 Deallocate (ijk)
00430
00431 Deallocate (sin_search, cos_search, z_search)
00432
00433
00434
00435 #ifdef VERBOSE
00436 print 9980, trim(ch_id), ierror
00437
00438 call psmile_flushstd
00439 #endif /* VERBOSE */
00440
00441
00442
00443 9990 format (1x, a, ': psmile_neigh_nearx_irreg2_dble')
00444 9980 format (1x, a, ': psmile_neigh_nearx_irreg2_dble: eof, ierror =', i3)
00445
00446 end subroutine PSMILe_Neigh_nearx_irreg2_dble