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 * Generated November 2015
32 *****************************************************************************/
33 
34 #include "lapacke_utils.h"
35 
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)36 lapack_int LAPACKE_dggevx_work( int matrix_layout, char balanc, char jobvl,
37                                 char jobvr, char sense, lapack_int n, double* a,
38                                 lapack_int lda, double* b, lapack_int ldb,
39                                 double* alphar, double* alphai, double* beta,
40                                 double* vl, lapack_int ldvl, double* vr,
41                                 lapack_int ldvr, lapack_int* ilo,
42                                 lapack_int* ihi, double* lscale, double* rscale,
43                                 double* abnrm, double* bbnrm, double* rconde,
44                                 double* rcondv, double* work, lapack_int lwork,
45                                 lapack_int* iwork, lapack_logical* bwork )
46 {
47     lapack_int info = 0;
48     if( matrix_layout == LAPACK_COL_MAJOR ) {
49         /* Call LAPACK function and adjust info */
50         LAPACK_dggevx( &balanc, &jobvl, &jobvr, &sense, &n, a, &lda, b, &ldb,
51                        alphar, alphai, beta, vl, &ldvl, vr, &ldvr, ilo, ihi,
52                        lscale, rscale, abnrm, bbnrm, rconde, rcondv, work,
53                        &lwork, iwork, bwork, &info );
54         if( info < 0 ) {
55             info = info - 1;
56         }
57     } else if( matrix_layout == LAPACK_ROW_MAJOR ) {
58         lapack_int lda_t = MAX(1,n);
59         lapack_int ldb_t = MAX(1,n);
60         lapack_int ldvl_t = MAX(1,n);
61         lapack_int ldvr_t = MAX(1,n);
62         double* a_t = NULL;
63         double* b_t = NULL;
64         double* vl_t = NULL;
65         double* vr_t = NULL;
66         /* Check leading dimension(s) */
67         if( lda < n ) {
68             info = -8;
69             LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
70             return info;
71         }
72         if( ldb < n ) {
73             info = -10;
74             LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
75             return info;
76         }
77         if( ldvl < n ) {
78             info = -15;
79             LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
80             return info;
81         }
82         if( ldvr < n ) {
83             info = -17;
84             LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
85             return info;
86         }
87         /* Query optimal working array(s) size if requested */
88         if( lwork == -1 ) {
89             LAPACK_dggevx( &balanc, &jobvl, &jobvr, &sense, &n, a, &lda_t, b,
90                            &ldb_t, alphar, alphai, beta, vl, &ldvl_t, vr,
91                            &ldvr_t, ilo, ihi, lscale, rscale, abnrm, bbnrm,
92                            rconde, rcondv, work, &lwork, iwork, bwork, &info );
93             return (info < 0) ? (info - 1) : info;
94         }
95         /* Allocate memory for temporary array(s) */
96         a_t = (double*)LAPACKE_malloc( sizeof(double) * lda_t * MAX(1,n) );
97         if( a_t == NULL ) {
98             info = LAPACK_TRANSPOSE_MEMORY_ERROR;
99             goto exit_level_0;
100         }
101         b_t = (double*)LAPACKE_malloc( sizeof(double) * ldb_t * MAX(1,n) );
102         if( b_t == NULL ) {
103             info = LAPACK_TRANSPOSE_MEMORY_ERROR;
104             goto exit_level_1;
105         }
106         if( LAPACKE_lsame( jobvl, 'v' ) ) {
107             vl_t = (double*)
108                 LAPACKE_malloc( sizeof(double) * ldvl_t * MAX(1,n) );
109             if( vl_t == NULL ) {
110                 info = LAPACK_TRANSPOSE_MEMORY_ERROR;
111                 goto exit_level_2;
112             }
113         }
114         if( LAPACKE_lsame( jobvr, 'v' ) ) {
115             vr_t = (double*)
116                 LAPACKE_malloc( sizeof(double) * ldvr_t * MAX(1,n) );
117             if( vr_t == NULL ) {
118                 info = LAPACK_TRANSPOSE_MEMORY_ERROR;
119                 goto exit_level_3;
120             }
121         }
122         /* Transpose input matrices */
123         LAPACKE_dge_trans( matrix_layout, n, n, a, lda, a_t, lda_t );
124         LAPACKE_dge_trans( matrix_layout, n, n, b, ldb, b_t, ldb_t );
125         /* Call LAPACK function and adjust info */
126         LAPACK_dggevx( &balanc, &jobvl, &jobvr, &sense, &n, a_t, &lda_t, b_t,
127                        &ldb_t, alphar, alphai, beta, vl_t, &ldvl_t, vr_t,
128                        &ldvr_t, ilo, ihi, lscale, rscale, abnrm, bbnrm, rconde,
129                        rcondv, work, &lwork, iwork, bwork, &info );
130         if( info < 0 ) {
131             info = info - 1;
132         }
133         /* Transpose output matrices */
134         LAPACKE_dge_trans( LAPACK_COL_MAJOR, n, n, a_t, lda_t, a, lda );
135         LAPACKE_dge_trans( LAPACK_COL_MAJOR, n, n, b_t, ldb_t, b, ldb );
136         if( LAPACKE_lsame( jobvl, 'v' ) ) {
137             LAPACKE_dge_trans( LAPACK_COL_MAJOR, n, n, vl_t, ldvl_t, vl, ldvl );
138         }
139         if( LAPACKE_lsame( jobvr, 'v' ) ) {
140             LAPACKE_dge_trans( LAPACK_COL_MAJOR, n, n, vr_t, ldvr_t, vr, ldvr );
141         }
142         /* Release memory and exit */
143         if( LAPACKE_lsame( jobvr, 'v' ) ) {
144             LAPACKE_free( vr_t );
145         }
146 exit_level_3:
147         if( LAPACKE_lsame( jobvl, 'v' ) ) {
148             LAPACKE_free( vl_t );
149         }
150 exit_level_2:
151         LAPACKE_free( b_t );
152 exit_level_1:
153         LAPACKE_free( a_t );
154 exit_level_0:
155         if( info == LAPACK_TRANSPOSE_MEMORY_ERROR ) {
156             LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
157         }
158     } else {
159         info = -1;
160         LAPACKE_xerbla( "LAPACKE_dggevx_work", info );
161     }
162     return info;
163 }
164