1c4762a1bSJed Brown static char help[] = "Test MatCreateRedundantMatrix for a BAIJ matrix.\n\ 2c4762a1bSJed Brown Contributed by Lawrence Mitchell, Feb. 21, 2017\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 59371c9d4SSatish Balay int main(int argc, char **args) { 6c4762a1bSJed Brown Mat A, B; 7c4762a1bSJed Brown Vec diag; 8c4762a1bSJed Brown PetscMPIInt size, rank; 9c4762a1bSJed Brown 10327415f7SBarry Smith PetscFunctionBeginUser; 119566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 14c4762a1bSJed Brown 159566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 169566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE)); 179566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A, 2)); 189566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATBAIJ)); 199566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 20c4762a1bSJed Brown 219566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &diag, NULL)); 229566063dSJacob Faibussowitsch PetscCall(VecSet(diag, 1.0)); 239566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A, diag, INSERT_VALUES)); 249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 269566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 27c4762a1bSJed Brown 289566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(A, size, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B)); 29*48a46eb9SPierre Jolivet if (rank == 0) PetscCall(MatView(B, PETSC_VIEWER_STDOUT_SELF)); 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 349566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 35b122ec5aSJacob Faibussowitsch return 0; 36c4762a1bSJed Brown } 37c4762a1bSJed Brown 38c4762a1bSJed Brown /*TEST 39c4762a1bSJed Brown 40c4762a1bSJed Brown test: 41c4762a1bSJed Brown 42c4762a1bSJed Brown test: 43c4762a1bSJed Brown suffix: 2 44c4762a1bSJed Brown nsize: 3 45c4762a1bSJed Brown 46c4762a1bSJed Brown TEST*/ 47