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_zgesv_incpiv.c
Go to the documentation of this file.
1 
15 #include <stdlib.h>
16 #include <stdio.h>
17 #include <string.h>
18 #include <math.h>
19 
20 #include <plasma.h>
21 #include <cblas.h>
22 #include <lapacke.h>
23 #include <core_blas.h>
24 #include "testing_zmain.h"
25 
28  blas_colmajor = 102 };
29 
31  blas_base = 151,
32  blas_t = 152,
33  blas_rnd = 153,
34  blas_ieee = 154,
35  blas_emin = 155,
36  blas_emax = 156,
37  blas_eps = 157,
38  blas_prec = 158,
41  blas_sfmin = 161};
42 
52 
53 static void
54 BLAS_error(char *rname, int err, int val, int x) {
55  fprintf( stderr, "%s %d %d %d\n", rname, err, val, x );
56  abort();
57 }
58 
59 static
60 void
61 BLAS_zge_norm(enum blas_order_type order, enum blas_norm_type norm,
62  int m, int n, const PLASMA_Complex64_t *a, int lda, double *res) {
63  int i, j; float anorm, v;
64  char rname[] = "BLAS_zge_norm";
65 
66  if (order != blas_colmajor) BLAS_error( rname, -1, order, 0 );
67 
68  if (norm == blas_frobenius_norm) {
69  anorm = 0.0f;
70  for (j = n; j; --j) {
71  for (i = m; i; --i) {
72  v = a[0];
73  anorm += v * v;
74  a++;
75  }
76  a += lda - m;
77  }
78  anorm = sqrt( anorm );
79  } else if (norm == blas_inf_norm) {
80  anorm = 0.0f;
81  for (i = 0; i < m; ++i) {
82  v = 0.0f;
83  for (j = 0; j < n; ++j) {
84  v += cabs( a[i + j * lda] );
85  }
86  if (v > anorm)
87  anorm = v;
88  }
89  } else {
90  BLAS_error( rname, -2, norm, 0 );
91  return;
92  }
93 
94  if (res) *res = anorm;
95 }
96 
97 static
98 double
99 BLAS_dpow_di(double x, int n) {
100  double rv = 1.0;
101 
102  if (n < 0) {
103  n = -n;
104  x = 1.0 / x;
105  }
106 
107  for (; n; n >>= 1, x *= x) {
108  if (n & 1)
109  rv *= x;
110  }
111 
112  return rv;
113 }
114 
115 static
116 double
117 BLAS_dfpinfo(enum blas_cmach_type cmach) {
118  double eps = 1.0, r = 1.0, o = 1.0, b = 2.0;
119  int t = 53, l = 1024, m = -1021;
120  char rname[] = "BLAS_dfpinfo";
121 
122  if ((sizeof eps) == sizeof(float)) {
123  t = 24;
124  l = 128;
125  m = -125;
126  } else {
127  t = 53;
128  l = 1024;
129  m = -1021;
130  }
131 
132  /* for (i = 0; i < t; ++i) eps *= half; */
133  eps = BLAS_dpow_di( b, -t );
134  /* for (i = 0; i >= m; --i) r *= half; */
135  r = BLAS_dpow_di( b, m-1 );
136 
137  o -= eps;
138  /* for (i = 0; i < l; ++i) o *= b; */
139  o = (o * BLAS_dpow_di( b, l-1 )) * b;
140 
141  switch (cmach) {
142  case blas_eps: return eps;
143  case blas_sfmin: return r;
144  default:
145  BLAS_error( rname, -1, cmach, 0 );
146  break;
147  }
148  return 0.0;
149 }
150 
151 static int check_solution(int, int , PLASMA_Complex64_t *, int, PLASMA_Complex64_t *, PLASMA_Complex64_t *, int, double);
152 
153 int testing_zgesv_incpiv(int argc, char **argv)
154 {
155  /* Check for valid arguments*/
156  if (argc != 4){
157  USAGE("GESV_INCPIV", "N LDA NRHS LDB",
158  " - N : the size of the matrix\n"
159  " - LDA : leading dimension of the matrix A\n"
160  " - NRHS : number of RHS\n"
161  " - LDB : leading dimension of the matrix B\n");
162  return -1;
163  }
164 
165  int N = atoi(argv[0]);
166  int LDA = atoi(argv[1]);
167  int NRHS = atoi(argv[2]);
168  int LDB = atoi(argv[3]);
169  double eps;
170  int info_solution;
171  int i,j;
172  int LDAxN = LDA*N;
173  int LDBxNRHS = LDB*NRHS;
174 
175  PLASMA_Complex64_t *A1 = (PLASMA_Complex64_t *)malloc(LDA*N*(sizeof*A1));
176  PLASMA_Complex64_t *A2 = (PLASMA_Complex64_t *)malloc(LDA*N*(sizeof*A2));
177  PLASMA_Complex64_t *B1 = (PLASMA_Complex64_t *)malloc(LDB*NRHS*(sizeof*B1));
178  PLASMA_Complex64_t *B2 = (PLASMA_Complex64_t *)malloc(LDB*NRHS*(sizeof*B2));
180  int *IPIV;
181 
182  /* Check if unable to allocate memory */
183  if ( (!A1) || (!A2)|| (!B1) || (!B2) ) {
184  printf("Out of Memory \n ");
185  return -2;
186  }
187 
188  eps = BLAS_dfpinfo(blas_eps);
189 
190  /*----------------------------------------------------------
191  * TESTING ZGESV
192  */
193 
194  /* Initialize A1 and A2 Matrix */
195  LAPACKE_zlarnv_work(IONE, ISEED, LDAxN, A1);
196  for ( i = 0; i < N; i++)
197  for ( j = 0; j < N; j++)
198  A2[LDA*j+i] = A1[LDA*j+i];
199 
200  /* Initialize B1 and B2 */
201  LAPACKE_zlarnv_work(IONE, ISEED, LDBxNRHS, B1);
202  for ( i = 0; i < N; i++)
203  for ( j = 0; j < NRHS; j++)
204  B2[LDB*j+i] = B1[LDB*j+i];
205 
206  /* PLASMA ZGESV */
208  PLASMA_zgesv_incpiv(N, NRHS, A2, LDA, L, IPIV, B2, LDB);
209 
210  printf("\n");
211  printf("------ TESTS FOR PLASMA INCPIV ZGESV ROUTINE ------- \n");
212  printf(" Size of the Matrix %d by %d\n", N, N);
213  printf("\n");
214  printf(" The matrix A is randomly generated for each test.\n");
215  printf("============\n");
216  printf(" The relative machine precision (eps) is to be %e \n", eps);
217  printf(" Computational tests pass if scaled residuals are less than 60.\n");
218 
219  /* Check the factorization and the solution */
220  info_solution = check_solution(N, NRHS, A1, LDA, B1, B2, LDB, eps);
221 
222  if ((info_solution == 0)){
223  printf("***************************************************\n");
224  printf(" ---- TESTING INCPIV ZGESV ............... PASSED !\n");
225  printf("***************************************************\n");
226  }
227  else{
228  printf("************************************************\n");
229  printf(" - TESTING INCPIV ZGESV ... FAILED !\n");
230  printf("************************************************\n");
231  }
232 
233  /*-------------------------------------------------------------
234  * TESTING ZGETRF + ZGETRS
235  */
236 
237  /* Initialize A1 and A2 */
238  LAPACKE_zlarnv_work(IONE, ISEED, LDAxN, A1);
239  for ( i = 0; i < N; i++)
240  for ( j = 0; j < N; j++)
241  A2[LDA*j+i] = A1[LDA*j+i];
242 
243  /* Initialize B1 and B2 */
244  LAPACKE_zlarnv_work(IONE, ISEED, LDBxNRHS, B1);
245  for ( i = 0; i < N; i++)
246  for ( j = 0; j < NRHS; j++)
247  B2[LDB*j+i] = B1[LDB*j+i];
248 
249  /* Plasma routines */
250  PLASMA_zgetrf_incpiv(N, N, A2, LDA, L, IPIV);
251  PLASMA_zgetrs_incpiv(PlasmaNoTrans, N, NRHS, A2, LDA, L, IPIV, B2, LDB);
252 
253  printf("\n");
254  printf("------ TESTS FOR PLASMA ZGETRF + ZGETRS ROUTINE ------- \n");
255  printf(" Size of the Matrix %d by %d\n", N, N);
256  printf("\n");
257  printf(" The matrix A is randomly generated for each test.\n");
258  printf("============\n");
259  printf(" The relative machine precision (eps) is to be %e \n", eps);
260  printf(" Computational tests pass if scaled residuals are less than 60.\n");
261 
262  /* Check the solution */
263  info_solution = check_solution(N, NRHS, A1, LDA, B1, B2, LDB, eps);
264 
265  if ((info_solution == 0)){
266  printf("***************************************************\n");
267  printf(" ---- TESTING INCPIV ZGETRF + ZGETRS ..... PASSED !\n");
268  printf("***************************************************\n");
269  }
270  else{
271  printf("***************************************************\n");
272  printf(" - TESTING INCPIV ZGETRF + ZGETRS ... FAILED !\n");
273  printf("***************************************************\n");
274  }
275 
276  /*-------------------------------------------------------------
277  * TESTING ZGETRF + ZTRSMPL + ZTRSM
278  */
279 
280  /* Initialize A1 and A2 */
281  LAPACKE_zlarnv_work(IONE, ISEED, LDAxN, A1);
282  for ( i = 0; i < N; i++)
283  for ( j = 0; j < N; j++)
284  A2[LDA*j+i] = A1[LDA*j+i];
285 
286  /* Initialize B1 and B2 */
287  LAPACKE_zlarnv_work(IONE, ISEED, LDBxNRHS, B1);
288  for ( i = 0; i < N; i++)
289  for ( j = 0; j < NRHS; j++)
290  B2[LDB*j+i] = B1[LDB*j+i];
291 
292  /* PLASMA routines */
293  PLASMA_zgetrf_incpiv(N, N, A2, LDA, L, IPIV);
294  PLASMA_ztrsmpl(N, NRHS, A2, LDA, L, IPIV, B2, LDB);
296  N, NRHS, 1.0, A2, LDA, B2, LDB);
297 
298  printf("\n");
299  printf("------ TESTS FOR PLASMA INCPIV ZGETRF + ZTRSMPL + ZTRSM ROUTINE ------- \n");
300  printf(" Size of the Matrix %d by %d\n", N, N);
301  printf("\n");
302  printf(" The matrix A is randomly generated for each test.\n");
303  printf("============\n");
304  printf(" The relative machine precision (eps) is to be %e \n", eps);
305  printf(" Computational tests pass if scaled residuals are less than 60.\n");
306 
307  /* Check the solution */
308  info_solution = check_solution(N, NRHS, A1, LDA, B1, B2, LDB, eps);
309 
310  if ((info_solution == 0)){
311  printf("***************************************************\n");
312  printf(" ---- TESTING INCPIV ZGETRF + ZTRSMPL + ZTRSM ... PASSED !\n");
313  printf("***************************************************\n");
314  }
315  else{
316  printf("**************************************************\n");
317  printf(" - TESTING INCPIV ZGETRF + ZTRSMPL + ZTRSM ... FAILED !\n");
318  printf("**************************************************\n");
319  }
320 
321  free(A1); free(A2); free(B1); free(B2); free(IPIV); free(L);
322 
323  return 0;
324 }
325 
326 /*------------------------------------------------------------------------
327  * Check the accuracy of the solution of the linear system
328  */
329 
330 static int check_solution(int N, int NRHS, PLASMA_Complex64_t *A1, int LDA, PLASMA_Complex64_t *B1, PLASMA_Complex64_t *B2, int LDB, double eps )
331 {
332  int info_solution;
333  double Rnorm, Anorm, Xnorm, Bnorm, result;
334  PLASMA_Complex64_t alpha, beta;
335  double *work = (double *)malloc(N*sizeof(double));
336 
337  alpha = 1.0;
338  beta = -1.0;
339 
340  BLAS_zge_norm( blas_colmajor, blas_inf_norm, N, NRHS, B2, LDB, &Xnorm );
341  BLAS_zge_norm( blas_colmajor, blas_inf_norm, N, N, A1, LDA, &Anorm );
342  BLAS_zge_norm( blas_colmajor, blas_inf_norm, N, NRHS, B1, LDB, &Bnorm );
343 
344  cblas_zgemm(CblasColMajor, CblasNoTrans, CblasNoTrans, N, NRHS, N, CBLAS_SADDR(alpha), A1, LDA, B2, LDB, CBLAS_SADDR(beta), B1, LDB);
345  BLAS_zge_norm( blas_colmajor, blas_inf_norm, N, NRHS, B1, LDB, &Rnorm );
346 
347  if (getenv("PLASMA_TESTING_VERBOSE"))
348  printf( "||A||_oo=%f\n||X||_oo=%f\n||B||_oo=%f\n||A X - B||_oo=%e\n", Anorm, Xnorm, Bnorm, Rnorm );
349 
350  result = Rnorm / ( (Anorm*Xnorm+Bnorm)*N*eps ) ;
351  printf("============\n");
352  printf("Checking the Residual of the solution \n");
353  printf("-- ||Ax-B||_oo/((||A||_oo||x||_oo+||B||_oo).N.eps) = %e \n", result);
354 
355  if ( isnan(Xnorm) || isinf(Xnorm) || isnan(result) || isinf(result) || (result > 60.0) ) {
356  printf("-- The solution is suspicious ! \n");
357  info_solution = 1;
358  }
359  else{
360  printf("-- The solution is CORRECT ! \n");
361  info_solution = 0;
362  }
363  free(work);
364  return info_solution;
365 }