Actual source code: ex35.c

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

  3: static char help[] = "Tests MatGetSubMatrices().nn";

 5:  #include petscmat.h

  7: int main(int argc,char **args)
  8: {
  9:   Mat    A,B,*Bsub;
 10:   int    i,j,m = 6,n = 6,N = 36,ierr,I,J;
 11:   PetscScalar v;
 12:   IS     isrow;

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

 16:   MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,PETSC_NULL,&A);
 17:   for (i=0; i<m; i++) {
 18:     for (j=0; j<n; j++) {
 19:       v = -1.0;  I = j + n*i;
 20:       if (i>0)   {J = I - n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 21:       if (i<m-1) {J = I + n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 22:       if (j>0)   {J = I - 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 23:       if (j<n-1) {J = I + 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 24:       v = 4.0; MatSetValues(A,1,&I,1,&I,&v,INSERT_VALUES);
 25:     }
 26:   }
 27:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 28:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 29:   MatView(A,PETSC_VIEWER_STDOUT_SELF);

 31:   /* take the first diagonal block */
 32:   ISCreateStride(PETSC_COMM_WORLD,m,0,1,&isrow);
 33:   MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub);
 34:   B = *Bsub;
 35:   PetscFree(Bsub);
 36:   ISDestroy(isrow);
 37:   MatView(B,PETSC_VIEWER_STDOUT_SELF);
 38:   MatDestroy(B);

 40:   /* take a strided block */
 41:   ISCreateStride(PETSC_COMM_WORLD,m,0,2,&isrow);
 42:   MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub);
 43:   B = *Bsub;
 44:   PetscFree(Bsub);
 45:   ISDestroy(isrow);
 46:   MatView(B,PETSC_VIEWER_STDOUT_SELF);
 47:   MatDestroy(B);

 49:   /* take the last block */
 50:   ISCreateStride(PETSC_COMM_WORLD,m,N-m-1,1,&isrow);
 51:   MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub);
 52:   B = *Bsub;
 53:   PetscFree(Bsub);
 54:   ISDestroy(isrow);
 55:   MatView(B,PETSC_VIEWER_STDOUT_SELF);
 56: 
 57:   MatDestroy(B);
 58:   MatDestroy(A);

 60:   PetscFinalize();
 61:   return 0;
 62: }