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