Actual source code: ex13.c

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

  3: static char help[] = "Tests copying and ordering uniprocessor row-based sparse matrices.nn";

 5:  #include petscmat.h

  7: int main(int argc,char **args)
  8: {
  9:   Mat     C,A;
 10:   int     i,j,m = 5,n = 5,I,J,ierr;
 11:   PetscScalar  v;
 12:   IS      perm,iperm;

 14:   PetscInitialize(&argc,&args,(char *)0,help);

 16:   MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,PETSC_NULL,&C);

 18:   /* create the matrix for the five point stencil, YET AGAIN*/
 19:   for (i=0; i<m; i++) {
 20:     for (j=0; j<n; j++) {
 21:       v = -1.0;  I = j + n*i;
 22:       if (i>0)   {J = I - n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 23:       if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 24:       if (j>0)   {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 25:       if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 26:       v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
 27:     }
 28:   }
 29:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 30:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 32:   MatConvert(C,MATSAME,&A);

 34:   MatGetOrdering(A,MATORDERING_ND,&perm,&iperm);
 35:   ISView(perm,PETSC_VIEWER_STDOUT_SELF);
 36:   ISView(iperm,PETSC_VIEWER_STDOUT_SELF);
 37:   MatView(A,PETSC_VIEWER_STDOUT_SELF);

 39:   ISDestroy(perm);
 40:   ISDestroy(iperm);
 41:   MatDestroy(C);
 42:   MatDestroy(A);
 43:   PetscFinalize();
 44:   return 0;
 45: }