00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_send_destinations (method_id, cpl_index, dir_index, &
00012 sender, tag, ierror)
00013
00014
00015
00016 use PRISM_constants
00017
00018 use PSMILe, dummy_interface => PSMILe_Send_destinations
00019
00020 implicit none
00021
00022
00023
00024 Integer, Intent (In) :: method_id
00025
00026
00027
00028 Integer, Intent (In) :: cpl_index
00029
00030
00031
00032
00033 Integer, Intent (In) :: dir_index
00034
00035
00036
00037
00038 Integer, Intent (In) :: sender
00039
00040
00041
00042 Integer, Intent (In) :: tag
00043
00044
00045
00046
00047
00048 integer, Intent (Out) :: ierror
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058 Type (Method), Pointer :: mp
00059 Type(Send_information), Pointer :: send_info
00060
00061
00062
00063 Integer :: lrequest
00064 Integer :: n, i
00065 #ifdef PRISM_ASSERTION
00066 Integer :: nar, ntot
00067 #endif
00068
00069
00070
00071 Integer, parameter :: nerrp = 3
00072 Integer :: ierrp (nerrp)
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093 Character(len=len_cvs_string), save :: mycvs =
00094 '$Id: psmile_send_destinations.F90 2887 2011-01-14 09:25:49Z redler $'
00095
00096
00097
00098
00099
00100 #ifdef VERBOSE
00101 print 9990, trim(ch_id), cpl_index, dir_index
00102
00103 call psmile_flushstd
00104 #endif /* VERBOSE */
00105
00106 ierror = 0
00107
00108 mp => Methods(method_id)
00109
00110
00111
00112
00113
00114
00115 if (cpl_index /= PRISM_UNDEFINED) then
00116 n = mp%send_infos_coupler(cpl_index)%nloc
00117
00118
00119 #ifdef DEBUGX
00120 print *, ' Preparing Isend to ', sender, 'with tag ', tag, &
00121 ' request ', lrequest, ' of size ', n*ndim_3d
00122 #endif /* DEBUGX */
00123
00124 call MPI_Isend (mp%send_infos_coupler(cpl_index)%dstijk, n*ndim_3d, &
00125 MPI_INTEGER, &
00126 sender, tag, comm_psmile, lrequest, ierror)
00127
00128 if (ierror /= MPI_SUCCESS) then
00129 ierrp (1) = ierror
00130 ierrp (2) = sender
00131 ierrp (3) = tag
00132 ierror = PRISM_Error_Send
00133
00134 call psmile_error (ierror, 'MPI_Isend', &
00135 ierrp, 3, __FILE__, __LINE__ )
00136 return
00137 endif
00138
00139 call MPI_request_free (lrequest, ierror)
00140 if (ierror /= MPI_SUCCESS) then
00141 ierrp (1) = ierror
00142 ierror = PRISM_Error_MPI
00143
00144 call psmile_error (ierror, 'MPI_Request_free', &
00145 ierrp, 1, __FILE__, __LINE__ )
00146 return
00147 endif
00148
00149 endif
00150
00151
00152
00153
00154
00155
00156 if (dir_index /= PRISM_UNDEFINED) then
00157 send_info => mp%send_infos_direct(dir_index)
00158
00159
00160
00161 n = 0
00162 do i = 1, send_info%nparts
00163 n = n + Product(send_info%npoints (1:send_info%nvec,i))
00164 enddo
00165
00166 #ifdef PRISM_ASSERTION
00167 ntot = n
00168 #endif
00169 if (n > 0) then
00170 #ifdef DEBUGX
00171 print *, ' Preparing Isend to ', sender, 'with tag ', tag, &
00172 ' request ', lrequest, ' of size ', n*ndim_3d
00173 #endif /* DEBUGX */
00174 call MPI_Isend (send_info%dstijk, n*ndim_3d, MPI_INTEGER, &
00175 sender, tag, comm_psmile, lrequest, ierror)
00176
00177 if (ierror /= MPI_SUCCESS) then
00178 ierrp (1) = ierror
00179 ierrp (2) = sender
00180 ierrp (3) = tag
00181 ierror = PRISM_Error_Send
00182
00183 call psmile_error (ierror, 'MPI_Isend(dstijk)', &
00184 ierrp, 3, __FILE__, __LINE__ )
00185 return
00186 endif
00187
00188 call MPI_request_free (lrequest, ierror)
00189 if (ierror /= MPI_SUCCESS) then
00190 ierrp (1) = ierror
00191 ierror = PRISM_Error_MPI
00192
00193 call psmile_error (ierror, 'MPI_Request_free', &
00194 ierrp, 1, __FILE__, __LINE__ )
00195 return
00196 endif
00197 endif
00198
00199
00200
00201
00202 n = send_info%nars(1,1)
00203
00204 if (n > 0) then
00205
00206 #ifdef PRISM_ASSERTION
00207
00208
00209
00210 if (send_info%nvec > 1 .or. send_info%nparts > 1) then
00211 print *, 'nvec nparts', send_info%nvec, send_info%nparts
00212 call psmile_assert ( __FILE__, __LINE__, &
00213 'Routine is designed only for nvec = nparts = 1 for clustered areas')
00214 endif
00215
00216
00217
00218 do nar = 1, n
00219 ntot = ntot + &
00220 (send_info%dstars(2,1,nar)-send_info%dstars(1,1,nar)+1) * &
00221 (send_info%dstars(2,2,nar)-send_info%dstars(1,2,nar)+1) * &
00222 (send_info%dstars(2,3,nar)-send_info%dstars(1,3,nar)+1)
00223 end do
00224
00225 if (ntot /= send_info%nloc) then
00226 print *, 'nloc, ntot, nars', send_info%nloc, ntot, n
00227 call psmile_assert ( __FILE__, __LINE__, &
00228 "nloc doesn't match npoints + SUM(clustered areas)")
00229 endif
00230 #endif
00231
00232
00233
00234
00235 #ifdef DEBUGX
00236 print *, ' Preparing Isend to ', sender, 'with tag ', tag, &
00237 ' request ', lrequest, ' of size ', n*ndim_3d
00238 #endif /* DEBUGX */
00239 call MPI_Isend (send_info%dstars, n*2*ndim_3d, &
00240 MPI_INTEGER, sender, tag, comm_psmile, &
00241 lrequest, ierror)
00242
00243 if (ierror /= MPI_SUCCESS) then
00244 ierrp (1) = ierror
00245 ierrp (2) = sender
00246 ierrp (3) = tag
00247 ierror = PRISM_Error_Send
00248
00249 call psmile_error (ierror, 'MPI_Isend(dstars)', &
00250 ierrp, 3, __FILE__, __LINE__ )
00251 return
00252 endif
00253
00254 call MPI_request_free (lrequest, ierror)
00255 if (ierror /= MPI_SUCCESS) then
00256 ierrp (1) = ierror
00257 ierror = PRISM_Error_MPI
00258
00259 call psmile_error (ierror, 'MPI_Request_free', &
00260 ierrp, 1, __FILE__, __LINE__ )
00261 return
00262 endif
00263 endif
00264
00265 endif
00266
00267
00268
00269
00270
00271 #ifdef VERBOSE
00272 print 9980, trim(ch_id), ierror
00273
00274 call psmile_flushstd
00275 #endif /* VERBOSE */
00276
00277
00278
00279
00280 #ifdef VERBOSE
00281
00282 9990 format (1x, a, ': psmile_send_destinations: cpl/dir index ', 2i10)
00283 9980 format (1x, a, ': psmile_send_destinations: eof ierror =', i3)
00284
00285 #endif /* VERBOSE */
00286
00287 end subroutine PSMILe_Send_destinations