xref: /petsc/src/mat/tests/ex207.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Test MatCreateRedundantMatrix for a BAIJ matrix.\n\
2*c4762a1bSJed Brown                       Contributed by Lawrence Mitchell, Feb. 21, 2017\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown int main(int argc,char **args)
6*c4762a1bSJed Brown {
7*c4762a1bSJed Brown   Mat               A,B;
8*c4762a1bSJed Brown   Vec               diag;
9*c4762a1bSJed Brown   PetscErrorCode    ierr;
10*c4762a1bSJed Brown   PetscMPIInt       size,rank;
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
13*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);
14*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD, &A); CHKERRQ(ierr);
17*c4762a1bSJed Brown   ierr = MatSetSizes(A, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
18*c4762a1bSJed Brown   ierr = MatSetBlockSize(A, 2);CHKERRQ(ierr);
19*c4762a1bSJed Brown   ierr = MatSetType(A, MATBAIJ);CHKERRQ(ierr);
20*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown   ierr = MatCreateVecs(A, &diag, NULL);CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = VecSet(diag, 1.0); CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = MatDiagonalSet(A, diag, INSERT_VALUES);CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
27*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown   ierr = MatCreateRedundantMatrix(A, size, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B);CHKERRQ(ierr);
30*c4762a1bSJed Brown   if (!rank) {
31*c4762a1bSJed Brown     ierr = MatView(B,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
32*c4762a1bSJed Brown   }
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = VecDestroy(&diag);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = PetscFinalize();
38*c4762a1bSJed Brown   return ierr;
39*c4762a1bSJed Brown }
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown /*TEST
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown    test:
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown    test:
47*c4762a1bSJed Brown       suffix: 2
48*c4762a1bSJed Brown       nsize: 3
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown TEST*/
51