1 2 static char help[] = "Passes a sparse matrix to MATLAB.\n\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc,char **args) 7 { 8 PetscInt m = 4,n = 5,i,j,II,J; 9 PetscScalar one = 1.0,v; 10 Vec x; 11 Mat A; 12 13 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 14 PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 15 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 16 17 PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 18 PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 19 PetscCall(MatSetFromOptions(A)); 20 PetscCall(MatSetUp(A)); 21 22 for (i=0; i<m; i++) { 23 for (j=0; j<n; j++) { 24 v = -1.0; II = j + n*i; 25 if (i>0) {J = II - n; PetscCall(MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES));} 26 if (i<m-1) {J = II + n; PetscCall(MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES));} 27 if (j>0) {J = II - 1; PetscCall(MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES));} 28 if (j<n-1) {J = II + 1; PetscCall(MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES));} 29 v = 4.0; PetscCall(MatSetValues(A,1,&II,1,&II,&v,INSERT_VALUES)); 30 } 31 } 32 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 33 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 34 #if defined(PETSC_USE_SOCKET_VIEWER) 35 PetscCall(MatView(A,PETSC_VIEWER_SOCKET_WORLD)); 36 #endif 37 PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&x)); 38 PetscCall(VecSet(x,one)); 39 #if defined(PETSC_USE_SOCKET_VIEWER) 40 PetscCall(VecView(x,PETSC_VIEWER_SOCKET_WORLD)); 41 #endif 42 43 PetscCall(PetscSleep(30)); 44 45 PetscCall(VecDestroy(&x)); 46 PetscCall(MatDestroy(&A)); 47 PetscCall(PetscFinalize()); 48 return 0; 49 } 50 51 /*TEST 52 53 test: 54 TODO: cannot test with socket viewer 55 56 TEST*/ 57