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