xref: /petsc/src/mat/tests/ex69.c (revision bbfc976bf0ef01cd099a3849dcc568d2489efeb0)
1*bbfc976bSStefano Zampini static char help[] = "Tests MatCreateDenseCUDA(), MatDenseCUDAPlaceArray(), MatDenseCUDAReplaceArray(), MatDenseCUDAResetArray()\n";
2*bbfc976bSStefano Zampini 
3*bbfc976bSStefano Zampini #include <petscmat.h>
4*bbfc976bSStefano Zampini 
5*bbfc976bSStefano Zampini static PetscErrorCode MatMult_S(Mat S,Vec x,Vec y)
6*bbfc976bSStefano Zampini {
7*bbfc976bSStefano Zampini   PetscErrorCode ierr;
8*bbfc976bSStefano Zampini   Mat            A;
9*bbfc976bSStefano Zampini 
10*bbfc976bSStefano Zampini   PetscFunctionBegin;
11*bbfc976bSStefano Zampini   ierr = MatShellGetContext(S,(void**)&A);CHKERRQ(ierr);
12*bbfc976bSStefano Zampini   ierr = MatMult(A,x,y);CHKERRQ(ierr);
13*bbfc976bSStefano Zampini   PetscFunctionReturn(0);
14*bbfc976bSStefano Zampini }
15*bbfc976bSStefano Zampini 
16*bbfc976bSStefano Zampini int main(int argc,char **argv)
17*bbfc976bSStefano Zampini {
18*bbfc976bSStefano Zampini   Mat            A,B,C,S;
19*bbfc976bSStefano Zampini   Vec            t,v;
20*bbfc976bSStefano Zampini   PetscScalar    *vv,*aa;
21*bbfc976bSStefano Zampini   PetscInt       n=30,k=6,l=0,i,Istart,Iend,nloc,bs,test=1;
22*bbfc976bSStefano Zampini   PetscBool      flg,reset,use_shell = PETSC_FALSE;
23*bbfc976bSStefano Zampini   PetscErrorCode ierr;
24*bbfc976bSStefano Zampini   VecType        vtype;
25*bbfc976bSStefano Zampini 
26*bbfc976bSStefano Zampini   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
27*bbfc976bSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
28*bbfc976bSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);CHKERRQ(ierr);
29*bbfc976bSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);CHKERRQ(ierr);
30*bbfc976bSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-test",&test,NULL);CHKERRQ(ierr);
31*bbfc976bSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-use_shell",&use_shell,NULL);CHKERRQ(ierr);
32*bbfc976bSStefano Zampini   if (k < 0) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_USER,"k %D must be positive",k);
33*bbfc976bSStefano Zampini   if (l < 0) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_USER,"l %D must be positive",l);
34*bbfc976bSStefano Zampini   if (l > k) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_USER,"l %D must be smaller or equal than k %D\n",l,k);
35*bbfc976bSStefano Zampini 
36*bbfc976bSStefano Zampini   /* sparse matrix */
37*bbfc976bSStefano Zampini   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
38*bbfc976bSStefano Zampini   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr);
39*bbfc976bSStefano Zampini   ierr = MatSetType(A,MATAIJCUSPARSE);CHKERRQ(ierr);
40*bbfc976bSStefano Zampini   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
41*bbfc976bSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
42*bbfc976bSStefano Zampini 
43*bbfc976bSStefano Zampini   ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
44*bbfc976bSStefano Zampini   for (i=Istart;i<Iend;i++) {
45*bbfc976bSStefano Zampini     if (i>0) { ierr = MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
46*bbfc976bSStefano Zampini     if (i<n-1) { ierr = MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
47*bbfc976bSStefano Zampini     ierr = MatSetValue(A,i,i,2.0,INSERT_VALUES);CHKERRQ(ierr);
48*bbfc976bSStefano Zampini   }
49*bbfc976bSStefano Zampini   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50*bbfc976bSStefano Zampini   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
51*bbfc976bSStefano Zampini 
52*bbfc976bSStefano Zampini   /* template vector */
53*bbfc976bSStefano Zampini   ierr = MatCreateVecs(A,NULL,&t);CHKERRQ(ierr);
54*bbfc976bSStefano Zampini   ierr = VecGetType(t,&vtype);CHKERRQ(ierr);
55*bbfc976bSStefano Zampini 
56*bbfc976bSStefano Zampini   /* long vector, contains the stacked columns of an nxk dense matrix */
57*bbfc976bSStefano Zampini   ierr = VecGetLocalSize(t,&nloc);CHKERRQ(ierr);
58*bbfc976bSStefano Zampini   ierr = VecGetBlockSize(t,&bs);CHKERRQ(ierr);
59*bbfc976bSStefano Zampini   ierr = VecCreate(PetscObjectComm((PetscObject)t),&v);CHKERRQ(ierr);
60*bbfc976bSStefano Zampini   ierr = VecSetType(v,vtype);CHKERRQ(ierr);
61*bbfc976bSStefano Zampini   ierr = VecSetSizes(v,k*nloc,k*n);CHKERRQ(ierr);
62*bbfc976bSStefano Zampini   ierr = VecSetBlockSize(v,bs);CHKERRQ(ierr);
63*bbfc976bSStefano Zampini   ierr = VecSetRandom(v,NULL);CHKERRQ(ierr);
64*bbfc976bSStefano Zampini 
65*bbfc976bSStefano Zampini   /* dense matrix that contains the columns of v */
66*bbfc976bSStefano Zampini   ierr = VecCUDAGetArray(v,&vv);CHKERRQ(ierr);
67*bbfc976bSStefano Zampini 
68*bbfc976bSStefano Zampini   /* test few cases for MatDenseCUDA handling pointers */
69*bbfc976bSStefano Zampini   switch (test) {
70*bbfc976bSStefano Zampini   case 1:
71*bbfc976bSStefano Zampini     ierr = MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,vv,&B);CHKERRQ(ierr); /* pass a pointer to avoid allocation of storage */
72*bbfc976bSStefano Zampini     ierr = MatDenseCUDAReplaceArray(B,NULL);CHKERRQ(ierr);  /* replace with a null pointer, the value after BVRestoreMat */
73*bbfc976bSStefano Zampini     ierr = MatDenseCUDAPlaceArray(B,vv+l*nloc);CHKERRQ(ierr);  /* set the actual pointer */
74*bbfc976bSStefano Zampini     reset = PETSC_TRUE;
75*bbfc976bSStefano Zampini     break;
76*bbfc976bSStefano Zampini   case 2:
77*bbfc976bSStefano Zampini     ierr = MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,NULL,&B);CHKERRQ(ierr);
78*bbfc976bSStefano Zampini     ierr = MatDenseCUDAPlaceArray(B,vv+l*nloc);CHKERRQ(ierr);  /* set the actual pointer */
79*bbfc976bSStefano Zampini     reset = PETSC_TRUE;
80*bbfc976bSStefano Zampini     break;
81*bbfc976bSStefano Zampini   default:
82*bbfc976bSStefano Zampini     ierr = MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,vv+l*nloc,&B);CHKERRQ(ierr);
83*bbfc976bSStefano Zampini     reset = PETSC_FALSE;
84*bbfc976bSStefano Zampini     break;
85*bbfc976bSStefano Zampini   }
86*bbfc976bSStefano Zampini   ierr = VecCUDARestoreArray(v,&vv);CHKERRQ(ierr);
87*bbfc976bSStefano Zampini 
88*bbfc976bSStefano Zampini   /* Test MatMatMult */
89*bbfc976bSStefano Zampini   if (use_shell) {
90*bbfc976bSStefano Zampini     ierr = MatCreateShell(PetscObjectComm((PetscObject)v),nloc,nloc,n,n,A,&S);CHKERRQ(ierr);
91*bbfc976bSStefano Zampini     ierr = MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_S);CHKERRQ(ierr);
92*bbfc976bSStefano Zampini     ierr = MatShellSetVecType(S,vtype);CHKERRQ(ierr);
93*bbfc976bSStefano Zampini     /* we could have called the general convertor also */
94*bbfc976bSStefano Zampini     /* ierr = MatConvert(A,MATSHELL,MAT_INITIAL_MATRIX,&S);CHKERRQ(ierr); */
95*bbfc976bSStefano Zampini   } else {
96*bbfc976bSStefano Zampini     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
97*bbfc976bSStefano Zampini     S    = A;
98*bbfc976bSStefano Zampini   }
99*bbfc976bSStefano Zampini 
100*bbfc976bSStefano Zampini   ierr = MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,NULL,&C);CHKERRQ(ierr);
101*bbfc976bSStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
102*bbfc976bSStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
103*bbfc976bSStefano Zampini   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
104*bbfc976bSStefano Zampini   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105*bbfc976bSStefano Zampini   ierr = MatProductCreateWithMat(S,B,NULL,C);CHKERRQ(ierr);
106*bbfc976bSStefano Zampini   ierr = MatProductSetType(C,MATPRODUCT_AB);CHKERRQ(ierr);
107*bbfc976bSStefano Zampini   ierr = MatProductSetFromOptions(C);CHKERRQ(ierr);
108*bbfc976bSStefano Zampini   ierr = MatProductSymbolic(C);CHKERRQ(ierr);
109*bbfc976bSStefano Zampini   ierr = MatProductNumeric(C);CHKERRQ(ierr);
110*bbfc976bSStefano Zampini 
111*bbfc976bSStefano Zampini   /* test MatMatMult */
112*bbfc976bSStefano Zampini   ierr = MatMatMultEqual(S,B,C,10,&flg);CHKERRQ(ierr);
113*bbfc976bSStefano Zampini   if (!flg) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error MatMatMult\n"); }
114*bbfc976bSStefano Zampini 
115*bbfc976bSStefano Zampini   ierr = MatDestroy(&C);CHKERRQ(ierr);
116*bbfc976bSStefano Zampini   ierr = MatDestroy(&S);CHKERRQ(ierr);
117*bbfc976bSStefano Zampini 
118*bbfc976bSStefano Zampini   /* finished using B */
119*bbfc976bSStefano Zampini   ierr = MatDenseCUDAGetArray(B,&aa);CHKERRQ(ierr);
120*bbfc976bSStefano Zampini   if (vv != aa-l*nloc) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong array");
121*bbfc976bSStefano Zampini   ierr = MatDenseCUDARestoreArray(B,&aa);CHKERRQ(ierr);
122*bbfc976bSStefano Zampini   if (reset) {
123*bbfc976bSStefano Zampini     ierr = MatDenseCUDAResetArray(B);CHKERRQ(ierr);
124*bbfc976bSStefano Zampini   }
125*bbfc976bSStefano Zampini   ierr = VecCUDARestoreArray(v,&vv);CHKERRQ(ierr);
126*bbfc976bSStefano Zampini 
127*bbfc976bSStefano Zampini   if (test == 1) {
128*bbfc976bSStefano Zampini     ierr = MatDenseCUDAGetArray(B,&aa);CHKERRQ(ierr);
129*bbfc976bSStefano Zampini     if (aa) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected a null pointer");
130*bbfc976bSStefano Zampini     ierr = MatDenseCUDARestoreArray(B,&aa);CHKERRQ(ierr);
131*bbfc976bSStefano Zampini   }
132*bbfc976bSStefano Zampini 
133*bbfc976bSStefano Zampini   /* free work space */
134*bbfc976bSStefano Zampini   ierr = MatDestroy(&B);CHKERRQ(ierr);
135*bbfc976bSStefano Zampini   ierr = MatDestroy(&A);CHKERRQ(ierr);
136*bbfc976bSStefano Zampini   ierr = VecDestroy(&t);CHKERRQ(ierr);
137*bbfc976bSStefano Zampini   ierr = VecDestroy(&v);CHKERRQ(ierr);
138*bbfc976bSStefano Zampini   ierr = PetscFinalize();
139*bbfc976bSStefano Zampini   return ierr;
140*bbfc976bSStefano Zampini }
141*bbfc976bSStefano Zampini 
142*bbfc976bSStefano Zampini /*TEST
143*bbfc976bSStefano Zampini 
144*bbfc976bSStefano Zampini   build:
145*bbfc976bSStefano Zampini     requires: cuda
146*bbfc976bSStefano Zampini 
147*bbfc976bSStefano Zampini   test:
148*bbfc976bSStefano Zampini     requires: cuda
149*bbfc976bSStefano Zampini     suffix: 1
150*bbfc976bSStefano Zampini     nsize: {{1 2}}
151*bbfc976bSStefano Zampini     args: -test {{0 1 2}} -k 6 -l {{0 5}} -use_shell {{0 1}}
152*bbfc976bSStefano Zampini 
153*bbfc976bSStefano Zampini TEST*/
154