psmile_global_search_nnx_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_nnx_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_global_search_nnx_real (comp_info, search,        &
00012                         var_id, tgt_coords_x, tgt_coords_y, tgt_coords_z, &
00013                         neighbors_3d, nloc, num_neigh, nb_extra,          &
00014                         extra_search, send_index,                         &
00015                         ierror)
00016 !
00017 ! !USES:
00018 !
00019       use PRISM
00020 !
00021       use psmile_grid, only : psmile_transform_extent_cyclic, &
00022                               psmile_transform_extent_back, &
00023                               max_num_trans_parts, &
00024                               code_no_trans, &
00025                               common_grid_range
00026 !
00027       use PSMILe, dummy_interface => PSMILe_global_search_nnx_real
00028 
00029       Implicit none
00030 !
00031 ! !INPUT PARAMETERS:
00032 
00033       Type (Enddef_comp), Intent (In) :: comp_info
00034 !
00035 !     Info on the component in which the donor cells
00036 !     should be searched.
00037 !
00038       Type (Enddef_search), Intent (InOut) :: search
00039 
00040 !     Info's on coordinates to be searched
00041 
00042       Integer,            Intent (In) :: var_id
00043 
00044 !     Handle to the grid function
00045 
00046       Integer,            Intent (In) :: send_index
00047 
00048 !     Send index of coupler send for the method and coordinates
00049 
00050       Integer,            Intent (In) :: nloc
00051 !
00052 !     Total  Number of locations to be transferred
00053 !
00054       Real,   Intent (In)             :: tgt_coords_x (nloc)
00055       Real,   Intent (In)             :: tgt_coords_y (nloc)
00056       Real,   Intent (In)             :: tgt_coords_z (nloc)
00057 
00058 !     Coordinates of the target points (which were searched and found)
00059 !
00060       Integer,            Intent (In) :: num_neigh
00061 !
00062 !     Last dimension of neighbors array "neighbors_3d" and
00063 !     number of neighbors to be searched.
00064 !     besser n_corners
00065 
00066       Integer,            Intent (In) :: nb_extra
00067 !
00068 !     Number of nearest neigbhbours to be searched per extra point
00069 !
00070 ! !INPUT/OUTPUT PARAMETERS:
00071 !
00072       Integer,            Intent (InOut) :: neighbors_3d (ndim_3d, nloc, 
00073                                                           num_neigh)
00074 
00075 !     Indices of neighbor locations (interpolation bases)
00076 !
00077       Type (Extra_search_info), Intent (InOut) :: extra_search
00078 !
00079 !     Number of locations where
00080 !       (*) required mask values were not true
00081 
00082 ! !OUTPUT PARAMETERS:
00083 !
00084       Integer,           Intent (Out) :: ierror
00085 
00086 !     Returns the error code of PSMILe_global_search_nnx_real;
00087 !             ierror = 0 : No error
00088 !             ierror > 0 : Severe error
00089 !
00090 ! !LOCAL PARAMETERS
00091 !
00092 !  DBLE_EARTH = Earth radius in meter
00093 !
00094 !  LEN_ITEM = Number of pieces of information
00095 !             in integer buffer "ibuf"
00096 !             1: index in indices/dist_real/send_mask
00097 !  INLP_LEN = Inner loop length for optimization
00098 !  lon   = Index of Longitudes in arrays "boxes"
00099 !  lat   = Index of Latitudes  in arrays "boxes"
00100 !  vrt   = Index of height     in arrays "boxes"
00101 !
00102 !
00103       Real, Parameter                 :: real_earth = 6400000.0d0
00104 !
00105       Integer, Parameter              :: len_item = 1
00106 !
00107       Integer, Parameter              :: inlp_len = 256
00108 !
00109       Integer, Parameter              :: lon = 1
00110       Integer, Parameter              :: lat = 2
00111       Integer, Parameter              :: vrt = 3
00112 !
00113 !
00114 ! !LOCAL VARIABLES
00115 !
00116 !     ... Method pointer
00117 !
00118       Integer                         :: comp_id, grid_id
00119       Integer                         :: global_grid_id
00120       Integer, Pointer                :: indices (:)
00121 !
00122 !     ... Info's on actual component
00123 !
00124       Integer                         :: igrid, rank
00125       Integer                         :: n_total
00126       Integer                         :: grid_type
00127 !
00128 !     ... lookup table
00129 !
00130       Integer, allocatable            :: igrid_to_dest_comp(:)
00131 !
00132 !     ... loop variables
00133 !
00134       Integer                         :: j, n
00135       Integer                         :: ibeg, iend
00136 !
00137 !     ... for data to other processes
00138 !
00139       Logical                         :: mask_changed
00140       Integer                         :: n_extra
00141 !
00142       Logical, Allocatable            :: send_mask (:)
00143       Real, Allocatable               :: boxes (:, :, :)
00144       Real (PSMILe_float_kind)        :: range_box (2, ndim_3d)
00145       Real (PSMILe_float_kind)        :: dist_max (2)
00146       Real                            :: real_huge
00147 !
00148 !     ... For locations controlled
00149 !
00150 !     Real, Pointer                   :: sin_search (:, :)
00151 !     Real, Pointer                   :: cos_search (:, :)
00152 !     Real, Pointer                   :: z_search (:)
00153       Real, Pointer                   :: dist_real (:, :)
00154       Real                            :: r_earth2, r_lat
00155       Real                            :: delta, sin_d, cos_lat
00156 !
00157 !     ... transformed faces
00158 !
00159       Integer                         :: n_trans, n_int
00160 !
00161       Integer                         :: tr_codes (max_num_trans_parts)
00162       Integer                         :: found (max_num_trans_parts)
00163       Real (PSMILe_float_kind)        :: dinter_trans (2, ndim_3d)
00164       Real (PSMILe_float_kind)        :: transformed (2, ndim_3d, 
00165                                                       max_num_trans_parts), 
00166                                          dinter (2, ndim_3d, max_num_trans_parts)
00167 !
00168       Real, Save                      :: period2 (ndim_3d)
00169 !
00170 !     ... for points found
00171 !
00172       Integer                         :: n_found, n_liste
00173       Integer                         :: ip_dist
00174 !
00175 !     ... for send to partners
00176 !
00177       Type (Send_information), Pointer :: send_info
00178       Type (enddef_msg_extra)         :: extra_msg
00179       Integer                         :: extra_msg_buf (msg_extra_size)
00180 !
00181       Integer                         :: i, ip, ipi, ireq
00182       Integer                         :: len_ibuf, ndibuf, n_send
00183       Integer                         :: len_rbuf, len_rtem, ndrbuf
00184       Integer, Allocatable            :: ibuf (:)
00185       Real, Pointer                   :: buf (:)
00186 !
00187 !     ... for recv from partners
00188 !
00189       Integer                         :: answer (msg_extra_size)
00190       Integer, Allocatable            :: selected (:, :)
00191       Type (Select_search_info), Pointer  :: sel_info (:)
00192 !
00193 !     ... for communication
00194 !
00195       Integer                         :: dest, index, dest_comp, sender
00196       Integer                         :: nrecv, n_recv_req
00197       Integer                         :: recv_req
00198       Integer                         :: save_lreq (2:3)
00199       Integer                         :: status (MPI_STATUS_SIZE)
00200 !
00201 !     ... for error handling
00202 !
00203       Integer, Parameter              :: nerrp = 3
00204       Integer                         :: ierrp (nerrp)
00205 !
00206 ! !DESCRIPTION:
00207 !
00208 ! Subroutine "PSMILe_global_search_nnx_real" performs a global 
00209 ! nearest neighbour search for the interpolation points 
00210 ! which require an extra seach.
00211 !
00212 ! TODO: Save buffer indices in indices/send_mask/dist_real (see len_item)
00213 !       locally instead of passing them !!!
00214 !
00215 ! !REVISION HISTORY:
00216 !
00217 !   Date      Programmer   Description
00218 ! ----------  ----------   -----------
00219 ! 09.10.06    H. Ritzdorf  created
00220 !
00221 !EOP
00222 !----------------------------------------------------------------------
00223 !
00224 !  $Id: psmile_global_search_nnx_real.F90 3130 2011-04-12 14:39:14Z hanke $
00225 !  $Author: hanke $
00226 !
00227    Character(len=len_cvs_string), save :: mycvs = 
00228        '$Id: psmile_global_search_nnx_real.F90 3130 2011-04-12 14:39:14Z hanke $'
00229 !
00230 !----------------------------------------------------------------------
00231 !
00232 !  Initialization
00233 !
00234 #ifdef VERBOSE
00235       print 9990, trim(ch_id), var_id, search%msg_intersections%field_info%tgt_method_id, search%sender
00236 #endif /* VERBOSE */
00237 
00238       call psmile_flushstd
00239 !
00240       ierror = 0
00241 !
00242       ndrbuf = 0
00243       ndibuf = 0
00244       n_recv_req = 0
00245 !
00246       r_earth2 = 1.0 / (real_earth * 2.0)
00247       r_lat = 1.0 / (real_earth * real_deg2rad)
00248 !
00249       period2 = (common_grid_range(2,:) - common_grid_range(1,:)) / 2.0d0
00250 !
00251       grid_id        = Methods(Fields(var_id)%method_id)%grid_id
00252       global_grid_id = Grids(grid_id)%global_grid_id
00253       grid_type      = Grids(grid_id)%grid_type
00254 !
00255       comp_id = comp_info%comp_id
00256 !
00257 ! Uebergeben oder ein-dimensional
00258 !
00259       send_info => Methods(Fields(var_id)%method_id)%send_infos_coupler(send_index)
00260 !
00261       indices    => extra_search%indices
00262       dist_real  => extra_search%dist_real
00263 !
00264       n_extra = extra_search%n_extra
00265 !
00266 #ifdef PRISM_ASSERTION
00267 #ifdef HUHU
00268       sin_search => extra_search%sin_search_real
00269       cos_search => extra_search%cos_search_real 
00270         z_search => extra_search%z_search_real
00271       if ( .not. Associated (extra_search%dist_real)       .or. &
00272            .not. Associated (extra_search%cos_search_real) .or. &
00273            .not. Associated (extra_search%sin_search_real) .or. &
00274            .not. Associated (extra_search%z_search_real)   ) then
00275 #else
00276       if ( .not. Associated (extra_search%dist_real) ) then
00277 #endif
00278 !
00279          call psmile_assert (__FILE__, __LINE__, &
00280                              "arrays should be allocated and set")
00281       endif
00282 #endif
00283 !
00284 !  Rank of actual process in Component communicator
00285 !
00286       rank    = Comps(comp_id)%rank
00287 !
00288 ! n_total   = Total number of grids/extents
00289 ! igrid_end = Last  grid/extent of actual process
00290 ! igrid     = Index of actual grid
00291 !
00292       n_total   = SUM (comp_info%Number_of_grids_vector(:))
00293 !
00294       igrid = search%msg_intersections%first_src_all_extents_grid_id
00295 !
00296 ! Can man die offset und partition information ausnutzen ?
00297 !
00298 !===> Allocate work vectors
00299 !
00300 ! boxes(i)     = Box around ith extra point with interior radius
00301 !                "dist_real (i)"
00302 ! send_mask(i) = Does the i-th extra point to be transferred 
00303 !                to the corresponding process ?
00304 !
00305       Allocate (boxes (2, ndim_3d, n_extra), STAT = ierror)
00306 
00307       if ( ierror > 0 ) then
00308          ierrp (1) = ierror
00309          ierrp (2) = n_extra * (2 * ndim_3d)
00310 
00311          ierror = PRISM_Error_Alloc
00312          call psmile_error ( ierror, 'boxes', &
00313                              ierrp, 2, __FILE__, __LINE__ )
00314          return
00315       endif
00316 !
00317       Allocate (send_mask (n_extra), STAT = ierror)
00318 
00319       if ( ierror > 0 ) then
00320          ierrp (1) = ierror
00321          ierrp (2) = n_extra
00322 
00323          ierror = PRISM_Error_Alloc
00324          call psmile_error ( ierror, 'send_mask', &
00325                              ierrp, 2, __FILE__, __LINE__ )
00326          return
00327       endif
00328 !
00329       mask_changed = .true.
00330 !
00331 !===> Compute range (bounding box) for global nearest neighbour
00332 !     search of extra points
00333 !
00334 !     ?dist_real = huge ? neighbors_3d (:, i, nb_extra) auswerten ?
00335 !     !!! TODO: Diese Punkte loesen eine Suche in allen Bloecken aus !
00336 !               Besser erst lokal in der Nachbarschaft suchen !
00337 !
00338 !==> Get maximal extension for 3rd direction (z-direction)
00339 !    ? einmal berechnen ?
00340 !
00341      real_huge = huge (dist_real(1,1))
00342      dist_max (1) = comp_info%all_extent_infos(igrid)%extent (1, vrt)
00343      dist_max (2) = comp_info%all_extent_infos(igrid)%extent (2, vrt)
00344 !
00345 #ifdef DEBUGX
00346      print *, "dist_max", dist_max
00347 #endif
00348 !
00349       do n = 1, n_total
00350         if (global_grid_id == comp_info%all_extent_infos(n)%global_grid_id) then
00351            dist_max (1) = min (dist_max(1), comp_info%all_extent_infos(n)%extent(1, vrt))
00352            dist_max (2) = max (dist_max(2), comp_info%all_extent_infos(n)%extent(2, vrt))
00353         endif
00354       end do
00355 #ifdef DEBUGX
00356      print *, "dist_max", dist_max
00357 #endif
00358 !
00359       do ibeg = 1, n_extra, inlp_len
00360          iend = min (n_extra, ibeg + inlp_len - 1)
00361 !
00362 !        Construct boxes of interest in lon direction
00363 !        Delta lon = 2 * asin (sin (dist_real/(2*real_earth)) / cos (lat))
00364 !                    / real_deg2rad
00365 !cdir vector
00366          do i = ibeg, iend
00367             if (dist_real(i, nb_extra) * r_earth2 >= real_pih) then
00368                sin_d = 1.0
00369             else
00370                sin_d = abs(sin(dist_real(i, nb_extra) * r_earth2))
00371             endif
00372 
00373             cos_lat = cos (tgt_coords_y(indices(i)) * real_deg2rad)
00374 
00375             if (sin_d >= cos_lat) then
00376                delta = period2(lon)
00377             else
00378                delta = 2.0d0 * asin (sin_d / cos_lat) / real_deg2rad
00379                delta = min (period2(lon), delta)
00380             endif
00381 
00382             boxes (1, lon, i) = tgt_coords_x(indices(i)) - delta
00383             boxes (2, lon, i) = tgt_coords_x(indices(i)) + delta
00384          end do
00385 
00386 !cdir vector
00387          do i = ibeg, iend
00388             if (boxes (2,lon,i) - boxes (1,lon,i) > period2(lon)*2) then
00389                boxes (1,lon,i) = -period2(lon)
00390                boxes (2,lon,i) =  period2(lon)
00391             else if (boxes (1,lon,i) < -period2(lon)*2) then
00392                boxes (1,lon,i) = boxes (1,lon,i) + period2(lon)*2
00393                boxes (2,lon,i) = boxes (2,lon,i) + period2(lon)*2
00394             else if (boxes (2,lon,i) >  period2(lon)*2) then
00395                boxes (1,lon,i) = boxes (1,lon,i) - period2(lon)*2
00396                boxes (2,lon,i) = boxes (2,lon,i) - period2(lon)*2
00397             endif
00398          end do ! i
00399 
00400 !        Construct boxes of interest in lat direction
00401 !        Delta lat = dist_real / (real_earth * real_deg2rad)
00402 
00403 !cdir vector
00404          do i = ibeg, iend
00405             boxes (1, lat, i) = max (tgt_coords_y(indices(i)) -                        &
00406                                      min (period2(lat), dist_real(i, nb_extra)*r_lat), &
00407                                      real(common_grid_range(1,lat)))
00408             boxes (2, lat, i) = min (tgt_coords_y(indices(i)) +                        &
00409                                      min (period2(lat), dist_real(i, nb_extra)*r_lat), &
00410                                      real(common_grid_range(2,lat)))
00411          end do
00412       end do ! j
00413 !
00414       range_box (1, lon) = minval (boxes (1, lon, :))
00415       range_box (2, lon) = maxval (boxes (2, lon, :))
00416 !
00417       range_box (1, lat) = minval (boxes (1, lat, :))
00418       range_box (2, lat) = maxval (boxes (2, lat, :))
00419 !
00420 !        Construct boxes of interest in vertical direction
00421 !cdir vector
00422       do i = 1, n_extra
00423          if (dist_real(i, nb_extra) == real_huge) then
00424             boxes (1, vrt, i) = dist_max (1)
00425             boxes (2, vrt, i) = dist_max (2)
00426          else
00427             boxes (1, vrt, i) = tgt_coords_z(indices(i)) -             &
00428                               dist_real(i, nb_extra)
00429             boxes (2, vrt, i) = tgt_coords_z(indices(i)) +             &
00430                               dist_real(i, nb_extra)
00431 #ifdef DEBUGX
00432      print *, "dist_real", i, indices(i), dist_real(i, nb_extra)
00433 #endif
00434          endif
00435       end do
00436 !
00437       range_box (1, vrt) = minval (boxes (1, vrt, :))
00438       range_box (2, vrt) = maxval (boxes (2, vrt, :))
00439 !
00440 #ifdef PRISM_ASSERTION
00441       if (range_box (2,lon) - range_box (1,lon) > period2(lon)*4) then
00442           print *, 'range in lon direction', range_box (1:2,lon)
00443          call psmile_assert ( __FILE__, __LINE__, &
00444                              'range_box too large in lon direction')
00445       endif
00446 !
00447       if (range_box (2,lat) - range_box (1,lat) > period2(lat)*4) then
00448           print *, 'range in lat direction', range_box (1:2,lat)
00449          call psmile_assert ( __FILE__, __LINE__, &
00450                              'range_box too large in lat direction')
00451       endif
00452 #endif
00453 !
00454 !===> Transform range "range_box" into cyclic ranges
00455 !
00456       call psmile_transform_extent_cyclic (grid_type,           &
00457                   range_box, transformed, tr_codes, n_trans, ierror)
00458       if (ierror > 0) return
00459 !
00460 !     create lookup table for translating igrid into a icomp
00461 !
00462       Allocate (igrid_to_dest_comp(SUM(comp_info%Number_of_grids_vector)), stat = ierror)
00463 
00464       if ( ierror > 0 ) then
00465          ierrp (1) = ierror
00466          ierrp (2) = SUM(comp_info%Number_of_grids_vector)
00467 
00468          ierror = PRISM_Error_Alloc
00469          call psmile_error ( ierror, 'igrid_to_dest_comp', &
00470                              ierrp, 2, __FILE__, __LINE__ )
00471          return
00472       endif
00473 
00474       igrid = 0
00475       do dest_comp = 1, comp_info%size
00476          do j = 1, comp_info%Number_of_grids_vector(dest_comp)
00477             igrid = igrid + 1
00478             igrid_to_dest_comp(igrid) = dest_comp
00479          enddo
00480       enddo
00481 !
00482 !===> Control extents of other domains in order to find partner
00483 !
00484 ! igrid = current grid
00485 !
00486 !     loop over all grids/blocks (local and remote ones)
00487       do igrid = 1, SUM(comp_info%Number_of_grids_vector)
00488 
00489          dest_comp = igrid_to_dest_comp(igrid)
00490 
00491          if (mask_changed) then
00492             mask_changed  = .false.
00493             send_mask (:) = .false.
00494          endif
00495 !
00496 ! TODO:
00497 ! Das ist mir zu teuer; erst mit den extents in einer einfachen
00498 ! schleife potentiell partner finden; siehe found (:)
00499 ! in psmile_find_intersect
00500 !
00501 !        if current block is not part of the same global grid
00502          if (global_grid_id /= comp_info%all_extent_infos(igrid)%global_grid_id) cycle
00503 !
00504 !===> ... Is this my info
00505 !         (note ranks start with 0, ...)
00506 !
00507          if (rank == dest_comp-1 .and. &
00508                comp_info%all_extent_infos(igrid)%local_grid_id == grid_id) cycle
00509 !
00510 !===> ... Does the "n"-th grid and "range_box" coincide ?
00511 !        Rausziehen und vektorisieren ?
00512 !
00513 !  n_int = Number of significant intersections
00514 !
00515 !cdir vector
00516          do i = 1, n_trans
00517             dinter (1, :, i) = max (transformed(1,:,i), &
00518                                     comp_info%all_extent_infos(igrid)%extent(1,:))
00519             dinter (2, :, i) = min (transformed(2,:,i), &
00520                                     comp_info%all_extent_infos(igrid)%extent(2,:))
00521          end do
00522 !
00523          n_int = 0
00524 !cdir vector
00525          do i = 1, n_trans
00526             if (minval(dinter (2,:,i) - dinter (1,:,i)) >= 0.0d0) then
00527                n_int = n_int + 1
00528                found (n_int) = i
00529             endif
00530          end do
00531 !
00532          if (n_int == 0) cycle
00533 !
00534 !===> ... Transform intersections back into original range
00535 !
00536          do i = 1, n_int
00537             if (tr_codes(found(i)) /= code_no_trans) then
00538                call psmile_transform_extent_back (tr_codes(found(i)), &
00539                   dinter(:, :, found(i)), dinter_trans, 1, ierror)
00540                if ( ierror /= 0 ) return
00541 
00542                dinter (:, :, found(i)) = dinter_trans
00543             endif
00544          end do
00545 !
00546 !===> ... Get the extra points whose radius has an intersection
00547 !         with the extent of the "n"-th grid
00548 !
00549 ! send_mask(i) = Does the i-th extra point to be transferred 
00550 !                to the corresponding process ?
00551 !???  What about cyclic boundaries ????
00552 !???  use global offsets here !
00553 !
00554 ! Bedingung: 1. Punkt is in dinter drin
00555 !            2. vereinfachte distanz kleiner als dist_real
00556 !
00557          do i = 1, n_int
00558 #ifdef DEBUGX
00559 print *, 'extent controlled', i, n_int, dinter(:,:,found(i))
00560 #endif
00561 !cdir vector
00562             do j = 1, n_extra
00563                send_mask (j) = send_mask (j) .or. &
00564                      (boxes(2,lon,j) >= dinter(1,lon,found(i)) .and. &
00565                       boxes(1,lon,j) <= dinter(2,lon,found(i)) .and. &
00566                       boxes(2,lat,j) >= dinter(1,lat,found(i)) .and. &
00567                       boxes(1,lat,j) <= dinter(2,lat,found(i)) .and. &
00568                       boxes(2,vrt,j) >= dinter(1,vrt,found(i)) .and. &
00569                       boxes(1,vrt,j) <= dinter(2,vrt,found(i)))
00570 #ifdef DEBUGX
00571 print *, j, boxes(:,:,j)
00572 #endif
00573             end do ! j
00574          end do ! i
00575 !
00576 !===> ... Get number of points to be sent
00577 !
00578          n_send = COUNT(send_mask)
00579          if (n_send == 0) cycle
00580 !
00581          mask_changed = .true.
00582 !
00583 !===> ... A grid of process "dest_comp" is a potential partner
00584 !
00585          dest = comp_info%psmile_ranks(dest_comp)
00586 !
00587          n_recv_req = n_recv_req + 1
00588 !
00589 #ifdef VERBOSE
00590          print 9970, trim(ch_id), rank, dest_comp, dest, n_send
00591 
00592          call psmile_flushstd
00593 #endif /* VERBOSE */
00594 !
00595 !===> ... Allocate send buffers
00596 !
00597 !  cf. subroutine PSMILe_Store_faces_3d_reg_real
00598 !
00599 !  len_rtem = Number of pieces of information
00600 !             per point in Real buffer "rbuf"
00601 !
00602 !  Real vector buf (1:ndrbuf):
00603 !  buf (:) = Coordinates of the point to be searched
00604 !            ndim_3d data items
00605 !  buf (:) = Maximal distance
00606 !            1 data item
00607 !
00608 !  Integer buffer ibuf (1:ndibuf):
00609 !
00610 !  ibuf (:) = Index in indices/dist_real/send_mask
00611 !
00612 !        len_rtem = ndim_3d + nb_extra ! coords + alle distancen
00613          len_rtem = ndim_3d + 1        ! coords + groesste Distance
00614 !
00615          len_rbuf = n_send * len_rtem
00616          len_ibuf = n_send * len_item
00617 !
00618          if (len_rbuf > ndrbuf) then
00619             if (ndrbuf > 0) then
00620                Deallocate (buf)
00621             endif
00622 !
00623             ndrbuf = len_rbuf
00624             Allocate (buf(ndrbuf), STAT = ierror)
00625 
00626             if ( ierror > 0 ) then
00627                ierrp (1) = ierror
00628                ierrp (2) = ndrbuf
00629 
00630                ierror = PRISM_Error_Alloc
00631                call psmile_error ( ierror, 'buf', &
00632                                    ierrp, 2, __FILE__, __LINE__ )
00633                return
00634             endif
00635          endif
00636 !
00637          if (len_ibuf > ndibuf) then
00638             if (ndibuf > 0) then
00639                Deallocate (ibuf)
00640             endif
00641 !
00642             ndibuf = len_ibuf
00643             Allocate (ibuf(ndibuf), STAT = ierror)
00644 
00645             if ( ierror > 0 ) then
00646                ierrp (1) = ierror
00647                ierrp (2) = ndibuf
00648 
00649                ierror = PRISM_Error_Alloc
00650                call psmile_error ( ierror, 'ibuf', &
00651                                    ierrp, 2, __FILE__, __LINE__ )
00652                return
00653             endif
00654          endif
00655 !
00656 !===> Store data in send buffers
00657 !     ??? compact machen !?!?
00658 !         Mehrfache Zellen Info's in ibuf und corner volumen
00659 !         eliminieren ?
00660 !
00661 #ifdef USE_PACK
00662          ibuf (1:n_send*len_item:len_item) = &
00663             PACK ((/ (i, i=1,n_extra) /), send_mask)
00664 #else
00665          ipi = 0
00666 
00667 !cdir vector
00668          do i = 1, n_extra
00669             if (send_mask (i)) then
00670                ibuf (ipi+1) = i
00671                ipi = ipi + len_item
00672             endif
00673          end do
00674 #endif /* USE_PACK */
00675 !
00676 !===> ... store coordinates and distance a single vectors 
00677 !         (to be consistent with storage of PRISM external coordinates)
00678 !         ? ist das wiklich besser, als die per Punkt zu speichern.
00679 !           Grund fuer die Feldweise Speicherung ist, die
00680 !           moegliche Wiederverwendung von Routinen.
00681 !
00682 ! ? besser direkt cos_search, sin_search, z versenden ?
00683 ! Das waeren 5 satt 3 items;
00684 ! Neu rechnen schneller als senden ?
00685 !              buf (ip+1) = cos_search (i, 1)
00686 !              buf (ip+2) = cos_search (i, 2)
00687 !              buf (ip+3) = sin_search (i, 1)
00688 !              buf (ip+4) = sin_search (i, 2)
00689 !              buf (ip+5) = z_search (i) = tgt_coords_z (indices(i))
00690 !
00691 #ifdef ALL_ITEM_PER_POINT
00692          ip = 0
00693 !cdir vector
00694          do i = 1, n_extra
00695             if (send_mask (i)) then
00696 !
00697                buf (ip+1) = tgt_coords_x (indices(i))
00698                buf (ip+2) = tgt_coords_y (indices(i))
00699                buf (ip+3) = tgt_coords_z (indices(i))
00700 !
00701                buf (ip+4) = dist_real (i, nb_extra)
00702                ip = ip + len_rtem
00703             endif
00704          end do
00705 #else
00706          ip = 0
00707 !cdir vector
00708          do i = 1, n_extra
00709             if (send_mask (i)) then
00710                ip = ip + 1
00711                buf (ip           ) = tgt_coords_x (indices(i))
00712                buf (ip+ n_send   ) = tgt_coords_y (indices(i))
00713                buf (ip+(n_send*2)) = tgt_coords_z (indices(i))
00714                buf (ip+(n_send*3)) = dist_real (i, nb_extra)
00715             endif
00716          end do
00717 #endif
00718 !
00719 #ifdef TODO
00720 ! TODO: If dest_comp == rank, direct call of psmile_search_donor_extra()
00721          if (dest_comp == rank) then
00722 !
00723             search_global%sender = dest
00724             search_global%msg_extra = extra_msg
00725 !     search_global%len_msg = nd_msg
00726 !
00727             call psmile_search_donor_extra (search_global, tol, ierror)
00728             if (ierror > 0) return
00729 
00730          endif
00731 #endif /* TODO */
00732 !
00733 !===> ... Send request to destination process
00734 !         extra_msg is received by routine PSMILe_Enddef_extra and
00735 !         the other buffer are received by routine
00736 !         PSMILe_Enddef_action_extra.
00737 !
00738          extra_msg%reqest_type         = PSMILe_nnghbr3D
00739          extra_msg%datatype            = PRISM_REAL
00740          extra_msg%len_int_data        = len_ibuf
00741          extra_msg%len_coord_data      = len_rbuf
00742          extra_msg%global_comp_id      = comp_info%global_comp_id
00743          extra_msg%transi_out_id       = search%msg_intersections%field_info%transient_out_id
00744 !
00745          extra_msg%num_volumes         = n_send
00746          extra_msg%num_int_per_vol     = len_item
00747          extra_msg%num_items_per_coord = len_rtem
00748 !
00749 !        ... No partition info
00750 !        Partition info senden; dann koennte man srcloc ausnutzen
00751 !
00752          extra_msg%partition_avail = .false.
00753 !
00754          extra_msg%idx_req       = n_recv_req
00755          extra_msg%num_neigh     = nb_extra
00756 
00757          extra_msg%local_grid_id = comp_info%all_extent_infos(igrid)%local_grid_id
00758 !
00759 #ifdef DEBUGX
00760          do i = 1, n_send
00761             print *, ibuf ((i-1)*len_item+1:i*len_item)
00762          end do
00763 #endif
00764 !
00765          call psmile_pack_msg_extra (extra_msg, extra_msg_buf)
00766 
00767          call psmile_bsend (extra_msg_buf, msg_extra_size, MPI_INTEGER, &
00768                             dest, exttag, comm_psmile, ierror)
00769          if (ierror /= MPI_SUCCESS) then
00770             ierrp (1) = ierror
00771             ierrp (2) = dest
00772             ierrp (3) = exttag
00773 
00774             ierror = PRISM_Error_Send
00775 
00776             call psmile_error (ierror, 'psmile_bsend(msg)', &
00777                                ierrp, 3, __FILE__, __LINE__ )
00778             return
00779          endif
00780 !
00781          call psmile_bsend (ibuf, len_ibuf, MPI_INTEGER, &
00782                             dest, exttag, comm_psmile, ierror)
00783          if (ierror /= MPI_SUCCESS) then
00784             ierrp (1) = ierror
00785             ierrp (2) = dest
00786             ierrp (3) = exttag
00787 
00788             ierror = PRISM_Error_Send
00789 
00790             call psmile_error (ierror, 'psmile_bsend(ibuf)', &
00791                                ierrp, 3, __FILE__, __LINE__ )
00792             return
00793          endif
00794 !
00795          call psmile_bsend (buf, len_rbuf, MPI_REAL, &
00796                             dest, exttag, comm_psmile, ierror)
00797          if (ierror /= MPI_SUCCESS) then
00798             ierrp (1) = ierror
00799             ierrp (2) = dest
00800             ierrp (3) = exttag
00801 
00802             ierror = PRISM_Error_Send
00803 
00804             call psmile_error (ierror, 'psmile_bsend(buf)', &
00805                                ierrp, 3, __FILE__, __LINE__ )
00806             return
00807          endif
00808 !
00809 !===> ... Setup request for answer of destination process
00810 !         ??? wirklich; sollte wohl von dem normalen abgehandelt werden
00811 !
00812          if (n_recv_req == 1) then
00813             call MPI_Irecv (answer, msg_extra_size, MPI_INTEGER, MPI_ANY_SOURCE, &
00814                             rexttag, comm_psmile, recv_req, ierror)
00815             if (ierror /= MPI_SUCCESS) then
00816 
00817                ierrp (1) = ierror
00818                ierrp (2) = dest
00819                ierrp (3) = rexttag
00820 
00821                ierror = PRISM_Error_Recv
00822 
00823                call psmile_error ( ierror, 'MPI_Irecv', &
00824                                    ierrp, 3, __FILE__, __LINE__ )
00825                return
00826 
00827             endif
00828          endif
00829 !
00830       end do ! igrid
00831 !
00832 !===> Deallocate
00833 !
00834       Deallocate (send_mask, boxes, igrid_to_dest_comp)
00835 !
00836       if (ndrbuf > 0) Deallocate (buf)
00837       if (ndibuf > 0) Deallocate (ibuf)
00838 !
00839 #ifdef __BUG__
00840       if (n_recv_req == 0) then
00841 
00842 !rr      send_info%nrecv and send_info%num2recv already may contain the results
00843 !rr      from the global search for filling the interpolation stencil. Here we
00844 !rr      perform the extra search for points for which no unmasked source points
00845 !rr      have been found. It must not be set to zero here as it is later used
00846 !rr      in psmile_compact_neighbours to compact points from the global search.
00847 !rr      send_info%nrecv and send_info%num2recv get updated in psmile_select_nn_found
00848 !rr      which is called from this routine.
00849 
00850          send_info%nrecv    = 0
00851          send_info%num2recv = 0
00852 
00853          go to 1000
00854       endif
00855 #endif
00856 
00857       if (n_recv_req == 0) go to 1000
00858 
00859 !
00860 !-------------------------------------------------------------------------
00861 !     Wait for answers
00862 !      Note: (2) Message from a partner with an intersection is disabled
00863 !                since request for grid receive (lrequest (3)) is setup
00864 !                if such a message is found.
00865 !            ??? Koennte man enablen wenn lrequest(3) == MPI_REQUEST_NULL
00866 !            (3) Old receive of grid data (paction%lrequest (3)) is disabled
00867 !-------------------------------------------------------------------------
00868 !
00869 ! nrecv      = Number of significant messages received
00870 !
00871       nrecv = 0
00872 !
00873       Allocate (sel_info (n_recv_req), stat = ierror)
00874       if ( ierror > 0 ) then
00875          ierrp (1) = ierror
00876          ierrp (2) = n_recv_req
00877 
00878          ierror = PRISM_Error_Alloc
00879          call psmile_error ( ierror, 'sel_info', ierrp, 2, &
00880                              __FILE__, __LINE__ )
00881          return
00882       endif
00883 !
00884 !===> Save old values of requests (2:3),
00885 !     ignore index 2 () and
00886 !     set index 3 to first receive of returned meta data
00887 !
00888       save_lreq (2:3) = paction%lrequest (2:3)
00889       paction%lrequest (2) = MPI_REQUEST_NULL
00890 !
00891       do n = 1, n_recv_req
00892          paction%lrequest (3) = recv_req
00893 !
00894          index = 0
00895          do while (index /= 3) 
00896 #ifdef DEBUG
00897             print *, trim(ch_id), paction%nreq, recv_req
00898             call psmile_flushstd
00899 #endif
00900 !
00901             call MPI_Waitany (paction%nreq, paction%lrequest, index, status, ierror)
00902 
00903             if ( ierror /= MPI_SUCCESS ) then
00904                ierrp (1) = ierror
00905                ierrp (2) = status (MPI_SOURCE)
00906                ierrp (3) = status (MPI_TAG)
00907 
00908                ierror = PRISM_Error_MPI
00909 
00910                call psmile_error ( ierror, 'MPI_Waitany', &
00911                                  ierrp, 3, __FILE__, __LINE__ )
00912                return
00913             endif
00914 
00915 #ifdef PRISM_ASSERTION
00916             if (index == MPI_UNDEFINED) then
00917                call psmile_assert ( __FILE__, __LINE__, &
00918                                     'request list is empty')
00919             endif
00920 #endif
00921 !
00922             if (index /= 3) then
00923                call psmile_enddef_action (search, index, status, ierror)
00924                if (ierror > 0) return
00925             endif
00926          end do ! while
00927 !
00928 !===> ... Get 
00929 !
00930          sender = status (MPI_SOURCE)
00931          len_ibuf = answer (3)
00932          len_rbuf = answer (4)
00933 !
00934 #ifdef VERBOSE
00935          print 9960, trim(ch_id), sender, len_ibuf, len_rbuf
00936 
00937          call psmile_flushstd
00938 #endif /* VERBOSE */
00939 !
00940 #ifdef PRISM_ASSERTION
00941          if (len_ibuf < 0) then
00942             print *, 'len_ibuf =', len_ibuf
00943             call psmile_assert (__FILE__, __LINE__, &
00944                  "len_ibuf should be >= 0")
00945          endif
00946 !
00947          if (len_rbuf < 0) then
00948             print *, 'len_rbuf =', len_rbuf
00949             call psmile_assert (__FILE__, __LINE__, &
00950                  "len_rbuf should be >= 0")
00951          endif
00952 #endif
00953 !
00954          if (len_ibuf > 0) then
00955             Allocate (ibuf (1:len_ibuf), stat = ierror)
00956 
00957             if ( ierror > 0 ) then
00958                ierrp (1) = ierror
00959                ierrp (2) = len_ibuf
00960 
00961                ierror = PRISM_Error_Alloc
00962                call psmile_error ( ierror, 'ibuf', ierrp, 2, &
00963                                    __FILE__, __LINE__ )
00964                return
00965             endif
00966 !
00967             call MPI_Recv (ibuf, len_ibuf, MPI_INTEGER, sender, &
00968                            rexttag, comm_psmile, status, ierror)
00969             if (ierror /= MPI_SUCCESS) then
00970 
00971                ierrp (1) = ierror
00972                ierrp (2) = sender
00973                ierrp (3) = rexttag
00974 
00975                ierror = PRISM_Error_Recv
00976 
00977                call psmile_error ( ierror, 'MPI_Recv(ibuf)', ierrp, 3, &
00978                                    __FILE__, __LINE__ )
00979                return
00980             endif
00981 !
00982             Allocate (buf (1:len_rbuf), stat = ierror)
00983 
00984             if ( ierror > 0 ) then
00985                ierrp (1) = ierror
00986                ierrp (2) = len_rbuf
00987 
00988                ierror = PRISM_Error_Alloc
00989                call psmile_error ( ierror, 'buf', ierrp, 2, &
00990                                    __FILE__, __LINE__ )
00991                return
00992             endif
00993 !
00994             call MPI_Recv (buf, len_rbuf, MPI_REAL, sender, &
00995                            rexttag, comm_psmile, status, ierror)
00996             if (ierror /= MPI_SUCCESS) then
00997 
00998                ierrp (1) = ierror
00999                ierrp (2) = sender
01000                ierrp (3) = rexttag
01001 
01002                ierror = PRISM_Error_Recv
01003 
01004                call psmile_error ( ierror, 'MPI_Recv(rbuf)', ierrp, 3, &
01005                                    __FILE__, __LINE__ )
01006                return
01007             endif
01008 !
01009 !===> ...
01010 !
01011 ! n_send  = Number of cells to be searched
01012 ! n_found = Number of interpolation bases found
01013 ! n_liste = Number of data items to be used for interpolation
01014 !           n_liste <= n_found since a data item can be used
01015 !           multiply for interpolation.
01016 !
01017             n_send  = answer (7)
01018             n_liste = answer (8)
01019             n_found = answer (9)
01020             ireq    = answer (15)
01021 !
01022 #ifdef PRISM_ASSERTION
01023             if (answer (1) /= PSMILe_nnghbr3D) then
01024                print *, 'answer(1)', answer(1), PSMILe_nnghbr3D
01025                call psmile_assert (__FILE__, __LINE__, &
01026                     "expected nearest neighbour interpolation")
01027             endif
01028 !
01029             if (ireq < 1 .or. ireq > n_recv_req) then
01030                print *, 'ireq, n_recv_req', ireq, n_recv_req
01031                call psmile_assert (__FILE__, __LINE__, &
01032                     "ireq must be in range of 1:n_recv_req")
01033             endif
01034 #endif
01035 !
01036             nrecv = nrecv + 1
01037 !
01038             sel_info(nrecv)%sender    = sender
01039             sel_info(nrecv)%n_liste   = n_liste
01040             sel_info(nrecv)%index     = answer (10)
01041             sel_info(nrecv)%method_id = answer (11)
01042             sel_info(nrecv)%msg_id    = answer (16)
01043 
01044             sel_info(nrecv)%real_buf => buf
01045 !
01046 !===> ... Allocate and initialize array "selected"
01047 !
01048 ! selected = Temporary vector for global nearest neighbour search
01049 !            for extra points with
01050 !            selected (1, *) = Number (index) of the receive
01051 !                              from which the nearest neigbour was
01052 !                              selected.
01053 !            selected (2, *) = Index in list of points 
01054 !                              (sent by neighbouring process)
01055 !            (cf. routine PSMILe_Add_nn_found_real)
01056 !            Dimension : selected (2, n_extra)
01057 !
01058 !
01059             if (.not. Allocated (selected)) then
01060                Allocate (selected(2, n_extra), stat = ierror)
01061 !
01062                if (ierror /= 0) then
01063                   ierrp (1) = n_extra * 2
01064 
01065                   ierror = PRISM_Error_Alloc
01066 
01067                   call psmile_error ( ierror, 'selected', ierrp, 1, &
01068                               __FILE__, __LINE__ )
01069                   return
01070                endif
01071 !
01072                selected (1, :) = 0
01073             endif
01074 !
01075 !===> ... Add information on points found by global search 
01076 !
01077             ip_dist = n_liste * ndim_3d
01078 !
01079             call psmile_add_nn_found_real (search, extra_search,     &
01080                         ibuf (1:n_send),                             &
01081                         ibuf   (n_send+1:2*n_send), n_send,          &
01082                         ibuf (2*n_send+1:2*n_send+n_found),          &
01083                         buf (ip_dist+1:ip_dist+n_found), n_found,    &
01084                         nb_extra, selected, sel_info, nrecv, ierror)
01085             if (ierror > 0) return
01086 !
01087             Deallocate (ibuf)
01088          endif
01089 !
01090 !===> ... Setup new request
01091 !
01092          if (n < n_recv_req) then
01093             call MPI_Irecv (answer, msg_extra_size, MPI_INTEGER, MPI_ANY_SOURCE, &
01094                               rexttag, comm_psmile, recv_req, &
01095                               ierror)
01096             if (ierror /= MPI_SUCCESS) then
01097 
01098                ierrp (1) = ierror
01099                ierrp (2) = dest
01100                ierrp (3) = rexttag
01101 
01102                ierror = PRISM_Error_Recv
01103 
01104                call psmile_error ( ierror, 'MPI_Irecv', &
01105                                     ierrp, 3, __FILE__, __LINE__ )
01106                return
01107             endif
01108          endif
01109       end do
01110 !
01111 !===> Determine which points are selected and return info to
01112 !     neighbouring processes
01113 !
01114       if (nrecv > 0) then
01115          call psmile_select_nn_found (search, extra_search,  &
01116                      send_info,                              &
01117                      selected, sel_info, nrecv, nb_extra,    &
01118                      neighbors_3d, nloc, num_neigh,          &
01119                      Grids(grid_id)%grid_shape, ierror)
01120          if (ierror > 0) return
01121       endif
01122 !
01123 !    ... Free memory
01124 !
01125       if (Allocated (selected)) then
01126          Deallocate (selected)
01127       endif
01128 !
01129       Deallocate (sel_info)
01130 !
01131 !    ... Restore original requests (2:3)
01132 !
01133 #ifdef PRISM_ASSERTION
01134       if (paction%lrequest (2) /= MPI_REQUEST_NULL .or. &
01135           paction%lrequest (3) /= MPI_REQUEST_NULL) then
01136          print *, 'request: ', paction%lrequest (2:3)
01137          call psmile_assert ( __FILE__, __LINE__, &
01138                              'Illegal request stored')
01139 
01140       endif
01141 #endif
01142 !
01143       paction%lrequest (2:3) = save_lreq (2:3)
01144 !
01145 !-------------------------------------------------------------------------
01146 !===> All done
01147 !-------------------------------------------------------------------------
01148 !
01149 1000  continue
01150 #ifdef VERBOSE
01151       print 9980, trim(ch_id), ierror
01152 
01153       call psmile_flushstd
01154 #endif /* VERBOSE */
01155 !
01156 !  Formats:
01157 !
01158 
01159 #ifdef VERBOSE
01160 
01161 9990 format (1x, a, ': psmile_global_search_nnx_real: var_id', i3, &
01162                     ' to ', i3, '(', i2, ')')
01163 9980 format (1x, a, ': psmile_global_search_nnx_real: eof ierror =', i3)
01164 9970 format (1x, a, ': psmile_global_search_nnx_real: send from', i3, &
01165                     ' to', i3, '[', i3, '], n_send =', i6)
01166 
01167 9960 format (1x, a, ': psmile_global_search_nnx_real: got rexttag message:', &
01168                     ' sender ', i4, ', len_ibuf, len_rbuf', 2i8)
01169 9950 format (1x, a, ': psmile_global_search_nnx_real: before waitany :', &
01170                     'nreq =', i4, ', recv_req ', i8)
01171 #endif /* VERBOSE */
01172 
01173 #ifdef DEBUG
01174 #endif
01175 
01176       end subroutine PSMILe_global_search_nnx_real

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1