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