File: testing_rotg_strided_batched.hpp

package info (click to toggle)
rocblas 5.5.1%2Bdfsg-7
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 565,372 kB
  • sloc: cpp: 198,491; python: 44,792; f90: 25,111; sh: 24,429; asm: 8,954; xml: 222; makefile: 147; ansic: 107; awk: 14
file content (315 lines) | stat: -rw-r--r-- 13,732 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
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
/* ************************************************************************
 * Copyright (C) 2018-2023 Advanced Micro Devices, Inc. All rights reserved.
 *
 * Permission is hereby granted, free of charge, to any person obtaining a copy
 * of this software and associated documentation files (the "Software"), to deal
 * in the Software without restriction, including without limitation the rights
 * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell cop-
 * ies of the Software, and to permit persons to whom the Software is furnished
 * to do so, subject to the following conditions:
 *
 * The above copyright notice and this permission notice shall be included in all
 * copies or substantial portions of the Software.
 *
 * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IM-
 * PLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS
 * FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR
 * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER
 * IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNE-
 * CTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
 *
 * ************************************************************************ */

#pragma once

#include "cblas_interface.hpp"
#include "norm.hpp"
#include "rocblas.hpp"
#include "rocblas_init.hpp"
#include "rocblas_math.hpp"
#include "rocblas_random.hpp"
#include "rocblas_test.hpp"
#include "rocblas_vector.hpp"
#include "unit.hpp"
#include "utility.hpp"

template <typename T, typename U = T>
void testing_rotg_strided_batched_bad_arg(const Arguments& arg)
{
    // clang-format off
    const bool FORTRAN                         = arg.fortran;
    auto       rocblas_rotg_strided_batched_fn = FORTRAN ? rocblas_rotg_strided_batched<T, U, true>
                                                         : rocblas_rotg_strided_batched<T, U, false>;
    // clang-format on

    rocblas_int    batch_count = 5;
    rocblas_stride stride_a    = 10;
    rocblas_stride stride_b    = 10;
    rocblas_stride stride_c    = 10;
    rocblas_stride stride_s    = 10;

    rocblas_local_handle handle{arg};

    // Allocate device memory
    device_strided_batch_vector<T> da(1, 1, stride_a, batch_count);
    device_strided_batch_vector<T> db(1, 1, stride_b, batch_count);
    device_strided_batch_vector<U> dc(1, 1, stride_c, batch_count);
    device_strided_batch_vector<T> ds(1, 1, stride_s, batch_count);

    // Check device memory allocation
    CHECK_DEVICE_ALLOCATION(da.memcheck());
    CHECK_DEVICE_ALLOCATION(db.memcheck());
    CHECK_DEVICE_ALLOCATION(dc.memcheck());
    CHECK_DEVICE_ALLOCATION(ds.memcheck());

    EXPECT_ROCBLAS_STATUS(
        (rocblas_rotg_strided_batched_fn(
            nullptr, da, stride_a, db, stride_b, dc, stride_c, ds, stride_s, batch_count)),
        rocblas_status_invalid_handle);
    EXPECT_ROCBLAS_STATUS(
        (rocblas_rotg_strided_batched_fn(
            handle, nullptr, stride_a, db, stride_b, dc, stride_c, ds, stride_s, batch_count)),
        rocblas_status_invalid_pointer);
    EXPECT_ROCBLAS_STATUS(
        (rocblas_rotg_strided_batched_fn(
            handle, da, stride_a, nullptr, stride_b, dc, stride_c, ds, stride_s, batch_count)),
        rocblas_status_invalid_pointer);
    EXPECT_ROCBLAS_STATUS(
        (rocblas_rotg_strided_batched_fn(
            handle, da, stride_a, db, stride_b, nullptr, stride_c, ds, stride_s, batch_count)),
        rocblas_status_invalid_pointer);
    EXPECT_ROCBLAS_STATUS(
        (rocblas_rotg_strided_batched_fn(
            handle, da, stride_a, db, stride_b, dc, stride_c, nullptr, stride_s, batch_count)),
        rocblas_status_invalid_pointer);
}

template <typename T, typename U = T>
void testing_rotg_strided_batched(const Arguments& arg)
{
    // clang-format off
    const bool FORTRAN                         = arg.fortran;
    auto       rocblas_rotg_strided_batched_fn = FORTRAN ? rocblas_rotg_strided_batched<T, U, true>
                                                         : rocblas_rotg_strided_batched<T, U, false>;
    // clang-format on

    rocblas_int stride_a    = arg.stride_a;
    rocblas_int stride_b    = arg.stride_b;
    rocblas_int stride_c    = arg.stride_c;
    rocblas_int stride_s    = arg.stride_d;
    rocblas_int batch_count = arg.batch_count;

    rocblas_local_handle handle{arg};
    double               gpu_time_used, cpu_time_used;
    double               norm_error_host = 0.0, norm_error_device = 0.0;
    const U              rel_error = std::numeric_limits<U>::epsilon() * 100;

    // check to prevent undefined memory allocation error
    if(batch_count <= 0)
    {
        CHECK_ROCBLAS_ERROR(rocblas_set_pointer_mode(handle, rocblas_pointer_mode_device));
        EXPECT_ROCBLAS_STATUS((rocblas_rotg_strided_batched_fn(handle,
                                                               nullptr,
                                                               stride_a,
                                                               nullptr,
                                                               stride_b,
                                                               nullptr,
                                                               stride_c,
                                                               nullptr,
                                                               stride_s,
                                                               batch_count)),
                              rocblas_status_success);
        return;
    }

    host_strided_batch_vector<T> ha(1, 1, stride_a, batch_count);
    host_strided_batch_vector<T> hb(1, 1, stride_b, batch_count);
    host_strided_batch_vector<U> hc(1, 1, stride_c, batch_count);
    host_strided_batch_vector<T> hs(1, 1, stride_s, batch_count);

    bool enable_near_check_general = true;

    // Initialize data on host memory
    rocblas_init_vector(ha, arg, rocblas_client_never_set_nan, true);
    rocblas_init_vector(hb, arg, rocblas_client_never_set_nan, false);
    rocblas_init_vector(hc, arg, rocblas_client_never_set_nan, false);
    rocblas_init_vector(hs, arg, rocblas_client_never_set_nan, false);

    ha[0][0] = arg.get_alpha<T>(); // reuse alpha in place of a to keep number of arguments small
    hb[0][0] = arg.get_beta<T>(); // reuse beta  in place of a to keep number of arguments small
    hc[0][0] = U(0);
    hs[0][0] = T(0);

    // CPU_BLAS
    host_strided_batch_vector<T> ha_gold(1, 1, stride_a, batch_count);
    host_strided_batch_vector<T> hb_gold(1, 1, stride_b, batch_count);
    host_strided_batch_vector<U> hc_gold(1, 1, stride_c, batch_count);
    host_strided_batch_vector<T> hs_gold(1, 1, stride_s, batch_count);

    ha_gold.copy_from(ha);
    hb_gold.copy_from(hb);
    hc_gold.copy_from(hc);
    hs_gold.copy_from(hs);

    cpu_time_used = get_time_us_no_sync();
    for(int b = 0; b < batch_count; b++)
    {
        cblas_rotg<T, U>(ha_gold[b], hb_gold[b], hc_gold[b], hs_gold[b]);
    }
    cpu_time_used = get_time_us_no_sync() - cpu_time_used;

    // Test rocblas_pointer_mode_host
    {
        host_strided_batch_vector<T> ra(1, 1, stride_a, batch_count);
        host_strided_batch_vector<T> rb(1, 1, stride_b, batch_count);
        host_strided_batch_vector<U> rc(1, 1, stride_c, batch_count);
        host_strided_batch_vector<T> rs(1, 1, stride_s, batch_count);

        ra.copy_from(ha);
        rb.copy_from(hb);
        rc.copy_from(hc);
        rs.copy_from(hs);

        CHECK_ROCBLAS_ERROR(rocblas_set_pointer_mode(handle, rocblas_pointer_mode_host));
        handle.pre_test(arg);
        CHECK_ROCBLAS_ERROR((rocblas_rotg_strided_batched_fn(
            handle, ra, stride_a, rb, stride_b, rc, stride_c, rs, stride_s, batch_count)));
        handle.post_test(arg);

        if(arg.unit_check)
        {
            if(enable_near_check_general)
            {
                near_check_general<T>(1, 1, 1, stride_a, ha_gold, ra, batch_count, rel_error);
                near_check_general<T>(1, 1, 1, stride_b, hb_gold, rb, batch_count, rel_error);
                near_check_general<U>(1, 1, 1, stride_c, hc_gold, rc, batch_count, rel_error);
                near_check_general<T>(1, 1, 1, stride_s, hs_gold, rs, batch_count, rel_error);
            }
        }

        if(arg.norm_check)
        {
            norm_error_host
                = norm_check_general<T>('F', 1, 1, 1, stride_a, ha_gold, ra, batch_count);
            norm_error_host
                += norm_check_general<T>('F', 1, 1, 1, stride_b, hb_gold, rb, batch_count);
            norm_error_host
                += norm_check_general<U>('F', 1, 1, 1, stride_c, hc_gold, rc, batch_count);
            norm_error_host
                += norm_check_general<T>('F', 1, 1, 1, stride_s, hs_gold, rs, batch_count);
        }
    }

    // Test rocblas_pointer_mode_device
    {
        // Allocate device memory
        device_strided_batch_vector<T> da(1, 1, stride_a, batch_count);
        device_strided_batch_vector<T> db(1, 1, stride_b, batch_count);
        device_strided_batch_vector<U> dc(1, 1, stride_c, batch_count);
        device_strided_batch_vector<T> ds(1, 1, stride_s, batch_count);

        // Check device memory allocation
        CHECK_DEVICE_ALLOCATION(da.memcheck());
        CHECK_DEVICE_ALLOCATION(db.memcheck());
        CHECK_DEVICE_ALLOCATION(dc.memcheck());
        CHECK_DEVICE_ALLOCATION(ds.memcheck());

        // Transfer from CPU to GPU
        CHECK_HIP_ERROR(da.transfer_from(ha));
        CHECK_HIP_ERROR(db.transfer_from(hb));
        CHECK_HIP_ERROR(dc.transfer_from(hc));
        CHECK_HIP_ERROR(ds.transfer_from(hs));

        CHECK_ROCBLAS_ERROR(rocblas_set_pointer_mode(handle, rocblas_pointer_mode_device));
        handle.pre_test(arg);
        CHECK_ROCBLAS_ERROR((rocblas_rotg_strided_batched_fn(
            handle, da, stride_a, db, stride_b, dc, stride_c, ds, stride_s, batch_count)));
        handle.post_test(arg);

        host_strided_batch_vector<T> ra(1, 1, stride_a, batch_count);
        host_strided_batch_vector<T> rb(1, 1, stride_b, batch_count);
        host_strided_batch_vector<U> rc(1, 1, stride_c, batch_count);
        host_strided_batch_vector<T> rs(1, 1, stride_s, batch_count);

        // Transfer from GPU to CPU
        CHECK_HIP_ERROR(ra.transfer_from(da));
        CHECK_HIP_ERROR(rb.transfer_from(db));
        CHECK_HIP_ERROR(rc.transfer_from(dc));
        CHECK_HIP_ERROR(rs.transfer_from(ds));

        if(arg.unit_check)
        {
            if(enable_near_check_general)
            {
                near_check_general<T>(1, 1, 1, stride_a, ha_gold, ra, batch_count, rel_error);
                near_check_general<T>(1, 1, 1, stride_b, hb_gold, rb, batch_count, rel_error);
                near_check_general<U>(1, 1, 1, stride_c, hc_gold, rc, batch_count, rel_error);
                near_check_general<T>(1, 1, 1, stride_s, hs_gold, rs, batch_count, rel_error);
            }
        }

        if(arg.norm_check)
        {
            norm_error_device
                = norm_check_general<T>('F', 1, 1, 1, stride_a, ha_gold, ra, batch_count);
            norm_error_device
                += norm_check_general<T>('F', 1, 1, 1, stride_b, hb_gold, rb, batch_count);
            norm_error_device
                += norm_check_general<U>('F', 1, 1, 1, stride_c, hc_gold, rc, batch_count);
            norm_error_device
                += norm_check_general<T>('F', 1, 1, 1, stride_s, hs_gold, rs, batch_count);
        }
    }

    if(arg.timing)
    {
        int number_cold_calls = arg.cold_iters;
        int number_hot_calls  = arg.iters;
        // Device mode will be quicker
        // (TODO: or is there another reason we are typically using host_mode for timing?)
        CHECK_ROCBLAS_ERROR(rocblas_set_pointer_mode(handle, rocblas_pointer_mode_device));

        // Allocate device memory
        device_strided_batch_vector<T> da(1, 1, stride_a, batch_count);
        device_strided_batch_vector<T> db(1, 1, stride_b, batch_count);
        device_strided_batch_vector<U> dc(1, 1, stride_c, batch_count);
        device_strided_batch_vector<T> ds(1, 1, stride_s, batch_count);

        // Check device memory allocation
        CHECK_DEVICE_ALLOCATION(da.memcheck());
        CHECK_DEVICE_ALLOCATION(db.memcheck());
        CHECK_DEVICE_ALLOCATION(dc.memcheck());
        CHECK_DEVICE_ALLOCATION(ds.memcheck());

        // Transfer from CPU to GPU
        CHECK_HIP_ERROR(da.transfer_from(ha));
        CHECK_HIP_ERROR(db.transfer_from(hb));
        CHECK_HIP_ERROR(dc.transfer_from(hc));
        CHECK_HIP_ERROR(ds.transfer_from(hs));

        for(int iter = 0; iter < number_cold_calls; iter++)
        {
            rocblas_rotg_strided_batched_fn(
                handle, da, stride_a, db, stride_b, dc, stride_c, ds, stride_s, batch_count);
        }
        hipStream_t stream;
        CHECK_ROCBLAS_ERROR(rocblas_get_stream(handle, &stream));
        gpu_time_used = get_time_us_sync(stream); // in microseconds
        for(int iter = 0; iter < number_hot_calls; iter++)
        {
            rocblas_rotg_strided_batched_fn(
                handle, da, stride_a, db, stride_b, dc, stride_c, ds, stride_s, batch_count);
        }
        gpu_time_used = get_time_us_sync(stream) - gpu_time_used;

        ArgumentModel<e_stride_a, e_stride_b, e_stride_c, e_stride_d, e_batch_count>{}.log_args<T>(
            rocblas_cout,
            arg,
            gpu_time_used,
            ArgumentLogging::NA_value,
            ArgumentLogging::NA_value,
            cpu_time_used,
            norm_error_host,
            norm_error_device);
    }
}