xref: /petsc/src/mat/tests/ex207.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
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   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr);
14   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr);
15 
16   ierr = MatCreate(PETSC_COMM_WORLD, &A);CHKERRQ(ierr);
17   ierr = MatSetSizes(A, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
18   ierr = MatSetBlockSize(A, 2);CHKERRQ(ierr);
19   ierr = MatSetType(A, MATBAIJ);CHKERRQ(ierr);
20   ierr = MatSetUp(A);CHKERRQ(ierr);
21 
22   ierr = MatCreateVecs(A, &diag, NULL);CHKERRQ(ierr);
23   ierr = VecSet(diag, 1.0);CHKERRQ(ierr);
24   ierr = MatDiagonalSet(A, diag, INSERT_VALUES);CHKERRQ(ierr);
25   ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
26   ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
27   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
28 
29   ierr = MatCreateRedundantMatrix(A, size, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B);CHKERRQ(ierr);
30   if (!rank) {
31     ierr = MatView(B,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
32   }
33 
34   ierr = MatDestroy(&A);CHKERRQ(ierr);
35   ierr = MatDestroy(&B);CHKERRQ(ierr);
36   ierr = VecDestroy(&diag);CHKERRQ(ierr);
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