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 { 7 Mat A,B; 8 Vec diag; 9 PetscMPIInt size,rank; 10 11 CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 12 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 13 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 14 15 CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 16 CHKERRQ(MatSetSizes(A, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE)); 17 CHKERRQ(MatSetBlockSize(A, 2)); 18 CHKERRQ(MatSetType(A, MATBAIJ)); 19 CHKERRQ(MatSetUp(A)); 20 21 CHKERRQ(MatCreateVecs(A, &diag, NULL)); 22 CHKERRQ(VecSet(diag, 1.0)); 23 CHKERRQ(MatDiagonalSet(A, diag, INSERT_VALUES)); 24 CHKERRQ(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 25 CHKERRQ(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 26 CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 27 28 CHKERRQ(MatCreateRedundantMatrix(A, size, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B)); 29 if (rank == 0) { 30 CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_SELF)); 31 } 32 33 CHKERRQ(MatDestroy(&A)); 34 CHKERRQ(MatDestroy(&B)); 35 CHKERRQ(VecDestroy(&diag)); 36 CHKERRQ(PetscFinalize()); 37 return 0; 38 } 39 40 /*TEST 41 42 test: 43 44 test: 45 suffix: 2 46 nsize: 3 47 48 TEST*/ 49