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