xref: /petsc/src/mat/tutorials/ex11.c (revision f97672e55eacc8688507b9471cd7ec2664d7f203)
1 
2 static char help[] = "Tests MatMeshToDual()\n\n";
3 
4 /*
5   Include "petscmat.h" so that we can use matrices.
6   automatically includes:
7      petscsys.h       - base PETSc routines   petscvec.h    - vectors
8      petscmat.h    - matrices
9      petscis.h     - index sets            petscviewer.h - viewers
10 */
11 #include <petscmat.h>
12 
13 int main(int argc,char **args)
14 {
15   Mat             mesh,dual;
16   PetscInt        Nvertices = 6;       /* total number of vertices */
17   PetscInt        ncells    = 2;       /* number cells on this process */
18   PetscInt        *ii,*jj;
19   PetscMPIInt     size,rank;
20   MatPartitioning part;
21   IS              is;
22 
23   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
24   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&size));
25   PetscCheck(size == 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This example is for exactly two processes");
26   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD,&rank));
27 
28   PetscCall(PetscMalloc1(3,&ii));
29   PetscCall(PetscMalloc1(6,&jj));
30   ii[0] = 0; ii[1] = 3; ii[2] = 6;
31   if (rank == 0) {
32     jj[0] = 0; jj[1] = 1; jj[2] = 2; jj[3] = 1; jj[4] = 2; jj[5] = 3;
33   } else {
34     jj[0] = 1; jj[1] = 4; jj[2] = 5; jj[3] = 1; jj[4] = 3; jj[5] = 5;
35   }
36   PetscCall(MatCreateMPIAdj(MPI_COMM_WORLD,ncells,Nvertices,ii,jj,NULL,&mesh));
37   PetscCall(MatMeshToCellGraph(mesh,2,&dual));
38   PetscCall(MatView(dual,PETSC_VIEWER_STDOUT_WORLD));
39 
40   PetscCall(MatPartitioningCreate(MPI_COMM_WORLD,&part));
41   PetscCall(MatPartitioningSetAdjacency(part,dual));
42   PetscCall(MatPartitioningSetFromOptions(part));
43   PetscCall(MatPartitioningApply(part,&is));
44   PetscCall(ISView(is,PETSC_VIEWER_STDOUT_WORLD));
45   PetscCall(ISDestroy(&is));
46   PetscCall(MatPartitioningDestroy(&part));
47 
48   PetscCall(MatDestroy(&mesh));
49   PetscCall(MatDestroy(&dual));
50   PetscCall(PetscFinalize());
51   return 0;
52 }
53 
54 /*TEST
55 
56    build:
57      requires: parmetis
58 
59    test:
60       nsize: 2
61       requires: parmetis
62 
63 TEST*/
64