File: testing_ztrsm.cpp

package info (click to toggle)
magma 2.9.0%2Bds-3
  • links: PTS, VCS
  • area: contrib
  • in suites: forky, sid
  • size: 83,556 kB
  • sloc: cpp: 709,115; fortran: 121,916; ansic: 32,343; python: 25,603; f90: 15,208; makefile: 945; xml: 253; csh: 232; sh: 203; perl: 104
file content (279 lines) | stat: -rw-r--r-- 12,296 bytes parent folder | download | duplicates (3)
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
/*
    -- MAGMA (version 2.9.0) --
       Univ. of Tennessee, Knoxville
       Univ. of California, Berkeley
       Univ. of Colorado, Denver
       @date January 2025

       @precisions normal z -> c d s
       @author Chongxiao Cao
*/
// includes, system
#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include <math.h>

// includes, project
#include "flops.h"
#include "magma_v2.h"
#include "magma_lapack.h"
#include "magma_operators.h"  // for MAGMA_Z_DIV
#include "testings.h"

/* ////////////////////////////////////////////////////////////////////////////
   -- Testing ztrsm
*/
int main( int argc, char** argv)
{
    #ifdef MAGMA_HAVE_OPENCL
    #define dA(i_, j_)  dA, ((i_) + (j_)*ldda)
    #define dB(i_, j_)  dB, ((i_) + (j_)*lddb)
    #else
    #define dA(i_, j_) (dA + (i_) + (j_)*ldda)
    #define dB(i_, j_) (dB + (i_) + (j_)*lddb)
    #endif

    #define hA(i_, j_) (hA + (i_) + (j_)*lda)

    TESTING_CHECK( magma_init() );
    magma_print_environment();

    real_Double_t   gflops, magma_perf=0, magma_time=0, dev_perf, dev_time, cpu_perf=0, cpu_time=0;
    double          magma_error=0, dev_error, lapack_error, *work;
    magma_int_t M, N, info;
    magma_int_t Ak;
    magma_int_t sizeB;
    magma_int_t lda, ldb, ldda, lddb;
    magma_int_t ione     = 1;
    magma_int_t ISEED[4] = {0,0,0,1};

    magmaDoubleComplex *hA, *hB, *hBdev, *hBmagma, *hBlapack, *hX;
    magmaDoubleComplex_ptr dA, dB;
    magmaDoubleComplex c_neg_one = MAGMA_Z_NEG_ONE;
    magmaDoubleComplex c_one = MAGMA_Z_ONE;
    magmaDoubleComplex alpha = MAGMA_Z_ONE; //MAGMA_Z_MAKE(  0.29, -0.86 );
    int status = 0;

    magma_opts opts;
    opts.matrix = "rand_dominant";  // default
    opts.tolerance = 100;           // default
    opts.parse_opts( argc, argv );

    double tol = opts.tolerance * lapackf77_dlamch("E");

    // pass ngpu = -1 to test multi-GPU code using 1 gpu
    magma_int_t abs_ngpu = abs( opts.ngpu );

    printf("%% side = %s, uplo = %s, transA = %s, diag = %s, ngpu = %lld\n",
           lapack_side_const(opts.side), lapack_uplo_const(opts.uplo),
           lapack_trans_const(opts.transA), lapack_diag_const(opts.diag), (long long) abs_ngpu);

    printf("%%   M     N  MAGMA Gflop/s (ms)  %s Gflop/s (ms)   CPU Gflop/s (ms)      MAGMA     %s   LAPACK error\n", g_platform_str, g_platform_str);
    printf("%%============================================================================================================\n");
    for( int itest = 0; itest < opts.ntest; ++itest ) {
        for( int iter = 0; iter < opts.niter; ++iter ) {
            M = opts.msize[itest];
            N = opts.nsize[itest];
            gflops = FLOPS_ZTRSM(opts.side, M, N) / 1e9;

            if ( opts.side == MagmaLeft ) {
                lda = M;
                Ak  = M;
            } else {
                lda = N;
                Ak  = N;
            }

            ldb = M;

            ldda = magma_roundup( lda, opts.align );  // multiple of 32 by default
            lddb = magma_roundup( ldb, opts.align );  // multiple of 32 by default

            //sizeA = lda*Ak;
            sizeB = ldb*N;

            TESTING_CHECK( magma_zmalloc_cpu( &hA,       lda*Ak    ));
            TESTING_CHECK( magma_zmalloc_cpu( &hB,       ldb*N     ));
            TESTING_CHECK( magma_zmalloc_cpu( &hX,       ldb*N     ));
            TESTING_CHECK( magma_zmalloc_cpu( &hBlapack, ldb*N     ));
            TESTING_CHECK( magma_zmalloc_cpu( &hBdev,    ldb*N     ));
            TESTING_CHECK( magma_zmalloc_cpu( &hBmagma,  ldb*N     ));
            TESTING_CHECK( magma_dmalloc_cpu( &work,     max(M, N) ));

            TESTING_CHECK( magma_zmalloc( &dA,       ldda*Ak ));
            TESTING_CHECK( magma_zmalloc( &dB,       lddb*N  ));

            /* Initialize the matrices */
            magma_generate_matrix( opts, Ak, Ak, hA, lda );

            // set unused data to nan
            magma_int_t Ak_1 = Ak - 1;
            if (opts.uplo == MagmaLower)
                lapackf77_zlaset( "upper", &Ak_1, &Ak_1, &MAGMA_Z_NAN, &MAGMA_Z_NAN, &hA[ 1*lda ], &lda );
            else
                lapackf77_zlaset( "lower", &Ak_1, &Ak_1, &MAGMA_Z_NAN, &MAGMA_Z_NAN, &hA[ 1     ], &lda );

            // Factor A into L L^H or U U^H to get a well-conditioned triangular matrix.
            // If diag == Unit, the diagonal is replaced; this is still well-conditioned.
            // First, brute force positive definiteness.
            for (int i = 0; i < Ak; ++i) {
                hA[ i + i*lda ] += Ak;
            }
            lapackf77_zpotrf( lapack_uplo_const(opts.uplo), &Ak, hA, &lda, &info );
            assert( info == 0 );

            lapackf77_zlarnv( &ione, ISEED, &sizeB, hB );
            lapackf77_zlacpy( MagmaFullStr, &M, &N, hB, &ldb, hBlapack, &ldb );
            lapackf77_zlacpy( MagmaFullStr, &M, &N, hB, &ldb, hBmagma,  &ldb );
            magma_zsetmatrix( Ak, Ak, hA, lda, dA(0,0), ldda, opts.queue );

            /* =====================================================================
               Performs operation using MAGMABLAS (only with CUDA)
               =================================================================== */
            #if defined(MAGMA_HAVE_CUDA) || defined(MAGMA_HAVE_HIP)
                magma_zsetmatrix( M, N, hB, ldb, dB(0,0), lddb, opts.queue );

                if (opts.ngpu == 1) {
                    magma_time = magma_sync_wtime( opts.queue );
                    magmablas_ztrsm( opts.side, opts.uplo, opts.transA, opts.diag,
                                     M, N,
                                     alpha, dA(0,0), ldda,
                                            dB(0,0), lddb, opts.queue );
                    magma_time = magma_sync_wtime( opts.queue ) - magma_time;
                    magma_zgetmatrix( M, N, dB(0,0), lddb, hBmagma, ldb, opts.queue );
                }
                else {
                    magma_time = magma_wtime();
                    magma_ztrsm_m( abs_ngpu, opts.side, opts.uplo, opts.transA, opts.diag,
                                   M, N,
                                   alpha, hA,      lda,
                                          hBmagma, ldb );
                    magma_time = magma_wtime() - magma_time;
                }
                magma_perf = gflops / magma_time;
            #endif

            /* =====================================================================
               Performs operation using cuBLAS / clBLAS
               =================================================================== */
            magma_zsetmatrix( M, N, hB, ldb, dB(0,0), lddb, opts.queue );

            dev_time = magma_sync_wtime( opts.queue );
            magma_ztrsm( opts.side, opts.uplo, opts.transA, opts.diag,
                         M, N,
                         alpha, dA(0,0), ldda,
                                dB(0,0), lddb, opts.queue );
            dev_time = magma_sync_wtime( opts.queue ) - dev_time;
            dev_perf = gflops / dev_time;

            magma_zgetmatrix( M, N, dB(0,0), lddb, hBdev, ldb, opts.queue );

            /* =====================================================================
               Performs operation using CPU BLAS
               =================================================================== */
            if ( opts.lapack ) {
                cpu_time = magma_wtime();
                blasf77_ztrsm( lapack_side_const(opts.side), lapack_uplo_const(opts.uplo),
                               lapack_trans_const(opts.transA), lapack_diag_const(opts.diag),
                               &M, &N,
                               &alpha, hA, &lda,
                                       hBlapack, &ldb );
                cpu_time = magma_wtime() - cpu_time;
                cpu_perf = gflops / cpu_time;
            }

            /* =====================================================================
               Check the result
               =================================================================== */
            // ||b - 1/alpha*A*x|| / (||A||*||x||)
            magmaDoubleComplex inv_alpha = MAGMA_Z_DIV( c_one, alpha );
            double normR, normX, normA;
            normA = lapackf77_zlantr( "I",
                                      lapack_uplo_const(opts.uplo),
                                      lapack_diag_const(opts.diag),
                                      &Ak, &Ak, hA, &lda, work );

            #if defined(MAGMA_HAVE_CUDA) || defined(MAGMA_HAVE_HIP)
                // check magma
                memcpy( hX, hBmagma, sizeB*sizeof(magmaDoubleComplex) );
                blasf77_ztrmm( lapack_side_const(opts.side), lapack_uplo_const(opts.uplo),
                               lapack_trans_const(opts.transA), lapack_diag_const(opts.diag),
                               &M, &N,
                               &inv_alpha, hA, &lda,
                                           hX, &ldb );

                blasf77_zaxpy( &sizeB, &c_neg_one, hB, &ione, hX, &ione );
                normR = lapackf77_zlange( "I", &M, &N, hX,      &ldb, work );
                normX = lapackf77_zlange( "I", &M, &N, hBmagma, &ldb, work );
                magma_error = normR/(normX*normA);
            #endif

            // check cuBLAS / clBLAS
            memcpy( hX, hBdev, sizeB*sizeof(magmaDoubleComplex) );
            blasf77_ztrmm( lapack_side_const(opts.side), lapack_uplo_const(opts.uplo),
                           lapack_trans_const(opts.transA), lapack_diag_const(opts.diag),
                           &M, &N,
                           &inv_alpha, hA, &lda,
                                       hX, &ldb );

            blasf77_zaxpy( &sizeB, &c_neg_one, hB, &ione, hX, &ione );
            normR = lapackf77_zlange( "I", &M, &N, hX,    &ldb, work );
            normX = lapackf77_zlange( "I", &M, &N, hBdev, &ldb, work );
            dev_error = normR/(normX*normA);

            bool okay = (magma_error < tol && dev_error < tol);
            status += ! okay;
            if ( opts.lapack ) {
                // check lapack
                // this verifies that the matrix wasn't so bad that it couldn't be solved accurately.
                memcpy( hX, hBlapack, sizeB*sizeof(magmaDoubleComplex) );
                blasf77_ztrmm( lapack_side_const(opts.side), lapack_uplo_const(opts.uplo),
                               lapack_trans_const(opts.transA), lapack_diag_const(opts.diag),
                               &M, &N,
                               &inv_alpha, hA, &lda,
                                           hX, &ldb );

                blasf77_zaxpy( &sizeB, &c_neg_one, hB, &ione, hX, &ione );
                normR = lapackf77_zlange( "I", &M, &N, hX,       &ldb, work );
                normX = lapackf77_zlange( "I", &M, &N, hBlapack, &ldb, work );
                lapack_error = normR/(normX*normA);

                printf("%5lld %5lld   %7.2f (%7.2f)   %7.2f (%7.2f)   %7.2f (%7.2f)   %8.2e   %8.2e   %8.2e   %s\n",
                        (long long) M, (long long) N,
                        magma_perf,  1000.*magma_time,
                        dev_perf,    1000.*dev_time,
                        cpu_perf,    1000.*cpu_time,
                        magma_error, dev_error, lapack_error,
                        (okay ? "ok" : "failed"));
            }
            else {
                printf("%5lld %5lld   %7.2f (%7.2f)   %7.2f (%7.2f)     ---   (  ---  )   %8.2e   %8.2e     ---      %s\n",
                        (long long) M, (long long) N,
                        magma_perf,  1000.*magma_time,
                        dev_perf,    1000.*dev_time,
                        magma_error, dev_error,
                        (okay ? "ok" : "failed"));
            }

            magma_free_cpu( hA );
            magma_free_cpu( hB );
            magma_free_cpu( hX );
            magma_free_cpu( hBlapack );
            magma_free_cpu( hBdev    );
            magma_free_cpu( hBmagma  );
            magma_free_cpu( work     );

            magma_free( dA );
            magma_free( dB );
            fflush( stdout );
        }
        if ( opts.niter > 1 ) {
            printf( "\n" );
        }
    }

    opts.cleanup();
    TESTING_CHECK( magma_finalize() );
    return status;
}