File: cgemm_vbatched.cpp

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

       @generated from magmablas/zgemm_vbatched.cpp, normal z -> c, Wed Jan 22 14:42:07 2025

       @author Ahmad Abdelfattah
*/
#include "magma_internal.h"
#include "commonblas_c.h"

#define PRECISION_c

/******************************************************************************/
extern "C" void
magmablas_cgemm_vbatched_max_nocheck(
    magma_trans_t transA, magma_trans_t transB,
    magma_int_t* m, magma_int_t* n, magma_int_t* k,
    magmaFloatComplex alpha,
    magmaFloatComplex const * const * dA_array, magma_int_t* ldda,
    magmaFloatComplex const * const * dB_array, magma_int_t* lddb,
    magmaFloatComplex beta,
    magmaFloatComplex **dC_array, magma_int_t* lddc,
    magma_int_t batchCount,
    magma_int_t max_m, magma_int_t max_n, magma_int_t max_k,
    magma_queue_t queue )
{
    magmablas_cgemm_vbatched_core(
            transA, transB,
            max_m, max_n, max_k,
            m, n, k,
            alpha, dA_array, 0, 0, ldda,
                   dB_array, 0, 0, lddb,
            beta,  dC_array, 0, 0, lddc,
            batchCount, queue );
}


/******************************************************************************/
extern "C" void
magmablas_cgemm_vbatched_max(
    magma_trans_t transA, magma_trans_t transB,
    magma_int_t* m, magma_int_t* n, magma_int_t* k,
    magmaFloatComplex alpha,
    magmaFloatComplex const * const * dA_array, magma_int_t* ldda,
    magmaFloatComplex const * const * dB_array, magma_int_t* lddb,
    magmaFloatComplex beta,
    magmaFloatComplex **dC_array, magma_int_t* lddc,
    magma_int_t batchCount,
    magma_int_t max_m, magma_int_t max_n, magma_int_t max_k,
    magma_queue_t queue )
{
    magma_int_t info = 0;

    info =  magma_gemm_vbatched_checker( transA, transB, m, n, k, ldda, lddb, lddc, batchCount, queue );

    if (info != 0) {
        magma_xerbla( __func__, -(info) );
        return;
    }

    magmablas_cgemm_vbatched_max_nocheck(
            transA, transB,
            m, n, k,
            alpha, dA_array, ldda,
                   dB_array, lddb,
            beta,  dC_array, lddc,
            batchCount,
            max_m, max_n, max_k,
            queue );
}


/******************************************************************************/
extern "C" void
magmablas_cgemm_vbatched_nocheck(
    magma_trans_t transA, magma_trans_t transB,
    magma_int_t* m, magma_int_t* n, magma_int_t* k,
    magmaFloatComplex alpha,
    magmaFloatComplex const * const * dA_array, magma_int_t* ldda,
    magmaFloatComplex const * const * dB_array, magma_int_t* lddb,
    magmaFloatComplex beta,
    magmaFloatComplex **dC_array, magma_int_t* lddc,
    magma_int_t batchCount, magma_queue_t queue )
{
    // compute the max. dimensions
    magma_imax_size_3(m, n, k, batchCount, queue);
    magma_int_t max_m, max_n, max_k;
    magma_igetvector_async(1, &m[batchCount], 1, &max_m, 1, queue);
    magma_igetvector_async(1, &n[batchCount], 1, &max_n, 1, queue);
    magma_igetvector_async(1, &k[batchCount], 1, &max_k, 1, queue);
    magma_queue_sync( queue );

    magmablas_cgemm_vbatched_max_nocheck(
            transA, transB,
            m, n, k,
            alpha, dA_array, ldda,
                   dB_array, lddb,
            beta,  dC_array, lddc,
            batchCount,
            max_m, max_n, max_k,
            queue );
}


/***************************************************************************//**
    Purpose
    -------
    CGEMM performs one of the matrix-matrix operations

        C = alpha*op( A )*op( B ) + beta*C,

    where op( X ) is one of

        op( X ) = X      or
        op( X ) = X**T   or
        op( X ) = X**H,

    alpha and beta are scalars, and A, B and C are matrices, with
    op( A ) an m by k matrix, op( B ) a k by n matrix and C an m by n matrix.

    Parameters
    ----------
    @param[in]
    transA  magma_trans_t.
            On entry, transA specifies the form of op( A ) to be used in
            the matrix multiplication as follows:
      -     = MagmaNoTrans:    op( A ) = A.
      -     = MagmaTrans:      op( A ) = A**T.
      -     = MagmaConjTrans:  op( A ) = A**H.

    @param[in]
    transB  magma_trans_t.
            On entry, transB specifies the form of op( B ) to be used in
            the matrix multiplication as follows:
      -     = MagmaNoTrans:    op( B ) = B.
      -     = MagmaTrans:      op( B ) = B**T.
      -     = MagmaConjTrans:  op( B ) = B**H.

    @param[in]
    m       Array of integers, dimension (batchCount + 1)
            On entry,  each INTEGER M  specifies  the number  of rows  of the
            corresponding matrices op( A ) and C.  M  must  be at least  zero.
            The last element of the array is used internally by the routine.

    @param[in]
    n       Array of integers, dimension (batchCount + 1).
            On entry,  each INTEGER N  specifies the number  of columns of the
            corresponding matrix op( B ) and the number of columns of the
            corresponding matrix C. N must be at least zero.
            The last element of the array is used internally by the routine.

    @param[in]
    k       Array of integers, dimension (batchCount + 1).
            On entry,  each INTEGER K  specifies  the number of columns of the
            corresponding matrix op( A ) and the number of rows of the
            corresponding matrix op( B ). K must be at least  zero.
            The last element of the array is used internally by the routine.

    @param[in]
    alpha   COMPLEX
            On entry, ALPHA specifies the scalar alpha.

    @param[in]
    dA_array      Array of pointers, dimension (batchCount).
             Each is a COMPLEX array A of DIMENSION ( LDDA, ka ), where ka is
             K  when  transA = MagmaNoTrans,  and is  M  otherwise.
             Before entry with  transA = MagmaNoTrans,  the leading  M by K
             part of the array A must contain the matrix A, otherwise
             the leading  K by M  part of the array A must contain  the
             matrix A.

    @param[in]
    ldda    Array of integers, dimension (batchCount + 1).
            On entry, each INTEGER LDDA specifies the first dimension of
            the corresponding matrix A as declared in the calling (sub) program.
            When  transA = MagmaNoTrans then LDDA must be at least  max( 1, M ),
            otherwise  ldda must be at least  max( 1, K ).
            The last element of the array is used internally by the routine.

    @param[in]
    dB_array      Array of pointers, dimension (batchCount).
             Each is a COMPLEX array B of DIMENSION ( LDDB, kb ), where kb is
             N  when  transB = MagmaNoTrans,  and is  K  otherwise.
             Before entry with  transB = MagmaNoTrans,  the leading  K by N
             part of the array B must contain the matrix B, otherwise
             the leading  N by K  part of the array B must contain  the
             matrix B.

    @param[in]
    lddb    Array of integers, dimension (batchCount + 1).
            On entry, each INTEGER LDDB specifies the first dimension of
            the corresponding matrix B as declared in the calling (sub) program.
            When  transB = MagmaNoTrans then LDDB must be at least  max( 1, K ),
            otherwise  LDDB must be at least  max( 1, N ).
            The last element of the array is used internally by the routine.

    @param[in]
    beta    COMPLEX.
            On entry,  BETA  specifies the scalar  beta.  When  BETA  is
            supplied as zero then C need not be set on input.

    @param[in,out]
    dC_array      Array of pointers, dimension (batchCount).
             Each is a COMPLEX array C of DIMENSION ( LDDC, N ).
             Before entry, the leading  M by N  part of the array  C must
             contain the matrix  C,  except when  beta  is zero, in which
             case C need not be set on entry.
             On exit, the array  C  is overwritten by the  M by N  matrix
             ( alpha*op( A )*op( B ) + beta*C ).

    @param[in]
    lddc    Array of integers, dimension (batchCount + 1).
            On entry, each INTEGER LDDC specifies the first dimension of
            the corresponding matrix C as declared in  the  calling  (sub)  program.
            LDDC  must  be  at  least max( 1, M ).
            The last element of the array is used internally by the routine.

    @param[in]
    batchCount  INTEGER
                The number of matrices to operate on.

    @param[in]
    queue   magma_queue_t
            Queue to execute in.

    @ingroup magma_gemm_batched
*******************************************************************************/
extern "C" void
magmablas_cgemm_vbatched(
    magma_trans_t transA, magma_trans_t transB,
    magma_int_t* m, magma_int_t* n, magma_int_t* k,
    magmaFloatComplex alpha,
    magmaFloatComplex const * const * dA_array, magma_int_t* ldda,
    magmaFloatComplex const * const * dB_array, magma_int_t* lddb,
    magmaFloatComplex beta,
    magmaFloatComplex **dC_array, magma_int_t* lddc,
    magma_int_t batchCount, magma_queue_t queue )
{
    magma_int_t info = 0;

    info =  magma_gemm_vbatched_checker( transA, transB, m, n, k, ldda, lddb, lddc, batchCount, queue );

    if (info != 0) {
        magma_xerbla( __func__, -(info) );
        return;
    }

    // compute the max. dimensions
    magma_imax_size_3(m, n, k, batchCount, queue);
    magma_int_t max_m, max_n, max_k;
    magma_igetvector_async(1, &m[batchCount], 1, &max_m, 1, queue);
    magma_igetvector_async(1, &n[batchCount], 1, &max_n, 1, queue);
    magma_igetvector_async(1, &k[batchCount], 1, &max_k, 1, queue);
    magma_queue_sync( queue );

    magmablas_cgemm_vbatched_max_nocheck(
            transA, transB,
            m, n, k,
            alpha, dA_array, ldda,
                   dB_array, lddb,
            beta,  dC_array, lddc,
            batchCount,
            max_m, max_n, max_k,
            queue );
}