Actual source code: ex15.c

petsc-3.3-p7 2013-05-11
  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>

 13: int main(int argc, char **args)
 14: {
 15:   Mat             A;
 16:   MatPartitioning part;
 17:   IS              is;
 18:   PetscInt        r,N = 10, start, end;
 19:   PetscErrorCode  ierr;
 20: 
 21:   PetscInitialize(&argc, &args, (char *) 0, help);
 22:   PetscOptionsGetInt(PETSC_NULL, "-N", &N, PETSC_NULL);
 23:   MatCreate(PETSC_COMM_WORLD, &A);
 24:   MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N);
 25:   MatSetFromOptions(A);
 26:   MatSeqAIJSetPreallocation(A, 3, PETSC_NULL);
 27:   MatMPIAIJSetPreallocation(A, 3, PETSC_NULL, 2, PETSC_NULL);

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

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

 43:       cols[0] = r-1; cols[1] = r;
 44:       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;
 52:       MatSetValues(A, 1, &r, 3, cols, vals, INSERT_VALUES);
 53:     }
 54:   }
 55:   MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
 56:   MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);

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

 68:   MatDestroy(&A);
 69:   PetscFinalize();
 70:   return 0;
 71: }