psmile_global_search_real.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PSMILe_global_search_real
00008 !
00009 ! !INTERFACE:
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 ! !USES:
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 ! !INPUT PARAMETERS:
00038 
00039       Type (Enddef_comp), Intent (In) :: comp_info
00040 !
00041 !     Info on the component in which the donor cells
00042 !     should be searched.
00043 !
00044       Type (Enddef_search), Intent (InOut) :: search
00045 
00046 !     Info's on coordinates to be searched
00047 !
00048       Type (real_vector), Intent (In) :: tgt_grid (ndim_3d)
00049 
00050 !     Coordinates to be searched
00051 
00052       Integer,            Intent (In) :: control (2, ndim_3d, search%npart)
00053 
00054 !     Index range found for partition "ipart".
00055 
00056       Integer,            Intent (In) :: len_cpl (search%npart)
00057 
00058 !     Number of points to be sent to the coupler for each partition.
00059 
00060       Integer,            Intent (In) :: var_id
00061 
00062 !     Handle to the grid function
00063 
00064       Integer,            Intent (In) :: grid_valid_shape (2, ndim_3d)
00065       
00066 !     Specifies the valid block shape for the "ndim_3d"-dimensional block
00067 
00068       Integer,            Intent (In) :: send_index
00069 
00070 !     Send index of coupler send for the method and coordinates
00071 !
00072       Integer,            Intent (In) :: nloc
00073 !
00074 !     Total  Number of locations to be transferred
00075 !
00076       Integer,            Intent (In) :: num_neigh
00077 !
00078 !     Last dimension of neighbors array "neighbors_3d" and
00079 !     number of neighbors to be searched.
00080 !     besser n_corners
00081 
00082       Integer,            Intent (In) :: n_intmethods
00083 !
00084 !     Number of interpolations methods used to interpolate the data.
00085 !     Dimension of vectors "interpolation_methods" and "interpolation_search".
00086 
00087       Integer,            Intent (In) :: interpolation_methods (n_intmethods)
00088 !
00089 !     Codes for the interpolation methods
00090 !
00091       Logical,            Intent (In) :: interpolation_search(n_intmethods)
00092 !
00093 !     Is global search required for corresponding interpolation method ?
00094 
00095       Logical,            Intent (In) :: mask_available
00096 !  !!!!! Fehlt da nicht noch ein "globales" mask_available !!!!
00097 !
00098       Logical,            Intent (In) :: use_mask
00099 
00100       Integer,            Intent (In) :: use_how
00101 
00102 !     Information about how to apply the mask
00103 
00104       Integer,            Intent (In) :: grid_type
00105 !
00106 !     Code for the type of grid generation
00107 !     PRISM_Reglonlatvrt            = Regular grid in all directions
00108 !     PRISM_Irrlonlat_regvrt        = Horizontal irregular, regular in the vertical
00109 !     PRISM_Irrlonlatvrt            = Irregular grid in all directions
00110 !     PRISM_Gaussreduced_regvrt     = Gaussian reduced grid with regular vertical
00111 !
00112 ! !INPUT/OUTPUT PARAMETERS:
00113 !
00114       Integer,            Intent (InOut) :: neighbors_3d (ndim_3d, nloc, 
00115                                                           num_neigh)
00116 
00117 !     Indices of neighbor locations (interpolation bases)
00118 !
00119       Type (Extra_search_info), Intent (InOut) :: extra_search
00120 !
00121 !     Number of locations where
00122 !       (*) required mask values were not true
00123 
00124 ! !OUTPUT PARAMETERS:
00125 !
00126       Integer,           Intent (Out) :: ierror
00127 
00128 !     Returns the error code of PSMILe_global_search_real;
00129 !             ierror = 0 : No error
00130 !             ierror > 0 : Severe error
00131 !
00132 ! !LOCAL PARAMETERS
00133 !
00134 !  NC_REG = Number of corners for regular directions
00135 !  masked_out = Hash value for a point which is masked out
00136 !
00137       Real, Parameter              :: eps = 0.0d0
00138 !
00139       Integer, Parameter           :: nc_reg = 2
00140       Integer, Parameter           :: masked_out = 0
00141 
00142 !
00143 ! !LOCAL VARIABLES
00144 !
00145 !     ... Method pointer
00146 !
00147       Integer                         :: comp_id, method_id
00148       Type (Method), Pointer          :: mp
00149 !
00150       Type (Corner_Block), Pointer    :: corner_pointer
00151 !
00152 !     ... for the grid
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 !     ... Info's on actual component
00162 !
00163       Integer                         :: rank, igrid
00164       Real (PSMILe_float_kind), Pointer :: extents (:, :, :)
00165 !
00166 !     ... lookup table
00167 !
00168       Integer, allocatable            :: igrid_to_dest_comp(:)
00169 !
00170 !     ... loop variables
00171 !
00172       Integer                         :: ipart, j, n
00173 !
00174 !     ... for data on faces to other processes
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 !     ... transformed faces
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 !     ... for points found
00196 !
00197       Integer                         :: n_found, n_liste
00198 !
00199 !     ... for send to partners
00200 !
00201       Type (Send_information), Pointer :: send_info
00202       Integer                         :: extra_msg (nd_msgextra)
00203       Integer                         :: answer (nd_msgextra)
00204       Logical                         :: virtual_cell_available
00205 !
00206       Integer                         :: interp
00207       Integer                         :: i, ip, ipi, ipi_prev, ireq
00208       Integer                         :: len_ibuf, len_item, ndibuf, n_send
00209       Integer                         :: len_rbuf, len_rtem, ndrbuf, nd_len
00210       Integer, Pointer                :: len_nsend (:, :)
00211       Integer, Pointer                :: old_len (:, :)
00212       Integer, Allocatable            :: ibuf (:)
00213       Integer, Allocatable            :: srcloc_ind (:, :)
00214       Integer, Pointer                :: virtual_ind (:)
00215       Integer, Target                 :: dummy_virtual_ind (1)
00216       Real, Pointer                   :: buf (:)
00217 !
00218 !     ... for communication
00219 !
00220       Integer                         :: dest, index, dest_comp, sender
00221       Integer                         :: nrecv, nprev_recv, n_recv_req 
00222       Integer                         :: recv_req
00223       Integer                         :: save_lreq (2:3)
00224       Integer                         :: status (MPI_STATUS_SIZE)
00225 !
00226 !     ... for error handling
00227 !
00228       Integer, parameter              :: nerrp = 3
00229       Integer                         :: ierrp (nerrp)
00230 !
00231 ! !DESCRIPTION:
00232 !
00233 ! Subroutine "PSMILe_global_search_real" performs a global search
00234 ! for the interpolation points not found in the local search.
00235 !
00236 ! TODO: tgt_grid statt coords uebergeben;
00237 !       kann man sich dann die unterschiedlichen Routine fuer
00238 !       Get_faces_..., store_faces__ sparen !?!?
00239 !
00240 ! !REVISION HISTORY:
00241 !
00242 !   Date      Programmer   Description
00243 ! ----------  ----------   -----------
00244 ! 11.11.04    H. Ritzdorf  created
00245 ! 06.04.10    M. Hanke     adjusted in order to support multiple
00246 !                          local grid id's for a single global grid id
00247 !
00248 !EOP
00249 !----------------------------------------------------------------------
00250 !
00251 !  $Id: psmile_global_search_real.F90 2988 2011-02-24 14:36:19Z hanke $
00252 !  $Author: hanke $
00253 !
00254    Character(len=len_cvs_string), save :: mycvs = 
00255        '$Id: psmile_global_search_real.F90 2988 2011-02-24 14:36:19Z hanke $'
00256 !
00257 !----------------------------------------------------------------------
00258 !
00259 !  Initialization
00260 !
00261 #ifdef VERBOSE
00262       print 9990, trim(ch_id), var_id, search%msg_intersections%first_tgt_method_id, search%sender
00263 
00264       call psmile_flushstd
00265 #endif /* VERBOSE */
00266 !
00267       ierror = 0
00268 !
00269       ndrbuf = 0
00270       ndibuf = 0
00271       n_recv_req = 0
00272 !
00273       method_id = Fields(var_id)%method_id
00274       mp => Methods(method_id)
00275 !
00276       grid_id   = mp%grid_id
00277       grid_info => Grids(grid_id)
00278 !
00279       global_grid_id = grid_info%global_grid_id
00280 !
00281 ! Uebergeben oder ein-dimensional
00282 !
00283       send_info => mp%send_infos_coupler(send_index)
00284       corner_pointer => grid_info%corner_pointer
00285 !
00286       comp_id = comp_info%comp_id
00287 !
00288       indices_req => extra_search%indices_req
00289       required    => extra_search%required
00290       len_req     => extra_search%len_req
00291       nreq        =  extra_search%n_req
00292 !
00293       virtual_cell_available = Associated (send_info%virtual)
00294 !
00295 #ifdef PRISM_ASSERTION
00296       if (grid_type /= grid_info%grid_type) then
00297 ! nicht uebergeben !
00298          print *, 'grid_type', grid_type, grid_info%grid_type
00299          call psmile_assert ( __FILE__, __LINE__, &
00300                  'inconsistent grid types')
00301       endif
00302 #endif
00303 !
00304 !  Rank of actual process in Component communicator
00305 !
00306       rank    = Comps(comp_id)%rank
00307 !
00308       extents => comp_info%all_extents
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 ! Can man die offset und partition information ausnutzen ?
00335 !
00336 !===> Allocate work vectors
00337 !
00338 !  n_faces = Number of block faces on which cell/control volume "i" is
00339 !            located.
00340 !  !!! n_faces wieder raus !!!
00341 !  faces   = Range on block faces which interfere with faces of
00342 !            cell/control volume "i".
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 !===> Get number of faces and intersection with block boundary
00383 !     for each cell to be searched
00384 !
00385 !  len_item = Number of pieces of information
00386 !             in integer buffer "ibuf"
00387 !             ibuf(1:3) = srcloc of extra point
00388 !                     (index of lower left corner)
00389 !             ibuf(4)   = Index of the extra point to be searched
00390 !                         This is a local index in the vectors
00391 !                         "required" and "indices_req" of the actual
00392 !                         partition "ipart".
00393 !             ibuf(5)   = Code for interpolation bases required
00394 !             ibuf(6)   = Virtual cell info (Only GaussReduced)
00395 !             TODO: ibuf (4, *) raus
00396 !  len_rtem = Number of pieces of information
00397 !             per point in Double Precison buffer "rbuf"
00398 !
00399 ! Ausserdem sollte man mehrfache srclocs rausfiltern
00400 !
00401       nx = nreq ! oder weniger, wenn nicht alle
00402                 ! cellen gesucht werden sollen
00403       len_item = ndim_3d + 2 ! default
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 !===> Get range of faces and transform range into cyclic range
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%npart, nd_len), stat = ierror)
00495 
00496       if ( ierror > 0 ) then
00497          ierrp (1) = ierror
00498          ierrp (2) = search%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 !     create lookup table for translating igrid into a icomp
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 !===> Control extents of other domains in order to find partner
00529 !
00530 ! igrid = current grid
00531 !
00532 !     loop over all grids/blocks (local and remote ones)
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 !        if current block is not part of the same global grid
00543          if (global_grid_id /= comp_info%all_extent_infos(4,igrid)) cycle
00544 !
00545 !        if the current block is the one we are doing the global search for
00546 !        (note ranks start with 0, ...)
00547          if (rank == dest_comp-1 .and. &
00548                comp_info%all_extent_infos(1,igrid) == grid_id) cycle
00549 !
00550 !===> ... Does the "igrid"-th grid and the range of faces coincide ?
00551 !
00552 !  n_int = Number of significant intersections
00553 !
00554 !cdir vector
00555          do i = 1, n_trans
00556 
00557             dinter (1, :, i) = max (face_range_transformed(1,:,i), &
00558                                     extents(1,:,igrid))
00559             dinter (2, :, i) = min (face_range_transformed(2,:,i), &
00560                                     extents(2,:,igrid))
00561 
00562             ! in addition to the intersections we need to check whether the grid
00563             ! extent and the transformed face range has a common face at the
00564             ! common_grid_range boundaries (normally at -180/180)
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                ! set dummy data
00580                dinter (1,:,i+n_trans) = 0
00581                dinter (2,:,i+n_trans) = -1
00582             endif
00583          end do ! i = 1, n_trans
00584 !
00585          n_int = 0
00586 !cdir vector
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 !===> ... Transform intersections back into original range
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 !===> ... Get control volumes which
00609 !         have intersection with the extent of the "n"-th grid
00610 !???  What about cyclic boundaries ????
00611 !???  use global offsets here !
00612 !
00613          do i = 1, n_int
00614 #ifdef DEBUGX
00615             print *, 'extent controlled', i, n_int, dinter(:,:,found(i))
00616 #endif
00617 !cdir vector
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 ! j
00630          end do ! i
00631 !
00632 !===> ... Get number of points to be sent
00633 !
00634          n_send = COUNT(send_mask)
00635          if (n_send == 0) cycle
00636 !
00637          mask_changed = .true.
00638 !
00639 !===> ... A grid of process "dest_comp" is a potential partner
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 !           Reallocate len_nsend
00648 
00649             old_len => len_nsend
00650             Allocate (len_nsend (search%npart, nd_len), stat = ierror)
00651 
00652             if ( ierror > 0 ) then
00653                ierrp (1) = ierror
00654                ierrp (2) = search%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 !===> ... Allocate send buffers
00674 !
00675 !  cf. subroutine PSMILe_Store_faces_3d_reg_real
00676 !
00677 !  Double precision vector buf (1:ndrbuf):
00678 !  buf (:) = Coordinates of the point to be searched
00679 !            ndim_3d data items
00680 !  buf (:) = Coordinates of corners of control volume/cell
00681 !            "len_rtem" data items
00682 !
00683 !  Integer buffer ibuf (1:ndibuf):
00684 !
00685 !  ibuf (:) = srcloc of extra point
00686 !             ndim_3d data items
00687 !  ibuf (:) = Index of the extra point to be searched
00688 !             1       data items
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 !===> Get Indices (in srclocs) of points to be sent
00732 !     TODO: Direkt in ibuf speichern !
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 ! das wurde fuer die alten Gaussreduced Gird (vor virtual_cell info)
00771 ! gebraucht; derzeit noch drin wg. conservative
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 !===> Store data in send buffers
00805 !     ??? compact machen !?!?
00806 !         Mehrfache Zellen Info's in ibuf und corner volumen
00807 !         eliminieren ?
00808 !     TODO: ipi und ip wieder vereinigen
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%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%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%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%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 ! TODO: If dest_comp == rank, direct call of psmile_search_donor_extra()
01000          if (dest_comp == rank) then
01001 !
01002             search_global%sender = dest
01003             search_global%msg_extra = msg_extra (1:nd_msg)
01004 !     search_global%len_msg = nd_msg
01005 !
01006             call psmile_search_donor_extra (search_global, tol, ierror)
01007             if (ierror > 0) return
01008 
01009          endif
01010 #endif /* TODO */
01011 !
01012 !===> ... Send request to destination process
01013 !         extra_msg is received by routine PSMILe_Enddef_extra and
01014 !         the other buffer are received by routine
01015 !         PSMILe_Enddef_action_extra.
01016 !
01017 !  extra_msg (1)  = Type of request
01018 !  extra_msg (2)  = Datatype
01019 !  extra_msg (3)  = Length of buffer "ibuf" containing the integer data.
01020 !  extra_msg (4)  = Length of buffer  "buf" containing the
01021 !                   real/real data.
01022 !  extra_msg (5)  = Global comp id of component
01023 !  extra_msg (6)  = Global var id "id_trans_out"
01024 !  extra_msg (7)  = Grid Type of sending grid
01025 !                   (PRISM_Reglonlatvrt, ...)
01026 !  extra_msg (8)  = Number of control volumes sent
01027 !  extra_msg (9)  = Number of integer data items per control
01028 !                   volume sent
01029 !  extra_msg (10) = Number of real/real data items per control
01030 !                   volume sent
01031 !  extra_msg (11) = Is data on global offset available ?
01032 !  extra_msg (12) = I-Index of global offset (partition)
01033 !  extra_msg (13) = J-Index of global offset (partition)
01034 !  extra_msg (14) = K-Index of global offset (partition)
01035 !
01036 !  extra_msg (15) = Index in len_nsend
01037 !  extra_msg (16) = num_neigh (Number of interpolation bases)
01038 !
01039 !  extra_msg (17) = local grid id on destination process
01040 !
01041          extra_msg (1)  = interpolation_methods (interp)
01042          extra_msg (2)  = PRISM_REAL
01043          extra_msg (3)  = len_ibuf
01044          extra_msg (4)  = len_rbuf
01045          extra_msg (5)  = comp_info%global_comp_id
01046          extra_msg (6)  = search%msg_intersections%transient_out_id
01047 !
01048          extra_msg (7)  = grid_type
01049          extra_msg (8)  = n_send
01050          extra_msg (9)  = len_item
01051          extra_msg (10) = len_rtem
01052 !
01053          extra_msg (15) = n_recv_req
01054          extra_msg (16) = num_neigh
01055 
01056          extra_msg (17) = comp_info%all_extent_infos(1,igrid) ! local grid id of destination process
01057 !
01058 #ifdef DEBUGX
01059          print *, 'extra_msg (1:4)', extra_msg (1:4)
01060          print *, 'grid shape', grid_valid_shape
01061          do i = 1, n_send
01062             print *, ibuf ((i-1)*len_item+1:i*len_item)
01063          end do
01064 #endif
01065 !
01066          if (Associated (grid_info%partition)) then
01067 !           First dimension: Blocks
01068 !           Last  dimension: Global indices of partition
01069 !
01070             extra_msg (11)  = 1
01071             extra_msg (12)  = grid_info%partition (1, 1)
01072             extra_msg (13)  = grid_info%partition (1, 2)
01073             extra_msg (14)  = grid_info%partition (1, 3)
01074 !
01075 !===> ... Transfer local coordinates into global coordinates
01076 !
01077             if (grid_info%grid_type == PRISM_Gaussreduced_regvrt) then
01078 
01079                call psmile_trans_loc2glob_gauss2 (grid_id, &
01080                          ibuf, len_item, n_send, ierror)
01081 
01082             else
01083 !
01084                call psmile_trans_loc2glob_3d (grid_id,     &
01085                          ibuf, len_item, n_send, ierror)
01086 
01087             endif
01088 !
01089             if (ierror > 0) return
01090 
01091          else
01092 !
01093 !           ... No partition info
01094 !
01095             extra_msg (11)   = 0
01096          endif
01097 !
01098 #ifdef DEBUG
01099          print *, ' Sending ', exttag, ' to destination ', dest
01100 #endif
01101          call psmile_bsend (extra_msg, nd_msgextra, MPI_INTEGER, &
01102                             dest, exttag, comm_psmile, ierror)
01103          if (ierror /= MPI_SUCCESS) then
01104             ierrp (1) = ierror
01105             ierrp (2) = dest
01106             ierrp (3) = exttag
01107 
01108             ierror = PRISM_Error_Send
01109 
01110             call psmile_error (ierror, 'psmile_bsend(msg)', &
01111                                ierrp, 3, __FILE__, __LINE__ )
01112             return
01113          endif
01114 !
01115 #ifdef DEBUG
01116          print *, ' Sending ', exttag, ' to destination ', dest
01117 #endif
01118          call psmile_bsend (ibuf, len_ibuf, MPI_INTEGER, &
01119                             dest, exttag, comm_psmile, ierror)
01120          if (ierror /= MPI_SUCCESS) then
01121             ierrp (1) = ierror
01122             ierrp (2) = dest
01123             ierrp (3) = exttag
01124 
01125             ierror = PRISM_Error_Send
01126 
01127             call psmile_error (ierror, 'psmile_bsend(ibuf)', &
01128                                ierrp, 3, __FILE__, __LINE__ )
01129             return
01130          endif
01131 !
01132 #ifdef DEBUG
01133          print *, ' Sending ', exttag, ' to destination ', dest
01134 #endif
01135          call psmile_bsend (buf, len_rbuf, MPI_REAL, &
01136                             dest, exttag, comm_psmile, ierror)
01137          if (ierror /= MPI_SUCCESS) then
01138             ierrp (1) = ierror
01139             ierrp (2) = dest
01140             ierrp (3) = exttag
01141 
01142             ierror = PRISM_Error_Send
01143 
01144             call psmile_error (ierror, 'psmile_bsend(buf)', &
01145                                ierrp, 3, __FILE__, __LINE__ )
01146             return
01147          endif
01148 !
01149 !===> ... Setup request for answer of destination process
01150 !         ??? wirklich; sollte wohl von dem normalen abgehandelt werden
01151 !
01152          if (n_recv_req == 1) then
01153             call MPI_Irecv (answer, nd_msgextra, MPI_INTEGER, MPI_ANY_SOURCE, &
01154                             rexttag, comm_psmile, recv_req, ierror)
01155             if (ierror /= MPI_SUCCESS) then
01156 
01157                ierrp (1) = ierror
01158                ierrp (2) = dest
01159                ierrp (3) = rexttag
01160 
01161                ierror = PRISM_Error_Recv
01162 
01163                call psmile_error ( ierror, 'MPI_Irecv', &
01164                                    ierrp, 3, __FILE__, __LINE__ )
01165                return
01166 
01167             endif
01168          endif
01169 !
01170       end do ! igrid
01171 !
01172 !===> Deallocate
01173 !
01174       Deallocate (send_mask, n_faces, faces, igrid_to_dest_comp)
01175 !
01176       if (ndrbuf > 0) Deallocate (buf)
01177       if (ndibuf > 0) Deallocate (ibuf)
01178 !
01179 !-------------------------------------------------------------------------
01180 !     Wait for answers
01181 !      Note: (2) Message from a partner with an intersection is disabled
01182 !                since request for grid receive (lrequest (3)) is setup
01183 !                if such a message is found.
01184 !            ??? Koennte man enablen wenn lrequest(3) == MPI_REQUEST_NULL
01185 !            (3) Old receive of grid data (paction%lrequest (3)) is disabled
01186 !-------------------------------------------------------------------------
01187 !
01188 ! nrecv      = Number of significant messages received
01189 ! nprev_recv = Number of data items (per coordinate) received
01190 !              in prevoius receives.
01191 !
01192       nprev_recv = 0
01193       nrecv = 0
01194 !
01195       if (n_recv_req > 0) then
01196          Allocate (extra_search%real_bufs(n_recv_req),  &
01197                    send_info%sender_global(n_recv_req), &
01198                    send_info%len_sent(n_recv_req),      &
01199                    send_info%msg_id(n_recv_req),        &
01200                    stat = ierror)
01201          if ( ierror > 0 ) then
01202             ierrp (1) = ierror
01203             ierrp (2) = n_recv_req
01204 
01205             ierror = PRISM_Error_Alloc
01206             call psmile_error ( ierror, 'extra_search%real_bufs', &
01207                                 ierrp, 2, __FILE__, __LINE__ )
01208             return
01209          endif
01210       endif
01211 !
01212       save_lreq (2:3) = paction%lrequest (2:3)
01213 !
01214       paction%lrequest (2) = MPI_REQUEST_NULL
01215 !
01216       do n = 1, n_recv_req
01217          paction%lrequest (3) = recv_req
01218 !
01219          index = 0
01220          do while (index /= 3) 
01221 #ifdef DEBUG
01222             print *, trim(ch_id), ': psmile_global_search_real: action%nreq, recv_req', &
01223                      paction%nreq, recv_req
01224             call psmile_flushstd
01225 #endif
01226 !
01227             call MPI_Waitany (paction%nreq, paction%lrequest, index, status, ierror)
01228 
01229             if ( ierror /= MPI_SUCCESS ) then
01230                ierrp (1) = ierror
01231                ierrp (2) = status (MPI_SOURCE)
01232                ierrp (3) = status (MPI_TAG)
01233 
01234                ierror = PRISM_Error_MPI
01235 
01236                call psmile_error ( ierror, 'MPI_Waitany', &
01237                                    ierrp, 3, __FILE__, __LINE__ )
01238                return
01239             endif
01240 
01241 #ifdef PRISM_ASSERTION
01242             if (index == MPI_UNDEFINED) then
01243                call psmile_assert ( __FILE__, __LINE__, &
01244                                     'request list is empty')
01245             endif
01246 #endif
01247 !
01248             if (index /= 3) then
01249                call psmile_enddef_action (search, index, status, ierror)
01250                if (ierror > 0) return
01251             endif
01252          end do ! while
01253 !
01254 !===> ... Get 
01255 !
01256          sender = status (MPI_SOURCE)
01257          len_ibuf = answer (3)
01258          len_rbuf = answer (4)
01259 !
01260 #ifdef VERBOSE
01261          print 9960, trim(ch_id), sender, len_ibuf, len_rbuf
01262 
01263          call psmile_flushstd
01264 #endif /* VERBOSE */
01265 !
01266 #ifdef PRISM_ASSERTION
01267          if (len_ibuf < 0) then
01268             print *, 'len_ibuf =', len_ibuf
01269             call psmile_assert (__FILE__, __LINE__, &
01270                  "len_ibuf should be >= 0")
01271          endif
01272 !
01273          if (len_rbuf < 0) then
01274             print *, 'len_rbuf =', len_rbuf
01275             call psmile_assert (__FILE__, __LINE__, &
01276                  "len_rbuf should be >= 0")
01277          endif
01278 #endif
01279 !
01280          if (len_ibuf > 0) then
01281             Allocate (ibuf (1:len_ibuf), stat = ierror)
01282 
01283             if ( ierror > 0 ) then
01284                ierrp (1) = ierror
01285                ierrp (2) = len_ibuf
01286 
01287                ierror = PRISM_Error_Alloc
01288                call psmile_error ( ierror, 'ibuf', &
01289                                    ierrp, 2, __FILE__, __LINE__ )
01290                return
01291             endif
01292 !
01293             call MPI_Recv (ibuf, len_ibuf, MPI_INTEGER, sender, &
01294                            rexttag, comm_psmile, status, ierror)
01295             if (ierror /= MPI_SUCCESS) then
01296 
01297                ierrp (1) = ierror
01298                ierrp (2) = sender
01299                ierrp (3) = rexttag
01300 
01301                ierror = PRISM_Error_Recv
01302 
01303                call psmile_error ( ierror, 'MPI_Recv(ibuf)', &
01304                                    ierrp, 3, __FILE__, __LINE__ )
01305                return
01306             endif
01307 !
01308             Allocate (buf (1:len_rbuf), stat = ierror)
01309 
01310             if ( ierror > 0 ) then
01311                ierrp (1) = ierror
01312                ierrp (2) = len_rbuf
01313 
01314                ierror = PRISM_Error_Alloc
01315                call psmile_error ( ierror, 'buf', &
01316                                    ierrp, 2, __FILE__, __LINE__ )
01317                return
01318             endif
01319 !
01320             call MPI_Recv (buf, len_rbuf, MPI_REAL, sender, &
01321                            rexttag, comm_psmile, status, ierror)
01322             if (ierror /= MPI_SUCCESS) then
01323 
01324                ierrp (1) = ierror
01325                ierrp (2) = sender
01326                ierrp (3) = rexttag
01327 
01328                ierror = PRISM_Error_Recv
01329 
01330                call psmile_error ( ierror, 'MPI_Recv(rbuf)', &
01331                                    ierrp, 3, __FILE__, __LINE__ )
01332                return
01333             endif
01334 !
01335 !===> ...
01336 !
01337 ! n_send  = Number of cells to be searched
01338 ! n_found = Number of interpolation bases found
01339 ! n_liste = Number of data items to be used for interpolation
01340 !           n_liste <= n_found since a data item can be used
01341 !           multiply for interpolation.
01342 !
01343             n_send  = answer (7)
01344             n_liste = answer (8)
01345             n_found = answer (9)
01346             ireq = answer (15)
01347 !
01348 #ifdef PRISM_ASSERTION
01349             if (answer (1) /= interpolation_methods (interp)) then
01350                 print *, 'answer(1)', answer (1), interpolation_methods (interp)
01351                call psmile_assert (__FILE__, __LINE__, &
01352                     "incorrect message received ? interpolations doesn't fit")
01353             endif
01354 !
01355             if (ireq < 1 .or. ireq > n_recv_req) then
01356                print *, 'ireq, n_recv_req', ireq, n_recv_req
01357                call psmile_assert (__FILE__, __LINE__, &
01358                     "ireq must be in range of 1:n_recv_req")
01359             endif
01360 #endif
01361 !
01362             nrecv = nrecv + 1
01363 !
01364             extra_search%real_bufs (nrecv)%vector  => buf
01365             send_info%sender_global (nrecv) = sender
01366             send_info%len_sent      (nrecv) = n_liste
01367             send_info%msg_id        (nrecv) = answer (16)
01368 !
01369             if (nprev_recv > 0) then
01370 !              Increment index_found (for collected list of buf's)
01371 !              except for point masked out
01372 !cdir vector
01373                do i = 2*n_send+1, 2*n_send+n_found
01374                   if (ibuf(i) /= masked_out) then
01375                      ibuf (i) = ibuf (i) + nprev_recv
01376                   endif
01377                end do
01378             endif
01379 !
01380             nprev_recv = nprev_recv + n_liste
01381 !
01382 !===> ... Add information on points found by global search 
01383 !
01384             call psmile_add_points_found (grid_id, search, extra_search, &
01385                         ibuf (1:n_send), ibuf (n_send+1:2*n_send),       &
01386                         n_send, len_nsend(1:search%npart,ireq),          &
01387                         ibuf (2*n_send+1:2*n_send+n_found), n_found,     &
01388                         neighbors_3d, nloc, num_neigh,                   &
01389                         grid_valid_shape, use_how, ierror)
01390                    if (ierror > 0) return
01391 !
01392             Deallocate (ibuf)
01393          endif
01394 !
01395 !===> ... Setup new request
01396 !
01397          if (n < n_recv_req) then
01398             call MPI_Irecv (answer, nd_msgextra, MPI_INTEGER, MPI_ANY_SOURCE, &
01399                               rexttag, comm_psmile, recv_req, &
01400                               ierror)
01401             if (ierror /= MPI_SUCCESS) then
01402 
01403                ierrp (1) = ierror
01404                ierrp (2) = dest
01405                ierrp (3) = rexttag
01406 
01407                ierror = PRISM_Error_Recv
01408 
01409                call psmile_error ( ierror, 'MPI_Irecv', &
01410                                     ierrp, 3, __FILE__, __LINE__ )
01411                return
01412             endif
01413          endif
01414       end do
01415 !
01416 !    ... Restore original requests (2:3)
01417 !
01418 #ifdef PRISM_ASSERTION
01419       if (paction%lrequest (2) /= MPI_REQUEST_NULL .or. &
01420           paction%lrequest (3) /= MPI_REQUEST_NULL) then
01421          print *, 'request: ', paction%lrequest (2:3)
01422          call psmile_assert ( __FILE__, __LINE__, &
01423                              'Illegal request stored')
01424 
01425       endif
01426 #endif
01427 !
01428       paction%lrequest (2:3) = save_lreq (2:3)
01429 !
01430       Deallocate (len_nsend)
01431 !
01432 !===> Save number of data to be received in "extra_search"
01433 !
01434       send_info%nrecv    = nrecv
01435       send_info%num2recv = nprev_recv
01436 !
01437 #ifdef PRISM_ASSERTION
01438       if (nrecv > 0) then
01439          if (send_info%num2recv /= SUM (send_info%len_sent(1:nrecv))) then
01440             print *, 'num2recv', send_info%num2recv, &
01441                      SUM (send_info%len_sent(1:nrecv)), nrecv
01442             call psmile_assert (__FILE__, __LINE__, &
01443                  "inconsistent num2recv")
01444          endif
01445       endif
01446 #endif
01447 !
01448 !-------------------------------------------------------------------------
01449 !===> All done
01450 !-------------------------------------------------------------------------
01451 !
01452 #ifdef VERBOSE
01453       print 9980, trim(ch_id), ierror
01454 
01455       call psmile_flushstd
01456 #endif /* VERBOSE */
01457 !
01458 !  Formats:
01459 !
01460 
01461 #ifdef VERBOSE
01462 
01463 9990 format (1x, a, ': psmile_global_search_real: var_id', i3, &
01464                     ' to ', i3, '(', i2, ')')
01465 9980 format (1x, a, ': psmile_global_search_real: eof ierror =', i3)
01466 9970 format (1x, a, ': psmile_global_search_real: send from', i3, &
01467                     ' to', i3, '[', i3, '], n_send =', i6)
01468 
01469 9960 format (1x, a, ': psmile_global_search_real: got rexttag message:', &
01470                     ' sender ', i4, ', len_ibuf, len_rbuf', 2i8)
01471 9950 format (1x, a, ': psmile_global_search_real: before waitany :', &
01472                     'nreq =', i4, ', recv_req ', i8)
01473 #endif /* VERBOSE */
01474 
01475 #ifdef DEBUG
01476 #endif
01477 
01478       end subroutine PSMILe_global_search_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1