Actual source code: ex15.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: static char help[] = "Example of using graph partitioning to segment an image\n\n";

  3: /*T
  4:    Concepts: Mat^mat partitioning
  5:    Concepts: Mat^image segmentation
  6:    Processors: n
  7: T*/

  9:  #include <petscmat.h>

 11: int main(int argc, char **args)
 12: {
 13:   Mat             A;
 14:   MatPartitioning part;
 15:   IS              is;
 16:   PetscInt        r,N = 10, start, end;
 17:   PetscErrorCode  ierr;

 19:   PetscInitialize(&argc, &args, (char*) 0, help);
 20:   PetscOptionsGetInt(NULL,NULL, "-N", &N, NULL);
 21:   MatCreate(PETSC_COMM_WORLD, &A);
 22:   MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N);
 23:   MatSetFromOptions(A);
 24:   MatSeqAIJSetPreallocation(A, 3, NULL);
 25:   MatMPIAIJSetPreallocation(A, 3, NULL, 2, NULL);

 27:   /* Create a linear mesh */
 28:   MatGetOwnershipRange(A, &start, &end);
 29:   for (r = start; r < end; ++r) {
 30:     if (r == 0) {
 31:       PetscInt    cols[2];
 32:       PetscScalar vals[2];

 34:       cols[0] = r;   cols[1] = r+1;
 35:       vals[0] = 1.0; vals[1] = 1.0;

 37:       MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES);
 38:     } else if (r == N-1) {
 39:       PetscInt    cols[2];
 40:       PetscScalar vals[2];

 42:       cols[0] = r-1; cols[1] = r;
 43:       vals[0] = 1.0; vals[1] = 1.0;

 45:       MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES);
 46:     } else {
 47:       PetscInt    cols[3];
 48:       PetscScalar vals[3];

 50:       cols[0] = r-1; cols[1] = r;   cols[2] = r+1;
 51:       vals[0] = 1.0; vals[1] = 1.0; vals[2] = 1.0;

 53:       MatSetValues(A, 1, &r, 3, cols, vals, INSERT_VALUES);
 54:     }
 55:   }
 56:   MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
 57:   MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);

 59:   MatPartitioningCreate(PETSC_COMM_WORLD, &part);
 60:   MatPartitioningSetAdjacency(part, A);
 61:   MatPartitioningSetFromOptions(part);
 62:   /*MatPartitioningSetVertexWeights(part, const PetscInt weights[]);*/
 63:   /*MatPartitioningSetPartitionWeights(part,const PetscReal weights[]);*/
 64:   MatPartitioningApply(part, &is);
 65:   ISView(is, PETSC_VIEWER_STDOUT_WORLD);
 66:   ISDestroy(&is);
 67:   MatPartitioningDestroy(&part);

 69:   MatDestroy(&A);
 70:   PetscFinalize();
 71:   return ierr;
 72: }


 75: /*TEST

 77:    test:
 78:       nsize: 3
 79:       requires: parmetis
 80:       args: -mat_partitioning_type parmetis

 82:    test:
 83:       suffix: 2
 84:       nsize: 3
 85:       requires: ptscotch
 86:       args: -mat_partitioning_type ptscotch

 88:    test:
 89:       suffix: 3
 90:       nsize: 4
 91:       requires: party
 92:       args: -mat_partitioning_type party

 94:    test:
 95:       suffix: 4
 96:       nsize: 3
 97:       requires: chaco
 98:       args: -mat_partitioning_type chaco

100: TEST*/