xref: /petsc/src/mat/tests/ex32.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] = "Tests MATSEQDENSECUDA\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **argv)
7 {
8   Mat            A,AC,B;
9   PetscErrorCode ierr;
10   PetscInt       m = 10,n = 10;
11   PetscReal      r,tol = 10*PETSC_SMALL;
12 
13   ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr;
14   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
15   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
16   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
17   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n));
18   CHKERRQ(MatSetType(A,MATSEQDENSE));
19   CHKERRQ(MatSetFromOptions(A));
20   CHKERRQ(MatSeqDenseSetPreallocation(A,NULL));
21   CHKERRQ(MatSetRandom(A,NULL));
22 #if 0
23   PetscInt       i,j;
24   PetscScalar    val;
25   for (i=0; i<m; i++) {
26     for (j=0; j<n; j++) {
27       val = (PetscScalar)(i+j);
28       CHKERRQ(MatSetValues(A,1,&i,1,&j,&val,INSERT_VALUES));
29     }
30   }
31   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
32   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
33 #endif
34 
35   /* Create a CUDA version of A */
36 #if defined(PETSC_HAVE_CUDA)
37   CHKERRQ(MatConvert(A,MATSEQDENSECUDA,MAT_INITIAL_MATRIX,&AC));
38 #else
39   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&AC));
40 #endif
41   CHKERRQ(MatDuplicate(AC,MAT_COPY_VALUES,&B));
42 
43   /* full CUDA AXPY */
44   CHKERRQ(MatAXPY(B,-1.0,AC,SAME_NONZERO_PATTERN));
45   CHKERRQ(MatNorm(B,NORM_INFINITY,&r));
46   PetscCheckFalse(r != 0.0,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatDuplicate + MatCopy + MatAXPY %g",(double)r);
47 
48   /* test Copy */
49   CHKERRQ(MatCopy(AC,B,SAME_NONZERO_PATTERN));
50 
51   /* call MatAXPY_Basic since B is CUDA, A is CPU,  */
52   CHKERRQ(MatAXPY(B,-1.0,A,SAME_NONZERO_PATTERN));
53   CHKERRQ(MatNorm(B,NORM_INFINITY,&r));
54   PetscCheckFalse(r != 0.0,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatDuplicate + MatCopy + MatAXPY_Basic %g",(double)r);
55 
56   if (m == n) {
57     Mat B1,B2;
58 
59     CHKERRQ(MatCopy(AC,B,SAME_NONZERO_PATTERN));
60     /* full CUDA PtAP */
61     CHKERRQ(MatPtAP(B,AC,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B1));
62 
63     /* CPU PtAP since A is on the CPU only */
64     CHKERRQ(MatPtAP(B,A,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B2));
65 
66     CHKERRQ(MatAXPY(B2,-1.0,B1,SAME_NONZERO_PATTERN));
67     CHKERRQ(MatNorm(B2,NORM_INFINITY,&r));
68     PetscCheckFalse(r > tol,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatPtAP %g",(double)r);
69 
70     /* test reuse */
71     CHKERRQ(MatPtAP(B,AC,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B1));
72     CHKERRQ(MatPtAP(B,A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B2));
73     CHKERRQ(MatAXPY(B2,-1.0,B1,SAME_NONZERO_PATTERN));
74     CHKERRQ(MatNorm(B2,NORM_INFINITY,&r));
75     PetscCheckFalse(r > tol,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatPtAP %g",(double)r);
76 
77     CHKERRQ(MatDestroy(&B1));
78     CHKERRQ(MatDestroy(&B2));
79   }
80 
81   CHKERRQ(MatDestroy(&B));
82   CHKERRQ(MatDestroy(&AC));
83   CHKERRQ(MatDestroy(&A));
84   ierr = PetscFinalize();
85   return ierr;
86 }
87 
88 /*TEST
89 
90    build:
91      requires: cuda
92 
93    test:
94      output_file: output/ex32_1.out
95      args: -m {{3 5 12}} -n {{3 5 12}}
96      suffix: seqdensecuda
97 
98 TEST*/
99