00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_locations_gauss2 ( &
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_gauss2
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%npart)
00031
00032
00033
00034
00035 Type (integer_vector) :: found (search%npart, 2)
00036
00037
00038
00039
00040
00041 Type (integer_vector) :: loc (search%npart, 2)
00042
00043
00044
00045
00046
00047 Integer, Intent (In) :: control (2, ndim_3d, search%npart)
00048
00049
00050
00051 Integer, Intent (In) :: method_id
00052
00053
00054
00055
00056
00057 Integer, Intent (Out) :: dir_index
00058
00059
00060
00061
00062
00063 Integer, Intent (Out) :: cpl_index
00064
00065
00066
00067
00068 Integer, Intent (Out) :: len_cpl (search%npart)
00069
00070
00071
00072 Integer, Intent (Out) :: ierror
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083 Integer, Parameter :: val_direct = 1
00084 Integer, Parameter :: val_coupler = -1
00085 Integer, Parameter :: val_both = 0
00086
00087 Integer, Parameter :: send_coupler_index = 1
00088 Integer, Parameter :: send_direct_index = 2
00089
00090 Integer, Parameter :: indl = 1
00091 Integer, Parameter :: indz = 2
00092
00093
00094
00095 Type (Method), Pointer :: mp
00096
00097 Real :: ratio
00098
00099 Integer :: ncpl_tot, ndir_tot
00100 Integer :: ncpl (search%npart)
00101 Integer :: ndir (search%npart)
00102 Integer :: ncplz (search%npart)
00103 Integer :: ndirz (search%npart)
00104 Integer :: index, n, val_cpl
00105
00106
00107
00108 Integer :: ipart, nadd, nprev
00109 Integer :: ibeg (search%npart, indz)
00110 Integer :: iend (search%npart, indz)
00111 Integer :: nprevi
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139 Character(len=len_cvs_string), save :: mycvs =
00140 '$Id: psmile_locations_gauss2.F90 2966 2011-02-18 09:47:30Z hanke $'
00141
00142
00143
00144
00145
00146 #ifdef VERBOSE
00147 print 9990, trim(ch_id), method_id, search%msg_intersections%first_tgt_method_id, &
00148 search%sender
00149
00150 call psmile_flushstd
00151 #endif /* VERBOSE */
00152
00153 ierror = 0
00154 cpl_index = PRISM_undefined
00155 dir_index = PRISM_undefined
00156
00157 len_cpl = 0
00158
00159 mp => Methods(method_id)
00160
00161 #ifdef PRISM_ASSERTION
00162 if (search%grid_type == PRISM_Irrlonlatvrt) then
00163 print *, "search%grid_type = ", search%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%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 ndirz(:) = 0
00198 ncplz(:) = 0
00199
00200 ibeg (:, :) = 1
00201
00202 do ipart = 1, search%npart
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218 iend (ipart, indl) = (range (2,1, ipart) - range(1,1, ipart) + 1) * &
00219 (range (2,2, ipart) - range(1,2, ipart) + 1)
00220 iend (ipart, indz) = range (2,3, ipart) - range(1,3, ipart) + 1
00221
00222
00223 do n = ibeg(ipart, indl), iend (ipart, indl)
00224 if (found(ipart,indl)%vector(n) == val_coupler) then
00225 ncpl(ipart) = ncpl(ipart) + 1
00226 else if (found(ipart,indl)%vector(n) == val_direct) then
00227 ndir(ipart) = ndir(ipart) + 1
00228 endif
00229 end do
00230
00231
00232
00233
00234 do n = ibeg(ipart, indz), iend (ipart, indz)
00235 if (found(ipart,indz)%vector(n) == val_coupler) then
00236 ncplz(ipart) = ncplz(ipart) + 1
00237 else if (found(ipart,indz)%vector(n) == val_direct) then
00238 ndirz(ipart) = ndirz(ipart) + 1
00239 endif
00240 end do
00241
00242 end do
00243
00244 #ifdef DEBUG
00245 print *, trim(ch_id), ': psmile_locations_gauss2:'
00246 print *, "range: ", range
00247 print *, "control:", control
00248 print *, "ibeg:", ibeg
00249 print *, "iend:", iend
00250 print *, "ncpl :", ncpl, ncplz
00251 print *, "ndir :", ndir, ndirz
00252
00253 call psmile_flushstd
00254 #endif
00255
00256
00257
00258
00259
00260 do ipart = 1, search%npart
00261 n = ndir(ipart) * ndirz(ipart)
00262 ndir_tot = ndir_tot + n
00263 ncpl_tot = ncpl_tot + ((ncpl (ipart)+ndir (ipart)) * &
00264 (ncplz(ipart)+ndirz(ipart)) - n)
00265 end do
00266
00267 if (ncpl_tot + ndir_tot > 0) then
00268 ratio = real(ndir_tot) / (ncpl_tot+ndir_tot)
00269 else
00270 ratio = 0.0
00271 endif
00272
00273 #ifdef ONLY_FOR_TESTING
00274 print *, '######## psmile_locations_3d_reg: ratio set to 0.0'
00275 ratio = 0.01
00276 #endif
00277 #ifdef DEBUG
00278 print *, 'ncpl_tot, ndir_tot, ratio:', ncpl_tot, ndir_tot, ratio
00279 #endif
00280
00281
00282
00283
00284
00285
00286 if (max(maxval (ndir), maxval (ndirz)) > 0 .and. ratio < 0.05) then
00287
00288
00289
00290
00291
00292 ncpl = ncpl + ndir
00293 ndir = 0
00294
00295 ncplz = ncplz + ndirz
00296 ndirz = 0
00297
00298 ncpl_tot = ncpl_tot + ndir_tot
00299 ndir_tot = 0
00300
00301 val_cpl = val_both
00302 else
00303 val_cpl = val_coupler
00304 endif
00305
00306
00307
00308 if (ncpl_tot > 0) then
00309
00310 call psmile_get_info_index (method_id, send_coupler_index, index, ierror)
00311
00312 if (ierror > 0) return
00313
00314 cpl_index = index
00315
00316
00317
00318
00319 mp%send_infos_coupler(index)%nvec = 2
00320 mp%send_infos_coupler(index)%nparts = search%npart
00321
00322 mp%send_infos_coupler(index)%remote_method_id = search%msg_intersections%first_tgt_method_id
00323
00324
00325
00326 mp%send_infos_coupler(index)%dest = 0
00327
00328
00329
00330 call psmile_locations_alloc (mp%send_infos_coupler(index), ierror)
00331 if (ierror > 0) return
00332
00333 nprev = 0
00334 nprevi = 0
00335
00336 do ipart = 1, search%npart
00337
00338 call psmile_store_dest_locs_21d ( &
00339 found(ipart,indl)%vector, range(1:2, 1:ndim_3d, ipart), &
00340 control (1:2, 1:ndim_3d, ipart), found(ipart,indz)%vector, &
00341 mp%send_infos_coupler(index), ncpl_tot, &
00342 val_cpl, nprev, nadd, ierror)
00343 if (ierror > 0) return
00344
00345 nprev = nprev + nadd
00346 len_cpl (ipart) = nadd
00347
00348 call psmile_store_source_locs_1d ( &
00349 found(ipart,indl)%vector, loc(ipart,indl)%vector, &
00350 ibeg (ipart,indl), iend(ipart,indl), &
00351 mp%send_infos_coupler(index), ncpl(ipart), &
00352 val_cpl, indl, ipart, nprevi, nadd, ierror)
00353 if (ierror > 0) return
00354
00355
00356 call psmile_store_source_locs_1d ( &
00357 found(ipart,indz)%vector, loc(ipart,indz)%vector, &
00358 ibeg (ipart,indz), iend(ipart,indz), &
00359 mp%send_infos_coupler(index), ncplz(ipart), &
00360 val_cpl, indz, ipart, nprevi, nadd, ierror)
00361 if (ierror > 0) return
00362
00363 end do
00364
00365 mp%send_infos_coupler(index)%nloc = nprev
00366
00367 endif
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
00381 mp%send_infos_direct(index)%nvec = 2
00382 mp%send_infos_direct(index)%nparts = search%npart
00383 mp%send_infos_direct(index)%remote_method_id = search%msg_intersections%first_tgt_method_id
00384
00385
00386
00387 call psmile_locations_alloc (mp%send_infos_direct(index), ierror)
00388 if (ierror > 0) return
00389
00390 nprev = 0
00391 nprevi = 0
00392
00393 do ipart = 1, search%npart
00394
00395 call psmile_store_dest_locs_21d ( &
00396 found(ipart,indl)%vector, range(1:2, 1:ndim_3d, ipart), &
00397 control (1:2, 1:ndim_3d, ipart), found(ipart,indz)%vector, &
00398 mp%send_infos_direct(index), ndir_tot, &
00399 val_direct, nprev, nadd, ierror)
00400 if (ierror > 0) return
00401
00402 nprev = nprev + nadd
00403
00404 call psmile_store_source_locs_1d ( &
00405 found(ipart,indl)%vector, loc(ipart,indl)%vector, &
00406 ibeg (ipart,indl), iend(ipart,indl), &
00407 mp%send_infos_direct(index), ndir(ipart), &
00408 val_direct, indl, ipart, nprevi, nadd, ierror)
00409 if (ierror > 0) return
00410
00411 call psmile_store_source_locs_1d ( &
00412 found(ipart,indz)%vector, loc(ipart,indz)%vector, &
00413 ibeg (ipart,indz), iend(ipart,indz), &
00414 mp%send_infos_direct(index), ndirz(ipart), &
00415 val_direct, indz, ipart, nprevi, nadd, ierror)
00416 if (ierror > 0) return
00417
00418 end do
00419
00420 mp%send_infos_direct(index)%nloc = nprev
00421
00422 endif
00423
00424
00425
00426 #ifdef VERBOSE
00427 print 9980, trim(ch_id), ierror, dir_index, cpl_index
00428
00429 call psmile_flushstd
00430 #endif /* VERBOSE */
00431
00432
00433
00434
00435 #ifdef VERBOSE
00436
00437 9990 format (1x, a, ': psmile_locations_gauss2: method_id', i3, &
00438 ' to ', i3, '(', i2, ')')
00439 9980 format (1x, a, ': psmile_locations_gauss2: eof ierror =', i3, &
00440 ', dir_index =', i10, ', cpl_index =', i10)
00441
00442 #endif /* VERBOSE */
00443
00444 end subroutine PSMILe_Locations_gauss2