File: ctrmm_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 (238 lines) | stat: -rw-r--r-- 8,922 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
/*
    -- MAGMA (version 2.9.0) --
       Univ. of Tennessee, Knoxville
       Univ. of California, Berkeley
       Univ. of Colorado, Denver
       @date January 2025

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

       @author Ahmad Abdelfattah
*/

#include "magma_internal.h"
#include "commonblas_c.h"

#define PRECISION_c

/******************************************************************************/
extern "C" void
magmablas_ctrmm_vbatched_max_nocheck(
        magma_side_t side, magma_uplo_t uplo, magma_trans_t transA, magma_diag_t diag,
        magma_int_t max_m, magma_int_t max_n, magma_int_t* m, magma_int_t* n,
        magmaFloatComplex alpha,
        magmaFloatComplex **dA_array, magma_int_t* ldda,
        magmaFloatComplex **dB_array, magma_int_t* lddb,
        magma_int_t batchCount, magma_queue_t queue )
{
    if ( max_m <= 0 || max_n <= 0 )
        return;

    magmablas_ctrmm_vbatched_core(
            side, uplo, transA, diag,
            max_m, max_n, m, n,
            alpha, dA_array, 0, 0, ldda,
                   dB_array, 0, 0, lddb,
            batchCount, queue );
}

/******************************************************************************/
extern "C" void
magmablas_ctrmm_vbatched_max(
        magma_side_t side, magma_uplo_t uplo, magma_trans_t transA, magma_diag_t diag,
        magma_int_t max_m, magma_int_t max_n, magma_int_t* m, magma_int_t* n,
        magmaFloatComplex alpha,
        magmaFloatComplex **dA_array, magma_int_t* ldda,
        magmaFloatComplex **dB_array, magma_int_t* lddb,
        magma_int_t batchCount, magma_queue_t queue )
{
    magma_int_t info = 0;
    info =  magma_trmm_vbatched_checker(side, uplo, transA, diag, m, n, ldda, lddb, batchCount, queue);

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

    magmablas_ctrmm_vbatched_max_nocheck(
            side, uplo, transA, diag,
            max_m, max_n, m, n,
            alpha, dA_array, ldda,
                   dB_array, lddb,
            batchCount,
            queue );
}

/******************************************************************************/
extern "C" void
magmablas_ctrmm_vbatched_nocheck(
        magma_side_t side, magma_uplo_t uplo, magma_trans_t transA, magma_diag_t diag,
        magma_int_t* m, magma_int_t* n,
        magmaFloatComplex alpha,
        magmaFloatComplex **dA_array, magma_int_t* ldda,
        magmaFloatComplex **dB_array, magma_int_t* lddb,
        magma_int_t batchCount, magma_queue_t queue )
{
    // compute the max. dimensions
    magma_imax_size_2(m, n, batchCount, queue);
    magma_int_t max_m, max_n;
    magma_igetvector_async(1, &m[batchCount], 1, &max_m, 1, queue);
    magma_igetvector_async(1, &n[batchCount], 1, &max_n, 1, queue);
    magma_queue_sync( queue );

    magmablas_ctrmm_vbatched_max_nocheck(
            side, uplo, transA, diag,
            max_m, max_n, m, n,
            alpha, dA_array, ldda,
                   dB_array, lddb,
            batchCount, queue );
}

/***************************************************************************//**
    Purpose
    =======

    CTRMM  performs one of the matrix-matrix operations

       B := alpha*op( A )*B,   or   B := alpha*B*op( A )

    where  alpha  is a scalar,  B  is an m by n matrix,  A  is a unit, or
    non-unit,  upper or lower triangular matrix  and  op( A )  is one  of

       op( A ) = A   or   op( A ) = A'   or   op( A ) = conjg( A' ).

    Parameters
    ==========

    @param[in]
    side     magma_side_t.
             On entry,  side specifies whether  op( A ) multiplies B from
             the left or right as follows:
                side = magmaLeft   B := alpha*op( A )*B.
                side = magmaRight  B := alpha*B*op( A ).
             Unchanged on exit.

    @param[in]
    uplo     magma_uplo_t.
             On entry, uplo specifies whether the matrix A is an upper or
             lower triangular matrix as follows:
                uplo = magmaUpper   A is an upper triangular matrix.
                uplo = magmaLower   A is a lower triangular matrix.
             Unchanged on exit.

    @param[in]
    transA   magma_trans_t.
             On entry, transA specifies the form of op( A ) to be used in
             the matrix multiplication as follows:
                transA = MagmaNoTrans     op( A ) = A.
                transA = MagmaTrans       op( A ) = A'.
                transA = MagmaConjTrans   op( A ) = conjg( A' ).
             Unchanged on exit.

    @param[in]
    diag     magma_diag_t.
             On entry, diag specifies whether or not A is unit triangular
             as follows:
                diag = MagmaUnit      A is assumed to be unit triangular.
                diag = MagmaNonUnit   A is not assumed to be unit
                                    triangular.
             Unchanged on exit.

    @param[in]
    m        INTEGER array, dimension(batchCount + 1).
             On entry, each integer M specifies the number of rows of the corresponding
             matrix B. M must be at least zero.
             Unchanged on exit.

    @param[in]
    n        INTEGER array, dimension(batchCount + 1).
             On entry, each integer n specifies the number of columns of the corresponding
             matrix B.  N must be at least zero.
             Unchanged on exit.

    @param[in]
    alpha    DOUBLE COMPLEX.
             On entry,  alpha specifies the scalar  alpha. When  alpha is
             zero then  A is not referenced and  B need not be set before
             entry.
             Unchanged on exit.

    @param[in]
    dA_array     Array of pointers, dimension(batchCount).
             Each is a DOUBLE COMPLEX array A of DIMENSION ( ldda, k ), where k is M
             when  side = magmaLeft  and is  N  when  side = magmaRight.
             Before entry  with  uplo = magmaUpper,  the  leading  k by k
             upper triangular part of the array  A must contain the upper
             triangular matrix  and the strictly lower triangular part of
             A is not referenced.
             Before entry  with  uplo = magmaLower,  the  leading  k by k
             lower triangular part of the array  A must contain the lower
             triangular matrix  and the strictly upper triangular part of
             A is not referenced.
             Note that when  diag = MagmaUnit,  the diagonal elements of
             A  are not referenced either,  but are assumed to be  unity.
             Unchanged on exit.

    @param[in]
    ldda     INTEGER array, dimension(batchCount + 1).
             On entry, ldda specifies the first dimension of A as declared
             in the calling (sub) program.  When  side = magmaLeft  then
             ldda  must be at least  max( 1, M ),  when  side = magmaRight
             then ldda must be at least max( 1, N ).
             Unchanged on exit.

    @param[in,out]
    dB_array     Array of pointers, dimension(batchCount).
             Each is a DOUBLE COMPLEX array B of DIMENSION ( lddb, N ).
             Before entry,  the leading  M by N part of the array  B must
             contain the matrix  B,  and  on exit  is overwritten  by the
             transformed matrix.

    @param[in]
    lddb     INTEGER array, dimension(batchCount + 1).
             On entry, lddb specifies the first dimension of B as declared
             in  the  calling  (sub)  program.   lddb  must  be  at  least
             max( 1, M ).
             Unchanged on exit.

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

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

    @ingroup magma_trmm_batched
*******************************************************************************/
extern "C" void
magmablas_ctrmm_vbatched(
        magma_side_t side, magma_uplo_t uplo, magma_trans_t transA, magma_diag_t diag,
        magma_int_t* m, magma_int_t* n,
        magmaFloatComplex alpha,
        magmaFloatComplex **dA_array, magma_int_t* ldda,
        magmaFloatComplex **dB_array, magma_int_t* lddb,
        magma_int_t batchCount, magma_queue_t queue )
{
    magma_int_t info = 0;
    info =  magma_trmm_vbatched_checker(side, uplo, transA, diag, m, n, ldda, lddb, batchCount, queue);

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

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

    magmablas_ctrmm_vbatched_max_nocheck(
            side, uplo, transA, diag,
            max_m, max_n, m, n,
            alpha, dA_array, ldda,
                   dB_array, lddb,
            batchCount, queue );
}