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