PLASMA  2.4.5
PLASMA - Parallel Linear Algebra for Scalable Multi-core Architectures
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups
example_sgetrs.c
Go to the documentation of this file.
1 
17 #include <stdlib.h>
18 #include <stdio.h>
19 #include <string.h>
20 #include <math.h>
21 
22 #include <plasma.h>
23 #include <cblas.h>
24 #include <lapacke.h>
25 #include <core_blas.h>
26 
27 int check_solution(int, int , float *, int, float *, float *, int);
28 
29 int IONE=1;
30 int ISEED[4] = {0,0,0,1}; /* initial seed for slarnv() */
31 
32 int main ()
33 {
34 
35  int cores = 2;
36  int N = 10;
37  int LDA = 10;
38  int NRHS = 5;
39  int LDB = 10;
40  int info;
41  int info_solution;
42  int i,j;
43  int LDAxN = LDA*N;
44  int LDBxNRHS = LDB*NRHS;
45 
46  float *A1 = (float *)malloc(LDA*N*(sizeof*A1));
47  float *A2 = (float *)malloc(LDA*N*(sizeof*A2));
48  float *B1 = (float *)malloc(LDB*NRHS*(sizeof*B1));
49  float *B2 = (float *)malloc(LDB*NRHS*(sizeof*B2));
50  float *L;
51  int *IPIV;
52 
53  /* Check if unable to allocate memory */
54  if ((!A1)||(!A2)||(!B1)||(!B2)){
55  printf("Out of Memory \n ");
56  return EXIT_SUCCESS;
57  }
58 
59  /*Plasma Initialize*/
60  PLASMA_Init(cores);
61  printf("-- PLASMA is initialized to run on %d cores. \n",cores);
62 
63  /* Initialize A1 and A2 Matrix */
64  LAPACKE_slarnv_work(IONE, ISEED, LDAxN, A1);
65  for ( i = 0; i < N; i++)
66  for ( j = 0; j < N; j++)
67  A2[LDA*j+i] = A1[LDA*j+i];
68 
69  /* Initialize B1 and B2 */
70  LAPACKE_slarnv_work(IONE, ISEED, LDBxNRHS, B1);
71  for ( i = 0; i < N; i++)
72  for ( j = 0; j < NRHS; j++)
73  B2[LDB*j+i] = B1[LDB*j+i];
74 
75 
76  /* Allocate L and IPIV */
77  info = PLASMA_Alloc_Workspace_sgetrf_incpiv(N, N, &L, &IPIV);
78 
79  /* LU factorization of the matrix A */
80  info = PLASMA_sgetrf_incpiv(N, N, A2, LDA, L, IPIV);
81 
82  /* Solve the problem */
83  info = PLASMA_sgetrs_incpiv(PlasmaNoTrans, N, NRHS, A2, LDA, L, IPIV, B2, LDB);
84 
85  /* Check the solution */
86  info_solution = check_solution(N, NRHS, A1, LDA, B1, B2, LDB);
87 
88  if ((info_solution != 0)|(info != 0))
89  printf("-- Error in SGETRS example ! \n");
90  else
91  printf("-- Run of SGETRS example successful ! \n");
92 
93  free(A1); free(A2); free(B1); free(B2); free(IPIV); free(L);
94 
96 
97  return EXIT_SUCCESS;
98 }
99 
100 /*------------------------------------------------------------------------
101  * Check the accuracy of the solution of the linear system
102  */
103 
104 int check_solution(int N, int NRHS, float *A1, int LDA, float *B1, float *B2, int LDB)
105 {
106  int info_solution;
107  float Rnorm, Anorm, Xnorm, Bnorm;
108  float alpha, beta;
109  float *work = (float *)malloc(N*sizeof(float));
110  float eps;
111 
112  eps = LAPACKE_slamch_work('e');
113 
114  alpha = 1.0;
115  beta = -1.0;
116 
117  Xnorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B2, LDB, work);
118  Anorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, N, A1, LDA, work);
119  Bnorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
120 
121  cblas_sgemm(CblasColMajor, CblasNoTrans, CblasNoTrans, N, NRHS, N, (alpha), A1, LDA, B2, LDB, (beta), B1, LDB);
122  Rnorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
123 
124  printf("============\n");
125  printf("Checking the Residual of the solution \n");
126  printf("-- ||Ax-B||_oo/((||A||_oo||x||_oo+||B||_oo).N.eps) = %e \n",Rnorm/((Anorm*Xnorm+Bnorm)*N*eps));
127 
128  if ( isnan(Rnorm/((Anorm*Xnorm+Bnorm)*N*eps)) || (Rnorm/((Anorm*Xnorm+Bnorm)*N*eps) > 10.0) ){
129  printf("-- The solution is suspicious ! \n");
130  info_solution = 1;
131  }
132  else{
133  printf("-- The solution is CORRECT ! \n");
134  info_solution = 0;
135  }
136 
137  free(work);
138 
139  return info_solution;
140 }