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