00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_ext_compact_irreg2_dble ( &
00012 send_info, array, shape, grid_valid_shape, &
00013 dest_vector, dest_size, ierror)
00014
00015
00016
00017 use PRISM_constants
00018
00019 use PSMILe, dummy_interface => PSMILe_Ext_compact_irreg2_dble
00020
00021 implicit none
00022
00023
00024
00025 Type (Send_information), Intent (InOut) :: send_info
00026
00027
00028
00029 Integer, Intent (In) :: shape (2, ndim_2d)
00030
00031
00032
00033 Double Precision, Intent (In) :: array (shape(1,1):shape(2,1),
00034 shape(1,2):shape(2,2))
00035
00036
00037
00038 Integer, Intent (In) :: grid_valid_shape (2, ndim_3d)
00039
00040
00041
00042 Integer, Intent (In) :: dest_size
00043
00044
00045
00046
00047
00048 Double Precision, Intent (Out) :: dest_vector (dest_size)
00049
00050
00051
00052 Integer, Intent (Out) :: ierror
00053
00054
00055
00056
00057
00058
00059
00060
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
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093 Character(len=len_cvs_string), save :: mycvs =
00094 '$Id: psmile_ext_compact_irreg2_dble.F90 2325 2010-04-21 15:00:07Z valcke $'
00095
00096
00097
00098
00099
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
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
00151
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
00160
00161 list_entries => send_info%list_entries
00162
00163 #ifdef PRISM_ASSERTION
00164
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
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
00189
00190 #ifdef VERBOSE
00191 print 9980, trim(ch_id)
00192
00193 call psmile_flushstd
00194 #endif /* VERBOSE */
00195
00196
00197
00198 9990 format (1x, a, ': psmile_ext_compact_irreg2_dble: send_entire_valid', l2)
00199 9980 format (1x, a, ': psmile_ext_compact_irreg2_dble: eof')
00200
00201 end subroutine psmile_Ext_compact_irreg2_dble