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
example_ztrsm.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 
29 
30 int IONE=1;
31 int ISEED[4] = {0,0,0,1}; /* initial seed for zlarnv() */
32 
33 int main ()
34 {
35 
36  int cores = 2;
37  int N = 10 ;
38  int LDA = 10 ;
39  int NRHS = 5 ;
40  int LDB = 10 ;
41  int info_solution;
42 
43  PLASMA_Complex64_t *A1 = (PLASMA_Complex64_t *)malloc(LDA*N*sizeof(PLASMA_Complex64_t));
44  PLASMA_Complex64_t *A2 = (PLASMA_Complex64_t *)malloc(LDA*N*sizeof(PLASMA_Complex64_t));
45  PLASMA_Complex64_t *B1 = (PLASMA_Complex64_t *)malloc(LDB*NRHS*sizeof(PLASMA_Complex64_t));
46  PLASMA_Complex64_t *B2 = (PLASMA_Complex64_t *)malloc(LDB*NRHS*sizeof(PLASMA_Complex64_t));
47 
48  /* Check if unable to allocate memory */
49  if ((!A1)||(!A2)||(!B1)||(!B2)){
50  printf("Out of Memory \n ");
51  return EXIT_SUCCESS;
52  }
53 
54  /* Plasma Initialize */
55  PLASMA_Init(cores);
56  printf("-- PLASMA is initialized to run on %d cores. \n",cores);
57 
58  /* Initialize A1 and A2 for Symmetric Positive Matrix */
59  PLASMA_zplghe( (double)N, N, A1, LDA, 51 );
60  PLASMA_zlacpy( PlasmaUpperLower, N, N, A1, LDA, A2, LDA );
61 
62  /* Initialize B1 and B2 */
63  PLASMA_zplrnt( N, NRHS, B1, LDB, 371 );
64  PLASMA_zlacpy( PlasmaUpperLower, N, NRHS, B1, LDB, B2, LDB );
65 
66  /* PLASMA routines */
67  PLASMA_zpotrf(PlasmaLower, N, A2, LDA);
69  N, NRHS, (PLASMA_Complex64_t)1.0, A2, LDA, B2, LDB);
71  N, NRHS, (PLASMA_Complex64_t)1.0, A2, LDA, B2, LDB);
72 
73  /* Check the solution */
74  info_solution = check_solution(N, NRHS, A1, LDA, B1, B2, LDB);
75 
76  if ( info_solution != 0 )
77  printf("-- Error in ZTRSM example ! \n");
78  else
79  printf("-- Run of ZTRSM example successful ! \n");
80 
81  free(A1); free(A2); free(B1); free(B2);
82 
84 
85  return EXIT_SUCCESS;
86 }
87 
88 
89 /*------------------------------------------------------------------------
90  * Check the accuracy of the solution of the linear system
91  */
92 
93 int check_solution(int N, int NRHS, PLASMA_Complex64_t *A1, int LDA, PLASMA_Complex64_t *B1, PLASMA_Complex64_t *B2, int LDB)
94 {
95  int info_solution;
96  double Rnorm, Anorm, Xnorm, Bnorm;
97  PLASMA_Complex64_t alpha, beta;
98  double *work = (double *)malloc(N*sizeof(double));
99  double eps;
100 
101  eps = LAPACKE_dlamch_work('e');
102 
103  alpha = 1.0;
104  beta = -1.0;
105 
106  Xnorm = LAPACKE_zlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B2, LDB, work);
107  Anorm = LAPACKE_zlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, N, A1, LDA, work);
108  Bnorm = LAPACKE_zlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
109 
110  cblas_zgemm(CblasColMajor, CblasNoTrans, CblasNoTrans, N, NRHS, N, CBLAS_SADDR(alpha), A1, LDA, B2, LDB, CBLAS_SADDR(beta), B1, LDB);
111  Rnorm = LAPACKE_zlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
112 
113  printf("============\n");
114  printf("Checking the Residual of the solution \n");
115  printf("-- ||Ax-B||_oo/((||A||_oo||x||_oo+||B||_oo).N.eps) = %e \n",Rnorm/((Anorm*Xnorm+Bnorm)*N*eps));
116 
117  if (Rnorm/((Anorm*Xnorm+Bnorm)*N*eps) > 10.0){
118  printf("-- The solution is suspicious ! \n");
119  info_solution = 1;
120  }
121  else{
122  printf("-- The solution is CORRECT ! \n");
123  info_solution = 0;
124  }
125 
126  free(work);
127 
128  return info_solution;
129 }