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