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_strsm.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_smain.h"
25 
26 #undef COMPLEX
27 #define REAL
28 
30  int M, int N, float alpha,
31  float *A, int LDA,
32  float *Bref, float *Bplasma, int LDB);
33 
34 int testing_strsm(int argc, char **argv)
35 {
36  /* Check for number of arguments*/
37  if ( argc != 5 ) {
38  USAGE("TRSM", "alpha M N LDA LDB",
39  " - alpha : alpha coefficient\n"
40  " - M : number of rows of matrices B\n"
41  " - N : number of columns of matrices B\n"
42  " - LDA : leading dimension of matrix A\n"
43  " - LDB : leading dimension of matrix B\n");
44  return -1;
45  }
46 
47  float alpha = (float) atol(argv[0]);
48  int M = atoi(argv[1]);
49  int N = atoi(argv[2]);
50  int LDA = atoi(argv[3]);
51  int LDB = atoi(argv[4]);
52 
53  float eps;
54  int info_solution;
55  int s, u, t, d, i;
56  int LDAxM = LDA*max(M,N);
57  int LDBxN = LDB*max(M,N);
58 
59  float *A = (float *)malloc(LDAxM*sizeof(float));
60  float *B = (float *)malloc(LDBxN*sizeof(float));
61  float *Binit = (float *)malloc(LDBxN*sizeof(float));
62  float *Bfinal = (float *)malloc(LDBxN*sizeof(float));
63 
64  /* Check if unable to allocate memory */
65  if ( (!A) || (!B) || (!Binit) || (!Bfinal)){
66  printf("Out of Memory \n ");
67  return -2;
68  }
69 
70  eps = LAPACKE_slamch_work('e');
71 
72  printf("\n");
73  printf("------ TESTS FOR PLASMA STRSM ROUTINE ------- \n");
74  printf(" Size of the Matrix B : %d by %d\n", M, N);
75  printf("\n");
76  printf(" The matrix A is randomly generated for each test.\n");
77  printf("============\n");
78  printf(" The relative machine precision (eps) is to be %e \n",eps);
79  printf(" Computational tests pass if scaled residuals are less than 10.\n");
80 
81  /*----------------------------------------------------------
82  * TESTING STRSM
83  */
84 
85  /* Initialize A, B, C */
86  LAPACKE_slarnv_work(IONE, ISEED, LDAxM, A);
87  LAPACKE_slarnv_work(IONE, ISEED, LDBxN, B);
88  for(i=0;i<max(M,N);i++)
89  A[LDA*i+i] = A[LDA*i+i] + 2.0;
90 
91  for (s=0; s<2; s++) {
92  for (u=0; u<2; u++) {
93 #ifdef COMPLEX
94  for (t=0; t<3; t++) {
95 #else
96  for (t=0; t<2; t++) {
97 #endif
98  for (d=0; d<2; d++) {
99 
100  memcpy(Binit, B, LDBxN*sizeof(float));
101  memcpy(Bfinal, B, LDBxN*sizeof(float));
102 
103  /* PLASMA STRSM */
104  PLASMA_strsm(side[s], uplo[u], trans[t], diag[d],
105  M, N, alpha, A, LDA, Bfinal, LDB);
106 
107  /* Check the solution */
108  info_solution = check_solution(side[s], uplo[u], trans[t], diag[d],
109  M, N, alpha, A, LDA, Binit, Bfinal, LDB);
110 
111  printf("***************************************************\n");
112  if (info_solution == 0) {
113  printf(" ---- TESTING STRSM (%s, %s, %s, %s) ...... PASSED !\n",
114  sidestr[s], uplostr[u], transstr[t], diagstr[d]);
115  }
116  else {
117  printf(" ---- TESTING STRSM (%s, %s, %s, %s) ... FAILED !\n",
118  sidestr[s], uplostr[u], transstr[t], diagstr[d]);
119  }
120  printf("***************************************************\n");
121  }
122  }
123  }
124  }
125 
126  free(A); free(B);
127  free(Binit); free(Bfinal);
128 
129  return 0;
130 }
131 
132 /*--------------------------------------------------------------
133  * Check the solution
134  */
136  int M, int N, float alpha,
137  float *A, int LDA,
138  float *Bref, float *Bplasma, int LDB)
139 {
140  int info_solution;
141  float Anorm, Binitnorm, Bplasmanorm, Blapacknorm, Rnorm, result;
142  float eps;
143  float mzone = (float)-1.0;
144 
145  float *work = (float *)malloc(max(M, N)* sizeof(float));
146  int Am, An;
147 
148  if (side == PlasmaLeft) {
149  Am = M; An = M;
150  } else {
151  Am = N; An = N;
152  }
153 
154  Anorm = LAPACKE_slantr_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), lapack_const(uplo), lapack_const(diag),
155  Am, An, A, LDA, work);
156  Binitnorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), M, N, Bref, LDB, work);
157  Bplasmanorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), M, N, Bplasma, LDB, work);
158 
160  (CBLAS_DIAG)diag, M, N, (alpha), A, LDA, Bref, LDB);
161 
162  Blapacknorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), M, N, Bref, LDB, work);
163 
164  cblas_saxpy(LDB * N, (mzone), Bplasma, 1, Bref, 1);
165 
166  Rnorm = LAPACKE_slange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), M, N, Bref, LDB, work);
167 
168  eps = LAPACKE_slamch_work('e');
169 
170  printf("Rnorm %e, Anorm %e, Binitnorm %e, Bplasmanorm %e, Blapacknorm %e\n",
171  Rnorm, Anorm, Binitnorm, Bplasmanorm, Blapacknorm);
172 
173  result = Rnorm / ((Anorm + Blapacknorm) * max(M,N) * eps);
174 
175  printf("============\n");
176  printf("Checking the norm of the difference against reference STRSM \n");
177  printf("-- ||Cplasma - Clapack||_oo/((||A||_oo+||B||_oo).N.eps) = %e \n", result);
178 
179  if ( isinf(Blapacknorm) || isinf(Bplasmanorm) || isnan(result) || isinf(result) || (result > 10.0) ) {
180  printf("-- The solution is suspicious ! \n");
181  info_solution = 1;
182  }
183  else {
184  printf("-- The solution is CORRECT ! \n");
185  info_solution= 0 ;
186  }
187  free(work);
188 
189  return info_solution;
190 }