1 static char help[] = "Test MatCreateRedundantMatrix for a BAIJ matrix.\n\ 2 Contributed by Lawrence Mitchell, Feb. 21, 2017\n\n"; 3 4 #include <petscmat.h> 5 int main(int argc, char **args) { 6 Mat A, B; 7 Vec diag; 8 PetscMPIInt size, rank; 9 10 PetscFunctionBeginUser; 11 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 12 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 13 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 14 15 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 16 PetscCall(MatSetSizes(A, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE)); 17 PetscCall(MatSetBlockSize(A, 2)); 18 PetscCall(MatSetType(A, MATBAIJ)); 19 PetscCall(MatSetUp(A)); 20 21 PetscCall(MatCreateVecs(A, &diag, NULL)); 22 PetscCall(VecSet(diag, 1.0)); 23 PetscCall(MatDiagonalSet(A, diag, INSERT_VALUES)); 24 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 25 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 26 PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 27 28 PetscCall(MatCreateRedundantMatrix(A, size, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B)); 29 if (rank == 0) { PetscCall(MatView(B, PETSC_VIEWER_STDOUT_SELF)); } 30 31 PetscCall(MatDestroy(&A)); 32 PetscCall(MatDestroy(&B)); 33 PetscCall(VecDestroy(&diag)); 34 PetscCall(PetscFinalize()); 35 return 0; 36 } 37 38 /*TEST 39 40 test: 41 42 test: 43 suffix: 2 44 nsize: 3 45 46 TEST*/ 47