LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
example_DGESV_colmajor.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 (col-major, high-level) Example Program Results
27
28 -- LAPACKE Example routine --
29 -- LAPACK is a software package provided by Univ. of Tennessee, --
30 -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
31*/
32/* Includes */
33#include <stdlib.h>
34#include <stdio.h>
35#include <string.h>
36#include "lapacke.h"
37#include "lapacke_example_aux.h"
38
39/* Main program */
40int main(int argc, char **argv) {
41
42 /* Locals */
43 lapack_int n, nrhs, lda, ldb, info;
44 int i, j;
45 /* Local arrays */
46 double *A, *b;
47 lapack_int *ipiv;
48
49 /* Default Value */
50 n = 5; nrhs = 1;
51
52 /* Arguments */
53 for( i = 1; i < argc; i++ ) {
54 if( strcmp( argv[i], "-n" ) == 0 ) {
55 n = atoi(argv[i+1]);
56 i++;
57 }
58 if( strcmp( argv[i], "-nrhs" ) == 0 ) {
59 nrhs = atoi(argv[i+1]);
60 i++;
61 }
62 }
63
64 /* Initialization */
65 lda=n, ldb=n;
66 A = (double *)malloc(n*n*sizeof(double)) ;
67 if (A==NULL){ printf("error of memory allocation\n"); exit(0); }
68 b = (double *)malloc(n*nrhs*sizeof(double)) ;
69 if (b==NULL){ printf("error of memory allocation\n"); exit(0); }
70 ipiv = (lapack_int *)malloc(n*sizeof(lapack_int)) ;
71 if (ipiv==NULL){ printf("error of memory allocation\n"); exit(0); }
72
73 for( i = 0; i < n; i++ ) {
74 for( j = 0; j < n; j++ ) A[i+j*lda] = ((double) rand()) / ((double) RAND_MAX) - 0.5;
75 }
76
77 for(i=0;i<n*nrhs;i++)
78 b[i] = ((double) rand()) / ((double) RAND_MAX) - 0.5;
79
80 /* Print Entry Matrix */
81 print_matrix_colmajor( "Entry Matrix A", n, n, A, lda );
82 /* Print Right Rand Side */
83 print_matrix_colmajor( "Right Rand Side b", n, nrhs, b, ldb );
84 printf( "\n" );
85
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_COL_MAJOR, n, nrhs, A, lda, ipiv,
90 b, ldb );
91
92 /* Check for the exact singularity */
93 if( info > 0 ) {
94 printf( "The diagonal element of the triangular factor of A,\n" );
95 printf( "U(%" LAPACK_IFMT ",%" LAPACK_IFMT ") is zero, so that A is singular;\n", info, info );
96 printf( "the solution could not be computed.\n" );
97 exit( 1 );
98 }
99 if (info <0) exit( 1 );
100 /* Print solution */
101 print_matrix_colmajor( "Solution", n, nrhs, b, ldb );
102 /* Print details of LU factorization */
103 print_matrix_colmajor( "Details of LU factorization", n, n, A, lda );
104 /* Print pivot indices */
105 print_vector( "Pivot indices", n, ipiv );
106 exit( 0 );
107} /* End of LAPACKE_dgesv Example */
108
int main()
Definition: cblas_example1.c:7
#define lapack_int
Definition: lapack.h:87
#define LAPACK_IFMT
Definition: lapack.h:98
#define LAPACK_COL_MAJOR
Definition: lapacke.h:53
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:35
void print_vector(char *desc, lapack_int n, lapack_int *vec)
void print_matrix_colmajor(char *desc, lapack_int m, lapack_int n, double *mat, lapack_int ldm)