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>
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown Mat A, B;
8c4762a1bSJed Brown Vec diag;
9c4762a1bSJed Brown PetscMPIInt size, rank;
10c4762a1bSJed Brown
11327415f7SBarry Smith PetscFunctionBeginUser;
12*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
15c4762a1bSJed Brown
169566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE));
189566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A, 2));
199566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATBAIJ));
209566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
21c4762a1bSJed Brown
229566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &diag, NULL));
239566063dSJacob Faibussowitsch PetscCall(VecSet(diag, 1.0));
249566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A, diag, INSERT_VALUES));
259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
279566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
28c4762a1bSJed Brown
299566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(A, size, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B));
3048a46eb9SPierre Jolivet if (rank == 0) PetscCall(MatView(B, PETSC_VIEWER_STDOUT_SELF));
31c4762a1bSJed Brown
329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag));
359566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
36b122ec5aSJacob Faibussowitsch return 0;
37c4762a1bSJed Brown }
38c4762a1bSJed Brown
39c4762a1bSJed Brown /*TEST
40c4762a1bSJed Brown
41c4762a1bSJed Brown test:
42c4762a1bSJed Brown
43c4762a1bSJed Brown test:
44c4762a1bSJed Brown suffix: 2
45c4762a1bSJed Brown nsize: 3
46c4762a1bSJed Brown
47c4762a1bSJed Brown TEST*/
48