libflame revision_anchor
Functions
bl1_trsm.c File Reference

(r)

Functions

void bl1_strsm (side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, float *alpha, float *a, int a_rs, int a_cs, float *b, int b_rs, int b_cs)
 
void bl1_dtrsm (side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, double *alpha, double *a, int a_rs, int a_cs, double *b, int b_rs, int b_cs)
 
void bl1_ctrsm (side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, scomplex *alpha, scomplex *a, int a_rs, int a_cs, scomplex *b, int b_rs, int b_cs)
 
void bl1_ztrsm (side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, dcomplex *alpha, dcomplex *a, int a_rs, int a_cs, dcomplex *b, int b_rs, int b_cs)
 
void bl1_strsm_blas (side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, float *alpha, float *a, int lda, float *b, int ldb)
 
void bl1_dtrsm_blas (side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, double *alpha, double *a, int lda, double *b, int ldb)
 
void bl1_ctrsm_blas (side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, scomplex *alpha, scomplex *a, int lda, scomplex *b, int ldb)
 
void bl1_ztrsm_blas (side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, dcomplex *alpha, dcomplex *a, int lda, dcomplex *b, int ldb)
 

Function Documentation

◆ bl1_ctrsm()

void bl1_ctrsm ( side1_t  side,
uplo1_t  uplo,
trans1_t  trans,
diag1_t  diag,
int  m,
int  n,
scomplex alpha,
scomplex a,
int  a_rs,
int  a_cs,
scomplex b,
int  b_rs,
int  b_cs 
)
220{
221 int m_save = m;
222 int n_save = n;
223 scomplex* a_save = a;
224 scomplex* b_save = b;
225 int a_rs_save = a_rs;
226 int a_cs_save = a_cs;
227 int b_rs_save = b_rs;
228 int b_cs_save = b_cs;
230 int dim_a;
231 int lda, inca;
232 int ldb, incb;
233 int lda_conj, inca_conj;
234 int a_was_copied;
235
236 // Return early if possible.
237 if ( bl1_zero_dim2( m, n ) ) return;
238
239 // If necessary, allocate, initialize, and use a temporary contiguous
240 // copy of each matrix rather than the original matrices.
243 dim_a,
244 dim_a,
246 &a, &a_rs, &a_cs );
247
249 n,
251 &b, &b_rs, &b_cs );
252
253 // Figure out whether A was copied to contiguous memory. This is used to
254 // prevent redundant copying.
255 a_was_copied = ( a != a_save );
256
257 // Initialize with values assuming column-major storage.
258 lda = a_cs;
259 inca = a_rs;
260 ldb = b_cs;
261 incb = b_rs;
262
263 // Adjust the parameters based on the storage of each matrix.
264 if ( bl1_is_col_storage( b_rs, b_cs ) )
265 {
266 if ( bl1_is_col_storage( a_rs, a_cs ) )
267 {
268 // requested operation: B_c := tr( uplo( A_c ) ) * B_c
269 // effective operation: B_c := tr( uplo( A_c ) ) * B_c
270 }
271 else // if ( bl1_is_row_storage( a_rs, a_cs ) )
272 {
273 // requested operation: B_c := tr( uplo( A_r ) ) * B_c
274 // effective operation: B_c := tr( ~uplo( A_c ) )^T * B_c
276
277 bl1_toggle_uplo( uplo );
279 }
280 }
281 else // if ( bl1_is_row_storage( b_rs, b_cs ) )
282 {
283 if ( bl1_is_col_storage( a_rs, a_cs ) )
284 {
285 // requested operation: B_r := tr( uplo( A_c ) ) * B_r
286 // effective operation: B_c := B_c * tr( uplo( A_c ) )^T
288
289 bl1_swap_ints( m, n );
290
293 }
294 else // if ( bl1_is_row_storage( a_rs, a_cs ) )
295 {
296 // requested operation: B_r := tr( uplo( A_r ) ) * B_r
297 // effective operation: B_c := B_c * tr( ~uplo( A_c ) )
300
301 bl1_swap_ints( m, n );
302
303 bl1_toggle_uplo( uplo );
305 }
306 }
307
308 // Initialize with values assuming that trans is not conjnotrans.
309 a_conj = a;
310 lda_conj = lda;
311 inca_conj = inca;
312
313 // We want to handle the conjnotrans case. The easiest way to do so is
314 // by making a conjugated copy of A.
316 {
317 int dim_a;
318
320
322 lda_conj = dim_a;
323 inca_conj = 1;
324
325 bl1_ccopymrt( uplo,
327 dim_a,
328 dim_a,
329 a, inca, lda,
331 }
332 else if ( bl1_is_conjnotrans( trans ) && a_was_copied )
333 {
334 int dim_a;
335
337
338 bl1_cconjmr( uplo,
339 dim_a,
340 dim_a,
342 }
343
344
346 uplo,
347 trans,
348 diag,
349 m,
350 n,
351 alpha,
353 b, ldb );
354
356 bl1_cfree( a_conj );
357
358 // Free any temporary contiguous matrices, copying the result back to
359 // the original matrix.
361 &a, &a_rs, &a_cs );
362
364 n_save,
366 &b, &b_rs, &b_cs );
367}
int i
Definition bl1_axmyv2.c:145
void bl1_cconjmr(uplo1_t uplo, int m, int n, scomplex *a, int a_rs, int a_cs)
Definition bl1_conjmr.c:23
void bl1_ccopymrt(uplo1_t uplo, trans1_t trans, int m, int n, scomplex *a, int a_rs, int a_cs, scomplex *b, int b_rs, int b_cs)
Definition bl1_copymrt.c:223
void bl1_ctrsm_blas(side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, scomplex *alpha, scomplex *a, int lda, scomplex *b, int ldb)
Definition bl1_trsm.c:614
int bl1_is_col_storage(int rs, int cs)
Definition bl1_is.c:90
int bl1_is_conjnotrans(trans1_t trans)
Definition bl1_is.c:25
int bl1_zero_dim2(int m, int n)
Definition bl1_is.c:118
void bl1_cfree_contigm(scomplex *a_save, int a_rs_save, int a_cs_save, scomplex **a, int *a_rs, int *a_cs)
Definition bl1_free_contigm.c:45
void bl1_cfree(scomplex *p)
Definition bl1_free.c:40
void bl1_ccreate_contigm(int m, int n, scomplex *a_save, int a_rs_save, int a_cs_save, scomplex **a, int *a_rs, int *a_cs)
Definition bl1_create_contigm.c:81
scomplex * bl1_callocm(unsigned int m, unsigned int n)
Definition bl1_allocm.c:40
void bl1_ccreate_contigmr(uplo1_t uplo, int m, int n, scomplex *a_save, int a_rs_save, int a_cs_save, scomplex **a, int *a_rs, int *a_cs)
Definition bl1_create_contigmr.c:77
void bl1_cfree_saved_contigm(int m, int n, scomplex *a_save, int a_rs_save, int a_cs_save, scomplex **a, int *a_rs, int *a_cs)
Definition bl1_free_saved_contigm.c:59
void bl1_set_dim_with_side(side1_t side, int m, int n, int *dim_new)
Definition bl1_set_dims.c:27
@ BLIS1_CONJ_NO_TRANSPOSE
Definition blis_type_defs.h:56
Definition blis_type_defs.h:133

References bl1_callocm(), bl1_cconjmr(), bl1_ccopymrt(), bl1_ccreate_contigm(), bl1_ccreate_contigmr(), bl1_cfree(), bl1_cfree_contigm(), bl1_cfree_saved_contigm(), bl1_ctrsm_blas(), bl1_is_col_storage(), bl1_is_conjnotrans(), bl1_set_dim_with_side(), bl1_zero_dim2(), and BLIS1_CONJ_NO_TRANSPOSE.

Referenced by bl1_ctrsmsx(), FLA_LU_nopiv_opc_var1(), FLA_LU_nopiv_opc_var2(), FLA_LU_nopiv_opc_var3(), FLA_LU_piv_opc_var3(), and FLA_Trsm_external().

◆ bl1_ctrsm_blas()

void bl1_ctrsm_blas ( side1_t  side,
uplo1_t  uplo,
trans1_t  trans,
diag1_t  diag,
int  m,
int  n,
scomplex alpha,
scomplex a,
int  lda,
scomplex b,
int  ldb 
)
615{
616#ifdef BLIS1_ENABLE_CBLAS_INTERFACES
622
627
633 m,
634 n,
635 alpha,
636 a, lda,
637 b, ldb );
638#else
639 char blas_side;
640 char blas_uplo;
641 char blas_trans;
642 char blas_diag;
643
648
650 &blas_uplo,
651 &blas_trans,
652 &blas_diag,
653 &m,
654 &n,
655 alpha,
656 a, &lda,
657 b, &ldb );
658#endif
659}
void F77_ctrsm(char *side, char *uplo, char *transa, char *diag, int *m, int *n, scomplex *alpha, scomplex *a, int *lda, scomplex *b, int *ldb)
CBLAS_ORDER
Definition blis_prototypes_cblas.h:17
@ CblasColMajor
Definition blis_prototypes_cblas.h:17
CBLAS_UPLO
Definition blis_prototypes_cblas.h:19
CBLAS_TRANSPOSE
Definition blis_prototypes_cblas.h:18
CBLAS_SIDE
Definition blis_prototypes_cblas.h:21
CBLAS_DIAG
Definition blis_prototypes_cblas.h:20
void cblas_ctrsm(const enum CBLAS_ORDER Order, const enum CBLAS_SIDE Side, const enum CBLAS_UPLO Uplo, const enum CBLAS_TRANSPOSE TransA, const enum CBLAS_DIAG Diag, const int M, const int N, const void *alpha, const void *A, const int lda, void *B, const int ldb)
void bl1_param_map_to_netlib_diag(diag1_t blis_diag, void *blas_diag)
Definition bl1_param_map.c:95
void bl1_param_map_to_netlib_trans(trans1_t blis_trans, void *blas_trans)
Definition bl1_param_map.c:15
void bl1_param_map_to_netlib_side(side1_t blis_side, void *blas_side)
Definition bl1_param_map.c:71
void bl1_param_map_to_netlib_uplo(uplo1_t blis_uplo, void *blas_uplo)
Definition bl1_param_map.c:47

References bl1_param_map_to_netlib_diag(), bl1_param_map_to_netlib_side(), bl1_param_map_to_netlib_trans(), bl1_param_map_to_netlib_uplo(), cblas_ctrsm(), CblasColMajor, and F77_ctrsm().

Referenced by bl1_ctrsm().

◆ bl1_dtrsm()

void bl1_dtrsm ( side1_t  side,
uplo1_t  uplo,
trans1_t  trans,
diag1_t  diag,
int  m,
int  n,
double alpha,
double a,
int  a_rs,
int  a_cs,
double b,
int  b_rs,
int  b_cs 
)
117{
118 int m_save = m;
119 int n_save = n;
120 double* a_save = a;
121 double* b_save = b;
122 int a_rs_save = a_rs;
123 int a_cs_save = a_cs;
124 int b_rs_save = b_rs;
125 int b_cs_save = b_cs;
126 int dim_a;
127 int lda, inca;
128 int ldb, incb;
129
130 // Return early if possible.
131 if ( bl1_zero_dim2( m, n ) ) return;
132
133 // If necessary, allocate, initialize, and use a temporary contiguous
134 // copy of each matrix rather than the original matrices.
137 dim_a,
138 dim_a,
140 &a, &a_rs, &a_cs );
141
143 n,
145 &b, &b_rs, &b_cs );
146
147 // Initialize with values assuming column-major storage.
148 lda = a_cs;
149 inca = a_rs;
150 ldb = b_cs;
151 incb = b_rs;
152
153 // Adjust the parameters based on the storage of each matrix.
154 if ( bl1_is_col_storage( b_rs, b_cs ) )
155 {
156 if ( bl1_is_col_storage( a_rs, a_cs ) )
157 {
158 // requested operation: B_c := tr( uplo( A_c ) ) * B_c
159 // effective operation: B_c := tr( uplo( A_c ) ) * B_c
160 }
161 else // if ( bl1_is_row_storage( a_rs, a_cs ) )
162 {
163 // requested operation: B_c := tr( uplo( A_r ) ) * B_c
164 // effective operation: B_c := tr( ~uplo( A_c ) )^T * B_c
166
167 bl1_toggle_uplo( uplo );
169 }
170 }
171 else // if ( bl1_is_row_storage( b_rs, b_cs ) )
172 {
173 if ( bl1_is_col_storage( a_rs, a_cs ) )
174 {
175 // requested operation: B_r := tr( uplo( A_c ) ) * B_r
176 // effective operation: B_c := B_c * tr( uplo( A_c ) )^T
178
179 bl1_swap_ints( m, n );
180
183 }
184 else // if ( bl1_is_row_storage( a_rs, a_cs ) )
185 {
186 // requested operation: B_r := tr( uplo( A_r ) ) * B_r
187 // effective operation: B_c := B_c * tr( ~uplo( A_c ) )
190
191 bl1_swap_ints( m, n );
192
193 bl1_toggle_uplo( uplo );
195 }
196 }
197
199 uplo,
200 trans,
201 diag,
202 m,
203 n,
204 alpha,
205 a, lda,
206 b, ldb );
207
208 // Free any temporary contiguous matrices, copying the result back to
209 // the original matrix.
211 &a, &a_rs, &a_cs );
212
214 n_save,
216 &b, &b_rs, &b_cs );
217}
void bl1_dtrsm_blas(side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, double *alpha, double *a, int lda, double *b, int ldb)
Definition bl1_trsm.c:567
void bl1_dfree_contigm(double *a_save, int a_rs_save, int a_cs_save, double **a, int *a_rs, int *a_cs)
Definition bl1_free_contigm.c:29
void bl1_dcreate_contigmr(uplo1_t uplo, int m, int n, double *a_save, int a_rs_save, int a_cs_save, double **a, int *a_rs, int *a_cs)
Definition bl1_create_contigmr.c:45
void bl1_dcreate_contigm(int m, int n, double *a_save, int a_rs_save, int a_cs_save, double **a, int *a_rs, int *a_cs)
Definition bl1_create_contigm.c:47
void bl1_dfree_saved_contigm(int m, int n, double *a_save, int a_rs_save, int a_cs_save, double **a, int *a_rs, int *a_cs)
Definition bl1_free_saved_contigm.c:36

References bl1_dcreate_contigm(), bl1_dcreate_contigmr(), bl1_dfree_contigm(), bl1_dfree_saved_contigm(), bl1_dtrsm_blas(), bl1_is_col_storage(), bl1_set_dim_with_side(), and bl1_zero_dim2().

Referenced by bl1_dtrsmsx(), FLA_LU_nopiv_opd_var1(), FLA_LU_nopiv_opd_var2(), FLA_LU_nopiv_opd_var3(), FLA_LU_piv_opd_var3(), and FLA_Trsm_external().

◆ bl1_dtrsm_blas()

void bl1_dtrsm_blas ( side1_t  side,
uplo1_t  uplo,
trans1_t  trans,
diag1_t  diag,
int  m,
int  n,
double alpha,
double a,
int  lda,
double b,
int  ldb 
)
568{
569#ifdef BLIS1_ENABLE_CBLAS_INTERFACES
575
580
586 m,
587 n,
588 *alpha,
589 a, lda,
590 b, ldb );
591#else
592 char blas_side;
593 char blas_uplo;
594 char blas_trans;
595 char blas_diag;
596
601
603 &blas_uplo,
604 &blas_trans,
605 &blas_diag,
606 &m,
607 &n,
608 alpha,
609 a, &lda,
610 b, &ldb );
611#endif
612}
void F77_dtrsm(char *side, char *uplo, char *transa, char *diag, int *m, int *n, double *alpha, double *a, int *lda, double *b, int *ldb)
void cblas_dtrsm(const enum CBLAS_ORDER Order, const enum CBLAS_SIDE Side, const enum CBLAS_UPLO Uplo, const enum CBLAS_TRANSPOSE TransA, const enum CBLAS_DIAG Diag, const int M, const int N, const double alpha, const double *A, const int lda, double *B, const int ldb)

References bl1_param_map_to_netlib_diag(), bl1_param_map_to_netlib_side(), bl1_param_map_to_netlib_trans(), bl1_param_map_to_netlib_uplo(), cblas_dtrsm(), CblasColMajor, and F77_dtrsm().

Referenced by bl1_dtrsm().

◆ bl1_strsm()

void bl1_strsm ( side1_t  side,
uplo1_t  uplo,
trans1_t  trans,
diag1_t  diag,
int  m,
int  n,
float alpha,
float a,
int  a_rs,
int  a_cs,
float b,
int  b_rs,
int  b_cs 
)
14{
15 int m_save = m;
16 int n_save = n;
17 float* a_save = a;
18 float* b_save = b;
19 int a_rs_save = a_rs;
20 int a_cs_save = a_cs;
21 int b_rs_save = b_rs;
22 int b_cs_save = b_cs;
23 int dim_a;
24 int lda, inca;
25 int ldb, incb;
26
27 // Return early if possible.
28 if ( bl1_zero_dim2( m, n ) ) return;
29
30 // If necessary, allocate, initialize, and use a temporary contiguous
31 // copy of each matrix rather than the original matrices.
34 dim_a,
35 dim_a,
37 &a, &a_rs, &a_cs );
38
40 n,
42 &b, &b_rs, &b_cs );
43
44 // Initialize with values assuming column-major storage.
45 lda = a_cs;
46 inca = a_rs;
47 ldb = b_cs;
48 incb = b_rs;
49
50 // Adjust the parameters based on the storage of each matrix.
52 {
54 {
55 // requested operation: B_c := tr( uplo( A_c ) ) * B_c
56 // effective operation: B_c := tr( uplo( A_c ) ) * B_c
57 }
58 else // if ( bl1_is_row_storage( a_rs, a_cs ) )
59 {
60 // requested operation: B_c := tr( uplo( A_r ) ) * B_c
61 // effective operation: B_c := tr( ~uplo( A_c ) )^T * B_c
63
64 bl1_toggle_uplo( uplo );
66 }
67 }
68 else // if ( bl1_is_row_storage( b_rs, b_cs ) )
69 {
71 {
72 // requested operation: B_r := tr( uplo( A_c ) ) * B_r
73 // effective operation: B_c := B_c * tr( uplo( A_c ) )^T
75
76 bl1_swap_ints( m, n );
77
80 }
81 else // if ( bl1_is_row_storage( a_rs, a_cs ) )
82 {
83 // requested operation: B_r := tr( uplo( A_r ) ) * B_r
84 // effective operation: B_c := B_c * tr( ~uplo( A_c ) )
87
88 bl1_swap_ints( m, n );
89
90 bl1_toggle_uplo( uplo );
92 }
93 }
94
96 uplo,
97 trans,
98 diag,
99 m,
100 n,
101 alpha,
102 a, lda,
103 b, ldb );
104
105 // Free any temporary contiguous matrices, copying the result back to
106 // the original matrix.
108 &a, &a_rs, &a_cs );
109
111 n_save,
113 &b, &b_rs, &b_cs );
114}
void bl1_strsm_blas(side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, float *alpha, float *a, int lda, float *b, int ldb)
Definition bl1_trsm.c:520
void bl1_sfree_contigm(float *a_save, int a_rs_save, int a_cs_save, float **a, int *a_rs, int *a_cs)
Definition bl1_free_contigm.c:13
void bl1_sfree_saved_contigm(int m, int n, float *a_save, int a_rs_save, int a_cs_save, float **a, int *a_rs, int *a_cs)
Definition bl1_free_saved_contigm.c:13
void bl1_screate_contigmr(uplo1_t uplo, int m, int n, float *a_save, int a_rs_save, int a_cs_save, float **a, int *a_rs, int *a_cs)
Definition bl1_create_contigmr.c:13
void bl1_screate_contigm(int m, int n, float *a_save, int a_rs_save, int a_cs_save, float **a, int *a_rs, int *a_cs)
Definition bl1_create_contigm.c:13

References bl1_is_col_storage(), bl1_screate_contigm(), bl1_screate_contigmr(), bl1_set_dim_with_side(), bl1_sfree_contigm(), bl1_sfree_saved_contigm(), bl1_strsm_blas(), and bl1_zero_dim2().

Referenced by bl1_strsmsx(), FLA_LU_nopiv_ops_var1(), FLA_LU_nopiv_ops_var2(), FLA_LU_nopiv_ops_var3(), FLA_LU_piv_ops_var3(), and FLA_Trsm_external().

◆ bl1_strsm_blas()

void bl1_strsm_blas ( side1_t  side,
uplo1_t  uplo,
trans1_t  trans,
diag1_t  diag,
int  m,
int  n,
float alpha,
float a,
int  lda,
float b,
int  ldb 
)
521{
522#ifdef BLIS1_ENABLE_CBLAS_INTERFACES
528
533
539 m,
540 n,
541 *alpha,
542 a, lda,
543 b, ldb );
544#else
545 char blas_side;
546 char blas_uplo;
547 char blas_trans;
548 char blas_diag;
549
554
556 &blas_uplo,
557 &blas_trans,
558 &blas_diag,
559 &m,
560 &n,
561 alpha,
562 a, &lda,
563 b, &ldb );
564#endif
565}
void F77_strsm(char *side, char *uplo, char *transa, char *diag, int *m, int *n, float *alpha, float *a, int *lda, float *b, int *ldb)
void cblas_strsm(const enum CBLAS_ORDER Order, const enum CBLAS_SIDE Side, const enum CBLAS_UPLO Uplo, const enum CBLAS_TRANSPOSE TransA, const enum CBLAS_DIAG Diag, const int M, const int N, const float alpha, const float *A, const int lda, float *B, const int ldb)

References bl1_param_map_to_netlib_diag(), bl1_param_map_to_netlib_side(), bl1_param_map_to_netlib_trans(), bl1_param_map_to_netlib_uplo(), cblas_strsm(), CblasColMajor, and F77_strsm().

Referenced by bl1_strsm().

◆ bl1_ztrsm()

void bl1_ztrsm ( side1_t  side,
uplo1_t  uplo,
trans1_t  trans,
diag1_t  diag,
int  m,
int  n,
dcomplex alpha,
dcomplex a,
int  a_rs,
int  a_cs,
dcomplex b,
int  b_rs,
int  b_cs 
)
370{
371 int m_save = m;
372 int n_save = n;
373 dcomplex* a_save = a;
374 dcomplex* b_save = b;
375 int a_rs_save = a_rs;
376 int a_cs_save = a_cs;
377 int b_rs_save = b_rs;
378 int b_cs_save = b_cs;
380 int dim_a;
381 int lda, inca;
382 int ldb, incb;
383 int lda_conj, inca_conj;
384 int a_was_copied;
385
386 // Return early if possible.
387 if ( bl1_zero_dim2( m, n ) ) return;
388
389 // If necessary, allocate, initialize, and use a temporary contiguous
390 // copy of each matrix rather than the original matrices.
393 dim_a,
394 dim_a,
396 &a, &a_rs, &a_cs );
397
399 n,
401 &b, &b_rs, &b_cs );
402
403 // Figure out whether A was copied to contiguous memory. This is used to
404 // prevent redundant copying.
405 a_was_copied = ( a != a_save );
406
407 // Initialize with values assuming column-major storage.
408 lda = a_cs;
409 inca = a_rs;
410 ldb = b_cs;
411 incb = b_rs;
412
413 // Adjust the parameters based on the storage of each matrix.
414 if ( bl1_is_col_storage( b_rs, b_cs ) )
415 {
416 if ( bl1_is_col_storage( a_rs, a_cs ) )
417 {
418 // requested operation: B_c := tr( uplo( A_c ) ) * B_c
419 // effective operation: B_c := tr( uplo( A_c ) ) * B_c
420 }
421 else // if ( bl1_is_row_storage( a_rs, a_cs ) )
422 {
423 // requested operation: B_c := tr( uplo( A_r ) ) * B_c
424 // effective operation: B_c := tr( ~uplo( A_c ) )^T * B_c
426
427 bl1_toggle_uplo( uplo );
429 }
430 }
431 else // if ( bl1_is_row_storage( b_rs, b_cs ) )
432 {
433 if ( bl1_is_col_storage( a_rs, a_cs ) )
434 {
435 // requested operation: B_r := tr( uplo( A_c ) ) * B_r
436 // effective operation: B_c := B_c * tr( uplo( A_c ) )^T
438
439 bl1_swap_ints( m, n );
440
443 }
444 else // if ( bl1_is_row_storage( a_rs, a_cs ) )
445 {
446 // requested operation: B_r := tr( uplo( A_r ) ) * B_r
447 // effective operation: B_c := B_c * tr( ~uplo( A_c ) )
450
451 bl1_swap_ints( m, n );
452
454 bl1_toggle_uplo( uplo );
455 }
456 }
457
458 // Initialize with values assuming that trans is not conjnotrans.
459 a_conj = a;
460 lda_conj = lda;
461 inca_conj = inca;
462
463 // We want to handle the conjnotrans case. The easiest way to do so is
464 // by making a conjugated copy of A.
466 {
467 int dim_a;
468
470
472 lda_conj = dim_a;
473 inca_conj = 1;
474
475 bl1_zcopymrt( uplo,
477 dim_a,
478 dim_a,
479 a, inca, lda,
481 }
482 else if ( bl1_is_conjnotrans( trans ) && a_was_copied )
483 {
484 int dim_a;
485
487
488 bl1_zconjmr( uplo,
489 dim_a,
490 dim_a,
492 }
493
495 uplo,
496 trans,
497 diag,
498 m,
499 n,
500 alpha,
502 b, ldb );
503
505 bl1_zfree( a_conj );
506
507 // Free any temporary contiguous matrices, copying the result back to
508 // the original matrix.
510 &a, &a_rs, &a_cs );
511
513 n_save,
515 &b, &b_rs, &b_cs );
516}
void bl1_zconjmr(uplo1_t uplo, int m, int n, dcomplex *a, int a_rs, int a_cs)
Definition bl1_conjmr.c:79
void bl1_zcopymrt(uplo1_t uplo, trans1_t trans, int m, int n, dcomplex *a, int a_rs, int a_cs, dcomplex *b, int b_rs, int b_cs)
Definition bl1_copymrt.c:328
void bl1_ztrsm_blas(side1_t side, uplo1_t uplo, trans1_t trans, diag1_t diag, int m, int n, dcomplex *alpha, dcomplex *a, int lda, dcomplex *b, int ldb)
Definition bl1_trsm.c:661
void bl1_zcreate_contigm(int m, int n, dcomplex *a_save, int a_rs_save, int a_cs_save, dcomplex **a, int *a_rs, int *a_cs)
Definition bl1_create_contigm.c:115
void bl1_zcreate_contigmr(uplo1_t uplo, int m, int n, dcomplex *a_save, int a_rs_save, int a_cs_save, dcomplex **a, int *a_rs, int *a_cs)
Definition bl1_create_contigmr.c:109
dcomplex * bl1_zallocm(unsigned int m, unsigned int n)
Definition bl1_allocm.c:45
void bl1_zfree(dcomplex *p)
Definition bl1_free.c:45
void bl1_zfree_contigm(dcomplex *a_save, int a_rs_save, int a_cs_save, dcomplex **a, int *a_rs, int *a_cs)
Definition bl1_free_contigm.c:61
void bl1_zfree_saved_contigm(int m, int n, dcomplex *a_save, int a_rs_save, int a_cs_save, dcomplex **a, int *a_rs, int *a_cs)
Definition bl1_free_saved_contigm.c:82
Definition blis_type_defs.h:138

References bl1_is_col_storage(), bl1_is_conjnotrans(), bl1_set_dim_with_side(), bl1_zallocm(), bl1_zconjmr(), bl1_zcopymrt(), bl1_zcreate_contigm(), bl1_zcreate_contigmr(), bl1_zero_dim2(), bl1_zfree(), bl1_zfree_contigm(), bl1_zfree_saved_contigm(), bl1_ztrsm_blas(), and BLIS1_CONJ_NO_TRANSPOSE.

Referenced by bl1_ztrsmsx(), FLA_LU_nopiv_opz_var1(), FLA_LU_nopiv_opz_var2(), FLA_LU_nopiv_opz_var3(), FLA_LU_piv_opz_var3(), and FLA_Trsm_external().

◆ bl1_ztrsm_blas()

void bl1_ztrsm_blas ( side1_t  side,
uplo1_t  uplo,
trans1_t  trans,
diag1_t  diag,
int  m,
int  n,
dcomplex alpha,
dcomplex a,
int  lda,
dcomplex b,
int  ldb 
)
662{
663#ifdef BLIS1_ENABLE_CBLAS_INTERFACES
669
674
680 m,
681 n,
682 alpha,
683 a, lda,
684 b, ldb );
685#else
686 char blas_side;
687 char blas_uplo;
688 char blas_trans;
689 char blas_diag;
690
695
697 &blas_uplo,
698 &blas_trans,
699 &blas_diag,
700 &m,
701 &n,
702 alpha,
703 a, &lda,
704 b, &ldb );
705#endif
706}
void F77_ztrsm(char *side, char *uplo, char *transa, char *diag, int *m, int *n, dcomplex *alpha, dcomplex *a, int *lda, dcomplex *b, int *ldb)
void cblas_ztrsm(const enum CBLAS_ORDER Order, const enum CBLAS_SIDE Side, const enum CBLAS_UPLO Uplo, const enum CBLAS_TRANSPOSE TransA, const enum CBLAS_DIAG Diag, const int M, const int N, const void *alpha, const void *A, const int lda, void *B, const int ldb)

References bl1_param_map_to_netlib_diag(), bl1_param_map_to_netlib_side(), bl1_param_map_to_netlib_trans(), bl1_param_map_to_netlib_uplo(), cblas_ztrsm(), CblasColMajor, and F77_ztrsm().

Referenced by bl1_ztrsm().