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 229fbee547SJacob Faibussowitsch static inline int petsc_geo_mg_compare(const void *a, const void *b) 23c8b0795cSMark F. Adams { 24c8b0795cSMark F. Adams return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree); 25c8b0795cSMark F. Adams } 26c8b0795cSMark F. Adams 272e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 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 */ 34302f38e8SMark F. Adams PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 352e68589bSMark F. Adams { 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 482c71b3e2SJacob Faibussowitsch PetscCheckFalse(nloc!=a_nloc && aloc!=a_nloc,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of local blocks %D must be %D or %D.",a_nloc,nloc,aloc); 492e68589bSMark F. Adams 50c8b0795cSMark F. Adams pc_gamg->data_cell_rows = 1; 51*7827d75bSBarry 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++) { 662e68589bSMark F. Adams for (ii = 0; ii < ndm; ii++) { 672e68589bSMark F. Adams pc_gamg->data[ii*nloc + kk] = coords[kk*ndm + ii]; 682e68589bSMark F. Adams } 692e68589bSMark F. Adams } 7090fbc344SStefano Zampini } else { /* assumes the coordinates are blocked */ 7190fbc344SStefano Zampini for (kk = 0; kk < nloc; kk++) { 7290fbc344SStefano Zampini for (ii = 0; ii < ndm; ii++) { 7390fbc344SStefano Zampini pc_gamg->data[ii*nloc + kk] = coords[bs*kk*ndm + ii]; 7490fbc344SStefano Zampini } 7590fbc344SStefano Zampini } 7690fbc344SStefano Zampini } 772c71b3e2SJacob Faibussowitsch PetscCheckFalse(pc_gamg->data[arrsz] != -99.,PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data[arrsz %D] %g != -99.",arrsz,pc_gamg->data[arrsz]); 782e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 792e68589bSMark F. Adams PetscFunctionReturn(0); 802e68589bSMark F. Adams } 812e68589bSMark F. Adams 822e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 832e68589bSMark F. Adams /* 842e68589bSMark F. Adams PCSetData_GEO 852e68589bSMark F. Adams 862e68589bSMark F. Adams Input Parameter: 872e68589bSMark F. Adams . pc - 882e68589bSMark F. Adams */ 89b8cd405aSMark F. Adams PetscErrorCode PCSetData_GEO(PC pc, Mat m) 902e68589bSMark F. Adams { 912e68589bSMark F. Adams PetscFunctionBegin; 92ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"GEO MG needs coordinates"); 932e68589bSMark F. Adams } 942e68589bSMark F. Adams 952e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 962e68589bSMark F. Adams /* 972e68589bSMark F. Adams PCSetFromOptions_GEO 982e68589bSMark F. Adams 992e68589bSMark F. Adams Input Parameter: 1002e68589bSMark F. Adams . pc - 1012e68589bSMark F. Adams */ 1024416b707SBarry Smith PetscErrorCode PCSetFromOptions_GEO(PetscOptionItems *PetscOptionsObject,PC pc) 1032e68589bSMark F. Adams { 1042e68589bSMark F. Adams PetscFunctionBegin; 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"GAMG-GEO options")); 1062e68589bSMark F. Adams { 1072e68589bSMark F. Adams /* -pc_gamg_sa_nsmooths */ 1082e68589bSMark F. Adams /* pc_gamg_sa->smooths = 0; */ 1092e68589bSMark F. Adams /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */ 1102e68589bSMark F. Adams /* "smoothing steps for smoothed aggregation, usually 1 (0)", */ 1112e68589bSMark F. Adams /* "PCGAMGSetNSmooths_AGG", */ 1122e68589bSMark F. Adams /* pc_gamg_sa->smooths, */ 1132e68589bSMark F. Adams /* &pc_gamg_sa->smooths, */ 1142e68589bSMark F. Adams /* &flag); */ 1152e68589bSMark F. Adams } 1169566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 1172e68589bSMark F. Adams PetscFunctionReturn(0); 1182e68589bSMark F. Adams } 1192e68589bSMark F. Adams 1202e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1212e68589bSMark F. Adams /* 1222e68589bSMark F. Adams triangulateAndFormProl 1232e68589bSMark F. Adams 1242e68589bSMark F. Adams Input Parameter: 1252e68589bSMark F. Adams . selected_2 - list of selected local ID, includes selected ghosts 126a2f3521dSMark F. Adams . data_stride - 127a2f3521dSMark F. Adams . coords[2*data_stride] - column vector of local coordinates w/ ghosts 1282adfe9d3SBarry Smith . nselected_1 - selected IDs that go with base (1) graph includes selected ghosts 1290cbbd2e1SMark F. Adams . clid_lid_1[nselected_1] - lids of selected (c) nodes ??????????? 1302adfe9d3SBarry Smith . agg_lists_1 - list of aggregates selected_1 vertices of aggregate unselected vertices 1312e68589bSMark F. Adams . crsGID[selected.size()] - global index for prolongation operator 1322e68589bSMark F. Adams . bs - block size 1332e68589bSMark F. Adams Output Parameter: 1342e68589bSMark F. Adams . a_Prol - prolongation operator 1352e68589bSMark F. Adams . a_worst_best - measure of worst missed fine vertex, 0 is no misses 1362e68589bSMark F. Adams */ 1372adfe9d3SBarry Smith static PetscErrorCode triangulateAndFormProl(IS selected_2,PetscInt data_stride,PetscReal coords[],PetscInt nselected_1,const PetscInt clid_lid_1[],const PetscCoarsenData *agg_lists_1, 1382adfe9d3SBarry Smith const PetscInt crsGID[],PetscInt bs,Mat a_Prol,PetscReal *a_worst_best) 1392e68589bSMark F. Adams { 1402e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE) 141b43b03e9SMark F. Adams PetscInt jj,tid,tt,idx,nselected_2; 1422e68589bSMark F. Adams struct triangulateio in,mid; 1430cbbd2e1SMark F. Adams const PetscInt *selected_idx_2; 14473911c69SBarry Smith PetscMPIInt rank; 1452e68589bSMark F. Adams PetscInt Istart,Iend,nFineLoc,myFine0; 1462e68589bSMark F. Adams int kk,nPlotPts,sid; 1473b4367a7SBarry Smith MPI_Comm comm; 148c8b0795cSMark F. Adams PetscReal tm; 149c8b0795cSMark F. Adams 1506e111a19SKarl Rupp PetscFunctionBegin; 1519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)a_Prol,&comm)); 1529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 1539566063dSJacob Faibussowitsch PetscCall(ISGetSize(selected_2, &nselected_2)); 1542e68589bSMark F. Adams if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */ 1552e68589bSMark F. Adams *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */ 156806fa848SBarry Smith } else *a_worst_best = 0.0; 1579566063dSJacob Faibussowitsch PetscCallMPI(MPIU_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm)); 158c8b0795cSMark F. Adams if (tm > 0.0) { 159c8b0795cSMark F. Adams *a_worst_best = 100.0; 1602e68589bSMark F. Adams PetscFunctionReturn(0); 1612e68589bSMark F. Adams } 1629566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend)); 1632e68589bSMark F. Adams nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs; 1642e68589bSMark F. Adams nPlotPts = nFineLoc; /* locals */ 1654c500f23SPierre Jolivet /* triangle */ 1662e68589bSMark F. Adams /* Define input points - in*/ 1672e68589bSMark F. Adams in.numberofpoints = nselected_2; 1682e68589bSMark F. Adams in.numberofpointattributes = 0; 1692e68589bSMark F. Adams /* get nselected points */ 1709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2*nselected_2, &in.pointlist)); 1719566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected_2, &selected_idx_2)); 1722e68589bSMark F. Adams 1732e68589bSMark F. Adams for (kk=0,sid=0; kk<nselected_2; kk++,sid += 2) { 1742e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 1752e68589bSMark F. Adams in.pointlist[sid] = coords[lid]; 176a2f3521dSMark F. Adams in.pointlist[sid+1] = coords[data_stride + lid]; 1772e68589bSMark F. Adams if (lid>=nFineLoc) nPlotPts++; 1782e68589bSMark F. Adams } 1792c71b3e2SJacob Faibussowitsch PetscCheckFalse(sid != 2*nselected_2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != 2*nselected_2 %D",sid,nselected_2); 1802e68589bSMark F. Adams 1812e68589bSMark F. Adams in.numberofsegments = 0; 1822e68589bSMark F. Adams in.numberofedges = 0; 1832e68589bSMark F. Adams in.numberofholes = 0; 1842e68589bSMark F. Adams in.numberofregions = 0; 1850a545947SLisandro Dalcin in.trianglelist = NULL; 1860a545947SLisandro Dalcin in.segmentmarkerlist = NULL; 1870a545947SLisandro Dalcin in.pointattributelist = NULL; 1880a545947SLisandro Dalcin in.pointmarkerlist = NULL; 1890a545947SLisandro Dalcin in.triangleattributelist = NULL; 1900a545947SLisandro Dalcin in.trianglearealist = NULL; 1910a545947SLisandro Dalcin in.segmentlist = NULL; 1920a545947SLisandro Dalcin in.holelist = NULL; 1930a545947SLisandro Dalcin in.regionlist = NULL; 1940a545947SLisandro Dalcin in.edgelist = NULL; 1950a545947SLisandro Dalcin in.edgemarkerlist = NULL; 1960a545947SLisandro Dalcin in.normlist = NULL; 1972fa5cd67SKarl Rupp 1982e68589bSMark F. Adams /* triangulate */ 1990a545947SLisandro Dalcin mid.pointlist = NULL; /* Not needed if -N switch used. */ 2002e68589bSMark F. Adams /* Not needed if -N switch used or number of point attributes is zero: */ 2010a545947SLisandro Dalcin mid.pointattributelist = NULL; 2020a545947SLisandro Dalcin mid.pointmarkerlist = NULL; /* Not needed if -N or -B switch used. */ 2030a545947SLisandro Dalcin mid.trianglelist = NULL; /* Not needed if -E switch used. */ 2042e68589bSMark F. Adams /* Not needed if -E switch used or number of triangle attributes is zero: */ 2050a545947SLisandro Dalcin mid.triangleattributelist = NULL; 2060a545947SLisandro Dalcin mid.neighborlist = NULL; /* Needed only if -n switch used. */ 2072e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P not used: */ 2080a545947SLisandro Dalcin mid.segmentlist = NULL; 2092e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P and -B not used: */ 2100a545947SLisandro Dalcin mid.segmentmarkerlist = NULL; 2110a545947SLisandro Dalcin mid.edgelist = NULL; /* Needed only if -e switch used. */ 2120a545947SLisandro Dalcin mid.edgemarkerlist = NULL; /* Needed if -e used and -B not used. */ 2132e68589bSMark F. Adams mid.numberoftriangles = 0; 2142e68589bSMark F. Adams 2152e68589bSMark F. Adams /* Triangulate the points. Switches are chosen to read and write a */ 2162e68589bSMark F. Adams /* PSLG (p), preserve the convex hull (c), number everything from */ 2172e68589bSMark F. Adams /* zero (z), assign a regional attribute to each element (A), and */ 2182e68589bSMark F. Adams /* produce an edge list (e), a Voronoi diagram (v), and a triangle */ 2192e68589bSMark F. Adams /* neighbor list (n). */ 2202e68589bSMark F. Adams if (nselected_2 != 0) { /* inactive processor */ 2212e68589bSMark F. Adams char args[] = "npczQ"; /* c is needed ? */ 2222e68589bSMark F. Adams triangulate(args, &in, &mid, (struct triangulateio*) NULL); 2232e68589bSMark F. Adams /* output .poly files for 'showme' */ 2242e68589bSMark F. Adams if (!PETSC_TRUE) { 2252e68589bSMark F. Adams static int level = 1; 2262e68589bSMark F. Adams FILE *file; char fname[32]; 2272e68589bSMark F. Adams 228c5df96a5SBarry Smith sprintf(fname,"C%d_%d.poly",level,rank); file = fopen(fname, "w"); 2292e68589bSMark F. Adams /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/ 2302e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); 2312e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2322e68589bSMark F. Adams for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid += 2) { 2332e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2342e68589bSMark F. Adams } 2352e68589bSMark F. Adams /*One line: <# of segments> <# of boundary markers (0 or 1)> */ 2362e68589bSMark F. Adams fprintf(file, "%d %d\n",0,0); 2372e68589bSMark F. Adams /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */ 2382e68589bSMark F. Adams /* One line: <# of holes> */ 2392e68589bSMark F. Adams fprintf(file, "%d\n",0); 2402e68589bSMark F. Adams /* Following lines: <hole #> <x> <y> */ 2412e68589bSMark F. Adams /* Optional line: <# of regional attributes and/or area constraints> */ 2422e68589bSMark F. Adams /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */ 2432e68589bSMark F. Adams fclose(file); 2442e68589bSMark F. Adams 2452e68589bSMark F. Adams /* elems */ 246c5df96a5SBarry Smith sprintf(fname,"C%d_%d.ele",level,rank); file = fopen(fname, "w"); 2472e68589bSMark F. Adams /* First line: <# of triangles> <nodes per triangle> <# of attributes> */ 2482e68589bSMark F. Adams fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0); 2492e68589bSMark F. Adams /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */ 2502e68589bSMark F. Adams for (kk=0,sid=0; kk<mid.numberoftriangles; kk++,sid += 3) { 2512e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]); 2522e68589bSMark F. Adams } 2532e68589bSMark F. Adams fclose(file); 2542e68589bSMark F. Adams 255c5df96a5SBarry Smith sprintf(fname,"C%d_%d.node",level,rank); file = fopen(fname, "w"); 2562e68589bSMark F. Adams /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */ 2572e68589bSMark F. Adams /* fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); */ 2582e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",nPlotPts,2,0,0); 2592e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2602e68589bSMark F. Adams for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid+=2) { 2612e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2622e68589bSMark F. Adams } 2632e68589bSMark F. Adams 2642e68589bSMark F. Adams sid /= 2; 2652e68589bSMark F. Adams for (jj=0; jj<nFineLoc; jj++) { 2662e68589bSMark F. Adams PetscBool sel = PETSC_TRUE; 2672e68589bSMark F. Adams for (kk=0; kk<nselected_2 && sel; kk++) { 2682e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 2692e68589bSMark F. Adams if (lid == jj) sel = PETSC_FALSE; 2702e68589bSMark F. Adams } 2712fa5cd67SKarl Rupp if (sel) fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]); 2722e68589bSMark F. Adams } 2732e68589bSMark F. Adams fclose(file); 2742c71b3e2SJacob Faibussowitsch PetscCheckFalse(sid != nPlotPts,PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != nPlotPts %D",sid,nPlotPts); 2752e68589bSMark F. Adams level++; 2762e68589bSMark F. Adams } 2772e68589bSMark F. Adams } 2789566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0)); 2792e68589bSMark F. Adams { /* form P - setup some maps */ 2800cbbd2e1SMark F. Adams PetscInt clid,mm,*nTri,*node_tri; 2812e68589bSMark F. Adams 2829566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nselected_2, &node_tri,nselected_2, &nTri)); 2832e68589bSMark F. Adams 2842e68589bSMark F. Adams /* need list of triangles on node */ 2852e68589bSMark F. Adams for (kk=0; kk<nselected_2; kk++) nTri[kk] = 0; 2862e68589bSMark F. Adams for (tid=0,kk=0; tid<mid.numberoftriangles; tid++) { 2872e68589bSMark F. Adams for (jj=0; jj<3; jj++) { 2882e68589bSMark F. Adams PetscInt cid = mid.trianglelist[kk++]; 2892e68589bSMark F. Adams if (nTri[cid] == 0) node_tri[cid] = tid; 2902e68589bSMark F. Adams nTri[cid]++; 2912e68589bSMark F. Adams } 2922e68589bSMark F. Adams } 2932e68589bSMark F. Adams #define EPS 1.e-12 2942e68589bSMark F. Adams /* find points and set prolongation */ 2950cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nFineLoc; mm++) { 296e78576d6SMark F. Adams PetscBool ise; 2979566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists_1,mm,&ise)); 298e78576d6SMark F. Adams if (!ise) { 2990cbbd2e1SMark F. Adams const PetscInt lid = mm; 3002e68589bSMark F. Adams PetscScalar AA[3][3]; 3012e68589bSMark F. Adams PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO; 302539c167fSBarry Smith PetscCDIntNd *pos; 303539c167fSBarry Smith 3049566063dSJacob Faibussowitsch PetscCall(PetscCDGetHeadPos(agg_lists_1,lid,&pos)); 305e78576d6SMark F. Adams while (pos) { 306e78576d6SMark F. Adams PetscInt flid; 3079566063dSJacob Faibussowitsch PetscCall(PetscCDIntNdGetID(pos, &flid)); 3089566063dSJacob Faibussowitsch PetscCall(PetscCDGetNextPos(agg_lists_1,lid,&pos)); 3090cbbd2e1SMark F. Adams 3102e68589bSMark F. Adams if (flid < nFineLoc) { /* could be a ghost */ 3112e68589bSMark F. Adams PetscInt bestTID = -1; PetscReal best_alpha = 1.e10; 3122e68589bSMark F. Adams const PetscInt fgid = flid + myFine0; 3132e68589bSMark F. Adams /* compute shape function for gid */ 314a2f3521dSMark F. Adams const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0}; 3152e68589bSMark F. Adams PetscBool haveit =PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3]; 3162fa5cd67SKarl Rupp 3172e68589bSMark F. Adams /* look for it */ 3180cbbd2e1SMark F. Adams for (tid = node_tri[clid], jj=0; 3192e68589bSMark F. Adams jj < 5 && !haveit && tid != -1; 3202e68589bSMark F. Adams jj++) { 3212e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3222e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3232e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 324a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3252e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3262e68589bSMark F. Adams } 3272e68589bSMark F. Adams 3282e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = (PetscScalar)fcoord[tt]; 3292e68589bSMark F. Adams 3302e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 3318b83055fSJed Brown PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3322e68589bSMark F. Adams { 3332e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal lowest=1.e10; 3342e68589bSMark F. Adams for (tt = 0, idx = 0; tt < 3; tt++) { 3352e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE; 3362e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) < lowest) { 3372e68589bSMark F. Adams lowest = PetscRealPart(alpha[tt]); 3382e68589bSMark F. Adams idx = tt; 3392e68589bSMark F. Adams } 3402e68589bSMark F. Adams } 3412e68589bSMark F. Adams haveit = have; 3422e68589bSMark F. Adams } 3432e68589bSMark F. Adams tid = mid.neighborlist[3*tid + idx]; 3442e68589bSMark F. Adams } 3452e68589bSMark F. Adams 3462e68589bSMark F. Adams if (!haveit) { 3472e68589bSMark F. Adams /* brute force */ 3482e68589bSMark F. Adams for (tid=0; tid<mid.numberoftriangles && !haveit; tid++) { 3492e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3502e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3512e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 352a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3532e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3542e68589bSMark F. Adams } 3552e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt]; 3562e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 3578b83055fSJed Brown PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3582e68589bSMark F. Adams { 3592e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal worst=0.0, v; 3602e68589bSMark F. Adams for (tt=0; tt<3 && have; tt++) { 3612e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE; 3622e68589bSMark F. Adams if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v; 3632e68589bSMark F. Adams } 3642e68589bSMark F. Adams if (worst < best_alpha) { 3652e68589bSMark F. Adams best_alpha = worst; bestTID = tid; 3662e68589bSMark F. Adams } 3672e68589bSMark F. Adams haveit = have; 3682e68589bSMark F. Adams } 3692e68589bSMark F. Adams } 3702e68589bSMark F. Adams } 3712e68589bSMark F. Adams if (!haveit) { 3722e68589bSMark F. Adams if (best_alpha > *a_worst_best) *a_worst_best = best_alpha; 3732e68589bSMark F. Adams /* use best one */ 3742e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3752e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*bestTID + tt]; 3762e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 377a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3782e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3792e68589bSMark F. Adams } 3802e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt]; 3812e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 3828b83055fSJed Brown PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3832e68589bSMark F. Adams } 3842e68589bSMark F. Adams 3852e68589bSMark F. Adams /* put in row of P */ 3862e68589bSMark F. Adams for (idx=0; idx<3; idx++) { 3872e68589bSMark F. Adams PetscScalar shp = alpha[idx]; 3882e68589bSMark F. Adams if (PetscAbs(PetscRealPart(shp)) > 1.e-6) { 3892e68589bSMark F. Adams PetscInt cgid = crsGID[clids[idx]]; 3902e68589bSMark F. Adams PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */ 3912e68589bSMark F. Adams for (tt=0; tt < bs; tt++, ii++, jj++) { 3929566063dSJacob Faibussowitsch PetscCall(MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES)); 3932e68589bSMark F. Adams } 3942e68589bSMark F. Adams } 3952e68589bSMark F. Adams } 3962e68589bSMark F. Adams } 3970cbbd2e1SMark F. Adams } /* aggregates iterations */ 3980cbbd2e1SMark F. Adams clid++; 3990cbbd2e1SMark F. Adams } /* a coarse agg */ 4000cbbd2e1SMark F. Adams } /* for all fine nodes */ 4010cbbd2e1SMark F. Adams 4029566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected_2, &selected_idx_2)); 4039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY)); 4049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY)); 4052e68589bSMark F. Adams 4069566063dSJacob Faibussowitsch PetscCall(PetscFree2(node_tri,nTri)); 4072e68589bSMark F. Adams } 4089566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0)); 4092e68589bSMark F. Adams free(mid.trianglelist); 4102e68589bSMark F. Adams free(mid.neighborlist); 411c6b50ea1SMark free(mid.segmentlist); 412c6b50ea1SMark free(mid.segmentmarkerlist); 413c6b50ea1SMark free(mid.pointlist); 414c6b50ea1SMark free(mid.pointmarkerlist); 4159566063dSJacob Faibussowitsch PetscCall(PetscFree(in.pointlist)); 4162e68589bSMark F. Adams PetscFunctionReturn(0); 4172e68589bSMark F. Adams #else 4183b4367a7SBarry Smith SETERRQ(PetscObjectComm((PetscObject)a_Prol),PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG"); 4192e68589bSMark F. Adams #endif 4202e68589bSMark F. Adams } 4212e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 4222e68589bSMark F. Adams /* 4232e68589bSMark F. Adams getGIDsOnSquareGraph - square graph, get 4242e68589bSMark F. Adams 4252e68589bSMark F. Adams Input Parameter: 4260cbbd2e1SMark F. Adams . nselected_1 - selected local indices (includes ghosts in input Gmat1) 427b43b03e9SMark F. Adams . clid_lid_1 - [nselected_1] lids of selected nodes 4282e68589bSMark F. Adams . Gmat1 - graph that goes with 'selected_1' 4292e68589bSMark F. Adams Output Parameter: 4302e68589bSMark F. Adams . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2) 4312e68589bSMark F. Adams . a_Gmat_2 - graph that is squared of 'Gmat_1' 4322e68589bSMark F. Adams . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes 4332e68589bSMark F. Adams */ 4344b1575e2SStefano Zampini 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) 4352e68589bSMark F. Adams { 43673911c69SBarry Smith PetscMPIInt size; 437b43b03e9SMark F. Adams PetscInt *crsGID, kk,my0,Iend,nloc; 4383b4367a7SBarry Smith MPI_Comm comm; 4392e68589bSMark F. Adams 4402e68589bSMark F. Adams PetscFunctionBegin; 4419566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Gmat1,&comm)); 4429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 4439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat1,&my0,&Iend)); /* AIJ */ 4442e68589bSMark F. Adams nloc = Iend - my0; /* this does not change */ 4452e68589bSMark F. Adams 446c5df96a5SBarry Smith if (size == 1) { /* not much to do in serial */ 4479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nselected_1, &crsGID)); 448b43b03e9SMark F. Adams for (kk=0; kk<nselected_1; kk++) crsGID[kk] = kk; 4490a545947SLisandro Dalcin *a_Gmat_2 = NULL; 4509566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2)); 451806fa848SBarry Smith } else { 452b43b03e9SMark F. Adams PetscInt idx,num_fine_ghosts,num_crs_ghost,myCrs0; 4532e68589bSMark F. Adams Mat_MPIAIJ *mpimat2; 4542e68589bSMark F. Adams Mat Gmat2; 4552e68589bSMark F. Adams Vec locState; 4562e68589bSMark F. Adams PetscScalar *cpcol_state; 4572e68589bSMark F. Adams 4582e68589bSMark F. Adams /* scan my coarse zero gid, set 'lid_state' with coarse GID */ 459b43b03e9SMark F. Adams kk = nselected_1; 4609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPI_SUM, comm)); 461b43b03e9SMark F. Adams myCrs0 -= nselected_1; 4622e68589bSMark F. Adams 463b43b03e9SMark F. Adams if (a_Gmat_2) { /* output */ 4642e68589bSMark F. Adams /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */ 4659566063dSJacob Faibussowitsch PetscCall(PCGAMGSquareGraph_GAMG(pc,Gmat1,&Gmat2)); 4662e68589bSMark F. Adams *a_Gmat_2 = Gmat2; /* output */ 467806fa848SBarry Smith } else Gmat2 = Gmat1; /* use local to get crsGIDs at least */ 4682e68589bSMark F. Adams /* get coarse grid GIDS for selected (locals and ghosts) */ 4692e68589bSMark F. Adams mpimat2 = (Mat_MPIAIJ*)Gmat2->data; 4709566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Gmat2, &locState, NULL)); 4719566063dSJacob Faibussowitsch PetscCall(VecSet(locState, (PetscScalar)(PetscReal)(-1))); /* set with UNKNOWN state */ 472b43b03e9SMark F. Adams for (kk=0; kk<nselected_1; kk++) { 473b43b03e9SMark F. Adams PetscInt fgid = clid_lid_1[kk] + my0; 4742e68589bSMark F. Adams PetscScalar v = (PetscScalar)(kk+myCrs0); 4759566063dSJacob Faibussowitsch PetscCall(VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES)); /* set with PID */ 4762e68589bSMark F. Adams } 4779566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(locState)); 4789566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(locState)); 4799566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD)); 4809566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD)); 4819566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts)); 4829566063dSJacob Faibussowitsch PetscCall(VecGetArray(mpimat2->lvec, &cpcol_state)); 4832e68589bSMark F. Adams for (kk=0,num_crs_ghost=0; kk<num_fine_ghosts; kk++) { 4842e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++; 4852e68589bSMark F. Adams } 4869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nselected_1+num_crs_ghost, &crsGID)); /* output */ 4872e68589bSMark F. Adams { 4882e68589bSMark F. Adams PetscInt *selected_set; 4899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nselected_1+num_crs_ghost, &selected_set)); 4902e68589bSMark F. Adams /* do ghost of 'crsGID' */ 491b43b03e9SMark F. Adams for (kk=0,idx=nselected_1; kk<num_fine_ghosts; kk++) { 4922e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) { 4932e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 4942e68589bSMark F. Adams selected_set[idx] = nloc + kk; 4952e68589bSMark F. Adams crsGID[idx++] = cgid; 4962e68589bSMark F. Adams } 4972e68589bSMark F. Adams } 4982c71b3e2SJacob Faibussowitsch PetscCheckFalse(idx != (nselected_1+num_crs_ghost),PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != (nselected_1 %D + num_crs_ghost %D)",idx,nselected_1,num_crs_ghost); 4999566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mpimat2->lvec, &cpcol_state)); 5002e68589bSMark F. Adams /* do locals in 'crsGID' */ 5019566063dSJacob Faibussowitsch PetscCall(VecGetArray(locState, &cpcol_state)); 5022e68589bSMark F. Adams for (kk=0,idx=0; kk<nloc; kk++) { 5032e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) { 5042e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5052e68589bSMark F. Adams selected_set[idx] = kk; 5062e68589bSMark F. Adams crsGID[idx++] = cgid; 5072e68589bSMark F. Adams } 5082e68589bSMark F. Adams } 5092c71b3e2SJacob Faibussowitsch PetscCheckFalse(idx != nselected_1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != nselected_1 %D",idx,nselected_1); 5109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(locState, &cpcol_state)); 5112e68589bSMark F. Adams 5120a545947SLisandro Dalcin if (a_selected_2 != NULL) { /* output */ 5139566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2)); 514806fa848SBarry Smith } else { 5159566063dSJacob Faibussowitsch PetscCall(PetscFree(selected_set)); 5162e68589bSMark F. Adams } 5170cbbd2e1SMark F. Adams } 5189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&locState)); 5192e68589bSMark F. Adams } 5202e68589bSMark F. Adams *a_crsGID = crsGID; /* output */ 5212e68589bSMark F. Adams PetscFunctionReturn(0); 5222e68589bSMark F. Adams } 5232e68589bSMark F. Adams 5242e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 5252e68589bSMark F. Adams /* 526fd1112cbSBarry Smith PCGAMGGraph_GEO 5272e68589bSMark F. Adams 5282e68589bSMark F. Adams Input Parameter: 5292e68589bSMark F. Adams . pc - this 5302e68589bSMark F. Adams . Amat - matrix on this fine level 5312e68589bSMark F. Adams Output Parameter: 532c8b0795cSMark F. Adams . a_Gmat 5332e68589bSMark F. Adams */ 5342adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_GEO(PC pc,Mat Amat,Mat *a_Gmat) 535c8b0795cSMark F. Adams { 536c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 537c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 538c1eae691SMark Adams const PetscReal vfilter = pc_gamg->threshold[0]; 5393b4367a7SBarry Smith MPI_Comm comm; 540c8b0795cSMark F. Adams Mat Gmat; 5410cbbd2e1SMark F. Adams PetscBool set,flg,symm; 5426e111a19SKarl Rupp 543c8b0795cSMark F. Adams PetscFunctionBegin; 5449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm)); 5459566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_GAMGGraph_GEO,0,0,0,0)); 546c8b0795cSMark F. Adams 5479566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(Amat, &set, &flg)); 548263489e9SJed Brown symm = (PetscBool)!(set && flg); 5490cbbd2e1SMark F. Adams 5509566063dSJacob Faibussowitsch PetscCall(PCGAMGCreateGraph(Amat, &Gmat)); 5519566063dSJacob Faibussowitsch PetscCall(PCGAMGFilterGraph(&Gmat, vfilter, symm)); 552c8b0795cSMark F. Adams 553c8b0795cSMark F. Adams *a_Gmat = Gmat; 5549566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_GAMGGraph_GEO,0,0,0,0)); 555c8b0795cSMark F. Adams PetscFunctionReturn(0); 556c8b0795cSMark F. Adams } 557c8b0795cSMark F. Adams 558c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 559c8b0795cSMark F. Adams /* 560fd1112cbSBarry Smith PCGAMGCoarsen_GEO 561c8b0795cSMark F. Adams 562c8b0795cSMark F. Adams Input Parameter: 563e0940f08SMark F. Adams . a_pc - this 564e0940f08SMark F. Adams . a_Gmat - graph 565c8b0795cSMark F. Adams Output Parameter: 566c8b0795cSMark F. Adams . a_llist_parent - linked list from selected indices for data locality only 567c8b0795cSMark F. Adams */ 568fd1112cbSBarry Smith PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc,Mat *a_Gmat,PetscCoarsenData **a_llist_parent) 569c8b0795cSMark F. Adams { 570c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,kk,Ii,ncols; 5710cbbd2e1SMark F. Adams IS perm; 572c8b0795cSMark F. Adams GAMGNode *gnodes; 573c8b0795cSMark F. Adams PetscInt *permute; 574e0940f08SMark F. Adams Mat Gmat = *a_Gmat; 5753b4367a7SBarry Smith MPI_Comm comm; 576b43b03e9SMark F. Adams MatCoarsen crs; 577c8b0795cSMark F. Adams 578c8b0795cSMark F. Adams PetscFunctionBegin; 5799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)a_pc,&comm)); 5809566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0)); 5819566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend)); 582c8b0795cSMark F. Adams nloc = (Iend-Istart); 583c8b0795cSMark F. Adams 584c8b0795cSMark F. Adams /* create random permutation with sort for geo-mg */ 5859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &gnodes)); 5869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &permute)); 587c8b0795cSMark F. Adams 588c8b0795cSMark F. Adams for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */ 5899566063dSJacob Faibussowitsch PetscCall(MatGetRow(Gmat,Ii,&ncols,NULL,NULL)); 590c8b0795cSMark F. Adams { 591c8b0795cSMark F. Adams PetscInt lid = Ii - Istart; 592c8b0795cSMark F. Adams gnodes[lid].lid = lid; 593c8b0795cSMark F. Adams gnodes[lid].degree = ncols; 594c8b0795cSMark F. Adams } 5959566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Gmat,Ii,&ncols,NULL,NULL)); 596c8b0795cSMark F. Adams } 597c8b0795cSMark F. Adams if (PETSC_TRUE) { 598e2c00dcbSBarry Smith PetscRandom rand; 599c8b0795cSMark F. Adams PetscBool *bIndexSet; 600e2c00dcbSBarry Smith PetscReal rr; 601e2c00dcbSBarry Smith PetscInt iSwapIndex; 602e2c00dcbSBarry Smith 6039566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm,&rand)); 6049566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nloc, &bIndexSet)); 6052fa5cd67SKarl Rupp for (Ii = 0; Ii < nloc; Ii++) { 6069566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand,&rr)); 607e2c00dcbSBarry Smith iSwapIndex = (PetscInt) (rr*nloc); 6082fa5cd67SKarl Rupp if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 609c8b0795cSMark F. Adams GAMGNode iTemp = gnodes[iSwapIndex]; 610c8b0795cSMark F. Adams gnodes[iSwapIndex] = gnodes[Ii]; 611c8b0795cSMark F. Adams gnodes[Ii] = iTemp; 612c8b0795cSMark F. Adams bIndexSet[Ii] = PETSC_TRUE; 613c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 614c8b0795cSMark F. Adams } 615c8b0795cSMark F. Adams } 6169566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 6179566063dSJacob Faibussowitsch PetscCall(PetscFree(bIndexSet)); 618c8b0795cSMark F. Adams } 619c8b0795cSMark F. Adams /* only sort locals */ 6200cbbd2e1SMark F. Adams qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare); 621c8b0795cSMark F. Adams /* create IS of permutation */ 6222fa5cd67SKarl Rupp for (kk=0; kk<nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */ 6239566063dSJacob Faibussowitsch PetscCall(PetscFree(gnodes)); 6249566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm)); 625c8b0795cSMark F. Adams 626c8b0795cSMark F. Adams /* get MIS aggs */ 627b43b03e9SMark F. Adams 6289566063dSJacob Faibussowitsch PetscCall(MatCoarsenCreate(comm, &crs)); 6299566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetType(crs, MATCOARSENMIS)); 6309566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetGreedyOrdering(crs, perm)); 6319566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetAdjacency(crs, Gmat)); 6329566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_FALSE)); 6339566063dSJacob Faibussowitsch PetscCall(MatCoarsenApply(crs)); 6349566063dSJacob Faibussowitsch PetscCall(MatCoarsenGetData(crs, a_llist_parent)); 6359566063dSJacob Faibussowitsch PetscCall(MatCoarsenDestroy(&crs)); 636c8b0795cSMark F. Adams 6379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 6389566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0)); 639c8b0795cSMark F. Adams PetscFunctionReturn(0); 640c8b0795cSMark F. Adams } 641c8b0795cSMark F. Adams 642c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 643c8b0795cSMark F. Adams /* 6440cbbd2e1SMark F. Adams PCGAMGProlongator_GEO 645c8b0795cSMark F. Adams 646c8b0795cSMark F. Adams Input Parameter: 647c8b0795cSMark F. Adams . pc - this 648c8b0795cSMark F. Adams . Amat - matrix on this fine level 649c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 6500cbbd2e1SMark F. Adams . selected_1 - [nselected] 6510cbbd2e1SMark F. Adams . agg_lists - [nselected] 652c8b0795cSMark F. Adams Output Parameter: 653c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 654c8b0795cSMark F. Adams */ 6552adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_GEO(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 6562e68589bSMark F. Adams { 6572e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 6582e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 659c8b0795cSMark F. Adams const PetscInt dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols; 660b43b03e9SMark F. Adams PetscInt Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid; 661c8b0795cSMark F. Adams Mat Prol; 662c5df96a5SBarry Smith PetscMPIInt rank, size; 6633b4367a7SBarry Smith MPI_Comm comm; 6640cbbd2e1SMark F. Adams IS selected_2,selected_1; 6652e68589bSMark F. Adams const PetscInt *selected_idx; 666d9558ea9SBarry Smith MatType mtype; 6672e68589bSMark F. Adams 6682e68589bSMark F. Adams PetscFunctionBegin; 6699566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm)); 6709566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0)); 6719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 6729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 6739566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend)); 6749566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat, &bs)); 67571959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 6762c71b3e2SJacob Faibussowitsch PetscCheckFalse((Iend-Istart) % bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) % bs %D",Iend,Istart,bs); 6772e68589bSMark F. Adams 6782e68589bSMark F. Adams /* get 'nLocalSelected' */ 6799566063dSJacob Faibussowitsch PetscCall(PetscCDGetMIS(agg_lists, &selected_1)); 6809566063dSJacob Faibussowitsch PetscCall(ISGetSize(selected_1, &jj)); 6819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jj, &clid_flid)); 6829566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected_1, &selected_idx)); 683c8b0795cSMark F. Adams for (kk=0,nLocalSelected=0; kk<jj; kk++) { 6842e68589bSMark F. Adams PetscInt lid = selected_idx[kk]; 685b43b03e9SMark F. Adams if (lid<nloc) { 6869566063dSJacob Faibussowitsch PetscCall(MatGetRow(Gmat,lid+my0,&ncols,NULL,NULL)); 6872fa5cd67SKarl Rupp if (ncols>1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */ 6889566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Gmat,lid+my0,&ncols,NULL,NULL)); 689b43b03e9SMark F. Adams } 6902e68589bSMark F. Adams } 6919566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected_1, &selected_idx)); 6929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected_1)); /* this is selected_1 in serial */ 6932e68589bSMark F. Adams 694d9558ea9SBarry Smith /* create prolongator matrix */ 6959566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat,&mtype)); 6969566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Prol)); 6979566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE)); 6989566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(Prol, bs, bs)); 6999566063dSJacob Faibussowitsch PetscCall(MatSetType(Prol, mtype)); 7009566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Prol,3*data_cols,NULL)); 7019566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Prol,3*data_cols,NULL,3*data_cols,NULL)); 7022e68589bSMark F. Adams 703c8b0795cSMark F. Adams /* can get all points "removed" - but not on geomg */ 7049566063dSJacob Faibussowitsch PetscCall(MatGetSize(Prol, &kk, &jj)); 7057f66b68fSMark Adams if (!jj) { 7069566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"ERROE: no selected points on coarse grid\n")); 7079566063dSJacob Faibussowitsch PetscCall(PetscFree(clid_flid)); 7089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 7090298fd71SBarry Smith *a_P_out = NULL; /* out */ 7102e68589bSMark F. Adams PetscFunctionReturn(0); 7112e68589bSMark F. Adams } 7122e68589bSMark F. Adams 7132e68589bSMark F. Adams { 7142e68589bSMark F. Adams PetscReal *coords; 715a2f3521dSMark F. Adams PetscInt data_stride; 7160298fd71SBarry Smith PetscInt *crsGID = NULL; 7172e68589bSMark F. Adams Mat Gmat2; 7182e68589bSMark F. Adams 7192c71b3e2SJacob Faibussowitsch PetscCheckFalse(dim != data_cols,PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim %D != data_cols %D",dim,data_cols); 7202e68589bSMark F. Adams /* grow ghost data for better coarse grid cover of fine grid */ 7219566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0)); 722a2f3521dSMark F. Adams /* messy method, squares graph and gets some data */ 7239566063dSJacob Faibussowitsch PetscCall(getGIDsOnSquareGraph(pc, nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID)); 7242e68589bSMark F. Adams /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */ 7259566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0)); 7262e68589bSMark F. Adams /* create global vector of coorindates in 'coords' */ 727c5df96a5SBarry Smith if (size > 1) { 7289566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords)); 729806fa848SBarry Smith } else { 730c8b0795cSMark F. Adams coords = (PetscReal*)pc_gamg->data; 731a2f3521dSMark F. Adams data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols; 7322e68589bSMark F. Adams } 7339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Gmat2)); 7342e68589bSMark F. Adams 7352e68589bSMark F. Adams /* triangulate */ 7362e68589bSMark F. Adams if (dim == 2) { 737c8b0795cSMark F. Adams PetscReal metric,tm; 7389566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0)); 7399566063dSJacob Faibussowitsch PetscCall(triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric)); 7409566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0)); 7419566063dSJacob Faibussowitsch PetscCall(PetscFree(crsGID)); 7422e68589bSMark F. Adams 7432e68589bSMark F. Adams /* clean up and create coordinates for coarse grid (output) */ 7449566063dSJacob Faibussowitsch if (size > 1) PetscCall(PetscFree(coords)); 7452e68589bSMark F. Adams 7469566063dSJacob Faibussowitsch PetscCallMPI(MPIU_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm)); 747c8b0795cSMark F. Adams if (tm > 1.) { /* needs to be globalized - should not happen */ 7489566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc," failed metric for coarse grid %e\n",(double)tm)); 7499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 750806fa848SBarry Smith } else if (metric > .0) { 7519566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"worst metric for coarse grid = %e\n",(double)metric)); 7522e68589bSMark F. Adams } 7533b4367a7SBarry Smith } else SETERRQ(comm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG"); 7542e68589bSMark F. Adams { /* create next coords - output */ 7552e68589bSMark F. Adams PetscReal *crs_crds; 7569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim*nLocalSelected, &crs_crds)); 7572e68589bSMark F. Adams for (kk=0; kk<nLocalSelected; kk++) { /* grab local select nodes to promote - output */ 758b43b03e9SMark F. Adams PetscInt lid = clid_flid[kk]; 759c8b0795cSMark F. Adams for (jj=0; jj<dim; jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid]; 7602e68589bSMark F. Adams } 761c8b0795cSMark F. Adams 7629566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 763c8b0795cSMark F. Adams pc_gamg->data = crs_crds; /* out */ 764c8b0795cSMark F. Adams pc_gamg->data_sz = dim*nLocalSelected; 7652e68589bSMark F. Adams } 7669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected_2)); 7672e68589bSMark F. Adams } 768a2f3521dSMark F. Adams 7692e68589bSMark F. Adams *a_P_out = Prol; /* out */ 7709566063dSJacob Faibussowitsch PetscCall(PetscFree(clid_flid)); 7719566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0)); 772c8b0795cSMark F. Adams PetscFunctionReturn(0); 773c8b0795cSMark F. Adams } 774c8b0795cSMark F. Adams 7759b8ffb57SJed Brown static PetscErrorCode PCDestroy_GAMG_GEO(PC pc) 7769b8ffb57SJed Brown { 7779b8ffb57SJed Brown PetscFunctionBegin; 7789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL)); 7799b8ffb57SJed Brown PetscFunctionReturn(0); 7809b8ffb57SJed Brown } 7819b8ffb57SJed Brown 782c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 783c8b0795cSMark F. Adams /* 784c8b0795cSMark F. Adams PCCreateGAMG_GEO 785c8b0795cSMark F. Adams 786c8b0795cSMark F. Adams Input Parameter: 787c8b0795cSMark F. Adams . pc - 788c8b0795cSMark F. Adams */ 789c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_GEO(PC pc) 790c8b0795cSMark F. Adams { 791c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 792c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 793c8b0795cSMark F. Adams 794c8b0795cSMark F. Adams PetscFunctionBegin; 7951ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO; 7969b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_GEO; 797c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 798c8b0795cSMark F. Adams 799c8b0795cSMark F. Adams /* set internal function pointers */ 800fd1112cbSBarry Smith pc_gamg->ops->graph = PCGAMGGraph_GEO; 801fd1112cbSBarry Smith pc_gamg->ops->coarsen = PCGAMGCoarsen_GEO; 8021ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_GEO; 803fd1112cbSBarry Smith pc_gamg->ops->optprolongator = NULL; 8041ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_GEO; 805c8b0795cSMark F. Adams 8069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_GEO)); 8072e68589bSMark F. Adams PetscFunctionReturn(0); 8082e68589bSMark F. Adams } 809