00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 subroutine psmile_put_field_21d_real (data_array, shape, nbr_fields, &
00013 srclocs, nparts, nloc, npoints, &
00014 dest, tag, comm, ierror)
00015
00016
00017
00018 use PRISM_constants
00019
00020 use PSMILe, dummy_interface => PSMILe_Put_field_21d_real
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 Real, 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 Integer, Intent (In) :: nparts
00044
00045
00046
00047 Integer, Intent (In) :: npoints (2, nparts)
00048
00049
00050
00051
00052
00053
00054 Type(integer_vector), Intent (In) :: srclocs (2, nparts)
00055
00056
00057
00058
00059 Integer, Intent (In) :: dest
00060
00061
00062
00063 Integer, Intent (In) :: tag
00064
00065
00066
00067 Integer, Intent (In) :: comm
00068
00069
00070
00071
00072
00073 Integer, Intent (Out) :: ierror
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083 Real, Pointer :: buffer (:)
00084 Integer :: i, j, n
00085
00086
00087
00088 Integer :: ipart, nprev
00089
00090
00091
00092 Integer, parameter :: nerrp = 3
00093 Integer :: ierrp (nerrp)
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118 Character(len=len_cvs_string), save :: mycvs =
00119 '$Id: psmile_put_field_21d_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00120
00121
00122
00123
00124
00125 #ifdef VERBOSE
00126 print 9990, trim(ch_id), dest
00127
00128 call psmile_flushstd
00129 #endif /* VERBOSE */
00130
00131 ierror = 0
00132
00133 #ifdef PRISM_ASSERTION
00134 nprev = 0
00135 do ipart = 1, nparts
00136 nprev = nprev + npoints(1,ipart)*npoints (2,ipart)
00137 end do
00138
00139 if (nloc /= nprev) then
00140 print *, "nloc, ntot, npoints", nloc, nprev, npoints
00141 call psmile_assert ( __FILE__, __LINE__, &
00142 "nloc /= Sum(npoints(1,i)*npoints(2,i))")
00143 endif
00144 #endif
00145
00146
00147
00148
00149 Allocate (buffer(nloc*nbr_fields), STAT = ierror)
00150 if ( ierror > 0 ) then
00151 ierrp (1) = ierror
00152 ierrp (2) = nloc * nbr_fields
00153
00154 ierror = PRISM_Error_Alloc
00155
00156 call psmile_error ( ierror, 'buffer', &
00157 ierrp, 2, __FILE__, __LINE__ )
00158 return
00159 endif
00160
00161
00162
00163 #ifdef DEBUG
00164 print 9960, nloc, nbr_fields
00165 9960 format (1x, 'psmile_put_field_21d_real: pack: nloc, nbr_fields', i7, i4)
00166 #endif
00167
00168 do n = 1, nbr_fields
00169
00170 nprev = (n-1)*nloc
00171 do ipart = 1, nparts
00172 do j = 1, npoints (2, ipart)
00173
00174 do i = 1, npoints (1, ipart)
00175 buffer (nprev+(j-1)*npoints(1, ipart)+i) = &
00176 data_array (srclocs(1,ipart)%vector (i*2-1), &
00177 srclocs(1,ipart)%vector (i*2 ), &
00178 srclocs(2,ipart)%vector (j), n)
00179 end do
00180 end do
00181
00182 nprev = nprev + npoints(1,ipart)*npoints(2,ipart)
00183 end do
00184 end do
00185
00186 #ifdef DEBUG
00187 ipart = 1
00188 print 9950, srclocs(1,ipart)%vector(1:2), srclocs(2,ipart)%vector(1), &
00189 data_array(srclocs(1,ipart)%vector (1), &
00190 srclocs(1,ipart)%vector (2), &
00191 srclocs(2,ipart)%vector (1), 1)
00192 9950 format (1x, 'psmile_put_field_21d_real: srcijk', 3i5, ' val', f10.5)
00193 #endif
00194
00195
00196
00197
00198
00199 #ifdef DEBUG
00200 print 9970, nloc, dest, tag, comm
00201 9970 format (1x, 'psmile_put_field_21d_real: nloc', i7, '->', i2, ', tag', i5, ', comm', i5)
00202 #endif /* DEBUG */
00203
00204 call psmile_bsend (buffer, nloc*nbr_fields, MPI_REAL, &
00205 dest, tag, comm, ierror)
00206
00207 if (ierror /= MPI_SUCCESS) then
00208 ierrp (1) = ierror
00209 ierrp (2) = dest
00210 ierrp (3) = tag
00211
00212 ierror = PRISM_Error_Send
00213
00214 call psmile_error (ierror, 'psmile_bsend', &
00215 ierrp, 3, __FILE__, __LINE__ )
00216 return
00217 endif
00218
00219
00220
00221 Deallocate (buffer)
00222
00223
00224
00225 #ifdef VERBOSE
00226 print 9980, trim(ch_id), ierror
00227
00228 call psmile_flushstd
00229 #endif /* VERBOSE */
00230
00231
00232
00233
00234
00235 #ifdef VERBOSE
00236
00237 9990 format (1x, a, ': psmile_put_field_21d_real: dest', i3)
00238 9980 format (1x, a, ': psmile_put_field_21d_real: eof ierror =', i3)
00239
00240 #endif /* VERBOSE */
00241
00242 end subroutine PSMILe_Put_field_21d_real