00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_enddef_action (search, index, status, ierror)
00012
00013
00014
00015 use PRISM_constants
00016 use PSMILe_common, only : num_req_types
00017 use PSMILe, dummy_interface => PSMILe_Enddef_action
00018
00019 implicit none
00020
00021
00022
00023
00024
00025 Integer, Intent (In) :: index
00026
00027
00028
00029
00030
00031 Integer, Intent (In) :: status (MPI_STATUS_SIZE)
00032
00033
00034
00035
00036
00037 Type (Enddef_search), Intent (InOut) :: search
00038
00039
00040
00041
00042
00043 Integer, Intent (Out) :: ierror
00044
00045
00046
00047
00048
00049
00050
00051
00052 Logical, Parameter :: new_search = .true.
00053
00054
00055
00056
00057
00058 logical :: new_intersection
00059 Type (enddef_msg_intersections) :: msg_intersections
00060
00061
00062
00063 Integer :: npart
00064
00065
00066
00067 Integer :: ind
00068
00069
00070
00071 Integer :: meslen, sender
00072 Integer :: lstatus (MPI_STATUS_SIZE)
00073
00074
00075
00076 Type (enddef_msg_locations) :: msg_locations
00077
00078
00079
00080 Integer, parameter :: nerrp = 3
00081 Integer :: ierrp (nerrp)
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100 Character(len=len_cvs_string), save :: mycvs =
00101 '$Id: psmile_enddef_action.F90 2887 2011-01-14 09:25:49Z redler $'
00102
00103
00104
00105
00106
00107 ierror = 0
00108
00109 new_intersection = .false.
00110 sender = status (MPI_SOURCE)
00111
00112 #ifdef VERBOSE
00113 print 9990, trim(ch_id), index, sender
00114
00115 call psmile_flushstd
00116 #endif /* VERBOSE */
00117
00118
00119
00120 if (index <= 2) then
00121 call MPI_Get_count (status, MPI_INTEGER, meslen, ierror)
00122 if (meslen > nd_msgint .or. meslen < ind_msgint_tag) then
00123 ierror = PRISM_Error_Wrong
00124 ierrp (1) = meslen
00125 ierrp (2) = status (MPI_SOURCE)
00126 ierrp (3) = status (MPI_TAG)
00127
00128 call psmile_error ( ierror, 'MPI_Get_count', &
00129 ierrp, 3, __FILE__, __LINE__ )
00130 return
00131 endif
00132 endif
00133
00134 select case (index)
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146 case ( 1 )
00147
00148 call psmile_init_enddef_msg_inters (msg_intersections)
00149 call psmile_unpack_msg_intersections (msg_intersections, &
00150 paction%msgreq)
00151
00152 ind = msg_intersections%relative_msg_tag
00153
00154 #ifdef PRISM_ASSERTION
00155 if (status(MPI_TAG) /= reqtag) then
00156 call psmile_assert ( __FILE__, __LINE__, &
00157 'wrong message tag')
00158 endif
00159
00160
00161
00162 if (ind < 1 .or. ind > paction%n_answer2recv) then
00163 call psmile_assert ( __FILE__, __LINE__, &
00164 'Illegal local tag')
00165 endif
00166 #endif
00167
00168 call psmile_send_req_subgrid (msg_intersections, &
00169 sender, grdtag, ierror)
00170 if (ierror > 0) return
00171
00172
00173
00174
00175
00176
00177 if (paction%lrequest (num_req_types-1+ind) /= MPI_REQUEST_NULL) then
00178
00179 call MPI_Wait (paction%lrequest (num_req_types-1+ind), lstatus, ierror)
00180
00181 if ( ierror /= MPI_SUCCESS ) then
00182 ierrp (1) = ierror
00183 ierrp (2) = lstatus (MPI_SOURCE)
00184 ierrp (3) = lstatus (MPI_TAG)
00185
00186 ierror = PRISM_Error_MPI
00187
00188 call psmile_error ( ierror, 'MPI_Wait', &
00189 ierrp, 3, __FILE__, __LINE__ )
00190 return
00191 endif
00192
00193 #ifdef PRISM_ASSERTION
00194 if (lstatus(MPI_TAG) /= loctag+ind) then
00195 call psmile_assert ( __FILE__, __LINE__, &
00196 'wrong message tag')
00197 endif
00198
00199 if (lstatus(MPI_SOURCE) /= sender) then
00200 call psmile_assert ( __FILE__, __LINE__, &
00201 'wrong sender')
00202 endif
00203 #endif
00204
00205 call psmile_unpack_msg_locations(msg_locations, &
00206 paction%loc_messages (:, ind))
00207
00208 if ( msg_locations%requires_conserv_remap == PSMILe_conserv2D ) then
00209
00210
00211
00212 #ifdef PRISM_ASSERTION
00213 if (msg_locations%src_rank /= sender) then
00214 call psmile_assert ( __FILE__, __LINE__, &
00215 'wrong sender')
00216 endif
00217 #endif
00218 call psmile_enddef_action_cell (msg_locations, ierror)
00219 if (ierror > 0) return
00220
00221
00222
00223 #ifdef DEBUGX
00224 print *, ' Receiving locations from', sender, 'with tag ', loctag+ind, &
00225 ' size ', msgloc_size
00226 #endif /* DEBUGX */
00227 call MPI_Recv (paction%loc_messages (1,ind), msgloc_size, &
00228 MPI_INTEGER, sender, loctag+ind, comm_psmile, lstatus, ierror)
00229
00230 if ( ierror /= MPI_SUCCESS ) then
00231 ierrp (1) = ierror
00232 ierrp (2) = sender
00233 ierrp (3) = loctag+num_req_types-1+ind
00234
00235 ierror = PRISM_Error_Recv
00236
00237 call psmile_error ( ierror, 'MPI_Recv', &
00238 ierrp, 3, __FILE__, __LINE__ )
00239 return
00240 endif
00241
00242 endif
00243
00244 call psmile_unpack_msg_locations(msg_locations, &
00245 paction%loc_messages (:, ind))
00246
00247 call psmile_enddef_action_loc (msg_locations, ierror)
00248 if (ierror > 0) return
00249
00250 endif
00251
00252
00253
00254 call MPI_Irecv (paction%loc_messages(1,ind), msgloc_size, MPI_INTEGER, &
00255 sender, loctag+ind, comm_psmile, &
00256 paction%lrequest (num_req_types-1+ind), ierror)
00257
00258 if ( ierror /= MPI_SUCCESS ) then
00259 ierrp (1) = ierror
00260 ierrp (2) = sender
00261 ierrp (3) = loctag+ind
00262
00263 ierror = PRISM_Error_Recv
00264
00265 call psmile_error ( ierror, 'MPI_Irecv', &
00266 ierrp, 3, __FILE__, __LINE__ )
00267 return
00268 endif
00269 #ifdef DEBUGX
00270 print *, ' Posting Irecv request(',num_req_types-1+ind,') ', &
00271 paction%lrequest(num_req_types-1+ind), 'with tag ', loctag+ind, &
00272 ' and size ', msgloc_size
00273 #endif /* DEBUGX */
00274
00275
00276
00277
00278
00279 paction%n_answer = paction%n_answer + 1
00280
00281 if (paction%n_answer < paction%n_answer2recv) then
00282 #ifdef PRISM_ASSERTION
00283 if (paction%lrequest(1) /= MPI_REQUEST_NULL) then
00284 call psmile_assert ( __FILE__, __LINE__, &
00285 'paction%lrequest(1) is not finished !')
00286 endif
00287 #endif
00288
00289 call MPI_Irecv (paction%msgreq, nd_msgint, MPI_INTEGER, &
00290 MPI_ANY_SOURCE, reqtag, comm_psmile, &
00291 paction%lrequest(1), ierror)
00292
00293 if ( ierror /= MPI_SUCCESS ) then
00294 ierrp (1) = ierror
00295 ierror = PRISM_Error_MPI
00296
00297 call psmile_error ( ierror, 'MPI_Irecv', &
00298 ierrp, 1, __FILE__, __LINE__ )
00299 return
00300 endif
00301 #ifdef DEBUGX
00302 print *, ' Posting Irecv request(1) ', paction%lrequest(1), 'with tag ', reqtag
00303 #endif /* DEBUGX */
00304 #ifdef PRISM_ASSERTION
00305 else if (paction%lrequest (1) /= MPI_REQUEST_NULL) then
00306 call psmile_assert ( __FILE__, __LINE__, &
00307 'bad request; should be MPI_REQUEST_NULL')
00308 #endif
00309 endif
00310
00311 case (2)
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321 call psmile_init_enddef_msg_inters (msg_intersections)
00322 call psmile_unpack_msg_intersections (msg_intersections, &
00323 paction%msgint)
00324
00325 #ifdef PRISM_ASSERTION
00326 if (status(MPI_TAG) /= paction%lastag) then
00327 call psmile_assert ( __FILE__, __LINE__, &
00328 'wrong message tag; must be lastag')
00329 endif
00330 #endif
00331 paction%n = paction%n + 1
00332
00333 npart = msg_intersections%num_parts
00334
00335
00336
00337 if (npart <= 0) then
00338 new_intersection = .true.
00339 else
00340
00341
00342
00343
00344
00345
00346 call psmile_bsend (paction%msgint, nd_msgint, MPI_INTEGER, &
00347 sender, reqtag, comm_psmile, ierror)
00348
00349 if (ierror /= MPI_SUCCESS) then
00350 call psmile_error (PRISM_Error_Send, 'MPI_Send', &
00351 (/ierror, sender, reqtag/), 3, &
00352 __FILE__, __LINE__ )
00353 ierror = PRISM_Error_Send
00354 return
00355 endif
00356
00357
00358
00359
00360 Allocate (paction%recv_req ((ndim_3d+2)*npart), stat = ierror)
00361 if ( ierror /= 0 ) then
00362 ierrp (1) = (ndim_3d+1) * npart
00363 ierror = PRISM_Error_Alloc
00364 call psmile_error ( ierror, 'paction%recv_req', &
00365 ierrp, 1, __FILE__, __LINE__ )
00366 return
00367 endif
00368
00369 paction%recv_req(:) = MPI_REQUEST_NULL
00370
00371 call psmile_recv_req_subgrid (msg_intersections, sender, grdtag, search, &
00372 paction%recv_req, new_search, ierror)
00373 if (ierror > 0) return
00374
00375 paction%lrequest (3) = paction%recv_req (1)
00376 paction%grid2receive = .true.
00377 new_intersection = .false.
00378 endif
00379
00380 case ( 3 )
00381
00382
00383
00384
00385
00386
00387
00388
00389 ierror = PRISM_Error_Internal
00390
00391 call psmile_error ( ierror, 'Receive of grid data is not supported in PSMILe_Enddef_action', &
00392 ierrp, 0, __FILE__, __LINE__ )
00393
00394 case ( 4 )
00395
00396
00397
00398
00399
00400
00401 if (paction%msg_extra(1) == PSMILe_Finalize_extra_search) then
00402 paction%n_fin = paction%n_fin + paction%msg_extra(2)
00403
00404 #ifdef PRISM_ASSERTION
00405 if (paction%n_fin > paction%n_fin2recv) then
00406 print *, 'n_fin, n_fin2recv', paction%n_fin, paction%n_fin2recv
00407 call psmile_assert (__FILE__, __LINE__, &
00408 "Too many Finalize messages for extra search received")
00409 endif
00410 #endif
00411 else
00412 call psmile_enddef_action_extra (paction%msg_extra, nd_msgextra, &
00413 sender, ierror)
00414 if (ierror > 0) return
00415 endif
00416
00417
00418
00419 if (paction%n_fin < paction%n_fin2recv) then
00420 call MPI_Irecv (paction%msg_extra, nd_msgextra, MPI_INTEGER, &
00421 MPI_ANY_SOURCE, exttag, comm_psmile, &
00422 paction%lrequest(4), ierror)
00423
00424 if ( ierror /= MPI_SUCCESS ) then
00425 ierrp (1) = ierror
00426 ierror = PRISM_Error_MPI
00427
00428 call psmile_error ( ierror, 'MPI_Irecv', &
00429 ierrp, 1, __FILE__, __LINE__ )
00430 return
00431 endif
00432 #ifdef DEBUGX
00433 print *, ' Posting Irecv request(4) ', paction%lrequest(4), 'with tag ', exttag
00434 #endif /* DEBUGX */
00435
00436 endif
00437
00438 case ( 5 )
00439
00440
00441
00442
00443
00444
00445 call psmile_enddef_action_sel (status(MPI_SOURCE), ierror)
00446 if (ierror > 0) return
00447
00448 case (num_req_types:)
00449
00450
00451
00452
00453
00454
00455 call psmile_unpack_msg_locations(msg_locations, &
00456 paction%loc_messages (:, index-num_req_types+1))
00457
00458
00459 if ( msg_locations%requires_conserv_remap == PSMILe_conserv2D ) then
00460
00461
00462
00463 #ifdef PRISM_ASSERTION
00464 if (msg_locations%src_rank /= sender) then
00465 call psmile_assert ( __FILE__, __LINE__, &
00466 'wrong sender')
00467 endif
00468 #endif
00469 call psmile_enddef_action_cell (msg_locations, ierror)
00470 if (ierror > 0) return
00471
00472
00473
00474 call MPI_Irecv (paction%loc_messages (1, index-num_req_types+1), msgloc_size, &
00475 MPI_INTEGER, sender, loctag+index-num_req_types+1, comm_psmile, &
00476 paction%lrequest(index), ierror)
00477
00478 if ( ierror /= MPI_SUCCESS ) then
00479 ierrp (1) = ierror
00480 ierrp (2) = sender
00481 ierrp (3) = loctag+index-num_req_types+1
00482
00483 ierror = PRISM_Error_Recv
00484
00485 call psmile_error ( ierror, 'MPI_Irecv', &
00486 ierrp, 3, __FILE__, __LINE__ )
00487 return
00488 endif
00489
00490 #ifdef DEBUGX
00491 print *, ' Posting Irecv request(', index, ')', &
00492 paction%lrequest(index), 'with tag ', loctag+index-num_req_types+1, &
00493 ' and size ', msgloc_size
00494 #endif /* DEBUGX */
00495 else
00496
00497 call psmile_enddef_action_loc (msg_locations, ierror)
00498 if (ierror > 0) return
00499
00500 endif
00501
00502 end select
00503
00504
00505
00506
00507 if (new_intersection .and. paction%n < paction%ninter) then
00508 #ifdef PRISM_ASSERTION
00509 if (paction%lrequest(2) /= MPI_REQUEST_NULL) then
00510 call psmile_assert ( __FILE__, __LINE__, &
00511 'paction%lrequest(2) is not finished !')
00512 endif
00513 #endif
00514
00515 #define CONS_REMAP_DEADLOCK_WORKAROUND
00516 #ifndef CONS_REMAP_DEADLOCK_WORKAROUND
00517 call MPI_Irecv (paction%msgint, nd_msgint, MPI_INTEGER, &
00518 MPI_ANY_SOURCE, paction%lastag, comm_psmile, &
00519 paction%lrequest(2), ierror)
00520 #else
00521 call MPI_Irecv (paction%msgint, nd_msgint, MPI_INTEGER, &
00522 maxval(paction%intersect_ranks), paction%lastag, comm_psmile, &
00523 paction%lrequest(2), ierror)
00524 paction%intersect_ranks(maxloc(paction%intersect_ranks)) = -1
00525 #endif
00526 if ( ierror /= MPI_SUCCESS ) then
00527 ierrp (1) = ierror
00528 ierror = PRISM_Error_MPI
00529
00530 call psmile_error ( ierror, 'MPI_Irecv', &
00531 ierrp, 1, __FILE__, __LINE__ )
00532 return
00533 endif
00534 #ifdef DEBUGX
00535 print *, ' Posting Irecv request(2)', paction%lrequest(2), 'with tag ', paction%lastag
00536 #endif /* DEBUGX */
00537 endif
00538
00539
00540
00541 #ifdef VERBOSE
00542 print 9980, trim(ch_id), ierror
00543 call psmile_flushstd
00544 #endif /* VERBOSE */
00545
00546
00547
00548
00549
00550 #ifdef VERBOSE
00551
00552 9990 format (1x, a, ': psmile_enddef_action: index', i4, ', sender', i4)
00553 9980 format (1x, a, ': psmile_enddef_action: eof ierror =', i3)
00554
00555 #endif /* VERBOSE */
00556
00557 end subroutine PSMILe_Enddef_action