psmile_compact_neighbors_3d.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_Compact_neighbors_3d
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_compact_neighbors_3d ( &
00012                neighbors_3d, nloc, num_neigh,  &
00013                grid_valid_shape, extra_search, &
00014                send_info, neighbors, ierror)
00015 !
00016 ! !USES:
00017 !
00018       use PRISM_constants
00019       use psmile_common, only : extra_search_info
00020       use PSMILe, dummy_interface => PSMILe_Compact_neighbors_3d
00021 
00022       implicit none
00023 !
00024 ! !INPUT PARAMETERS:
00025 
00026       Integer, Intent (In)            :: nloc
00027 !
00028 !     Number of locations to be transferred
00029 !
00030       Integer, Intent (In)            :: num_neigh
00031 !
00032 !     Last dimension of neighbors array "neighbors_3d"
00033 !
00034       Integer, Intent (In)            :: neighbors_3d (ndim_3d, nloc, num_neigh)
00035 !
00036 !     Indices of neighbor locations (interpolation bases)
00037 !
00038       Integer, Intent (In)            :: grid_valid_shape (2, ndim_3d)
00039 
00040 !     Specifies the valid block shape for the "ndim_3d"-dimensional block
00041 
00042       Type(Extra_search_info), Intent(In) :: extra_search
00043 
00044 !     Contains inforation on how to identify points/cells which are found by global search
00045 !
00046 ! !INPUT/OUTPUT PARAMETERS:
00047 !
00048       Type (Send_information), Intent (InOut)  :: send_info
00049 !
00050 ! !OUTPUT PARAMETERS:
00051 !
00052       Integer, Intent (Out)           :: neighbors (nloc, num_neigh)
00053 !
00054 !     1D-Indices of neighbor locations (interpolation bases)
00055 !
00056       Integer, Intent (Out)           :: ierror
00057 
00058 !     Returns the error code of PSMILE_Compact_neighbors_3d;
00059 !             ierror = 0 : No error
00060 !             ierror > 0 : Severe error
00061 !
00062 ! !LOCAL VARIABLES
00063 !
00064 !     ... for compact list
00065 !
00066       Integer                         :: i, j, k, n, n_list
00067       Integer                         :: size (ndim_3d)
00068 !
00069       Integer, Allocatable            :: indices (:, :, :)
00070       Integer, Pointer                :: list_entries (:, :)
00071 !
00072 !     ... for error handling
00073 !
00074       Integer, parameter              :: nerrp = 2
00075       Integer                         :: ierrp (nerrp)
00076 !
00077 !
00078 ! !DESCRIPTION:
00079 !
00080 ! Subroutine "PSMILe_Compact_neighbors_3d" transforms the 3d-indices
00081 ! stored in "neighbours_3d" into a 1D-index of a compact list of 
00082 ! the interpolation bases to be required for interpolation.
00083 !
00084 !
00085 ! !REVISION HISTORY:
00086 !
00087 !   Date      Programmer   Description
00088 ! ----------  ----------   -----------
00089 ! 03.07.21    H. Ritzdorf  created
00090 !
00091 !EOP
00092 !----------------------------------------------------------------------
00093 !
00094 !  $Id: psmile_compact_neighbors_3d.F90 2695 2010-10-29 15:48:58Z hanke $
00095 !  $Author: hanke $
00096 !
00097    Character(len=len_cvs_string), save :: mycvs = 
00098        '$Id: psmile_compact_neighbors_3d.F90 2695 2010-10-29 15:48:58Z hanke $'
00099 !
00100 !----------------------------------------------------------------------
00101 !
00102 !  Initialization
00103 !
00104 #ifdef VERBOSE
00105       print 9990, trim(ch_id)
00106 
00107       call psmile_flushstd
00108 #endif /* VERBOSE */
00109 !
00110       ierror = 0
00111 !
00112       size (:) = grid_valid_shape(2,:) - grid_valid_shape(1,:) + 1
00113 !
00114 !===> Generate indices within list of interpolations bases
00115 !
00116       Allocate (indices (grid_valid_shape(1,1):grid_valid_shape(2,1), &
00117                          grid_valid_shape(1,2):grid_valid_shape(2,2), &
00118                          grid_valid_shape(1,3):grid_valid_shape(2,3)), &
00119                 STAT = ierror)
00120 !
00121       if (ierror > 0) then
00122          ierrp (1) = ierror
00123          ierrp (2) = PRODUCT (size)
00124 
00125          ierror = PRISM_Error_Alloc
00126 
00127          call psmile_error ( ierror, 'indices', &
00128                              ierrp, 2, __FILE__, __LINE__ )
00129          return
00130       endif
00131 !
00132       indices = 0
00133 !
00134 !===> Generate mask of valid points
00135 !
00136       do n = 1, num_neigh
00137 
00138 !cdir vector
00139          do i = 1, nloc
00140          if (grid_valid_shape(1,1) <= neighbors_3d (1,i,n) .and. &
00141              grid_valid_shape(2,1) >= neighbors_3d (1,i,n) .and. &
00142              grid_valid_shape(1,2) <= neighbors_3d (2,i,n) .and. &
00143              grid_valid_shape(2,2) >= neighbors_3d (2,i,n) .and. &
00144              grid_valid_shape(1,3) <= neighbors_3d (3,i,n) .and. &
00145              grid_valid_shape(2,3) >= neighbors_3d (3,i,n)) then
00146 
00147             neighbors (i, n) = 1
00148 
00149          else
00150 !
00151 !           ... no interpolation base available
00152 !
00153             neighbors (i, n) = 0
00154          endif
00155          end do ! i = 1, nloc
00156       end do ! n = 1, num_neigh
00157 !
00158 !===> Generate list
00159 !
00160       n_list = 0
00161 
00162       do n = 1, num_neigh
00163 !cdir vector
00164          do i = 1, nloc
00165          if (neighbors (i, n) > 0) then
00166             if (indices (neighbors_3d (1,i,n), &
00167                          neighbors_3d (2,i,n), &
00168                          neighbors_3d (3,i,n)) == 0) then
00169                n_list = n_list + 1
00170                indices (neighbors_3d (1,i,n), &
00171                         neighbors_3d (2,i,n), &
00172                         neighbors_3d (3,i,n)) = n_list
00173             endif
00174          endif
00175          end do ! i = 1, nloc
00176       end do ! n = 1, num_neigh
00177 !
00178 #ifdef PRISM_ASSERTION
00179       if (n_list > PRODUCT (size)) then
00180          print *, trim(ch_id), ": n_list =", n_list, nloc, size
00181          call psmile_assert (__FILE__, __LINE__, &
00182               "Invalid number of entries n_list generated")
00183       endif
00184 #endif
00185 !
00186 !===> Special case: all points
00187 !
00188 #ifdef ONLY_FOR_TESTING
00189       print *, '#### psmile_compact_neighbors_3d: ONLY FOR TESTING ACTIVE'
00190       if ( .false. ) then
00191 #else
00192 !
00193 ! Do not jump into send_entire_valid_shape for corner neighbours
00194 ! Currently this is avieved by requiring num_neigh > 1. nneighbour
00195 ! search with only 1 neighbour required would also fall into
00196 ! this category. For the rather rare case we would also run
00197 ! into the less optimal part and establish a list.
00198 !
00199       if (n_list == PRODUCT (size) .and. num_neigh > 1 ) then
00200 #endif
00201          do n = 1, num_neigh
00202 !cdir vector
00203             do i = 1, nloc
00204 
00205             neighbors (i, n) = neighbors (i, n) *                    &
00206            (((neighbors_3d (3,i,n)-grid_valid_shape(1,3))*size(2)  + &
00207              (neighbors_3d (2,i,n)-grid_valid_shape(1,2)))*size(1) + &
00208              (neighbors_3d (1,i,n)-grid_valid_shape(1,1)) + 1)
00209 
00210             end do ! i = 1, nloc
00211          end do ! n = 1, num_neigh
00212 !
00213          Nullify (send_info%list_entries)
00214          send_info%send_entire_valid_shape = .true.
00215 
00216       else
00217 
00218          send_info%send_entire_valid_shape = .false.
00219 !
00220          Allocate (send_info%list_entries (ndim_3d, n_list), &
00221                    STAT = ierror)
00222 !
00223          if (ierror > 0) then
00224             ierrp (1) = ierror
00225             ierrp (2) = ndim_3d * n_list
00226 
00227             ierror = PRISM_Error_Alloc
00228 
00229             call psmile_error ( ierror, 'send_info%list_entries', &
00230                                 ierrp, 2, __FILE__, __LINE__ )
00231             return
00232          endif
00233 !
00234          list_entries => send_info%list_entries
00235 !
00236 !===> Generate compact list of entries to be sent
00237 !
00238          do n = 1, num_neigh
00239 !cdir vector
00240             do i = 1, nloc
00241             if (neighbors (i,n) > 0) then
00242                 neighbors (i,n) = indices (neighbors_3d (1,i,n), &
00243                                            neighbors_3d (2,i,n), &
00244                                            neighbors_3d (3,i,n))
00245             endif
00246             end do ! i = 1, nloc
00247          end do ! n = 1, num_neigh
00248 !
00249          do k = grid_valid_shape(1,3), grid_valid_shape(2,3)
00250             do j = grid_valid_shape(1,2), grid_valid_shape(2,2)
00251                do i = grid_valid_shape(1,1), grid_valid_shape(2,1)
00252                if (indices (i,j,k) > 0) then
00253                   list_entries (1, indices (i,j,k)) = i
00254                   list_entries (2, indices (i,j,k)) = j
00255                   list_entries (3, indices (i,j,k)) = k
00256                endif
00257                end do
00258             end do
00259          end do
00260 
00261       endif
00262 !
00263       Deallocate (indices)
00264 !
00265 #ifdef PRISM_ASSERTION
00266 !
00267 !===> Control list generated
00268 !
00269       do n = 1, num_neigh
00270          do i = 1, nloc
00271          if (neighbors (i, n) < 0 .or. neighbors (i, n) > n_list) exit
00272          end do
00273 !
00274          if (i < nloc) then
00275              print *, trim(ch_id), i, n, neighbors (i, n), n_list
00276              call psmile_assert (__FILE__, __LINE__, &
00277                   "Invalid neighbour list generated ")
00278          endif
00279       end do
00280 #endif
00281 !
00282 !-----------------------------------------------------------------------
00283 !     Add points found in global search.
00284 !     The points are added at the end of the buffer.
00285 !-----------------------------------------------------------------------
00286 !
00287       if ( send_info%nrecv > 0 .and. send_info%num2recv > 0 ) then
00288 !
00289 !
00290 #ifdef PRISM_ASSERTION
00291 !
00292 !===> Control list generated
00293 !
00294             do n = 1, num_neigh
00295                do i = 1, nloc
00296                if (neighbors_3d (1, i, n) == extra_search%global_marker .and. &
00297                   (neighbors_3d (2, i, n) < 1 .or.               &
00298                    neighbors_3d (2, i, n) > send_info%num2recv)) exit
00299                end do
00300 
00301                if (i < nloc) then
00302                    print *, trim(ch_id), ": i, n, neighbors_3d", &
00303                          i, n, neighbors_3d (1:ndim_3d, i, n)
00304                    print *, trim(ch_id), ": nrecv, num2recv", &
00305                          send_info%nrecv, send_info%num2recv
00306 !
00307                    call psmile_assert (__FILE__, __LINE__, &
00308                         "Invalid index of global search ")
00309                endif
00310             end do
00311 #endif
00312 !
00313             do n = 1, num_neigh
00314 !cdir vector
00315                do i = 1, nloc
00316                if (neighbors_3d (1, i, n) == extra_search%global_marker) then
00317                    neighbors (i, n) = n_list + neighbors_3d (2, i, n)
00318                endif
00319                end do ! i = 1, nloc
00320             end do ! n = 1, num_neigh
00321 !
00322          n_list = n_list + send_info%num2recv
00323       endif
00324 !
00325 !===> All done
00326 !
00327       send_info%n_list = n_list
00328 !
00329 #ifdef VERBOSE
00330       print 9980, trim(ch_id), n_list
00331 
00332       call psmile_flushstd
00333 #endif /* VERBOSE */
00334 !
00335 !  Formats:
00336 !
00337 9990 format (1x, a, ': psmile_compact_neighbors_3d')
00338 9980 format (1x, a, ': psmile_compact_neighbors_3d: eof, n_list', i9)
00339 
00340       end subroutine psmile_compact_neighbors_3d

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1