xref: /petsc/src/mat/tests/ex69.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
1bbfc976bSStefano Zampini static char help[] = "Tests MatCreateDenseCUDA(), MatDenseCUDAPlaceArray(), MatDenseCUDAReplaceArray(), MatDenseCUDAResetArray()\n";
2bbfc976bSStefano Zampini 
3bbfc976bSStefano Zampini #include <petscmat.h>
4bbfc976bSStefano Zampini 
5bbfc976bSStefano Zampini static PetscErrorCode MatMult_S(Mat S,Vec x,Vec y)
6bbfc976bSStefano Zampini {
7bbfc976bSStefano Zampini   PetscErrorCode ierr;
8bbfc976bSStefano Zampini   Mat            A;
9bbfc976bSStefano Zampini 
10d9b9dddeSStefano Zampini   PetscFunctionBeginUser;
11bbfc976bSStefano Zampini   ierr = MatShellGetContext(S,(void**)&A);CHKERRQ(ierr);
12bbfc976bSStefano Zampini   ierr = MatMult(A,x,y);CHKERRQ(ierr);
13bbfc976bSStefano Zampini   PetscFunctionReturn(0);
14bbfc976bSStefano Zampini }
15bbfc976bSStefano Zampini 
1654da937aSStefano Zampini static PetscBool test_cusparse_transgen = PETSC_FALSE;
1754da937aSStefano Zampini 
18d9b9dddeSStefano Zampini static PetscErrorCode MatMultTranspose_S(Mat S,Vec x,Vec y)
19d9b9dddeSStefano Zampini {
20d9b9dddeSStefano Zampini   PetscErrorCode ierr;
21d9b9dddeSStefano Zampini   Mat            A;
22d9b9dddeSStefano Zampini 
23d9b9dddeSStefano Zampini   PetscFunctionBeginUser;
24d9b9dddeSStefano Zampini   ierr = MatShellGetContext(S,(void**)&A);CHKERRQ(ierr);
25d9b9dddeSStefano Zampini   ierr = MatMultTranspose(A,x,y);CHKERRQ(ierr);
2654da937aSStefano Zampini 
2754da937aSStefano Zampini   /* alternate transgen true and false to test code logic */
281a2c6b5cSJunchao Zhang   ierr = MatSetOption(A,MAT_FORM_EXPLICIT_TRANSPOSE,test_cusparse_transgen);CHKERRQ(ierr);
2954da937aSStefano Zampini   test_cusparse_transgen = (PetscBool)!test_cusparse_transgen;
30d9b9dddeSStefano Zampini   PetscFunctionReturn(0);
31d9b9dddeSStefano Zampini }
32d9b9dddeSStefano Zampini 
33bbfc976bSStefano Zampini int main(int argc,char **argv)
34bbfc976bSStefano Zampini {
35bbfc976bSStefano Zampini   Mat            A,B,C,S;
36bbfc976bSStefano Zampini   Vec            t,v;
37bbfc976bSStefano Zampini   PetscScalar    *vv,*aa;
38bbfc976bSStefano Zampini   PetscInt       n=30,k=6,l=0,i,Istart,Iend,nloc,bs,test=1;
39bbfc976bSStefano Zampini   PetscBool      flg,reset,use_shell = PETSC_FALSE;
40bbfc976bSStefano Zampini   PetscErrorCode ierr;
41bbfc976bSStefano Zampini   VecType        vtype;
42bbfc976bSStefano Zampini 
43bbfc976bSStefano Zampini   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
44bbfc976bSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
45bbfc976bSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);CHKERRQ(ierr);
46bbfc976bSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);CHKERRQ(ierr);
47bbfc976bSStefano Zampini   ierr = PetscOptionsGetInt(NULL,NULL,"-test",&test,NULL);CHKERRQ(ierr);
48bbfc976bSStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-use_shell",&use_shell,NULL);CHKERRQ(ierr);
49bbfc976bSStefano Zampini   if (k < 0) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_USER,"k %D must be positive",k);
50bbfc976bSStefano Zampini   if (l < 0) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_USER,"l %D must be positive",l);
51bbfc976bSStefano Zampini   if (l > k) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_USER,"l %D must be smaller or equal than k %D\n",l,k);
52bbfc976bSStefano Zampini 
53bbfc976bSStefano Zampini   /* sparse matrix */
54bbfc976bSStefano Zampini   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
55bbfc976bSStefano Zampini   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr);
56bbfc976bSStefano Zampini   ierr = MatSetType(A,MATAIJCUSPARSE);CHKERRQ(ierr);
570a7ee96eSStefano Zampini   ierr = MatSetOptionsPrefix(A,"A_");CHKERRQ(ierr);
58bbfc976bSStefano Zampini   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
59bbfc976bSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
60bbfc976bSStefano Zampini 
61d9b9dddeSStefano Zampini   /* test special case for SeqAIJCUSPARSE to generate explicit transpose (not default) */
621a2c6b5cSJunchao Zhang   ierr = MatSetOption(A,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);CHKERRQ(ierr);
63d9b9dddeSStefano Zampini 
64bbfc976bSStefano Zampini   ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
65bbfc976bSStefano Zampini   for (i=Istart;i<Iend;i++) {
66bbfc976bSStefano Zampini     if (i>0) { ierr = MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
67bbfc976bSStefano Zampini     if (i<n-1) { ierr = MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
68bbfc976bSStefano Zampini     ierr = MatSetValue(A,i,i,2.0,INSERT_VALUES);CHKERRQ(ierr);
69bbfc976bSStefano Zampini   }
70bbfc976bSStefano Zampini   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71bbfc976bSStefano Zampini   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
72bbfc976bSStefano Zampini 
73bbfc976bSStefano Zampini   /* template vector */
74bbfc976bSStefano Zampini   ierr = MatCreateVecs(A,NULL,&t);CHKERRQ(ierr);
75bbfc976bSStefano Zampini   ierr = VecGetType(t,&vtype);CHKERRQ(ierr);
76bbfc976bSStefano Zampini 
77bbfc976bSStefano Zampini   /* long vector, contains the stacked columns of an nxk dense matrix */
78bbfc976bSStefano Zampini   ierr = VecGetLocalSize(t,&nloc);CHKERRQ(ierr);
79bbfc976bSStefano Zampini   ierr = VecGetBlockSize(t,&bs);CHKERRQ(ierr);
80bbfc976bSStefano Zampini   ierr = VecCreate(PetscObjectComm((PetscObject)t),&v);CHKERRQ(ierr);
81bbfc976bSStefano Zampini   ierr = VecSetType(v,vtype);CHKERRQ(ierr);
82bbfc976bSStefano Zampini   ierr = VecSetSizes(v,k*nloc,k*n);CHKERRQ(ierr);
83bbfc976bSStefano Zampini   ierr = VecSetBlockSize(v,bs);CHKERRQ(ierr);
84bbfc976bSStefano Zampini   ierr = VecSetRandom(v,NULL);CHKERRQ(ierr);
85bbfc976bSStefano Zampini 
86bbfc976bSStefano Zampini   /* dense matrix that contains the columns of v */
87bbfc976bSStefano Zampini   ierr = VecCUDAGetArray(v,&vv);CHKERRQ(ierr);
88bbfc976bSStefano Zampini 
89bbfc976bSStefano Zampini   /* test few cases for MatDenseCUDA handling pointers */
90bbfc976bSStefano Zampini   switch (test) {
91bbfc976bSStefano Zampini   case 1:
92bbfc976bSStefano Zampini     ierr = MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,vv,&B);CHKERRQ(ierr); /* pass a pointer to avoid allocation of storage */
93bbfc976bSStefano Zampini     ierr = MatDenseCUDAReplaceArray(B,NULL);CHKERRQ(ierr);  /* replace with a null pointer, the value after BVRestoreMat */
94bbfc976bSStefano Zampini     ierr = MatDenseCUDAPlaceArray(B,vv+l*nloc);CHKERRQ(ierr);  /* set the actual pointer */
95bbfc976bSStefano Zampini     reset = PETSC_TRUE;
96bbfc976bSStefano Zampini     break;
97bbfc976bSStefano Zampini   case 2:
98bbfc976bSStefano Zampini     ierr = MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,NULL,&B);CHKERRQ(ierr);
99bbfc976bSStefano Zampini     ierr = MatDenseCUDAPlaceArray(B,vv+l*nloc);CHKERRQ(ierr);  /* set the actual pointer */
100bbfc976bSStefano Zampini     reset = PETSC_TRUE;
101bbfc976bSStefano Zampini     break;
102bbfc976bSStefano Zampini   default:
103bbfc976bSStefano Zampini     ierr = MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,vv+l*nloc,&B);CHKERRQ(ierr);
104bbfc976bSStefano Zampini     reset = PETSC_FALSE;
105bbfc976bSStefano Zampini     break;
106bbfc976bSStefano Zampini   }
107bbfc976bSStefano Zampini   ierr = VecCUDARestoreArray(v,&vv);CHKERRQ(ierr);
108bbfc976bSStefano Zampini 
109bbfc976bSStefano Zampini   /* Test MatMatMult */
110bbfc976bSStefano Zampini   if (use_shell) {
11154da937aSStefano Zampini     /* we could have called the general convertor below, but we explicit set the operations
11254da937aSStefano Zampini        ourselves to test MatProductSymbolic_X_Dense, MatProductNumeric_X_Dense code */
11354da937aSStefano Zampini     /* ierr = MatConvert(A,MATSHELL,MAT_INITIAL_MATRIX,&S);CHKERRQ(ierr); */
114bbfc976bSStefano Zampini     ierr = MatCreateShell(PetscObjectComm((PetscObject)v),nloc,nloc,n,n,A,&S);CHKERRQ(ierr);
115bbfc976bSStefano Zampini     ierr = MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_S);CHKERRQ(ierr);
116d9b9dddeSStefano Zampini     ierr = MatShellSetOperation(S,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_S);CHKERRQ(ierr);
117bbfc976bSStefano Zampini     ierr = MatShellSetVecType(S,vtype);CHKERRQ(ierr);
118bbfc976bSStefano Zampini   } else {
119bbfc976bSStefano Zampini     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
120bbfc976bSStefano Zampini     S    = A;
121bbfc976bSStefano Zampini   }
122bbfc976bSStefano Zampini 
123bbfc976bSStefano Zampini   ierr = MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,NULL,&C);CHKERRQ(ierr);
124bbfc976bSStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
125bbfc976bSStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
126bbfc976bSStefano Zampini   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
127bbfc976bSStefano Zampini   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
128d9b9dddeSStefano Zampini 
129d9b9dddeSStefano Zampini   /* test MatMatMult */
130bbfc976bSStefano Zampini   ierr = MatProductCreateWithMat(S,B,NULL,C);CHKERRQ(ierr);
131bbfc976bSStefano Zampini   ierr = MatProductSetType(C,MATPRODUCT_AB);CHKERRQ(ierr);
132bbfc976bSStefano Zampini   ierr = MatProductSetFromOptions(C);CHKERRQ(ierr);
133bbfc976bSStefano Zampini   ierr = MatProductSymbolic(C);CHKERRQ(ierr);
134bbfc976bSStefano Zampini   ierr = MatProductNumeric(C);CHKERRQ(ierr);
135bbfc976bSStefano Zampini   ierr = MatMatMultEqual(S,B,C,10,&flg);CHKERRQ(ierr);
136*1e1ea65dSPierre Jolivet   if (!flg) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error MatMatMult\n");CHKERRQ(ierr); }
137bbfc976bSStefano Zampini 
138d9b9dddeSStefano Zampini   /* test MatTransposeMatMult */
139d9b9dddeSStefano Zampini   ierr = MatProductCreateWithMat(S,B,NULL,C);CHKERRQ(ierr);
140d9b9dddeSStefano Zampini   ierr = MatProductSetType(C,MATPRODUCT_AtB);CHKERRQ(ierr);
141d9b9dddeSStefano Zampini   ierr = MatProductSetFromOptions(C);CHKERRQ(ierr);
142d9b9dddeSStefano Zampini   ierr = MatProductSymbolic(C);CHKERRQ(ierr);
143d9b9dddeSStefano Zampini   ierr = MatProductNumeric(C);CHKERRQ(ierr);
144d9b9dddeSStefano Zampini   ierr = MatTransposeMatMultEqual(S,B,C,10,&flg);CHKERRQ(ierr);
145*1e1ea65dSPierre Jolivet   if (!flg) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error MatTransposeMatMult\n");CHKERRQ(ierr); }
146d9b9dddeSStefano Zampini 
147bbfc976bSStefano Zampini   ierr = MatDestroy(&C);CHKERRQ(ierr);
148bbfc976bSStefano Zampini   ierr = MatDestroy(&S);CHKERRQ(ierr);
149bbfc976bSStefano Zampini 
150bbfc976bSStefano Zampini   /* finished using B */
151bbfc976bSStefano Zampini   ierr = MatDenseCUDAGetArray(B,&aa);CHKERRQ(ierr);
152bbfc976bSStefano Zampini   if (vv != aa-l*nloc) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong array");
153bbfc976bSStefano Zampini   ierr = MatDenseCUDARestoreArray(B,&aa);CHKERRQ(ierr);
154bbfc976bSStefano Zampini   if (reset) {
155bbfc976bSStefano Zampini     ierr = MatDenseCUDAResetArray(B);CHKERRQ(ierr);
156bbfc976bSStefano Zampini   }
157bbfc976bSStefano Zampini   ierr = VecCUDARestoreArray(v,&vv);CHKERRQ(ierr);
158bbfc976bSStefano Zampini 
159bbfc976bSStefano Zampini   if (test == 1) {
160bbfc976bSStefano Zampini     ierr = MatDenseCUDAGetArray(B,&aa);CHKERRQ(ierr);
161bbfc976bSStefano Zampini     if (aa) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected a null pointer");
162bbfc976bSStefano Zampini     ierr = MatDenseCUDARestoreArray(B,&aa);CHKERRQ(ierr);
163bbfc976bSStefano Zampini   }
164bbfc976bSStefano Zampini 
165bbfc976bSStefano Zampini   /* free work space */
166bbfc976bSStefano Zampini   ierr = MatDestroy(&B);CHKERRQ(ierr);
167bbfc976bSStefano Zampini   ierr = MatDestroy(&A);CHKERRQ(ierr);
168bbfc976bSStefano Zampini   ierr = VecDestroy(&t);CHKERRQ(ierr);
169bbfc976bSStefano Zampini   ierr = VecDestroy(&v);CHKERRQ(ierr);
170bbfc976bSStefano Zampini   ierr = PetscFinalize();
171bbfc976bSStefano Zampini   return ierr;
172bbfc976bSStefano Zampini }
173bbfc976bSStefano Zampini 
174bbfc976bSStefano Zampini /*TEST
175bbfc976bSStefano Zampini 
176bbfc976bSStefano Zampini   build:
177bbfc976bSStefano Zampini     requires: cuda
178bbfc976bSStefano Zampini 
179bbfc976bSStefano Zampini   test:
180bbfc976bSStefano Zampini     requires: cuda
181bbfc976bSStefano Zampini     suffix: 1
182bbfc976bSStefano Zampini     nsize: {{1 2}}
1830a7ee96eSStefano Zampini     args: -A_mat_type {{aij aijcusparse}} -test {{0 1 2}} -k 6 -l {{0 5}} -use_shell {{0 1}}
184bbfc976bSStefano Zampini 
185bbfc976bSStefano Zampini TEST*/
186