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