12e68589bSMark F. Adams /* 22e68589bSMark F. Adams GAMG geometric-algebric multiogrid PC - Mark Adams 2011 32e68589bSMark F. Adams */ 42e68589bSMark F. Adams 52e68589bSMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 62e68589bSMark F. Adams 72e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE) 85a449d36SBarry Smith #if !defined(ANSI_DECLARATORS) 95a449d36SBarry Smith #define ANSI_DECLARATORS 105a449d36SBarry Smith #endif 112e68589bSMark F. Adams #include <triangle.h> 122e68589bSMark F. Adams #endif 132e68589bSMark F. Adams 142e68589bSMark F. Adams #include <petscblaslapack.h> 152e68589bSMark F. Adams 16c8b0795cSMark F. Adams /* Private context for the GAMG preconditioner */ 17c8b0795cSMark F. Adams typedef struct { 18c8b0795cSMark F. Adams PetscInt lid; /* local vertex index */ 19c8b0795cSMark F. Adams PetscInt degree; /* vertex degree */ 20c8b0795cSMark F. Adams } GAMGNode; 212fa5cd67SKarl Rupp 22d71ae5a4SJacob Faibussowitsch static inline int petsc_geo_mg_compare(const void *a, const void *b) 23d71ae5a4SJacob Faibussowitsch { 24c8b0795cSMark F. Adams return (((GAMGNode *)a)->degree - ((GAMGNode *)b)->degree); 25c8b0795cSMark F. Adams } 26c8b0795cSMark F. Adams 2704c3f3b8SBarry Smith // PetscClangLinter pragma disable: -fdoc-sowing-chars 282e68589bSMark F. Adams /* 292e68589bSMark F. Adams PCSetCoordinates_GEO 302e68589bSMark F. Adams 312e68589bSMark F. Adams Input Parameter: 322e68589bSMark F. Adams . pc - the preconditioner context 332e68589bSMark F. Adams */ 3466976f2fSJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 35d71ae5a4SJacob Faibussowitsch { 362e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 372e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 3890fbc344SStefano Zampini PetscInt arrsz, bs, my0, kk, ii, nloc, Iend, aloc; 392e68589bSMark F. Adams Mat Amat = pc->pmat; 402e68589bSMark F. Adams 412e68589bSMark F. Adams PetscFunctionBegin; 422e68589bSMark F. Adams PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1); 439566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat, &bs)); 449566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Amat, &my0, &Iend)); 4590fbc344SStefano Zampini aloc = (Iend - my0); 462e68589bSMark F. Adams nloc = (Iend - my0) / bs; 47a2f3521dSMark F. Adams 482472a847SBarry Smith PetscCheck(nloc == a_nloc || aloc == a_nloc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of local blocks %" PetscInt_FMT " must be %" PetscInt_FMT " or %" PetscInt_FMT ".", a_nloc, nloc, aloc); 492e68589bSMark F. Adams 50c8b0795cSMark F. Adams pc_gamg->data_cell_rows = 1; 517827d75bSBarry Smith PetscCheck(coords || (nloc <= 0), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'."); 52c8b0795cSMark F. Adams pc_gamg->data_cell_cols = ndm; /* coordinates */ 532e68589bSMark F. Adams 54c8b0795cSMark F. Adams arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols; 552e68589bSMark F. Adams 562e68589bSMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 577f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 589566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data)); 602e68589bSMark F. Adams } 612e68589bSMark F. Adams for (kk = 0; kk < arrsz; kk++) pc_gamg->data[kk] = -999.; 622e68589bSMark F. Adams pc_gamg->data[arrsz] = -99.; 632e68589bSMark F. Adams /* copy data in - column oriented */ 6490fbc344SStefano Zampini if (nloc == a_nloc) { 652e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) { 66ad540459SPierre Jolivet for (ii = 0; ii < ndm; ii++) pc_gamg->data[ii * nloc + kk] = coords[kk * ndm + ii]; 672e68589bSMark F. Adams } 6890fbc344SStefano Zampini } else { /* assumes the coordinates are blocked */ 6990fbc344SStefano Zampini for (kk = 0; kk < nloc; kk++) { 70ad540459SPierre Jolivet for (ii = 0; ii < ndm; ii++) pc_gamg->data[ii * nloc + kk] = coords[bs * kk * ndm + ii]; 7190fbc344SStefano Zampini } 7290fbc344SStefano Zampini } 7363a3b9bcSJacob Faibussowitsch PetscCheck(pc_gamg->data[arrsz] == -99., PETSC_COMM_SELF, PETSC_ERR_PLIB, "pc_gamg->data[arrsz %" PetscInt_FMT "] %g != -99.", arrsz, (double)pc_gamg->data[arrsz]); 742e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 762e68589bSMark F. Adams } 772e68589bSMark F. Adams 7804c3f3b8SBarry Smith // PetscClangLinter pragma disable: -fdoc-sowing-chars 792e68589bSMark F. Adams /* 802e68589bSMark F. Adams PCSetData_GEO 812e68589bSMark F. Adams 822e68589bSMark F. Adams Input Parameter: 832e68589bSMark F. Adams . pc - 842e68589bSMark F. Adams */ 8566976f2fSJacob Faibussowitsch static PetscErrorCode PCSetData_GEO(PC pc, Mat m) 86d71ae5a4SJacob Faibussowitsch { 872e68589bSMark F. Adams PetscFunctionBegin; 88ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "GEO MG needs coordinates"); 892e68589bSMark F. Adams } 902e68589bSMark F. Adams 9166976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GEO(PC pc, PetscOptionItems *PetscOptionsObject) 92d71ae5a4SJacob Faibussowitsch { 932e68589bSMark F. Adams PetscFunctionBegin; 94d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-GEO options"); 952e68589bSMark F. Adams { 962e68589bSMark F. Adams /* -pc_gamg_sa_nsmooths */ 972e68589bSMark F. Adams /* pc_gamg_sa->smooths = 0; */ 982e68589bSMark F. Adams /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */ 992e68589bSMark F. Adams /* "smoothing steps for smoothed aggregation, usually 1 (0)", */ 1002e68589bSMark F. Adams /* "PCGAMGSetNSmooths_AGG", */ 1012e68589bSMark F. Adams /* pc_gamg_sa->smooths, */ 1022e68589bSMark F. Adams /* &pc_gamg_sa->smooths, */ 1032e68589bSMark F. Adams /* &flag); */ 1042e68589bSMark F. Adams } 105d0609cedSBarry Smith PetscOptionsHeadEnd(); 1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1072e68589bSMark F. Adams } 1082e68589bSMark F. Adams 10904c3f3b8SBarry Smith // PetscClangLinter pragma disable: -fdoc-sowing-chars 1102e68589bSMark F. Adams /* 1112e68589bSMark F. Adams triangulateAndFormProl 1122e68589bSMark F. Adams 1132e68589bSMark F. Adams Input Parameter: 1142e68589bSMark F. Adams . selected_2 - list of selected local ID, includes selected ghosts 115a2f3521dSMark F. Adams . data_stride - 116a2f3521dSMark F. Adams . coords[2*data_stride] - column vector of local coordinates w/ ghosts 1172adfe9d3SBarry Smith . nselected_1 - selected IDs that go with base (1) graph includes selected ghosts 1180cbbd2e1SMark F. Adams . clid_lid_1[nselected_1] - lids of selected (c) nodes ??????????? 1192adfe9d3SBarry Smith . agg_lists_1 - list of aggregates selected_1 vertices of aggregate unselected vertices 1202e68589bSMark F. Adams . crsGID[selected.size()] - global index for prolongation operator 1212e68589bSMark F. Adams . bs - block size 1222e68589bSMark F. Adams Output Parameter: 1232e68589bSMark F. Adams . a_Prol - prolongation operator 1242e68589bSMark F. Adams . a_worst_best - measure of worst missed fine vertex, 0 is no misses 1252e68589bSMark F. Adams */ 126d71ae5a4SJacob Faibussowitsch static PetscErrorCode triangulateAndFormProl(IS selected_2, PetscInt data_stride, PetscReal coords[], PetscInt nselected_1, const PetscInt clid_lid_1[], const PetscCoarsenData *agg_lists_1, const PetscInt crsGID[], PetscInt bs, Mat a_Prol, PetscReal *a_worst_best) 127d71ae5a4SJacob Faibussowitsch { 1282e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE) 129b43b03e9SMark F. Adams PetscInt jj, tid, tt, idx, nselected_2; 1302e68589bSMark F. Adams struct triangulateio in, mid; 1310cbbd2e1SMark F. Adams const PetscInt *selected_idx_2; 13273911c69SBarry Smith PetscMPIInt rank; 1332e68589bSMark F. Adams PetscInt Istart, Iend, nFineLoc, myFine0; 1342e68589bSMark F. Adams int kk, nPlotPts, sid; 1353b4367a7SBarry Smith MPI_Comm comm; 136c8b0795cSMark F. Adams PetscReal tm; 137c8b0795cSMark F. Adams 1386e111a19SKarl Rupp PetscFunctionBegin; 1399566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm)); 1409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1419566063dSJacob Faibussowitsch PetscCall(ISGetSize(selected_2, &nselected_2)); 1422e68589bSMark F. Adams if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */ 1432e68589bSMark F. Adams *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */ 144806fa848SBarry Smith } else *a_worst_best = 0.0; 1451c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm)); 146c8b0795cSMark F. Adams if (tm > 0.0) { 147c8b0795cSMark F. Adams *a_worst_best = 100.0; 1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1492e68589bSMark F. Adams } 1509566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend)); 1519371c9d4SSatish Balay nFineLoc = (Iend - Istart) / bs; 1529371c9d4SSatish Balay myFine0 = Istart / bs; 1532e68589bSMark F. Adams nPlotPts = nFineLoc; /* locals */ 1544c500f23SPierre Jolivet /* triangle */ 1552e68589bSMark F. Adams /* Define input points - in*/ 1562e68589bSMark F. Adams in.numberofpoints = nselected_2; 1572e68589bSMark F. Adams in.numberofpointattributes = 0; 1582e68589bSMark F. Adams /* get nselected points */ 1599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * nselected_2, &in.pointlist)); 1609566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected_2, &selected_idx_2)); 1612e68589bSMark F. Adams 1622e68589bSMark F. Adams for (kk = 0, sid = 0; kk < nselected_2; kk++, sid += 2) { 1632e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 1642e68589bSMark F. Adams in.pointlist[sid] = coords[lid]; 165a2f3521dSMark F. Adams in.pointlist[sid + 1] = coords[data_stride + lid]; 1662e68589bSMark F. Adams if (lid >= nFineLoc) nPlotPts++; 1672e68589bSMark F. Adams } 168c05e6a5bSJed Brown PetscCheck(sid == 2 * nselected_2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "sid %d != 2*nselected_2 %" PetscInt_FMT, sid, nselected_2); 1692e68589bSMark F. Adams 1702e68589bSMark F. Adams in.numberofsegments = 0; 1712e68589bSMark F. Adams in.numberofedges = 0; 1722e68589bSMark F. Adams in.numberofholes = 0; 1732e68589bSMark F. Adams in.numberofregions = 0; 1740a545947SLisandro Dalcin in.trianglelist = NULL; 1750a545947SLisandro Dalcin in.segmentmarkerlist = NULL; 1760a545947SLisandro Dalcin in.pointattributelist = NULL; 1770a545947SLisandro Dalcin in.pointmarkerlist = NULL; 1780a545947SLisandro Dalcin in.triangleattributelist = NULL; 1790a545947SLisandro Dalcin in.trianglearealist = NULL; 1800a545947SLisandro Dalcin in.segmentlist = NULL; 1810a545947SLisandro Dalcin in.holelist = NULL; 1820a545947SLisandro Dalcin in.regionlist = NULL; 1830a545947SLisandro Dalcin in.edgelist = NULL; 1840a545947SLisandro Dalcin in.edgemarkerlist = NULL; 1850a545947SLisandro Dalcin in.normlist = NULL; 1862fa5cd67SKarl Rupp 1872e68589bSMark F. Adams /* triangulate */ 1880a545947SLisandro Dalcin mid.pointlist = NULL; /* Not needed if -N switch used. */ 1892e68589bSMark F. Adams /* Not needed if -N switch used or number of point attributes is zero: */ 1900a545947SLisandro Dalcin mid.pointattributelist = NULL; 1910a545947SLisandro Dalcin mid.pointmarkerlist = NULL; /* Not needed if -N or -B switch used. */ 1920a545947SLisandro Dalcin mid.trianglelist = NULL; /* Not needed if -E switch used. */ 1932e68589bSMark F. Adams /* Not needed if -E switch used or number of triangle attributes is zero: */ 1940a545947SLisandro Dalcin mid.triangleattributelist = NULL; 1950a545947SLisandro Dalcin mid.neighborlist = NULL; /* Needed only if -n switch used. */ 1962e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P not used: */ 1970a545947SLisandro Dalcin mid.segmentlist = NULL; 1982e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P and -B not used: */ 1990a545947SLisandro Dalcin mid.segmentmarkerlist = NULL; 2000a545947SLisandro Dalcin mid.edgelist = NULL; /* Needed only if -e switch used. */ 2010a545947SLisandro Dalcin mid.edgemarkerlist = NULL; /* Needed if -e used and -B not used. */ 2022e68589bSMark F. Adams mid.numberoftriangles = 0; 2032e68589bSMark F. Adams 2042e68589bSMark F. Adams /* Triangulate the points. Switches are chosen to read and write a */ 2052e68589bSMark F. Adams /* PSLG (p), preserve the convex hull (c), number everything from */ 2062e68589bSMark F. Adams /* zero (z), assign a regional attribute to each element (A), and */ 2072e68589bSMark F. Adams /* produce an edge list (e), a Voronoi diagram (v), and a triangle */ 2082e68589bSMark F. Adams /* neighbor list (n). */ 2092e68589bSMark F. Adams if (nselected_2 != 0) { /* inactive processor */ 2102e68589bSMark F. Adams char args[] = "npczQ"; /* c is needed ? */ 2112e68589bSMark F. Adams triangulate(args, &in, &mid, (struct triangulateio *)NULL); 2122e68589bSMark F. Adams /* output .poly files for 'showme' */ 2132e68589bSMark F. Adams if (!PETSC_TRUE) { 2142e68589bSMark F. Adams static int level = 1; 2159371c9d4SSatish Balay FILE *file; 2169371c9d4SSatish Balay char fname[32]; 2172e68589bSMark F. Adams 218a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(fname, PETSC_STATIC_ARRAY_LENGTH(fname), "C%d_%d.poly", level, rank)); 2199371c9d4SSatish Balay file = fopen(fname, "w"); 2202e68589bSMark F. Adams /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/ 2212e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n", in.numberofpoints, 2, 0, 0); 2222e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 223ad540459SPierre Jolivet for (kk = 0, sid = 0; kk < in.numberofpoints; kk++, sid += 2) fprintf(file, "%d %e %e\n", kk, in.pointlist[sid], in.pointlist[sid + 1]); 2242e68589bSMark F. Adams /*One line: <# of segments> <# of boundary markers (0 or 1)> */ 2252e68589bSMark F. Adams fprintf(file, "%d %d\n", 0, 0); 2262e68589bSMark F. Adams /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */ 2272e68589bSMark F. Adams /* One line: <# of holes> */ 2282e68589bSMark F. Adams fprintf(file, "%d\n", 0); 2292e68589bSMark F. Adams /* Following lines: <hole #> <x> <y> */ 2302e68589bSMark F. Adams /* Optional line: <# of regional attributes and/or area constraints> */ 2312e68589bSMark F. Adams /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */ 2322e68589bSMark F. Adams fclose(file); 2332e68589bSMark F. Adams 2342e68589bSMark F. Adams /* elems */ 235a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(fname, PETSC_STATIC_ARRAY_LENGTH(fname), "C%d_%d.ele", level, rank)); 2369371c9d4SSatish Balay file = fopen(fname, "w"); 2372e68589bSMark F. Adams /* First line: <# of triangles> <nodes per triangle> <# of attributes> */ 2382e68589bSMark F. Adams fprintf(file, "%d %d %d\n", mid.numberoftriangles, 3, 0); 2392e68589bSMark F. Adams /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */ 240ad540459SPierre Jolivet for (kk = 0, sid = 0; kk < mid.numberoftriangles; kk++, sid += 3) fprintf(file, "%d %d %d %d\n", kk, mid.trianglelist[sid], mid.trianglelist[sid + 1], mid.trianglelist[sid + 2]); 2412e68589bSMark F. Adams fclose(file); 2422e68589bSMark F. Adams 243a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(fname, PETSC_STATIC_ARRAY_LENGTH(fname), "C%d_%d.node", level, rank)); 2449371c9d4SSatish Balay file = fopen(fname, "w"); 2452e68589bSMark F. Adams /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */ 2462e68589bSMark F. Adams /* fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); */ 2472e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n", nPlotPts, 2, 0, 0); 2482e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 249ad540459SPierre Jolivet for (kk = 0, sid = 0; kk < in.numberofpoints; kk++, sid += 2) fprintf(file, "%d %e %e\n", kk, in.pointlist[sid], in.pointlist[sid + 1]); 2502e68589bSMark F. Adams 2512e68589bSMark F. Adams sid /= 2; 2522e68589bSMark F. Adams for (jj = 0; jj < nFineLoc; jj++) { 2532e68589bSMark F. Adams PetscBool sel = PETSC_TRUE; 2542e68589bSMark F. Adams for (kk = 0; kk < nselected_2 && sel; kk++) { 2552e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 2562e68589bSMark F. Adams if (lid == jj) sel = PETSC_FALSE; 2572e68589bSMark F. Adams } 2582fa5cd67SKarl Rupp if (sel) fprintf(file, "%d %e %e\n", sid++, coords[jj], coords[data_stride + jj]); 2592e68589bSMark F. Adams } 2602e68589bSMark F. Adams fclose(file); 261c05e6a5bSJed Brown PetscCheck(sid == nPlotPts, PETSC_COMM_SELF, PETSC_ERR_PLIB, "sid %d != nPlotPts %d", sid, nPlotPts); 2622e68589bSMark F. Adams level++; 2632e68589bSMark F. Adams } 2642e68589bSMark F. Adams } 2652e68589bSMark F. Adams { /* form P - setup some maps */ 2660cbbd2e1SMark F. Adams PetscInt clid, mm, *nTri, *node_tri; 2672e68589bSMark F. Adams 2689566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nselected_2, &node_tri, nselected_2, &nTri)); 2692e68589bSMark F. Adams 2702e68589bSMark F. Adams /* need list of triangles on node */ 2712e68589bSMark F. Adams for (kk = 0; kk < nselected_2; kk++) nTri[kk] = 0; 2722e68589bSMark F. Adams for (tid = 0, kk = 0; tid < mid.numberoftriangles; tid++) { 2732e68589bSMark F. Adams for (jj = 0; jj < 3; jj++) { 2742e68589bSMark F. Adams PetscInt cid = mid.trianglelist[kk++]; 2752e68589bSMark F. Adams if (nTri[cid] == 0) node_tri[cid] = tid; 2762e68589bSMark F. Adams nTri[cid]++; 2772e68589bSMark F. Adams } 2782e68589bSMark F. Adams } 2792e68589bSMark F. Adams #define EPS 1.e-12 2802e68589bSMark F. Adams /* find points and set prolongation */ 2810cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nFineLoc; mm++) { 282e78576d6SMark F. Adams PetscBool ise; 2835e62d202SMark Adams PetscCall(PetscCDIsEmptyAt(agg_lists_1, mm, &ise)); 284e78576d6SMark F. Adams if (!ise) { 2850cbbd2e1SMark F. Adams const PetscInt lid = mm; 2862e68589bSMark F. Adams PetscScalar AA[3][3]; 2872e68589bSMark F. Adams PetscBLASInt N = 3, NRHS = 1, LDA = 3, IPIV[3], LDB = 3, INFO; 288539c167fSBarry Smith PetscCDIntNd *pos; 289539c167fSBarry Smith 2909566063dSJacob Faibussowitsch PetscCall(PetscCDGetHeadPos(agg_lists_1, lid, &pos)); 291e78576d6SMark F. Adams while (pos) { 292e78576d6SMark F. Adams PetscInt flid; 2939566063dSJacob Faibussowitsch PetscCall(PetscCDIntNdGetID(pos, &flid)); 2949566063dSJacob Faibussowitsch PetscCall(PetscCDGetNextPos(agg_lists_1, lid, &pos)); 2950cbbd2e1SMark F. Adams 2962e68589bSMark F. Adams if (flid < nFineLoc) { /* could be a ghost */ 2979371c9d4SSatish Balay PetscInt bestTID = -1; 2989371c9d4SSatish Balay PetscReal best_alpha = 1.e10; 2992e68589bSMark F. Adams const PetscInt fgid = flid + myFine0; 3002e68589bSMark F. Adams /* compute shape function for gid */ 301a2f3521dSMark F. Adams const PetscReal fcoord[3] = {coords[flid], coords[data_stride + flid], 1.0}; 3029371c9d4SSatish Balay PetscBool haveit = PETSC_FALSE; 3039371c9d4SSatish Balay PetscScalar alpha[3]; 3049371c9d4SSatish Balay PetscInt clids[3]; 3052fa5cd67SKarl Rupp 3062e68589bSMark F. Adams /* look for it */ 3079371c9d4SSatish Balay for (tid = node_tri[clid], jj = 0; jj < 5 && !haveit && tid != -1; jj++) { 3082e68589bSMark F. Adams for (tt = 0; tt < 3; tt++) { 3092e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3 * tid + tt]; 3102e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 3119371c9d4SSatish Balay AA[tt][0] = coords[lid2]; 3129371c9d4SSatish Balay AA[tt][1] = coords[data_stride + lid2]; 3139371c9d4SSatish Balay AA[tt][2] = 1.0; 3142e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3152e68589bSMark F. Adams } 3162e68589bSMark F. Adams 3172e68589bSMark F. Adams for (tt = 0; tt < 3; tt++) alpha[tt] = (PetscScalar)fcoord[tt]; 3182e68589bSMark F. Adams 3192e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 320792fecdfSBarry Smith PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&N, &NRHS, (PetscScalar *)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3212e68589bSMark F. Adams { 3229371c9d4SSatish Balay PetscBool have = PETSC_TRUE; 3239371c9d4SSatish Balay PetscReal lowest = 1.e10; 3242e68589bSMark F. Adams for (tt = 0, idx = 0; tt < 3; tt++) { 3252e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) > (1.0 + EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE; 3262e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) < lowest) { 3272e68589bSMark F. Adams lowest = PetscRealPart(alpha[tt]); 3282e68589bSMark F. Adams idx = tt; 3292e68589bSMark F. Adams } 3302e68589bSMark F. Adams } 3312e68589bSMark F. Adams haveit = have; 3322e68589bSMark F. Adams } 3332e68589bSMark F. Adams tid = mid.neighborlist[3 * tid + idx]; 3342e68589bSMark F. Adams } 3352e68589bSMark F. Adams 3362e68589bSMark F. Adams if (!haveit) { 3372e68589bSMark F. Adams /* brute force */ 3382e68589bSMark F. Adams for (tid = 0; tid < mid.numberoftriangles && !haveit; tid++) { 3392e68589bSMark F. Adams for (tt = 0; tt < 3; tt++) { 3402e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3 * tid + tt]; 3412e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 3429371c9d4SSatish Balay AA[tt][0] = coords[lid2]; 3439371c9d4SSatish Balay AA[tt][1] = coords[data_stride + lid2]; 3449371c9d4SSatish Balay AA[tt][2] = 1.0; 3452e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3462e68589bSMark F. Adams } 3472e68589bSMark F. Adams for (tt = 0; tt < 3; tt++) alpha[tt] = fcoord[tt]; 3482e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 349792fecdfSBarry Smith PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&N, &NRHS, (PetscScalar *)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3502e68589bSMark F. Adams { 3519371c9d4SSatish Balay PetscBool have = PETSC_TRUE; 3529371c9d4SSatish Balay PetscReal worst = 0.0, v; 3532e68589bSMark F. Adams for (tt = 0; tt < 3 && have; tt++) { 3542e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) > 1.0 + EPS || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE; 3552e68589bSMark F. Adams if ((v = PetscAbs(PetscRealPart(alpha[tt]) - 0.5)) > worst) worst = v; 3562e68589bSMark F. Adams } 3572e68589bSMark F. Adams if (worst < best_alpha) { 3589371c9d4SSatish Balay best_alpha = worst; 3599371c9d4SSatish Balay bestTID = tid; 3602e68589bSMark F. Adams } 3612e68589bSMark F. Adams haveit = have; 3622e68589bSMark F. Adams } 3632e68589bSMark F. Adams } 3642e68589bSMark F. Adams } 3652e68589bSMark F. Adams if (!haveit) { 3662e68589bSMark F. Adams if (best_alpha > *a_worst_best) *a_worst_best = best_alpha; 3672e68589bSMark F. Adams /* use best one */ 3682e68589bSMark F. Adams for (tt = 0; tt < 3; tt++) { 3692e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3 * bestTID + tt]; 3702e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 3719371c9d4SSatish Balay AA[tt][0] = coords[lid2]; 3729371c9d4SSatish Balay AA[tt][1] = coords[data_stride + lid2]; 3739371c9d4SSatish Balay AA[tt][2] = 1.0; 3742e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3752e68589bSMark F. Adams } 3762e68589bSMark F. Adams for (tt = 0; tt < 3; tt++) alpha[tt] = fcoord[tt]; 3772e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 378792fecdfSBarry Smith PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&N, &NRHS, (PetscScalar *)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3792e68589bSMark F. Adams } 3802e68589bSMark F. Adams 3812e68589bSMark F. Adams /* put in row of P */ 3822e68589bSMark F. Adams for (idx = 0; idx < 3; idx++) { 3832e68589bSMark F. Adams PetscScalar shp = alpha[idx]; 3842e68589bSMark F. Adams if (PetscAbs(PetscRealPart(shp)) > 1.e-6) { 3852e68589bSMark F. Adams PetscInt cgid = crsGID[clids[idx]]; 3862e68589bSMark F. Adams PetscInt jj = cgid * bs, ii = fgid * bs; /* need to gloalize */ 38748a46eb9SPierre Jolivet for (tt = 0; tt < bs; tt++, ii++, jj++) PetscCall(MatSetValues(a_Prol, 1, &ii, 1, &jj, &shp, INSERT_VALUES)); 3882e68589bSMark F. Adams } 3892e68589bSMark F. Adams } 3902e68589bSMark F. Adams } 3910cbbd2e1SMark F. Adams } /* aggregates iterations */ 3920cbbd2e1SMark F. Adams clid++; 3930cbbd2e1SMark F. Adams } /* a coarse agg */ 3940cbbd2e1SMark F. Adams } /* for all fine nodes */ 3950cbbd2e1SMark F. Adams 3969566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected_2, &selected_idx_2)); 3979566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY)); 3989566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY)); 3992e68589bSMark F. Adams 4009566063dSJacob Faibussowitsch PetscCall(PetscFree2(node_tri, nTri)); 4012e68589bSMark F. Adams } 4022e68589bSMark F. Adams free(mid.trianglelist); 4032e68589bSMark F. Adams free(mid.neighborlist); 404c6b50ea1SMark free(mid.segmentlist); 405c6b50ea1SMark free(mid.segmentmarkerlist); 406c6b50ea1SMark free(mid.pointlist); 407c6b50ea1SMark free(mid.pointmarkerlist); 4089566063dSJacob Faibussowitsch PetscCall(PetscFree(in.pointlist)); 4093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4102e68589bSMark F. Adams #else 4113b4367a7SBarry Smith SETERRQ(PetscObjectComm((PetscObject)a_Prol), PETSC_ERR_PLIB, "configure with TRIANGLE to use geometric MG"); 4122e68589bSMark F. Adams #endif 4132e68589bSMark F. Adams } 414f1580f4eSBarry Smith 41504c3f3b8SBarry Smith // PetscClangLinter pragma disable: -fdoc-sowing-chars 4162e68589bSMark F. Adams /* 4172e68589bSMark F. Adams getGIDsOnSquareGraph - square graph, get 4182e68589bSMark F. Adams 4192e68589bSMark F. Adams Input Parameter: 4200cbbd2e1SMark F. Adams . nselected_1 - selected local indices (includes ghosts in input Gmat1) 421b43b03e9SMark F. Adams . clid_lid_1 - [nselected_1] lids of selected nodes 4222e68589bSMark F. Adams . Gmat1 - graph that goes with 'selected_1' 4232e68589bSMark F. Adams Output Parameter: 4242e68589bSMark F. Adams . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2) 4252e68589bSMark F. Adams . a_Gmat_2 - graph that is squared of 'Gmat_1' 4262e68589bSMark F. Adams . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes 4272e68589bSMark F. Adams */ 428d71ae5a4SJacob Faibussowitsch static PetscErrorCode getGIDsOnSquareGraph(PC pc, PetscInt nselected_1, const PetscInt clid_lid_1[], const Mat Gmat1, IS *a_selected_2, Mat *a_Gmat_2, PetscInt **a_crsGID) 429d71ae5a4SJacob Faibussowitsch { 43073911c69SBarry Smith PetscMPIInt size; 431b43b03e9SMark F. Adams PetscInt *crsGID, kk, my0, Iend, nloc; 4323b4367a7SBarry Smith MPI_Comm comm; 4332e68589bSMark F. Adams 4342e68589bSMark F. Adams PetscFunctionBegin; 4359566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm)); 4369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 4379566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat1, &my0, &Iend)); /* AIJ */ 4382e68589bSMark F. Adams nloc = Iend - my0; /* this does not change */ 4392e68589bSMark F. Adams 440c5df96a5SBarry Smith if (size == 1) { /* not much to do in serial */ 4419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nselected_1, &crsGID)); 442b43b03e9SMark F. Adams for (kk = 0; kk < nselected_1; kk++) crsGID[kk] = kk; 4430a545947SLisandro Dalcin *a_Gmat_2 = NULL; 4449566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nselected_1, clid_lid_1, PETSC_COPY_VALUES, a_selected_2)); 445806fa848SBarry Smith } else { 446b43b03e9SMark F. Adams PetscInt idx, num_fine_ghosts, num_crs_ghost, myCrs0; 4472e68589bSMark F. Adams Mat_MPIAIJ *mpimat2; 4482e68589bSMark F. Adams Mat Gmat2; 4492e68589bSMark F. Adams Vec locState; 4502e68589bSMark F. Adams PetscScalar *cpcol_state; 4512e68589bSMark F. Adams 4522e68589bSMark F. Adams /* scan my coarse zero gid, set 'lid_state' with coarse GID */ 453b43b03e9SMark F. Adams kk = nselected_1; 4549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPI_SUM, comm)); 455b43b03e9SMark F. Adams myCrs0 -= nselected_1; 4562e68589bSMark F. Adams 457b43b03e9SMark F. Adams if (a_Gmat_2) { /* output */ 4582e68589bSMark F. Adams /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */ 4599566063dSJacob Faibussowitsch PetscCall(PCGAMGSquareGraph_GAMG(pc, Gmat1, &Gmat2)); 4602e68589bSMark F. Adams *a_Gmat_2 = Gmat2; /* output */ 461806fa848SBarry Smith } else Gmat2 = Gmat1; /* use local to get crsGIDs at least */ 4622e68589bSMark F. Adams /* get coarse grid GIDS for selected (locals and ghosts) */ 4632e68589bSMark F. Adams mpimat2 = (Mat_MPIAIJ *)Gmat2->data; 4649566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Gmat2, &locState, NULL)); 4659566063dSJacob Faibussowitsch PetscCall(VecSet(locState, (PetscScalar)(PetscReal)(-1))); /* set with UNKNOWN state */ 466b43b03e9SMark F. Adams for (kk = 0; kk < nselected_1; kk++) { 467b43b03e9SMark F. Adams PetscInt fgid = clid_lid_1[kk] + my0; 4682e68589bSMark F. Adams PetscScalar v = (PetscScalar)(kk + myCrs0); 4699566063dSJacob Faibussowitsch PetscCall(VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES)); /* set with PID */ 4702e68589bSMark F. Adams } 4719566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(locState)); 4729566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(locState)); 4739566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mpimat2->Mvctx, locState, mpimat2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 4749566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mpimat2->Mvctx, locState, mpimat2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 4759566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts)); 4769566063dSJacob Faibussowitsch PetscCall(VecGetArray(mpimat2->lvec, &cpcol_state)); 4772e68589bSMark F. Adams for (kk = 0, num_crs_ghost = 0; kk < num_fine_ghosts; kk++) { 4782e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++; 4792e68589bSMark F. Adams } 4809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nselected_1 + num_crs_ghost, &crsGID)); /* output */ 4812e68589bSMark F. Adams { 4822e68589bSMark F. Adams PetscInt *selected_set; 4839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nselected_1 + num_crs_ghost, &selected_set)); 4842e68589bSMark F. Adams /* do ghost of 'crsGID' */ 485b43b03e9SMark F. Adams for (kk = 0, idx = nselected_1; kk < num_fine_ghosts; kk++) { 4862e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) { 4872e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 4882e68589bSMark F. Adams selected_set[idx] = nloc + kk; 4892e68589bSMark F. Adams crsGID[idx++] = cgid; 4902e68589bSMark F. Adams } 4912e68589bSMark F. Adams } 49263a3b9bcSJacob Faibussowitsch PetscCheck(idx == (nselected_1 + num_crs_ghost), PETSC_COMM_SELF, PETSC_ERR_PLIB, "idx %" PetscInt_FMT " != (nselected_1 %" PetscInt_FMT " + num_crs_ghost %" PetscInt_FMT ")", idx, nselected_1, num_crs_ghost); 4939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mpimat2->lvec, &cpcol_state)); 4942e68589bSMark F. Adams /* do locals in 'crsGID' */ 4959566063dSJacob Faibussowitsch PetscCall(VecGetArray(locState, &cpcol_state)); 4962e68589bSMark F. Adams for (kk = 0, idx = 0; kk < nloc; kk++) { 4972e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) { 4982e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 4992e68589bSMark F. Adams selected_set[idx] = kk; 5002e68589bSMark F. Adams crsGID[idx++] = cgid; 5012e68589bSMark F. Adams } 5022e68589bSMark F. Adams } 50363a3b9bcSJacob Faibussowitsch PetscCheck(idx == nselected_1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "idx %" PetscInt_FMT " != nselected_1 %" PetscInt_FMT, idx, nselected_1); 5049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(locState, &cpcol_state)); 5052e68589bSMark F. Adams 5060a545947SLisandro Dalcin if (a_selected_2 != NULL) { /* output */ 5079566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, (nselected_1 + num_crs_ghost), selected_set, PETSC_OWN_POINTER, a_selected_2)); 508806fa848SBarry Smith } else { 5099566063dSJacob Faibussowitsch PetscCall(PetscFree(selected_set)); 5102e68589bSMark F. Adams } 5110cbbd2e1SMark F. Adams } 5129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&locState)); 5132e68589bSMark F. Adams } 5142e68589bSMark F. Adams *a_crsGID = crsGID; /* output */ 5153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5162e68589bSMark F. Adams } 5172e68589bSMark F. Adams 51866976f2fSJacob Faibussowitsch static PetscErrorCode PCGAMGCreateGraph_GEO(PC pc, Mat Amat, Mat *a_Gmat) 519d71ae5a4SJacob Faibussowitsch { 520c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 521c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 522c1eae691SMark Adams const PetscReal vfilter = pc_gamg->threshold[0]; 5236e111a19SKarl Rupp 524c8b0795cSMark F. Adams PetscFunctionBegin; 525*e02fb3cdSMark Adams PetscCall(MatCreateGraph(Amat, PETSC_TRUE, PETSC_TRUE, vfilter, 0, NULL, a_Gmat)); 5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 527c8b0795cSMark F. Adams } 528c8b0795cSMark F. Adams 52966976f2fSJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc, Mat *a_Gmat, PetscCoarsenData **a_llist_parent) 530d71ae5a4SJacob Faibussowitsch { 531c8b0795cSMark F. Adams PetscInt Istart, Iend, nloc, kk, Ii, ncols; 5320cbbd2e1SMark F. Adams IS perm; 533c8b0795cSMark F. Adams GAMGNode *gnodes; 534c8b0795cSMark F. Adams PetscInt *permute; 535e0940f08SMark F. Adams Mat Gmat = *a_Gmat; 5363b4367a7SBarry Smith MPI_Comm comm; 537b43b03e9SMark F. Adams MatCoarsen crs; 538c8b0795cSMark F. Adams 539c8b0795cSMark F. Adams PetscFunctionBegin; 5409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)a_pc, &comm)); 541849bee69SMark Adams 5429566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend)); 543c8b0795cSMark F. Adams nloc = (Iend - Istart); 544c8b0795cSMark F. Adams 545c8b0795cSMark F. Adams /* create random permutation with sort for geo-mg */ 5469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &gnodes)); 5479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &permute)); 548c8b0795cSMark F. Adams 549c8b0795cSMark F. Adams for (Ii = Istart; Ii < Iend; Ii++) { /* locals only? */ 5509566063dSJacob Faibussowitsch PetscCall(MatGetRow(Gmat, Ii, &ncols, NULL, NULL)); 551c8b0795cSMark F. Adams { 552c8b0795cSMark F. Adams PetscInt lid = Ii - Istart; 553c8b0795cSMark F. Adams gnodes[lid].lid = lid; 554c8b0795cSMark F. Adams gnodes[lid].degree = ncols; 555c8b0795cSMark F. Adams } 5569566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Gmat, Ii, &ncols, NULL, NULL)); 557c8b0795cSMark F. Adams } 558c8b0795cSMark F. Adams if (PETSC_TRUE) { 559e2c00dcbSBarry Smith PetscRandom rand; 560c8b0795cSMark F. Adams PetscBool *bIndexSet; 561e2c00dcbSBarry Smith PetscReal rr; 562e2c00dcbSBarry Smith PetscInt iSwapIndex; 563e2c00dcbSBarry Smith 5649566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rand)); 5659566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nloc, &bIndexSet)); 5662fa5cd67SKarl Rupp for (Ii = 0; Ii < nloc; Ii++) { 5679566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &rr)); 568e2c00dcbSBarry Smith iSwapIndex = (PetscInt)(rr * nloc); 5692fa5cd67SKarl Rupp if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 570c8b0795cSMark F. Adams GAMGNode iTemp = gnodes[iSwapIndex]; 571c8b0795cSMark F. Adams gnodes[iSwapIndex] = gnodes[Ii]; 572c8b0795cSMark F. Adams gnodes[Ii] = iTemp; 573c8b0795cSMark F. Adams bIndexSet[Ii] = PETSC_TRUE; 574c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 575c8b0795cSMark F. Adams } 576c8b0795cSMark F. Adams } 5779566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 5789566063dSJacob Faibussowitsch PetscCall(PetscFree(bIndexSet)); 579c8b0795cSMark F. Adams } 580c8b0795cSMark F. Adams /* only sort locals */ 5810cbbd2e1SMark F. Adams qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare); 582c8b0795cSMark F. Adams /* create IS of permutation */ 5832fa5cd67SKarl Rupp for (kk = 0; kk < nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */ 5849566063dSJacob Faibussowitsch PetscCall(PetscFree(gnodes)); 5859566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm)); 586c8b0795cSMark F. Adams 587c8b0795cSMark F. Adams /* get MIS aggs */ 588b43b03e9SMark F. Adams 5899566063dSJacob Faibussowitsch PetscCall(MatCoarsenCreate(comm, &crs)); 5909566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetType(crs, MATCOARSENMIS)); 5919566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetGreedyOrdering(crs, perm)); 5929566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetAdjacency(crs, Gmat)); 5939566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_FALSE)); 5949566063dSJacob Faibussowitsch PetscCall(MatCoarsenApply(crs)); 5959566063dSJacob Faibussowitsch PetscCall(MatCoarsenGetData(crs, a_llist_parent)); 5969566063dSJacob Faibussowitsch PetscCall(MatCoarsenDestroy(&crs)); 597c8b0795cSMark F. Adams 5989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 599849bee69SMark Adams 6003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 601c8b0795cSMark F. Adams } 602c8b0795cSMark F. Adams 60366976f2fSJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_GEO(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out) 604d71ae5a4SJacob Faibussowitsch { 6052e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 6062e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 607c8b0795cSMark F. Adams const PetscInt dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols; 608b43b03e9SMark F. Adams PetscInt Istart, Iend, nloc, my0, jj, kk, ncols, nLocalSelected, bs, *clid_flid; 609c8b0795cSMark F. Adams Mat Prol; 610c5df96a5SBarry Smith PetscMPIInt rank, size; 6113b4367a7SBarry Smith MPI_Comm comm; 6120cbbd2e1SMark F. Adams IS selected_2, selected_1; 6132e68589bSMark F. Adams const PetscInt *selected_idx; 614d9558ea9SBarry Smith MatType mtype; 6152e68589bSMark F. Adams 6162e68589bSMark F. Adams PetscFunctionBegin; 6179566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 618849bee69SMark Adams 6199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 6209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 6219566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend)); 6229566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat, &bs)); 6239371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 6249371c9d4SSatish Balay my0 = Istart / bs; 62563a3b9bcSJacob Faibussowitsch PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT ") %% bs %" PetscInt_FMT, Iend, Istart, bs); 6262e68589bSMark F. Adams 6272e68589bSMark F. Adams /* get 'nLocalSelected' */ 6285e62d202SMark Adams PetscCall(PetscCDGetNonemptyIS(agg_lists, &selected_1)); 6299566063dSJacob Faibussowitsch PetscCall(ISGetSize(selected_1, &jj)); 6309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jj, &clid_flid)); 6319566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected_1, &selected_idx)); 632c8b0795cSMark F. Adams for (kk = 0, nLocalSelected = 0; kk < jj; kk++) { 6332e68589bSMark F. Adams PetscInt lid = selected_idx[kk]; 634b43b03e9SMark F. Adams if (lid < nloc) { 6359566063dSJacob Faibussowitsch PetscCall(MatGetRow(Gmat, lid + my0, &ncols, NULL, NULL)); 636aaa8cc7dSPierre Jolivet if (ncols > 1) clid_flid[nLocalSelected++] = lid; /* filter out singletons */ 6379566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Gmat, lid + my0, &ncols, NULL, NULL)); 638b43b03e9SMark F. Adams } 6392e68589bSMark F. Adams } 6409566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected_1, &selected_idx)); 6419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected_1)); /* this is selected_1 in serial */ 6422e68589bSMark F. Adams 643d9558ea9SBarry Smith /* create prolongator matrix */ 6449566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat, &mtype)); 6459566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Prol)); 6469566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * bs, PETSC_DETERMINE, PETSC_DETERMINE)); 6479566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(Prol, bs, bs)); 6489566063dSJacob Faibussowitsch PetscCall(MatSetType(Prol, mtype)); 6499566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Prol, 3 * data_cols, NULL)); 6509566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Prol, 3 * data_cols, NULL, 3 * data_cols, NULL)); 6512e68589bSMark F. Adams 652c8b0795cSMark F. Adams /* can get all points "removed" - but not on geomg */ 6539566063dSJacob Faibussowitsch PetscCall(MatGetSize(Prol, &kk, &jj)); 6547f66b68fSMark Adams if (!jj) { 6559566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "ERROE: no selected points on coarse grid\n")); 6569566063dSJacob Faibussowitsch PetscCall(PetscFree(clid_flid)); 6579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 6580298fd71SBarry Smith *a_P_out = NULL; /* out */ 6593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6602e68589bSMark F. Adams } 6612e68589bSMark F. Adams 6622e68589bSMark F. Adams { 6632e68589bSMark F. Adams PetscReal *coords; 664a2f3521dSMark F. Adams PetscInt data_stride; 6650298fd71SBarry Smith PetscInt *crsGID = NULL; 6662e68589bSMark F. Adams Mat Gmat2; 6672e68589bSMark F. Adams 66863a3b9bcSJacob Faibussowitsch PetscCheck(dim == data_cols, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " != data_cols %" PetscInt_FMT, dim, data_cols); 6692e68589bSMark F. Adams /* grow ghost data for better coarse grid cover of fine grid */ 670a2f3521dSMark F. Adams /* messy method, squares graph and gets some data */ 6719566063dSJacob Faibussowitsch PetscCall(getGIDsOnSquareGraph(pc, nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID)); 6722e68589bSMark F. Adams /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */ 6732e68589bSMark F. Adams /* create global vector of coorindates in 'coords' */ 674c5df96a5SBarry Smith if (size > 1) { 6759566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords)); 676806fa848SBarry Smith } else { 677c8b0795cSMark F. Adams coords = (PetscReal *)pc_gamg->data; 678a2f3521dSMark F. Adams data_stride = pc_gamg->data_sz / pc_gamg->data_cell_cols; 6792e68589bSMark F. Adams } 6809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Gmat2)); 6812e68589bSMark F. Adams 6822e68589bSMark F. Adams /* triangulate */ 6830fdf79fbSJacob Faibussowitsch { 684c8b0795cSMark F. Adams PetscReal metric, tm; 6850fdf79fbSJacob Faibussowitsch 6860fdf79fbSJacob Faibussowitsch PetscCheck(dim == 2, comm, PETSC_ERR_PLIB, "3D not implemented for 'geo' AMG"); 6879566063dSJacob Faibussowitsch PetscCall(triangulateAndFormProl(selected_2, data_stride, coords, nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric)); 6889566063dSJacob Faibussowitsch PetscCall(PetscFree(crsGID)); 6892e68589bSMark F. Adams 6902e68589bSMark F. Adams /* clean up and create coordinates for coarse grid (output) */ 6919566063dSJacob Faibussowitsch if (size > 1) PetscCall(PetscFree(coords)); 6922e68589bSMark F. Adams 6931c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm)); 694c8b0795cSMark F. Adams if (tm > 1.) { /* needs to be globalized - should not happen */ 6959566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, " failed metric for coarse grid %e\n", (double)tm)); 6969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 697806fa848SBarry Smith } else if (metric > .0) { 6989566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "worst metric for coarse grid = %e\n", (double)metric)); 6992e68589bSMark F. Adams } 7000fdf79fbSJacob Faibussowitsch } 7012e68589bSMark F. Adams { /* create next coords - output */ 7022e68589bSMark F. Adams PetscReal *crs_crds; 7039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * nLocalSelected, &crs_crds)); 7042e68589bSMark F. Adams for (kk = 0; kk < nLocalSelected; kk++) { /* grab local select nodes to promote - output */ 705b43b03e9SMark F. Adams PetscInt lid = clid_flid[kk]; 706c8b0795cSMark F. Adams for (jj = 0; jj < dim; jj++) crs_crds[jj * nLocalSelected + kk] = pc_gamg->data[jj * nloc + lid]; 7072e68589bSMark F. Adams } 708c8b0795cSMark F. Adams 7099566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 710c8b0795cSMark F. Adams pc_gamg->data = crs_crds; /* out */ 711c8b0795cSMark F. Adams pc_gamg->data_sz = dim * nLocalSelected; 7122e68589bSMark F. Adams } 7139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected_2)); 7142e68589bSMark F. Adams } 715a2f3521dSMark F. Adams 7162e68589bSMark F. Adams *a_P_out = Prol; /* out */ 7179566063dSJacob Faibussowitsch PetscCall(PetscFree(clid_flid)); 718849bee69SMark Adams 7193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 720c8b0795cSMark F. Adams } 721c8b0795cSMark F. Adams 722d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_GEO(PC pc) 723d71ae5a4SJacob Faibussowitsch { 7249b8ffb57SJed Brown PetscFunctionBegin; 7259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 7263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7279b8ffb57SJed Brown } 7289b8ffb57SJed Brown 729d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_GEO(PC pc) 730d71ae5a4SJacob Faibussowitsch { 731c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 732c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 733c8b0795cSMark F. Adams 734c8b0795cSMark F. Adams PetscFunctionBegin; 7351ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO; 7369b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_GEO; 737c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 738c8b0795cSMark F. Adams 739c8b0795cSMark F. Adams /* set internal function pointers */ 7402d776b49SBarry Smith pc_gamg->ops->creategraph = PCGAMGCreateGraph_GEO; 741fd1112cbSBarry Smith pc_gamg->ops->coarsen = PCGAMGCoarsen_GEO; 7421ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_GEO; 743fd1112cbSBarry Smith pc_gamg->ops->optprolongator = NULL; 7441ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_GEO; 745c8b0795cSMark F. Adams 7469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_GEO)); 7473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7482e68589bSMark F. Adams } 749