00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 subroutine psmile_search_donor_extra_nn (comp_info, search, &
00013 var_id, tol, ierror)
00014
00015
00016
00017 use PRISM
00018
00019 use PSMILe, dummy_interface => PSMILe_Search_donor_extra_nn
00020
00021 implicit none
00022
00023
00024
00025 Type (Enddef_comp), Intent (In) :: comp_info
00026
00027
00028
00029
00030 Integer, Intent (In) :: var_id
00031
00032
00033
00034 Double Precision, Intent (In) :: tol
00035
00036
00037
00038
00039
00040
00041
00042 Type (Enddef_global_search), Intent (InOut) :: search
00043
00044
00045
00046
00047
00048 Integer, Intent (Out) :: ierror
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058 Type (Gridfunction), Pointer :: field
00059
00060 Integer :: method_id
00061
00062
00063
00064 Type (Grid), Pointer :: grid_info
00065 Integer :: datatype, grid_id
00066
00067
00068
00069
00070
00071 Integer :: len_item, n_send
00072 Integer :: ip_dist, nb_extra
00073
00074 Integer :: ip (ndim_3d)
00075 Real :: rtol
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085 Integer :: mask_id
00086 Logical :: src_mask_available
00087
00088 Logical, Target :: dummy_mask_array (1)
00089
00090 Integer :: dummy_mask_shape (2, ndim_3d)
00091
00092 Logical, Pointer :: mask_array (:)
00093
00094 Integer :: mask_shape (2, ndim_3d)
00095
00096
00097
00098 Integer :: i, n, nlocs
00099
00100 Integer, Allocatable :: nfound (:)
00101 Integer, Pointer :: locations (:, :), locs (:, :)
00102 Integer, Allocatable :: hash (:)
00103
00104
00105
00106 Integer, Parameter :: nerrp = 2
00107 Integer :: ierrp (nerrp)
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127 Character(len=len_cvs_string), save :: mycvs =
00128 '$Id: psmile_search_donor_extra_nn.F90 3112 2011-04-07 15:03:18Z hanke $'
00129
00130
00131
00132
00133 #ifdef VERBOSE
00134 print 9990, trim(ch_id), comp_info%comp_id, search%sender
00135
00136 call psmile_flushstd
00137 #endif /* VERBOSE */
00138
00139 ierror = 0
00140
00141 field => Fields (var_id)
00142
00143 method_id = field%method_id
00144 mask_id = field%mask_id
00145
00146 grid_id = Methods(method_id)%grid_id
00147 grid_info => Grids (grid_id)
00148
00149 src_mask_available = mask_id /= PRISM_UNDEFINED
00150
00151 if (src_mask_available) then
00152 mask_array => Masks(mask_id)%mask_array
00153
00154 mask_shape = Masks(mask_id)%mask_shape
00155 else
00156 mask_array => dummy_mask_array
00157
00158 mask_shape = dummy_mask_shape
00159 endif
00160
00161
00162
00163
00164
00165 n_send = search%msg_extra%num_volumes
00166 len_item = search%msg_extra%num_int_per_vol
00167 nb_extra = search%msg_extra%num_neigh
00168
00169 #ifdef PRISM_ASSERTION
00170
00171
00172
00173 if (search%msg_extra%reqest_type /= PSMILe_nnghbr3D) then
00174 call psmile_assert ( __FILE__, __LINE__, &
00175 'Interpolation method must be Nearest Neighbour')
00176 endif
00177
00178 if (search%msg_extra%num_items_per_coord < ndim_3d + 1 ) then
00179 call psmile_assert ( __FILE__, __LINE__, &
00180 'expected at least len_rtem = 4')
00181 endif
00182
00183 #endif
00184
00185 #ifdef TODO_GLOBAL_COORDS
00186
00187
00188
00189
00190 shift (1:ndim_3d) = grid_info%grid_shape(1, 1:ndim_3d) - &
00191 grid_info%partition (1, 1:ndim_3d) - 1
00192
00193
00194 do i = 1, n_send
00195 search%ibuf ((i-1)*len_item+1) = &
00196 search%ibuf ((i-1)*len_item+1) + shift (1)
00197 search%ibuf ((i-1)*len_item+2) = &
00198 search%ibuf ((i-1)*len_item+2) + shift (2)
00199 search%ibuf ((i-1)*len_item+3) = &
00200 search%ibuf ((i-1)*len_item+3) + shift (3)
00201 end do
00202
00203 #ifdef DEBUGX
00204 print *, 'grid shape', grid_info%grid_shape
00205 print *, 'srcloc (1:3), indices (j), len_item', len_item
00206 do i = 1, n_send
00207 print 8900, search%ibuf ((i-1)*len_item+1:i*len_item)
00208 end do
00209 #endif
00210 #endif
00211
00212
00213
00214 Allocate (nfound(n_send), locations (ndim_3d, n_send), STAT = ierror)
00215
00216
00217
00218 if ( ierror > 0 ) then
00219 ierrp (1) = ierror
00220 ierrp (2) = (ndim_3d + 1) * n_send
00221
00222 ierror = PRISM_Error_Alloc
00223 call psmile_error ( ierror, 'locations, nfound', &
00224 ierrp, 2, __FILE__, __LINE__ )
00225 return
00226 endif
00227
00228 nfound = 0
00229
00230
00231
00232
00233
00234
00235
00236 ip(1) = 1
00237 ip(2) = ip(1) + n_send
00238 ip(3) = ip(2) + n_send
00239 ip_dist = ip(3) + n_send
00240
00241 datatype = grid_info%corner_pointer%corner_datatype
00242
00243 if (datatype == MPI_REAL) then
00244
00245 rtol = tol
00246
00247 call psmile_search_donor_nnx_real (comp_info, &
00248 search, var_id, &
00249 search%rbuf(ip(1):ip(1)+n_send-1), &
00250 search%rbuf(ip(2):ip(2)+n_send-1), &
00251 search%rbuf(ip(3):ip(3)+n_send-1), &
00252 search%rbuf(ip_dist:ip_dist+n_send-1), &
00253 nfound, locations, n_send, nb_extra, rtol, ierror)
00254
00255 else if (datatype == MPI_DOUBLE_PRECISION) then
00256
00257 call psmile_search_donor_nnx_dble (comp_info, &
00258 search, var_id, &
00259 search%dbuf(ip(1):ip(1)+n_send-1), &
00260 search%dbuf(ip(2):ip(2)+n_send-1), &
00261 search%dbuf(ip(3):ip(3)+n_send-1), &
00262 search%dbuf(ip_dist:ip_dist+n_send-1), &
00263 nfound, locations, n_send, nb_extra, tol, ierror)
00264
00265 #if defined ( PRISM_QUAD_TYPE )
00266 else if (datatype == MPI_REAL16) the
00267 #endif
00268 endif
00269
00270 if (ierror > 0) return
00271
00272
00273
00274
00275
00276 nlocs = SUM (nfound)
00277
00278 #ifdef DEBUGX
00279 print *, 'psmile_search_donor_extra_nn: nlocs, n_send', nlocs, n_send
00280 do i = 1, n_send
00281 if (nfound(i) < 0 .or. nfound (i) > 1) then
00282 print *, 'i, nfound', i, nfound(i)
00283 endif
00284 end do
00285 #endif
00286
00287 if (nlocs > 0) then
00288 Allocate (hash (nlocs), stat = ierror)
00289 if (ierror /= 0) then
00290 ierrp (1) = nlocs
00291
00292 ierror = PRISM_Error_Alloc
00293
00294 call psmile_error ( ierror, 'hash', ierrp, 1, &
00295 __FILE__, __LINE__ )
00296 return
00297 endif
00298
00299 if (nlocs < n_send) then
00300 Allocate (locs (ndim_3d, nlocs), stat = ierror)
00301 if (ierror /= 0) then
00302 ierrp (1) = ndim_3d * nlocs
00303
00304 ierror = PRISM_Error_Alloc
00305
00306 call psmile_error ( ierror, 'locs', ierrp, 1, &
00307 __FILE__, __LINE__ )
00308 return
00309 endif
00310
00311 n = 0
00312
00313 do i = 1, n_send
00314 if (nfound (i) > 0) then
00315 n = n + 1
00316 locs (1, n) = locations (1, i)
00317 locs (2, n) = locations (2, i)
00318 locs (3, n) = locations (3, i)
00319 endif
00320 end do
00321
00322 #ifdef PRISM_ASSERTION
00323 if (n /= nlocs) then
00324 print *, 'nlocs, n', nlocs, n, n_send
00325 call psmile_assert ( __FILE__, __LINE__, &
00326 'nlocs != n')
00327 endif
00328 #endif
00329
00330 Deallocate (locations)
00331 else
00332
00333 locs => locations
00334
00335 endif
00336
00337
00338
00339
00340 call psmile_hash_extra (search, locs, hash, nlocs, &
00341 mask_array, mask_shape, src_mask_available, &
00342 grid_info%grid_shape, ierror)
00343 if (ierror > 0) return
00344
00345 Deallocate (locs, hash)
00346
00347 else
00348
00349
00350
00351 search%n_liste = 0
00352 search%n_found = 0
00353
00354 NULLIFY (search%neigh_3d)
00355 NULLIFY (search%index_found)
00356
00357 Deallocate (locations)
00358
00359 endif
00360
00361
00362
00363
00364
00365 if (datatype == MPI_REAL) then
00366 call psmile_return_extra_off_real (comp_info, search, var_id, &
00367 nfound, search%rbuf(ip_dist:ip_dist+n_send-1), n_send, &
00368 nb_extra, ierror)
00369
00370 else if (datatype == MPI_DOUBLE_PRECISION) then
00371 call psmile_return_extra_off_dble (comp_info, search, var_id, &
00372 nfound, search%dbuf(ip_dist:ip_dist+n_send-1), n_send, &
00373 nb_extra, ierror)
00374
00375 #if defined ( PRISM_QUAD_TYPE )
00376 else if (datatype == MPI_REAL16) the
00377 call psmile_return_extra_off_quad (comp_info, search, var_id, &
00378 nfound, search%qbuf(ip_dist:ip_dist+n_send-1), n_send, &
00379 nb_extra, ierror)
00380 #endif
00381 endif
00382
00383 Deallocate (nfound)
00384
00385
00386
00387 #ifdef VERBOSE
00388 print 9980, trim(ch_id), grid_id, search%sender, ierror
00389
00390 call psmile_flushstd
00391 #endif /* VERBOSE */
00392
00393 return
00394
00395
00396
00397 #ifdef VERBOSE
00398 9990 format (1x, a, ': psmile_search_donor_extra_nn: comp_id =', i3, &
00399 '; sender =', i4)
00400 9980 format (1x, a, ': psmile_search_donor_extra_nn: grid id =', i3, &
00401 '; eof sender =', i3, ', ierror =', i4)
00402
00403 9970 format (1x, a, ': psmile_search_donor_extra_nn:', &
00404 1x, 'offset local ', 3i7, &
00405 /1x, 'extent ', 3i7, &
00406 /1x, 'offset remote ', 3i7)
00407 #endif /* VERBOSE */
00408 8900 format (1x, 3i6, 5(';', i6) )
00409
00410 end subroutine PSMILe_Search_donor_extra_nn