1*524fbf86SStefano Zampini static char help[] = "Test MatDenseGetSubMatrix() on a CUDA matrix.\n"; 2*524fbf86SStefano Zampini 3*524fbf86SStefano Zampini #include <petscmat.h> 4*524fbf86SStefano Zampini 5*524fbf86SStefano Zampini int main(int argc,char **argv) 6*524fbf86SStefano Zampini { 7*524fbf86SStefano Zampini Mat A,B; 8*524fbf86SStefano Zampini PetscScalar *b; 9*524fbf86SStefano Zampini PetscInt n=4,lda=5,i,k; 10*524fbf86SStefano Zampini PetscBool cuda=PETSC_FALSE; 11*524fbf86SStefano Zampini 12*524fbf86SStefano Zampini PetscCall(PetscInitialize(&argc,&argv,0,help)); 13*524fbf86SStefano Zampini PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 14*524fbf86SStefano Zampini PetscCall(PetscOptionsGetInt(NULL,NULL,"-lda",&lda,NULL)); 15*524fbf86SStefano Zampini PetscCall(PetscOptionsGetBool(NULL,NULL,"-cuda",&cuda,NULL)); 16*524fbf86SStefano Zampini PetscCheck(lda>=n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"lda %" PetscInt_FMT " < n %" PetscInt_FMT,lda,n); 17*524fbf86SStefano Zampini 18*524fbf86SStefano Zampini #if defined(PETSC_HAVE_CUDA) 19*524fbf86SStefano Zampini if (cuda) PetscCall(MatCreateSeqDenseCUDA(PETSC_COMM_SELF,lda,n,NULL,&A)); 20*524fbf86SStefano Zampini else 21*524fbf86SStefano Zampini #endif 22*524fbf86SStefano Zampini PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,lda,n,NULL,&A)); 23*524fbf86SStefano Zampini 24*524fbf86SStefano Zampini for (k=0;k<3;k++) { 25*524fbf86SStefano Zampini PetscCall(MatDenseGetSubMatrix(A,0,n,0,n,&B)); 26*524fbf86SStefano Zampini PetscCall(MatDenseGetArray(B,&b)); 27*524fbf86SStefano Zampini for (i=0; i<n; i++) { 28*524fbf86SStefano Zampini b[i+i*lda] = 2.0*(i+1); 29*524fbf86SStefano Zampini if (i>0) b[i+(i-1)*lda] = (PetscReal)(k+1); 30*524fbf86SStefano Zampini } 31*524fbf86SStefano Zampini PetscCall(MatDenseRestoreArray(B,&b)); 32*524fbf86SStefano Zampini PetscCall(MatDenseRestoreSubMatrix(A,&B)); 33*524fbf86SStefano Zampini PetscCall(MatView(A,NULL)); 34*524fbf86SStefano Zampini } 35*524fbf86SStefano Zampini 36*524fbf86SStefano Zampini PetscCall(MatDestroy(&A)); 37*524fbf86SStefano Zampini PetscCall(PetscFinalize()); 38*524fbf86SStefano Zampini return 0; 39*524fbf86SStefano Zampini } 40*524fbf86SStefano Zampini 41*524fbf86SStefano Zampini /*TEST 42*524fbf86SStefano Zampini 43*524fbf86SStefano Zampini testset: 44*524fbf86SStefano Zampini output_file: output/ex257_1.out 45*524fbf86SStefano Zampini diff_args: -j 46*524fbf86SStefano Zampini test: 47*524fbf86SStefano Zampini suffix: 1 48*524fbf86SStefano Zampini test: 49*524fbf86SStefano Zampini suffix: 1_cuda 50*524fbf86SStefano Zampini args: -cuda 51*524fbf86SStefano Zampini requires: cuda 52*524fbf86SStefano Zampini filter: sed -e "s/seqdensecuda/seqdense/" 53*524fbf86SStefano Zampini 54*524fbf86SStefano Zampini TEST*/ 55