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_dsgesv.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 #include "testing_zmain.h"
27 
28 static int check_solution(int, int , double *, int, double *, double *, int, double);
29 
30 int testing_dsgesv(int argc, char **argv)
31 {
32  /* Check for valid arguments*/
33  if (argc != 4){
34  USAGE("CGESV", "N LDA NRHS LDB",
35  " - N : the size of the matrix\n"
36  " - LDA : leading dimension of the matrix A\n"
37  " - NRHS : number of RHS\n"
38  " - LDB : leading dimension of the matrix B\n");
39  return -1;
40  }
41 
42  int N = atoi(argv[0]);
43  int LDA = atoi(argv[1]);
44  int NRHS = atoi(argv[2]);
45  int LDB = atoi(argv[3]);
46  int ITER;
47  double eps;
48  int info, info_solution;
49  int i,j;
50  int LDAxN = LDA*N;
51  int LDBxNRHS = LDB*NRHS;
52 
53  double *A1 = (double *)malloc(LDA*N*sizeof(double));
54  double *A2 = (double *)malloc(LDA*N*sizeof(double));
55  double *B1 = (double *)malloc(LDB*NRHS*sizeof(double));
56  double *B2 = (double *)malloc(LDB*NRHS*sizeof(double));
57  int *ipiv = (int *)malloc(N*sizeof(int));
58 
59  /* Check if unable to allocate memory */
60  if ( (!A1) || (!A2) || (!B1) || (!B2) ) {
61  printf("Out of Memory \n ");
62  exit(0);
63  }
64 
65  eps = LAPACKE_dlamch_work('e');
66 
67  /*----------------------------------------------------------
68  * TESTING DSGESV
69  */
70 
71  /* Initialize A1 and A2 Matrix */
72  LAPACKE_dlarnv_work(IONE, ISEED, LDAxN, A1);
73  for ( i = 0; i < N; i++)
74  for ( j = 0; j < N; j++)
75  A2[LDA*j+i] = A1[LDA*j+i];
76 
77  /* Initialize B1 and B2 */
78  LAPACKE_dlarnv_work(IONE, ISEED, LDBxNRHS, B1);
79  for ( i = 0; i < N; i++)
80  for ( j = 0; j < NRHS; j++)
81  B2[LDB*j+i] = B1[LDB*j+i];
82 
83  printf("\n");
84  printf("------ TESTS FOR PLASMA DSGESV ROUTINE ------- \n");
85  printf(" Size of the Matrix %d by %d\n", N, N);
86  printf("\n");
87  printf(" The matrix A is randomly generated for each test.\n");
88  printf("============\n");
89  printf(" The relative machine precision (eps) is to be %e \n", eps);
90  printf(" Computational tests pass if scaled residuals are less than 60.\n");
91 
92  /* PLASMA DSGESV */
93  info = PLASMA_dsgesv(N, NRHS, A2, LDA, ipiv, B1, LDB, B2, LDB, &ITER);
94  if (info != PLASMA_SUCCESS ) {
95  printf("PLASMA_dsgesv is not completed: info = %d\n", info);
96  info_solution = 1;
97  } else {
98  printf(" Solution obtained with %d iterations\n", ITER);
99 
100  /* Check the factorization and the solution */
101  info_solution = check_solution(N, NRHS, A1, LDA, B1, B2, LDB, eps);
102  }
103 
104  if ((info_solution == 0)){
105  printf("***************************************************\n");
106  printf(" ---- TESTING DSGESV ..................... PASSED !\n");
107  printf("***************************************************\n");
108  }
109  else{
110  printf("************************************************\n");
111  printf(" ---- TESTING DSGESV ... FAILED !\n");
112  printf("************************************************\n");
113  }
114 
115  free(A1); free(A2); free(B1); free(B2); free(ipiv);
116 
117  return 0;
118 }
119 
120 /*------------------------------------------------------------------------
121  * Check the accuracy of the solution of the linear system
122  */
123 static int check_solution(int N, int NRHS, double *A1, int LDA, double *B1, double *B2, int LDB, double eps )
124 {
125  int info_solution;
126  double Rnorm, Anorm, Xnorm, Bnorm;
127  double alpha, beta;
128  double *work = (double *)malloc(N*(sizeof *work));
129  double result;
130 
131  alpha = 1.0;
132  beta = -1.0;
133 
134  Xnorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B2, LDB, work);
135  Anorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, N, A1, LDA, work);
136  Bnorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
137 
138  cblas_dgemm(CblasColMajor, CblasNoTrans, CblasNoTrans, N, NRHS, N, (alpha), A1, LDA, B2, LDB, (beta), B1, LDB);
139  Rnorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
140 
141  result = Rnorm/((Anorm*Xnorm+Bnorm)*N*eps);
142 
143  if (getenv("PLASMA_TESTING_VERBOSE"))
144  printf( "||A||_oo=%f\n||X||_oo=%f\n||B||_oo=%f\n||A X - B||_oo=%e\n", Anorm, Xnorm, Bnorm, Rnorm );
145 
146  printf("============\n");
147  printf("Checking the Residual of the solution \n");
148  printf("-- ||Ax-B||_oo/((||A||_oo||x||_oo+||B||_oo).N.eps) = %e \n", result);
149 
150  if ( isnan(Xnorm) || isinf(Xnorm) || isnan(result) || isinf(result) || (result > 60.0) ){
151  printf("-- The solution is suspicious ! \n");
152  info_solution = 1;
153  }
154  else{
155  printf("-- The solution is CORRECT ! \n");
156  info_solution = 0;
157  }
158 
159  free(work);
160 
161  return info_solution;
162 }