00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 subroutine psmile_enddef_action_sel (sender, ierror)
00013
00014
00015
00016 use PRISM_constants
00017
00018 use PSMILe, dummy_interface => PSMILe_Enddef_action_sel
00019
00020 Implicit none
00021
00022
00023
00024 Integer, Intent (In) :: sender
00025
00026
00027
00028
00029
00030 Integer, Intent (Out) :: ierror
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040 Type (Method), Pointer :: mp
00041
00042
00043
00044 Integer :: i
00045
00046
00047
00048 Type (Send_appl_information), Pointer :: send_info
00049 Integer :: index, n_liste, nsel
00050 Integer, Pointer :: neigh_3d (:, :)
00051
00052
00053
00054 Integer, Allocatable :: ibuf (:)
00055 Integer :: status (MPI_STATUS_SIZE)
00056
00057
00058
00059 Integer, Parameter :: nerrp = 1
00060 Integer :: ierrp (nerrp)
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082 Character(len=len_cvs_string), save :: mycvs =
00083 '$Id: psmile_enddef_action_sel.F90 2082 2009-10-21 13:31:19Z hanke $'
00084
00085
00086
00087
00088
00089 #ifdef VERBOSE
00090 print 9990, trim(ch_id), sender
00091
00092 call psmile_flushstd
00093 #endif /* VERBOSE */
00094
00095 #ifdef PRISM_ASSERTION
00096 if (paction%lrequest(5) /= MPI_REQUEST_NULL) then
00097 call psmile_assert ( __FILE__, __LINE__, &
00098 'paction%lrequest(5) should be MPI_REQUEST_NULL !')
00099 endif
00100
00101 if (paction%msg_sel (1) /= PSMILe_nnghbr3D) then
00102 print *, 'msg_sel(1)', paction%msg_sel (1), PSMILe_nnghbr3D
00103 call psmile_assert ( __FILE__, __LINE__, &
00104 'Incorrect message received !')
00105 endif
00106
00107 if (paction%msg_sel (2) < paction%msg_sel (3)) then
00108 print *, 'msg_sel(2:3)', paction%msg_sel (2:3)
00109 call psmile_assert ( __FILE__, __LINE__, &
00110 'n_liste must be >= nsel !')
00111 endif
00112
00113 if (paction%msg_sel (5) < 1 .or. &
00114 paction%msg_sel (5) > Number_of_Methods_allocated) then
00115 print *, trim(ch_id), "method id =", &
00116 paction%msg_sel (5), Number_of_Methods_allocated
00117 call psmile_assert ( __FILE__, __LINE__, &
00118 "invalid method id")
00119 endif
00120
00121 if (paction%n_selected < 1) then
00122 print *, "n_selected", paction%n_selected
00123 call psmile_assert ( __FILE__, __LINE__, &
00124 'No request outstanding !')
00125 endif
00126 #endif
00127
00128
00129
00130 paction%n_selected = paction%n_selected - 1
00131
00132 n_liste = paction%msg_sel (2)
00133 nsel = paction%msg_sel (3)
00134 index = paction%msg_sel (4)
00135 #ifdef DEBUG
00136 print *, 'nsel, nliste', nsel, n_liste
00137 #endif
00138
00139 mp => Methods (paction%msg_sel (5))
00140 send_info => mp%send_infos_appl(index)
00141
00142 #ifdef PRISM_ASSERTION
00143 if (index < 1 .or. index > mp%n_send_info_appl) then
00144 print *, "index", index, mp%n_send_info_appl
00145 call psmile_assert ( __FILE__, __LINE__, &
00146 'invalid send_info_appl index !')
00147 endif
00148
00149 if (send_info%dest /= sender) then
00150 call psmile_assert ( __FILE__, __LINE__, &
00151 "Sender doesn't fit !")
00152 endif
00153 #endif
00154
00155 if (nsel == 0) then
00156
00157
00158
00159 send_info%nloc = 0
00160
00161 Deallocate (send_info%neigh_3d)
00162 Nullify (send_info%neigh_3d)
00163
00164 #ifdef TODO
00165
00166
00167 if (mp%n_send_info_appl == index) then
00168
00169
00170
00171
00172
00173 endif
00174 #endif
00175
00176 else if (nsel < n_liste) then
00177
00178
00179
00180
00181 send_info%nloc = nsel
00182
00183 Allocate (neigh_3d(ndim_3d, nsel), ibuf (nsel), stat = ierror)
00184
00185 if (ierror /= 0) then
00186 ierrp (1) = ndim_3d * nsel + nsel
00187
00188 ierror = PRISM_Error_Alloc
00189
00190 call psmile_error ( ierror, 'neigh_3d', ierrp, 1, &
00191 __FILE__, __LINE__ )
00192 return
00193 endif
00194
00195 call MPI_Recv (ibuf, nsel, MPI_INTEGER, &
00196 sender, seltag, comm_psmile, &
00197 status, ierror)
00198 if ( ierror /= MPI_SUCCESS ) then
00199 ierrp (1) = ierror
00200 ierror = PRISM_Error_MPI
00201
00202 call psmile_error ( ierror, 'MPI_Recv(ibuf)', ierrp, 1, &
00203 __FILE__, __LINE__ )
00204 return
00205 endif
00206
00207
00208 do i = 1, nsel
00209 neigh_3d (:, i) = send_info%neigh_3d (:, ibuf(i))
00210 end do
00211
00212 Deallocate (send_info%neigh_3d)
00213 Deallocate (ibuf)
00214
00215 send_info%neigh_3d => neigh_3d
00216 send_info%nloc = nsel
00217 endif
00218
00219
00220
00221 if (paction%n_selected > 0) then
00222 call MPI_Irecv (paction%msg_sel, nd_msgsel, MPI_INTEGER, &
00223 MPI_ANY_SOURCE, seltag, comm_psmile, &
00224 paction%lrequest(5), ierror)
00225
00226 if ( ierror /= MPI_SUCCESS ) then
00227 ierrp (1) = ierror
00228 ierror = PRISM_Error_MPI
00229
00230 call psmile_error ( ierror, 'MPI_Irecv(seltag)', ierrp, 1, &
00231 __FILE__, __LINE__ )
00232 return
00233 endif
00234 endif
00235
00236
00237
00238 #ifdef VERBOSE
00239 print 9980, trim(ch_id), ierror
00240
00241 call psmile_flushstd
00242 #endif /* VERBOSE */
00243
00244
00245
00246
00247 #ifdef VERBOSE
00248
00249 9990 format (1x, a, ': psmile_enddef_action_sel: sender ', i6)
00250 9980 format (1x, a, ': psmile_enddef_action_sel: eof ierror =', i3)
00251
00252 #endif /* VERBOSE */
00253
00254 end subroutine PSMILe_Enddef_action_sel