00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 subroutine psmile_put_irr_field_dble (data_array, shape, nbr_fields, &
00013 srcloc, npoints, srcars, nars, nloc, &
00014 dest, tag, comm, ierror)
00015
00016
00017
00018 use PRISM_constants
00019
00020 use PSMILe, dummy_interface => PSMILe_Put_irr_field_dble
00021
00022 implicit none
00023
00024
00025
00026 Integer, Intent (In) :: shape (2, ndim_3d)
00027
00028
00029
00030 Integer, Intent (In) :: nbr_fields
00031
00032
00033
00034 Double Precision, Intent (In) :: data_array (shape(1,1):shape(2,1),
00035 shape(1,2):shape(2,2),
00036 shape(1,3):shape(2,3),
00037 nbr_fields)
00038
00039 Integer, Intent (In) :: nloc
00040
00041
00042
00043
00044 Integer, Intent (In) :: npoints
00045
00046
00047
00048 Integer, Intent (In) :: srcloc (ndim_3d, npoints)
00049
00050
00051
00052 Integer, Intent (In) :: nars
00053
00054
00055
00056 Integer, Intent (In) :: srcars (2, ndim_3d, nars)
00057
00058
00059
00060 Integer, Intent (In) :: dest
00061
00062
00063
00064 Integer, Intent (In) :: tag
00065
00066
00067
00068 Integer, Intent (In) :: comm
00069
00070
00071
00072
00073
00074 Integer, Intent (Out) :: ierror
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084 Double Precision, Pointer :: buffer (:)
00085
00086
00087
00088 Integer :: i, n, nar, np, nprev
00089 #define USE_PACK
00090 #ifndef USE_PACK
00091 Integer :: j, k
00092 #endif
00093
00094
00095
00096 Integer, parameter :: nerrp = 3
00097 Integer :: ierrp (nerrp)
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121 Character(len=len_cvs_string), save :: mycvs =
00122 '$Id: psmile_put_irr_field_dble.F90 2325 2010-04-21 15:00:07Z valcke $'
00123
00124
00125
00126
00127
00128 #ifdef VERBOSE
00129 print 9990, trim(ch_id), dest
00130
00131 call psmile_flushstd
00132 #endif /* VERBOSE */
00133
00134 ierror = 0
00135
00136 #ifdef PRISM_ASSERTION
00137
00138
00139
00140 nprev = npoints
00141 do nar = 1, nars
00142 nprev = nprev + (srcars(2,1,nar)-srcars(1,1,nar)+1) * &
00143 (srcars(2,2,nar)-srcars(1,2,nar)+1) * &
00144 (srcars(2,3,nar)-srcars(1,3,nar)+1)
00145 enddo
00146
00147 if (nloc /= nprev) then
00148 print *, 'nloc, nprev, npoints, nars', nloc, nprev, npoints, nars
00149 call psmile_assert ( __FILE__, __LINE__, &
00150 "nloc /= npoints + SUM(points in ars)")
00151 endif
00152 #endif
00153
00154
00155
00156
00157 Allocate (buffer(nloc*nbr_fields), STAT = ierror)
00158 if ( ierror > 0 ) then
00159 ierrp (1) = ierror
00160 ierrp (2) = nloc * nbr_fields
00161
00162 ierror = PRISM_Error_Alloc
00163
00164 call psmile_error ( ierror, 'buffer', &
00165 ierrp, 2, __FILE__, __LINE__ )
00166 return
00167 endif
00168
00169
00170
00171
00172 #ifdef DEBUG
00173 print 9960, nloc, nbr_fields
00174 9960 format (1x, 'psmile_put_irr_field_dble: unpack: nloc, nbr_fields', i7, i4)
00175 print 9961, shape
00176 9961 format (1x, 'psmile_put_irr_field_dble: shape', 6i5)
00177 #endif
00178
00179 do n = 1, nbr_fields
00180
00181 do i = 1, npoints
00182 buffer ((n-1)*nloc + i) = data_array (srcloc (1, i), &
00183 srcloc (2, i), &
00184 srcloc (3, i), n)
00185 end do
00186
00187 nprev = (n-1)*nloc + npoints
00188 do nar = 1, nars
00189 #ifdef USE_PACK
00190 np = (srcars(2,1,nar)-srcars(1,1,nar)+1) * &
00191 (srcars(2,2,nar)-srcars(1,2,nar)+1) * &
00192 (srcars(2,3,nar)-srcars(1,3,nar)+1)
00193
00194 buffer (nprev+1:nprev+np) = &
00195 pack (data_array(srcars(1,1,nar):srcars(2,1,nar), &
00196 srcars(1,2,nar):srcars(2,2,nar), &
00197 srcars(1,3,nar):srcars(2,3,nar), n), .true.)
00198 nprev = nprev + np
00199 #else
00200 np = srcars(2,1,nar)-srcars(1,1,nar)+1
00201
00202 do k = srcars(1,3,nar), srcars(2,3,nar)
00203 do j = srcars(1,2,nar), srcars(2,2,nar)
00204 buffer (nprev+1:nprev+np) = &
00205 data_array (srcars(1,1,nar):srcars(2,1,nar), j, k, n)
00206 nprev = nprev + np
00207 end do
00208 end do
00209 #endif /* HUHU */
00210 end do
00211 end do
00212
00213 #ifdef PRISM_ASSERTION
00214 if (nloc*nbr_fields /= nprev) then
00215 print *, 'nloc*nbr_fields, nprev', nloc*nbr_fields, nprev
00216 call psmile_assert ( __FILE__, __LINE__, &
00217 "Error in packing: nprev /= nloc*nbr_fields")
00218 endif
00219 #endif
00220
00221 #ifdef DEBUG
00222 if (npoints > 0) then
00223 print 9950, srcloc (:,1), &
00224 data_array(srcloc (1,1), srcloc (2,1), srcloc (3,1), 1)
00225 9950 format (1x, 'psmile_put_irr_field_dble: srcijk', 3i5, ' val', f10.5)
00226 endif
00227 #endif
00228
00229
00230
00231
00232
00233 call psmile_bsend (buffer, nloc*nbr_fields, MPI_DOUBLE_PRECISION, &
00234 dest, tag, comm, ierror)
00235
00236 if (ierror /= MPI_SUCCESS) then
00237 ierrp (1) = ierror
00238 ierrp (2) = dest
00239 ierrp (3) = tag
00240
00241 ierror = PRISM_Error_Send
00242
00243 call psmile_error (ierror, 'psmile_bsend', &
00244 ierrp, 3, __FILE__, __LINE__ )
00245 return
00246 endif
00247
00248
00249
00250 Deallocate (buffer)
00251
00252
00253
00254 #ifdef VERBOSE
00255 print 9980, trim(ch_id), ierror
00256
00257 call psmile_flushstd
00258 #endif /* VERBOSE */
00259
00260
00261
00262
00263
00264 #ifdef VERBOSE
00265
00266 9990 format (1x, a, ': psmile_put_irr_field_dble: dest', i3)
00267 9980 format (1x, a, ': psmile_put_irr_field_dble: eof ierror =', i3)
00268
00269 #endif /* VERBOSE */
00270
00271 end subroutine PSMILe_Put_irr_field_dble