1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245
|
subroutine m4_func_NGA_ADD_PATCH(m4_test_type, m4_ndim)
implicit none
#include "mafdecls.fh"
#include "global.fh"
c
integer n,m
integer ndim
parameter (n = m4_n)
parameter (m = (m4_n**m4_ndim)/100)
parameter (ndim = m4_ndim)
m4_data_type a(substr(m4_array, 1, eval(m4_ndim*2-1)))
m4_data_type b(substr(m4_array, 1, eval(m4_ndim*2-1)))
integer dims(ndim)
integer g_a, g_b,g_c
integer chunk(ndim)
integer i, total
integer elems, count_elems
integer loop
double precision drand
m4_data_type alpha, beta
integer lop(ndim), hip(ndim), hipl(ndim)
integer alo(ndim), ahi(ndim)
integer blo(ndim), bhi(ndim)
integer clo(ndim), chi(ndim)
c for different array dimensions
ifelse(m4_ndim,1,`',`
m4_data_type d(substr(m4_array, 1, eval((m4_ndim-1)*2-1)))
integer dndim
parameter (dndim = m4_ndim-1)
integer ddims(dndim),dlo(dndim),dhi(dndim),dtotal
')
c
integer nproc, me
logical status
integer repeat
c
nproc = ga_nnodes()
me = ga_nodeid()
c
c---------------------- initialize the GA -----------------------
c initialize the chunk, dims, ld, and calculate the number
c of elements
total=1
do i = 1,ndim
chunk(i) = 0
dims(i) = n
total = total * dims(i)
enddo
c
c*** Create global arrays
if (.not. nga_create(m4_MT, ndim, dims, 'a', chunk, g_a))
$ call ga_error(' ga_create failed ',1)
c
c test the same distribution and different distribution seperately
do repeat=1,2
if(repeat.eq.1) then
status = ga_duplicate(g_a, g_b, 'a_duplicated')
if(.not.ga_compare_distr(g_a, g_b))
$ call ga_error("g_b distribution different",0)
c
status = ga_duplicate(g_a, g_c, 'a_duplicated_again')
if(.not.ga_compare_distr(g_a, g_c))
$ call ga_error("g_c distribution different",0)
c
else
do i = 1,ndim
if(mod(i,2).eq.0) chunk(i) = n
enddo
if (.not. nga_create(m4_MT, ndim, dims, 'b', chunk, g_b))
$ call ga_error(' ga_create failed ',1)
do i = 1,ndim
if(mod(i,2).eq.0) chunk(i) = 0
if(mod(i,2).eq.1) chunk(i) = n
enddo
if (.not. nga_create(m4_MT, ndim, dims, 'c', chunk, g_c))
$ call ga_error(' ga_create failed ',1)
endif
c
call ga_sync()
c
c---------------------------NGA_ADD_PATCH -------------------------
c
if(repeat.eq.1) then
m4_print_info(nga_add_patch)
if(me.eq.0) print *, 'Testing with the same distributions'
else
if(me.eq.0) print *, 'Testing with different distributions'
endif
c
c initialize GA
call m4_util_init_array(m4_test_type)(a,total)
call nga_distribution(g_a, me, lop, hip)
elems = count_elems(lop, hip, ndim)
if(elems.gt.0) call nga_put(g_a,lop,hip,
$ a(substr(m4_lop_all, 1, eval(m4_ndim*7-1))),dims)
call m4_util_init_array(m4_test_type)(b,total)
call nga_distribution(g_b, me, lop, hip)
elems = count_elems(lop, hip, ndim)
if(elems.gt.0) call nga_put(g_b,lop,hip,
$ b(substr(m4_lop_all, 1, eval(m4_ndim*7-1))),dims)
c
call ga_sync()
do i = 1,ndim
lop(i) = 1
hipl(i) = n-2
hip(i) = n
enddo
c
c---
do loop=1, 10
call random_range(lop,hipl,alo,ahi,ndim)
do i=1, ndim
blo(i) = alo(i) + 1
bhi(i) = ahi(i) + 1
clo(i) = alo(i) + 2
chi(i) = ahi(i) + 2
enddo
if(me.eq.0)then
call add_range(loop,alo,ahi,ndim,blo,bhi,ndim)
c$$$ print *, loop,': [',(alo(i),':',ahi(i), i=1,ndim),']',
c$$$ $ '+', '[',(blo(i),':',bhi(i), i=1,ndim),']'
endif
c
alpha = m4_rand(1)
beta = m4_rand(1)
c
c keep copies of the origian arrays
call nga_get(g_a,alo,ahi,
$ a(substr(m4_alo_all, 1, eval(m4_ndim*7-1))),dims)
call nga_get(g_b,blo,bhi,
$ b(substr(m4_blo_all, 1, eval(m4_ndim*7-1))),dims)
c the result should be (in a)
call m4_util_scale_patch(m4_test_type)(total,
$ alpha,a,alo,ahi,ndim,dims,
$ beta,b,blo,bhi,ndim,dims)
c
call nga_add_patch(alpha, g_a, alo, ahi, beta, g_b, blo, bhi,
$ g_c, clo,chi)
c
call nga_get(g_c,clo,chi,
$ b(substr(m4_clo_all, 1, eval(m4_ndim*7-1))),dims)
c
call m4_util_compare_patches(m4_test_type)(1d-5,total,
$ a,alo,ahi,ndim,dims,total,b,clo,chi,ndim,dims)
enddo
c
call ga_sync()
if(me.eq.0)then
print *, 'OK'
print *, ' '
call ffflush(6)
endif
c---------------------------
c
status = ga_destroy(g_b)
enddo
c-----------------------------------------------------------------
changequote({,})
ifelse(m4_ndim,1,{},{
c testing copy on differet dimensions
dtotal = 1
do i = 1,dndim
ddims(i) = n
dtotal = dtotal * ddims(i)
enddo
c
if (.not. nga_create(m4_MT, dndim, ddims, 'd', chunk, g_b))
$ call ga_error(' ga_create failed ',1)
c
if(me.eq.0)
$ print *, 'Testing adding patch on different dimensions'
c
call ga_sync()
c
c initialize g_b
call m4_util_init_array(m4_test_type)(d,dtotal)
call nga_distribution(g_b, me, dlo, dhi)
elems = count_elems(dlo, dhi, dndim)
if(elems.gt.0) call nga_put(g_b,dlo,dhi,
$ d(substr(m4_dlo_all, 1, eval((m4_ndim-1)*7-1))),ddims)
c
c calculate the maximum range of g_a that can fit into g_b
do i = 1,ndim
lop(i) = 1
hip(i) = n
enddo
hip(dndim) = 1
c
call ga_sync()
c
do loop=1, 10
call random_range(lop,hip,alo,ahi,ndim)
c
do i=1, dndim
dlo(i) = alo(dndim-i+1)
dhi(i) = ahi(dndim-i+1)
enddo
dlo(1) = alo(ndim)
dhi(1) = ahi(ndim)
c
if(me.eq.0) then
call add_range(loop,alo,ahi,ndim,dlo,dhi,dndim)
c$$$ print *, loop,': [',(alo(i),':',ahi(i), i=1,ndim),']',
c$$$ $ '+', '[',(dlo(i),':',dhi(i), i=1,dndim),']'
endif
c
alpha = m4_rand(1)
beta = m4_rand(1)
c
c keep copies of the origian arrays
call nga_get(g_a,alo,ahi,
$ a(substr(m4_alo_all, 1, eval(m4_ndim*7-1))),dims)
call nga_get(g_b,dlo,dhi,
$ d(substr(m4_dlo_all, 1, eval((m4_ndim-1)*7-1))),ddims)
c
c the result should be (in a)
call m4_util_scale_patch(m4_test_type)(total,
$ alpha,a,alo,ahi,ndim,dims,
$ beta,d,dlo,dhi,dndim,ddims)
c
call nga_add_patch(alpha,g_a,alo,ahi,beta,g_b,dlo,dhi,
$ g_c,alo,ahi)
c
call nga_get(g_c,alo,ahi,
$ b(substr(m4_alo_all, 1, eval(m4_ndim*7-1))),dims)
c
call m4_util_compare_patches(m4_test_type)(1d-5,total,
$ a,alo,ahi,ndim,dims,total,b,alo,ahi,ndim,dims)
enddo
c
call ga_sync()
if(me.eq.0)then
print *, ' add patches on different dimensions: OK'
print *, ' '
call ffflush(6)
endif
c
status = ga_destroy(g_b)
})
changequote(`,')
c---
c
status = ga_destroy(g_c)
status = ga_destroy(g_a)
end
|