File: test_memset.jf90

package info (click to toggle)
devicexlib 0.8.6-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 10,364 kB
  • sloc: f90: 77,678; sh: 3,701; fortran: 773; makefile: 268; python: 246; ansic: 69; awk: 36
file content (258 lines) | stat: -rw-r--r-- 8,726 bytes parent folder | download
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
246
247
248
249
250
251
252
253
254
255
256
257
258
!
! Copyright (C) 2023, MaX CoE
! Distributed under the MIT License 
! (license terms are at http://opensource.org/licenses/MIT).
!
!--
!
! Testing: devxlib_memset, devxlib_memcpy
!
!==================================================================
!==================================================================
! *DO NOT EDIT*: automatically generated from test_memcpy.jf90
!==================================================================
!==================================================================
!
#include<devxlib_macros.h>
#include<devxlib_defs.h>
!
program test_memset
  !
  ! This program tests the routines and interfaces related to
  ! devxlib_memset and devxlib_memcpy
  !
  use iso_fortran_env, only : int32, int64, real32, real64
  use devxlib_memset
  use devxlib_memcpy
  use devxlib_mapping

  implicit none

  integer(int64) :: cnt, count_max
  real(real64), parameter :: thr=1.0d-6
  integer,      parameter :: nranks={{dimensions}}
  integer :: ndim(nranks)
  integer :: vrange(2,nranks)
  integer :: vlbound(nranks)
  real(real64):: t0, t1, count_rate
  logical :: is_alloc, lfail_malloc
  integer :: npass,nfail
  integer :: nfail_malloc
  !
  integer(int32)  :: i_const_int32
  integer(int64)  :: i_const_int64
  real(real32)    :: r_const_real32
  real(real64)    :: r_const_real64
  complex(real32) :: c_const_real32
  complex(real64) :: c_const_real64
  logical         :: l_const


{%- for t in types %}
{%- for p in kinds[t] %}
{%- for d in range(1,dimensions+1) %}
  {{t}}({{p.val}}), allocatable :: A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d({% for dd in range(d) %}:{% if not loop.last %}, {%- endif %}{% endfor %})
  {{t}}({{p.val}}), allocatable :: A_hst1__{{p.name}}_{{t[0]|lower}}{{d}}d({% for dd in range(d) %}:{% if not loop.last %}, {%- endif %}{% endfor %})
  {{t}}({{p.val}}), allocatable :: A_hst2__{{p.name}}_{{t[0]|lower}}{{d}}d({% for dd in range(d) %}:{% if not loop.last %}, {%- endif %}{% endfor %})
  {{t}}({{p.val}}), allocatable DEV_ATTR :: A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d({% for dd in range(d) %}:{% if not loop.last %}, {%- endif %}{% endfor %})
{% endfor %}
{% endfor %}
{% endfor %}
  
  integer :: i,ierr
  integer :: {% for dd in range(1,dimensions+1) %}ndim{{dd}}{% if not loop.last %}, {%- endif %} {% endfor %}
  integer :: {% for dd in range(1,dimensions+1) %}lbound{{dd}}{% if not loop.last %}, {%- endif %} {% endfor %}
  integer :: {% for dd in range(1,dimensions+1) %}range{{dd}}(2){% if not loop.last %}, {%- endif %} {% endfor %}
  integer :: {% for dd in range(1,dimensions+1) %}bound{{dd}}(2){% if not loop.last %}, {%- endif %} {% endfor %}
  character(256) :: arg, str

!
!============================
! get dims
!============================
!
  ! defaults
  ndim(:)=100
  vrange(1,:)=1
  vrange(2,:)=ndim
  vlbound(:)=1

  i=0
  do
    call get_command_argument(i, arg)
    if (len_trim(arg) == 0) exit
    !
    select case (trim(arg))
    case("-h","--help")
      write(6,"(a)") "Usage: "
      write(6,"(a)") "   ./test_memset.x [--dims <vals>] [--range <vals>] [--lbound <vals>]"
      stop
    end select
    !
    i = i+1
    call get_command_argument(i, str)
    if (len_trim(str) == 0) exit
    !
    select case (trim(arg))
    case("-dims","--dims")
      read(str,*,iostat=ierr) ndim(:)
      if (ierr/=0) STOP "reading cmd-line args: dims"
    case("-range","--range")
      read(str,*,iostat=ierr) vrange(:,:)
      if (ierr/=0) STOP "reading cmd-line args: range"
    case("-lbound","--lbound")
      read(str,*,iostat=ierr) vlbound(:)
      if (ierr/=0) STOP "reading cmd-line args: lbound"
    end select
  enddo
  !
  write(6,"(/,a,/)") "Running test with params: "
  write(6,"(3x,a,10i5)") "  ndim: ", ndim(:)
  write(6,"(3x,a,10i5)") "lbound: ", vlbound(:)
  do i = 1, nranks
     write(6,"(3x,a,i2,3x,10i5)") " range: ", i, vrange(:,i)
  enddo
  write(6,"()")
  !
  npass=0
  nfail=0
  nfail_malloc=0

{% for d in range(1,dimensions+1) %}
  ndim{{d}}=ndim({{d}})
  lbound{{d}}=vlbound({{d}})
  range{{d}}=vrange(:,{{d}})
  bound{{d}}(1)=lbound{{d}}
  bound{{d}}(2)=lbound{{d}}+ndim{{d}}-1
{%- endfor %}

!
!============================
! check memcpy
!============================
!

{%- for t in types %}
{%- for p in kinds[t] %}
{%- for d in range(1,dimensions+1) %}
  !
  !=====================
  write(6,"(/,3x,a)") "checking {{p.name}}_{{t[0]|lower}}{{d}}d ..." 
  call system_clock(cnt, count_rate, count_max)
  t0=cnt/count_rate
  lfail_malloc=.false.
  !
  is_alloc = devxlib_mapped( A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d )
  !
  if (is_alloc) then
     nfail_malloc=nfail_malloc+1
     lfail_malloc=.true.
     write(6,"(3x,a)") "mem allocation check-pre FAILED" 
  endif
  !
  ! allocations
  allocate( A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d({% for dd in range(d)%}bound{{dd+1}}(1):bound{{dd+1}}(2){% if not loop.last %},{%- endif %}{%endfor%}) )
  allocate( A_hst1__{{p.name}}_{{t[0]|lower}}{{d}}d({% for dd in range(d)%}bound{{dd+1}}(1):bound{{dd+1}}(2){% if not loop.last %},{%- endif %}{%endfor%}) )
  allocate( A_hst2__{{p.name}}_{{t[0]|lower}}{{d}}d({% for dd in range(d)%}bound{{dd+1}}(1):bound{{dd+1}}(2){% if not loop.last %},{%- endif %}{%endfor%}) )
  allocate( A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d({% for dd in range(d)%}bound{{dd+1}}(1):bound{{dd+1}}(2){% if not loop.last %},{%- endif %}{%endfor%}) )
#if defined __DXL_OPENACC || defined __DXL_OPENMP_GPU
  call devxlib_map(A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d)
#endif
  !
  is_alloc = devxlib_mapped( A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d )
  !
  if (.not.is_alloc) then
    nfail_malloc=nfail_malloc+1
    lfail_malloc=.true.
    write(6,"(3x,a)") "mem allocation check-post FAILED" 
  endif
  !
  ! init
  {%- if t=="complex" %}
    c_const_{{p.val}} = (1.0,2.0)
  {% elif t=="integer" %}
    i_const_{{p.val}} = 3
  {% elif t=="logical" %}
    l_const = .true.
  {%else%}
    r_const_{{p.val}} = 4.0
  {%endif%}
  !
  ! memset
  !
  {%- if t=="complex" %}
    A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d=c_const_{{p.val}}
    call devxlib_memset_h(A_hst1__{{p.name}}_{{t[0]|lower}}{{d}}d,c_const_{{p.val}})
    call devxlib_memset_d(A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d,c_const_{{p.val}})
  {% elif t=="integer" %}
    A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d=i_const_{{p.val}}
    call devxlib_memset_h(A_hst1__{{p.name}}_{{t[0]|lower}}{{d}}d,i_const_{{p.val}})
    call devxlib_memset_d(A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d,i_const_{{p.val}})
  {% elif t=="logical" %}
    A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d=l_const
    call devxlib_memset_h(A_hst1__{{p.name}}_{{t[0]|lower}}{{d}}d,l_const)
    call devxlib_memset_d(A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d,l_const)
  {%else%}
    A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d=r_const_{{p.val}}
    call devxlib_memset_h(A_hst1__{{p.name}}_{{t[0]|lower}}{{d}}d,r_const_{{p.val}})
    call devxlib_memset_d(A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d,r_const_{{p.val}})
  {%endif%}
  ! 
  ! mem copy d2h
  call devxlib_memcpy_d2h(A_hst2__{{p.name}}_{{t[0]|lower}}{{d}}d, A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d)
  !
  ! check
  if (lfail_malloc) write(6,"(3x,a)") "Malloc FAILED" 
  !
  {%- if t=="logical" %}
  if ( any( A_hst1__{{p.name}}_{{t[0]|lower}}{{d}}d .neqv. A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d ) .or. &
       any( A_hst2__{{p.name}}_{{t[0]|lower}}{{d}}d .neqv. A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d )    ) then
  {%- else %}
  if ( any(abs(A_hst1__{{p.name}}_{{t[0]|lower}}{{d}}d -A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d )> thr) .or. &
       any(abs(A_hst2__{{p.name}}_{{t[0]|lower}}{{d}}d -A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d )> thr) ) then
  {%- endif %}
     !
     write(6,"(3x,a)") "FAILED" 
     nfail=nfail+1
  else
     write(6,"(3x,a)") "passed" 
     npass=npass+1
  endif
  !
  deallocate(A_hst0__{{p.name}}_{{t[0]|lower}}{{d}}d)
  deallocate(A_hst1__{{p.name}}_{{t[0]|lower}}{{d}}d)
  deallocate(A_hst2__{{p.name}}_{{t[0]|lower}}{{d}}d)
  !
#if defined __DXL_OPENACC || defined __DXL_OPENMP_GPU
  call devxlib_unmap(A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d)
  is_alloc = devxlib_mapped( A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d )
  !
  if (is_alloc) then
    nfail_malloc=nfail_malloc+1
    lfail_malloc=.true.
    write(6,"(3x,a)") "mem allocation check-post2 FAILED" 
  endif
  !
#endif
  !
  deallocate(A_dev__{{p.name}}_{{t[0]|lower}}{{d}}d)
  !
  call system_clock(cnt, count_rate, count_max)
  t1=cnt/count_rate
  write(6,"(3x,a,f12.6,' sec')") "Timing: ", t1-t0
  !
{% endfor %}
{% endfor %}
{% endfor %}

  !
  ! summary
  !
  write(6,"(/,a)") "# Test SUMMARY:"
  write(6,"(3x,a,i5)") "# passed: ", npass
  write(6,"(3x,a,i5)") "# failed: ", nfail
  write(6,"(3x,a,i5)") "# failed malloc: ", nfail_malloc
  write(6,"()")

end program test_memset