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