00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_select_nn_found (search, extra_search, &
00012 send_info, &
00013 selected, sel_info, nrecv, nb_extra, &
00014 neighbors_3d, nloc, num_neigh, &
00015 grid_valid_shape, ierror)
00016
00017
00018
00019 use PRISM
00020
00021 use PSMILe, dummy_interface => PSMILe_Select_nn_found
00022
00023 Implicit none
00024
00025
00026
00027 Type (Enddef_search), Intent (In) :: search
00028
00029
00030
00031 Integer, Intent (In) :: nb_extra
00032
00033
00034
00035 Integer, Intent (In) :: nloc
00036
00037
00038
00039
00040
00041 Integer, Intent (In) :: num_neigh
00042
00043
00044
00045
00046 Integer, Intent (In) :: grid_valid_shape (2, ndim_3d)
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056 Integer, Intent (InOut) :: nrecv
00057
00058
00059
00060
00061 Type (Extra_search_info), Intent (InOut) :: extra_search
00062
00063
00064
00065
00066 Type (Send_information), Intent (InOut) :: send_info
00067
00068
00069
00070
00071 Integer, Intent (InOut) :: selected (2, extra_search%n_extra)
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081 Type (Select_search_info), Intent (InOut) :: sel_info (nrecv)
00082
00083
00084 Integer, Intent (InOut) :: neighbors_3d (ndim_3d, nloc,
00085 num_neigh)
00086
00087
00088
00089 Integer, Intent (Out) :: ierror
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099 Integer, Parameter :: masked_out = 0
00100
00101
00102
00103
00104
00105 Integer :: i, j, k, n
00106
00107
00108
00109 Integer :: ndibuf
00110 Integer :: n_liste, nprev_recv
00111 Integer, Pointer :: used (:), ibuf (:)
00112 Integer, Pointer :: indices (:)
00113
00114
00115
00116 Integer :: nold
00117 Integer, Pointer :: len_sent (:)
00118 Integer, Pointer :: sender_global (:)
00119 Integer, Pointer :: msg_id (:)
00120 Type (dble_vector), Pointer :: dble_bufs (:)
00121 Type (real_vector), Pointer :: real_bufs (:)
00122
00123
00124
00125 Integer :: nreq, sender
00126 Integer :: msg_sel (nd_msgsel), req (2)
00127 #ifndef PRISM_with_MPI2
00128 Integer :: statuses (MPI_STATUS_SIZE, 2)
00129 #endif
00130
00131
00132
00133 Integer, Parameter :: nerrp = 3
00134 Integer :: ierrp (nerrp)
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156 Character(len=len_cvs_string), save :: mycvs =
00157 '$Id: psmile_select_nn_found.F90 2799 2010-12-03 14:33:35Z hanke $'
00158
00159
00160
00161
00162
00163 #ifdef VERBOSE
00164 print 9990, trim(ch_id), nrecv
00165
00166 call psmile_flushstd
00167 #endif /* VERBOSE */
00168
00169 ierror = 0
00170 ndibuf = 0
00171 nprev_recv = send_info%num2recv
00172
00173 indices => extra_search%indices
00174
00175 sel_info(1:nrecv)%num_req = (/ (i, i=1,nrecv) /)
00176
00177 #ifdef PRISM_ASSERTION
00178
00179
00180
00181 if (nb_extra > 1) then
00182 call psmile_assert (__FILE__, __LINE__, &
00183 'nb_extra > 1 currently not supported')
00184 endif
00185 #endif
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196 msg_sel (1) = PSMILe_nnghbr3D
00197
00198
00199
00200 i = 1
00201 do while ( i <= nrecv)
00202 used => sel_info(i)%used
00203 n_liste = sel_info(i)%n_liste
00204 sender = sel_info(i)%sender
00205
00206 n = Count (used > 0)
00207 #ifdef DEBUGX
00208 print *, 'psmile_select_nn_found: i, sender, n, n_liste', &
00209 i, sender, n, n_liste
00210 #endif
00211
00212 msg_sel (2) = n_liste
00213 msg_sel (3) = n
00214 msg_sel (4) = sel_info(i)%index
00215 msg_sel (5) = sel_info(i)%method_id
00216
00217 nreq = 1
00218
00219
00220
00221 call MPI_Isend (msg_sel, nd_msgsel, MPI_INTEGER, &
00222 sender, seltag, comm_psmile, req(1), ierror)
00223
00224 if (ierror /= MPI_SUCCESS) then
00225 ierrp (1) = ierror
00226 ierrp (2) = sender
00227 ierrp (3) = seltag
00228
00229 ierror = PRISM_Error_Send
00230
00231 call psmile_error (ierror, 'MPI_Isend (msg_sel)', ierrp, 3, &
00232 __FILE__, __LINE__ )
00233 return
00234 endif
00235
00236 if (n == 0) then
00237
00238
00239
00240
00241 if (search%datatype == MPI_REAL) then
00242 Deallocate (sel_info(i)%real_buf)
00243
00244 else if (search%datatype == MPI_DOUBLE_PRECISION) then
00245 Deallocate (sel_info(i)%dble_buf)
00246
00247 #if defined ( PRISM_QUAD_TYPE )
00248 else if (search%datatype == MPI_REAL16) then
00249 Deallocate (sel_info(i)%quad_buf)
00250 #endif
00251 endif
00252
00253 Deallocate (sel_info(i)%used)
00254
00255
00256
00257 if (i /= nrecv) then
00258 sel_info(i:nrecv-1) = sel_info(i+1:nrecv)
00259 endif
00260
00261
00262
00263 nrecv = nrecv - 1
00264
00265 cycle
00266 endif
00267
00268
00269
00270
00271 if (n < n_liste) then
00272 if (n > ndibuf) then
00273 if (ndibuf > 0) then
00274 Deallocate (ibuf)
00275 endif
00276
00277 Allocate (ibuf (n), stat = ierror)
00278 if (ierror /= 0) then
00279 ierrp (1) = ierror
00280 ierrp (2) = n
00281
00282 ierror = PRISM_Error_Alloc
00283
00284 call psmile_error ( ierror, 'ibuf', ierrp, 2, &
00285 __FILE__, __LINE__ )
00286 return
00287 endif
00288
00289 ndibuf = n
00290 endif
00291
00292
00293
00294 #ifdef USE_PACK
00295 ibuf (1:n) = PACK ((/ (j, j=1,n_liste) /), used(:) > 0)
00296 #else
00297 k = 0
00298
00299 do j = 1, n_liste
00300 if (used(j) > 0) then
00301 k = k + 1
00302 ibuf (k) = j
00303 endif
00304 end do
00305
00306 #ifdef PRISM_ASSERTION
00307 if (k /= n) then
00308 print *, trim(ch_id), ": k, n", k, n
00309 call psmile_assert (__FILE__, __LINE__, &
00310 'Inconsistent number of used points')
00311 endif
00312
00313 #endif
00314 #endif /* USE_PACK */
00315
00316
00317
00318 call MPI_Isend (ibuf, n, MPI_INTEGER, &
00319 sender, seltag, comm_psmile, req(2), ierror)
00320
00321 if (ierror /= MPI_SUCCESS) then
00322 ierrp (1) = ierror
00323 ierrp (2) = sender
00324 ierrp (3) = seltag
00325
00326 ierror = PRISM_Error_Send
00327
00328 call psmile_error (ierror, 'MPI_Isend (ibuf)', ierrp, 3, &
00329 __FILE__, __LINE__ )
00330 return
00331 endif
00332
00333 nreq = 2
00334
00335
00336
00337
00338
00339 k = 0
00340
00341 do j = 1, n_liste
00342 if (used (j) > 0) then
00343 k = k + 1
00344 used (j) = k
00345 endif
00346 end do
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357 do j = 1, extra_search%n_extra
00358 if (selected(1,j) == sel_info(i)%num_req) then
00359 neighbors_3d (1, indices(j), 1) = extra_search%global_marker
00360 neighbors_3d (2, indices(j), 1) = used(selected (2,j)) + nprev_recv
00361 endif
00362 enddo
00363
00364
00365
00366
00367 sel_info(i)%n_liste = n
00368
00369 if (search%datatype == MPI_REAL) then
00370
00371 call psmile_assert (__FILE__, __LINE__, "TODO")
00372
00373 else if (search%datatype == MPI_DOUBLE_PRECISION) then
00374
00375
00376
00377 #ifdef DEBUG
00378 print *, 'psmile_select_nn_found: #### real implementation is needed !!!'
00379 #endif
00380
00381 do j = 1, n
00382 if (search%datatype == MPI_REAL) then
00383 sel_info(i)%real_buf(j) = sel_info(i)%real_buf(ibuf(j))
00384 sel_info(i)%real_buf(j+n) = sel_info(i)%real_buf(ibuf(j)+n_liste)
00385 sel_info(i)%real_buf(j+n*2) = sel_info(i)%real_buf(ibuf(j)+n_liste*2)
00386 else if (search%datatype == MPI_DOUBLE_PRECISION) then
00387 sel_info(i)%dble_buf(j) = sel_info(i)%dble_buf(ibuf(j))
00388 sel_info(i)%dble_buf(j+n) = sel_info(i)%dble_buf(ibuf(j)+n_liste)
00389 sel_info(i)%dble_buf(j+n*2) = sel_info(i)%dble_buf(ibuf(j)+n_liste*2)
00390 endif
00391 end do
00392 #if defined ( PRISM_QUAD_TYPE )
00393 else if (search%datatype == MPI_REAL16) then
00394
00395 call psmile_assert (__FILE__, __LINE__, "TODO")
00396
00397 #endif
00398 endif
00399
00400 else
00401
00402 do j = 1, extra_search%n_extra
00403 if (selected(1,j) == sel_info(i)%num_req) then
00404 neighbors_3d (1, indices(j), 1) = extra_search%global_marker
00405 neighbors_3d (2, indices(j), 1) = selected (2,j) + nprev_recv
00406 endif
00407 enddo
00408 endif
00409
00410 Deallocate (sel_info(i)%used)
00411
00412
00413
00414 nprev_recv = nprev_recv + n
00415 i = i + 1
00416
00417
00418
00419 #ifdef PRISM_with_MPI2
00420 call MPI_Waitall (nreq, req, MPI_STATUSES_IGNORE, ierror)
00421 #else
00422 call MPI_Waitall (nreq, req, statuses, ierror)
00423 #endif
00424
00425 if (ierror /= MPI_SUCCESS) then
00426 ierrp (1) = ierror
00427
00428 ierror = PRISM_Error_MPI
00429
00430 call psmile_error (ierror, 'MPI_Waitall', ierrp, 1, &
00431 __FILE__, __LINE__ )
00432 return
00433 endif
00434
00435 enddo
00436
00437
00438
00439 if (ndibuf > 0) then
00440 Deallocate (ibuf)
00441 endif
00442
00443
00444
00445 if (nrecv > 0) then
00446
00447
00448
00449 nold = send_info%nrecv
00450 n = send_info%nrecv + nrecv
00451
00452 Allocate (len_sent(n), sender_global(n), dble_bufs(n), &
00453 real_bufs(n), msg_id(n), stat = ierror)
00454 if (ierror /= 0) then
00455 ierrp (1) = ierror
00456 ierrp (2) = n * 3
00457
00458 ierror = PRISM_Error_Alloc
00459
00460 call psmile_error ( ierror, 'len_sent, ...', ierrp, 2, &
00461 __FILE__, __LINE__ )
00462 return
00463 endif
00464
00465 if (nold > 0) then
00466 len_sent(1:nold) = send_info%len_sent (1:nold)
00467 sender_global(1:nold) = send_info%sender_global (1:nold)
00468 msg_id(1:nold) = send_info%msg_id (1:nold)
00469
00470 if (search%datatype == MPI_REAL) then
00471
00472 do i = 1, nold
00473 real_bufs(i)%vector => extra_search%real_bufs(i)%vector
00474 enddo
00475
00476 Deallocate (extra_search%real_bufs)
00477
00478 else if (search%datatype == MPI_DOUBLE_PRECISION) then
00479
00480 do i = 1, nold
00481 dble_bufs(i)%vector => extra_search%dble_bufs(i)%vector
00482 enddo
00483
00484 Deallocate (extra_search%dble_bufs)
00485 #if defined ( PRISM_QUAD_TYPE )
00486 else if (search%datatype == MPI_REAL16) then
00487 #endif
00488 endif
00489
00490
00491 Deallocate (send_info%len_sent, send_info%sender_global, send_info%msg_id)
00492 endif
00493
00494
00495 do i = 1, nrecv
00496 len_sent (nold+i) = sel_info(i)%n_liste
00497 sender_global (nold+i) = sel_info(i)%sender
00498 msg_id(nold+i) = sel_info(i)%msg_id
00499 if (search%datatype == MPI_REAL) then
00500 real_bufs(nold+i)%vector => sel_info(i)%real_buf
00501 else if (search%datatype == MPI_DOUBLE_PRECISION) then
00502 dble_bufs(nold+i)%vector => sel_info(i)%dble_buf
00503 endif
00504 end do
00505
00506 send_info%len_sent => len_sent
00507 send_info%sender_global => sender_global
00508 send_info%msg_id => msg_id
00509
00510 if (search%datatype == MPI_REAL) then
00511 extra_search%real_bufs => real_bufs
00512 else if (search%datatype == MPI_DOUBLE_PRECISION) then
00513 extra_search%dble_bufs => dble_bufs
00514 endif
00515
00516 #ifdef DEBUGX
00517 print *, 'psmile_select_nn_found: nold, nrecv, num2recv, nprec_recv', &
00518 nold, nrecv, send_info%num2recv, nprev_recv
00519 print *, 'psmile_select_nn_found: len_sent', len_sent(1), len_sent(2)
00520 #endif
00521
00522 send_info%nrecv = nold + nrecv
00523 send_info%num2recv = nprev_recv
00524
00525 #ifdef PRISM_ASSERTION
00526 if (send_info%num2recv /= &
00527 SUM (send_info%len_sent(1:send_info%nrecv))) then
00528 print *, 'num2recv', send_info%num2recv, &
00529 SUM (send_info%len_sent(1:send_info%nrecv)), &
00530 send_info%nrecv
00531 call psmile_assert (__FILE__, __LINE__, &
00532 "inconsistent num2recv")
00533 endif
00534 #endif
00535 endif
00536
00537
00538
00539 #ifdef VERBOSE
00540 print 9980, trim(ch_id), ierror, nrecv
00541
00542 call psmile_flushstd
00543 #endif /* VERBOSE */
00544
00545
00546
00547 #ifdef VERBOSE
00548
00549 9990 format (1x, a, ': psmile_select_nn_found: nrecv', i7)
00550 9980 format (1x, a, ': psmile_select_nn_found: eof ierror =', i3, &
00551 '; nrecv', i7)
00552
00553 #endif /* VERBOSE */
00554
00555 #ifdef DEBUG
00556 #endif
00557
00558 end subroutine PSMILe_Select_nn_found