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; 109566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(S,&A)); 119566063dSJacob Faibussowitsch PetscCall(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; 229566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(S,&A)); 239566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,x,y)); 2454da937aSStefano Zampini 2554da937aSStefano Zampini /* alternate transgen true and false to test code logic */ 269566063dSJacob Faibussowitsch PetscCall(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 VecType vtype; 39bbfc976bSStefano Zampini 409566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL)); 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL)); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-test",&test,NULL)); 459566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_shell",&use_shell,NULL)); 46*08401ef6SPierre Jolivet PetscCheck(k >= 0,PETSC_COMM_WORLD,PETSC_ERR_USER,"k %" PetscInt_FMT " must be positive",k); 47*08401ef6SPierre Jolivet PetscCheck(l >= 0,PETSC_COMM_WORLD,PETSC_ERR_USER,"l %" PetscInt_FMT " must be positive",l); 48*08401ef6SPierre Jolivet PetscCheck(l <= k,PETSC_COMM_WORLD,PETSC_ERR_USER,"l %" PetscInt_FMT " must be smaller or equal than k %" PetscInt_FMT,l,k); 49bbfc976bSStefano Zampini 50bbfc976bSStefano Zampini /* sparse matrix */ 519566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n)); 539566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATAIJCUSPARSE)); 549566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(A,"A_")); 559566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 569566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 57bbfc976bSStefano Zampini 58d9b9dddeSStefano Zampini /* test special case for SeqAIJCUSPARSE to generate explicit transpose (not default) */ 599566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE)); 60d9b9dddeSStefano Zampini 619566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 62bbfc976bSStefano Zampini for (i=Istart;i<Iend;i++) { 639566063dSJacob Faibussowitsch if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES)); 649566063dSJacob Faibussowitsch if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES)); 659566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES)); 66bbfc976bSStefano Zampini } 679566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 689566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 69bbfc976bSStefano Zampini 70bbfc976bSStefano Zampini /* template vector */ 719566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,NULL,&t)); 729566063dSJacob Faibussowitsch PetscCall(VecGetType(t,&vtype)); 73bbfc976bSStefano Zampini 74bbfc976bSStefano Zampini /* long vector, contains the stacked columns of an nxk dense matrix */ 759566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(t,&nloc)); 769566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(t,&bs)); 779566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)t),&v)); 789566063dSJacob Faibussowitsch PetscCall(VecSetType(v,vtype)); 799566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v,k*nloc,k*n)); 809566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(v,bs)); 819566063dSJacob Faibussowitsch PetscCall(VecSetRandom(v,NULL)); 82bbfc976bSStefano Zampini 83bbfc976bSStefano Zampini /* dense matrix that contains the columns of v */ 849566063dSJacob Faibussowitsch PetscCall(VecCUDAGetArray(v,&vv)); 85bbfc976bSStefano Zampini 86bbfc976bSStefano Zampini /* test few cases for MatDenseCUDA handling pointers */ 87bbfc976bSStefano Zampini switch (test) { 88bbfc976bSStefano Zampini case 1: 899566063dSJacob Faibussowitsch PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,vv,&B)); /* pass a pointer to avoid allocation of storage */ 909566063dSJacob Faibussowitsch PetscCall(MatDenseCUDAReplaceArray(B,NULL)); /* replace with a null pointer, the value after BVRestoreMat */ 919566063dSJacob Faibussowitsch PetscCall(MatDenseCUDAPlaceArray(B,vv+l*nloc)); /* set the actual pointer */ 92bbfc976bSStefano Zampini reset = PETSC_TRUE; 93bbfc976bSStefano Zampini break; 94bbfc976bSStefano Zampini case 2: 959566063dSJacob Faibussowitsch PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,NULL,&B)); 969566063dSJacob Faibussowitsch PetscCall(MatDenseCUDAPlaceArray(B,vv+l*nloc)); /* set the actual pointer */ 97bbfc976bSStefano Zampini reset = PETSC_TRUE; 98bbfc976bSStefano Zampini break; 99bbfc976bSStefano Zampini default: 1009566063dSJacob Faibussowitsch PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,vv+l*nloc,&B)); 101bbfc976bSStefano Zampini reset = PETSC_FALSE; 102bbfc976bSStefano Zampini break; 103bbfc976bSStefano Zampini } 1049566063dSJacob Faibussowitsch PetscCall(VecCUDARestoreArray(v,&vv)); 105bbfc976bSStefano Zampini 106bbfc976bSStefano Zampini /* Test MatMatMult */ 107bbfc976bSStefano Zampini if (use_shell) { 10854da937aSStefano Zampini /* we could have called the general convertor below, but we explicit set the operations 10954da937aSStefano Zampini ourselves to test MatProductSymbolic_X_Dense, MatProductNumeric_X_Dense code */ 1109566063dSJacob Faibussowitsch /* PetscCall(MatConvert(A,MATSHELL,MAT_INITIAL_MATRIX,&S)); */ 1119566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)v),nloc,nloc,n,n,A,&S)); 1129566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_S)); 1139566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_S)); 1149566063dSJacob Faibussowitsch PetscCall(MatShellSetVecType(S,vtype)); 115bbfc976bSStefano Zampini } else { 1169566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 117bbfc976bSStefano Zampini S = A; 118bbfc976bSStefano Zampini } 119bbfc976bSStefano Zampini 1209566063dSJacob Faibussowitsch PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v),nloc,PETSC_DECIDE,n,k-l,NULL,&C)); 1219566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1229566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 1239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 1249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 125d9b9dddeSStefano Zampini 126d9b9dddeSStefano Zampini /* test MatMatMult */ 1279566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(S,B,NULL,C)); 1289566063dSJacob Faibussowitsch PetscCall(MatProductSetType(C,MATPRODUCT_AB)); 1299566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(C)); 1309566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(C)); 1319566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(C)); 1329566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(S,B,C,10,&flg)); 1339566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error MatMatMult\n")); 134bbfc976bSStefano Zampini 135d9b9dddeSStefano Zampini /* test MatTransposeMatMult */ 1369566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(S,B,NULL,C)); 1379566063dSJacob Faibussowitsch PetscCall(MatProductSetType(C,MATPRODUCT_AtB)); 1389566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(C)); 1399566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(C)); 1409566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(C)); 1419566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMultEqual(S,B,C,10,&flg)); 1429566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error MatTransposeMatMult\n")); 143d9b9dddeSStefano Zampini 1449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 146bbfc976bSStefano Zampini 147bbfc976bSStefano Zampini /* finished using B */ 1489566063dSJacob Faibussowitsch PetscCall(MatDenseCUDAGetArray(B,&aa)); 149*08401ef6SPierre Jolivet PetscCheck(vv == aa-l*nloc,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong array"); 1509566063dSJacob Faibussowitsch PetscCall(MatDenseCUDARestoreArray(B,&aa)); 151bbfc976bSStefano Zampini if (reset) { 1529566063dSJacob Faibussowitsch PetscCall(MatDenseCUDAResetArray(B)); 153bbfc976bSStefano Zampini } 1549566063dSJacob Faibussowitsch PetscCall(VecCUDARestoreArray(v,&vv)); 155bbfc976bSStefano Zampini 156bbfc976bSStefano Zampini if (test == 1) { 1579566063dSJacob Faibussowitsch PetscCall(MatDenseCUDAGetArray(B,&aa)); 15828b400f6SJacob Faibussowitsch PetscCheck(!aa,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected a null pointer"); 1599566063dSJacob Faibussowitsch PetscCall(MatDenseCUDARestoreArray(B,&aa)); 160bbfc976bSStefano Zampini } 161bbfc976bSStefano Zampini 162bbfc976bSStefano Zampini /* free work space */ 1639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&t)); 1669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1679566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 168b122ec5aSJacob Faibussowitsch return 0; 169bbfc976bSStefano Zampini } 170bbfc976bSStefano Zampini 171bbfc976bSStefano Zampini /*TEST 172bbfc976bSStefano Zampini 173bbfc976bSStefano Zampini build: 174bbfc976bSStefano Zampini requires: cuda 175bbfc976bSStefano Zampini 176bbfc976bSStefano Zampini test: 177bbfc976bSStefano Zampini requires: cuda 178bbfc976bSStefano Zampini suffix: 1 179bbfc976bSStefano Zampini nsize: {{1 2}} 1800a7ee96eSStefano Zampini args: -A_mat_type {{aij aijcusparse}} -test {{0 1 2}} -k 6 -l {{0 5}} -use_shell {{0 1}} 181bbfc976bSStefano Zampini 182bbfc976bSStefano Zampini TEST*/ 183