xref: /petsc/src/mat/tutorials/ex11.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
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   PetscErrorCode  ierr;
22   PetscInt        Nvertices = 6;       /* total number of vertices */
23   PetscInt        ncells    = 2;       /* number cells on this process */
24   PetscInt        *ii,*jj;
25   PetscMPIInt     size,rank;
26   MatPartitioning part;
27   IS              is;
28 
29   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
30   ierr = MPI_Comm_size(MPI_COMM_WORLD,&size);CHKERRQ(ierr);
31   if (size != 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example is for exactly two processes");
32   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&rank);CHKERRQ(ierr);
33 
34   ierr  = PetscMalloc1(3,&ii);CHKERRQ(ierr);
35   ierr  = PetscMalloc1(6,&jj);CHKERRQ(ierr);
36   ii[0] = 0; ii[1] = 3; ii[2] = 6;
37   if (!rank) {
38     jj[0] = 0; jj[1] = 1; jj[2] = 2; jj[3] = 1; jj[4] = 2; jj[5] = 3;
39   } else {
40     jj[0] = 1; jj[1] = 4; jj[2] = 5; jj[3] = 1; jj[4] = 3; jj[5] = 5;
41   }
42   ierr = MatCreateMPIAdj(MPI_COMM_WORLD,ncells,Nvertices,ii,jj,NULL,&mesh);CHKERRQ(ierr);
43   ierr = MatMeshToCellGraph(mesh,2,&dual);CHKERRQ(ierr);
44   ierr = MatView(dual,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
45 
46   ierr = MatPartitioningCreate(MPI_COMM_WORLD,&part);CHKERRQ(ierr);
47   ierr = MatPartitioningSetAdjacency(part,dual);CHKERRQ(ierr);
48   ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr);
49   ierr = MatPartitioningApply(part,&is);CHKERRQ(ierr);
50   ierr = ISView(is,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
51   ierr = ISDestroy(&is);CHKERRQ(ierr);
52   ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr);
53 
54   ierr = MatDestroy(&mesh);CHKERRQ(ierr);
55   ierr = MatDestroy(&dual);CHKERRQ(ierr);
56   ierr = PetscFinalize();
57   return ierr;
58 }
59 
60 
61 
62 
63 
64 /*TEST
65 
66    build:
67      requires: parmetis
68 
69    test:
70       nsize: 2
71       requires: parmetis
72 
73 TEST*/
74