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