Actual source code: ex71.c

  1: /*$Id: ex71.c,v 1.45 2001/08/07 03:03:07 balay Exp $*/

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

 5:  #include petscmat.h

  9: int main(int argc,char **args)
 10: {
 11:   int         ierr,m = 4,n = 5,i,j,I,J;
 12:   PetscScalar one = 1.0,v;
 13:   Vec         x;
 14:   Mat         A;
 15:   PetscViewer viewer;

 17:   PetscInitialize(&argc,&args,(char *)0,help);
 18:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 19:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

 21:   PetscViewerSocketOpen(PETSC_COMM_WORLD,"eagle",-1,&viewer);
 22:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&A);
 23:   MatSetFromOptions(A);

 25:   for (i=0; i<m; i++) {
 26:     for (j=0; j<n; j++) {
 27:       v = -1.0;  I = j + n*i;
 28:       if (i>0)   {J = I - n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 29:       if (i<m-1) {J = I + n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 30:       if (j>0)   {J = I - 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 31:       if (j<n-1) {J = I + 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 32:       v = 4.0; MatSetValues(A,1,&I,1,&I,&v,INSERT_VALUES);
 33:     }
 34:   }
 35:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 36:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 38:   MatView(A,viewer);

 40:   VecCreateSeq(PETSC_COMM_SELF,m,&x);
 41:   VecSet(&one,x);
 42:   VecView(x,viewer);
 43: 
 44:   PetscSleep(30);
 45:   PetscObjectDestroy((PetscObject)viewer);

 47:   VecDestroy(x);
 48:   MatDestroy(A);
 49:   PetscFinalize();
 50:   return 0;
 51: }
 52: