00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_add_nn_found_dble (search, extra_search, &
00012 index_sent, found, n_send, &
00013 index_found, distance, n_found, nb_extra, &
00014 selected, sel_info, nrecv, ierror)
00015
00016
00017
00018 use PRISM
00019
00020 use PSMILe, dummy_interface => PSMILe_Add_nn_found_dble
00021
00022 Implicit none
00023
00024
00025
00026 Type (Enddef_search), Intent (In) :: search
00027
00028
00029
00030 Integer, Intent (In) :: n_send
00031
00032
00033
00034
00035 Integer, Intent (In) :: index_sent (n_send)
00036
00037
00038
00039
00040
00041 Integer, Intent (In) :: found (n_send)
00042
00043
00044
00045 Integer, Intent (In) :: n_found
00046
00047
00048
00049
00050 Integer, Intent (In) :: index_found (n_found)
00051
00052
00053
00054
00055 Double Precision, Intent (In) :: distance (n_found)
00056
00057
00058
00059 Integer, Intent (In) :: nb_extra
00060
00061
00062
00063
00064
00065
00066
00067 Integer, Intent (InOut) :: nrecv
00068
00069
00070
00071
00072
00073 Type (Extra_search_info), Intent (InOut) :: extra_search
00074
00075
00076
00077
00078 Integer, Intent (InOut) :: selected (2, extra_search%n_extra)
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088 Type (Select_search_info), Intent (InOut) :: sel_info (nrecv)
00089
00090
00091
00092
00093
00094
00095
00096 Integer, Intent (Out) :: ierror
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106 Integer, Parameter :: masked_out = 0
00107
00108
00109
00110
00111
00112 Integer :: i, j, n
00113
00114
00115
00116 Integer :: n_liste, nsel
00117 Integer :: prev
00118 Integer, Pointer :: used (:)
00119 Double Precision, Pointer :: dist_dble (:, :)
00120
00121
00122
00123 Integer :: msg_sel (nd_msgsel)
00124
00125
00126
00127 Integer, Parameter :: nerrp = 3
00128 Integer :: ierrp (nerrp)
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147 Character(len=len_cvs_string), save :: mycvs =
00148 '$Id: psmile_add_nn_found_dble.F90 2082 2009-10-21 13:31:19Z hanke $'
00149
00150
00151
00152
00153
00154 #ifdef VERBOSE
00155 print 9990, trim(ch_id), n_send, n_found, sel_info(nrecv)%n_liste
00156
00157 call psmile_flushstd
00158 #endif /* VERBOSE */
00159
00160 ierror = 0
00161
00162 dist_dble => extra_search%dist_dble
00163
00164
00165
00166 n_liste = sel_info(nrecv)%n_liste
00167
00168 #ifdef PRISM_ASSERTION
00169
00170
00171
00172
00173
00174
00175 do i = 1, n_found
00176 if (index_found (i) <= 0 .and. &
00177 index_found (i) /= masked_out) exit
00178 end do
00179
00180 if (i <= n_found) then
00181 print *, 'i, n_found, index_found', i, n_found, index_found (i)
00182 call psmile_assert ( __FILE__, __LINE__, &
00183 'invalid index in found list')
00184 endif
00185
00186
00187 do i = 1, n_found
00188 if (index_found (i) > n_liste) exit
00189 end do
00190
00191 if (i <= n_found) then
00192 print *, 'i, n_found, index_found, n_liste', &
00193 i, n_found, index_found (i), n_liste
00194 call psmile_assert ( __FILE__, __LINE__, &
00195 'invalid index (> n_liste) in found list')
00196 endif
00197
00198 #endif
00199
00200
00201
00202
00203 Allocate (sel_info(nrecv)%used(n_liste), stat = ierror)
00204 if (ierror /= 0) then
00205 ierrp (1) = n_liste
00206
00207 ierror = PRISM_Error_Alloc
00208
00209 call psmile_error ( ierror, 'sel_info(nrecv)%used', ierrp, 1, &
00210 __FILE__, __LINE__ )
00211 return
00212 endif
00213
00214 used => sel_info(nrecv)%used
00215
00216 used (:) = 0
00217
00218
00219 do n = 1, n_found
00220 if (index_found(n) > 0) then
00221 used (index_found(n)) = used (index_found(n)) + 1
00222 endif
00223 end do
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236 nsel = 0
00237 n = 0
00238 do j = 1, n_send
00239 if (found (j) > 0) then
00240 n = n + 1
00241 if (index_found(n) > 0) then
00242 i = index_sent (j)
00243
00244 #ifdef DEBUGX
00245 print *, 'psmile_add_nn_found_dble: n, j, i, i_liste, dists', &
00246 n, j, i, index_found(n), dist_dble (i,1), distance (n)
00247 #endif
00248
00249 if (dist_dble (i,1) <= distance (n)) then
00250
00251
00252
00253
00254
00255 used (index_found(n)) = used (index_found(n)) - 1
00256
00257 else
00258
00259
00260
00261
00262
00263 if (selected (1, i) > 0) then
00264 prev = selected (1, i)
00265
00266 #ifdef PRISM_ASSERTION
00267 if (prev > nrecv) then
00268 print *, trim(ch_id), "prev, nrecv", prev, nrecv
00269 call psmile_assert ( __FILE__, __LINE__, &
00270 'invalid recv request')
00271 endif
00272 #endif
00273
00274 sel_info(prev)%used (selected (2, i)) = &
00275 sel_info(prev)%used (selected (2, i)) - 1
00276 endif
00277
00278 nsel = nsel + 1
00279
00280 dist_dble (i,1) = distance (n)
00281 selected (1, i) = nrecv
00282 selected (2, i) = index_found(n)
00283 endif
00284 endif
00285 endif
00286 end do
00287
00288
00289
00290
00291
00292
00293
00294
00295 if (nsel == 0) then
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305 msg_sel (1) = PSMILe_nnghbr3D
00306 msg_sel (2) = n_liste
00307 msg_sel (3) = 0
00308 msg_sel (4) = sel_info(nrecv)%index
00309 msg_sel (5) = sel_info(nrecv)%method_id
00310
00311
00312
00313 call MPI_Send (msg_sel, nd_msgsel, MPI_INTEGER, &
00314 sel_info(nrecv)%sender, seltag, comm_psmile, ierror)
00315
00316 if (ierror /= MPI_SUCCESS) then
00317 ierrp (1) = ierror
00318 ierrp (2) = sel_info(nrecv)%sender
00319 ierrp (3) = seltag
00320
00321 ierror = PRISM_Error_Send
00322
00323 call psmile_error (ierror, 'MPI_send (msg_sel)', ierrp, 3, &
00324 __FILE__, __LINE__ )
00325 return
00326 endif
00327
00328
00329
00330 Deallocate (sel_info(nrecv)%used)
00331 Deallocate (sel_info(nrecv)%dble_buf)
00332
00333
00334
00335 nrecv = nrecv - 1
00336 endif
00337
00338
00339
00340 #ifdef VERBOSE
00341 print 9980, trim(ch_id), ierror, nsel
00342
00343 call psmile_flushstd
00344 #endif /* VERBOSE */
00345
00346
00347
00348
00349 #ifdef VERBOSE
00350
00351 9990 format (1x, a, ': psmile_add_nn_found_dble: n_send', i7, &
00352 '; n_found', i7, '; n_liste ', i7)
00353 9980 format (1x, a, ': psmile_add_nn_found_dble: eof ierror =', i3, &
00354 '; nsel', i7)
00355
00356 #endif /* VERBOSE */
00357
00358 #ifdef DEBUG
00359 #endif
00360
00361 end subroutine PSMILe_Add_nn_found_dble