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