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