LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
example_DGESV_rowmajor.c
Go to the documentation of this file.
1 /*
2  LAPACKE_dgesv Example
3  =====================
4 
5  The program computes the solution to the system of linear
6  equations with a square matrix A and multiple
7  right-hand sides B, where A is the coefficient matrix
8  and b is the right-hand side matrix:
9 
10  Description
11  ===========
12 
13  The routine solves for X the system of linear equations A*X = B,
14  where A is an n-by-n matrix, the columns of matrix B are individual
15  right-hand sides, and the columns of X are the corresponding
16  solutions.
17 
18  The LU decomposition with partial pivoting and row interchanges is
19  used to factor A as A = P*L*U, where P is a permutation matrix, L
20  is unit lower triangular, and U is upper triangular. The factored
21  form of A is then used to solve the system of equations A*X = B.
22 
23  LAPACKE Interface
24  =================
25 
26  LAPACKE_dgesv (row-major, high-level) Example Program Results
27 
28  -- LAPACKE Example routine (version 3.6.0) --
29  -- LAPACK is a software package provided by Univ. of Tennessee, --
30  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
31  November 2015
32 
33 */
34 #include <stdlib.h>
35 #include <stdio.h>
36 #include <string.h>
37 #include <lapacke.h>
38 #include "lapacke_example_aux.h"
39 
40 /* Main program */
41 int main(int argc, char **argv) {
42 
43  /* Locals */
44  lapack_int n, nrhs, lda, ldb, info;
45  int i, j;
46  /* Local arrays */
47  double *A, *b;
48  lapack_int *ipiv;
49 
50  /* Default Value */
51  n = 5; nrhs = 1;
52 
53  /* Arguments */
54  for( i = 1; i < argc; i++ ) {
55  if( strcmp( argv[i], "-n" ) == 0 ) {
56  n = atoi(argv[i+1]);
57  i++;
58  }
59  if( strcmp( argv[i], "-nrhs" ) == 0 ) {
60  nrhs = atoi(argv[i+1]);
61  i++;
62  }
63  }
64 
65  /* Initialization */
66  lda=n, ldb=nrhs;
67  A = (double *)malloc(n*n*sizeof(double)) ;
68  if (A==NULL){ printf("error of memory allocation\n"); exit(0); }
69  b = (double *)malloc(n*nrhs*sizeof(double)) ;
70  if (b==NULL){ printf("error of memory allocation\n"); exit(0); }
71  ipiv = (lapack_int *)malloc(n*sizeof(lapack_int)) ;
72  if (ipiv==NULL){ printf("error of memory allocation\n"); exit(0); }
73 
74  for( i = 0; i < n; i++ ) {
75  for( j = 0; j < n; j++ ) A[i*lda+j] = ((double) rand()) / ((double) RAND_MAX) - 0.5;
76  }
77 
78  for(i=0;i<n*nrhs;i++)
79  b[i] = ((double) rand()) / ((double) RAND_MAX) - 0.5;
80 
81  /* Print Entry Matrix */
82  print_matrix_rowmajor( "Entry Matrix A", n, n, A, lda );
83  /* Print Right Rand Side */
84  print_matrix_rowmajor( "Right Rand Side b", n, nrhs, b, ldb );
85  printf( "\n" );
86  /* Executable statements */
87  printf( "LAPACKE_dgesv (row-major, high-level) Example Program Results\n" );
88  /* Solve the equations A*X = B */
89  info = LAPACKE_dgesv( LAPACK_ROW_MAJOR, n, nrhs, A, lda, ipiv,
90  b, ldb );
91  /* Check for the exact singularity */
92  if( info > 0 ) {
93  printf( "The diagonal element of the triangular factor of A,\n" );
94  printf( "U(%i,%i) is zero, so that A is singular;\n", info, info );
95  printf( "the solution could not be computed.\n" );
96  exit( 1 );
97  }
98  if (info <0) exit( 1 );
99  /* Print solution */
100  print_matrix_rowmajor( "Solution", n, nrhs, b, ldb );
101  /* Print details of LU factorization */
102  print_matrix_rowmajor( "Details of LU factorization", n, n, A, lda );
103  /* Print pivot indices */
104  print_vector( "Pivot indices", n, ipiv );
105  exit( 0 );
106 } /* End of LAPACKE_dgesv Example */
107 
#define LAPACK_ROW_MAJOR
Definition: lapacke.h:119
int main(int argc, char **argv)
void print_vector(char *desc, lapack_int n, lapack_int *vec)
lapack_int LAPACKE_dgesv(int matrix_layout, lapack_int n, lapack_int nrhs, double *a, lapack_int lda, lapack_int *ipiv, double *b, lapack_int ldb)
Definition: lapacke_dgesv.c:36
#define lapack_int
Definition: lapacke.h:47
void print_matrix_rowmajor(char *desc, lapack_int m, lapack_int n, double *mat, lapack_int ldm)