aboutsummaryrefslogtreecommitdiff
path: root/LAPACKE/example/example_DGESV_colmajor.c
blob: 938bbfe36a811a2ee1734a1c77f2e5362a8ccecd (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
/*
   LAPACKE_dgesv Example
   =====================

   The program computes the solution to the system of linear
   equations with a square matrix A and multiple
   right-hand sides B, where A is the coefficient matrix
   and b is the right-hand side matrix:

   Description
   ===========

   The routine solves for X the system of linear equations A*X = B,
   where A is an n-by-n matrix, the columns of matrix B are individual
   right-hand sides, and the columns of X are the corresponding
   solutions.

   The LU decomposition with partial pivoting and row interchanges is
   used to factor A as A = P*L*U, where P is a permutation matrix, L
   is unit lower triangular, and U is upper triangular. The factored
   form of A is then used to solve the system of equations A*X = B.

   LAPACKE Interface
   =================

   LAPACKE_dgesv (col-major, high-level) Example Program Results

  -- LAPACKE Example routine (version 3.6.0) --
  -- LAPACK is a software package provided by Univ. of Tennessee,    --
  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
     November 2015

*/
/* Includes */
#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include "lapacke.h"
#include "lapacke_example_aux.h"

/* Main program */
int main(int argc, char **argv) {

        /* Locals */
        lapack_int n, nrhs, lda, ldb, info;
		int i, j;
        /* Local arrays */
		double *A, *b;
		lapack_int *ipiv;

        /* Default Value */
	    n = 5; nrhs = 1;

        /* Arguments */
	    for( i = 1; i < argc; i++ ) {
	    	if( strcmp( argv[i], "-n" ) == 0 ) {
		    	n  = atoi(argv[i+1]);
			    i++;
		    }
			if( strcmp( argv[i], "-nrhs" ) == 0 ) {
				nrhs  = atoi(argv[i+1]);
				i++;
			}
		}

        /* Initialization */
        lda=n, ldb=n;
		A = (double *)malloc(n*n*sizeof(double)) ;
		if (A==NULL){ printf("error of memory allocation\n"); exit(0); }
		b = (double *)malloc(n*nrhs*sizeof(double)) ;
		if (b==NULL){ printf("error of memory allocation\n"); exit(0); }
		ipiv = (lapack_int *)malloc(n*sizeof(lapack_int)) ;
		if (ipiv==NULL){ printf("error of memory allocation\n"); exit(0); }

        for( i = 0; i < n; i++ ) {
                for( j = 0; j < n; j++ ) A[i+j*lda] = ((double) rand()) / ((double) RAND_MAX) - 0.5;
		}

		for(i=0;i<n*nrhs;i++)
			b[i] = ((double) rand()) / ((double) RAND_MAX) - 0.5;

        /* Print Entry Matrix */
        print_matrix_colmajor( "Entry Matrix A", n, n, A, lda );
        /* Print Right Rand Side */
        print_matrix_colmajor( "Right Rand Side b", n, nrhs, b, ldb );
        printf( "\n" );

        /* Executable statements */
        printf( "LAPACKE_dgesv (row-major, high-level) Example Program Results\n" );
        /* Solve the equations A*X = B */
        info = LAPACKE_dgesv( LAPACK_COL_MAJOR, n, nrhs, A, lda, ipiv,
                        b, ldb );

        /* Check for the exact singularity */
        if( info > 0 ) {
                printf( "The diagonal element of the triangular factor of A,\n" );
                printf( "U(%i,%i) is zero, so that A is singular;\n", info, info );
                printf( "the solution could not be computed.\n" );
                exit( 1 );
        }
        if (info <0) exit( 1 );
        /* Print solution */
        print_matrix_colmajor( "Solution", n, nrhs, b, ldb );
        /* Print details of LU factorization */
        print_matrix_colmajor( "Details of LU factorization", n, n, A, lda );
        /* Print pivot indices */
        print_vector( "Pivot indices", n, ipiv );
        exit( 0 );
} /* End of LAPACKE_dgesv Example */