static char help[] = "Tests MatCreateDenseCUDA(), MatDenseCUDAPlaceArray(), MatDenseCUDAReplaceArray(), MatDenseCUDAResetArray()\n"; #include static PetscErrorCode MatMult_S(Mat S,Vec x,Vec y) { PetscErrorCode ierr; Mat A; PetscFunctionBeginUser; ierr = MatShellGetContext(S,&A);CHKERRQ(ierr); ierr = MatMult(A,x,y);CHKERRQ(ierr); PetscFunctionReturn(0); } static PetscBool test_cusparse_transgen = PETSC_FALSE; static PetscErrorCode MatMultTranspose_S(Mat S,Vec x,Vec y) { PetscErrorCode ierr; Mat A; PetscFunctionBeginUser; ierr = MatShellGetContext(S,&A);CHKERRQ(ierr); ierr = MatMultTranspose(A,x,y);CHKERRQ(ierr); /* alternate transgen true and false to test code logic */ ierr = MatSetOption(A,MAT_FORM_EXPLICIT_TRANSPOSE,test_cusparse_transgen);CHKERRQ(ierr); test_cusparse_transgen = (PetscBool)!test_cusparse_transgen; PetscFunctionReturn(0); } int main(int argc,char **argv) { Mat A,B,C,S; Vec t,v; PetscScalar *vv,*aa; PetscInt n=30,k=6,l=0,i,Istart,Iend,nloc,bs,test=1; PetscBool flg,reset,use_shell = PETSC_FALSE; PetscErrorCode ierr; VecType vtype; ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-test",&test,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,NULL,"-use_shell",&use_shell,NULL);CHKERRQ(ierr); if (k < 0) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_USER,"k %D must be positive",k); if (l < 0) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_USER,"l %D must be positive",l); if (l > k) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_USER,"l %D must be smaller or equal than k %D\n",l,k); /* sparse matrix */ ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr); ierr = MatSetType(A,MATAIJCUSPARSE);CHKERRQ(ierr); ierr = MatSetOptionsPrefix(A,"A_");CHKERRQ(ierr); ierr = MatSetFromOptions(A);CHKERRQ(ierr); ierr = MatSetUp(A);CHKERRQ(ierr); /* test special case for SeqAIJCUSPARSE to generate explicit transpose (not default) */ ierr = MatSetOption(A,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);CHKERRQ(ierr); ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr); for (i=Istart;i0) { ierr = MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);CHKERRQ(ierr); } if (i