xref: /petsc/src/mat/tutorials/ex17.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Example of using graph partitioning with a matrix in which some procs have empty ownership\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5*9371c9d4SSatish Balay int main(int argc, char **args) {
6c4762a1bSJed Brown   Mat             A;
7c4762a1bSJed Brown   MatPartitioning part;
8c4762a1bSJed Brown   IS              is;
9c4762a1bSJed Brown   PetscInt        i, m, N, rstart, rend, nemptyranks, *emptyranks, nbigranks, *bigranks;
10c4762a1bSJed Brown   PetscMPIInt     rank, size;
11c4762a1bSJed Brown 
12327415f7SBarry Smith   PetscFunctionBeginUser;
139566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   nemptyranks = 10;
18c4762a1bSJed Brown   nbigranks   = 10;
199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nemptyranks, &emptyranks, nbigranks, &bigranks));
20c4762a1bSJed Brown 
21d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Partitioning example options", NULL);
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-emptyranks", "Ranks to be skipped by partition", "", emptyranks, &nemptyranks, NULL));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-bigranks", "Ranks to be overloaded", "", bigranks, &nbigranks, NULL));
24d0609cedSBarry Smith   PetscOptionsEnd();
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   m = 1;
27*9371c9d4SSatish Balay   for (i = 0; i < nemptyranks; i++)
28*9371c9d4SSatish Balay     if (rank == emptyranks[i]) m = 0;
29*9371c9d4SSatish Balay   for (i = 0; i < nbigranks; i++)
30*9371c9d4SSatish Balay     if (rank == bigranks[i]) m = 5;
31c4762a1bSJed Brown 
329566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
339566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, m, m, PETSC_DECIDE, PETSC_DECIDE));
349566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
359566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A, 3, NULL));
369566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A, 3, NULL, 2, NULL));
379566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(A, 1, 3, NULL));
389566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(A, 1, 3, NULL, 2, NULL));
399566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(A, 1, 2, NULL));
409566063dSJacob Faibussowitsch   PetscCall(MatMPISBAIJSetPreallocation(A, 1, 2, NULL, 1, NULL));
41c4762a1bSJed Brown 
429566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &N));
439566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
44c4762a1bSJed Brown   for (i = rstart; i < rend; i++) {
45c4762a1bSJed Brown     const PetscInt    cols[] = {(i + N - 1) % N, i, (i + 1) % N};
46c4762a1bSJed Brown     const PetscScalar vals[] = {1, 1, 1};
479566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, cols, vals, INSERT_VALUES));
48c4762a1bSJed Brown   }
499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
519566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(PETSC_COMM_WORLD, &part));
549566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(part, A));
559566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(part));
569566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(part, &is));
579566063dSJacob Faibussowitsch   PetscCall(ISView(is, PETSC_VIEWER_STDOUT_WORLD));
589566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
599566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&part));
609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
619566063dSJacob Faibussowitsch   PetscCall(PetscFree2(emptyranks, bigranks));
629566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
63b122ec5aSJacob Faibussowitsch   return 0;
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /*TEST
67c4762a1bSJed Brown 
68c4762a1bSJed Brown    test:
69c4762a1bSJed Brown       nsize: 8
70c4762a1bSJed Brown       args: -emptyranks 0,2,4 -bigranks 1,3,7 -mat_partitioning_type average
71c4762a1bSJed Brown       # cannot test with external package partitioners since they produce different results on different systems
72c4762a1bSJed Brown 
73c4762a1bSJed Brown TEST*/
74