1 2 static char help[] = "Tests MatGetColumnVector()."; 3 4 #include <petscmat.h> 5 6 int main(int argc,char **args) 7 { 8 Mat C; 9 PetscInt i,j,m = 3,n = 2,Ii,J,col = 0; 10 PetscMPIInt size,rank; 11 PetscScalar v; 12 Vec yy; 13 14 CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 15 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-col",&col,NULL)); 16 17 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 18 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 19 n = 2*size; 20 21 /* create the matrix for the five point stencil, YET AGAIN*/ 22 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 23 CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 24 CHKERRQ(MatSetFromOptions(C)); 25 CHKERRQ(MatSeqAIJSetPreallocation(C,5,NULL)); 26 CHKERRQ(MatMPIAIJSetPreallocation(C,5,NULL,5,NULL)); 27 CHKERRQ(MatSetUp(C)); 28 29 for (i=0; i<m; i++) { 30 for (j=2*rank; j<2*rank+2; j++) { 31 v = -1.0; Ii = j + n*i; 32 if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 33 if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 34 if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 35 if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 36 v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 37 } 38 } 39 CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 40 CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 41 CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 42 43 CHKERRQ(MatCreateVecs(C,NULL,&yy)); 44 CHKERRQ(VecSetFromOptions(yy)); 45 46 CHKERRQ(MatGetColumnVector(C,yy,col)); 47 48 CHKERRQ(VecView(yy,PETSC_VIEWER_STDOUT_WORLD)); 49 50 CHKERRQ(VecDestroy(&yy)); 51 CHKERRQ(MatDestroy(&C)); 52 CHKERRQ(PetscFinalize()); 53 return 0; 54 } 55 56 /*TEST 57 58 test: 59 nsize: 3 60 args: -col 7 61 62 test: 63 suffix: dense 64 nsize: 3 65 args: -col 7 -mat_type dense -vec_type {{mpi standard}} 66 filter: grep -v type 67 68 test: 69 requires: cuda 70 suffix: dense_cuda 71 nsize: 3 72 output_file: output/ex60_dense.out 73 args: -col 7 -mat_type {{mpidense mpidensecuda}} -vec_type {{mpi standard cuda mpicuda}} 74 filter: grep -v type 75 76 TEST*/ 77