/*MC PETSCVIEWERDRAW - A viewer that generates graphics, either to the screen or a file .seealso: PetscViewerDrawOpen(), PetscViewerDrawGetDraw(), PETSC_VIEWER_DRAW_(),PETSC_VIEWER_DRAW_SELF, PETSC_VIEWER_DRAW_WORLD, PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() Level: beginner M*/ PETSC_EXTERN PetscErrorCode PetscViewerCreate_Draw(PetscViewer viewer) { PetscErrorCode ierr; PetscViewer_Draw *vdraw; PetscFunctionBegin; ierr = PetscNewLog(viewer,&vdraw);CHKERRQ(ierr); viewer->data = (void*)vdraw; viewer->ops->flush = PetscViewerFlush_Draw; viewer->ops->view = PetscViewerView_Draw; viewer->ops->destroy = PetscViewerDestroy_Draw; viewer->ops->setfromoptions = PetscViewerSetFromOptions_Draw; viewer->ops->getsubviewer = PetscViewerGetSubViewer_Draw; viewer->ops->restoresubviewer = PetscViewerRestoreSubViewer_Draw; /* these are created on the fly if requested */ vdraw->draw_max = 5; vdraw->draw_base = 0; vdraw->w = PETSC_DECIDE; vdraw->h = PETSC_DECIDE; ierr = PetscCalloc3(vdraw->draw_max,&vdraw->draw,vdraw->draw_max,&vdraw->drawlg,vdraw->draw_max,&vdraw->drawaxis);CHKERRQ(ierr); vdraw->singleton_made = PETSC_FALSE; PetscFunctionReturn(0); }
/*@C PetscViewerDrawGetDraw - Returns PetscDraw object from PetscViewer object. This PetscDraw object may then be used to perform graphics using PetscDrawXXX() commands. Collective on PetscViewer Input Parameters: + viewer - the PetscViewer (created with PetscViewerDrawOpen()) - windownumber - indicates which subwindow (usually 0) Ouput Parameter: . draw - the draw object Level: intermediate Concepts: drawing^accessing PetscDraw context from PetscViewer Concepts: graphics .seealso: PetscViewerDrawGetLG(), PetscViewerDrawGetAxis(), PetscViewerDrawOpen() @*/ PetscErrorCode PetscViewerDrawGetDraw(PetscViewer viewer,PetscInt windownumber,PetscDraw *draw) { PetscViewer_Draw *vdraw; PetscErrorCode ierr; PetscBool isdraw; PetscFunctionBegin; PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); PetscValidLogicalCollectiveInt(viewer,windownumber,2); if (draw) PetscValidPointer(draw,3); ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); if (!isdraw) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be draw type PetscViewer"); if (windownumber < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Window number cannot be negative"); vdraw = (PetscViewer_Draw*)viewer->data; windownumber += vdraw->draw_base; if (windownumber >= vdraw->draw_max) { /* allocate twice as many slots as needed */ PetscInt draw_max = vdraw->draw_max; PetscDraw *tdraw = vdraw->draw; PetscDrawLG *drawlg = vdraw->drawlg; PetscDrawAxis *drawaxis = vdraw->drawaxis; vdraw->draw_max = 2*windownumber; ierr = PetscCalloc3(vdraw->draw_max,&vdraw->draw,vdraw->draw_max,&vdraw->drawlg,vdraw->draw_max,&vdraw->drawaxis);CHKERRQ(ierr); ierr = PetscMemcpy(vdraw->draw,tdraw,draw_max*sizeof(PetscDraw));CHKERRQ(ierr); ierr = PetscMemcpy(vdraw->drawlg,drawlg,draw_max*sizeof(PetscDrawLG));CHKERRQ(ierr); ierr = PetscMemcpy(vdraw->drawaxis,drawaxis,draw_max*sizeof(PetscDrawAxis));CHKERRQ(ierr); ierr = PetscFree3(tdraw,drawlg,drawaxis);CHKERRQ(ierr); } if (!vdraw->draw[windownumber]) { char *title = vdraw->title, tmp_str[128]; if (windownumber) { ierr = PetscSNPrintf(tmp_str,sizeof(tmp_str),"%s:%d",vdraw->title?vdraw->title:"",windownumber);CHKERRQ(ierr); title = tmp_str; } ierr = PetscDrawCreate(PetscObjectComm((PetscObject)viewer),vdraw->display,title,PETSC_DECIDE,PETSC_DECIDE,vdraw->w,vdraw->h,&vdraw->draw[windownumber]);CHKERRQ(ierr); ierr = PetscLogObjectParent((PetscObject)viewer,(PetscObject)vdraw->draw[windownumber]);CHKERRQ(ierr); if (vdraw->drawtype) { ierr = PetscDrawSetType(vdraw->draw[windownumber],vdraw->drawtype);CHKERRQ(ierr); } ierr = PetscDrawSetPause(vdraw->draw[windownumber],vdraw->pause);CHKERRQ(ierr); ierr = PetscDrawSetOptionsPrefix(vdraw->draw[windownumber],((PetscObject)viewer)->prefix);CHKERRQ(ierr); ierr = PetscDrawSetFromOptions(vdraw->draw[windownumber]);CHKERRQ(ierr); } if (draw) *draw = vdraw->draw[windownumber]; if (draw) PetscValidHeaderSpecific(*draw,PETSC_DRAW_CLASSID,-1); PetscFunctionReturn(0); }
PetscErrorCode PetscFEGeomCreate(PetscQuadrature quad, PetscInt numCells, PetscInt dimEmbed, PetscBool faceData, PetscFEGeom **geom) { PetscFEGeom *g; PetscInt dim, Nq, N; const PetscReal *p; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscQuadratureGetData(quad,&dim,NULL,&Nq,&p,NULL);CHKERRQ(ierr); ierr = PetscNew(&g);CHKERRQ(ierr); g->xi = p; g->numCells = numCells; g->numPoints = Nq; g->dim = dim; g->dimEmbed = dimEmbed; N = numCells * Nq; ierr = PetscCalloc3(N * dimEmbed, &g->v, N * dimEmbed * dimEmbed, &g->J, N, &g->detJ);CHKERRQ(ierr); if (faceData) { ierr = PetscCalloc4(numCells, &g->face, N * dimEmbed, &g->n, N * dimEmbed * dimEmbed, &(g->suppInvJ[0]), N * dimEmbed * dimEmbed, &(g->suppInvJ[1]));CHKERRQ(ierr); } ierr = PetscCalloc1(N * dimEmbed * dimEmbed, &g->invJ);CHKERRQ(ierr); *geom = g; PetscFunctionReturn(0); }
/*@ DMPlexGetOrdering - Calculate a reordering of the mesh Collective on DM Input Parameter: + dm - The DMPlex object . otype - type of reordering, one of the following: $ MATORDERINGNATURAL - Natural $ MATORDERINGND - Nested Dissection $ MATORDERING1WD - One-way Dissection $ MATORDERINGRCM - Reverse Cuthill-McKee $ MATORDERINGQMD - Quotient Minimum Degree - label - [Optional] Label used to segregate ordering into sets, or NULL Output Parameter: . perm - The point permutation as an IS, perm[old point number] = new point number Note: The label is used to group sets of points together by label value. This makes it easy to reorder a mesh which has different types of cells, and then loop over each set of reordered cells for assembly. Level: intermediate .keywords: mesh .seealso: MatGetOrdering() @*/ PetscErrorCode DMPlexGetOrdering(DM dm, MatOrderingType otype, DMLabel label, IS *perm) { PetscInt numCells = 0; PetscInt *start = NULL, *adjacency = NULL, *cperm, *clperm, *invclperm, *mask, *xls, pStart, pEnd, c, i; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(dm, DM_CLASSID, 1); PetscValidPointer(perm, 3); ierr = DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency);CHKERRQ(ierr); ierr = PetscMalloc3(numCells,&cperm,numCells,&mask,numCells*2,&xls);CHKERRQ(ierr); if (numCells) { /* Shift for Fortran numbering */ for (i = 0; i < start[numCells]; ++i) ++adjacency[i]; for (i = 0; i <= numCells; ++i) ++start[i]; ierr = SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls);CHKERRQ(ierr); } ierr = PetscFree(start);CHKERRQ(ierr); ierr = PetscFree(adjacency);CHKERRQ(ierr); /* Shift for Fortran numbering */ for (c = 0; c < numCells; ++c) --cperm[c]; /* Segregate */ if (label) { IS valueIS; const PetscInt *values; PetscInt numValues, numPoints = 0; PetscInt *sperm, *vsize, *voff, v; ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr); ierr = ISSort(valueIS);CHKERRQ(ierr); ierr = ISGetLocalSize(valueIS, &numValues);CHKERRQ(ierr); ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); ierr = PetscCalloc3(numCells,&sperm,numValues,&vsize,numValues+1,&voff);CHKERRQ(ierr); for (v = 0; v < numValues; ++v) { ierr = DMLabelGetStratumSize(label, values[v], &vsize[v]);CHKERRQ(ierr); if (v < numValues-1) voff[v+2] += vsize[v] + voff[v+1]; numPoints += vsize[v]; } if (numPoints != numCells) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %D cells < %D total", numPoints, numCells); for (c = 0; c < numCells; ++c) { const PetscInt oldc = cperm[c]; PetscInt val, vloc; ierr = DMLabelGetValue(label, oldc, &val);CHKERRQ(ierr); if (val == -1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %D not present in label", oldc); ierr = PetscFindInt(val, numValues, values, &vloc);CHKERRQ(ierr); if (vloc < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %D not present label", val); sperm[voff[vloc+1]++] = oldc; } for (v = 0; v < numValues; ++v) { if (voff[v+1] - voff[v] != vsize[v]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of %D values found is %D != %D", values[v], voff[v+1] - voff[v], vsize[v]); } ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); ierr = ISDestroy(&valueIS);CHKERRQ(ierr); ierr = PetscMemcpy(cperm, sperm, numCells * sizeof(PetscInt));CHKERRQ(ierr); ierr = PetscFree3(sperm, vsize, voff);CHKERRQ(ierr); } /* Construct closure */ ierr = DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm);CHKERRQ(ierr); ierr = PetscFree3(cperm,mask,xls);CHKERRQ(ierr); ierr = PetscFree(clperm);CHKERRQ(ierr); /* Invert permutation */ ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, invclperm, PETSC_OWN_POINTER, perm);CHKERRQ(ierr); PetscFunctionReturn(0); }
PETSC_EXTERN PetscErrorCode MatGetOrdering_AWBM(Mat A, MatOrderingType type, IS *permR, IS *permC) { Vec *scalR, *scalC, scalRVec, scalCVec; scalR = &scalRVec; scalC = &scalCVec; /* EVERYTHING IS WRITTEN AS IF THE MATRIX WERE COLUMN-MAJOR */ Mat_SeqAIJ *aij = (Mat_SeqAIJ *) A->data; PetscInt n = A->rmap->n; /* Number of local columns */ PetscInt m = A->cmap->n; /* Number of local rows */ PetscInt *match; /* The row matched to each column, and inverse column permutation */ PetscInt *matchR; /* The column matched to each row */ PetscInt *p; /* The column permutation */ const PetscInt *ia = aij->i; const PetscInt *ja = aij->j; const MatScalar *a = aij->a; Vec colMax; PetscScalar *a_j, *sr, *sc; PetscReal *weights /* c_ij */, *u /* u_i */, *v /* v_j */, eps = PETSC_SQRT_MACHINE_EPSILON; PetscInt debug = 0, r, c, r1, c1; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscOptionsGetInt(NULL, "-debug", &debug, NULL);CHKERRQ(ierr); ierr = MatGetVecs(A, NULL, &colMax);CHKERRQ(ierr); ierr = MatGetRowMaxAbs(A, colMax, NULL);CHKERRQ(ierr); ierr = PetscMalloc2(n, &match, m, &matchR);CHKERRQ(ierr); ierr = PetscMalloc1(n, &p);CHKERRQ(ierr); ierr = PetscCalloc3(m, &u, n, &v, ia[n], &weights);CHKERRQ(ierr); for (c = 0; c < n; ++c) match[c] = -1; /* Compute weights */ ierr = VecGetArray(colMax, &a_j);CHKERRQ(ierr); for (c = 0; c < n; ++c) { for (r = ia[c]; r < ia[c+1]; ++r) { PetscReal ar = PetscAbsScalar(a[r]); if (ar == 0.0) weights[r] = PETSC_MAX_REAL; else weights[r] = log(a_j[c]/ar); } } /* Compute local row weights */ for (r = 0; r < m; ++r) u[r] = PETSC_MAX_REAL; for (c = 0; c < n; ++c) { for (r = ia[c]; r < ia[c+1]; ++r) { u[ja[r]] = PetscMin(u[ja[r]], weights[r]); } } /* Compute local column weights */ for (c = 0; c < n; ++c) { v[c] = PETSC_MAX_REAL; for (r = ia[c]; r < ia[c+1]; ++r) { v[c] = PetscMin(v[c], weights[r] - u[ja[r]]); } } for (r = 0; r < m; ++r) matchR[r] = -1; /* Match columns */ ierr = CheckUnmatched(n, match, matchR);CHKERRQ(ierr); for (c = 0; c < n; ++c) { /* if (match[c] >= 0) continue; */ if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Row %d\n Weights:", c);CHKERRQ(ierr);} for (r = ia[c]; r < ia[c+1]; ++r) { PetscReal weight = weights[r] - u[ja[r]] - v[c]; if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " %g", weight);CHKERRQ(ierr);} if ((weight <= eps) && (matchR[ja[r]] < 0)) { if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Matched %d -- %d\n", c, ja[r]);CHKERRQ(ierr);} match[c] = ja[r]; matchR[ja[r]] = c; break; } } if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);} } /* Deal with unmatched columns */ ierr = CheckUnmatched(n, match, matchR);CHKERRQ(ierr); for (c = 0; c < n; ++c) { if (match[c] >= 0) continue; for (r = ia[c]; r < ia[c+1]; ++r) { PetscReal weight = weights[r] - u[ja[r]] - v[c]; if (weight > eps) continue; /* \bar c_ij = 0 and (r, j1) \in M */ c1 = matchR[ja[r]]; for (r1 = ia[c1]; r1 < ia[c1+1]; ++r1) { PetscReal weight1 = weights[r1] - u[ja[r1]] - v[c1]; if ((matchR[ja[r1]] < 0) && (weight1 <= eps)) { /* (r, c1) in M is replaced by (r, c) and (r1, c1) */ if (debug) { ierr = PetscPrintf(PETSC_COMM_SELF, "Replaced match %d -- %d\n", c1, ja[r]);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_SELF, " Added match %d -- %d\n", c, ja[r]);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_SELF, " Added match %d -- %d\n", c1, ja[r1]);CHKERRQ(ierr); } match[c] = ja[r]; matchR[ja[r]] = c; match[c1] = ja[r1]; matchR[ja[r1]] = c1; break; } } if (match[c] >= 0) break; } } /* Allow matching with non-optimal rows */ ierr = CheckUnmatched(n, match, matchR);CHKERRQ(ierr); for (c = 0; c < n; ++c) { if (match[c] >= 0) continue; for (r = ia[c]; r < ia[c+1]; ++r) { if (matchR[ja[r]] < 0) { if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Matched non-opt %d -- %d\n", c, ja[r]);CHKERRQ(ierr);} match[c] = ja[r]; matchR[ja[r]] = c; break; } } } /* Deal with non-optimal unmatched columns */ ierr = CheckUnmatched(n, match, matchR);CHKERRQ(ierr); for (c = 0; c < n; ++c) { if (match[c] >= 0) continue; for (r = ia[c]; r < ia[c+1]; ++r) { /* \bar c_ij = 0 and (r, j1) \in M */ c1 = matchR[ja[r]]; for (r1 = ia[c1]; r1 < ia[c1+1]; ++r1) { if (matchR[ja[r1]] < 0) { /* (r, c1) in M is replaced by (r, c) and (r1, c1) */ if (debug) { ierr = PetscPrintf(PETSC_COMM_SELF, "Replaced match %d -- %d\n", c1, ja[r]);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_SELF, " Added match %d -- %d\n", c, ja[r]);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_SELF, " Added match %d -- %d\n", c1, ja[r1]);CHKERRQ(ierr); } match[c] = ja[r]; matchR[ja[r]] = c; match[c1] = ja[r1]; matchR[ja[r1]] = c1; break; } } if (match[c] >= 0) break; } } /* Complete matching */ ierr = CheckUnmatched(n, match, matchR);CHKERRQ(ierr); for (c = 0, r = 0; c < n; ++c) { if (match[c] >= n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Column %d matched to invalid row %d", c, match[c]); if (match[c] < 0) { for (; r < n; ++r) { if (matchR[r] < 0) { if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Matched default %d -- %d\n", c, r);CHKERRQ(ierr);} match[c] = r; matchR[r] = c; break; } } } } /* Check matching */ ierr = CheckUnmatched(n, match, matchR);CHKERRQ(ierr); for (c = 0; c < n; ++c) { if (match[c] < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Column %d unmatched", c); if (match[c] >= n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Column %d matched to invalid row %d", c, match[c]); } /* Make permutation */ for (c = 0; c < n; ++c) {p[match[c]] = c;} ierr = ISCreateGeneral(PETSC_COMM_SELF, n, p, PETSC_OWN_POINTER, permR);CHKERRQ(ierr); ierr = ISSetPermutation(*permR);CHKERRQ(ierr); ierr = ISCreateStride(PETSC_COMM_SELF, n, 0, 1, permC);CHKERRQ(ierr); ierr = ISSetPermutation(*permC);CHKERRQ(ierr); ierr = PetscFree2(match, matchR);CHKERRQ(ierr); /* Make scaling */ ierr = VecCreateSeq(PETSC_COMM_SELF, n, scalR);CHKERRQ(ierr); ierr = VecCreateSeq(PETSC_COMM_SELF, n, scalC);CHKERRQ(ierr); ierr = VecGetArray(*scalR, &sr);CHKERRQ(ierr); ierr = VecGetArray(*scalC, &sc);CHKERRQ(ierr); for (c = 0; c < n; ++c) { sr[c] = PetscExpReal(v[c])/a_j[c]; sc[c] = PetscExpReal(u[c]); } ierr = VecRestoreArray(*scalR, &sr);CHKERRQ(ierr); ierr = VecRestoreArray(*scalC, &sc);CHKERRQ(ierr); ierr = VecRestoreArray(colMax, &a_j);CHKERRQ(ierr); ierr = VecDestroy(&colMax);CHKERRQ(ierr); ierr = PetscFree3(u,v,weights);CHKERRQ(ierr); ierr = VecDestroy(scalR);CHKERRQ(ierr); ierr = VecDestroy(scalC);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@ DMPlexOrient - Give a consistent orientation to the input mesh Input Parameters: . dm - The DM Note: The orientation data for the DM are change in-place. $ This routine will fail for non-orientable surfaces, such as the Moebius strip. Level: advanced .seealso: DMCreate(), DMPLEX @*/ PetscErrorCode DMPlexOrient(DM dm) { MPI_Comm comm; PetscSF sf; const PetscInt *lpoints; const PetscSFNode *rpoints; PetscSFNode *rorntComp = NULL, *lorntComp = NULL; PetscInt *numNeighbors, **neighbors; PetscSFNode *nrankComp; PetscBool *match, *flipped; PetscBT seenCells, flippedCells, seenFaces; PetscInt *faceFIFO, fTop, fBottom, *cellComp, *faceComp; PetscInt numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0; PetscMPIInt rank, size, numComponents, comp = 0; PetscBool flg, flg2; PetscViewer viewer = NULL, selfviewer = NULL; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr); ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2);CHKERRQ(ierr); ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr); /* Truth Table mismatch flips do action mismatch flipA ^ flipB action F 0 flips no F F F F 1 flip yes F T T F 2 flips no T F T T 0 flips yes T T F T 1 flip no T 2 flips yes */ ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr); ierr = DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);CHKERRQ(ierr); ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr); ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr); ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr); ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr); ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr); ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr); ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr); ierr = PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp);CHKERRQ(ierr); /* OLD STYLE - Add an integer array over cells and faces (component) for connected component number Foreach component - Mark the initial cell as seen - Process component as usual - Set component for all seenCells - Wipe seenCells and seenFaces (flippedCells can stay) - Generate parallel adjacency for component using SF and seenFaces - Collect numComponents adj data from each proc to 0 - Build same serial graph - Use same solver - Use Scatterv to to send back flipped flags for each component - Negate flippedCells by component NEW STYLE - Create the adj on each process - Bootstrap to complete graph on proc 0 */ /* Loop over components */ for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1; do { /* Look for first unmarked cell */ for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break; if (cell >= cEnd) break; /* Initialize FIFO with first cell in component */ { const PetscInt *cone; PetscInt coneSize; fTop = fBottom = 0; ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); for (c = 0; c < coneSize; ++c) { faceFIFO[fBottom++] = cone[c]; ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr); } ierr = PetscBTSet(seenCells, cell-cStart);CHKERRQ(ierr); } /* Consider each face in FIFO */ while (fTop < fBottom) { ierr = DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr); } /* Set component for cells and faces */ for (cell = 0; cell < cEnd-cStart; ++cell) { if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp; } for (face = 0; face < fEnd-fStart; ++face) { if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp; } /* Wipe seenCells and seenFaces for next component */ ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr); ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr); ++comp; } while (1); numComponents = comp; if (flg) { PetscViewer v; ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr); ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr); ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr); ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); ierr = PetscViewerFlush(v);CHKERRQ(ierr); ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr); } /* Now all subdomains are oriented, but we need a consistent parallel orientation */ if (numLeaves >= 0) { /* Store orientations of boundary faces*/ ierr = PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);CHKERRQ(ierr); for (face = fStart; face < fEnd; ++face) { const PetscInt *cone, *support, *ornt; PetscInt coneSize, supportSize; ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); if (supportSize != 1) continue; ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr); ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr); ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr); for (c = 0; c < coneSize; ++c) if (cone[c] == face) break; if (dim == 1) { /* Use cone position instead, shifted to -1 or 1 */ if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2; else rorntComp[face].rank = c*2-1; } else { if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1; else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1; } rorntComp[face].index = faceComp[face-fStart]; } /* Communicate boundary edge orientations */ ierr = PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr); ierr = PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr); } /* Get process adjacency */ ierr = PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);CHKERRQ(ierr); viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm)); if (flg2) {ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);} ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr); for (comp = 0; comp < numComponents; ++comp) { PetscInt l, n; numNeighbors[comp] = 0; ierr = PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);CHKERRQ(ierr); /* I know this is p^2 time in general, but for bounded degree its alright */ for (l = 0; l < numLeaves; ++l) { const PetscInt face = lpoints[l]; /* Find a representative face (edge) separating pairs of procs */ if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) { const PetscInt rrank = rpoints[l].rank; const PetscInt rcomp = lorntComp[face].index; for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break; if (n >= numNeighbors[comp]) { PetscInt supportSize; ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize); if (flg) {ierr = PetscViewerASCIIPrintf(selfviewer, "[%d]: component %d, Found representative leaf %d (face %d) connecting to face %d on (%d, %d) with orientation %d\n", rank, comp, l, face, rpoints[l].index, rrank, rcomp, lorntComp[face].rank);CHKERRQ(ierr);} neighbors[comp][numNeighbors[comp]++] = l; } } } totNeighbors += numNeighbors[comp]; } ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr); ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); if (flg2) {ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);} ierr = PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);CHKERRQ(ierr); for (comp = 0, off = 0; comp < numComponents; ++comp) { PetscInt n; for (n = 0; n < numNeighbors[comp]; ++n, ++off) { const PetscInt face = lpoints[neighbors[comp][n]]; const PetscInt o = rorntComp[face].rank*lorntComp[face].rank; if (o < 0) match[off] = PETSC_TRUE; else if (o > 0) match[off] = PETSC_FALSE; else SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %d (%d, %d) neighbor: %d comp: %d", face, rorntComp[face], lorntComp[face], neighbors[comp][n], comp); nrankComp[off].rank = rpoints[neighbors[comp][n]].rank; nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index; } ierr = PetscFree(neighbors[comp]);CHKERRQ(ierr); } /* Collect the graph on 0 */ if (numLeaves >= 0) { Mat G; PetscBT seenProcs, flippedProcs; PetscInt *procFIFO, pTop, pBottom; PetscInt *N = NULL, *Noff; PetscSFNode *adj = NULL; PetscBool *val = NULL; PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o; PetscMPIInt size = 0; ierr = PetscCalloc1(numComponents, &flipped);CHKERRQ(ierr); if (!rank) {ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);} ierr = PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff);CHKERRQ(ierr); ierr = MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);CHKERRQ(ierr); for (p = 0; p < size; ++p) { displs[p+1] = displs[p] + Nc[p]; } if (!rank) {ierr = PetscMalloc1(displs[size],&N);CHKERRQ(ierr);} ierr = MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);CHKERRQ(ierr); for (p = 0, o = 0; p < size; ++p) { recvcounts[p] = 0; for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o]; displs[p+1] = displs[p] + recvcounts[p]; } if (!rank) {ierr = PetscMalloc2(displs[size], &adj, displs[size], &val);CHKERRQ(ierr);} ierr = MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);CHKERRQ(ierr); ierr = MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr); ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr); if (!rank) { for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];} if (flg) { PetscInt n; for (p = 0, off = 0; p < size; ++p) { for (c = 0; c < Nc[p]; ++c) { ierr = PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);CHKERRQ(ierr); for (n = 0; n < N[Noff[p]+c]; ++n, ++off) { ierr = PetscPrintf(PETSC_COMM_SELF, " edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);CHKERRQ(ierr); } } } } /* Symmetrize the graph */ ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr); ierr = MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);CHKERRQ(ierr); ierr = MatSetUp(G);CHKERRQ(ierr); for (p = 0, off = 0; p < size; ++p) { for (c = 0; c < Nc[p]; ++c) { const PetscInt r = Noff[p]+c; PetscInt n; for (n = 0; n < N[r]; ++n, ++off) { const PetscInt q = Noff[adj[off].rank] + adj[off].index; const PetscScalar o = val[off] ? 1.0 : 0.0; ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr); ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr); } } } ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = PetscBTCreate(Noff[size], &seenProcs);CHKERRQ(ierr); ierr = PetscBTMemzero(Noff[size], seenProcs);CHKERRQ(ierr); ierr = PetscBTCreate(Noff[size], &flippedProcs);CHKERRQ(ierr); ierr = PetscBTMemzero(Noff[size], flippedProcs);CHKERRQ(ierr); ierr = PetscMalloc1(Noff[size], &procFIFO);CHKERRQ(ierr); pTop = pBottom = 0; for (p = 0; p < Noff[size]; ++p) { if (PetscBTLookup(seenProcs, p)) continue; /* Initialize FIFO with next proc */ procFIFO[pBottom++] = p; ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr); /* Consider each proc in FIFO */ while (pTop < pBottom) { const PetscScalar *ornt; const PetscInt *neighbors; PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n; proc = procFIFO[pTop++]; flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr); /* Loop over neighboring procs */ for (n = 0; n < numNeighbors; ++n) { nproc = neighbors[n]; mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; seen = PetscBTLookup(seenProcs, nproc); flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; if (mismatch ^ (flippedA ^ flippedB)) { if (seen) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %d and %d do not match: Fault mesh is non-orientable", proc, nproc); if (!flippedB) { ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr); } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); if (!seen) { procFIFO[pBottom++] = nproc; ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr); } } } } ierr = PetscFree(procFIFO);CHKERRQ(ierr); ierr = MatDestroy(&G);CHKERRQ(ierr); ierr = PetscFree2(adj, val);CHKERRQ(ierr); ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr); } /* Scatter flip flags */ { PetscBool *flips = NULL; if (!rank) { ierr = PetscMalloc1(Noff[size], &flips);CHKERRQ(ierr); for (p = 0; p < Noff[size]; ++p) { flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; if (flg && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);CHKERRQ(ierr);} } for (p = 0; p < size; ++p) { displs[p+1] = displs[p] + Nc[p]; } } ierr = MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);CHKERRQ(ierr); ierr = PetscFree(flips);CHKERRQ(ierr); } if (!rank) {ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);} ierr = PetscFree(N);CHKERRQ(ierr); ierr = PetscFree4(recvcounts, displs, Nc, Noff);CHKERRQ(ierr); ierr = PetscFree2(nrankComp, match);CHKERRQ(ierr); /* Decide whether to flip cells in each component */ for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {ierr = PetscBTNegate(flippedCells, c);CHKERRQ(ierr);}} ierr = PetscFree(flipped);CHKERRQ(ierr); } if (flg) { PetscViewer v; ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr); ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr); ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr); ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); ierr = PetscViewerFlush(v);CHKERRQ(ierr); ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr); } /* Reverse flipped cells in the mesh */ for (c = cStart; c < cEnd; ++c) { if (PetscBTLookup(flippedCells, c-cStart)) { ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr); } } ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr); ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr); ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr); ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr); ierr = PetscFree2(rorntComp, lorntComp);CHKERRQ(ierr); ierr = PetscFree3(faceFIFO, cellComp, faceComp);CHKERRQ(ierr); PetscFunctionReturn(0); }