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
testing_sgels.c
Go to the documentation of this file.
1 
16 #include <stdlib.h>
17 #include <stdio.h>
18 #include <string.h>
19 #include <math.h>
20 
21 #include <plasma.h>
22 #include <cblas.h>
23 #include <lapacke.h>
24 #include <core_blas.h>
25 #include "testing_smain.h"
26 
27 #undef COMPLEX
28 #define REAL
29 
32  blas_colmajor = 102 };
33 
35  blas_upper = 121,
36  blas_lower = 122 };
37 
39  blas_base = 151,
40  blas_t = 152,
41  blas_rnd = 153,
42  blas_ieee = 154,
43  blas_emin = 155,
44  blas_emax = 156,
45  blas_eps = 157,
46  blas_prec = 158,
49  blas_sfmin = 161};
50 
60 
61 static void
62 BLAS_error(char *rname, int err, int val, int x) {
63  fprintf( stderr, "%s %d %d %d\n", rname, err, val, x );
64  abort();
65 }
66 
67 static
68 void
69 BLAS_ssy_norm(enum blas_order_type order, enum blas_norm_type norm,
70  enum blas_uplo_type uplo, int n, const float *a, int lda, float *res) {
71  int i, j; float anorm, v;
72  char rname[] = "BLAS_ssy_norm";
73 
74  if (order != blas_colmajor) BLAS_error( rname, -1, order, 0 );
75 
76  if (norm == blas_inf_norm) {
77  anorm = 0.0;
78  if (blas_upper == uplo) {
79  for (i = 0; i < n; ++i) {
80  v = 0.0;
81  for (j = 0; j < i; ++j) {
82  v += fabsf( a[j + i * lda] );
83  }
84  for (j = i; j < n; ++j) {
85  v += fabsf( a[i + j * lda] );
86  }
87  if (v > anorm)
88  anorm = v;
89  }
90  } else {
91  BLAS_error( rname, -3, norm, 0 );
92  return;
93  }
94  } else {
95  BLAS_error( rname, -2, norm, 0 );
96  return;
97  }
98 
99  if (res) *res = anorm;
100 }
101 
102 static
103 void
104 BLAS_sge_norm(enum blas_order_type order, enum blas_norm_type norm,
105  int m, int n, const float *a, int lda, float *res) {
106  int i, j; float anorm, v;
107  char rname[] = "BLAS_sge_norm";
108 
109  if (order != blas_colmajor) BLAS_error( rname, -1, order, 0 );
110 
111  if (norm == blas_frobenius_norm) {
112  anorm = 0.0f;
113  for (j = n; j; --j) {
114  for (i = m; i; --i) {
115  v = a[0];
116  anorm += v * v;
117  a++;
118  }
119  a += lda - m;
120  }
121  anorm = sqrt( anorm );
122  } else if (norm == blas_inf_norm) {
123  anorm = 0.0f;
124  for (i = 0; i < m; ++i) {
125  v = 0.0f;
126  for (j = 0; j < n; ++j) {
127  v += fabsf( a[i + j * lda] );
128  }
129  if (v > anorm)
130  anorm = v;
131  }
132  } else {
133  BLAS_error( rname, -2, norm, 0 );
134  return;
135  }
136 
137  if (res) *res = anorm;
138 }
139 
140 static
141 float
142 BLAS_spow_di(float x, int n) {
143  float rv = 1.0;
144 
145  if (n < 0) {
146  n = -n;
147  x = 1.0 / x;
148  }
149 
150  for (; n; n >>= 1, x *= x) {
151  if (n & 1)
152  rv *= x;
153  }
154 
155  return rv;
156 }
157 
158 static
159 float
160 BLAS_sfpinfo(enum blas_cmach_type cmach) {
161  float eps = 1.0, r = 1.0, o = 1.0, b = 2.0;
162  int t = 53, l = 1024, m = -1021;
163  char rname[] = "BLAS_sfpinfo";
164 
165  if ((sizeof eps) == sizeof(float)) {
166  t = 24;
167  l = 128;
168  m = -125;
169  } else {
170  t = 53;
171  l = 1024;
172  m = -1021;
173  }
174 
175  /* for (i = 0; i < t; ++i) eps *= half; */
176  eps = BLAS_spow_di( b, -t );
177  /* for (i = 0; i >= m; --i) r *= half; */
178  r = BLAS_spow_di( b, m-1 );
179 
180  o -= eps;
181  /* for (i = 0; i < l; ++i) o *= b; */
182  o = (o * BLAS_spow_di( b, l-1 )) * b;
183 
184  switch (cmach) {
185  case blas_eps: return eps;
186  case blas_sfmin: return r;
187  default:
188  BLAS_error( rname, -1, cmach, 0 );
189  break;
190  }
191  return 0.0;
192 }
193 
194 static int check_orthogonality(int, int, int, float*, float);
195 static int check_factorization(int, int, float*, float*, int, float*, float);
196 static int check_solution(int, int, int, float*, int, float*, float*, int, float);
197 
198 int testing_sgels(int argc, char **argv)
199 {
200  int mode = 0;
201 
202  if ( argc < 1 ){
203  goto usage;
204  } else {
205  mode = atoi(argv[0]);
206  }
207 
208  /* Check for number of arguments*/
209  if ( ((mode == 0) && (argc != 6)) ||
210  ((mode != 0) && (argc != 7)) ){
211  usage:
212  USAGE("GELS", "MODE M N LDA NRHS LDB [RH]",
213  " - MODE : 0: flat, 1: tree (RH needed)\n"
214  " - M : number of rows of the matrix A\n"
215  " - N : number of columns of the matrix A\n"
216  " - LDA : leading dimension of the matrix A\n"
217  " - NRHS : number of RHS\n"
218  " - LDB : leading dimension of the matrix B\n"
219  " - RH : Size of each subdomains\n");
220  return -1;
221  }
222 
223  int M = atoi(argv[1]);
224  int N = atoi(argv[2]);
225  int LDA = atoi(argv[3]);
226  int NRHS = atoi(argv[4]);
227  int LDB = atoi(argv[5]);
228  int rh;
229 
230  int K = min(M, N);
231  float eps;
232  int info_ortho, info_solution, info_factorization;
233  int i,j;
234  int LDAxN = LDA*N;
235  int LDBxNRHS = LDB*NRHS;
236 
237  float *A1 = (float *)malloc(LDA*N*sizeof(float));
238  float *A2 = (float *)malloc(LDA*N*sizeof(float));
239  float *B1 = (float *)malloc(LDB*NRHS*sizeof(float));
240  float *B2 = (float *)malloc(LDB*NRHS*sizeof(float));
241  float *Q = (float *)malloc(LDA*N*sizeof(float));
242  float *T;
243 
244  /* Check if unable to allocate memory */
245  if ((!A1)||(!A2)||(!B1)||(!B2)||(!Q)){
246  printf("Out of Memory \n ");
247  return -2;
248  }
249 
250  if ( mode ) {
251  rh = atoi(argv[6]);
252 
255  }
256 
258  eps = BLAS_sfpinfo( blas_eps );
259 
260  /*----------------------------------------------------------
261  * TESTING SGELS
262  */
263 
264  /* Initialize A1 and A2 */
265  LAPACKE_slarnv_work(IONE, ISEED, LDAxN, A1);
266  for (i = 0; i < M; i++)
267  for (j = 0; j < N; j++)
268  A2[LDA*j+i] = A1[LDA*j+i] ;
269 
270  /* Initialize B1 and B2 */
271  LAPACKE_slarnv_work(IONE, ISEED, LDBxNRHS, B1);
272  for (i = 0; i < M; i++)
273  for (j = 0; j < NRHS; j++)
274  B2[LDB*j+i] = B1[LDB*j+i] ;
275 
276  memset((void*)Q, 0, LDA*N*sizeof(float));
277  for (i = 0; i < K; i++)
278  Q[LDA*i+i] = 1.0;
279 
280  /* PLASMA SGELS */
281  PLASMA_sgels(PlasmaNoTrans, M, N, NRHS, A2, LDA, T, B2, LDB);
282 
283  /* PLASMA SGELS */
284  if (M >= N)
285  /* Building the economy-size Q */
286  PLASMA_sorgqr(M, N, K, A2, LDA, T, Q, LDA);
287  else
288  /* Building the economy-size Q */
289  PLASMA_sorglq(M, N, K, A2, LDA, T, Q, LDA);
290 
291  printf("\n");
292  printf("------ TESTS FOR PLASMA SGELS ROUTINE ------- \n");
293  printf(" Size of the Matrix %d by %d\n", M, N);
294  printf("\n");
295  printf(" The matrix A is randomly generated for each test.\n");
296  printf("============\n");
297  printf(" The relative machine precision (eps) is to be %e \n",eps);
298  printf(" Computational tests pass if scaled residuals are less than 60.\n");
299 
300  /* Check the orthogonality, factorization and the solution */
301  info_ortho = check_orthogonality(M, N, LDA, Q, eps);
302  info_factorization = check_factorization(M, N, A1, A2, LDA, Q, eps);
303  info_solution = check_solution(M, N, NRHS, A1, LDA, B1, B2, LDB, eps);
304 
305  if ((info_solution == 0)&(info_factorization == 0)&(info_ortho == 0)) {
306  printf("***************************************************\n");
307  printf(" ---- TESTING SGELS ...................... PASSED !\n");
308  printf("***************************************************\n");
309  }
310  else {
311  printf("************************************************\n");
312  printf(" - TESTING SGELS ... FAILED !\n");
313  printf("************************************************\n");
314  }
315 
316  /*-------------------------------------------------------------
317  * TESTING SGEQRF + SGEQRS or SGELQF + SGELQS
318  */
319 
320  /* Initialize A1 and A2 */
321  LAPACKE_slarnv_work(IONE, ISEED, LDAxN, A1);
322  for (i = 0; i < M; i++)
323  for (j = 0; j < N; j++)
324  A2[LDA*j+i] = A1[LDA*j+i];
325 
326  /* Initialize B1 and B2 */
327  LAPACKE_slarnv_work(IONE, ISEED, LDBxNRHS, B1);
328  for (i = 0; i < M; i++)
329  for (j = 0; j < NRHS; j++)
330  B2[LDB*j+i] = B1[LDB*j+i];
331 
332  memset((void*)Q, 0, LDA*N*sizeof(float));
333  for (i = 0; i < K; i++)
334  Q[LDA*i+i] = 1.0;
335 
336  if (M >= N) {
337  printf("\n");
338  printf("------ TESTS FOR PLASMA SGEQRF + SGEQRS ROUTINE ------- \n");
339  printf(" Size of the Matrix %d by %d\n", M, N);
340  printf("\n");
341  printf(" The matrix A is randomly generated for each test.\n");
342  printf("============\n");
343  printf(" The relative machine precision (eps) is to be %e \n", eps);
344  printf(" Computational tests pass if scaled residuals are less than 60.\n");
345 
346  /* Plasma routines */
347  PLASMA_sgeqrf(M, N, A2, LDA, T);
348  PLASMA_sorgqr(M, N, K, A2, LDA, T, Q, LDA);
349  PLASMA_sgeqrs(M, N, NRHS, A2, LDA, T, B2, LDB);
350 
351  /* Check the orthogonality, factorization and the solution */
352  info_ortho = check_orthogonality(M, N, LDA, Q, eps);
353  info_factorization = check_factorization(M, N, A1, A2, LDA, Q, eps);
354  info_solution = check_solution(M, N, NRHS, A1, LDA, B1, B2, LDB, eps);
355 
356  if ((info_solution == 0)&(info_factorization == 0)&(info_ortho == 0)) {
357  printf("***************************************************\n");
358  printf(" ---- TESTING SGEQRF + SGEQRS ............ PASSED !\n");
359  printf("***************************************************\n");
360  }
361  else{
362  printf("***************************************************\n");
363  printf(" - TESTING SGEQRF + SGEQRS ... FAILED !\n");
364  printf("***************************************************\n");
365  }
366  }
367  else {
368  printf("\n");
369  printf("------ TESTS FOR PLASMA SGELQF + SGELQS ROUTINE ------- \n");
370  printf(" Size of the Matrix %d by %d\n", M, N);
371  printf("\n");
372  printf(" The matrix A is randomly generated for each test.\n");
373  printf("============\n");
374  printf(" The relative machine precision (eps) is to be %e \n", eps);
375  printf(" Computational tests pass if scaled residuals are less than 60.\n");
376 
377  /* Plasma routines */
378  PLASMA_sgelqf(M, N, A2, LDA, T);
379  PLASMA_sorglq(M, N, K, A2, LDA, T, Q, LDA);
380  PLASMA_sgelqs(M, N, NRHS, A2, LDA, T, B2, LDB);
381 
382  /* Check the orthogonality, factorization and the solution */
383  info_ortho = check_orthogonality(M, N, LDA, Q, eps);
384  info_factorization = check_factorization(M, N, A1, A2, LDA, Q, eps);
385  info_solution = check_solution(M, N, NRHS, A1, LDA, B1, B2, LDB, eps);
386 
387  if ( (info_solution == 0) & (info_factorization == 0) & (info_ortho == 0) ) {
388  printf("***************************************************\n");
389  printf(" ---- TESTING SGELQF + SGELQS ............ PASSED !\n");
390  printf("***************************************************\n");
391  }
392  else {
393  printf("***************************************************\n");
394  printf(" - TESTING SGELQF + SGELQS ... FAILED !\n");
395  printf("***************************************************\n");
396  }
397  }
398 
399  /*----------------------------------------------------------
400  * TESTING SGEQRF + ZORMQR + STRSM
401  */
402 
403  /* Initialize A1 and A2 */
404  LAPACKE_slarnv_work(IONE, ISEED, LDAxN, A1);
405  for (i = 0; i < M; i++)
406  for (j = 0; j < N; j++)
407  A2[LDA*j+i] = A1[LDA*j+i];
408 
409  /* Initialize B1 and B2 */
410  memset(B2, 0, LDB*NRHS*sizeof(float));
411  LAPACKE_slarnv_work(IONE, ISEED, LDBxNRHS, B1);
412  for (i = 0; i < M; i++)
413  for (j = 0; j < NRHS; j++)
414  B2[LDB*j+i] = B1[LDB*j+i];
415 
416  /* PLASMA SGEQRF+ SORMQR + STRSM */
417  memset((void*)Q, 0, LDA*N*sizeof(float));
418  for (i = 0; i < K; i++)
419  Q[LDA*i+i] = 1.0;
420 
421  if (M >= N) {
422  printf("\n");
423  printf("------ TESTS FOR PLASMA SGEQRF + SORMQR + STRSM ROUTINE ------- \n");
424  printf(" Size of the Matrix %d by %d\n", M, N);
425  printf("\n");
426  printf(" The matrix A is randomly generated for each test.\n");
427  printf("============\n");
428  printf(" The relative machine precision (eps) is to be %e \n",eps);
429  printf(" Computational tests pass if scaled residuals are less than 60.\n");
430 
431  PLASMA_sgeqrf(M, N, A2, LDA, T);
432  PLASMA_sorgqr(M, N, K, A2, LDA, T, Q, LDA);
433  PLASMA_sormqr(PlasmaLeft, PlasmaTrans, M, NRHS, N, A2, LDA, T, B2, LDB);
434  PLASMA_strsm(PlasmaLeft, PlasmaUpper, PlasmaNoTrans, PlasmaNonUnit, N, NRHS, 1.0, A2, LDA, B2, LDB);
435  }
436  else {
437  printf("\n");
438  printf("------ TESTS FOR PLASMA SGELQF + SORMLQ + STRSM ROUTINE ------- \n");
439  printf(" Size of the Matrix %d by %d\n", M, N);
440  printf("\n");
441  printf(" The matrix A is randomly generated for each test.\n");
442  printf("============\n");
443  printf(" The relative machine precision (eps) is to be %e \n",eps);
444  printf(" Computational tests pass if scaled residuals are less than 60.\n");
445 
446  PLASMA_sgelqf(M, N, A2, LDA, T);
447  PLASMA_strsm(PlasmaLeft, PlasmaLower, PlasmaNoTrans, PlasmaNonUnit, M, NRHS, 1.0, A2, LDA, B2, LDB);
448  PLASMA_sorglq(M, N, K, A2, LDA, T, Q, LDA);
449  PLASMA_sormlq(PlasmaLeft, PlasmaTrans, N, NRHS, M, A2, LDA, T, B2, LDB);
450  }
451 
452  /* Check the orthogonality, factorization and the solution */
453  info_ortho = check_orthogonality(M, N, LDA, Q, eps);
454  info_factorization = check_factorization(M, N, A1, A2, LDA, Q, eps);
455  info_solution = check_solution(M, N, NRHS, A1, LDA, B1, B2, LDB, eps);
456 
457  if ( (info_solution == 0) & (info_factorization == 0) & (info_ortho == 0) ) {
458  if (M >= N) {
459  printf("***************************************************\n");
460  printf(" ---- TESTING SGEQRF + SORMQR + STRSM .... PASSED !\n");
461  printf("***************************************************\n");
462  }
463  else {
464  printf("***************************************************\n");
465  printf(" ---- TESTING SGELQF + STRSM + SORMLQ .... PASSED !\n");
466  printf("***************************************************\n");
467  }
468  }
469  else {
470  if (M >= N) {
471  printf("***************************************************\n");
472  printf(" - TESTING SGEQRF + SORMQR + STRSM ... FAILED !\n");
473  printf("***************************************************\n");
474  }
475  else {
476  printf("***************************************************\n");
477  printf(" - TESTING SGELQF + STRSM + SORMLQ ... FAILED !\n");
478  printf("***************************************************\n");
479  }
480  }
481 
482  free(A1); free(A2); free(B1); free(B2); free(Q); free(T);
483 
484  return 0;
485 }
486 
487 /*-------------------------------------------------------------------
488  * Check the orthogonality of Q
489  */
490 
491 static int check_orthogonality(int M, int N, int LDQ, float *Q, float eps)
492 {
493  float alpha, beta;
494  float normQ;
495  int info_ortho;
496  int i;
497  int minMN = min(M, N);
498 
499  float *work = (float *)malloc(minMN*sizeof(float));
500 
501  alpha = 1.0;
502  beta = -1.0;
503 
504  /* Build the idendity matrix USE DLASET?*/
505  float *Id = (float *) malloc(minMN*minMN*sizeof(float));
506  memset((void*)Id, 0, minMN*minMN*sizeof(float));
507  for (i = 0; i < minMN; i++)
508  Id[i*minMN+i] = (float)1.0;
509 
510  /* Perform Id - Q'Q */
511  if (M >= N)
512  cblas_ssyrk(CblasColMajor, CblasUpper, CblasTrans, N, M, alpha, Q, LDQ, beta, Id, N);
513  else
514  cblas_ssyrk(CblasColMajor, CblasUpper, CblasNoTrans, M, N, alpha, Q, LDQ, beta, Id, M);
515 
516  BLAS_ssy_norm( blas_colmajor, blas_inf_norm, blas_upper, minMN, Id, minMN, &normQ );
517 
518  printf("============\n");
519  printf("Checking the orthogonality of Q \n");
520  printf("||Id-Q'*Q||_oo / (N*eps) = %e \n", normQ/(minMN*eps));
521 
522  if ( isnan(normQ / (minMN * eps)) || isinf(normQ / (minMN * eps)) || (normQ / (minMN * eps) > 60.0) ) {
523  printf("-- Orthogonality is suspicious ! \n");
524  info_ortho=1;
525  }
526  else {
527  printf("-- Orthogonality is CORRECT ! \n");
528  info_ortho=0;
529  }
530 
531  free(work); free(Id);
532 
533  return info_ortho;
534 }
535 
536 /*------------------------------------------------------------
537  * Check the factorization QR
538  */
539 
540 static int check_factorization(int M, int N, float *A1, float *A2, int LDA, float *Q, float eps )
541 {
542  float Anorm, Rnorm;
543  float alpha, beta;
544  int info_factorization;
545  int i,j;
546 
547  float *Ql = (float *)malloc(M*N*sizeof(float));
548  float *Residual = (float *)malloc(M*N*sizeof(float));
549  float *work = (float *)malloc(max(M,N)*sizeof(float));
550 
551  alpha=1.0;
552  beta=0.0;
553 
554  if (M >= N) {
555  /* Extract the R */
556  float *R = (float *)malloc(N*N*sizeof(float));
557  memset((void*)R, 0, N*N*sizeof(float));
558  LAPACKE_slacpy_work(LAPACK_COL_MAJOR,'u', M, N, A2, LDA, R, N);
559 
560  /* Perform Ql=Q*R */
561  memset((void*)Ql, 0, M*N*sizeof(float));
562  cblas_sgemm(CblasColMajor, CblasNoTrans, CblasNoTrans, M, N, N, (alpha), Q, LDA, R, N, (beta), Ql, M);
563  free(R);
564  }
565  else {
566  /* Extract the L */
567  float *L = (float *)malloc(M*M*sizeof(float));
568  memset((void*)L, 0, M*M*sizeof(float));
569  LAPACKE_slacpy_work(LAPACK_COL_MAJOR,'l', M, N, A2, LDA, L, M);
570 
571  /* Perform Ql=LQ */
572  memset((void*)Ql, 0, M*N*sizeof(float));
573  cblas_sgemm(CblasColMajor, CblasNoTrans, CblasNoTrans, M, N, M, (alpha), L, M, Q, LDA, (beta), Ql, M);
574  free(L);
575  }
576 
577  /* Compute the Residual */
578  for (i = 0; i < M; i++)
579  for (j = 0 ; j < N; j++)
580  Residual[j*M+i] = A1[j*LDA+i]-Ql[j*M+i];
581 
582  BLAS_sge_norm( blas_colmajor, blas_inf_norm, M, N, Residual, M, &Rnorm );
583  BLAS_sge_norm( blas_colmajor, blas_inf_norm, M, N, A2, LDA, &Anorm );
584 
585  if (M >= N) {
586  printf("============\n");
587  printf("Checking the QR Factorization \n");
588  printf("-- ||A-QR||_oo/(||A||_oo.N.eps) = %e \n",Rnorm/(Anorm*N*eps));
589  }
590  else {
591  printf("============\n");
592  printf("Checking the LQ Factorization \n");
593  printf("-- ||A-LQ||_oo/(||A||_oo.N.eps) = %e \n",Rnorm/(Anorm*N*eps));
594  }
595 
596  if (isnan(Rnorm / (Anorm * N *eps)) || isinf(Rnorm / (Anorm * N *eps)) || (Rnorm / (Anorm * N * eps) > 60.0) ) {
597  printf("-- Factorization is suspicious ! \n");
598  info_factorization = 1;
599  }
600  else {
601  printf("-- Factorization is CORRECT ! \n");
602  info_factorization = 0;
603  }
604 
605  free(work); free(Ql); free(Residual);
606 
607  return info_factorization;
608 }
609 
610 /*--------------------------------------------------------------
611  * Check the solution
612  */
613 
614 static int check_solution(int M, int N, int NRHS, float *A1, int LDA, float *B1, float *B2, int LDB, float eps)
615 {
616  int info_solution;
617  float Rnorm, Anorm, Xnorm, Bnorm;
618  float alpha, beta;
619  float result;
620  float *work = (float *)malloc(max(M, N)* sizeof(float));
621 
622  alpha = 1.0;
623  beta = -1.0;
624 
625  BLAS_sge_norm( blas_colmajor, blas_inf_norm, M, N, A1, LDA, &Anorm );
626  BLAS_sge_norm( blas_colmajor, blas_inf_norm, M, NRHS, B2, LDB, &Xnorm );
627  BLAS_sge_norm( blas_colmajor, blas_inf_norm, N, NRHS, B1, LDB, &Bnorm );
628 
629  cblas_sgemm(CblasColMajor, CblasNoTrans, CblasNoTrans, M, NRHS, N, (alpha), A1, LDA, B2, LDB, (beta), B1, LDB);
630 
631  if (M >= N) {
632  float *Residual = (float *)malloc(M*NRHS*sizeof(float));
633  memset((void*)Residual, 0, M*NRHS*sizeof(float));
634  cblas_sgemm(CblasColMajor, CblasTrans, CblasNoTrans, N, NRHS, M, (alpha), A1, LDA, B1, LDB, (beta), Residual, M);
635  BLAS_sge_norm( blas_colmajor, blas_inf_norm, M, NRHS, Residual, M, &Rnorm );
636  free(Residual);
637  }
638  else {
639  float *Residual = (float *)malloc(N*NRHS*sizeof(float));
640  memset((void*)Residual, 0, N*NRHS*sizeof(float));
641  cblas_sgemm(CblasColMajor, CblasTrans, CblasNoTrans, N, NRHS, M, (alpha), A1, LDA, B1, LDB, (beta), Residual, N);
642  BLAS_sge_norm( blas_colmajor, blas_inf_norm, N, NRHS, Residual, N, &Rnorm );
643  free(Residual);
644  }
645 
646  if (getenv("PLASMA_TESTING_VERBOSE"))
647  printf( "||A||_oo=%f\n||X||_oo=%f\n||B||_oo=%f\n||A X - B||_oo=%e\n", Anorm, Xnorm, Bnorm, Rnorm );
648 
649  result = Rnorm / ( (Anorm*Xnorm+Bnorm)*N*eps ) ;
650  printf("============\n");
651  printf("Checking the Residual of the solution \n");
652  printf("-- ||Ax-B||_oo/((||A||_oo||x||_oo+||B||_oo).N.eps) = %e \n", result);
653 
654  if ( isnan(Xnorm) || isinf(Xnorm) || isnan(result) || isinf(result) || (result > 60.0) ) {
655  printf("-- The solution is suspicious ! \n");
656  info_solution = 1;
657  }
658  else{
659  printf("-- The solution is CORRECT ! \n");
660  info_solution = 0;
661  }
662  free(work);
663  return info_solution;
664 }