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