1abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/ 2abe9303eSLisandro Dalcin 3abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_PTSCOTCH) 4abe9303eSLisandro Dalcin EXTERN_C_BEGIN 5abe9303eSLisandro Dalcin #include <ptscotch.h> 6abe9303eSLisandro Dalcin EXTERN_C_END 7abe9303eSLisandro Dalcin #endif 8abe9303eSLisandro Dalcin 9abe9303eSLisandro Dalcin PetscBool PTScotchPartitionerCite = PETSC_FALSE; 10abe9303eSLisandro Dalcin const char PTScotchPartitionerCitation[] = 11abe9303eSLisandro Dalcin "@article{PTSCOTCH,\n" 12abe9303eSLisandro Dalcin " author = {C. Chevalier and F. Pellegrini},\n" 13abe9303eSLisandro Dalcin " title = {{PT-SCOTCH}: a tool for efficient parallel graph ordering},\n" 14abe9303eSLisandro Dalcin " journal = {Parallel Computing},\n" 15abe9303eSLisandro Dalcin " volume = {34},\n" 16abe9303eSLisandro Dalcin " number = {6},\n" 17abe9303eSLisandro Dalcin " pages = {318--331},\n" 18abe9303eSLisandro Dalcin " year = {2008},\n" 19abe9303eSLisandro Dalcin " doi = {https://doi.org/10.1016/j.parco.2007.12.001}\n" 20abe9303eSLisandro Dalcin "}\n"; 21abe9303eSLisandro Dalcin 22abe9303eSLisandro Dalcin typedef struct { 23abe9303eSLisandro Dalcin MPI_Comm pcomm; 24abe9303eSLisandro Dalcin PetscInt strategy; 25abe9303eSLisandro Dalcin PetscReal imbalance; 26abe9303eSLisandro Dalcin } PetscPartitioner_PTScotch; 27abe9303eSLisandro Dalcin 28abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_PTSCOTCH) 29abe9303eSLisandro Dalcin 30abe9303eSLisandro Dalcin #define CHKERRPTSCOTCH(ierr) do { if (ierr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error calling PT-Scotch library"); } while (0) 31abe9303eSLisandro Dalcin 32abe9303eSLisandro Dalcin static int PTScotch_Strategy(PetscInt strategy) 33abe9303eSLisandro Dalcin { 34abe9303eSLisandro Dalcin switch (strategy) { 35abe9303eSLisandro Dalcin case 0: return SCOTCH_STRATDEFAULT; 36abe9303eSLisandro Dalcin case 1: return SCOTCH_STRATQUALITY; 37abe9303eSLisandro Dalcin case 2: return SCOTCH_STRATSPEED; 38abe9303eSLisandro Dalcin case 3: return SCOTCH_STRATBALANCE; 39abe9303eSLisandro Dalcin case 4: return SCOTCH_STRATSAFETY; 40abe9303eSLisandro Dalcin case 5: return SCOTCH_STRATSCALABILITY; 41abe9303eSLisandro Dalcin case 6: return SCOTCH_STRATRECURSIVE; 42abe9303eSLisandro Dalcin case 7: return SCOTCH_STRATREMAP; 43abe9303eSLisandro Dalcin default: return SCOTCH_STRATDEFAULT; 44abe9303eSLisandro Dalcin } 45abe9303eSLisandro Dalcin } 46abe9303eSLisandro Dalcin 47abe9303eSLisandro Dalcin static PetscErrorCode PTScotch_PartGraph_Seq(SCOTCH_Num strategy, double imbalance, SCOTCH_Num n, SCOTCH_Num xadj[], SCOTCH_Num adjncy[], 48abe9303eSLisandro Dalcin SCOTCH_Num vtxwgt[], SCOTCH_Num adjwgt[], SCOTCH_Num nparts, SCOTCH_Num tpart[], SCOTCH_Num part[]) 49abe9303eSLisandro Dalcin { 50*d7cc930eSLisandro Dalcin SCOTCH_Arch archdat; 51abe9303eSLisandro Dalcin SCOTCH_Graph grafdat; 52abe9303eSLisandro Dalcin SCOTCH_Strat stradat; 53abe9303eSLisandro Dalcin SCOTCH_Num vertnbr = n; 54abe9303eSLisandro Dalcin SCOTCH_Num edgenbr = xadj[n]; 55abe9303eSLisandro Dalcin SCOTCH_Num* velotab = vtxwgt; 56abe9303eSLisandro Dalcin SCOTCH_Num* edlotab = adjwgt; 57abe9303eSLisandro Dalcin SCOTCH_Num flagval = strategy; 58abe9303eSLisandro Dalcin double kbalval = imbalance; 59abe9303eSLisandro Dalcin PetscErrorCode ierr; 60abe9303eSLisandro Dalcin 61abe9303eSLisandro Dalcin PetscFunctionBegin; 62abe9303eSLisandro Dalcin { 63abe9303eSLisandro Dalcin PetscBool flg = PETSC_TRUE; 64abe9303eSLisandro Dalcin ierr = PetscOptionsDeprecatedNoObject("-petscpartititoner_ptscotch_vertex_weight",NULL,"3.13","Use -petscpartitioner_use_vertex_weights");CHKERRQ(ierr); 65abe9303eSLisandro Dalcin ierr = PetscOptionsGetBool(NULL, NULL, "-petscpartititoner_ptscotch_vertex_weight", &flg, NULL);CHKERRQ(ierr); 66abe9303eSLisandro Dalcin if (!flg) velotab = NULL; 67abe9303eSLisandro Dalcin } 68abe9303eSLisandro Dalcin ierr = SCOTCH_graphInit(&grafdat);CHKERRPTSCOTCH(ierr); 69abe9303eSLisandro Dalcin ierr = SCOTCH_graphBuild(&grafdat, 0, vertnbr, xadj, xadj + 1, velotab, NULL, edgenbr, adjncy, edlotab);CHKERRPTSCOTCH(ierr); 70abe9303eSLisandro Dalcin ierr = SCOTCH_stratInit(&stradat);CHKERRPTSCOTCH(ierr); 71abe9303eSLisandro Dalcin ierr = SCOTCH_stratGraphMapBuild(&stradat, flagval, nparts, kbalval);CHKERRPTSCOTCH(ierr); 72abe9303eSLisandro Dalcin ierr = SCOTCH_archInit(&archdat);CHKERRPTSCOTCH(ierr); 73*d7cc930eSLisandro Dalcin if (tpart) { 74abe9303eSLisandro Dalcin ierr = SCOTCH_archCmpltw(&archdat, nparts, tpart);CHKERRPTSCOTCH(ierr); 75*d7cc930eSLisandro Dalcin } else { 76*d7cc930eSLisandro Dalcin ierr = SCOTCH_archCmplt(&archdat, nparts);CHKERRPTSCOTCH(ierr); 77*d7cc930eSLisandro Dalcin } 78abe9303eSLisandro Dalcin ierr = SCOTCH_graphMap(&grafdat, &archdat, &stradat, part);CHKERRPTSCOTCH(ierr); 79abe9303eSLisandro Dalcin SCOTCH_archExit(&archdat); 80abe9303eSLisandro Dalcin SCOTCH_stratExit(&stradat); 81abe9303eSLisandro Dalcin SCOTCH_graphExit(&grafdat); 82abe9303eSLisandro Dalcin PetscFunctionReturn(0); 83abe9303eSLisandro Dalcin } 84abe9303eSLisandro Dalcin 85abe9303eSLisandro Dalcin static PetscErrorCode PTScotch_PartGraph_MPI(SCOTCH_Num strategy, double imbalance, SCOTCH_Num vtxdist[], SCOTCH_Num xadj[], SCOTCH_Num adjncy[], 86abe9303eSLisandro Dalcin SCOTCH_Num vtxwgt[], SCOTCH_Num adjwgt[], SCOTCH_Num nparts, SCOTCH_Num tpart[], SCOTCH_Num part[], MPI_Comm comm) 87abe9303eSLisandro Dalcin { 88abe9303eSLisandro Dalcin PetscMPIInt procglbnbr; 89abe9303eSLisandro Dalcin PetscMPIInt proclocnum; 90abe9303eSLisandro Dalcin SCOTCH_Arch archdat; 91abe9303eSLisandro Dalcin SCOTCH_Dgraph grafdat; 92abe9303eSLisandro Dalcin SCOTCH_Dmapping mappdat; 93abe9303eSLisandro Dalcin SCOTCH_Strat stradat; 94abe9303eSLisandro Dalcin SCOTCH_Num vertlocnbr; 95abe9303eSLisandro Dalcin SCOTCH_Num edgelocnbr; 96abe9303eSLisandro Dalcin SCOTCH_Num* veloloctab = vtxwgt; 97abe9303eSLisandro Dalcin SCOTCH_Num* edloloctab = adjwgt; 98abe9303eSLisandro Dalcin SCOTCH_Num flagval = strategy; 99abe9303eSLisandro Dalcin double kbalval = imbalance; 100abe9303eSLisandro Dalcin PetscErrorCode ierr; 101abe9303eSLisandro Dalcin 102abe9303eSLisandro Dalcin PetscFunctionBegin; 103abe9303eSLisandro Dalcin { 104abe9303eSLisandro Dalcin PetscBool flg = PETSC_TRUE; 105abe9303eSLisandro Dalcin ierr = PetscOptionsDeprecatedNoObject("-petscpartititoner_ptscotch_vertex_weight",NULL,"3.13","Use -petscpartitioner_use_vertex_weights");CHKERRQ(ierr); 106abe9303eSLisandro Dalcin ierr = PetscOptionsGetBool(NULL, NULL, "-petscpartititoner_ptscotch_vertex_weight", &flg, NULL);CHKERRQ(ierr); 107abe9303eSLisandro Dalcin if (!flg) veloloctab = NULL; 108abe9303eSLisandro Dalcin } 109abe9303eSLisandro Dalcin ierr = MPI_Comm_size(comm, &procglbnbr);CHKERRQ(ierr); 110abe9303eSLisandro Dalcin ierr = MPI_Comm_rank(comm, &proclocnum);CHKERRQ(ierr); 111abe9303eSLisandro Dalcin vertlocnbr = vtxdist[proclocnum + 1] - vtxdist[proclocnum]; 112abe9303eSLisandro Dalcin edgelocnbr = xadj[vertlocnbr]; 113abe9303eSLisandro Dalcin 114abe9303eSLisandro Dalcin ierr = SCOTCH_dgraphInit(&grafdat, comm);CHKERRPTSCOTCH(ierr); 115abe9303eSLisandro Dalcin ierr = SCOTCH_dgraphBuild(&grafdat, 0, vertlocnbr, vertlocnbr, xadj, xadj + 1, veloloctab, NULL, edgelocnbr, edgelocnbr, adjncy, NULL, edloloctab);CHKERRPTSCOTCH(ierr); 116abe9303eSLisandro Dalcin ierr = SCOTCH_stratInit(&stradat);CHKERRPTSCOTCH(ierr); 117abe9303eSLisandro Dalcin ierr = SCOTCH_stratDgraphMapBuild(&stradat, flagval, procglbnbr, nparts, kbalval);CHKERRQ(ierr); 118abe9303eSLisandro Dalcin ierr = SCOTCH_archInit(&archdat);CHKERRPTSCOTCH(ierr); 119abe9303eSLisandro Dalcin if (tpart) { /* target partition weights */ 120abe9303eSLisandro Dalcin ierr = SCOTCH_archCmpltw(&archdat, nparts, tpart);CHKERRPTSCOTCH(ierr); 121abe9303eSLisandro Dalcin } else { 122abe9303eSLisandro Dalcin ierr = SCOTCH_archCmplt(&archdat, nparts);CHKERRPTSCOTCH(ierr); 123abe9303eSLisandro Dalcin } 124abe9303eSLisandro Dalcin ierr = SCOTCH_dgraphMapInit(&grafdat, &mappdat, &archdat, part);CHKERRPTSCOTCH(ierr); 125abe9303eSLisandro Dalcin ierr = SCOTCH_dgraphMapCompute(&grafdat, &mappdat, &stradat);CHKERRPTSCOTCH(ierr); 126abe9303eSLisandro Dalcin SCOTCH_dgraphMapExit(&grafdat, &mappdat); 127abe9303eSLisandro Dalcin SCOTCH_archExit(&archdat); 128abe9303eSLisandro Dalcin SCOTCH_stratExit(&stradat); 129abe9303eSLisandro Dalcin SCOTCH_dgraphExit(&grafdat); 130abe9303eSLisandro Dalcin PetscFunctionReturn(0); 131abe9303eSLisandro Dalcin } 132abe9303eSLisandro Dalcin 133abe9303eSLisandro Dalcin #endif /* PETSC_HAVE_PTSCOTCH */ 134abe9303eSLisandro Dalcin 135abe9303eSLisandro Dalcin static const char *const 136abe9303eSLisandro Dalcin PTScotchStrategyList[] = { 137abe9303eSLisandro Dalcin "DEFAULT", 138abe9303eSLisandro Dalcin "QUALITY", 139abe9303eSLisandro Dalcin "SPEED", 140abe9303eSLisandro Dalcin "BALANCE", 141abe9303eSLisandro Dalcin "SAFETY", 142abe9303eSLisandro Dalcin "SCALABILITY", 143abe9303eSLisandro Dalcin "RECURSIVE", 144abe9303eSLisandro Dalcin "REMAP" 145abe9303eSLisandro Dalcin }; 146abe9303eSLisandro Dalcin 147abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerDestroy_PTScotch(PetscPartitioner part) 148abe9303eSLisandro Dalcin { 149abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *) part->data; 150abe9303eSLisandro Dalcin PetscErrorCode ierr; 151abe9303eSLisandro Dalcin 152abe9303eSLisandro Dalcin PetscFunctionBegin; 153abe9303eSLisandro Dalcin ierr = MPI_Comm_free(&p->pcomm);CHKERRQ(ierr); 154abe9303eSLisandro Dalcin ierr = PetscFree(part->data);CHKERRQ(ierr); 155abe9303eSLisandro Dalcin PetscFunctionReturn(0); 156abe9303eSLisandro Dalcin } 157abe9303eSLisandro Dalcin 158abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_PTScotch_ASCII(PetscPartitioner part, PetscViewer viewer) 159abe9303eSLisandro Dalcin { 160abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *) part->data; 161abe9303eSLisandro Dalcin PetscErrorCode ierr; 162abe9303eSLisandro Dalcin 163abe9303eSLisandro Dalcin PetscFunctionBegin; 164abe9303eSLisandro Dalcin ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 165abe9303eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"using partitioning strategy %s\n",PTScotchStrategyList[p->strategy]);CHKERRQ(ierr); 166abe9303eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"using load imbalance ratio %g\n",(double)p->imbalance);CHKERRQ(ierr); 167abe9303eSLisandro Dalcin ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 168abe9303eSLisandro Dalcin PetscFunctionReturn(0); 169abe9303eSLisandro Dalcin } 170abe9303eSLisandro Dalcin 171abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_PTScotch(PetscPartitioner part, PetscViewer viewer) 172abe9303eSLisandro Dalcin { 173abe9303eSLisandro Dalcin PetscBool iascii; 174abe9303eSLisandro Dalcin PetscErrorCode ierr; 175abe9303eSLisandro Dalcin 176abe9303eSLisandro Dalcin PetscFunctionBegin; 177abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 178abe9303eSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 179abe9303eSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 180abe9303eSLisandro Dalcin if (iascii) {ierr = PetscPartitionerView_PTScotch_ASCII(part, viewer);CHKERRQ(ierr);} 181abe9303eSLisandro Dalcin PetscFunctionReturn(0); 182abe9303eSLisandro Dalcin } 183abe9303eSLisandro Dalcin 184abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerSetFromOptions_PTScotch(PetscOptionItems *PetscOptionsObject, PetscPartitioner part) 185abe9303eSLisandro Dalcin { 186abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *) part->data; 187abe9303eSLisandro Dalcin const char *const *slist = PTScotchStrategyList; 188abe9303eSLisandro Dalcin PetscInt nlist = (PetscInt)(sizeof(PTScotchStrategyList)/sizeof(PTScotchStrategyList[0])); 189abe9303eSLisandro Dalcin PetscBool flag; 190abe9303eSLisandro Dalcin PetscErrorCode ierr; 191abe9303eSLisandro Dalcin 192abe9303eSLisandro Dalcin PetscFunctionBegin; 193abe9303eSLisandro Dalcin ierr = PetscOptionsHead(PetscOptionsObject, "PetscPartitioner PTScotch Options");CHKERRQ(ierr); 194abe9303eSLisandro Dalcin ierr = PetscOptionsEList("-petscpartitioner_ptscotch_strategy","Partitioning strategy","",slist,nlist,slist[p->strategy],&p->strategy,&flag);CHKERRQ(ierr); 195abe9303eSLisandro Dalcin ierr = PetscOptionsReal("-petscpartitioner_ptscotch_imbalance","Load imbalance ratio","",p->imbalance,&p->imbalance,&flag);CHKERRQ(ierr); 196abe9303eSLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 197abe9303eSLisandro Dalcin PetscFunctionReturn(0); 198abe9303eSLisandro Dalcin } 199abe9303eSLisandro Dalcin 200abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerPartition_PTScotch(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) 201abe9303eSLisandro Dalcin { 202abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_PTSCOTCH) 203abe9303eSLisandro Dalcin MPI_Comm comm; 204abe9303eSLisandro Dalcin PetscInt nvtxs = numVertices; /* The number of vertices in full graph */ 205abe9303eSLisandro Dalcin PetscInt *vtxdist; /* Distribution of vertices across processes */ 206abe9303eSLisandro Dalcin PetscInt *xadj = start; /* Start of edge list for each vertex */ 207abe9303eSLisandro Dalcin PetscInt *adjncy = adjacency; /* Edge lists for all vertices */ 208abe9303eSLisandro Dalcin PetscInt *vwgt = NULL; /* Vertex weights */ 209abe9303eSLisandro Dalcin PetscInt *adjwgt = NULL; /* Edge weights */ 210abe9303eSLisandro Dalcin PetscInt v, i, *assignment, *points; 211abe9303eSLisandro Dalcin PetscMPIInt size, rank, p; 212abe9303eSLisandro Dalcin PetscBool hasempty = PETSC_FALSE; 213abe9303eSLisandro Dalcin PetscInt *tpwgts = NULL; 214abe9303eSLisandro Dalcin PetscErrorCode ierr; 215abe9303eSLisandro Dalcin 216abe9303eSLisandro Dalcin PetscFunctionBegin; 217abe9303eSLisandro Dalcin ierr = PetscObjectGetComm((PetscObject)part,&comm);CHKERRQ(ierr); 218abe9303eSLisandro Dalcin ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 219abe9303eSLisandro Dalcin ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 220abe9303eSLisandro Dalcin ierr = PetscMalloc2(size+1,&vtxdist,PetscMax(nvtxs,1),&assignment);CHKERRQ(ierr); 221abe9303eSLisandro Dalcin /* Calculate vertex distribution */ 222abe9303eSLisandro Dalcin vtxdist[0] = 0; 223abe9303eSLisandro Dalcin ierr = MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm);CHKERRQ(ierr); 224abe9303eSLisandro Dalcin for (p = 2; p <= size; ++p) { 225abe9303eSLisandro Dalcin hasempty = (PetscBool)(hasempty || !vtxdist[p-1] || !vtxdist[p]); 226abe9303eSLisandro Dalcin vtxdist[p] += vtxdist[p-1]; 227abe9303eSLisandro Dalcin } 228abe9303eSLisandro Dalcin /* null graph */ 229abe9303eSLisandro Dalcin if (vtxdist[size] == 0) { 230abe9303eSLisandro Dalcin ierr = PetscFree2(vtxdist, assignment);CHKERRQ(ierr); 231abe9303eSLisandro Dalcin ierr = ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition);CHKERRQ(ierr); 232abe9303eSLisandro Dalcin PetscFunctionReturn(0); 233abe9303eSLisandro Dalcin } 234abe9303eSLisandro Dalcin 235abe9303eSLisandro Dalcin /* Calculate vertex weights */ 236abe9303eSLisandro Dalcin if (vertSection) { 237abe9303eSLisandro Dalcin ierr = PetscMalloc1(nvtxs,&vwgt);CHKERRQ(ierr); 238abe9303eSLisandro Dalcin for (v = 0; v < nvtxs; ++v) { 239abe9303eSLisandro Dalcin ierr = PetscSectionGetDof(vertSection, v, &vwgt[v]);CHKERRQ(ierr); 240abe9303eSLisandro Dalcin } 241abe9303eSLisandro Dalcin } 242abe9303eSLisandro Dalcin 243abe9303eSLisandro Dalcin /* Calculate partition weights */ 244abe9303eSLisandro Dalcin if (targetSection) { 245abe9303eSLisandro Dalcin PetscInt sumw; 246abe9303eSLisandro Dalcin 247abe9303eSLisandro Dalcin ierr = PetscCalloc1(nparts,&tpwgts);CHKERRQ(ierr); 248abe9303eSLisandro Dalcin for (p = 0, sumw = 0; p < nparts; ++p) { 249abe9303eSLisandro Dalcin ierr = PetscSectionGetDof(targetSection,p,&tpwgts[p]);CHKERRQ(ierr); 250abe9303eSLisandro Dalcin sumw += tpwgts[p]; 251abe9303eSLisandro Dalcin } 252*d7cc930eSLisandro Dalcin if (!sumw) {ierr = PetscFree(tpwgts);CHKERRQ(ierr);} 253abe9303eSLisandro Dalcin } 254abe9303eSLisandro Dalcin 255abe9303eSLisandro Dalcin { 256abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *pts = (PetscPartitioner_PTScotch *) part->data; 257abe9303eSLisandro Dalcin int strat = PTScotch_Strategy(pts->strategy); 258abe9303eSLisandro Dalcin double imbal = (double)pts->imbalance; 259abe9303eSLisandro Dalcin 260abe9303eSLisandro Dalcin for (p = 0; !vtxdist[p+1] && p < size; ++p); 261abe9303eSLisandro Dalcin if (vtxdist[p+1] == vtxdist[size]) { 262abe9303eSLisandro Dalcin if (rank == p) { 263abe9303eSLisandro Dalcin ierr = PTScotch_PartGraph_Seq(strat, imbal, nvtxs, xadj, adjncy, vwgt, adjwgt, nparts, tpwgts, assignment);CHKERRQ(ierr); 264abe9303eSLisandro Dalcin } 265abe9303eSLisandro Dalcin } else { 266abe9303eSLisandro Dalcin MPI_Comm pcomm = pts->pcomm; 267abe9303eSLisandro Dalcin 268abe9303eSLisandro Dalcin if (hasempty) { 269abe9303eSLisandro Dalcin PetscInt cnt; 270abe9303eSLisandro Dalcin 271abe9303eSLisandro Dalcin ierr = MPI_Comm_split(pts->pcomm,!!nvtxs,rank,&pcomm);CHKERRQ(ierr); 272abe9303eSLisandro Dalcin for (p=0,cnt=0;p<size;p++) { 273abe9303eSLisandro Dalcin if (vtxdist[p+1] != vtxdist[p]) { 274abe9303eSLisandro Dalcin vtxdist[cnt+1] = vtxdist[p+1]; 275abe9303eSLisandro Dalcin cnt++; 276abe9303eSLisandro Dalcin } 277abe9303eSLisandro Dalcin } 278abe9303eSLisandro Dalcin }; 279abe9303eSLisandro Dalcin if (nvtxs) { 280abe9303eSLisandro Dalcin ierr = PTScotch_PartGraph_MPI(strat, imbal, vtxdist, xadj, adjncy, vwgt, adjwgt, nparts, tpwgts, assignment, pcomm);CHKERRQ(ierr); 281abe9303eSLisandro Dalcin } 282abe9303eSLisandro Dalcin if (hasempty) { 283abe9303eSLisandro Dalcin ierr = MPI_Comm_free(&pcomm);CHKERRQ(ierr); 284abe9303eSLisandro Dalcin } 285abe9303eSLisandro Dalcin } 286abe9303eSLisandro Dalcin } 287abe9303eSLisandro Dalcin ierr = PetscFree(vwgt);CHKERRQ(ierr); 288abe9303eSLisandro Dalcin ierr = PetscFree(tpwgts);CHKERRQ(ierr); 289abe9303eSLisandro Dalcin 290abe9303eSLisandro Dalcin /* Convert to PetscSection+IS */ 291abe9303eSLisandro Dalcin for (v = 0; v < nvtxs; ++v) {ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);} 292abe9303eSLisandro Dalcin ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr); 293abe9303eSLisandro Dalcin for (p = 0, i = 0; p < nparts; ++p) { 294abe9303eSLisandro Dalcin for (v = 0; v < nvtxs; ++v) { 295abe9303eSLisandro Dalcin if (assignment[v] == p) points[i++] = v; 296abe9303eSLisandro Dalcin } 297abe9303eSLisandro Dalcin } 298abe9303eSLisandro Dalcin if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs); 299abe9303eSLisandro Dalcin ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr); 300abe9303eSLisandro Dalcin 301abe9303eSLisandro Dalcin ierr = PetscFree2(vtxdist,assignment);CHKERRQ(ierr); 302abe9303eSLisandro Dalcin PetscFunctionReturn(0); 303abe9303eSLisandro Dalcin #else 304abe9303eSLisandro Dalcin SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-ptscotch."); 305abe9303eSLisandro Dalcin #endif 306abe9303eSLisandro Dalcin } 307abe9303eSLisandro Dalcin 308abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerInitialize_PTScotch(PetscPartitioner part) 309abe9303eSLisandro Dalcin { 310abe9303eSLisandro Dalcin PetscFunctionBegin; 311abe9303eSLisandro Dalcin part->noGraph = PETSC_FALSE; 312abe9303eSLisandro Dalcin part->ops->view = PetscPartitionerView_PTScotch; 313abe9303eSLisandro Dalcin part->ops->destroy = PetscPartitionerDestroy_PTScotch; 314abe9303eSLisandro Dalcin part->ops->partition = PetscPartitionerPartition_PTScotch; 315abe9303eSLisandro Dalcin part->ops->setfromoptions = PetscPartitionerSetFromOptions_PTScotch; 316abe9303eSLisandro Dalcin PetscFunctionReturn(0); 317abe9303eSLisandro Dalcin } 318abe9303eSLisandro Dalcin 319abe9303eSLisandro Dalcin /*MC 320abe9303eSLisandro Dalcin PETSCPARTITIONERPTSCOTCH = "ptscotch" - A PetscPartitioner object using the PT-Scotch library 321abe9303eSLisandro Dalcin 322abe9303eSLisandro Dalcin Level: intermediate 323abe9303eSLisandro Dalcin 324abe9303eSLisandro Dalcin Options Database Keys: 325abe9303eSLisandro Dalcin + -petscpartitioner_ptscotch_strategy <string> - PT-Scotch strategy. Choose one of default quality speed balance safety scalability recursive remap 326abe9303eSLisandro Dalcin - -petscpartitioner_ptscotch_imbalance <val> - Load imbalance ratio 327abe9303eSLisandro Dalcin 328abe9303eSLisandro Dalcin Notes: when the graph is on a single process, this partitioner actually uses Scotch and not PT-Scotch 329abe9303eSLisandro Dalcin 330abe9303eSLisandro Dalcin .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 331abe9303eSLisandro Dalcin M*/ 332abe9303eSLisandro Dalcin 333abe9303eSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_PTScotch(PetscPartitioner part) 334abe9303eSLisandro Dalcin { 335abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p; 336abe9303eSLisandro Dalcin PetscErrorCode ierr; 337abe9303eSLisandro Dalcin 338abe9303eSLisandro Dalcin PetscFunctionBegin; 339abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 340abe9303eSLisandro Dalcin ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 341abe9303eSLisandro Dalcin part->data = p; 342abe9303eSLisandro Dalcin 343abe9303eSLisandro Dalcin ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)part),&p->pcomm);CHKERRQ(ierr); 344abe9303eSLisandro Dalcin p->strategy = 0; 345abe9303eSLisandro Dalcin p->imbalance = 0.01; 346abe9303eSLisandro Dalcin 347abe9303eSLisandro Dalcin ierr = PetscPartitionerInitialize_PTScotch(part);CHKERRQ(ierr); 348abe9303eSLisandro Dalcin ierr = PetscCitationsRegister(PTScotchPartitionerCitation, &PTScotchPartitionerCite);CHKERRQ(ierr); 349abe9303eSLisandro Dalcin PetscFunctionReturn(0); 350abe9303eSLisandro Dalcin } 351