psmile_search_donor_extra_nn.F90

Go to the documentation of this file.
00001 ! #include "psmile.prep"
00002 !-----------------------------------------------------------------------
00003 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00004 ! All rights reserved. Use is subject to OASIS4 license terms.
00005 !-----------------------------------------------------------------------
00006 !BOP
00007 !
00008 ! !ROUTINE: PSMILe_Search_donor_extra_nn
00009 !
00010 ! !INTERFACE:
00011 
00012       subroutine psmile_search_donor_extra_nn (comp_info, search, &
00013                         var_id, tol, ierror)
00014 !
00015 ! !USES:
00016 !
00017       use PRISM
00018 !
00019       use PSMILe, dummy_interface => PSMILe_Search_donor_extra_nn
00020 
00021       implicit none
00022 !
00023 ! !INPUT PARAMETERS:
00024 !
00025       Type (Enddef_comp),          Intent (In)    :: comp_info
00026 !
00027 !     Info on the component in which the donor cells
00028 !     should be searched.
00029 
00030       Integer,                     Intent (In)    :: var_id
00031 
00032 !     Handle to the grid function
00033 
00034       Double Precision,            Intent (In)    :: tol
00035 
00036 !     Absolute tolerance for search of "identical" points
00037 !     TOL >= 0.0
00038 
00039 !
00040 ! !INPUT/OUTPUT PARAMETERS:
00041 !
00042       Type (Enddef_global_search), Intent (InOut) :: search
00043 
00044 !     Data on the points to be searched
00045 !
00046 ! !OUTPUT PARAMETERS:
00047 !
00048       Integer,                     Intent (Out)   :: ierror
00049 
00050 !     Returns the error code of PSMILe_Search_donor_extra_nn;
00051 !             ierror = 0 : No error
00052 !             ierror > 0 : Severe error
00053 !
00054 ! !LOCAL VARIABLES
00055 !
00056 !     ... field pointer
00057 !
00058       Type (Gridfunction), Pointer  :: field
00059 !
00060       Integer                       :: method_id
00061 !
00062 !     ... Grid pointer
00063 !
00064       Type (Grid), Pointer          :: grid_info
00065       Integer                       :: datatype, grid_id
00066 !
00067 !     ... for points to be searched
00068 !
00069 ! len_item = Length of a data item per point in search%ibuf
00070 !
00071       Integer                       :: len_item, n_send
00072       Integer                       :: ip_dist, nb_extra
00073 !
00074       Integer                       :: ip (ndim_3d)
00075       Real                          :: rtol
00076 
00077 !     Integer                       :: shift (ndim_3d)
00078 !
00079 !     ... for masks
00080 !         dummy_mask_array is a dummy array for to be
00081 !         transferred to interpolation routines
00082 !         Note: The target attributes (see !rr) were removed
00083 !               since there problems with a compiler and array bound checking
00084 !
00085       Integer                       :: mask_id
00086       Logical                       :: src_mask_available
00087 !
00088       Logical, Target               :: dummy_mask_array (1)
00089 !rr   Integer, target               :: dummy_mask_shape (2, ndim_3d)
00090       Integer                       :: dummy_mask_shape (2, ndim_3d)
00091 !
00092       Logical, Pointer              :: mask_array (:) 
00093 !rr   Integer, Pointer              :: mask_shape (:, :)
00094       Integer                       :: mask_shape (2, ndim_3d)
00095 !
00096 !     ... for locations found
00097 !
00098       Integer                       :: i, n, nlocs
00099 !
00100       Integer, Allocatable          :: nfound (:)
00101       Integer, Pointer              :: locations (:, :), locs (:, :)
00102       Integer, Allocatable          :: hash (:)
00103 !
00104 !     ... for error parameters
00105 !
00106       Integer, Parameter            :: nerrp = 2
00107       Integer                       :: ierrp (nerrp)
00108 !
00109 ! !DESCRIPTION:
00110 !
00111 ! Subroutine "PSMILe_Search_donor_extra_nn" performs the additional (global)
00112 ! search for coordinates sent by the requesting process if nearest neigbour
00113 ! search is required.
00114 !
00115 ! !REVISION HISTORY:
00116 !
00117 !   Date      Programmer   Description
00118 ! ----------  ----------   -----------
00119 ! 23.10.06    H. Ritzdorf  created
00120 !
00121 !EOP
00122 !----------------------------------------------------------------------
00123 !
00124 ! $Id: psmile_search_donor_extra_nn.F90 2082 2009-10-21 13:31:19Z hanke $
00125 ! $Author: hanke $
00126 !
00127    Character(len=len_cvs_string), save :: mycvs = 
00128        '$Id: psmile_search_donor_extra_nn.F90 2082 2009-10-21 13:31:19Z hanke $'
00129 !----------------------------------------------------------------------
00130 !
00131 !  Initialization
00132 !
00133 #ifdef VERBOSE
00134       print 9990, trim(ch_id), comp_info%comp_id, search%sender
00135 
00136       call psmile_flushstd
00137 #endif /* VERBOSE */
00138 !
00139       ierror = 0
00140 !
00141       field     => Fields (var_id)
00142 !
00143       method_id = field%method_id
00144       mask_id   = field%mask_id
00145 !
00146       grid_id   = Methods(method_id)%grid_id
00147       grid_info => Grids (grid_id)
00148 !
00149       src_mask_available = mask_id /= PRISM_UNDEFINED
00150 !
00151       if (src_mask_available) then
00152          mask_array => Masks(mask_id)%mask_array
00153 !rr      mask_shape => Masks(mask_id)%mask_shape
00154          mask_shape =  Masks(mask_id)%mask_shape
00155       else
00156          mask_array => dummy_mask_array
00157 !rr      mask_shape => dummy_mask_shape
00158          mask_shape =  dummy_mask_shape
00159       endif
00160 !
00161 !  n_send   = Number of points to be searched
00162 !  len_item = Number of pieces of information in integer buffer "ibuf".
00163 !             #1  : Index in indices/dist_dble/send_mask of sending process
00164 !
00165       n_send   = search%msg_extra (8)
00166       len_item = search%msg_extra (9)
00167       nb_extra = search%msg_extra (16)
00168 !
00169 #ifdef PRISM_ASSERTION
00170 !
00171 !===> Internal control
00172 !
00173       if (search%msg_extra(1) /= PSMILe_nnghbr3D) then
00174           call psmile_assert ( __FILE__, __LINE__, &
00175                               'Interpolation method must be Nearest Neighbour')
00176       endif
00177 !
00178       if (search%msg_extra(10) < ndim_3d + 1 ) then
00179           call psmile_assert ( __FILE__, __LINE__, &
00180                               'expected at least len_rtem = 4')
00181       endif
00182 !
00183 #endif
00184 
00185 #ifdef VERBOSE
00186       print *, trim(ch_id), ': psmile_search_donor_extra_nn search%msg_extra(1:4)', &
00187                search%msg_extra(1:4)
00188 !
00189       call psmile_flushstd
00190 #endif /* VERBOSE */
00191 !
00192 #ifdef TODO_GLOBAL_COORDS
00193 !  Transfer global grid indices from remote grid
00194 !      into local  grid indices of  current grid
00195 !  Weiss nicht, ob sich das wirklich fuer nearest neigbour lohnt ?
00196 !
00197       shift (1:ndim_3d) = grid_info%grid_shape(1, 1:ndim_3d) - &
00198                           grid_info%partition (1, 1:ndim_3d) - 1
00199 !
00200 !cdir vector
00201           do i = 1, n_send
00202           search%ibuf ((i-1)*len_item+1) = &
00203           search%ibuf ((i-1)*len_item+1) + shift (1)
00204           search%ibuf ((i-1)*len_item+2) = &
00205           search%ibuf ((i-1)*len_item+2) + shift (2)
00206           search%ibuf ((i-1)*len_item+3) = &
00207           search%ibuf ((i-1)*len_item+3) + shift (3)
00208           end do
00209 !
00210 #ifdef DEBUGX
00211       print *, 'grid shape', grid_info%grid_shape
00212       print *, 'srcloc (1:3), indices (j), len_item', len_item
00213          do i = 1, n_send
00214          print 8900, search%ibuf ((i-1)*len_item+1:i*len_item)
00215          end do
00216 #endif
00217 #endif
00218 !
00219 !===> Allocations locations array
00220 !
00221       Allocate (nfound(n_send), locations (ndim_3d, n_send), STAT = ierror)
00222 !     TRACE_ALLOCATE ("nfound", nfound,    n_send)
00223 !     TRACE_ALLOCATE ("locations", locations, n_send*ndim_3d)
00224 
00225       if ( ierror > 0 ) then
00226          ierrp (1) = ierror
00227          ierrp (2) = (ndim_3d + 1) * n_send
00228 
00229          ierror = PRISM_Error_Alloc
00230          call psmile_error ( ierror, 'locations, nfound', &
00231                  ierrp, 2, __FILE__, __LINE__ )
00232          return
00233       endif
00234 !
00235       nfound = 0
00236 !
00237 !----------------------------------------------------------------------
00238 !    Get nearest neighbour
00239 !----------------------------------------------------------------------
00240 !
00241 !===> Pointer to (lon, lat, z, distance) in buffer
00242 !
00243       ip(1)   = 1
00244       ip(2)   = ip(1) + n_send
00245       ip(3)   = ip(2) + n_send
00246       ip_dist = ip(3) + n_send
00247 !
00248       datatype  = grid_info%corner_pointer%corner_datatype
00249 !
00250       if (datatype == MPI_REAL) then
00251 
00252          rtol = tol
00253 
00254          call psmile_search_donor_nnx_real (comp_info,  &
00255                  search, var_id,                        &
00256                  search%rbuf(ip(1):ip(1)+n_send-1),     &
00257                  search%rbuf(ip(2):ip(2)+n_send-1),     &
00258                  search%rbuf(ip(3):ip(3)+n_send-1),     &
00259                  search%rbuf(ip_dist:ip_dist+n_send-1), &
00260                  nfound, locations, n_send, nb_extra, rtol, ierror)
00261 
00262       else if (datatype == MPI_DOUBLE_PRECISION) then
00263 
00264          call psmile_search_donor_nnx_dble (comp_info,  &
00265                  search, var_id,                        &
00266                  search%dbuf(ip(1):ip(1)+n_send-1),     &
00267                  search%dbuf(ip(2):ip(2)+n_send-1),     &
00268                  search%dbuf(ip(3):ip(3)+n_send-1),     &
00269                  search%dbuf(ip_dist:ip_dist+n_send-1), &
00270                  nfound, locations, n_send, nb_extra, tol, ierror)
00271 
00272 #if defined ( PRISM_QUAD_TYPE )
00273       else if (datatype == MPI_REAL16) the
00274 #endif
00275       endif
00276 
00277       if (ierror > 0) return
00278 !
00279 !----------------------------------------------------------------------
00280 !    Get number of found locations
00281 !----------------------------------------------------------------------
00282 !
00283       nlocs = SUM (nfound)
00284 !
00285 #ifdef DEBUGX
00286       print *, 'psmile_search_donor_extra_nn: nlocs, n_send', nlocs, n_send
00287       do i = 1, n_send
00288       if (nfound(i) < 0 .or. nfound (i) > 1) then
00289          print *, 'i, nfound', i, nfound(i)
00290       endif
00291       end do
00292 #endif
00293 !
00294       if (nlocs > 0) then
00295          Allocate (hash (nlocs), stat = ierror)
00296          if (ierror /= 0) then
00297             ierrp (1) = nlocs
00298 
00299             ierror = PRISM_Error_Alloc
00300 
00301             call psmile_error ( ierror, 'hash', ierrp, 1, &
00302                     __FILE__, __LINE__ )
00303             return
00304          endif
00305 !
00306          if (nlocs < n_send) then
00307             Allocate (locs (ndim_3d, nlocs), stat = ierror)
00308             if (ierror /= 0) then
00309                ierrp (1) = ndim_3d * nlocs
00310 
00311                ierror = PRISM_Error_Alloc
00312 
00313                call psmile_error ( ierror, 'locs', ierrp, 1, &
00314                        __FILE__, __LINE__ )
00315                return
00316             endif
00317 !
00318             n = 0
00319 !cdir vector
00320                do i = 1, n_send
00321                if (nfound (i) > 0) then
00322                   n = n + 1
00323                   locs (1, n) = locations (1, i)
00324                   locs (2, n) = locations (2, i)
00325                   locs (3, n) = locations (3, i)
00326                endif
00327                end do
00328 
00329 #ifdef PRISM_ASSERTION
00330             if (n /= nlocs) then
00331                print *, 'nlocs, n', nlocs, n, n_send
00332                call psmile_assert ( __FILE__, __LINE__, &
00333                        'nlocs != n')
00334             endif
00335 #endif
00336 
00337             Deallocate (locations)
00338          else
00339 
00340             locs => locations
00341 
00342          endif
00343 !
00344 !===> Create compact list of points to be sent and 
00345 !     return hash value "hash" (index in fortran order) for each location
00346 !
00347          call psmile_hash_extra (search, locs, hash, nlocs,  &
00348                  mask_array, mask_shape, src_mask_available, &
00349                  grid_info%grid_shape, ierror)
00350          if (ierror > 0) return
00351 !
00352          Deallocate (locs, hash)
00353 !
00354       else
00355 !
00356 !===> ... No point found. Set values in "search"
00357 !
00358          search%n_liste = 0
00359          search%n_found = 0
00360 
00361          NULLIFY (search%neigh_3d)
00362          NULLIFY (search%index_found)
00363 
00364          Deallocate (locations)
00365 
00366       endif
00367 !
00368 !----------------------------------------------------------------------
00369 !    Return info to the destination process
00370 !----------------------------------------------------------------------
00371 !
00372       if (datatype == MPI_REAL) then
00373          call psmile_return_extra_off_real (comp_info, search, var_id,      &
00374                      nfound, search%rbuf(ip_dist:ip_dist+n_send-1), n_send, &
00375                      nb_extra, ierror)
00376 
00377       else if (datatype == MPI_DOUBLE_PRECISION) then
00378          call psmile_return_extra_off_dble (comp_info, search, var_id,      &
00379                      nfound, search%dbuf(ip_dist:ip_dist+n_send-1), n_send, &
00380                      nb_extra, ierror)
00381 
00382 #if defined ( PRISM_QUAD_TYPE )
00383       else if (datatype == MPI_REAL16) the
00384          call psmile_return_extra_off_quad (comp_info, search, var_id,      &
00385                      nfound, search%qbuf(ip_dist:ip_dist+n_send-1), n_send, &
00386                      nb_extra, ierror)
00387 #endif
00388       endif
00389 !
00390       Deallocate (nfound)
00391 !
00392 !===> All done
00393 !
00394 #ifdef VERBOSE
00395       print 9980, trim(ch_id), grid_id, search%sender, ierror
00396 
00397       call psmile_flushstd
00398 #endif /* VERBOSE */
00399 !
00400       return
00401 !
00402 !  Formats:
00403 !
00404 #ifdef VERBOSE
00405 9990 format (1x, a, ': psmile_search_donor_extra_nn: comp_id =', i3, &
00406                     '; sender =', i4)
00407 9980 format (1x, a, ': psmile_search_donor_extra_nn: grid id =', i3, &
00408                     '; eof sender =', i3, ', ierror =', i4)
00409 
00410 9970 format (1x, a, ': psmile_search_donor_extra_nn:', &
00411              1x, 'offset local  ', 3i7, &
00412             /1x, 'extent        ', 3i7, &
00413             /1x, 'offset remote ', 3i7)
00414 #endif /* VERBOSE */
00415 8900 format (1x, 3i6, 5(';', i6) )
00416 
00417       end subroutine PSMILe_Search_donor_extra_nn

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1