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_sgelqs.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 #ifndef max
28 #define max(a, b) ((a) > (b) ? (a) : (b))
29 #endif
30 #ifndef min
31 #define min(a, b) ((a) < (b) ? (a) : (b))
32 #endif
33 
34 int check_solution(int, int, int, float*, int, float*, float*, int);
35 
36 int IONE=1;
37 int ISEED[4] = {0,0,0,1}; /* initial seed for slarnv() */
38 
39 int main ()
40 {
41 
42  int cores = 2;
43  int M = 10;
44  int N = 15;
45  int LDA = 10;
46  int NRHS = 5;
47  int LDB = 15;
48 
49  int info;
50  int info_solution;
51  int i,j;
52  int LDAxN = LDA*N;
53  int LDBxNRHS = LDB*NRHS;
54 
55  float *A1 = (float *)malloc(LDA*N*sizeof(float));
56  float *A2 = (float *)malloc(LDA*N*sizeof(float));
57  float *B1 = (float *)malloc(LDB*NRHS*sizeof(float));
58  float *B2 = (float *)malloc(LDB*NRHS*sizeof(float));
59  float *T;
60 
61  /* Check if unable to allocate memory */
62  if ((!A1)||(!A2)||(!B1)||(!B2)){
63  printf("Out of Memory \n ");
64  return EXIT_SUCCESS;
65  }
66 
67  /* Plasma Initialization */
68  PLASMA_Init(cores);
69  printf("-- PLASMA is initialized to run on %d cores. \n",cores);
70 
71  /* Allocate T */
73 
74  /* Initialize A1 and A2 */
75  LAPACKE_slarnv_work(IONE, ISEED, LDAxN, A1);
76  for (i = 0; i < M; i++)
77  for (j = 0; j < N; j++)
78  A2[LDA*j+i] = A1[LDA*j+i] ;
79 
80  /* Initialize B1 and B2 */
81  LAPACKE_slarnv_work(IONE, ISEED, LDBxNRHS, B1);
82  for (i = 0; i < M; i++)
83  for (j = 0; j < NRHS; j++)
84  B2[LDB*j+i] = B1[LDB*j+i] ;
85 
86  /* Factorization QR of the matrix A2 */
87  info = PLASMA_sgelqf(M, N, A2, LDA, T);
88 
89  /* Solve the problem */
90  info = PLASMA_sgelqs(M, N, NRHS, A2, LDA, T, B2, LDB);
91 
92  /* Check the solution */
93  info_solution = check_solution(M, N, NRHS, A1, LDA, B1, B2, LDB);
94 
95  if ((info_solution != 0)|(info != 0))
96  printf("-- Error in SGELQS example ! \n");
97  else
98  printf("-- Run of SGELQS example successful ! \n");
99 
100  free(A1); free(A2); free(B1); free(B2); free(T);
101 
102  PLASMA_Finalize();
103 
104  return EXIT_SUCCESS;
105 }
106 
107 /*--------------------------------------------------------------
108  * Check the solution
109  */
110 
111 int check_solution(int M, int N, int NRHS, float *A1, int LDA, float *B1, float *B2, int LDB)
112 {
113  int info_solution;
114  float Rnorm, Anorm, Xnorm, Bnorm;
115  float alpha, beta;
116  float *work = (float *)malloc(max(M, N)* sizeof(float));
117  float eps;
118 
119  eps = LAPACKE_slamch_work('e');
120 
121  alpha = 1.0;
122  beta = -1.0;
123 
124  Anorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), M, N, A1, LDA, work);
125  Xnorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), M, NRHS, B2, LDB, work);
126  Bnorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
127 
128  cblas_sgemm(CblasColMajor, CblasNoTrans, CblasNoTrans, M, NRHS, N, (alpha), A1, LDA, B2, LDB, (beta), B1, LDB);
129 
130  if (M >= N) {
131  float *Residual = (float *)malloc(M*NRHS*sizeof(float));
132  memset((void*)Residual, 0, M*NRHS*sizeof(float));
133  cblas_sgemm(CblasColMajor, CblasTrans, CblasNoTrans, N, NRHS, M, (alpha), A1, LDA, B1, LDB, (beta), Residual, M);
134  Rnorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), M, NRHS, Residual, M, work);
135  free(Residual);
136  }
137  else {
138  float *Residual = (float *)malloc(N*NRHS*sizeof(float));
139  memset((void*)Residual, 0, N*NRHS*sizeof(float));
140  cblas_sgemm(CblasColMajor, CblasTrans, CblasNoTrans, N, NRHS, M, (alpha), A1, LDA, B1, LDB, (beta), Residual, N);
141  Rnorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, Residual, N, work);
142  free(Residual);
143  }
144 
145  printf("============\n");
146  printf("Checking the Residual of the solution \n");
147  printf("-- ||Ax-B||_oo/((||A||_oo||x||_oo+||B||)_oo.N.eps) = %e \n",Rnorm/((Anorm*Xnorm+Bnorm)*N*eps));
148 
149  if (isnan(Rnorm / ((Anorm * Xnorm + Bnorm) * N * eps)) || (Rnorm / ((Anorm * Xnorm + Bnorm) * N * eps) > 10.0) ) {
150  printf("-- The solution is suspicious ! \n");
151  info_solution = 1;
152  }
153  else {
154  printf("-- The solution is CORRECT ! \n");
155  info_solution= 0 ;
156  }
157 
158  free(work);
159 
160  return info_solution;
161 }