Actual source code: ex207.c

  1: static char help[] = "Test MatCreateRedundantMatrix for a BAIJ matrix.\n\
  2:                       Contributed by Lawrence Mitchell, Feb. 21, 2017\n\n";

  4: #include <petscmat.h>
  5: int main(int argc,char **args)
  6: {
  7:   Mat               A,B;
  8:   Vec               diag;
  9:   PetscMPIInt       size,rank;

 11:   PetscInitialize(&argc,&args,(char*)0,help);
 12:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 13:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);

 15:   MatCreate(PETSC_COMM_WORLD, &A);
 16:   MatSetSizes(A, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE);
 17:   MatSetBlockSize(A, 2);
 18:   MatSetType(A, MATBAIJ);
 19:   MatSetUp(A);

 21:   MatCreateVecs(A, &diag, NULL);
 22:   VecSet(diag, 1.0);
 23:   MatDiagonalSet(A, diag, INSERT_VALUES);
 24:   MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
 25:   MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
 26:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);

 28:   MatCreateRedundantMatrix(A, size, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B);
 29:   if (rank == 0) {
 30:     MatView(B,PETSC_VIEWER_STDOUT_SELF);
 31:   }

 33:   MatDestroy(&A);
 34:   MatDestroy(&B);
 35:   VecDestroy(&diag);
 36:   PetscFinalize();
 37:   return 0;
 38: }

 40: /*TEST

 42:    test:

 44:    test:
 45:       suffix: 2
 46:       nsize: 3

 48: TEST*/