Actual source code: ex71.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2: static char help[] = "Passes a sparse matrix to MATLAB.\n\n";

  4:  #include <petscmat.h>

  6: int main(int argc,char **args)
  7: {
  9:   PetscInt       m   = 4,n = 5,i,j,II,J;
 10:   PetscScalar    one = 1.0,v;
 11:   Vec            x;
 12:   Mat            A;

 14:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 15:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 16:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);

 18:   MatCreate(PETSC_COMM_WORLD,&A);
 19:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 20:   MatSetFromOptions(A);

 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; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES);}
 26:       if (i<m-1) {J = II + n; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES);}
 27:       if (j>0)   {J = II - 1; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES);}
 28:       if (j<n-1) {J = II + 1; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES);}
 29:       v = 4.0; MatSetValues(A,1,&II,1,&II,&v,INSERT_VALUES);
 30:     }
 31:   }
 32:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 33:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 34: #if defined(PETSC_USE_SOCKET_VIEWER)
 35:   MatView(A,PETSC_VIEWER_SOCKET_WORLD);
 36: #endif
 37:   VecCreateSeq(PETSC_COMM_SELF,m,&x);
 38:   VecSet(x,one);
 39: #if defined(PETSC_USE_SOCKET_VIEWER)
 40:   VecView(x,PETSC_VIEWER_SOCKET_WORLD);
 41: #endif

 43:   PetscSleep(30);

 45:   VecDestroy(&x);
 46:   MatDestroy(&A);
 47:   PetscFinalize();
 48:   return ierr;
 49: }