1 /*****************************************************************************
2 Copyright (c) 2014, Intel Corp.
3 All rights reserved.
4
5 Redistribution and use in source and binary forms, with or without
6 modification, are permitted provided that the following conditions are met:
7
8 * Redistributions of source code must retain the above copyright notice,
9 this list of conditions and the following disclaimer.
10 * Redistributions in binary form must reproduce the above copyright
11 notice, this list of conditions and the following disclaimer in the
12 documentation and/or other materials provided with the distribution.
13 * Neither the name of Intel Corporation nor the names of its contributors
14 may be used to endorse or promote products derived from this software
15 without specific prior written permission.
16
17 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23 SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24 INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25 CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26 ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF
27 THE POSSIBILITY OF SUCH DAMAGE.
28 *****************************************************************************
29 * Contents: Native middle-level C interface to LAPACK function dggevx
30 * Author: Intel Corporation
31 *****************************************************************************/
32
33 #include "lapacke_utils.h"
34
LAPACKE_dggevx_work(int matrix_layout,char balanc,char jobvl,char jobvr,char sense,lapack_int n,double * a,lapack_int lda,double * b,lapack_int ldb,double * alphar,double * alphai,double * beta,double * vl,lapack_int ldvl,double * vr,lapack_int ldvr,lapack_int * ilo,lapack_int * ihi,double * lscale,double * rscale,double * abnrm,double * bbnrm,double * rconde,double * rcondv,double * work,lapack_int lwork,lapack_int * iwork,lapack_logical * bwork)35 lapack_int LAPACKE_dggevx_work( int matrix_layout, char balanc, char jobvl,
36 char jobvr, char sense, lapack_int n, double* a,
37 lapack_int lda, double* b, lapack_int ldb,
38 double* alphar, double* alphai, double* beta,
39 double* vl, lapack_int ldvl, double* vr,
40 lapack_int ldvr, lapack_int* ilo,
41 lapack_int* ihi, double* lscale, double* rscale,
42 double* abnrm, double* bbnrm, double* rconde,
43 double* rcondv, double* work, lapack_int lwork,
44 lapack_int* iwork, lapack_logical* bwork )
45 {
46 lapack_int info = 0;
47 if( matrix_layout == LAPACK_COL_MAJOR ) {
48 /* Call LAPACK function and adjust info */
49 LAPACK_dggevx( &balanc, &jobvl, &jobvr, &sense, &n, a, &lda, b, &ldb,
50 alphar, alphai, beta, vl, &ldvl, vr, &ldvr, ilo, ihi,
51 lscale, rscale, abnrm, bbnrm, rconde, rcondv, work,
52 &lwork, iwork, bwork, &info );
53 if( info < 0 ) {
54 info = info - 1;
55 }
56 } else if( matrix_layout == LAPACK_ROW_MAJOR ) {
57 lapack_int lda_t = MAX(1,n);
58 lapack_int ldb_t = MAX(1,n);
59 lapack_int ldvl_t = MAX(1,n);
60 lapack_int ldvr_t = MAX(1,n);
61 double* a_t = NULL;
62 double* b_t = NULL;
63 double* vl_t = NULL;
64 double* vr_t = NULL;
65 /* Check leading dimension(s) */
66 if( lda < n ) {
67 info = -8;
68 LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
69 return info;
70 }
71 if( ldb < n ) {
72 info = -10;
73 LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
74 return info;
75 }
76 if( ldvl < n ) {
77 info = -15;
78 LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
79 return info;
80 }
81 if( ldvr < n ) {
82 info = -17;
83 LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
84 return info;
85 }
86 /* Query optimal working array(s) size if requested */
87 if( lwork == -1 ) {
88 LAPACK_dggevx( &balanc, &jobvl, &jobvr, &sense, &n, a, &lda_t, b,
89 &ldb_t, alphar, alphai, beta, vl, &ldvl_t, vr,
90 &ldvr_t, ilo, ihi, lscale, rscale, abnrm, bbnrm,
91 rconde, rcondv, work, &lwork, iwork, bwork, &info );
92 return (info < 0) ? (info - 1) : info;
93 }
94 /* Allocate memory for temporary array(s) */
95 a_t = (double*)LAPACKE_malloc( sizeof(double) * lda_t * MAX(1,n) );
96 if( a_t == NULL ) {
97 info = LAPACK_TRANSPOSE_MEMORY_ERROR;
98 goto exit_level_0;
99 }
100 b_t = (double*)LAPACKE_malloc( sizeof(double) * ldb_t * MAX(1,n) );
101 if( b_t == NULL ) {
102 info = LAPACK_TRANSPOSE_MEMORY_ERROR;
103 goto exit_level_1;
104 }
105 if( LAPACKE_lsame( jobvl, 'v' ) ) {
106 vl_t = (double*)
107 LAPACKE_malloc( sizeof(double) * ldvl_t * MAX(1,n) );
108 if( vl_t == NULL ) {
109 info = LAPACK_TRANSPOSE_MEMORY_ERROR;
110 goto exit_level_2;
111 }
112 }
113 if( LAPACKE_lsame( jobvr, 'v' ) ) {
114 vr_t = (double*)
115 LAPACKE_malloc( sizeof(double) * ldvr_t * MAX(1,n) );
116 if( vr_t == NULL ) {
117 info = LAPACK_TRANSPOSE_MEMORY_ERROR;
118 goto exit_level_3;
119 }
120 }
121 /* Transpose input matrices */
122 LAPACKE_dge_trans( matrix_layout, n, n, a, lda, a_t, lda_t );
123 LAPACKE_dge_trans( matrix_layout, n, n, b, ldb, b_t, ldb_t );
124 /* Call LAPACK function and adjust info */
125 LAPACK_dggevx( &balanc, &jobvl, &jobvr, &sense, &n, a_t, &lda_t, b_t,
126 &ldb_t, alphar, alphai, beta, vl_t, &ldvl_t, vr_t,
127 &ldvr_t, ilo, ihi, lscale, rscale, abnrm, bbnrm, rconde,
128 rcondv, work, &lwork, iwork, bwork, &info );
129 if( info < 0 ) {
130 info = info - 1;
131 }
132 /* Transpose output matrices */
133 LAPACKE_dge_trans( LAPACK_COL_MAJOR, n, n, a_t, lda_t, a, lda );
134 LAPACKE_dge_trans( LAPACK_COL_MAJOR, n, n, b_t, ldb_t, b, ldb );
135 if( LAPACKE_lsame( jobvl, 'v' ) ) {
136 LAPACKE_dge_trans( LAPACK_COL_MAJOR, n, n, vl_t, ldvl_t, vl, ldvl );
137 }
138 if( LAPACKE_lsame( jobvr, 'v' ) ) {
139 LAPACKE_dge_trans( LAPACK_COL_MAJOR, n, n, vr_t, ldvr_t, vr, ldvr );
140 }
141 /* Release memory and exit */
142 if( LAPACKE_lsame( jobvr, 'v' ) ) {
143 LAPACKE_free( vr_t );
144 }
145 exit_level_3:
146 if( LAPACKE_lsame( jobvl, 'v' ) ) {
147 LAPACKE_free( vl_t );
148 }
149 exit_level_2:
150 LAPACKE_free( b_t );
151 exit_level_1:
152 LAPACKE_free( a_t );
153 exit_level_0:
154 if( info == LAPACK_TRANSPOSE_MEMORY_ERROR ) {
155 LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
156 }
157 } else {
158 info = -1;
159 LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
160 }
161 return info;
162 }
163