Actual source code: ex80.c
1: /*$Id: ex80.c,v 1.10 2001/04/10 19:35:44 bsmith Exp $*/
3: static char help[] = "Partition tiny grid.nn";
5: /*T
6: Concepts: partitioning
7: Processors: 4
8: T*/
10: /*
11: Include "petscmat.h" so that we can use matrices. Note that this file
12: automatically includes:
13: petsc.h - base PETSc routines petscvec.h - vectors
14: petscsys.h - system routines petscmat.h - matrices
15: petscis.h - index sets
16: petscviewer.h - viewers
17: */
18: #include petscsles.h
20: int main(int argc,char **args)
21: {
22: Mat A;
23: int ierr,rank,*ia,*ja,size;
24: MatPartitioning part;
25: IS is,isn;
27: PetscInitialize(&argc,&args,(char *)0,help);
28: MPI_Comm_size(PETSC_COMM_WORLD,&size);
29: if (size != 4) SETERRQ(1,"Must run with 4 processors");
30: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
32: PetscMalloc(5*sizeof(int),&ia);
33: PetscMalloc(16*sizeof(int),&ja);
34: if (rank == 0) {
35: ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6;
36: ja[8] = 2; ja[9] = 7;
37: ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
38: } else if (rank == 1) {
39: ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2;
40: ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11;
41: ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
42: } else if (rank == 2) {
43: ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6;
44: ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15;
45: ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
46: } else {
47: ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15;
48: ja[8] = 11; ja[9] = 14;
49: ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
50: }
52: MatCreateMPIAdj(PETSC_COMM_WORLD,4,16,ia,ja,PETSC_NULL,&A);
53: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
55: /*
56: Partition the graph of the matrix
57: */
58: MatPartitioningCreate(PETSC_COMM_WORLD,&part);
59: MatPartitioningSetAdjacency(part,A);
60: MatPartitioningSetFromOptions(part);
61: /* get new processor owner number of each vertex */
62: MatPartitioningApply(part,&is);
63: /* get new global number of each old global number */
64: ISPartitioningToNumbering(is,&isn);
65: ISView(isn,PETSC_VIEWER_STDOUT_WORLD);
66: ISDestroy(is);
68: ISDestroy(isn);
69: MatPartitioningDestroy(part);
71: /*
72: Free work space. All PETSc objects should be destroyed when they
73: are no longer needed.
74: */
75: MatDestroy(A);
78: PetscFinalize();
79: return 0;
80: }