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_dposv.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 
27 int check_solution(int, int, double*, int, double*, double*, int);
28 
29 int IONE=1;
30 int ISEED[4] = {0,0,0,1}; /* initial seed for dlarnv() */
31 
32 int main ()
33 {
34 
35  int cores = 2;
36  int N = 10;
37  int LDA = 10;
38  int NRHS = 5;
39  int LDB = 10;
40  int info;
41  int info_solution;
42 
43  double *A1 = (double *)malloc(LDA*N*sizeof(double));
44  double *A2 = (double *)malloc(LDA*N*sizeof(double));
45  double *B1 = (double *)malloc(LDB*NRHS*sizeof(double));
46  double *B2 = (double *)malloc(LDB*NRHS*sizeof(double));
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  /*-------------------------------------------------------------
59  * TESTING DPOSV
60  */
61 
62  /* Initialize A1 and A2 for Symmetric Positif Matrix (Hessenberg in the complex case) */
63  PLASMA_dplgsy( (double)N, N, A1, LDA, 51 );
64  PLASMA_dlacpy( PlasmaUpperLower, N, N, A1, LDA, A2, LDA );
65 
66  /* Initialize B1 and B2 */
67  PLASMA_dplrnt( N, NRHS, B1, LDB, 371 );
68  PLASMA_dlacpy( PlasmaUpperLower, N, NRHS, B1, LDB, B2, LDB );
69 
70  /* PLASMA DPOSV */
71  info = PLASMA_dposv(PlasmaUpper, N, NRHS, 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)|(info != 0))
77  printf("-- Error in DPOSV example ! \n");
78  else
79  printf("-- Run of DPOSV example successful ! \n");
80 
81  free(A1); free(A2); free(B1); free(B2);
82 
84 
85  return EXIT_SUCCESS;
86 }
87 
88 /*------------------------------------------------------------------------
89  * Check the accuracy of the solution of the linear system
90  */
91 
92 int check_solution(int N, int NRHS, double *A1, int LDA, double *B1, double *B2, int LDB)
93 {
94  int info_solution;
95  double Rnorm, Anorm, Xnorm, Bnorm;
96  double alpha, beta;
97  double *work = (double *)malloc(N*sizeof(double));
98  double eps;
99 
100  eps = LAPACKE_dlamch_work('e');
101 
102  alpha = 1.0;
103  beta = -1.0;
104 
105  Xnorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B2, LDB, work);
106  Anorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, N, A1, LDA, work);
107  Bnorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
108 
109  cblas_dgemm(CblasColMajor, CblasNoTrans, CblasNoTrans, N, NRHS, N, (alpha), A1, LDA, B2, LDB, (beta), B1, LDB);
110  Rnorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
111 
112  printf("============\n");
113  printf("Checking the Residual of the solution \n");
114  printf("-- ||Ax-B||_oo/((||A||_oo||x||_oo+||B||_oo).N.eps) = %e \n",Rnorm/((Anorm*Xnorm+Bnorm)*N*eps));
115 
116  if (Rnorm/((Anorm*Xnorm+Bnorm)*N*eps) > 10.0){
117  printf("-- The solution is suspicious ! \n");
118  info_solution = 1;
119  }
120  else{
121  printf("-- The solution is CORRECT ! \n");
122  info_solution = 0;
123  }
124 
125  free(work);
126 
127  return info_solution;
128 }