You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

example_DGESV_colmajor.c 3.6 kB

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