00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 subroutine psmile_trili_gauss2_extra (search, grid_id, mask_array, &
00013 mask_shape, mask_available, ibuf, len_item, &
00014 n_send, num_neigh, ierror)
00015
00016
00017
00018 use PRISM
00019 use PSMILe, dummy_interface => psmile_trili_gauss2_extra
00020 use psmile_grid_reduced_gauss
00021
00022 Implicit none
00023
00024
00025
00026 Integer, Intent (In) :: grid_id
00027
00028 Integer, Intent (In) :: mask_shape (2, ndim_3d)
00029
00030
00031
00032 Logical, Intent (In) ::
00033 mask_array ( mask_shape (1,1):mask_shape (2,1),
00034 mask_shape (1,2):mask_shape (2,2),
00035 mask_shape (1,3):mask_shape (2,3) )
00036
00037
00038 Logical, Intent (In) :: mask_available
00039
00040
00041
00042
00043
00044 Integer, Intent (In) :: len_item
00045
00046
00047
00048
00049
00050
00051 Integer, Intent (In) :: n_send
00052
00053
00054
00055 Integer, Intent (In) :: num_neigh
00056
00057
00058
00059
00060
00061 Type (Enddef_global_search), Intent (InOut) :: search
00062
00063
00064
00065 Integer, Intent (InOut) :: ibuf (len_item, n_send)
00066
00067
00068
00069
00070
00071
00072 Integer, Intent (Out) :: ierror
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084 Integer, Parameter :: n_corners_3d = 2 * 2 * 2
00085
00086 Integer, Parameter :: masked_out = 0
00087
00088
00089
00090
00091
00092 Integer :: j, n, i, k
00093
00094
00095
00096 Type (Grid), Pointer :: grid_info
00097 Integer, Pointer :: grid_valid_shape (:, :)
00098
00099
00100
00101 Logical :: virtual_cell_available
00102
00103
00104
00105 Integer :: code, n_corners, outside
00106
00107 Integer, Allocatable :: locs (:, :)
00108 Integer, Allocatable :: neighbours_3d_extra (:, :, :)
00109 Integer, Allocatable :: hash (:)
00110
00111
00112
00113 Logical :: missing_points(num_neigh, n_send)
00114 Integer :: decode_mask(num_neigh)
00115 Integer :: num_missing_points
00116
00117 integer, parameter :: shift_3d(ndim_3d, 9) =
00118 reshape ((/-1,-1,0, 0,-1,0, +1,-1,0,
00119 -1, 0,0, 0, 0,0, +1, 0,0,
00120 -1,+1,0, 0,+1,0, +1,+1,0/),
00121 (/ndim_3d,9/))
00122
00123
00124
00125 Integer, parameter :: nerrp = 1
00126 Integer :: ierrp (nerrp)
00127
00128 #ifdef DEBUG_TRACE
00129 Integer :: ictl_req, ind_loc (n_corners_3d)
00130 #endif
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151 Character(len=len_cvs_string), save :: mycvs =
00152 '$Id: psmile_trili_gauss2_extra.F90 3032 2011-03-18 17:29:36Z hanke $'
00153
00154
00155
00156
00157
00158 #ifdef VERBOSE
00159 print 9990, trim(ch_id), grids(grid_id)%comp_id
00160
00161 call psmile_flushstd
00162 #endif /* VERBOSE */
00163
00164 ierror = 0
00165 n_corners = num_neigh
00166
00167
00168
00169
00170
00171
00172
00173 grid_info => Grids(grid_id)
00174
00175 grid_valid_shape => grid_info%grid_shape
00176
00177 virtual_cell_available = len_item > ndim_3d + 2
00178
00179 if (grid_info%grid_type /= PRISM_Gaussreduced_Regvrt) then
00180 print *, "grid_type", grid_info%grid_type
00181 call psmile_assert (__FILE__, __LINE__, &
00182 "This routine is not designed for this grid type")
00183 endif
00184
00185 outside = grid_valid_shape (1, 1) - 11
00186
00187 #ifdef PRISM_ASSERTION
00188
00189 if (num_neigh < 1 .or. num_neigh > n_corners_3d) then
00190 print *, 'num_neigh', num_neigh, n_corners_3d
00191 call psmile_assert ( __FILE__, __LINE__, &
00192 'Invalid value for num_neigh (number of interpolation bases')
00193 endif
00194
00195
00196
00197
00198 do n = 1, n_send
00199 if (ibuf (2, n) < 1 .or. &
00200 ibuf (2, n) > size(grid_info%reduced_gauss_data%nbr_points_per_lat)) exit
00201 end do
00202
00203 if (n <= n_send) then
00204 print *, "i,j,k, gnbr_lats", ibuf(1:3,n), &
00205 size(grid_info%reduced_gauss_data%nbr_points_per_lat)
00206 call psmile_assert (__FILE__, __LINE__, &
00207 "second (virtual) index out of range")
00208 endif
00209
00210
00211 do n = 1, n_send
00212 if (ibuf (1, n) < 1 .or. &
00213 ibuf (1, n) > grid_info%reduced_gauss_data%nbr_points_per_lat(ibuf(2,n))) exit
00214 end do
00215
00216 if (n <= n_send) then
00217 print *, "i,j,k, points_per_lat", &
00218 ibuf(1:3,n), grid_info%reduced_gauss_data%nbr_points_per_lat(ibuf(2,n))
00219 call psmile_assert (__FILE__, __LINE__, &
00220 "first (virtual) index out of range")
00221 endif
00222 #endif
00223
00224 #ifdef DEBUG_TRACE
00225
00226 do ictl_req = 1, n_send
00227 if (ibuf (4, ictl_req) == 1) exit
00228 end do
00229 #endif
00230
00231
00232
00233
00234 allocate (neighbours_3d_extra (ndim_3d, n_send, 4), stat=ierror)
00235 if (ierror /= 0) then
00236 call psmile_error ( PRISM_Error_Alloc, 'locs', (/ndim_3d * n_send * 4/), 1, &
00237 __FILE__, __LINE__ )
00238 return
00239 endif
00240
00241 if (virtual_cell_available) then
00242
00243 do n = 1, n_send
00244
00245 if (ibuf(6, n) == 0) then
00246
00247 neighbours_3d_extra(:,n,1:4) = &
00248 psmile_gauss_get_bili_stencil (grid_id, ibuf(1:3,n))
00249
00250 else
00251
00252 neighbours_3d_extra(:,n,1:4) = &
00253 psmile_gauss_shift_bili_stencil (grid_id, ibuf(1:3,n), shift_3d(:,ibuf(6, n)))
00254
00255 endif
00256 enddo
00257 else
00258 neighbours_3d_extra(:,n,1:4) = psmile_gauss_get_bili_stencil (grid_id, ibuf(1:3,n))
00259 endif
00260
00261
00262
00263 do i = 1, num_neigh
00264 decode_mask(i) = ishft (1, i - 1)
00265 enddo
00266
00267 do n = 1, n_send
00268 missing_points(:, n) = iand (decode_mask(:), ibuf(5, n)) /= 0
00269 enddo
00270
00271
00272
00273 num_missing_points = count (missing_points)
00274 #ifdef PRISM_ASSERTION
00275 if (num_missing_points <= 0) then
00276 call psmile_assert ( __FILE__, __LINE__, &
00277 'Number of points to be searched == 0 ?!?')
00278 endif
00279 #endif
00280
00281 #ifdef DEBUGX
00282 print *, 'code for locations to be searched (n_send =', num_missing_points, ')'
00283 do n = 1, n_send
00284 print *, n, 'srcloc =', ibuf(1:3,n), 'code =', ibuf (5, n)
00285 print *, ' decoded =', missing_points(:,n)
00286 end do
00287 #endif
00288
00289 #ifdef DEBUG_TRACE
00290 if (ictl_req <= n_send) then
00291 print *, 'code for indices_req', ibuf (4, ictl_req)
00292 print *, 'ictl_req', ictl_req, 'srcloc', ibuf(1:3,ictl_req), &
00293 'code', ibuf (5, ictl_req), missing_points(:, ictl_req)
00294 endif
00295 #endif
00296
00297
00298
00299 Allocate (locs (ndim_3d, num_missing_points), stat = ierror)
00300 if (ierror /= 0) then
00301 call psmile_error ( PRISM_Error_Alloc, 'locs', (/ndim_3d * num_missing_points/), 1, &
00302 __FILE__, __LINE__ )
00303 return
00304 endif
00305
00306 k = 0
00307
00308 do i = 1, n_corners
00309
00310 do n = 1, n_send
00311 if (missing_points (i,n)) then
00312
00313 k = k + 1
00314
00315 locs (:, k) = neighbours_3d_extra(:, n, iand (i-1, 3)+1)
00316
00317 if (grid_valid_shape(2,3) - grid_valid_shape(1,3) >= ishft(n_corners, -2) - 1 .and. &
00318 n_corners > 4) then
00319 locs (3, k) = locs (3, k) + ishft (i-1, -2)
00320 endif
00321 end if
00322 end do
00323 end do
00324
00325 #ifdef PRISM_ASSERTION
00326 if (k /= num_missing_points) then
00327 print *, 'k, num_missing_points', k, num_missing_points
00328 call psmile_assert ( __FILE__, __LINE__, &
00329 'Inconsistent values of k and num_missing_points; must be same value')
00330 endif
00331 #endif
00332
00333 #ifdef DEBUG_TRACE
00334 if (ictl_req <= n_send) then
00335 ind_loc (:) = 0
00336 n = 0
00337 do j = 1, n_corners
00338
00339 do i = 1, n_send
00340 if (missing_points (j,i)) then
00341 n = n + 1
00342 if (i == ictl_req) ind_loc (j) = n
00343 end if
00344 end do
00345 end do
00346
00347 print *, 'locs searched for indices_req', ibuf (4, ictl_req)
00348 do j = 1, n_corners
00349 n = ind_loc (j)
00350 if (n > 0) &
00351 print *, 'ictl_req', ictl_req, ', n ', n, &
00352 ': locs', locs (1:ndim_3d, n)
00353 end do
00354 endif
00355 #endif
00356
00357
00358
00359
00360
00361 locs(1,:) = psmile_gauss_3d_to_local_1d (grid_id, locs(:,:), &
00362 num_missing_points, outside)
00363 locs(2,:) = 1
00364
00365
00366
00367
00368
00369 Allocate (hash (num_missing_points), stat = ierror)
00370 if (ierror /= 0) then
00371 ierrp (1) = num_missing_points
00372
00373 ierror = PRISM_Error_Alloc
00374
00375 call psmile_error ( ierror, 'hash', ierrp, 1, &
00376 __FILE__, __LINE__ )
00377 return
00378 endif
00379
00380 call psmile_hash_extra (search, locs, hash, num_missing_points, &
00381 mask_array, mask_shape, mask_available, grid_valid_shape, &
00382 ierror)
00383 if (ierror > 0) return
00384
00385 #ifdef DEBUGX
00386 print *, 'hash values and locs ', grid_valid_shape
00387 n = 0
00388 do j = 1, n_corners
00389
00390 do i = 1, n_send
00391 if (missing_points (j,i)) then
00392 n = n + 1
00393 print *, 'i, j, n, hash(n), locs(:,n) ', &
00394 i, j, n, hash(n), locs(:,n)
00395 end if
00396 end do
00397 end do
00398 #endif
00399
00400
00401
00402
00403 if (search%n_liste > 0) then
00404 ibuf (5, 1:n_send) = 0
00405
00406 n = 0
00407 code = 1
00408 do j = 1, n_corners
00409
00410 do i = 1, n_send
00411 if (missing_points (j,i)) then
00412
00413 n = n + 1
00414
00415 if (hash (n) >= 0) then
00416
00417
00418
00419
00420
00421 ibuf (5, i) = ibuf (5, i) + code
00422 endif
00423
00424 endif
00425 end do
00426
00427 code = code * 2
00428 end do
00429
00430 #ifdef PRISM_ASSERTION
00431 if (n /= num_missing_points) then
00432 print *, 'num_missing_points, n', num_missing_points, n
00433 call psmile_assert ( __FILE__, __LINE__, &
00434 'n /= num_missing_points')
00435 endif
00436 #endif
00437
00438 #ifdef DEBUGX
00439 print *, 'code for locations found', n_send
00440 do i = 1, n_send
00441 print *, 'code', i, ibuf (5, i), missing_points(:, i)
00442 end do
00443 #endif
00444
00445 #ifdef DEBUG_TRACE
00446 if (ictl_req <= n_send) then
00447 print *, 'code for locations found for indices_req', ibuf (4, ictl_req)
00448 print *, 'ictl_req', ictl_req, ': code', ibuf (5, ictl_req), missing_points(:, ictl_req)
00449 endif
00450 #endif
00451
00452 endif
00453
00454
00455
00456 Deallocate (locs, hash, neighbours_3d_extra)
00457
00458
00459
00460 #ifdef VERBOSE
00461 print 9980, trim(ch_id), ierror, search%n_found, n_send
00462
00463 call psmile_flushstd
00464 #endif /* VERBOSE */
00465
00466 return
00467
00468
00469
00470 #ifdef VERBOSE
00471 9990 format (1x, a, ': psmile_trili_gauss2_extra: comp_id =', i3)
00472 9980 format (1x, a, ': psmile_trili_gauss2_extra: eof', &
00473 ', ierror =', i4, ', n_found =', i8, i8)
00474 #endif /* VERBOSE */
00475
00476 end subroutine psmile_trili_gauss2_extra