xref: /petsc/src/mat/tutorials/ex17.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 static char help[] = "Example of using graph partitioning with a matrix in which some procs have empty ownership\n\n";
2 
3 /*T
4    Concepts: Mat^mat partitioning
5    Concepts: Mat^image segmentation
6    Processors: n
7 T*/
8 
9 #include <petscmat.h>
10 
11 int main(int argc, char **args)
12 {
13   Mat             A;
14   MatPartitioning part;
15   IS              is;
16   PetscInt        i,m,N,rstart,rend,nemptyranks,*emptyranks,nbigranks,*bigranks;
17   PetscMPIInt     rank,size;
18   PetscErrorCode  ierr;
19 
20   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
21   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
23 
24   nemptyranks = 10;
25   nbigranks   = 10;
26   CHKERRQ(PetscMalloc2(nemptyranks,&emptyranks,nbigranks,&bigranks));
27 
28   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Partitioning example options",NULL);CHKERRQ(ierr);
29   CHKERRQ(PetscOptionsIntArray("-emptyranks","Ranks to be skipped by partition","",emptyranks,&nemptyranks,NULL));
30   CHKERRQ(PetscOptionsIntArray("-bigranks","Ranks to be overloaded","",bigranks,&nbigranks,NULL));
31   ierr = PetscOptionsEnd();CHKERRQ(ierr);
32 
33   m = 1;
34   for (i=0; i<nemptyranks; i++) if (rank == emptyranks[i]) m = 0;
35   for (i=0; i<nbigranks; i++) if (rank == bigranks[i]) m = 5;
36 
37   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
38   CHKERRQ(MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE));
39   CHKERRQ(MatSetFromOptions(A));
40   CHKERRQ(MatSeqAIJSetPreallocation(A,3,NULL));
41   CHKERRQ(MatMPIAIJSetPreallocation(A,3,NULL,2,NULL));
42   CHKERRQ(MatSeqBAIJSetPreallocation(A,1,3,NULL));
43   CHKERRQ(MatMPIBAIJSetPreallocation(A,1,3,NULL,2,NULL));
44   CHKERRQ(MatSeqSBAIJSetPreallocation(A,1,2,NULL));
45   CHKERRQ(MatMPISBAIJSetPreallocation(A,1,2,NULL,1,NULL));
46 
47   CHKERRQ(MatGetSize(A,NULL,&N));
48   CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
49   for (i=rstart; i<rend; i++) {
50     const PetscInt    cols[] = {(i+N-1)%N,i,(i+1)%N};
51     const PetscScalar vals[] = {1,1,1};
52     CHKERRQ(MatSetValues(A,1,&i,3,cols,vals,INSERT_VALUES));
53   }
54   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
55   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
56   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
57 
58   CHKERRQ(MatPartitioningCreate(PETSC_COMM_WORLD,&part));
59   CHKERRQ(MatPartitioningSetAdjacency(part,A));
60   CHKERRQ(MatPartitioningSetFromOptions(part));
61   CHKERRQ(MatPartitioningApply(part,&is));
62   CHKERRQ(ISView(is,PETSC_VIEWER_STDOUT_WORLD));
63   CHKERRQ(ISDestroy(&is));
64   CHKERRQ(MatPartitioningDestroy(&part));
65   CHKERRQ(MatDestroy(&A));
66   CHKERRQ(PetscFree2(emptyranks,bigranks));
67   ierr = PetscFinalize();
68   return ierr;
69 }
70 
71 /*TEST
72 
73    test:
74       nsize: 8
75       args: -emptyranks 0,2,4 -bigranks 1,3,7 -mat_partitioning_type average
76       # cannot test with external package partitioners since they produce different results on different systems
77 
78 TEST*/
79