psmile_ccompact_3d_reg_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_ccompact_3d_reg_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_ccompact_3d_reg_real ( 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_3d_reg_real
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       Real, Intent (In)               :: array_x ( shape(1,1):shape(2,1), 
00045                                                    nb_corners )
00046       Real, Intent (In)               :: array_y ( shape(1,2):shape(2,2), 
00047                                                    nb_corners )
00048       Real, Intent (In)               :: array_z ( shape(1,3):shape(2,3) )
00049 !
00050 !     Arrays from which the compact list should be extracted
00051 !
00052       Type (Extra_search_info), Intent (InOut) :: extra_search
00053 !
00054 !     Info's on global coordinates to be searched
00055 !
00056       Integer, Intent (In)            :: dest_size
00057 !
00058 !     Size of destination arrays
00059 !
00060       Integer, Intent (In)            :: nbr_cells_tot
00061 !
00062 !     Total number of source cells
00063 !
00064       Integer, Intent (InOut)         :: source_cell_index (nbr_cells_tot)
00065 !
00066 !     Index of source cells
00067 !
00068       Integer, Intent (InOut)         :: neighcells (nbr_cells_tot, 4)
00069 !
00070 !     Dimension of output vector "dest_vector"
00071 !
00072 ! !OUTPUT PARAMETERS:
00073 !
00074       Real, Intent (Out)              :: dest_x (dest_size)
00075       Real, Intent (Out)              :: dest_y (dest_size)
00076       Real, Intent (Out)              :: dest_z (dest_size)
00077 !
00078 !     Output vectors of compact list
00079 !
00080       Integer, Intent (Out)           :: ierror
00081 
00082 !     Returns the error code of PSMILE_ccompact_3d_reg_real;
00083 !             ierror = 0 : No error
00084 !             ierror > 0 : Severe error
00085 !
00086 ! !LOCAL VARIABLES
00087 !
00088 !     ... for compact list
00089 !
00090       Integer, Pointer                :: list_entries (:, :)
00091       Integer                         :: nbr_cells_loc
00092 !
00093 !     ... for extracting the compact list
00094 
00095       Integer                         :: l, n, ni, nk
00096       Integer                         :: nij, nij_loc
00097 !
00098 !     ... for globally searched data
00099 !
00100 
00101       Real, Parameter                 :: tol = 1.0d-12 
00102 
00103       Integer                         :: i, j, k
00104       Integer                         :: dest_idx_off
00105       Integer                         :: idx1, idx2, jdx1, jdx2, stride
00106 !
00107 ! !DESCRIPTION:
00108 !
00109 ! Subroutine "PSMILe_ccompact_3d_reg_real" extracts the data
00110 ! required from the 2-dimensional corner array "array" and stores the data in
00111 ! destination vector "dest_vector".
00112 !
00113 ! Work on array_x and array_y in 1 role
00114 ! detect sense
00115 ! determine index
00116 ! we have to build two lists, one for the corners of the cells and another one for the
00117 ! cells (data) itself!!! 
00118 !
00119 !
00120 ! !REVISION HISTORY:
00121 !
00122 !   Date      Programmer   Description
00123 ! ----------  ----------   -----------
00124 ! 03.07.21    H. Ritzdorf  created
00125 !
00126 !EOP
00127 !----------------------------------------------------------------------
00128 !
00129 !  $Id: psmile_ccompact_3d_reg_real.F90 2590 2010-09-23 15:56:33Z hanke $
00130 !  $Author: hanke $
00131 !
00132    Character(len=len_cvs_string), save :: mycvs = 
00133        '$Id: psmile_ccompact_3d_reg_real.F90 2590 2010-09-23 15:56:33Z hanke $'
00134 !
00135 !----------------------------------------------------------------------
00136 !
00137 !  Initialization
00138 !
00139 #ifdef VERBOSE
00140       print 9990, trim(ch_id), send_info%send_entire_valid_shape
00141 
00142       call psmile_flushstd
00143 #endif /* VERBOSE */
00144 !
00145 #ifdef PRISM_ASSERTION
00146       if (dest_size < send_info%n_list) then
00147          print *, trim(ch_id), 'dest_size, n_list', dest_size, send_info%n_list
00148          call psmile_assert (__FILE__, __LINE__, &
00149                              "dest_size is not sufficient")
00150       endif
00151 #endif
00152 !
00153       ierror = 0
00154 
00155       nbr_cells_loc = nbr_cells_tot - send_info%num2recv
00156 
00157 !------------------------------------------------------------------------
00158 !     Store corners
00159 !------------------------------------------------------------------------
00160 !
00161 !===> Get compact list
00162 !
00163       ! Initialization of output variables
00164       dest_x = -1
00165       dest_y = -1
00166       dest_z = -1
00167 
00168       if ( send_info%send_entire_valid_shape ) then
00169 
00170          ni  =      (grid_valid_shape (2,1) - grid_valid_shape (1,1) + 1)
00171          nij = ni * (grid_valid_shape (2,2) - grid_valid_shape (1,2) + 1) * 4
00172          nk =       (grid_valid_shape (2,3) - grid_valid_shape (1,3) + 1)
00173 
00174 #ifdef PRISM_ASSERTION
00175          if (dest_size < nij * nk ) then
00176             print *, trim(ch_id), 'dest_size, nijk', dest_size, nij*nk
00177             call psmile_assert ( __FILE__, __LINE__, &
00178                                 "dest_size is not sufficient")
00179          endif
00180 #endif
00181          n = 0
00182 
00183          do k = grid_valid_shape(1,3), grid_valid_shape(2,3)
00184             do j = grid_valid_shape(1,2), grid_valid_shape(2,2)
00185                do i = grid_valid_shape(1,1), grid_valid_shape(2,1)
00186                   n = n + 1
00187                   dest_x(n) = min(array_x(i,1), array_x(i,2))
00188                   dest_y(n) = min(array_y(j,1), array_y(j,2))
00189                   dest_z(n) = array_z(k)
00190                end do
00191             end do
00192          end do
00193 
00194          do k = grid_valid_shape(1,3), grid_valid_shape(2,3)
00195             do j = grid_valid_shape(1,2), grid_valid_shape(2,2)
00196                do i = grid_valid_shape(1,1), grid_valid_shape(2,1)
00197                   n = n + 1
00198                   dest_x(n) = max(array_x(i,1), array_x(i,2))
00199                   dest_y(n) = min(array_y(j,1), array_y(j,2))
00200                   dest_z(n) = array_z(k)
00201                end do
00202             end do
00203          end do
00204 
00205          do k = grid_valid_shape(1,3), grid_valid_shape(2,3)
00206             do j = grid_valid_shape(1,2), grid_valid_shape(2,2)
00207                do i = grid_valid_shape(1,1), grid_valid_shape(2,1)
00208                   n = n + 1
00209                   dest_x(n) = max(array_x(i,1), array_x(i,2))
00210                   dest_y(n) = max(array_y(j,1), array_y(j,2))
00211                   dest_z(n) = array_z(k)
00212                end do
00213             end do
00214          end do
00215 
00216          do k = grid_valid_shape(1,3), grid_valid_shape(2,3)
00217             do j = grid_valid_shape(1,2), grid_valid_shape(2,2)
00218                do i = grid_valid_shape(1,1), grid_valid_shape(2,1)
00219                   n = n + 1
00220                   dest_x(n) = min(array_x(i,1), array_x(i,2))
00221                   dest_y(n) = max(array_y(j,1), array_y(j,2))
00222                   dest_z(n) = array_z(k)
00223                end do
00224             end do
00225          end do
00226 
00227       else
00228 !
00229 !===> Get compact list
00230 !
00231          list_entries => send_info%list_entries
00232 
00233          nij     = send_info%n_list
00234          nij_loc = nij - send_info%num2recv
00235 
00236 #ifdef PRISM_ASSERTION
00237          if (.not. Associated (list_entries)) then
00238             call psmile_assert (__FILE__, __LINE__, &
00239                  "list_entries are not available")
00240          endif
00241 
00242 !cdir vector
00243          do n = 1, nij_loc
00244             if ( list_entries(1, n) < grid_valid_shape (1,1) .or. &
00245                  list_entries(1, n) > grid_valid_shape (2,1) .or. &
00246                  list_entries(2, n) < grid_valid_shape (1,2) .or. &
00247                  list_entries(2, n) > grid_valid_shape (2,2)) exit
00248          end do
00249 
00250          if ( n < nij_loc ) then
00251             print *, 'n, list_entry, grid_valid_shape', n, list_entries(1, n), &
00252                  list_entries(2, n), grid_valid_shape
00253             call psmile_assert (__FILE__, __LINE__, &
00254                  "incorrect index in list_entries")
00255          endif
00256 #endif
00257          do n = 1, nij_loc
00258             dest_x(n) = min(array_x(list_entries(1,n),1), &
00259                             array_x(list_entries(1,n),2))
00260             dest_y(n) = min(array_y(list_entries(2,n),1), &
00261                             array_y(list_entries(2,n),2))
00262             dest_z(n) = array_z(list_entries(3,n))
00263          end do
00264 
00265          do n = 1, nij_loc
00266             dest_x(nij+n) = max(array_x(list_entries(1,n),1), &
00267                                 array_x(list_entries(1,n),2))
00268             dest_y(nij+n) = min(array_y(list_entries(2,n),1), &
00269                                 array_y(list_entries(2,n),2))
00270             dest_z(nij+n) = array_z(list_entries(3,n))
00271          end do
00272 
00273          do n = 1, nij_loc
00274             dest_x(2*nij+n) = max(array_x(list_entries(1,n),1), &
00275                                   array_x(list_entries(1,n),2))
00276             dest_y(2*nij+n) = max(array_y(list_entries(2,n),1), &
00277                                   array_y(list_entries(2,n),2))
00278             dest_z(2*nij+n) = array_z(list_entries(3,n))
00279 
00280          end do
00281 
00282          do n = 1, nij_loc
00283             dest_x(3*nij+n) = min(array_x(list_entries(1,n),1), &
00284                                   array_x(list_entries(1,n),2))
00285             dest_y(3*nij+n) = max(array_y(list_entries(2,n),1), &
00286                                   array_y(list_entries(2,n),2))
00287             dest_z(3*nij+n) = array_z(list_entries(3,n))
00288          end do
00289 
00290       endif
00291 
00292 !-----------------------------------------------------------------------
00293 !     Create source cell index list
00294 !-----------------------------------------------------------------------
00295 
00296       do n = 2, 4
00297          do i = 1, nbr_cells_loc
00298             neighcells (i,n) = neighcells (i,1) + (n-1) * nij
00299          enddo
00300       enddo
00301 
00302       source_cell_index(:) =  neighcells(:,1)
00303 !
00304 !-----------------------------------------------------------------------
00305 !     Add cells found in global search.
00306 !     The cells are added at the end of the buffer.
00307 !-----------------------------------------------------------------------
00308 !
00309       if ( send_info%num2recv > 0 ) then
00310 
00311          do n = 2, 4
00312             do i = 1, send_info%num2recv
00313                neighcells (nbr_cells_loc+i,n) = &
00314                neighcells (nbr_cells_loc+i,1) + (n-1) * nij
00315             end do
00316          enddo
00317 
00318          ! calucalate offset for dest_*
00319          ! we write behind the data of the local process
00320          dest_idx_off = nij_loc
00321 
00322          do l = 1, send_info%nrecv
00323 
00324             stride = nb_corners*send_info%len_sent(l)
00325 
00326             do i = 1, send_info%len_sent(l)
00327                idx1 = i
00328                idx2 = send_info%len_sent(l)+i
00329 
00330                jdx1 = stride + idx1
00331                jdx2 = stride + idx2
00332 
00333                n = 0
00334                dest_x (n*nij+dest_idx_off+i) = extra_search%real_bufs(l)%vector(idx1)
00335                dest_y (n*nij+dest_idx_off+i) = extra_search%real_bufs(l)%vector(jdx1)
00336                ! Rene : vertical coordinate not clean yet !!
00337                dest_z (n*nij+dest_idx_off+i) = array_z(1)
00338 
00339                n = 1
00340                dest_x (n*nij+dest_idx_off+i) = extra_search%real_bufs(l)%vector(idx2)
00341                dest_y (n*nij+dest_idx_off+i) = extra_search%real_bufs(l)%vector(jdx1)
00342                ! Rene : vertical coordinate not clean yet !!
00343                dest_z (n*nij+dest_idx_off+i) = array_z(1)
00344 
00345                n = 2
00346                dest_x (n*nij+dest_idx_off+i) = extra_search%real_bufs(l)%vector(idx2)
00347                dest_y (n*nij+dest_idx_off+i) = extra_search%real_bufs(l)%vector(jdx2)
00348                ! Rene : vertical coordinate not clean yet !!
00349                dest_z (n*nij+dest_idx_off+i) = array_z(1)
00350 
00351                n = 3
00352                dest_x (n*nij+dest_idx_off+i) = extra_search%real_bufs(l)%vector(idx1)
00353                dest_y (n*nij+dest_idx_off+i) = extra_search%real_bufs(l)%vector(jdx2)
00354                ! Rene : vertical coordinate not clean yet !!
00355                dest_z (n*nij+dest_idx_off+i) = array_z(1)
00356 
00357             end do
00358 
00359             ! we do not want to overwrite the data just written
00360             dest_idx_off = dest_idx_off + send_info%len_sent(l)
00361 
00362             Deallocate (extra_search%real_bufs(l)%vector)
00363 
00364          enddo
00365 
00366       endif
00367 !
00368 !===> All done
00369 !
00370 #ifdef VERBOSE
00371       print 9980, trim(ch_id)
00372 
00373       call psmile_flushstd
00374 #endif /* VERBOSE */
00375 !
00376 !  Formats:
00377 !
00378 9990 format (1x, a, ': psmile_ccompact_3d_reg_real: send_entire_valid', l2)
00379 9980 format (1x, a, ': psmile_ccompact_3d_reg_real: eof')
00380 
00381       end subroutine psmile_ccompact_3d_reg_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1