00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_locations_3d_reg ( &
00012 found, loc, range, control, &
00013 search, method_id, &
00014 dir_index, cpl_index, len_cpl, ierror)
00015
00016
00017
00018 use PRISM_constants
00019
00020 use PSMILe, dummy_interface => PSMILe_Locations_3d_reg
00021
00022 implicit none
00023
00024
00025
00026 Type (Enddef_search), Intent (InOut) :: search
00027
00028
00029
00030 Integer, Intent (In) :: range (2, ndim_3d, search%search_data%npart)
00031
00032
00033
00034
00035 Type (integer_vector) :: found (search%search_data%npart, ndim_3d)
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047 Type (integer_vector) :: loc (search%search_data%npart, ndim_3d)
00048
00049
00050
00051
00052 Integer, Intent (In) :: control (2, ndim_3d, search%search_data%npart)
00053
00054
00055
00056 Integer, Intent (In) :: method_id
00057
00058
00059
00060
00061
00062 Integer, Intent (Out) :: dir_index
00063
00064
00065
00066
00067
00068 Integer, Intent (Out) :: cpl_index
00069
00070
00071
00072
00073 Integer, Intent (Out) :: len_cpl (search%search_data%npart)
00074
00075
00076
00077 Integer, Intent (Out) :: ierror
00078
00079
00080
00081
00082
00083
00084
00085 Integer, Parameter :: val_direct = 1
00086 Integer, Parameter :: val_coupler = -1
00087 Integer, Parameter :: val_both = 0
00088
00089 Integer, Parameter :: send_coupler_index = 1
00090 Integer, Parameter :: send_direct_index = 2
00091
00092
00093
00094
00095
00096 Type (Method), Pointer :: mp
00097
00098
00099
00100 Real :: ratio
00101
00102 Integer :: n, ncpl_tot, ndir_tot
00103 Integer :: i, index, val_cpl
00104 Integer :: ncpl(ndim_3d, search%search_data%npart)
00105 Integer :: ndir(ndim_3d, search%search_data%npart)
00106
00107
00108
00109 Integer :: ipart, nadd, nprev
00110 Integer :: ibeg (ndim_3d, search%search_data%npart)
00111 Integer :: iend (ndim_3d, search%search_data%npart)
00112 Integer :: nprevi
00113
00114
00115
00116 Integer, parameter :: nerrp = 1
00117 Integer :: ierrp (nerrp)
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137 Character(len=len_cvs_string), save :: mycvs =
00138 '$Id: psmile_locations_3d_reg.F90 3248 2011-06-23 13:03:19Z coquart $'
00139
00140
00141
00142
00143
00144 #ifdef VERBOSE
00145 print 9990, trim(ch_id), method_id, search%msg_intersections%field_info%tgt_method_id, &
00146 search%sender
00147
00148 call psmile_flushstd
00149 #endif /* VERBOSE */
00150
00151 ierror = 0
00152
00153 cpl_index = PRISM_undefined
00154 dir_index = PRISM_undefined
00155
00156 len_cpl = 0
00157
00158 mp => Methods(method_id)
00159
00160 #ifdef PRISM_ASSERTION
00161 if (search%search_data%grid_type == PRISM_Irrlonlatvrt .or. &
00162 search%search_data%grid_type == PRISM_Irrlonlat_Regvrt) then
00163 print *, 'search%search_data%grid_type = ', search%search_data%grid_type
00164 call psmile_assert (__FILE__, __LINE__, &
00165 "This routine is not designed for this grid type")
00166 endif
00167
00168 if ( mp%status == PSMILe_Status_free .or. &
00169 mp%status == PSMILe_Status_undefined ) then
00170 call psmile_assert (__FILE__, __LINE__, "Incorrect method")
00171 endif
00172
00173 do ipart = 1, search%search_data%npart
00174 if (control (1,1,ipart) < range (1,1,ipart) .or. &
00175 range (2,1,ipart) < control (2,1,ipart) .or. &
00176 control (1,2,ipart) < range (1,2,ipart) .or. &
00177 range (2,2,ipart) < control (2,2,ipart) .or. &
00178 control (1,3,ipart) < range (1,3,ipart) .or. &
00179 range (2,3,ipart) < control (2,3,ipart)) then
00180 print *, ipart, control (:, :, ipart), range (:, :, ipart)
00181 call psmile_assert (__FILE__, __LINE__, "control is out of range")
00182 endif
00183 end do
00184 #endif
00185
00186
00187
00188
00189
00190
00191 ndir_tot = 0
00192 ncpl_tot = 0
00193
00194 ndir(:, :) = 0
00195 ncpl(:, :) = 0
00196
00197 ibeg (:, :) = 1
00198
00199 do ipart = 1, search%search_data%npart
00200 iend (:, ipart) = range (2,:, ipart) - range(1,:, ipart) + 1
00201
00202 do i = 1, ndim_3d
00203
00204
00205
00206
00207
00208
00209 do n = ibeg (i,ipart), iend (i, ipart)
00210 if (found(ipart, i)%vector(n) == val_coupler) then
00211 ncpl(i, ipart) = ncpl(i, ipart) + 1
00212 else if (found (ipart, i)%vector(n) == val_direct) then
00213 ndir(i, ipart) = ndir(i, ipart) + 1
00214 endif
00215 end do
00216 end do
00217
00218 end do
00219
00220
00221
00222
00223
00224 if (search%search_data%grid_type == PRISM_Reglonlatvrt) then
00225
00226 do ipart = 1, search%search_data%npart
00227
00228 n = ndir(1,ipart) * ndir (2,ipart) * ndir (3,ipart)
00229 ndir_tot = ndir_tot + n
00230 ncpl_tot = ncpl_tot + &
00231 ((ncpl(1, ipart)+ndir(1, ipart)) * &
00232 (ncpl(2, ipart)+ndir(2, ipart)) * &
00233 (ncpl(3, ipart)+ndir(3, ipart)) - n)
00234 end do
00235
00236 else
00237 ierror = PRISM_Error_Internal
00238 ierrp (1) = search%search_data%grid_type
00239 call psmile_error ( ierror, 'Hubert Error: search%search_data%grid_type', &
00240 ierrp, 1, __FILE__, __LINE__)
00241 return
00242 endif
00243
00244
00245 #ifdef DEBUG
00246 print *, "psmile_locations_3d_reg: ndir_tot, ndir", ndir_tot, ndir
00247 print *, "psmile_locations_3d_reg: ncpl_tot, ncpl", ncpl_tot, ncpl
00248 #endif
00249
00250 if (ncpl_tot + ndir_tot > 0) then
00251 ratio = real(ndir_tot) / (ncpl_tot+ndir_tot)
00252 else
00253 ratio = 0.0
00254 endif
00255
00256 #ifdef ONLY_FOR_TESTING
00257 print *, '######## psmile_locations_3d_reg: ratio set to 0.0'
00258 ratio = 0.01
00259 #endif
00260
00261
00262
00263
00264
00265
00266 if (maxval (ndir) > 0 .and. ratio < 0.05) then
00267
00268
00269
00270
00271
00272 ncpl = ncpl + ndir
00273 ndir = 0
00274
00275 ncpl_tot = ncpl_tot + ndir_tot
00276 ndir_tot = 0
00277
00278 val_cpl = val_both
00279 else
00280 val_cpl = val_coupler
00281
00282 endif
00283
00284
00285
00286
00287
00288 if (ncpl_tot > 0) then
00289 call psmile_get_info_index (method_id, send_coupler_index, index, ierror)
00290 if (ierror > 0) return
00291
00292 cpl_index = index
00293
00294
00295
00296 if (ndir_tot == 0) then
00297 mp%send_infos_coupler(index)%nvec = ndim_3d
00298 mp%send_infos_coupler(index)%nparts = search%search_data%npart
00299 else
00300
00301
00302
00303
00304 mp%send_infos_coupler(index)%nvec = 1
00305 mp%send_infos_coupler(index)%nparts = 1
00306 endif
00307
00308 mp%send_infos_coupler(index)%remote_method_id = &
00309 search%msg_intersections%field_info%tgt_method_id
00310
00311 mp%send_infos_coupler(index)%nrecv = 0
00312 mp%send_infos_coupler(index)%num2recv = 0
00313
00314
00315
00316 mp%send_infos_coupler(index)%dest = 0
00317
00318
00319
00320 call psmile_locations_alloc (mp%send_infos_coupler(index), ierror)
00321 if (ierror > 0) return
00322
00323 nprev = 0
00324 nprevi = 0
00325
00326 do ipart = 1, search%search_data%npart
00327
00328 call psmile_store_dest_locs_3d_reg (found (ipart, 1:ndim_3d), &
00329 loc (ipart, 1:ndim_3d), &
00330 range (:, :, ipart), &
00331 control (:, :, ipart), &
00332 mp%send_infos_coupler(index), ncpl_tot, &
00333 val_cpl, nprev, nadd, ierror)
00334 if (ierror > 0) return
00335
00336 len_cpl (ipart) = nadd
00337
00338 if (ndir_tot == 0) then
00339
00340 nprev = nprev + nadd
00341
00342 do i = 1, ndim_3d
00343 call psmile_store_source_locs_1d (found(ipart, i)%vector, &
00344 loc(ipart, i)%vector, &
00345 ibeg (i, ipart), iend (i, ipart), &
00346 mp%send_infos_coupler(index), &
00347 ncpl(i, ipart), &
00348 val_cpl, i, ipart, nprevi, nadd, ierror)
00349 if (ierror > 0) return
00350
00351 end do
00352 else
00353 call psmile_store_source_locs_3d_reg ( &
00354 found (ipart, 1:ndim_3d), loc (ipart, 1:ndim_3d), &
00355 range (:, :, ipart), control (:, :, ipart), &
00356 mp%send_infos_coupler(index), ncpl_tot, &
00357 val_cpl, nprev, nadd, ierror)
00358 if (ierror > 0) return
00359
00360 nprev = nprev + nadd
00361 endif
00362 end do
00363
00364 mp%send_infos_coupler(index)%nloc = nprev
00365 endif
00366
00367
00368
00369
00370
00371 if (ndir_tot > 0) then
00372
00373 call psmile_get_info_index (method_id, send_direct_index, index, ierror)
00374 if (ierror > 0) return
00375
00376 dir_index = index
00377
00378 mp%send_infos_direct(index)%dest = search%sender
00379
00380 mp%send_infos_direct(index)%nvec = ndim_3d
00381 mp%send_infos_direct(index)%nparts = search%search_data%npart
00382 mp%send_infos_direct(index)%remote_method_id = search%msg_intersections%field_info%tgt_method_id
00383
00384 mp%send_infos_direct(index)%nrecv = 0
00385 mp%send_infos_direct(index)%num2recv = 0
00386
00387
00388
00389 call psmile_locations_alloc (mp%send_infos_direct(index), ierror)
00390 if (ierror > 0) return
00391
00392 nprev = 0
00393 nprevi = 0
00394
00395 do ipart = 1, search%search_data%npart
00396
00397 call psmile_store_dest_locs_3d_reg (found(ipart, 1:ndim_3d), &
00398 loc (ipart, 1:ndim_3d), &
00399 range (:, :, ipart), &
00400 control (:, :, ipart), &
00401 mp%send_infos_direct(index), ndir_tot, &
00402 val_direct, nprev, nadd, ierror)
00403 if (ierror > 0) return
00404
00405 nprev = nprev + nadd
00406
00407 do i = 1, ndim_3d
00408 call psmile_store_source_locs_1d (found(ipart, i)%vector, &
00409 loc(ipart, i)%vector, &
00410 ibeg (i, ipart), iend (i, ipart), &
00411 mp%send_infos_direct(index), &
00412 ndir(i, ipart), &
00413 val_direct, i, ipart, nprevi, nadd, ierror)
00414 if (ierror > 0) return
00415
00416 end do
00417 end do
00418
00419 mp%send_infos_direct(index)%nloc = nprev
00420 endif
00421
00422
00423
00424 #ifdef VERBOSE
00425 print 9980, trim(ch_id), ierror
00426
00427 call psmile_flushstd
00428 #endif /* VERBOSE */
00429
00430
00431
00432
00433 #ifdef VERBOSE
00434
00435 9990 format (1x, a, ': psmile_locations_3d_reg: method_id', i3, &
00436 ' to ', i3, '(', i2, ')')
00437 9980 format (1x, a, ': psmile_locations_3d_reg: eof ierror =', i3)
00438
00439 #endif /* VERBOSE */
00440
00441 end subroutine PSMILe_Locations_3d_reg