LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
lapacke_zhprfs_work.c
Go to the documentation of this file.
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 zhprfs
30 * Author: Intel Corporation
31 * Generated November 2015
32 *****************************************************************************/
33 
34 #include "lapacke_utils.h"
35 
36 lapack_int LAPACKE_zhprfs_work( int matrix_layout, char uplo, lapack_int n,
37  lapack_int nrhs,
38  const lapack_complex_double* ap,
39  const lapack_complex_double* afp,
40  const lapack_int* ipiv,
41  const lapack_complex_double* b, lapack_int ldb,
43  double* ferr, double* berr,
44  lapack_complex_double* work, double* rwork )
45 {
46  lapack_int info = 0;
47  if( matrix_layout == LAPACK_COL_MAJOR ) {
48  /* Call LAPACK function and adjust info */
49  LAPACK_zhprfs( &uplo, &n, &nrhs, ap, afp, ipiv, b, &ldb, x, &ldx, ferr,
50  berr, work, rwork, &info );
51  if( info < 0 ) {
52  info = info - 1;
53  }
54  } else if( matrix_layout == LAPACK_ROW_MAJOR ) {
55  lapack_int ldb_t = MAX(1,n);
56  lapack_int ldx_t = MAX(1,n);
57  lapack_complex_double* b_t = NULL;
58  lapack_complex_double* x_t = NULL;
59  lapack_complex_double* ap_t = NULL;
60  lapack_complex_double* afp_t = NULL;
61  /* Check leading dimension(s) */
62  if( ldb < nrhs ) {
63  info = -9;
64  LAPACKE_xerbla( "LAPACKE_zhprfs_work", info );
65  return info;
66  }
67  if( ldx < nrhs ) {
68  info = -11;
69  LAPACKE_xerbla( "LAPACKE_zhprfs_work", info );
70  return info;
71  }
72  /* Allocate memory for temporary array(s) */
73  b_t = (lapack_complex_double*)
75  ldb_t * MAX(1,nrhs) );
76  if( b_t == NULL ) {
78  goto exit_level_0;
79  }
80  x_t = (lapack_complex_double*)
82  ldx_t * MAX(1,nrhs) );
83  if( x_t == NULL ) {
85  goto exit_level_1;
86  }
87  ap_t = (lapack_complex_double*)
89  ( MAX(1,n) * MAX(2,n+1) ) / 2 );
90  if( ap_t == NULL ) {
92  goto exit_level_2;
93  }
94  afp_t = (lapack_complex_double*)
96  ( MAX(1,n) * MAX(2,n+1) ) / 2 );
97  if( afp_t == NULL ) {
99  goto exit_level_3;
100  }
101  /* Transpose input matrices */
102  LAPACKE_zge_trans( matrix_layout, n, nrhs, b, ldb, b_t, ldb_t );
103  LAPACKE_zge_trans( matrix_layout, n, nrhs, x, ldx, x_t, ldx_t );
104  LAPACKE_zhp_trans( matrix_layout, uplo, n, ap, ap_t );
105  LAPACKE_zhp_trans( matrix_layout, uplo, n, afp, afp_t );
106  /* Call LAPACK function and adjust info */
107  LAPACK_zhprfs( &uplo, &n, &nrhs, ap_t, afp_t, ipiv, b_t, &ldb_t, x_t,
108  &ldx_t, ferr, berr, work, rwork, &info );
109  if( info < 0 ) {
110  info = info - 1;
111  }
112  /* Transpose output matrices */
113  LAPACKE_zge_trans( LAPACK_COL_MAJOR, n, nrhs, x_t, ldx_t, x, ldx );
114  /* Release memory and exit */
115  LAPACKE_free( afp_t );
116 exit_level_3:
117  LAPACKE_free( ap_t );
118 exit_level_2:
119  LAPACKE_free( x_t );
120 exit_level_1:
121  LAPACKE_free( b_t );
122 exit_level_0:
123  if( info == LAPACK_TRANSPOSE_MEMORY_ERROR ) {
124  LAPACKE_xerbla( "LAPACKE_zhprfs_work", info );
125  }
126  } else {
127  info = -1;
128  LAPACKE_xerbla( "LAPACKE_zhprfs_work", info );
129  }
130  return info;
131 }
void LAPACK_zhprfs(char *uplo, lapack_int *n, lapack_int *nrhs, const lapack_complex_double *ap, const lapack_complex_double *afp, const lapack_int *ipiv, const lapack_complex_double *b, lapack_int *ldb, lapack_complex_double *x, lapack_int *ldx, double *ferr, double *berr, lapack_complex_double *work, double *rwork, lapack_int *info)
#define LAPACK_ROW_MAJOR
Definition: lapacke.h:119
#define lapack_complex_double
Definition: lapacke.h:90
void LAPACKE_zhp_trans(int matrix_layout, char uplo, lapack_int n, const lapack_complex_double *in, lapack_complex_double *out)
#define MAX(x, y)
Definition: lapacke_utils.h:47
#define LAPACKE_free(p)
Definition: lapacke.h:113
#define LAPACKE_malloc(size)
Definition: lapacke.h:110
lapack_int LAPACKE_zhprfs_work(int matrix_layout, char uplo, lapack_int n, lapack_int nrhs, const lapack_complex_double *ap, const lapack_complex_double *afp, const lapack_int *ipiv, const lapack_complex_double *b, lapack_int ldb, lapack_complex_double *x, lapack_int ldx, double *ferr, double *berr, lapack_complex_double *work, double *rwork)
#define LAPACK_COL_MAJOR
Definition: lapacke.h:120
void LAPACKE_xerbla(const char *name, lapack_int info)
#define lapack_int
Definition: lapacke.h:47
#define LAPACK_TRANSPOSE_MEMORY_ERROR
Definition: lapacke.h:123
void LAPACKE_zge_trans(int matrix_layout, lapack_int m, lapack_int n, const lapack_complex_double *in, lapack_int ldin, lapack_complex_double *out, lapack_int ldout)