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