Actual source code: ex35.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: static char help[] = "Tests MatCreateSubMatrices().\n\n";

  4:  #include <petscmat.h>

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

 14:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 15:   MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,NULL,&A);
 16:   for (i=0; i<m; i++) {
 17:     for (j=0; j<n; j++) {
 18:       v = -1.0;  Ii = j + n*i;
 19:       if (i>0)   {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 20:       if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 21:       if (j>0)   {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 22:       if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 23:       v = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 24:     }
 25:   }
 26:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 27:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 28:   MatView(A,PETSC_VIEWER_STDOUT_SELF);

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

 38:   /* take a strided block */
 39:   ISCreateStride(PETSC_COMM_WORLD,m,0,2,&isrow);
 40:   MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub);
 41:   B    = *Bsub;
 42:   ISDestroy(&isrow);
 43:   MatView(B,PETSC_VIEWER_STDOUT_SELF);
 44:   MatDestroySubMatrices(1,&Bsub);

 46:   /* take the last block */
 47:   ISCreateStride(PETSC_COMM_WORLD,m,N-m-1,1,&isrow);
 48:   MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub);
 49:   B    = *Bsub;
 50:   ISDestroy(&isrow);
 51:   MatView(B,PETSC_VIEWER_STDOUT_SELF);

 53:   MatDestroySubMatrices(1,&Bsub);
 54:   MatDestroy(&A);

 56:   PetscFinalize();
 57:   return ierr;
 58: }



 62: /*TEST

 64:    test:

 66: TEST*/