psmile_ccompact_irreg2_dble.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_ccompact_irreg2_dble
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_ccompact_irreg2_dble ( send_info,   &
00012                grid_valid_shape, shape, nb_corners,         &
00013                array_x, array_y, array_z,                   &
00014                extra_search, dest_size, nbr_cells_tot,      &
00015                source_cell_index,                           &
00016                neighcells, dest_x, dest_y, dest_z, ierror )
00017 !
00018 ! !USES:
00019 !
00020       use PRISM_constants
00021 !
00022       use PSMILe, dummy_interface => PSMILe_ccompact_irreg2_dble
00023 
00024       implicit none
00025 !
00026 ! !INPUT PARAMETERS:
00027 !
00028       Type (Send_information), Intent (InOut)  :: send_info
00029 !
00030 !     Structure containing the send information
00031 !
00032       Integer, Intent (In)            :: grid_valid_shape (2, ndim_3d)
00033 
00034 !     Specifies the valid block shape for the "ndim_3d"-dimensional block
00035 !
00036       Integer, Intent (In)            :: shape (2, ndim_3d)
00037 !
00038 !     Shape of the array "array"
00039 !
00040       Integer, Intent (In)            :: nb_corners
00041 !
00042 !     Number of corners in "array_x/y"
00043 !
00044       Double Precision, Intent (In)   :: array_x ( shape(1,1):shape(2,1), 
00045                                                    shape(1,2):shape(2,2), 
00046                                                    nb_corners )
00047       Double Precision, Intent (In)   :: array_y ( shape(1,1):shape(2,1), 
00048                                                    shape(1,2):shape(2,2), 
00049                                                    nb_corners )
00050       Double Precision, Intent (In)   :: array_z ( shape(1,3):shape(2,3) )
00051 !
00052 !     Arrays from which the compact list should be extracted
00053 !
00054       Type (Extra_search_info), Intent (InOut) :: extra_search
00055 !
00056 !     Info's on global coordinates to be searched
00057 !
00058       Integer, Intent (In)            :: dest_size
00059 !
00060 !     Size of destination arrays
00061 !
00062       Integer, Intent (In)            :: nbr_cells_tot
00063 !
00064 !     Total number of source cells
00065 !
00066       Integer, Intent (InOut)         :: source_cell_index (nbr_cells_tot)
00067 !
00068 !     Index of source cells
00069 !
00070       Integer, Intent (InOut)         :: neighcells (nbr_cells_tot, nb_corners)
00071 !
00072 !
00073 ! !OUTPUT PARAMETERS:
00074 !
00075       Double Precision, Intent (Out)  :: dest_x (dest_size)
00076       Double Precision, Intent (Out)  :: dest_y (dest_size)
00077       Double Precision, Intent (Out)  :: dest_z (dest_size)
00078 !
00079 !     Output vectors of compact list
00080 !
00081       Integer, Intent (Out)           :: ierror
00082 
00083 !     Returns the error code of PSMILE_ccompact_irreg2_dble;
00084 !             ierror = 0 : No error
00085 !             ierror > 0 : Severe error
00086 !
00087 ! !LOCAL VARIABLES
00088 !
00089 !     ... for compact list
00090 !
00091       Integer, Pointer                :: list_entries (:, :)
00092       Integer                         :: nbr_cells_loc
00093 !
00094 !     ... for extracting the compact list
00095 !
00096       Integer                         :: i, k, l, n
00097       Integer                         :: ni, nij, nk, nij_loc
00098 #if !defined DONT_USE_RESHAPE
00099       Integer                         :: dest_shape (1)
00100 #endif
00101 !
00102 !     ... for globally searched data
00103 !
00104       Integer                         :: dest_idx_off
00105       Integer                         :: stride, idx1, idx2
00106 !
00107 ! !DESCRIPTION:
00108 !
00109 ! Subroutine "PSMILe_ccompact_irreg2_dble" extracts the data
00110 ! required from the 2-dimensional corner array "array" and stores the data in
00111 ! destination vector "dest_vector".
00112 !
00113 ! !REVISION HISTORY:
00114 !
00115 !   Date      Programmer   Description
00116 ! ----------  ----------   -----------
00117 ! 03.07.21    R. Redler    created
00118 !
00119 !EOP
00120 !----------------------------------------------------------------------
00121 !
00122 !  $Id: psmile_ccompact_irreg2_dble.F90 2590 2010-09-23 15:56:33Z hanke $
00123 !  $Author: hanke $
00124 !
00125    Character(len=len_cvs_string), save :: mycvs = 
00126        '$Id: psmile_ccompact_irreg2_dble.F90 2590 2010-09-23 15:56:33Z hanke $'
00127 !
00128 !----------------------------------------------------------------------
00129 !
00130 !  Initialization
00131 !
00132 #ifdef VERBOSE
00133       print 9990, trim(ch_id), send_info%send_entire_valid_shape
00134 
00135       call psmile_flushstd
00136 #endif /* VERBOSE */
00137 !
00138 #ifdef PRISM_ASSERTION
00139       if (dest_size < send_info%n_list) then
00140          print *, trim(ch_id), 'dest_size, n_list', dest_size
00141          call psmile_assert (__FILE__, __LINE__, &
00142                              "dest_size is not sufficient")
00143       endif
00144 #endif
00145 !
00146       ierror = 0
00147 
00148       nbr_cells_loc = nbr_cells_tot - send_info%num2recv
00149 
00150       ! Initialization of output variables
00151       dest_x = -1
00152       dest_y = -1
00153       dest_z = -1
00154 !
00155 !------------------------------------------------------------------------
00156 !     Store corners
00157 !------------------------------------------------------------------------
00158 !
00159 !===> Special case: All points of grid shape
00160 !
00161       if (send_info%send_entire_valid_shape) then
00162 !
00163          ni  =      (grid_valid_shape (2, 1) - grid_valid_shape (1, 1) + 1)
00164          nij = ni * (grid_valid_shape (2, 2) - grid_valid_shape (1, 2) + 1) * nb_corners
00165          nk  =      (grid_valid_shape (2, 3) - grid_valid_shape (1, 3) + 1)
00166 !
00167 #ifdef PRISM_ASSERTION
00168          if (dest_size < nij * nk ) then
00169             print *, trim(ch_id), 'dest_size, nijk', dest_size, nij*nk
00170      
00171             call psmile_assert (__FILE__, __LINE__, &
00172                                 "dest_size is not sufficient")
00173          endif
00174 #endif
00175 
00176 #ifdef DONT_USE_RESHAPE
00177          n = 0
00178 !
00179          do l = 1, nb_corners
00180             do j = grid_valid_shape (1, 2), grid_valid_shape (2, 2)
00181                dest_x (n+1:n+ni) = &
00182                     array_x (grid_valid_shape (1,1):grid_valid_shape (2,1), j, l)
00183                dest_y (n+1:n+ni) = &
00184                     array_y (grid_valid_shape (1,1):grid_valid_shape (2,1), j, l)
00185                n = n + ni
00186             end do
00187          enddo
00188 #else
00189          dest_shape (1) = nij
00190 
00191          dest_x (1:nij) = RESHAPE ( &
00192                            array_x(grid_valid_shape(1,1):grid_valid_shape(2,1),  &
00193                                    grid_valid_shape(1,2):grid_valid_shape(2,2),  &
00194                                    nb_corners), dest_shape)
00195          dest_y (1:nij) = RESHAPE ( &
00196                            array_y(grid_valid_shape(1,1):grid_valid_shape(2,1),  &
00197                                    grid_valid_shape(1,2):grid_valid_shape(2,2),  &
00198                                    nb_corners), dest_shape)
00199 #endif
00200          dest_z (1:nij) = array_z(grid_valid_shape(1,3))
00201 !
00202 !===> ... Copy values for 3rd index
00203 !         ! Das kann man auch noch besser machen !
00204 !
00205          do k = 2, grid_valid_shape (2, 3) - grid_valid_shape (1, 3) + 1
00206             dest_x((k-1)*nij+1:k*nij) = dest_x(1:nij)
00207             dest_y((k-1)*nij+1:k*nij) = dest_y(1:nij)
00208             dest_z((k-1)*nij+1:k*nij) = array_z(k)
00209          end do
00210 
00211       else
00212 !
00213 !===> Get compact list
00214 !
00215          list_entries => send_info%list_entries
00216 
00217          nij     = send_info%n_list
00218          nij_loc = nij - send_info%num2recv
00219 
00220 #ifdef PRISM_ASSERTION
00221          if (.not. Associated (list_entries)) then
00222             call psmile_assert (__FILE__, __LINE__, &
00223                  "list_entries are not available")
00224          endif
00225 
00226 !cdir vector
00227          do n = 1, nij_loc
00228             if ( list_entries(1, n) < shape (1,1) .or. &
00229                  list_entries(1, n) > shape (2,1) .or. &
00230                  list_entries(2, n) < shape (1,2) .or. &
00231                  list_entries(2, n) > shape (2,2)) exit
00232          end do
00233 
00234          if ( n < nij_loc ) then
00235             print *, 'n, list_entry, shape', n, list_entries(1, n), &
00236                  list_entries(2, n), shape
00237             call psmile_assert (__FILE__, __LINE__, &
00238                  "incorrect index in list_entries")
00239          endif
00240 #endif
00241 
00242 !cdir vector
00243          do l = 1, nb_corners
00244             do n = 1, nij_loc
00245                dest_x ((l-1)*nij+n) = array_x(list_entries(1, n), &
00246                                               list_entries(2, n),l)
00247                dest_y ((l-1)*nij+n) = array_y(list_entries(1, n), &
00248                                               list_entries(2, n),l)
00249                dest_z ((l-1)*nij+n) = array_z(list_entries(3, n))
00250             end do
00251          enddo
00252 
00253       endif
00254 
00255 !------------------------------------------------------------------------
00256 !     Create source cell index list
00257 !------------------------------------------------------------------------
00258 
00259       do n = 2, nb_corners
00260          do i = 1, nbr_cells_loc
00261             neighcells (i,n) = neighcells (i,1) + (n-1) * nij
00262          end do
00263       enddo
00264 
00265       source_cell_index(:) = neighcells(:,1)
00266 !
00267 !-----------------------------------------------------------------------
00268 !     Add cells found in global search.
00269 !     The cells are added at the end of the buffer.
00270 !-----------------------------------------------------------------------
00271 !
00272       if ( send_info%num2recv > 0 ) then
00273 
00274          do n = 2, nb_corners
00275             do i = 1, send_info%num2recv
00276                neighcells (nbr_cells_loc+i,n) = &
00277                neighcells (nbr_cells_loc+i,1) + (n-1) * nij
00278             end do
00279          enddo
00280 
00281          ! calucalate offset for dest_*
00282          ! we write behind the data of the local process
00283          dest_idx_off = nij_loc
00284 
00285          do l = 1, extra_search%nrecv
00286 
00287             stride = nb_corners*send_info%len_sent(l)
00288 
00289             do k = 1, nb_corners
00290                do n = 1, send_info%len_sent(l)
00291                   idx1 = (k-1)*send_info%len_sent(l)+n
00292                   idx2 = stride + idx1
00293 
00294                   dest_x ((k-1)*nij+dest_idx_off+n) = extra_search%dble_bufs(l)%vector(idx1)
00295                   dest_y ((k-1)*nij+dest_idx_off+n) = extra_search%dble_bufs(l)%vector(idx2)
00296 ! Rene : vertical coordinate not clean yet !!
00297                   dest_z ((k-1)*nij+dest_idx_off+n) = array_z(1)
00298 
00299                end do
00300             enddo
00301 
00302             ! we do not want to overwrite the data just written
00303             dest_idx_off = dest_idx_off + send_info%len_sent(l)
00304 
00305             Deallocate (extra_search%dble_bufs(l)%vector)
00306 
00307          enddo
00308       endif
00309 !
00310 !===> All done
00311 !
00312 #ifdef VERBOSE
00313       print 9980, trim(ch_id), ierror
00314 
00315       call psmile_flushstd
00316 #endif /* VERBOSE */
00317 !
00318 !  Formats:
00319 !
00320 9990 format (1x, a, ': psmile_ccompact_irreg2_dble: send_entire_valid', l2)
00321 9980 format (1x, a, ': psmile_ccompact_irreg2_dble: eof ierror = ', i8)
00322 
00323       end subroutine psmile_ccompact_irreg2_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1