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_dpotrs.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  /* Initialize A1 and A2 for Symmetric Positif Matrix (Hessenberg in the complex case) */
59  PLASMA_dplgsy( (double)N, N, A1, LDA, 51 );
60  PLASMA_dlacpy( PlasmaUpperLower, N, N, A1, LDA, A2, LDA );
61 
62  /* Initialize B1 and B2 */
63  PLASMA_dplrnt( N, NRHS, B1, LDB, 371 );
64  PLASMA_dlacpy( PlasmaUpperLower, N, NRHS, B1, LDB, B2, LDB );
65 
66  /* Plasma routines */
67  info = PLASMA_dpotrf(PlasmaUpper, N, A2, LDA);
68  info = PLASMA_dpotrs(PlasmaUpper, N, NRHS, A2, LDA, B2, LDB);
69 
70  /* Check the solution */
71  info_solution = check_solution(N, NRHS, A1, LDA, B1, B2, LDB);
72 
73  if ((info_solution != 0)|(info != 0))
74  printf("-- Error in DPOTRS example ! \n");
75  else
76  printf("-- Run of DPOTRS example successful ! \n");
77 
78  free(A1); free(A2); free(B1); free(B2);
79 
81 
82  return EXIT_SUCCESS;
83 }
84 
85 /*------------------------------------------------------------------------
86  * Check the accuracy of the solution of the linear system
87  */
88 
89 int check_solution(int N, int NRHS, double *A1, int LDA, double *B1, double *B2, int LDB )
90 {
91  int info_solution;
92  double Rnorm, Anorm, Xnorm, Bnorm;
93  double alpha, beta;
94  double *work = (double *)malloc(N*sizeof(double));
95  double eps;
96 
97  eps = LAPACKE_dlamch_work('e');
98 
99  /* Initialize A1 and A2 for Symmetric Positive Matrix */
100  alpha = 1.0;
101  beta = -1.0;
102 
103  Xnorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B2, LDB, work);
104  Anorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, N, A1, LDA, work);
105  Bnorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
106 
107  cblas_dgemm(CblasColMajor, CblasNoTrans, CblasNoTrans, N, NRHS, N, (alpha), A1, LDA, B2, LDB, (beta), B1, LDB);
108  Rnorm = LAPACKE_dlange_work(LAPACK_COL_MAJOR, lapack_const(PlasmaInfNorm), N, NRHS, B1, LDB, work);
109 
110  printf("============\n");
111  printf("Checking the Residual of the solution \n");
112  printf("-- ||Ax-B||_oo/((||A||_oo||x||_oo+||B||_oo).N.eps) = %e \n",Rnorm/((Anorm*Xnorm+Bnorm)*N*eps));
113 
114  if (Rnorm/((Anorm*Xnorm+Bnorm)*N*eps) > 10.0){
115  printf("-- The solution is suspicious ! \n");
116  info_solution = 1;
117  }
118  else{
119  printf("-- The solution is CORRECT ! \n");
120  info_solution = 0;
121  }
122 
123  free(work);
124 
125  return info_solution;
126 }