Actual source code: ex83.c
petsc-3.7.7 2017-09-25
2: static char help[] = "Partition tiny grid using hierarchical partitioning and increase overlap using MatIncreaseOverlapSplit.\n\n";
4: /*T
5: Concepts: partitioning
6: Processors: 4
7: T*/
9: /*
10: Include "petscmat.h" so that we can use matrices. Note that this file
11: automatically includes:
12: petscsys.h - base PETSc routines petscvec.h - vectors
13: petscmat.h - matrices
14: petscis.h - index sets
15: petscviewer.h - viewers
16: */
17: #include <petscmat.h>
21: int main(int argc,char **args)
22: {
23: Mat A,B;
24: PetscErrorCode ierr;
25: PetscMPIInt rank,size,membershipKey;
26: PetscInt *ia,*ja,*indices_sc,isrows_localsize;
27: const PetscInt *indices;
28: MatPartitioning part;
29: IS is,isrows,isrows_sc;
30: IS coarseparts,fineparts;
31: MPI_Comm comm,scomm;
33: PetscInitialize(&argc,&args,(char*)0,help);
34: comm = PETSC_COMM_WORLD;
35: MPI_Comm_size(comm,&size);
36: if (size != 4) SETERRQ(comm,1,"Must run with 4 processors \n");
37: MPI_Comm_rank(comm,&rank);
38: /*set a small matrix */
39: PetscMalloc1(5,&ia);
40: PetscMalloc1(16,&ja);
41: if (!rank) {
42: ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6;
43: ja[8] = 2; ja[9] = 7;
44: ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
45: membershipKey = 0;
46: } else if (rank == 1) {
47: ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2;
48: ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11;
49: ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
50: membershipKey = 0;
51: } else if (rank == 2) {
52: ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6;
53: ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15;
54: ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
55: membershipKey = 1;
56: } else {
57: ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15;
58: ja[8] = 11; ja[9] = 14;
59: ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
60: membershipKey = 1;
61: }
62: MatCreateMPIAdj(comm,4,16,ia,ja,NULL,&A);
63: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
64: /*
65: Partition the graph of the matrix
66: */
67: MatPartitioningCreate(comm,&part);
68: MatPartitioningSetAdjacency(part,A);
69: MatPartitioningSetType(part,MATPARTITIONINGHIERARCH);
70: MatPartitioningHierarchicalSetNcoarseparts(part,2);
71: MatPartitioningHierarchicalSetNfineparts(part,2);
72: MatPartitioningSetFromOptions(part);
73: /* get new processor owner number of each vertex */
74: MatPartitioningApply(part,&is);
75: /* coarse parts */
76: MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts);
77: ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD);
78: /* fine parts */
79: MatPartitioningHierarchicalGetFineparts(part,&fineparts);
80: ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD);
81: /* partitioning */
82: ISView(is,PETSC_VIEWER_STDOUT_WORLD);
83: /* compute coming rows */
84: ISBuildTwoSided(is,NULL,&isrows);
85: ISView(isrows,PETSC_VIEWER_STDOUT_WORLD);
86: /*create a sub-communicator */
87: MPI_Comm_split(comm, membershipKey,rank,&scomm);
88: ISGetLocalSize(isrows,&isrows_localsize);
89: PetscCalloc1(isrows_localsize,&indices_sc);
90: ISGetIndices(isrows,&indices);
91: PetscMemcpy(indices_sc,indices,sizeof(PetscInt)*isrows_localsize);
92: ISRestoreIndices(isrows,&indices);
93: ISDestroy(&is);
94: ISDestroy(&coarseparts);
95: ISDestroy(&fineparts);
96: ISDestroy(&isrows);
97: MatPartitioningDestroy(&part);
98: /*create a sub-IS on the sub communicator */
99: ISCreateGeneral(scomm,isrows_localsize,indices_sc,PETSC_OWN_POINTER,&isrows_sc);
100: MatConvert(A,MATMPIAIJ,MAT_INITIAL_MATRIX,&B);
101: #if 1
102: MatView(B,PETSC_VIEWER_STDOUT_WORLD);
103: #endif
104: /*increase overlap */
105: MatIncreaseOverlapSplit(B,1,&isrows_sc,1);
106: ISView(isrows_sc,PETSC_NULL);
107: ISDestroy(&isrows_sc);
108: /*
109: Free work space. All PETSc objects should be destroyed when they
110: are no longer needed.
111: */
112: MatDestroy(&A);
113: MatDestroy(&B);
114: PetscFinalize();
115: return 0;
116: }