00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_global_search_real (comp_info, &
00012 control, len_cpl, &
00013 var_id, grid_valid_shape, &
00014 search, tgt_grid, &
00015 neighbors_3d, nloc, num_neigh, &
00016 extra_search, &
00017 interpolation_methods, &
00018 interpolation_search, n_intmethods, &
00019 send_index, &
00020 mask_available, use_mask, use_how, &
00021 grid_type, ierror)
00022
00023
00024
00025 use PRISM
00026
00027 use psmile_grid, only : psmile_transform_extent_cyclic, &
00028 psmile_transform_extent_back, &
00029 max_num_trans_parts, &
00030 code_no_trans, &
00031 common_grid_range
00032
00033 use PSMILe, dummy_interface => PSMILe_global_search_real
00034
00035 implicit none
00036
00037
00038
00039 Type (Enddef_comp), Intent (In) :: comp_info
00040
00041
00042
00043
00044 Type (Enddef_search), Intent (InOut) :: search
00045
00046
00047
00048 Type (real_vector), Intent (In) :: tgt_grid (ndim_3d)
00049
00050
00051
00052 Integer, Intent (In) :: control (2, ndim_3d, search%search_data%npart)
00053
00054
00055
00056 Integer, Intent (In) :: len_cpl (search%search_data%npart)
00057
00058
00059
00060 Integer, Intent (In) :: var_id
00061
00062
00063
00064 Integer, Intent (In) :: grid_valid_shape (2, ndim_3d)
00065
00066
00067
00068 Integer, Intent (In) :: send_index
00069
00070
00071
00072 Integer, Intent (In) :: nloc
00073
00074
00075
00076 Integer, Intent (In) :: num_neigh
00077
00078
00079
00080
00081
00082 Integer, Intent (In) :: n_intmethods
00083
00084
00085
00086
00087 Integer, Intent (In) :: interpolation_methods (n_intmethods)
00088
00089
00090
00091 Logical, Intent (In) :: interpolation_search(n_intmethods)
00092
00093
00094
00095 Logical, Intent (In) :: mask_available
00096
00097
00098 Logical, Intent (In) :: use_mask
00099
00100 Integer, Intent (In) :: use_how
00101
00102
00103
00104 Integer, Intent (In) :: grid_type
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114 Integer, Intent (InOut) :: neighbors_3d (ndim_3d, nloc,
00115 num_neigh)
00116
00117
00118
00119 Type (Extra_search_info), Intent (InOut) :: extra_search
00120
00121
00122
00123
00124
00125
00126 Integer, Intent (Out) :: ierror
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137 Real, Parameter :: eps = 0.0d0
00138
00139 Integer, Parameter :: nc_reg = 2
00140 Integer, Parameter :: masked_out = 0
00141
00142
00143
00144
00145
00146
00147 Integer :: comp_id, method_id
00148 Type (Method), Pointer :: mp
00149
00150 Type (Corner_Block), Pointer :: corner_pointer
00151
00152
00153
00154 Integer :: grid_id, global_grid_id
00155 Type(Grid), Pointer :: grid_info
00156
00157 Type (integer_vector), Pointer :: indices_req (:)
00158 Type (integer_vector), Pointer :: required (:)
00159 Integer, Pointer :: len_req (:)
00160
00161
00162
00163 Integer :: rank, igrid
00164
00165
00166
00167
00168 Integer, allocatable :: igrid_to_dest_comp(:)
00169
00170
00171
00172 Integer :: ipart, j, n
00173
00174
00175
00176 Logical :: mask_changed
00177 Integer :: nextra_prev, nreq, nx
00178
00179 Integer, Allocatable :: n_faces (:)
00180 Logical, Allocatable :: send_mask (:)
00181 Real, Allocatable :: faces (:, :, :)
00182
00183
00184
00185 Integer :: n_trans, n_int
00186
00187 Integer :: tr_codes (max_num_trans_parts)
00188 Integer :: found (2*max_num_trans_parts)
00189 Real (PSMILe_float_kind) :: face_range (2, ndim_3d)
00190 Real (PSMILe_float_kind) :: dinter_trans (2, ndim_3d)
00191 Real (PSMILe_float_kind) :: face_range_transformed (2, ndim_3d,
00192 max_num_trans_parts),
00193 dinter (2, ndim_3d, 2*max_num_trans_parts)
00194
00195
00196
00197 Integer :: n_found, n_liste
00198
00199
00200
00201 Type (Send_information), Pointer :: send_info
00202 Type (enddef_msg_extra) :: extra_msg
00203 Integer :: extra_msg_buf (msg_extra_size)
00204 Integer :: answer (msg_extra_size)
00205 Logical :: virtual_cell_available
00206
00207 Integer :: interp
00208 Integer :: i, ip, ipi, ipi_prev, ireq
00209 Integer :: len_ibuf, len_item, ndibuf, n_send
00210 Integer :: len_rbuf, len_rtem, ndrbuf, nd_len
00211 Integer, Pointer :: len_nsend (:, :)
00212 Integer, Pointer :: old_len (:, :)
00213 Integer, Allocatable :: ibuf (:)
00214 Integer, Allocatable :: srcloc_ind (:, :)
00215 Integer, Pointer :: virtual_ind (:)
00216 Integer, Target :: dummy_virtual_ind (1)
00217 Real, Pointer :: buf (:)
00218
00219
00220
00221 Integer :: dest, index, dest_comp, sender
00222 Integer :: nrecv, nprev_recv, n_recv_req
00223 Integer :: recv_req
00224 Integer :: save_lreq (2:3)
00225 Integer :: status (MPI_STATUS_SIZE)
00226
00227
00228
00229 Integer, parameter :: nerrp = 3
00230 Integer :: ierrp (nerrp)
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255 Character(len=len_cvs_string), save :: mycvs =
00256 '$Id: psmile_global_search_real.F90 3174 2011-05-04 09:41:43Z hanke $'
00257
00258
00259
00260
00261
00262 #ifdef VERBOSE
00263 print 9990, trim(ch_id), var_id, search%msg_intersections%field_info%tgt_method_id, &
00264 search%sender
00265
00266 call psmile_flushstd
00267 #endif /* VERBOSE */
00268
00269 ierror = 0
00270
00271 ndrbuf = 0
00272 ndibuf = 0
00273 n_recv_req = 0
00274
00275 method_id = Fields(var_id)%method_id
00276 mp => Methods(method_id)
00277
00278 grid_id = mp%grid_id
00279 grid_info => Grids(grid_id)
00280
00281 global_grid_id = grid_info%global_grid_id
00282
00283
00284
00285 send_info => mp%send_infos_coupler(send_index)
00286 corner_pointer => grid_info%corner_pointer
00287
00288 comp_id = comp_info%comp_id
00289
00290 indices_req => extra_search%indices_req
00291 required => extra_search%required
00292 len_req => extra_search%len_req
00293 nreq = extra_search%n_req
00294
00295 virtual_cell_available = Associated (send_info%virtual)
00296
00297 #ifdef PRISM_ASSERTION
00298 if (grid_type /= grid_info%grid_type) then
00299
00300 print *, 'grid_type', grid_type, grid_info%grid_type
00301 call psmile_assert ( __FILE__, __LINE__, &
00302 'inconsistent grid types')
00303 endif
00304 #endif
00305
00306
00307
00308 rank = Comps(comp_id)%rank
00309
00310 if (count(interpolation_search(:)) > 1) then
00311 ierrp (1) = count(interpolation_search(:))
00312
00313 ierror = PRISM_Error_Internal
00314
00315 call psmile_error ( ierror, 'Global search currently supported only for one interpolation method', &
00316 ierrp, 1, __FILE__, __LINE__ )
00317
00318 return
00319 endif
00320
00321 do interp = 1, n_intmethods
00322 if (interpolation_search(interp) ) exit
00323 end do
00324
00325 if (interp > n_intmethods) then
00326 ierror = PRISM_Error_Internal
00327 ierrp (1) = n_intmethods
00328
00329 call psmile_error ( ierror, 'No Global search for interpolation methods', &
00330 ierrp, 1, __FILE__, __LINE__ )
00331 return
00332 endif
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344 Allocate (faces (2, ndim_3d, nreq), STAT = ierror)
00345
00346 if ( ierror > 0 ) then
00347 ierrp (1) = ierror
00348 ierrp (2) = 2 * ndim_3d * nreq
00349
00350 ierror = PRISM_Error_Alloc
00351 call psmile_error ( ierror, 'faces', &
00352 ierrp, 2, __FILE__, __LINE__ )
00353 return
00354 endif
00355
00356 Allocate (n_faces (nreq), STAT = ierror)
00357
00358 if ( ierror > 0 ) then
00359 ierrp (1) = ierror
00360 ierrp (2) = nreq
00361
00362 ierror = PRISM_Error_Alloc
00363 call psmile_error ( ierror, 'n_faces', &
00364 ierrp, 2, __FILE__, __LINE__ )
00365 return
00366 endif
00367
00368 Allocate (send_mask (nreq), STAT = ierror)
00369
00370 if ( ierror > 0 ) then
00371 ierrp (1) = ierror
00372 ierrp (2) = nreq
00373
00374 ierror = PRISM_Error_Alloc
00375 call psmile_error ( ierror, 'send_mask', &
00376 ierrp, 2, __FILE__, __LINE__ )
00377 return
00378 endif
00379
00380 mask_changed = .true.
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401 nx = nreq
00402
00403 len_item = ndim_3d + 2
00404
00405 select case (grid_type)
00406
00407 case (PRISM_Reglonlatvrt)
00408 len_rtem = ndim_3d * (nc_reg + 1)
00409
00410 call psmile_get_faces_3d_reg_real ( &
00411 search, extra_search, &
00412 corner_pointer%corners_real(1)%vector, &
00413 corner_pointer%corners_real(2)%vector, &
00414 corner_pointer%corners_real(3)%vector, &
00415 corner_pointer%corner_shape, &
00416 corner_pointer%nbr_corners, &
00417 grid_valid_shape, &
00418 neighbors_3d, nloc, num_neigh, &
00419 faces, n_faces, nreq, ierror)
00420
00421 case (PRISM_Irrlonlat_regvrt)
00422 len_rtem = ndim_3d + &
00423 (corner_pointer%nbr_corners/nc_reg)*ndim_2d + ndim_2d
00424
00425 call psmile_get_faces_irreg2_real ( &
00426 search, extra_search, &
00427 corner_pointer%corners_real(1)%vector, &
00428 corner_pointer%corners_real(2)%vector, &
00429 corner_pointer%corners_real(3)%vector, &
00430 corner_pointer%corner_shape, &
00431 corner_pointer%nbr_corners, &
00432 grid_valid_shape, &
00433 neighbors_3d, nloc, num_neigh, &
00434 faces, n_faces, nreq, ierror)
00435
00436 case (PRISM_Irrlonlatvrt)
00437 len_rtem = ndim_3d + corner_pointer%nbr_corners*ndim_3d
00438
00439 call psmile_get_faces_3d_real ( &
00440 search, extra_search, &
00441 corner_pointer%corners_real(1)%vector, &
00442 corner_pointer%corners_real(2)%vector, &
00443 corner_pointer%corners_real(3)%vector, &
00444 corner_pointer%corner_shape, &
00445 corner_pointer%nbr_corners, &
00446 grid_valid_shape, &
00447 neighbors_3d, nloc, num_neigh, &
00448 faces, n_faces, nreq, ierror)
00449
00450 case (PRISM_Gaussreduced_regvrt)
00451 if (virtual_cell_available) len_item = ndim_3d + 3
00452 len_rtem = ndim_3d * (nc_reg + 1)
00453
00454 call psmile_get_faces_gauss2_real ( &
00455 search, extra_search, &
00456 corner_pointer%corners_real(1)%vector, &
00457 corner_pointer%corners_real(2)%vector, &
00458 corner_pointer%corners_real(3)%vector, &
00459 corner_pointer%corner_shape, &
00460 corner_pointer%nbr_corners, &
00461 grid_valid_shape, &
00462 neighbors_3d, nloc, num_neigh, &
00463 faces, nreq, ierror)
00464
00465 case default
00466 ierrp (1) = grid_type
00467 ierror = PRISM_Error_Internal
00468
00469 call psmile_error ( ierror, 'unsupported grid type', &
00470 ierrp, 1, __FILE__, __LINE__ )
00471 return
00472
00473 end select
00474
00475 if (ierror > 0) return
00476
00477
00478
00479 face_range (1, 1) = minval (faces (1, 1, 1:nreq))
00480 face_range (2, 1) = maxval (faces (2, 1, 1:nreq))
00481 face_range (1, 2) = minval (faces (1, 2, 1:nreq))
00482 face_range (2, 2) = maxval (faces (2, 2, 1:nreq))
00483 face_range (1, 3) = minval (faces (1, 3, 1:nreq))
00484 face_range (2, 3) = maxval (faces (2, 3, 1:nreq))
00485
00486 call psmile_transform_extent_cyclic (grid_info%grid_type, face_range, &
00487 face_range_transformed, tr_codes, &
00488 n_trans, ierror)
00489 if (ierror > 0) return
00490
00491
00492
00493 nd_len = 8
00494 Allocate (len_nsend (search%search_data%npart, nd_len), stat = ierror)
00495
00496 if ( ierror > 0 ) then
00497 ierrp (1) = ierror
00498 ierrp (2) = search%search_data%npart * nd_len
00499
00500 ierror = PRISM_Error_Alloc
00501 call psmile_error ( ierror, 'len_nsend', &
00502 ierrp, 2, __FILE__, __LINE__ )
00503 return
00504 endif
00505
00506
00507
00508 Allocate (igrid_to_dest_comp(SUM(comp_info%Number_of_grids_vector)), stat = ierror)
00509
00510 if ( ierror > 0 ) then
00511 ierrp (1) = ierror
00512 ierrp (2) = SUM(comp_info%Number_of_grids_vector)
00513
00514 ierror = PRISM_Error_Alloc
00515 call psmile_error ( ierror, 'igrid_to_dest_comp', &
00516 ierrp, 2, __FILE__, __LINE__ )
00517 return
00518 endif
00519
00520 igrid = 0
00521 do dest_comp = 1, comp_info%size
00522 do j = 1, comp_info%Number_of_grids_vector(dest_comp)
00523 igrid = igrid + 1
00524 igrid_to_dest_comp(igrid) = dest_comp
00525 enddo
00526 enddo
00527
00528
00529
00530
00531
00532
00533 do igrid = 1, SUM(comp_info%Number_of_grids_vector)
00534
00535 dest_comp = igrid_to_dest_comp(igrid)
00536
00537 if (mask_changed) then
00538 send_mask (:) = .false.
00539 mask_changed = .false.
00540 endif
00541
00542
00543 if (global_grid_id /= comp_info%all_extent_infos(igrid)%global_grid_id) cycle
00544
00545
00546
00547 if (rank == dest_comp-1 .and. &
00548 comp_info%all_extent_infos(igrid)%local_grid_id == grid_id) cycle
00549
00550
00551
00552
00553
00554
00555 do i = 1, n_trans
00556
00557 dinter (1, :, i) = max (face_range_transformed(1,:,i), &
00558 comp_info%all_extent_infos(igrid)%extent(1,:))
00559 dinter (2, :, i) = min (face_range_transformed(2,:,i), &
00560 comp_info%all_extent_infos(igrid)%extent(2,:))
00561
00562
00563
00564
00565
00566 if ( face_range_transformed(1,1,i) == common_grid_range(1,1) .and. &
00567 .not. face_range_transformed(2,1,i) == common_grid_range(2,1)) then
00568
00569 dinter (:, :, i+n_trans) = dinter (:, :, i)
00570 dinter (:, 1, i+n_trans) = common_grid_range(1,1)
00571
00572 else if (.not. face_range_transformed(1,1,i) == common_grid_range(1,1) .and. &
00573 face_range_transformed(2,1,i) == common_grid_range(2,1)) then
00574
00575 dinter (:, :, i+n_trans) = dinter (:, :, i)
00576 dinter (:, 1, i+n_trans) = common_grid_range(2,1)
00577
00578 else
00579
00580 dinter (1,:,i+n_trans) = 0
00581 dinter (2,:,i+n_trans) = -1
00582 endif
00583 end do
00584
00585 n_int = 0
00586
00587 do i = 1, 2 * n_trans
00588 if (minval(dinter (2,:,i) - dinter (1,:,i)) >= 0.0d0) then
00589 n_int = n_int + 1
00590 found (n_int) = i
00591 endif
00592 end do
00593
00594 if (n_int == 0) cycle
00595
00596
00597
00598 do i = 1, n_int
00599 if (tr_codes(found(i)) /= code_no_trans) then
00600 call psmile_transform_extent_back (tr_codes(found(i)), &
00601 dinter(:, :, found(i)), dinter_trans, 1, ierror)
00602 if ( ierror /= 0 ) return
00603
00604 dinter (:, :, found(i)) = dinter_trans
00605 endif
00606 end do
00607
00608
00609
00610
00611
00612
00613 do i = 1, n_int
00614 #ifdef DEBUGX
00615 print *, 'extent controlled', i, n_int, dinter(:,:,found(i))
00616 #endif
00617
00618 do j = 1, nx
00619 send_mask (j) = send_mask (j) .or. &
00620 (faces(2,1,j) >= dinter(1,1,found(i)) .and. &
00621 faces(1,1,j) <= dinter(2,1,found(i)) .and. &
00622 faces(2,2,j) >= dinter(1,2,found(i)) .and. &
00623 faces(1,2,j) <= dinter(2,2,found(i)) .and. &
00624 faces(2,3,j) >= dinter(1,3,found(i)) .and. &
00625 faces(1,3,j) <= dinter(2,3,found(i)))
00626 #ifdef DEBUGX
00627 print *, j, faces(:,:,j)
00628 #endif
00629 end do
00630 end do
00631
00632
00633
00634 n_send = COUNT(send_mask)
00635 if (n_send == 0) cycle
00636
00637 mask_changed = .true.
00638
00639
00640
00641 dest = comp_info%psmile_ranks(dest_comp)
00642
00643 n_recv_req = n_recv_req + 1
00644 if (n_recv_req > nd_len) then
00645 nd_len = nd_len + 16
00646
00647
00648
00649 old_len => len_nsend
00650 Allocate (len_nsend (search%search_data%npart, nd_len), stat = ierror)
00651
00652 if ( ierror > 0 ) then
00653 ierrp (1) = ierror
00654 ierrp (2) = search%search_data%npart * nd_len
00655
00656 ierror = PRISM_Error_Alloc
00657 call psmile_error ( ierror, 'len_nsend', &
00658 ierrp, 2, __FILE__, __LINE__ )
00659 return
00660 endif
00661
00662 len_nsend (:, 1:n_recv_req-1) = old_len (:, 1:n_recv_req-1)
00663
00664 Deallocate (old_len)
00665 endif
00666
00667 #ifdef VERBOSE
00668 print 9970, trim(ch_id), rank, dest_comp, dest, n_send
00669
00670 call psmile_flushstd
00671 #endif /* VERBOSE */
00672
00673
00674
00675
00676
00677
00678
00679
00680
00681
00682
00683
00684
00685
00686
00687
00688
00689
00690 len_rbuf = n_send * len_rtem
00691 len_ibuf = n_send * len_item
00692
00693 if (len_rbuf > ndrbuf) then
00694 if (ndrbuf > 0) then
00695 Deallocate (buf)
00696 endif
00697
00698 ndrbuf = len_rbuf
00699 Allocate (buf(ndrbuf), STAT = ierror)
00700
00701 if ( ierror > 0 ) then
00702 ierrp (1) = ierror
00703 ierrp (2) = ndrbuf
00704
00705 ierror = PRISM_Error_Alloc
00706 call psmile_error ( ierror, 'buf', &
00707 ierrp, 2, __FILE__, __LINE__ )
00708 return
00709 endif
00710 endif
00711
00712 if (len_ibuf > ndibuf) then
00713 if (ndibuf > 0) then
00714 Deallocate (ibuf)
00715 endif
00716
00717 ndibuf = len_ibuf
00718 Allocate (ibuf(ndibuf), STAT = ierror)
00719
00720 if ( ierror > 0 ) then
00721 ierrp (1) = ierror
00722 ierrp (2) = ndibuf
00723
00724 ierror = PRISM_Error_Alloc
00725 call psmile_error ( ierror, 'ibuf', &
00726 ierrp, 2, __FILE__, __LINE__ )
00727 return
00728 endif
00729 endif
00730
00731
00732
00733
00734 Allocate (srcloc_ind(ndim_3d, n_send), STAT = ierror)
00735
00736 if ( ierror > 0 ) then
00737 ierrp (1) = ierror
00738 ierrp (2) = ndim_3d * n_send
00739
00740 ierror = PRISM_Error_Alloc
00741 call psmile_error ( ierror, 'srcloc_ind', &
00742 ierrp, 2, __FILE__, __LINE__ )
00743 return
00744 endif
00745
00746 select case (send_info%nvec)
00747 case (1)
00748 call psmile_get_face_ind_3d (search, extra_search, &
00749 send_info, len_cpl, &
00750 send_mask, nreq, srcloc_ind, n_send, &
00751 ierror)
00752
00753 case (2)
00754 call psmile_get_face_ind_21d (search, extra_search, &
00755 send_info, len_cpl, &
00756 send_mask, nreq, srcloc_ind, n_send, &
00757 ierror)
00758
00759 case (ndim_3d)
00760 call psmile_get_face_ind_reg (search, extra_search, &
00761 send_info, len_cpl, &
00762 send_mask, nreq, srcloc_ind, n_send, &
00763 ierror)
00764 end select
00765
00766 if (ierror > 0) return
00767
00768 if (grid_type == PRISM_Gaussreduced_regvrt) then
00769 #ifdef DEBUG
00770
00771
00772 print *, '### psmile_global_search_real.F90 raus damit !!', send_info%nvec
00773 #endif
00774 do i = 1, n_send
00775 srcloc_ind (1, i) = abs (srcloc_ind(1,i))
00776 end do
00777 endif
00778
00779 if (virtual_cell_available) then
00780 Allocate (virtual_ind(n_send), STAT = ierror)
00781
00782 if ( ierror > 0 ) then
00783 ierrp (1) = ierror
00784 ierrp (2) = n_send
00785
00786 ierror = PRISM_Error_Alloc
00787 call psmile_error ( ierror, 'virtual_ind', &
00788 ierrp, 2, __FILE__, __LINE__ )
00789 return
00790 endif
00791
00792 call psmile_get_faces_virtual_ind (search, extra_search, &
00793 send_info, len_cpl, &
00794 send_mask, nreq, virtual_ind, n_send, &
00795 ierror)
00796 if (ierror > 0) return
00797
00798 else
00799
00800 virtual_ind => dummy_virtual_ind
00801
00802 endif
00803
00804
00805
00806
00807
00808
00809
00810 nextra_prev = 0
00811 ipi = 0
00812 ip = 0
00813
00814 select case (grid_type)
00815
00816 case (PRISM_Reglonlatvrt)
00817 do ipart = 1, search%search_data%npart
00818
00819 ipi_prev = ipi
00820
00821 if (len_req (ipart) > 0) then
00822
00823 call psmile_store_faces_3d_reg_real ( &
00824 indices_req(ipart)%vector, &
00825 required (ipart)%vector, &
00826 len_req (ipart), &
00827 tgt_grid (1)%vector, &
00828 tgt_grid (2)%vector, &
00829 tgt_grid (3)%vector, nloc, &
00830 corner_pointer%corners_real(1)%vector, &
00831 corner_pointer%corners_real(2)%vector, &
00832 corner_pointer%corners_real(3)%vector, &
00833 corner_pointer%corner_shape, &
00834 corner_pointer%nbr_corners, &
00835 grid_valid_shape, &
00836 send_mask(nextra_prev+1:), &
00837 srcloc_ind, &
00838 ibuf, len_item, n_send, ipi, &
00839 buf, len_rtem, n_send, ip, ierror)
00840
00841 nextra_prev = nextra_prev + len_req (ipart)
00842 endif
00843
00844 len_nsend (ipart, n_recv_req) = ipi - ipi_prev
00845 end do
00846
00847 #ifdef PRISM_ASSERTION
00848 if (ipi /= n_send) then
00849 print *, 'ipi, n_send', ipi, n_send
00850 call psmile_assert ( __FILE__, __LINE__, &
00851 'inconsistent length for pack buffer ibuf')
00852 endif
00853
00854 if (ip /= n_send) then
00855 print *, 'ip, n_send', ip, n_send
00856 call psmile_assert ( __FILE__, __LINE__, &
00857 'inconsistent length for pack buffer buf')
00858 endif
00859 #endif
00860
00861 case (PRISM_Irrlonlat_regvrt)
00862
00863 do ipart = 1, search%search_data%npart
00864
00865 ipi_prev = ipi
00866
00867 if (len_req (ipart) > 0) then
00868
00869 call psmile_store_faces_irreg2_real ( &
00870 indices_req(ipart)%vector, &
00871 required (ipart)%vector, &
00872 len_req (ipart), &
00873 tgt_grid (1)%vector, &
00874 tgt_grid (2)%vector, &
00875 tgt_grid (3)%vector, nloc, &
00876 corner_pointer%corners_real(1)%vector, &
00877 corner_pointer%corners_real(2)%vector, &
00878 corner_pointer%corners_real(3)%vector, &
00879 corner_pointer%corner_shape, &
00880 corner_pointer%nbr_corners, &
00881 grid_valid_shape, &
00882 send_mask(nextra_prev+1:), &
00883 srcloc_ind, &
00884 ibuf, len_item, n_send, ipi, &
00885 buf, len_rtem, n_send, ip, ierror)
00886
00887 nextra_prev = nextra_prev + len_req (ipart)
00888 endif
00889
00890 len_nsend (ipart, n_recv_req) = ipi - ipi_prev
00891 end do
00892
00893 #ifdef PRISM_ASSERTION
00894 if (ipi /= n_send) then
00895 print *, 'ipi, n_send', ipi, n_send
00896 call psmile_assert ( __FILE__, __LINE__, &
00897 'inconsistent length for pack buffer ibuf')
00898 endif
00899
00900 if (ip /= n_send) then
00901 print *, 'ip, n_send', ip, n_send
00902 call psmile_assert ( __FILE__, __LINE__, &
00903 'inconsistent length for pack buffer buf')
00904 endif
00905 #endif
00906
00907 case (PRISM_Gaussreduced_regvrt)
00908
00909 do ipart = 1, search%search_data%npart
00910
00911 ipi_prev = ipi
00912
00913 if (len_req (ipart) > 0) then
00914
00915 call psmile_store_faces_gauss2_real ( &
00916 indices_req(ipart)%vector, &
00917 required (ipart)%vector, &
00918 len_req (ipart), &
00919 tgt_grid (1)%vector, &
00920 tgt_grid (2)%vector, &
00921 tgt_grid (3)%vector, nloc, &
00922 corner_pointer%corners_real(1)%vector, &
00923 corner_pointer%corners_real(2)%vector, &
00924 corner_pointer%corners_real(3)%vector, &
00925 corner_pointer%corner_shape, &
00926 corner_pointer%nbr_corners, &
00927 grid_id, grid_valid_shape, &
00928 send_mask(nextra_prev+1:), &
00929 srcloc_ind, virtual_ind, &
00930 virtual_cell_available, &
00931 ibuf, len_item, n_send, ipi, &
00932 buf, len_rtem, n_send, ip, ierror)
00933
00934 nextra_prev = nextra_prev + len_req (ipart)
00935 endif
00936
00937 len_nsend (ipart, n_recv_req) = ipi - ipi_prev
00938 end do
00939
00940 #ifdef PRISM_ASSERTION
00941 if (ipi /= n_send) then
00942 print *, 'ipi, n_send', ipi, n_send
00943 call psmile_assert ( __FILE__, __LINE__, &
00944 'inconsistent length for pack buffer ibuf')
00945 endif
00946
00947 if (ip /= n_send) then
00948 print *, 'ip, n_send', ip, n_send
00949 call psmile_assert ( __FILE__, __LINE__, &
00950 'inconsistent length for pack buffer buf')
00951 endif
00952 #endif
00953
00954 case (PRISM_Irrlonlatvrt)
00955 do ipart = 1, search%search_data%npart
00956
00957 ipi_prev = ipi
00958
00959 if (len_req (ipart) > 0) then
00960
00961 call psmile_store_faces_3d_real ( &
00962 indices_req(ipart)%vector, &
00963 required (ipart)%vector, &
00964 len_req (ipart), &
00965 tgt_grid (1)%vector, &
00966 tgt_grid (2)%vector, &
00967 tgt_grid (3)%vector, nloc, &
00968 corner_pointer%corners_real(1)%vector, &
00969 corner_pointer%corners_real(2)%vector, &
00970 corner_pointer%corners_real(3)%vector, &
00971 corner_pointer%corner_shape, &
00972 corner_pointer%nbr_corners, &
00973 grid_valid_shape, &
00974 send_mask(nextra_prev+1:), &
00975 srcloc_ind, &
00976 ibuf, len_item, n_send, ipi, &
00977 buf, len_rtem, n_send, ip, ierror)
00978
00979 nextra_prev = nextra_prev + len_req (ipart)
00980 endif
00981
00982 len_nsend (ipart, n_recv_req) = ipi - ipi_prev
00983 end do
00984
00985 case default
00986 ierrp (1) = grid_type
00987 ierror = PRISM_Error_Internal
00988
00989 call psmile_error ( ierror, 'unsupported grid type', &
00990 ierrp, 1, __FILE__, __LINE__ )
00991 return
00992
00993 end select
00994
00995 Deallocate (srcloc_ind)
00996 if (virtual_cell_available) Deallocate (virtual_ind)
00997
00998 #ifdef TODO
00999
01000 if (dest_comp == rank) then
01001
01002 search_global%sender = dest
01003 search_global%msg_extra = extra_msg
01004
01005 call psmile_search_donor_extra (search_global, tol, ierror)
01006 if (ierror > 0) return
01007
01008 endif
01009 #endif /* TODO */
01010
01011
01012
01013
01014
01015
01016 call psmile_init_enddef_msg_extra(extra_msg)
01017
01018 extra_msg%reqest_type = interpolation_methods (interp)
01019 extra_msg%datatype = PRISM_REAL
01020 extra_msg%len_int_data = len_ibuf
01021 extra_msg%len_coord_data = len_rbuf
01022 extra_msg%global_comp_id = comp_info%global_comp_id
01023 extra_msg%transi_out_id = search%msg_intersections%field_info%transient_out_id
01024
01025 extra_msg%num_volumes = n_send
01026 extra_msg%num_int_per_vol = len_item
01027 extra_msg%num_items_per_coord = len_rtem
01028
01029 extra_msg%idx_req = n_recv_req
01030 extra_msg%num_neigh = num_neigh
01031
01032 extra_msg%local_grid_id = comp_info%all_extent_infos(igrid)%local_grid_id
01033
01034 #ifdef DEBUGX
01035 print *, 'grid shape', grid_valid_shape
01036 do i = 1, n_send
01037 print *, ibuf ((i-1)*len_item+1:i*len_item)
01038 end do
01039 #endif
01040
01041 if (Associated (grid_info%partition)) then
01042
01043
01044
01045 extra_msg%partition_avail = .true.
01046 extra_msg%partition(:) = grid_info%partition (1, 1:3)
01047
01048
01049
01050 if (grid_info%grid_type == PRISM_Gaussreduced_regvrt) then
01051
01052 call psmile_trans_loc2glob_gauss2 (grid_id, &
01053 ibuf, len_item, n_send, ierror)
01054
01055 else
01056
01057 call psmile_trans_loc2glob_3d (grid_id, &
01058 ibuf, len_item, n_send, ierror)
01059
01060 endif
01061
01062 if (ierror > 0) return
01063
01064 else
01065
01066
01067
01068 extra_msg%partition_avail = .false.
01069 endif
01070
01071 #ifdef DEBUG
01072 print *, ' Sending ', exttag, ' to destination ', dest
01073 #endif
01074 call psmile_pack_msg_extra (extra_msg, extra_msg_buf)
01075
01076 call psmile_bsend (extra_msg_buf, msg_extra_size, MPI_INTEGER, &
01077 dest, exttag, comm_psmile, ierror)
01078 if (ierror /= MPI_SUCCESS) then
01079 ierrp (1) = ierror
01080 ierrp (2) = dest
01081 ierrp (3) = exttag
01082
01083 ierror = PRISM_Error_Send
01084
01085 call psmile_error (ierror, 'psmile_bsend(msg)', &
01086 ierrp, 3, __FILE__, __LINE__ )
01087 return
01088 endif
01089
01090 #ifdef DEBUG
01091 print *, ' Sending ', exttag, ' to destination ', dest
01092 #endif
01093 call psmile_bsend (ibuf, len_ibuf, MPI_INTEGER, &
01094 dest, exttag, comm_psmile, ierror)
01095 if (ierror /= MPI_SUCCESS) then
01096 ierrp (1) = ierror
01097 ierrp (2) = dest
01098 ierrp (3) = exttag
01099
01100 ierror = PRISM_Error_Send
01101
01102 call psmile_error (ierror, 'psmile_bsend(ibuf)', &
01103 ierrp, 3, __FILE__, __LINE__ )
01104 return
01105 endif
01106
01107 #ifdef DEBUG
01108 print *, ' Sending ', exttag, ' to destination ', dest
01109 #endif
01110 call psmile_bsend (buf, len_rbuf, MPI_REAL, &
01111 dest, exttag, comm_psmile, ierror)
01112 if (ierror /= MPI_SUCCESS) then
01113 ierrp (1) = ierror
01114 ierrp (2) = dest
01115 ierrp (3) = exttag
01116
01117 ierror = PRISM_Error_Send
01118
01119 call psmile_error (ierror, 'psmile_bsend(buf)', &
01120 ierrp, 3, __FILE__, __LINE__ )
01121 return
01122 endif
01123
01124
01125
01126
01127 if (n_recv_req == 1) then
01128
01129
01130
01131
01132
01133 call MPI_Irecv (answer, msg_extra_size, MPI_INTEGER, MPI_ANY_SOURCE, &
01134 rexttag, comm_psmile, recv_req, ierror)
01135 if (ierror /= MPI_SUCCESS) then
01136
01137 ierrp (1) = ierror
01138 ierrp (2) = dest
01139 ierrp (3) = rexttag
01140
01141 ierror = PRISM_Error_Recv
01142
01143 call psmile_error ( ierror, 'MPI_Irecv', &
01144 ierrp, 3, __FILE__, __LINE__ )
01145 return
01146
01147 endif
01148 endif
01149
01150 end do
01151
01152
01153
01154 Deallocate (send_mask, n_faces, faces, igrid_to_dest_comp)
01155
01156 if (ndrbuf > 0) Deallocate (buf)
01157 if (ndibuf > 0) Deallocate (ibuf)
01158
01159
01160
01161
01162
01163
01164
01165
01166
01167
01168
01169
01170
01171
01172 nprev_recv = 0
01173 nrecv = 0
01174
01175 if (n_recv_req > 0) then
01176 Allocate (extra_search%real_bufs(n_recv_req), &
01177 send_info%sender_global(n_recv_req), &
01178 send_info%len_sent(n_recv_req), &
01179 send_info%msg_id(n_recv_req), &
01180 stat = ierror)
01181 if ( ierror > 0 ) then
01182 ierrp (1) = ierror
01183 ierrp (2) = n_recv_req
01184
01185 ierror = PRISM_Error_Alloc
01186 call psmile_error ( ierror, 'extra_search%real_bufs', &
01187 ierrp, 2, __FILE__, __LINE__ )
01188 return
01189 endif
01190 endif
01191
01192 save_lreq (2:3) = paction%lrequest (2:3)
01193
01194 paction%lrequest (2) = MPI_REQUEST_NULL
01195
01196 do n = 1, n_recv_req
01197 paction%lrequest (3) = recv_req
01198
01199 index = 0
01200 do while (index /= 3)
01201 #ifdef DEBUG
01202 print *, trim(ch_id), ': psmile_global_search_real: action%nreq, recv_req', &
01203 paction%nreq, recv_req
01204 call psmile_flushstd
01205 #endif
01206
01207 call MPI_Waitany (paction%nreq, paction%lrequest, index, status, ierror)
01208
01209 if ( ierror /= MPI_SUCCESS ) then
01210 ierrp (1) = ierror
01211 ierrp (2) = status (MPI_SOURCE)
01212 ierrp (3) = status (MPI_TAG)
01213
01214 ierror = PRISM_Error_MPI
01215
01216 call psmile_error ( ierror, 'MPI_Waitany', &
01217 ierrp, 3, __FILE__, __LINE__ )
01218 return
01219 endif
01220
01221 #ifdef PRISM_ASSERTION
01222 if (index == MPI_UNDEFINED) then
01223 call psmile_assert ( __FILE__, __LINE__, &
01224 'request list is empty')
01225 endif
01226 #endif
01227
01228 if (index /= 3) then
01229 call psmile_enddef_action (search, index, status, ierror)
01230 if (ierror > 0) return
01231 endif
01232 end do
01233
01234
01235
01236 sender = status (MPI_SOURCE)
01237 len_ibuf = answer (3)
01238 len_rbuf = answer (4)
01239
01240 #ifdef VERBOSE
01241 print 9960, trim(ch_id), sender, len_ibuf, len_rbuf
01242
01243 call psmile_flushstd
01244 #endif /* VERBOSE */
01245
01246 #ifdef PRISM_ASSERTION
01247 if (len_ibuf < 0) then
01248 print *, 'len_ibuf =', len_ibuf
01249 call psmile_assert (__FILE__, __LINE__, &
01250 "len_ibuf should be >= 0")
01251 endif
01252
01253 if (len_rbuf < 0) then
01254 print *, 'len_rbuf =', len_rbuf
01255 call psmile_assert (__FILE__, __LINE__, &
01256 "len_rbuf should be >= 0")
01257 endif
01258 #endif
01259
01260 if (len_ibuf > 0) then
01261 Allocate (ibuf (1:len_ibuf), stat = ierror)
01262
01263 if ( ierror > 0 ) then
01264 ierrp (1) = ierror
01265 ierrp (2) = len_ibuf
01266
01267 ierror = PRISM_Error_Alloc
01268 call psmile_error ( ierror, 'ibuf', &
01269 ierrp, 2, __FILE__, __LINE__ )
01270 return
01271 endif
01272
01273 call MPI_Recv (ibuf, len_ibuf, MPI_INTEGER, sender, &
01274 rexttag, comm_psmile, status, ierror)
01275 if (ierror /= MPI_SUCCESS) then
01276
01277 ierrp (1) = ierror
01278 ierrp (2) = sender
01279 ierrp (3) = rexttag
01280
01281 ierror = PRISM_Error_Recv
01282
01283 call psmile_error ( ierror, 'MPI_Recv(ibuf)', &
01284 ierrp, 3, __FILE__, __LINE__ )
01285 return
01286 endif
01287
01288 Allocate (buf (1:len_rbuf), stat = ierror)
01289
01290 if ( ierror > 0 ) then
01291 ierrp (1) = ierror
01292 ierrp (2) = len_rbuf
01293
01294 ierror = PRISM_Error_Alloc
01295 call psmile_error ( ierror, 'buf', &
01296 ierrp, 2, __FILE__, __LINE__ )
01297 return
01298 endif
01299
01300 call MPI_Recv (buf, len_rbuf, MPI_REAL, sender, &
01301 rexttag, comm_psmile, status, ierror)
01302 if (ierror /= MPI_SUCCESS) then
01303
01304 ierrp (1) = ierror
01305 ierrp (2) = sender
01306 ierrp (3) = rexttag
01307
01308 ierror = PRISM_Error_Recv
01309
01310 call psmile_error ( ierror, 'MPI_Recv(rbuf)', &
01311 ierrp, 3, __FILE__, __LINE__ )
01312 return
01313 endif
01314
01315
01316
01317
01318
01319
01320
01321
01322
01323 n_send = answer (7)
01324 n_liste = answer (8)
01325 n_found = answer (9)
01326 ireq = answer (15)
01327
01328 #ifdef PRISM_ASSERTION
01329 if (answer (1) /= interpolation_methods (interp)) then
01330 print *, 'answer(1)', answer (1), interpolation_methods (interp)
01331 call psmile_assert (__FILE__, __LINE__, &
01332 "incorrect message received ? interpolations doesn't fit")
01333 endif
01334
01335 if (ireq < 1 .or. ireq > n_recv_req) then
01336 print *, 'ireq, n_recv_req', ireq, n_recv_req
01337 call psmile_assert (__FILE__, __LINE__, &
01338 "ireq must be in range of 1:n_recv_req")
01339 endif
01340 #endif
01341
01342 nrecv = nrecv + 1
01343
01344 extra_search%real_bufs (nrecv)%vector => buf
01345 send_info%sender_global (nrecv) = sender
01346 send_info%len_sent (nrecv) = n_liste
01347 send_info%msg_id (nrecv) = answer (16)
01348
01349 if (nprev_recv > 0) then
01350
01351
01352
01353 do i = 2*n_send+1, 2*n_send+n_found
01354 if (ibuf(i) /= masked_out) then
01355 ibuf (i) = ibuf (i) + nprev_recv
01356 endif
01357 end do
01358 endif
01359
01360 nprev_recv = nprev_recv + n_liste
01361
01362
01363
01364 call psmile_add_points_found (grid_id, search, extra_search, &
01365 ibuf (1:n_send), ibuf (n_send+1:2*n_send), &
01366 n_send, len_nsend(1:search%search_data%npart,ireq), &
01367 ibuf (2*n_send+1:2*n_send+n_found), n_found, &
01368 neighbors_3d, nloc, num_neigh, &
01369 grid_valid_shape, use_how, ierror)
01370 if (ierror > 0) return
01371
01372 Deallocate (ibuf)
01373 endif
01374
01375
01376
01377 if (n < n_recv_req) then
01378
01379
01380
01381
01382
01383 call MPI_Irecv (answer, msg_extra_size, MPI_INTEGER, MPI_ANY_SOURCE, &
01384 rexttag, comm_psmile, recv_req, &
01385 ierror)
01386 if (ierror /= MPI_SUCCESS) then
01387
01388 ierrp (1) = ierror
01389 ierrp (2) = dest
01390 ierrp (3) = rexttag
01391
01392 ierror = PRISM_Error_Recv
01393
01394 call psmile_error ( ierror, 'MPI_Irecv', &
01395 ierrp, 3, __FILE__, __LINE__ )
01396 return
01397 endif
01398 endif
01399 end do
01400
01401
01402
01403 #ifdef PRISM_ASSERTION
01404 if (paction%lrequest (2) /= MPI_REQUEST_NULL .or. &
01405 paction%lrequest (3) /= MPI_REQUEST_NULL) then
01406 print *, 'request: ', paction%lrequest (2:3)
01407 call psmile_assert ( __FILE__, __LINE__, &
01408 'Illegal request stored')
01409
01410 endif
01411 #endif
01412
01413 paction%lrequest (2:3) = save_lreq (2:3)
01414
01415 Deallocate (len_nsend)
01416
01417
01418
01419 send_info%nrecv = nrecv
01420 send_info%num2recv = nprev_recv
01421
01422 #ifdef PRISM_ASSERTION
01423 if (nrecv > 0) then
01424 if (send_info%num2recv /= SUM (send_info%len_sent(1:nrecv))) then
01425 print *, 'num2recv', send_info%num2recv, &
01426 SUM (send_info%len_sent(1:nrecv)), nrecv
01427 call psmile_assert (__FILE__, __LINE__, &
01428 "inconsistent num2recv")
01429 endif
01430 endif
01431 #endif
01432
01433
01434
01435
01436
01437 #ifdef VERBOSE
01438 print 9980, trim(ch_id), ierror
01439
01440 call psmile_flushstd
01441 #endif /* VERBOSE */
01442
01443
01444
01445
01446 #ifdef VERBOSE
01447
01448 9990 format (1x, a, ': psmile_global_search_real: var_id', i3, &
01449 ' to ', i3, '(', i2, ')')
01450 9980 format (1x, a, ': psmile_global_search_real: eof ierror =', i3)
01451 9970 format (1x, a, ': psmile_global_search_real: send from', i3, &
01452 ' to', i3, '[', i3, '], n_send =', i6)
01453
01454 9960 format (1x, a, ': psmile_global_search_real: got rexttag message:', &
01455 ' sender ', i4, ', len_ibuf, len_rbuf', 2i8)
01456 9950 format (1x, a, ': psmile_global_search_real: before waitany :', &
01457 'nreq =', i4, ', recv_req ', i8)
01458 #endif /* VERBOSE */
01459
01460 #ifdef DEBUG
01461 #endif
01462
01463 end subroutine PSMILe_global_search_real