psmile_ext_compact_irreg2_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_ext_compact_irreg2_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_ext_compact_irreg2_real (          &
00012                send_info, array, shape, grid_valid_shape,  &
00013                dest_vector, dest_size, ierror)
00014 !
00015 ! !USES:
00016 !
00017       use PRISM_constants
00018 !
00019       use PSMILe, dummy_interface => PSMILe_Ext_compact_irreg2_real
00020 
00021       implicit none
00022 !
00023 ! !INPUT PARAMETERS:
00024 !
00025       Type (Send_information), Intent (InOut)  :: send_info
00026 !
00027 !     Structure containing the send information
00028 !
00029       Integer, Intent (In)            :: shape (2, ndim_2d)
00030 !
00031 !     Shape of the array "array"
00032 !
00033       Real, Intent (In)               :: array (shape(1,1):shape(2,1), 
00034                                                 shape(1,2):shape(2,2))
00035 !
00036 !     Array from which the compact list should be extracted
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       Integer, Intent (In)            :: dest_size
00043 !
00044 !     Dimension of output vector "dest_vector"
00045 !
00046 ! !OUTPUT PARAMETERS:
00047 !
00048       Real, Intent (Out)              :: dest_vector (dest_size)
00049 !
00050 !     Output vector of compact list
00051 !
00052       Integer, Intent (Out)           :: ierror
00053 
00054 !     Returns the error code of PSMILE_Ext_compact_irreg2_real;
00055 !             ierror = 0 : No error
00056 !             ierror > 0 : Severe error
00057 !
00058 ! !LOCAL VARIABLES
00059 !
00060 !     ... for extracting the compact list
00061 !
00062       Integer                         :: k, n, ni, nij
00063 #ifdef DONT_USE_RESHAPE
00064       Integer                         :: j
00065 #else
00066       Integer                         :: dest_shape (1)
00067 #endif
00068       Integer, Pointer                :: list_entries (:, :)
00069 !
00070 !     ... for error handling
00071 !
00072 !     Integer, parameter              :: nerrp = 2
00073 !     Integer                         :: ierrp (nerrp)
00074 !
00075 ! !DESCRIPTION:
00076 !
00077 ! Subroutine "PSMILe_Ext_compact_irreg2_real" extracts the data
00078 ! required from the 2-dimensional array "array" and stores the data in
00079 ! destination vector "dest_vector".
00080 !
00081 ! !REVISION HISTORY:
00082 !
00083 !   Date      Programmer   Description
00084 ! ----------  ----------   -----------
00085 ! 03.07.21    H. Ritzdorf  created
00086 !
00087 !EOP
00088 !----------------------------------------------------------------------
00089 !
00090 !  $Id: psmile_ext_compact_irreg2_real.F90 2325 2010-04-21 15:00:07Z valcke $
00091 !  $Author: valcke $
00092 !
00093    Character(len=len_cvs_string), save :: mycvs = 
00094        '$Id: psmile_ext_compact_irreg2_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00095 !
00096 !----------------------------------------------------------------------
00097 !----------------------------------------------------------------------
00098 !
00099 !  Initialization
00100 !
00101 #ifdef VERBOSE
00102       print 9990, trim(ch_id), send_info%send_entire_valid_shape
00103 
00104       call psmile_flushstd
00105 #endif /* VERBOSE */
00106 !
00107 #ifdef PRISM_ASSERTION
00108       if (dest_size < send_info%n_list) then
00109          print *, trim(ch_id), 'dest_size, n_list', dest_size, send_info%n_list
00110          call psmile_assert (__FILE__, __LINE__, &
00111                              "dest_size is not sufficient")
00112       endif
00113 #endif
00114 !
00115       ierror = 0
00116 !
00117 !===> Special case: All points of grid shape
00118 !
00119       if (send_info%send_entire_valid_shape) then
00120 !
00121          ni  =       grid_valid_shape (2, 1) - grid_valid_shape (1, 1) + 1
00122          nij = ni * (grid_valid_shape (2, 2) - grid_valid_shape (1, 2) + 1)
00123 !
00124 #ifdef PRISM_ASSERTION
00125          if (dest_size < nij*(grid_valid_shape(2,3)-grid_valid_shape(1,3)+1)) then
00126             print *, trim(ch_id), 'dest_size, nijk', dest_size, &
00127                      nij*(grid_valid_shape(2,3)-grid_valid_shape(1,3)+1)
00128      
00129             call psmile_assert (__FILE__, __LINE__, &
00130                                 "dest_size is not sufficient")
00131          endif
00132 #endif
00133 
00134 #ifdef DONT_USE_RESHAPE
00135          n = 0
00136 !
00137             do j = grid_valid_shape (1, 2), grid_valid_shape (2, 2)
00138             dest_vector (n+1:n+ni) = &
00139               array (grid_valid_shape (1,1):grid_valid_shape (2,1), j)
00140             n = n + ni
00141             end do
00142 #else
00143          dest_shape (1) = nij
00144          dest_vector (1:nij) = RESHAPE ( &
00145                        array(grid_valid_shape(1,1):grid_valid_shape(2,1),  &
00146                              grid_valid_shape(1,2):grid_valid_shape(2,2)), &
00147                        dest_shape)
00148 #endif
00149 !
00150 !===> ... Copy values for 3rd index
00151 !         ! Das kann man auch noch besser machen !
00152 !
00153          do k = 2, grid_valid_shape (2, 3) - grid_valid_shape (1, 3) + 1
00154             dest_vector((k-1)*nij+1:k*nij) = dest_vector(1:nij)
00155          end do
00156 
00157       else
00158 !
00159 !===> Get compact list
00160 !
00161          list_entries => send_info%list_entries
00162 
00163 #ifdef PRISM_ASSERTION
00164 !cdir vector
00165             do n = 1, send_info%n_list - send_info%num2recv
00166             if (list_entries(1, n) < shape (1,1) .or. &
00167                 list_entries(1, n) > shape (2,1) .or. &
00168                 list_entries(2, n) < shape (1,2) .or. &
00169                 list_entries(2, n) > shape (2,2)) exit
00170             end do
00171 
00172          if (n < send_info%n_list - send_info%num2recv) then
00173             print *, 'n, list_entry, shape', n, list_entries(1, n), &
00174                      list_entries(2, n), shape
00175             call psmile_assert (__FILE__, __LINE__, &
00176                                 "incorrect index in list_entries")
00177          endif
00178 #endif
00179 
00180 !cdir vector
00181             do n = 1, send_info%n_list - send_info%num2recv
00182             dest_vector (n) = array (list_entries(1, n), &
00183                                      list_entries(2, n))
00184             end do
00185 
00186       endif
00187 !
00188 !===> All done
00189 !
00190 #ifdef VERBOSE
00191       print 9980, trim(ch_id)
00192 
00193       call psmile_flushstd
00194 #endif /* VERBOSE */
00195 !
00196 !  Formats:
00197 !
00198 9990 format (1x, a, ': psmile_ext_compact_irreg2_real: send_entire_valid', l2)
00199 9980 format (1x, a, ': psmile_ext_compact_irreg2_real: eof')
00200 
00201       end subroutine psmile_Ext_compact_irreg2_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1