psmile_add_points_found.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_Add_points_found
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_add_points_found (grid_id, search,   &
00012                            extra_search, indices_returned,   &
00013                            found, n_send, len_nsend,         &
00014                            index_found, n_found,             &
00015                            neighbors_3d, nloc, num_neigh,    &
00016                            grid_valid_shape, use_how, ierror)
00017 !
00018 ! !USES:
00019 !
00020       use PRISM
00021 !
00022       use PSMILe, dummy_interface => PSMILe_Add_points_found
00023 #ifdef DEBUG_TRACE
00024       use psmile_debug_trace
00025 #endif
00026 
00027       Implicit none
00028 !
00029 ! !INPUT PARAMETERS:
00030 !
00031       Integer,           Intent (In)     :: grid_id
00032 !
00033 !     Handle to the pointer to the receive list
00034 !
00035       Type (Enddef_search),     Intent (InOut) :: search
00036 
00037 !     Info's on coordinates to be searched
00038 !
00039       Integer,           Intent (In)     :: n_send
00040 !
00041 !     Number of points for which additional interpolation data was found.
00042 !     Dimension of "found" and "indices_returned".
00043 
00044       Integer,            Intent (In)    :: indices_returned (n_send)
00045 !
00046 !     Indices of the extra points in vectors "required" and
00047 !     "indices_req" for which interpolation data was found.
00048 !     This is a local index of the partition "ipart".
00049 !     (cf. routine PSMILe_Store_faces_3d_reg_dble)
00050 
00051       Integer,            Intent (In)    :: n_found
00052 !
00053 !     Total number of points found which are used for interpolation.
00054 !     Dimension of "index_found".
00055 !
00056       Integer,            Intent (In)    :: index_found (n_found)
00057 !
00058 !     Index of point in compact list "liste"
00059 !
00060       Integer,            Intent (In)    :: len_nsend (search%search_data%npart)
00061 !
00062 !     Number of points for which additional interpolation data was required
00063 !     per partition.
00064 !
00065       Integer,            Intent (In)    :: nloc
00066 !
00067 !     Total  Number of locations to be transferred
00068 !     Last dimension of neighbors array "neighbors_3d" and
00069 !     number of neighbors to be searched.
00070 
00071       Integer,            Intent (In)    :: num_neigh
00072 !
00073 !     Last dimension of neighbors array "neighbors_3d" and
00074 !     number of neighbors to be searched.
00075 
00076       Integer,            Intent (In)    :: grid_valid_shape (2, ndim_3d)
00077 
00078 !     Specifies the valid block shape for the "ndim_3d"-dimensional block
00079 
00080       Integer,            Intent (In)    :: use_how
00081 
00082 !     Information about how to apply the mask
00083 
00084 ! !INPUT/OUTPUT PARAMETERS:
00085 !
00086       Type (Extra_search_info), Intent (InOut) :: extra_search
00087 !
00088 !     Number of locations where
00089 !       (*) required mask values were not true
00090 
00091       Integer,            Intent (InOut) :: found    (n_send)
00092 !
00093 !     Input: Code of the interpolation bases found for corresponding point.
00094 !            (cf. "required")
00095 !     Output: Code cleaned by points masked out.
00096 
00097       Integer,            Intent (InOut) :: neighbors_3d (ndim_3d, nloc, 
00098                                                           num_neigh)
00099 !
00100 ! !OUTPUT PARAMETERS:
00101 !
00102       Integer,            Intent (Out)   :: ierror
00103 
00104 !     Returns the error code of PSMILe_Add_points_found;
00105 !             ierror = 0 : No error
00106 !             ierror > 0 : Severe error
00107 !
00108 ! !LOCAL PARAMETERS:
00109 !
00110 ! masked_out   = Hash value for a point which is masked out
00111 !
00112       Integer, Parameter           :: masked_out = 0
00113       Integer, Parameter           :: undef      = -12
00114 !
00115 ! !LOCAL VARIABLES
00116 !
00117 !     ... loop variables
00118 !
00119       Integer                         :: i, k, n
00120 !
00121 !     ... Bit value for found points
00122 !
00123       Integer                         :: iband (n_send)
00124       Integer                         :: code
00125 !
00126 !     ... 
00127 !
00128       Integer, Allocatable            :: corner2ind (:, :)
00129       Integer                         :: n_corners
00130       Integer                         :: outside
00131 !
00132 !     ... for partitions
00133 !
00134       Integer                         :: ipart, ip
00135       Integer, Pointer                :: indices_req (:)
00136       Integer, Pointer                :: required (:)
00137       Integer, Pointer                :: len_req (:)
00138 !
00139 !     ... for error handling
00140 !
00141       Integer, Parameter              :: nerrp = 2
00142       Integer                         :: ierrp (nerrp)
00143 !
00144 ! !DESCRIPTION:
00145 !
00146 ! Subroutine "PSMILe_Add_points_found" add information on
00147 ! points found by global search to the already existing information.
00148 !
00149 ! !REVISION HISTORY:
00150 !
00151 !   Date      Programmer   Description
00152 ! ----------  ----------   -----------
00153 ! 02.02.06    H. Ritzdorf  created
00154 !
00155 !EOP
00156 !----------------------------------------------------------------------
00157 !
00158 !  $Id: psmile_add_points_found.F90 3174 2011-05-04 09:41:43Z hanke $
00159 !  $Author: hanke $
00160 !
00161    Character(len=len_cvs_string), save :: mycvs = 
00162        '$Id: psmile_add_points_found.F90 3174 2011-05-04 09:41:43Z hanke $'
00163 !
00164 !----------------------------------------------------------------------
00165 !
00166 !  Initialization
00167 !
00168 #ifdef VERBOSE
00169       print 9990, trim(ch_id), n_send, n_found
00170 
00171       call psmile_flushstd
00172 #endif /* VERBOSE */
00173 !
00174       ierror = 0
00175 !
00176       if (associated(Grids(grid_id)%reduced_gauss_data%remote_index)) then
00177           !
00178           ! in particular for Gauss-reduced grids
00179           !
00180           outside   = grid_valid_shape(2,1) + &
00181                       size(Grids(grid_id)%reduced_gauss_data%remote_index) + 11
00182       else
00183           outside   = grid_valid_shape(2,1) + 11
00184       endif
00185 
00186       n_corners = num_neigh
00187 !
00188       len_req     => extra_search%len_req
00189 !
00190 #ifdef PRISM_ASSERTION
00191       do i = 1, n_found
00192          if (index_found (i) <= 0 .and. &
00193              index_found (i) /= masked_out) exit
00194       end do
00195 !
00196       if (i <= n_found) then
00197          print *, 'i, n_found, index_found', i, n_found, index_found (i)
00198          call psmile_assert ( __FILE__, __LINE__, &
00199                  'invalid index in found list')
00200       endif
00201 
00202       ip = 0
00203       do ipart = 1, search%search_data%npart
00204          do i = 1, len_nsend (ipart)
00205             if (indices_returned (ip+i) < 1 .or. &
00206                 indices_returned (ip+i) > len_req(ipart)) exit
00207          end do
00208 !
00209          if (i <= len_nsend (ipart)) then
00210             print *, 'ipart, ip, i, nreq, indices_returned(ip+i)', &
00211                      ipart, ip, i, len_req(ipart),                 &
00212                      indices_returned (ip+i)
00213             call psmile_assert ( __FILE__, __LINE__, &
00214                      'invalid index in indices_returned')
00215          endif
00216 !
00217          ip = ip + len_nsend (ipart)
00218       end do
00219 #endif
00220 !
00221 #ifdef DEBUGX
00222       do i = 1, n_found
00223          print *, 'i index_found', i, index_found(i)
00224       end do
00225 
00226       code = 1
00227       k = 0
00228       do i = 1, n_corners
00229 
00230          do n = 1, n_send
00231             if (IAND (found(n), code) /= 0) then
00232                k = k + 1
00233                print *, 'found: n', n, code, index_found (k)
00234             endif
00235          end do
00236          code = code * 2
00237       end do ! i
00238 
00239       ip = 0
00240       do ipart = 1, search%search_data%npart
00241          do i = 1, len_nsend (ipart)
00242                print *, 'ip i indices_returned', ip, i, indices_returned (ip+i), found (ip+i)
00243 ! indices_req(indices_returned (ip+i))
00244          end do
00245             ip = ip + len_nsend (ipart)
00246       end do
00247 #endif
00248 !
00249 !===> Allocate array in order to store the indices of
00250 !     of the compact list of "interpolation bases found"
00251 !     (cf. PSMILe_Trili_3d_extra_off)
00252 !
00253       Allocate (corner2ind (n_send, n_corners), stat = ierror)
00254       if ( ierror > 0 ) then
00255          ierrp (1) = ierror
00256          ierrp (2) = n_send * n_corners
00257 
00258          ierror = PRISM_Error_Alloc
00259          call psmile_error ( ierror, 'corner2ind', &
00260                              ierrp, 2, __FILE__, __LINE__ )
00261          return
00262       endif
00263 !
00264       corner2ind (:, :) = undef
00265 !
00266 !===> Decode mask "found" of points found
00267 !     (cf. PSMILe_Trili_3d_extra_off)
00268 !
00269       code = 1
00270       k = 0
00271       do i = 1, n_corners
00272 
00273          do n = 1, n_send
00274             if (IAND (found(n), code) /= 0) then
00275 !
00276 !===> ...... Required "i-th" interpolation base was found.
00277 !            Corresponding value is stored in position "index_found(k)".
00278 !            index_found (k) == masked_out : Point is masked out
00279 !
00280                k = k + 1
00281                corner2ind (n, i) = index_found (k)
00282 !
00283             endif
00284          end do
00285 !
00286 !===> ... next interpolation base
00287 !
00288          code = code * 2
00289       end do ! i
00290 !
00291 #ifdef PRISM_ASSERTION
00292       if (k /= n_found) then
00293          print *, 'k, n_found', k, n_found
00294          call psmile_assert ( __FILE__, __LINE__, &
00295                  'k /= n_found')
00296       endif
00297 #endif
00298 !
00299 !===> Remove points from code "found" which are masked out
00300 !
00301       if (use_how /= PSMILE_novalue) then
00302 !
00303          code = 1
00304          do i = 1, n_corners
00305             do n = 1, n_send
00306                if (corner2ind (n,i) == masked_out) then
00307                   found (n) = found (n) - code
00308                endif
00309             end do
00310             code = code * 2
00311          end do ! i
00312       endif
00313 !
00314 !===> For each part
00315 !
00316 !     indices_req = Indices of the extra points in entire (!) list of
00317 !                   all points to be searched.
00318 !     required    = Code for the extra points
00319 !                   which are additionally required.
00320 !
00321       ip = 0
00322       do ipart = 1, search%search_data%npart
00323          if (len_req (ipart) > 0 .and. &
00324              len_nsend (ipart) > 0) then
00325 !
00326             indices_req => extra_search%indices_req(ipart)%vector
00327             required    => extra_search%required   (ipart)%vector
00328 !
00329 !===> ... Apply rules for points masked out
00330 !
00331             if (use_how == PSMILE_novalue) then
00332 !
00333 !                 If one neighbor is invalid reject all 
00334 !                 by moving all neighbors out of the valid range.
00335 
00336                do i = 1, n_corners
00337 !cdir vector
00338                   do n = 1, len_nsend (ipart)
00339                      if (corner2ind (ip+n,i) == masked_out) then
00340 !                       This point is not required any more
00341                         required (indices_returned (ip+n)) = 0
00342                         neighbors_3d (1, indices_req(indices_returned (ip+n)), &
00343                                       1:n_corners) = outside
00344                      endif
00345                   end do
00346                end do
00347             else 
00348 !
00349 ! HR: This should be work also for "use_how == PSMILE_novalue"
00350 !     "psmile_neigh_extra_points" should do the correct things in this case.
00351 !     TODO: TEST whether this statement is really correct.
00352 !
00353                do i = 1, n_corners
00354 !cdir vector
00355                   do n = 1, len_nsend (ipart)
00356 #ifdef DEBUG
00357                         print *, 'indices_returned', ip+n, indices_returned (ip+n), indices_req(indices_returned (ip+n))
00358                         print *, 'corner info', i, corner2ind (ip+n,i)
00359 #endif
00360                      if (corner2ind (ip+n,i) == masked_out) then
00361                         neighbors_3d (1, indices_req(indices_returned (ip+n)), &
00362                                       i) = outside
00363 #ifdef DEBUG_TRACE
00364                         if (indices_req(indices_returned (ip+n)) == ictl) then
00365                            print *, '#### psmile_add_points_found outside: ', outside
00366                            print *, '#### psmile_add_points_found ictl: set outside for ictl', ictl
00367                            print *, '#### psmile_add_points_found corner ', i
00368                            print *, ' -> neigh 1 ', neighbors_3d (1,indices_req(indices_returned (ip+n)),:)
00369                            print *, ' -> neigh 2 ', neighbors_3d (2,indices_req(indices_returned (ip+n)),:)
00370                            print *, ' -> neigh 3 ', neighbors_3d (3,indices_req(indices_returned (ip+n)),:)
00371                         endif
00372 #endif
00373                      endif
00374                   end do
00375                end do
00376             endif
00377 !
00378 !===> ... Get common points required by "found" and "required" and
00379 !         remove them from the coded list.
00380 !
00381 !cdir vector
00382             do n = 1, len_nsend (ipart)
00383                iband (n) = IAND (found (ip+n), required(indices_returned (ip+n)))
00384             end do
00385 !
00386 !cdir vector
00387             do n = 1, len_nsend (ipart)
00388                required(indices_returned (ip+n)) = &
00389                required(indices_returned (ip+n)) - iband (n)
00390             end do
00391 !
00392 !===> ... Set
00393 !         neighbors_3d (1, i) = "extra_search%global_marker" in order to mark
00394 !         that this is a point out of list of globally searched bases.
00395 !         neighbors_3d (2, i) = Index in list of globally searched bases.
00396 !
00397             code = 1
00398             do i = 1, n_corners
00399 !
00400 !cdir vector
00401                do n = 1, len_nsend (ipart)
00402                   if (IAND (iband(n), code) /= 0) then
00403 #ifdef DEBUG_TRACE
00404                      if (indices_req(indices_returned (ip+n)) == ictl) then
00405                         print *, '#### ictl:', ictl, ' set global marker and', corner2ind(ip+n, i)
00406                         print *, 'i, ip, indices_returned', &
00407                                   i, ip, indices_returned (ip+n), iband (n)
00408                      endif
00409 #endif
00410                      neighbors_3d (1,                              & 
00411                         indices_req(indices_returned (ip+n)), i) = &
00412                         extra_search%global_marker
00413 !
00414                      neighbors_3d (2,                              & 
00415                         indices_req(indices_returned (ip+n)), i) = &
00416                         corner2ind(ip+n, i)
00417                   else
00418 #ifdef DEBUG_TRACE
00419                      if (indices_req(indices_returned (ip+n)) == ictl) then
00420                         print *, '#### ictl:', ictl, ' global marker not set'
00421                         print *, 'i, ip, indices_returned', &
00422                                   i, ip, indices_returned (ip+n), iband (n), code
00423                      endif
00424 #endif
00425                   endif
00426                end do
00427 !
00428 !===> ... next interpolation base
00429 !
00430                code = code * 2
00431             end do
00432 !
00433 #ifdef PRISM_ASSERTION
00434 assert1:    do i = 1, n_corners
00435 !cdir vector
00436                do n = 1, len_nsend (ipart)
00437                   if (neighbors_3d (1, indices_req(indices_returned (ip+n)), i) &
00438                       == extra_search%global_marker .and. &
00439                       neighbors_3d (2, indices_req(indices_returned (ip+n)), i) &
00440                       <= 0)  exit assert1
00441                end do
00442             end do assert1
00443 !
00444             if (i <= n_corners) then
00445                print *, 'ipart, i, n', ipart, i, n, corner2ind(ip+n, i), &
00446                                        found(ip+n), iband(n)
00447                call psmile_assert ( __FILE__, __LINE__, &
00448                         'invalid index for interpolation base generated')
00449             endif
00450 #endif
00451          endif
00452 !
00453          ip = ip + len_nsend (ipart)
00454       end do ! ipart 
00455 !
00456 !===> All done
00457 !
00458 #ifdef VERBOSE
00459       print 9980, trim(ch_id), ierror
00460 
00461       call psmile_flushstd
00462 #endif /* VERBOSE */
00463 !
00464 !  Formats:
00465 !
00466 #ifdef VERBOSE
00467 
00468 9990 format (1x, a, ': psmile_add_points_found: n_send ', i6, &
00469              '; n_found ', i6)
00470 9980 format (1x, a, ': psmile_add_points_found: eof ierror =', i3)
00471 
00472 #endif /* VERBOSE */
00473 
00474 #ifdef DEBUG
00475 #endif
00476 
00477       end subroutine PSMILe_Add_points_found

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1