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