psmile_enddef_action.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PSMILe_Enddef_action
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_enddef_action (search, index, status, ierror)
00012 !
00013 ! !USES:
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 ! !INPUT PARAMETERS:
00022 !
00023 !     Infos on the components (located on the actual process)
00024 
00025       Integer,            Intent (In)      :: index
00026 
00027 !     Index in paction%lrequest of requests received.
00028 !     Note: Index = 3 which corresponds to the receive of grid data
00029 !           (tag = grdtag) is not allowed for this check routine.
00030 
00031       Integer,            Intent (In)      :: status (MPI_STATUS_SIZE)
00032 
00033 !     Status of the receive request
00034 !
00035 ! !INPUT/OUTPUT PARAMETERS:
00036 !
00037       Type (Enddef_search), Intent (InOut) :: search
00038 
00039 !     Data on the points to be searched
00040 !
00041 ! !OUTPUT PARAMETERS:
00042 !
00043       Integer, Intent (Out)                :: ierror
00044 
00045 !     Returns the error code of PSMILe_Enddef_action;
00046 !             ierror = 0 : No error
00047 !             ierror > 0 : Severe error
00048 !
00049 !
00050 ! !DEFINED PARAMETERS:
00051 !
00052       Logical, Parameter              :: new_search = .true.
00053 !
00054 ! !LOCAL VARIABLES
00055 !
00056 !  ... for messages contaning data on intersections and requests
00057 !
00058       logical                         :: new_intersection
00059       Type (enddef_msg_intersections) :: msg_intersections
00060 !
00061 !  ... for intersections sent
00062 !
00063       Integer                         :: npart
00064 !
00065 !  ... for grid data sent by a partner process
00066 !
00067       Integer                         :: ind
00068 !
00069 !  ... for communication
00070 !
00071       Integer                         :: meslen, sender
00072       Integer                         :: lstatus (MPI_STATUS_SIZE)
00073 !
00074 !  ... for handling extra search messages
00075 !
00076       Type (enddef_msg_extra)         :: msg_extra
00077 !
00078 !  ... for receiving locations
00079 !
00080       Type (enddef_msg_locations)     :: msg_locations
00081 !
00082 !  ... for error parameters
00083 !
00084       Integer, parameter              :: nerrp = 3
00085       Integer                         :: ierrp (nerrp)
00086 !
00087 ! !DESCRIPTION:
00088 !
00089 ! Subroutine "PSMILe_Enddef_action" performs the action required after
00090 ! a request/message from another process was received.
00091 !
00092 ! !REVISION HISTORY:
00093 !
00094 !   Date      Programmer   Description
00095 ! ----------  ----------   -----------
00096 ! 03.07.06    H. Ritzdorf  created
00097 !
00098 !EOP
00099 !----------------------------------------------------------------------
00100 !
00101 ! $Id: psmile_enddef_action.F90 3248 2011-06-23 13:03:19Z coquart $
00102 ! $Author: coquart $
00103 !
00104    Character(len=len_cvs_string), save :: mycvs = 
00105        '$Id: psmile_enddef_action.F90 3248 2011-06-23 13:03:19Z coquart $'
00106 !
00107 !----------------------------------------------------------------------
00108 !
00109 !  Initialization
00110 !
00111       ierror = 0
00112 !
00113       new_intersection = .false.
00114       sender = status (MPI_SOURCE)
00115 !
00116 #ifdef VERBOSE
00117       print 9990, trim(ch_id), index, sender
00118 
00119       call psmile_flushstd
00120 #endif /* VERBOSE */
00121 !
00122 !===> A request was completed
00123 !
00124       if (index <= 2) then
00125          call MPI_Get_count (status, MPI_INTEGER, meslen, ierror)
00126          if (meslen > nd_msgint .or. meslen < ind_msgint_tag) then
00127             ierror = PRISM_Error_Wrong
00128             ierrp (1) = meslen
00129             ierrp (2) = status (MPI_SOURCE)
00130             ierrp (3) = status (MPI_TAG)
00131 
00132             call psmile_error ( ierror, 'MPI_Get_count', &
00133                                 ierrp, 3, __FILE__, __LINE__ )
00134             return
00135          endif
00136       endif
00137 !
00138       select case (index)
00139 !
00140 !-----------------------------------------------------------------------
00141 !
00142 !        Perform action required
00143 !
00144 !        index = 1: Request for sending grid data found.
00145 !                   Send subgrid to partner which performs search of
00146 !                   donor cells.
00147 !
00148 !-----------------------------------------------------------------------
00149 !
00150       case ( 1 )
00151 
00152          call psmile_init_enddef_msg_inters (msg_intersections)
00153          call psmile_unpack_msg_intersections (msg_intersections, &
00154                                                paction%msgreq)
00155 !
00156          ind = msg_intersections%relative_msg_tag
00157 !
00158 #ifdef PRISM_ASSERTION
00159          if (status(MPI_TAG) /= reqtag) then
00160             call psmile_assert ( __FILE__, __LINE__, &
00161                                 'wrong message tag')
00162          endif
00163 !
00164 !    ... Control Range of loc_messages
00165 !
00166          if (ind < 1 .or. ind > paction%n_answer2recv) then
00167             call psmile_assert ( __FILE__, __LINE__, &
00168                                 'Illegal local tag')
00169          endif
00170 #endif
00171 !
00172          call psmile_send_req_subgrid (msg_intersections, &
00173                                        sender, grdtag, ierror)
00174          if (ierror > 0) return
00175 !
00176 !===> ... Control whether the locations of the previous grid request.
00177 !         The locations are sent shortly before the request for a new
00178 !         grid. Therefore, this routine may find the request for 
00179 !         a new grid before the locations sent.
00180 !
00181          if (paction%lrequest (num_req_types-1+ind) /= MPI_REQUEST_NULL) then
00182 
00183             call MPI_Wait (paction%lrequest (num_req_types-1+ind), lstatus, ierror)
00184 
00185             if ( ierror /= MPI_SUCCESS ) then
00186                ierrp (1) = ierror
00187                ierrp (2) = lstatus (MPI_SOURCE)
00188                ierrp (3) = lstatus (MPI_TAG)
00189 
00190                ierror = PRISM_Error_MPI
00191 
00192                call psmile_error ( ierror, 'MPI_Wait', &
00193                                    ierrp, 3, __FILE__, __LINE__ )
00194                return
00195             endif
00196 
00197 #ifdef PRISM_ASSERTION
00198             if (lstatus(MPI_TAG) /= loctag+ind) then
00199                call psmile_assert ( __FILE__, __LINE__, &
00200                                    'wrong message tag')
00201             endif
00202 !
00203             if (lstatus(MPI_SOURCE) /= sender) then
00204                call psmile_assert ( __FILE__, __LINE__, &
00205                                    'wrong sender')
00206             endif
00207 #endif
00208 
00209             call psmile_unpack_msg_locations(msg_locations, &
00210                         paction%loc_messages (:, ind))
00211 
00212             if ( msg_locations%requires_conserv_remap == PSMILe_conserv2D ) then
00213 !
00214 !===> Message indicates that we have to work on the cell-information first.
00215 !
00216 #ifdef PRISM_ASSERTION
00217                if (msg_locations%src_rank /= sender) then
00218                   call psmile_assert ( __FILE__, __LINE__, &
00219                        'wrong sender')
00220                endif
00221 #endif
00222                call psmile_enddef_action_cell (msg_locations, ierror)
00223                if (ierror > 0) return
00224 !
00225 !===> ... Finally receive the actual locations
00226 !
00227 #ifdef DEBUGX
00228                print *, ' Receiving locations from', sender, 'with tag ', loctag+ind, &
00229                         ' size ', msgloc_size
00230 #endif /* DEBUGX */
00231                call MPI_Recv (paction%loc_messages (1,ind), msgloc_size, &
00232                     MPI_INTEGER, sender, loctag+ind, comm_psmile, lstatus, ierror)
00233 
00234                if ( ierror /= MPI_SUCCESS ) then
00235                   ierrp (1) = ierror
00236                   ierrp (2) = sender
00237                   ierrp (3) = loctag+num_req_types-1+ind
00238 
00239                   ierror = PRISM_Error_Recv
00240 
00241                   call psmile_error ( ierror, 'MPI_Recv', &
00242                        ierrp, 3, __FILE__, __LINE__ )
00243                   return
00244                endif
00245 
00246             endif ! msg_locations%requires_conserv_remap == PSMILe_conserv2D
00247 
00248             call psmile_unpack_msg_locations(msg_locations, &
00249                         paction%loc_messages (:, ind))
00250 
00251             call psmile_enddef_action_loc (msg_locations, ierror)
00252             if (ierror > 0) return
00253 
00254          endif ! paction%lrequest (num_req_types-1+ind) /= MPI_REQUEST_NULL
00255 !
00256 !===> ... Set up request for receive of locations
00257 !
00258          call MPI_Irecv (paction%loc_messages(1,ind), msgloc_size, MPI_INTEGER, &
00259                          sender, loctag+ind, comm_psmile, &
00260                          paction%lrequest (num_req_types-1+ind), ierror)
00261 
00262          if ( ierror /= MPI_SUCCESS ) then
00263             ierrp (1) = ierror
00264             ierrp (2) = sender
00265             ierrp (3) = loctag+ind
00266 
00267             ierror = PRISM_Error_Recv
00268 
00269             call psmile_error ( ierror, 'MPI_Irecv', &
00270                                 ierrp, 3, __FILE__, __LINE__ )
00271             return
00272          endif
00273 #ifdef DEBUGX
00274          print *, ' Posting Irecv request(',num_req_types-1+ind,') ', &
00275              paction%lrequest(num_req_types-1+ind), 'with tag ', loctag+ind, &
00276                   ' and size ', msgloc_size
00277 #endif /* DEBUGX */
00278 !
00279 !===> ... Set up request for the next grid transfer
00280 !         Note: paction%n_answer is incremented here since it
00281 !               is also used by PSMILe_Enddef_action_loc.
00282 !
00283          paction%n_answer = paction%n_answer + 1
00284 !
00285          if (paction%n_answer < paction%n_answer2recv) then
00286 #ifdef PRISM_ASSERTION
00287             if (paction%lrequest(1) /= MPI_REQUEST_NULL) then
00288                call psmile_assert ( __FILE__, __LINE__, &
00289                                    'paction%lrequest(1) is not finished !')
00290             endif
00291 #endif
00292 !
00293             call MPI_Irecv (paction%msgreq, nd_msgint, MPI_INTEGER,      &
00294                             MPI_ANY_SOURCE, reqtag, comm_psmile, &
00295                             paction%lrequest(1), ierror)
00296 
00297             if ( ierror /= MPI_SUCCESS ) then
00298                ierrp (1) = ierror
00299                ierror = PRISM_Error_MPI
00300 
00301                call psmile_error ( ierror, 'MPI_Irecv', &
00302                                    ierrp, 1, __FILE__, __LINE__ )
00303                return
00304             endif
00305 #ifdef DEBUGX
00306             print *, ' Posting Irecv request(1) ', paction%lrequest(1), 'with tag ', reqtag
00307 #endif /* DEBUGX */
00308 #ifdef PRISM_ASSERTION
00309          else if (paction%lrequest (1) /= MPI_REQUEST_NULL) then
00310             call psmile_assert ( __FILE__, __LINE__, &
00311                                 'bad request; should be MPI_REQUEST_NULL')
00312 #endif
00313          endif
00314 !
00315       case (2)
00316 !
00317 !---------------------------------------------------------------------------
00318 !        index = 2: Message from a partner with an intersection was found
00319 !                   (cf. routine psmile_find_intersect)
00320 !        (1) send acknowledge (using tag reqtag)
00321 !        (2) allocate space
00322 !        (3) receive grid function
00323 !---------------------------------------------------------------------------
00324 !
00325          call psmile_init_enddef_msg_inters (msg_intersections)
00326          call psmile_unpack_msg_intersections (msg_intersections, &
00327                                                paction%msgint)
00328 
00329 #ifdef PRISM_ASSERTION
00330          if (status(MPI_TAG) /= paction%lastag) then
00331             call psmile_assert ( __FILE__, __LINE__, &
00332                                 'wrong message tag; must be lastag')
00333          endif
00334 #endif
00335         paction%n = paction%n + 1
00336 !
00337         npart = msg_intersections%num_parts
00338 !
00339 !===> ... Ignore dummy message (empty intersection)
00340 !
00341          if (npart <= 0) then
00342             new_intersection = .true.
00343          else
00344 !
00345 !===> ...   Send message back as acknowledge to sending process and
00346 !           as request to send the grid data.
00347 !           Note: The sending process may be the same process
00348 !               (i.e. the send is not allowed to block).
00349 !
00350             call psmile_bsend (paction%msgint, nd_msgint, MPI_INTEGER, &
00351                                sender, reqtag, comm_psmile, ierror)
00352 
00353             if (ierror /= MPI_SUCCESS) then
00354                call psmile_error (PRISM_Error_Send, 'MPI_Send', &
00355                                   (/ierror, sender, reqtag/), 3, &
00356                                   __FILE__, __LINE__ )
00357                ierror = PRISM_Error_Send
00358                return
00359             endif
00360 !
00361 !===> ... Set up receive the subgrid sent by process sender
00362 !         Note: search%search_data%npart = npart is set in PSMILe_Recv_req_subgrid
00363 !
00364             Allocate (paction%recv_req ((ndim_3d+2)*npart), stat = ierror)
00365             if ( ierror /= 0 ) then
00366                ierrp (1) = (ndim_3d+1) * npart
00367                ierror = PRISM_Error_Alloc
00368                call psmile_error ( ierror, 'paction%recv_req', &
00369                                    ierrp, 1, __FILE__, __LINE__ )
00370                return
00371             endif
00372 
00373             paction%recv_req(:) = MPI_REQUEST_NULL
00374 
00375             call psmile_recv_req_subgrid (msg_intersections, sender, grdtag, search, &
00376                                           paction%recv_req, new_search, ierror)
00377             if (ierror > 0) return
00378 
00379             paction%lrequest (3) = paction%recv_req (1)
00380             paction%grid2receive = .true.
00381             new_intersection = .false.
00382          endif ! if (npart <= 0)
00383 
00384       case ( 3 )
00385 !
00386 !-----------------------------------------------------------------------
00387 !        index = 3: Grid data was received (tag ``grdtag'')
00388 !                   Note: The first request of "paction%recv_req" was already
00389 !                         fulfilled.
00390 !
00391 !-----------------------------------------------------------------------
00392 !
00393          ierror = PRISM_Error_Internal
00394 
00395          call psmile_error ( ierror, 'Receive of grid data is not supported in PSMILe_Enddef_action', &
00396                              ierrp, 0, __FILE__, __LINE__ )
00397 !
00398       case ( 4 )
00399 !
00400 !-----------------------------------------------------------------------
00401 !        index = 4: Extra search request (tag ``exttag'') was found
00402 !
00403 !-----------------------------------------------------------------------
00404 !
00405          if (paction%msg_extra(1) == PSMILe_Finalize_extra_search) then
00406             paction%n_fin = paction%n_fin + paction%msg_extra(2)
00407 
00408 #ifdef PRISM_ASSERTION
00409             if (paction%n_fin > paction%n_fin2recv) then
00410                print *, 'n_fin, n_fin2recv', paction%n_fin, paction%n_fin2recv
00411                call psmile_assert (__FILE__, __LINE__, &
00412                   "Too many Finalize messages for extra search received")
00413             endif
00414 #endif
00415          else
00416             call psmile_unpack_msg_extra (msg_extra, paction%msg_extra)
00417          
00418             call psmile_enddef_action_extra (msg_extra, sender, ierror)
00419             if (ierror > 0) return
00420          endif
00421 !
00422 !===> ... Set up new request
00423 !
00424          if (paction%n_fin < paction%n_fin2recv) then
00425             call MPI_Irecv (paction%msg_extra, msg_extra_size, MPI_INTEGER, &
00426                             MPI_ANY_SOURCE, exttag, comm_psmile,            &
00427                             paction%lrequest(4), ierror)
00428 
00429             if ( ierror /= MPI_SUCCESS ) then
00430                ierrp (1) = ierror
00431                ierror = PRISM_Error_MPI
00432 
00433                call psmile_error ( ierror, 'MPI_Irecv', &
00434                                    ierrp, 1, __FILE__, __LINE__ )
00435                return
00436             endif
00437 #ifdef DEBUGX
00438             print *, ' Posting Irecv request(4) ', paction%lrequest(4), 'with tag ', exttag
00439 #endif /* DEBUGX */
00440 
00441          endif
00442 !
00443       case ( 5 )
00444 !
00445 !-----------------------------------------------------------------------
00446 !        index = 5: Info on seleced point of nearest neighbour search
00447 !                    (tag ``seltag'') was found
00448 !-----------------------------------------------------------------------
00449 !
00450          call psmile_enddef_action_sel (status(MPI_SOURCE), ierror)
00451          if (ierror > 0) return
00452 !
00453       case (num_req_types:)
00454 !
00455 !-----------------------------------------------------------------------
00456 !     index in (num_req_types:nreq): Locations (tag "loctag" were received)
00457 !     Request were set up with tag "loctag"
00458 !-----------------------------------------------------------------------
00459 !
00460          call psmile_unpack_msg_locations(msg_locations, &
00461                      paction%loc_messages (:, index-num_req_types+1))
00462 
00463 
00464          if ( msg_locations%requires_conserv_remap == PSMILe_conserv2D ) then
00465 !
00466 !===> Message indicaes that we have to work on the cell-information first.
00467 !
00468 #ifdef PRISM_ASSERTION
00469             if (msg_locations%src_rank /= sender) then
00470                call psmile_assert ( __FILE__, __LINE__, &
00471                                    'wrong sender')
00472             endif
00473 #endif
00474             call psmile_enddef_action_cell (msg_locations, ierror)
00475             if (ierror > 0) return
00476 !
00477 !===> ... Finally set up a new request for receive of locations
00478 !
00479             call MPI_Irecv (paction%loc_messages (1, index-num_req_types+1), msgloc_size, &
00480                             MPI_INTEGER, sender, loctag+index-num_req_types+1, comm_psmile, &
00481                             paction%lrequest(index), ierror)
00482 
00483             if ( ierror /= MPI_SUCCESS ) then
00484                ierrp (1) = ierror
00485                ierrp (2) = sender
00486                ierrp (3) = loctag+index-num_req_types+1
00487 
00488                ierror = PRISM_Error_Recv
00489 
00490                call psmile_error ( ierror, 'MPI_Irecv', &
00491                     ierrp, 3, __FILE__, __LINE__ )
00492                return
00493             endif
00494 
00495 #ifdef DEBUGX
00496             print *, ' Posting Irecv request(', index, ')', &
00497                  paction%lrequest(index), 'with tag ',  loctag+index-num_req_types+1, &
00498                      ' and size ', msgloc_size
00499 #endif /* DEBUGX */
00500          else
00501 
00502             call psmile_enddef_action_loc (msg_locations, ierror)
00503             if (ierror > 0) return
00504 
00505          endif
00506 
00507       end select 
00508 !
00509 !===> ... Set up request for receive of an intersection
00510 !
00511 !
00512       if (new_intersection .and. paction%n < paction%ninter) then
00513 #ifdef PRISM_ASSERTION
00514          if (paction%lrequest(2) /= MPI_REQUEST_NULL) then
00515             call psmile_assert ( __FILE__, __LINE__, &
00516                                 'paction%lrequest(2) is not finished !')
00517          endif
00518 #endif
00519 
00520 #define CONS_REMAP_DEADLOCK_WORKAROUND
00521 #ifndef CONS_REMAP_DEADLOCK_WORKAROUND
00522          call MPI_Irecv (paction%msgint, nd_msgint, MPI_INTEGER,      &
00523                          MPI_ANY_SOURCE, paction%lastag, comm_psmile, &
00524                          paction%lrequest(2), ierror)
00525 #else
00526          call MPI_Irecv (paction%msgint, nd_msgint, MPI_INTEGER,      &
00527                          maxval(paction%intersect_ranks), paction%lastag, comm_psmile, &
00528                          paction%lrequest(2), ierror)
00529          paction%intersect_ranks(maxloc(paction%intersect_ranks)) = -1
00530 #endif
00531          if ( ierror /= MPI_SUCCESS ) then
00532             ierrp (1) = ierror
00533             ierror = PRISM_Error_MPI
00534 
00535             call psmile_error ( ierror, 'MPI_Irecv', &
00536                                 ierrp, 1, __FILE__, __LINE__ )
00537             return
00538          endif
00539 #ifdef DEBUGX
00540          print *, ' Posting Irecv request(2)', paction%lrequest(2), 'with tag ', paction%lastag
00541 #endif /* DEBUGX */
00542       endif
00543 !
00544 !===> All done
00545 !
00546 #ifdef VERBOSE
00547       print 9980, trim(ch_id), ierror
00548       call psmile_flushstd
00549 #endif /* VERBOSE */
00550 
00551 !
00552 !  Formats:
00553 !
00554 
00555 #ifdef VERBOSE
00556 
00557 9990 format (1x, a, ': psmile_enddef_action: index', i4, ', sender', i4)
00558 9980 format (1x, a, ': psmile_enddef_action: eof ierror =', i3)
00559 
00560 #endif /* VERBOSE */
00561 
00562       end subroutine PSMILe_Enddef_action

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1