Actual source code: ex17.c

  1: static char help[] = "Example of using graph partitioning with a matrix in which some procs have empty ownership\n\n";

  3: #include <petscmat.h>

  5: int main(int argc, char **args)
  6: {
  7:   Mat             A;
  8:   MatPartitioning part;
  9:   IS              is;
 10:   PetscInt        i, m, N, rstart, rend, nemptyranks, *emptyranks, nbigranks, *bigranks;
 11:   PetscMPIInt     rank, size;

 14:   PetscInitialize(&argc, &args, (char *)0, help);
 15:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
 16:   MPI_Comm_size(PETSC_COMM_WORLD, &size);

 18:   nemptyranks = 10;
 19:   nbigranks   = 10;
 20:   PetscMalloc2(nemptyranks, &emptyranks, nbigranks, &bigranks);

 22:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Partitioning example options", NULL);
 23:   PetscOptionsIntArray("-emptyranks", "Ranks to be skipped by partition", "", emptyranks, &nemptyranks, NULL);
 24:   PetscOptionsIntArray("-bigranks", "Ranks to be overloaded", "", bigranks, &nbigranks, NULL);
 25:   PetscOptionsEnd();

 27:   m = 1;
 28:   for (i = 0; i < nemptyranks; i++)
 29:     if (rank == emptyranks[i]) m = 0;
 30:   for (i = 0; i < nbigranks; i++)
 31:     if (rank == bigranks[i]) m = 5;

 33:   MatCreate(PETSC_COMM_WORLD, &A);
 34:   MatSetSizes(A, m, m, PETSC_DECIDE, PETSC_DECIDE);
 35:   MatSetFromOptions(A);
 36:   MatSeqAIJSetPreallocation(A, 3, NULL);
 37:   MatMPIAIJSetPreallocation(A, 3, NULL, 2, NULL);
 38:   MatSeqBAIJSetPreallocation(A, 1, 3, NULL);
 39:   MatMPIBAIJSetPreallocation(A, 1, 3, NULL, 2, NULL);
 40:   MatSeqSBAIJSetPreallocation(A, 1, 2, NULL);
 41:   MatMPISBAIJSetPreallocation(A, 1, 2, NULL, 1, NULL);

 43:   MatGetSize(A, NULL, &N);
 44:   MatGetOwnershipRange(A, &rstart, &rend);
 45:   for (i = rstart; i < rend; i++) {
 46:     const PetscInt    cols[] = {(i + N - 1) % N, i, (i + 1) % N};
 47:     const PetscScalar vals[] = {1, 1, 1};
 48:     MatSetValues(A, 1, &i, 3, cols, vals, INSERT_VALUES);
 49:   }
 50:   MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
 51:   MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
 52:   MatView(A, PETSC_VIEWER_STDOUT_WORLD);

 54:   MatPartitioningCreate(PETSC_COMM_WORLD, &part);
 55:   MatPartitioningSetAdjacency(part, A);
 56:   MatPartitioningSetFromOptions(part);
 57:   MatPartitioningApply(part, &is);
 58:   ISView(is, PETSC_VIEWER_STDOUT_WORLD);
 59:   ISDestroy(&is);
 60:   MatPartitioningDestroy(&part);
 61:   MatDestroy(&A);
 62:   PetscFree2(emptyranks, bigranks);
 63:   PetscFinalize();
 64:   return 0;
 65: }

 67: /*TEST

 69:    test:
 70:       nsize: 8
 71:       args: -emptyranks 0,2,4 -bigranks 1,3,7 -mat_partitioning_type average
 72:       # cannot test with external package partitioners since they produce different results on different systems

 74: TEST*/