1524fbf86SStefano Zampini static char help[] = "Test MatDenseGetSubMatrix() on a CUDA matrix.\n"; 2524fbf86SStefano Zampini 3524fbf86SStefano Zampini #include <petscmat.h> 4524fbf86SStefano Zampini 5524fbf86SStefano Zampini int main(int argc,char **argv) 6524fbf86SStefano Zampini { 7524fbf86SStefano Zampini Mat A,B; 8524fbf86SStefano Zampini PetscScalar *b; 9524fbf86SStefano Zampini PetscInt n=4,lda=5,i,k; 10524fbf86SStefano Zampini PetscBool cuda=PETSC_FALSE; 11524fbf86SStefano Zampini 12*327415f7SBarry Smith PetscFunctionBeginUser; 13524fbf86SStefano Zampini PetscCall(PetscInitialize(&argc,&argv,0,help)); 14524fbf86SStefano Zampini PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 15524fbf86SStefano Zampini PetscCall(PetscOptionsGetInt(NULL,NULL,"-lda",&lda,NULL)); 16524fbf86SStefano Zampini PetscCall(PetscOptionsGetBool(NULL,NULL,"-cuda",&cuda,NULL)); 17524fbf86SStefano Zampini PetscCheck(lda>=n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"lda %" PetscInt_FMT " < n %" PetscInt_FMT,lda,n); 18524fbf86SStefano Zampini 19524fbf86SStefano Zampini #if defined(PETSC_HAVE_CUDA) 20524fbf86SStefano Zampini if (cuda) PetscCall(MatCreateSeqDenseCUDA(PETSC_COMM_SELF,lda,n,NULL,&A)); 21524fbf86SStefano Zampini else 22524fbf86SStefano Zampini #endif 23524fbf86SStefano Zampini PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,lda,n,NULL,&A)); 24524fbf86SStefano Zampini 25524fbf86SStefano Zampini for (k=0;k<3;k++) { 26524fbf86SStefano Zampini PetscCall(MatDenseGetSubMatrix(A,0,n,0,n,&B)); 27524fbf86SStefano Zampini PetscCall(MatDenseGetArray(B,&b)); 28524fbf86SStefano Zampini for (i=0; i<n; i++) { 29524fbf86SStefano Zampini b[i+i*lda] = 2.0*(i+1); 30524fbf86SStefano Zampini if (i>0) b[i+(i-1)*lda] = (PetscReal)(k+1); 31524fbf86SStefano Zampini } 32524fbf86SStefano Zampini PetscCall(MatDenseRestoreArray(B,&b)); 33524fbf86SStefano Zampini PetscCall(MatDenseRestoreSubMatrix(A,&B)); 34524fbf86SStefano Zampini PetscCall(MatView(A,NULL)); 35524fbf86SStefano Zampini } 36524fbf86SStefano Zampini 37524fbf86SStefano Zampini PetscCall(MatDestroy(&A)); 38524fbf86SStefano Zampini PetscCall(PetscFinalize()); 39524fbf86SStefano Zampini return 0; 40524fbf86SStefano Zampini } 41524fbf86SStefano Zampini 42524fbf86SStefano Zampini /*TEST 43524fbf86SStefano Zampini 44524fbf86SStefano Zampini testset: 45524fbf86SStefano Zampini output_file: output/ex257_1.out 46524fbf86SStefano Zampini diff_args: -j 47524fbf86SStefano Zampini test: 48524fbf86SStefano Zampini suffix: 1 49524fbf86SStefano Zampini test: 50524fbf86SStefano Zampini suffix: 1_cuda 51524fbf86SStefano Zampini args: -cuda 52524fbf86SStefano Zampini requires: cuda 53524fbf86SStefano Zampini filter: sed -e "s/seqdensecuda/seqdense/" 54524fbf86SStefano Zampini 55524fbf86SStefano Zampini TEST*/ 56