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) 8*5a449d36SBarry Smith #if !defined(ANSI_DECLARATORS) 9*5a449d36SBarry Smith #define ANSI_DECLARATORS 10*5a449d36SBarry 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 22b817416eSBarry Smith PETSC_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; 382e68589bSMark F. Adams PetscErrorCode ierr; 3990fbc344SStefano Zampini PetscInt arrsz,bs,my0,kk,ii,nloc,Iend,aloc; 402e68589bSMark F. Adams Mat Amat = pc->pmat; 412e68589bSMark F. Adams 422e68589bSMark F. Adams PetscFunctionBegin; 432e68589bSMark F. Adams PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1); 442e68589bSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 452e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr); 4690fbc344SStefano Zampini aloc = (Iend-my0); 472e68589bSMark F. Adams nloc = (Iend-my0)/bs; 48a2f3521dSMark F. Adams 4990fbc344SStefano Zampini if (nloc!=a_nloc && aloc!=a_nloc) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of local blocks %D must be %D or %D.",a_nloc,nloc,aloc); 502e68589bSMark F. Adams 51c8b0795cSMark F. Adams pc_gamg->data_cell_rows = 1; 5290fbc344SStefano Zampini if (!coords && nloc > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'."); 53c8b0795cSMark F. Adams pc_gamg->data_cell_cols = ndm; /* coordinates */ 542e68589bSMark F. Adams 55c8b0795cSMark F. Adams arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 562e68589bSMark F. Adams 572e68589bSMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 587f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 592e68589bSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 60854ce69bSBarry Smith ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr); 612e68589bSMark F. Adams } 622e68589bSMark F. Adams for (kk=0; kk<arrsz; kk++) pc_gamg->data[kk] = -999.; 632e68589bSMark F. Adams pc_gamg->data[arrsz] = -99.; 642e68589bSMark F. Adams /* copy data in - column oriented */ 6590fbc344SStefano Zampini if (nloc == a_nloc) { 662e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) { 672e68589bSMark F. Adams for (ii = 0; ii < ndm; ii++) { 682e68589bSMark F. Adams pc_gamg->data[ii*nloc + kk] = coords[kk*ndm + ii]; 692e68589bSMark F. Adams } 702e68589bSMark F. Adams } 7190fbc344SStefano Zampini } else { /* assumes the coordinates are blocked */ 7290fbc344SStefano Zampini for (kk = 0; kk < nloc; kk++) { 7390fbc344SStefano Zampini for (ii = 0; ii < ndm; ii++) { 7490fbc344SStefano Zampini pc_gamg->data[ii*nloc + kk] = coords[bs*kk*ndm + ii]; 7590fbc344SStefano Zampini } 7690fbc344SStefano Zampini } 7790fbc344SStefano Zampini } 7871959b99SBarry Smith if (pc_gamg->data[arrsz] != -99.) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data[arrsz %D] %g != -99.",arrsz,pc_gamg->data[arrsz]); 792e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 802e68589bSMark F. Adams PetscFunctionReturn(0); 812e68589bSMark F. Adams } 822e68589bSMark F. Adams 832e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 842e68589bSMark F. Adams /* 852e68589bSMark F. Adams PCSetData_GEO 862e68589bSMark F. Adams 872e68589bSMark F. Adams Input Parameter: 882e68589bSMark F. Adams . pc - 892e68589bSMark F. Adams */ 90b8cd405aSMark F. Adams PetscErrorCode PCSetData_GEO(PC pc, Mat m) 912e68589bSMark F. Adams { 922e68589bSMark F. Adams PetscFunctionBegin; 93ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"GEO MG needs coordinates"); 942e68589bSMark F. Adams } 952e68589bSMark F. Adams 962e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 972e68589bSMark F. Adams /* 982e68589bSMark F. Adams PCSetFromOptions_GEO 992e68589bSMark F. Adams 1002e68589bSMark F. Adams Input Parameter: 1012e68589bSMark F. Adams . pc - 1022e68589bSMark F. Adams */ 1034416b707SBarry Smith PetscErrorCode PCSetFromOptions_GEO(PetscOptionItems *PetscOptionsObject,PC pc) 1042e68589bSMark F. Adams { 1052e68589bSMark F. Adams PetscErrorCode ierr; 1062e68589bSMark F. Adams 1072e68589bSMark F. Adams PetscFunctionBegin; 108e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-GEO options");CHKERRQ(ierr); 1092e68589bSMark F. Adams { 1102e68589bSMark F. Adams /* -pc_gamg_sa_nsmooths */ 1112e68589bSMark F. Adams /* pc_gamg_sa->smooths = 0; */ 1122e68589bSMark F. Adams /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */ 1132e68589bSMark F. Adams /* "smoothing steps for smoothed aggregation, usually 1 (0)", */ 1142e68589bSMark F. Adams /* "PCGAMGSetNSmooths_AGG", */ 1152e68589bSMark F. Adams /* pc_gamg_sa->smooths, */ 1162e68589bSMark F. Adams /* &pc_gamg_sa->smooths, */ 1172e68589bSMark F. Adams /* &flag); */ 1182e68589bSMark F. Adams /* CHKERRQ(ierr); */ 1192e68589bSMark F. Adams } 1202e68589bSMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1212e68589bSMark F. Adams PetscFunctionReturn(0); 1222e68589bSMark F. Adams } 1232e68589bSMark F. Adams 1242e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1252e68589bSMark F. Adams /* 1262e68589bSMark F. Adams triangulateAndFormProl 1272e68589bSMark F. Adams 1282e68589bSMark F. Adams Input Parameter: 1292e68589bSMark F. Adams . selected_2 - list of selected local ID, includes selected ghosts 130a2f3521dSMark F. Adams . data_stride - 131a2f3521dSMark F. Adams . coords[2*data_stride] - column vector of local coordinates w/ ghosts 1322adfe9d3SBarry Smith . nselected_1 - selected IDs that go with base (1) graph includes selected ghosts 1330cbbd2e1SMark F. Adams . clid_lid_1[nselected_1] - lids of selected (c) nodes ??????????? 1342adfe9d3SBarry Smith . agg_lists_1 - list of aggregates selected_1 vertices of aggregate unselected vertices 1352e68589bSMark F. Adams . crsGID[selected.size()] - global index for prolongation operator 1362e68589bSMark F. Adams . bs - block size 1372e68589bSMark F. Adams Output Parameter: 1382e68589bSMark F. Adams . a_Prol - prolongation operator 1392e68589bSMark F. Adams . a_worst_best - measure of worst missed fine vertex, 0 is no misses 1402e68589bSMark F. Adams */ 1412adfe9d3SBarry 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, 1422adfe9d3SBarry Smith const PetscInt crsGID[],PetscInt bs,Mat a_Prol,PetscReal *a_worst_best) 1432e68589bSMark F. Adams { 1442e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE) 1452e68589bSMark F. Adams PetscErrorCode ierr; 146b43b03e9SMark F. Adams PetscInt jj,tid,tt,idx,nselected_2; 1472e68589bSMark F. Adams struct triangulateio in,mid; 1480cbbd2e1SMark F. Adams const PetscInt *selected_idx_2; 14973911c69SBarry Smith PetscMPIInt rank; 1502e68589bSMark F. Adams PetscInt Istart,Iend,nFineLoc,myFine0; 1512e68589bSMark F. Adams int kk,nPlotPts,sid; 1523b4367a7SBarry Smith MPI_Comm comm; 153c8b0795cSMark F. Adams PetscReal tm; 154c8b0795cSMark F. Adams 1556e111a19SKarl Rupp PetscFunctionBegin; 1563b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr); 1573b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 158c8b0795cSMark F. Adams ierr = ISGetSize(selected_2, &nselected_2);CHKERRQ(ierr); 1592e68589bSMark F. Adams if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */ 1602e68589bSMark F. Adams *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */ 161806fa848SBarry Smith } else *a_worst_best = 0.0; 162b2566f29SBarry Smith ierr = MPIU_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm);CHKERRQ(ierr); 163c8b0795cSMark F. Adams if (tm > 0.0) { 164c8b0795cSMark F. Adams *a_worst_best = 100.0; 1652e68589bSMark F. Adams PetscFunctionReturn(0); 1662e68589bSMark F. Adams } 1672e68589bSMark F. Adams ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr); 1682e68589bSMark F. Adams nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs; 1692e68589bSMark F. Adams nPlotPts = nFineLoc; /* locals */ 1702e68589bSMark F. Adams /* traingle */ 1712e68589bSMark F. Adams /* Define input points - in*/ 1722e68589bSMark F. Adams in.numberofpoints = nselected_2; 1732e68589bSMark F. Adams in.numberofpointattributes = 0; 1742e68589bSMark F. Adams /* get nselected points */ 175e632b94dSBarry Smith ierr = PetscMalloc1(2*nselected_2, &in.pointlist);CHKERRQ(ierr); 1762e68589bSMark F. Adams ierr = ISGetIndices(selected_2, &selected_idx_2);CHKERRQ(ierr); 1772e68589bSMark F. Adams 1782e68589bSMark F. Adams for (kk=0,sid=0; kk<nselected_2; kk++,sid += 2) { 1792e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 1802e68589bSMark F. Adams in.pointlist[sid] = coords[lid]; 181a2f3521dSMark F. Adams in.pointlist[sid+1] = coords[data_stride + lid]; 1822e68589bSMark F. Adams if (lid>=nFineLoc) nPlotPts++; 1832e68589bSMark F. Adams } 18471959b99SBarry Smith if (sid != 2*nselected_2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != 2*nselected_2 %D",sid,nselected_2); 1852e68589bSMark F. Adams 1862e68589bSMark F. Adams in.numberofsegments = 0; 1872e68589bSMark F. Adams in.numberofedges = 0; 1882e68589bSMark F. Adams in.numberofholes = 0; 1892e68589bSMark F. Adams in.numberofregions = 0; 1902e68589bSMark F. Adams in.trianglelist = 0; 1912e68589bSMark F. Adams in.segmentmarkerlist = 0; 1922e68589bSMark F. Adams in.pointattributelist = 0; 1932e68589bSMark F. Adams in.pointmarkerlist = 0; 1942e68589bSMark F. Adams in.triangleattributelist = 0; 1952e68589bSMark F. Adams in.trianglearealist = 0; 1962e68589bSMark F. Adams in.segmentlist = 0; 1972e68589bSMark F. Adams in.holelist = 0; 1982e68589bSMark F. Adams in.regionlist = 0; 1992e68589bSMark F. Adams in.edgelist = 0; 2002e68589bSMark F. Adams in.edgemarkerlist = 0; 2012e68589bSMark F. Adams in.normlist = 0; 2022fa5cd67SKarl Rupp 2032e68589bSMark F. Adams /* triangulate */ 2042e68589bSMark F. Adams mid.pointlist = 0; /* Not needed if -N switch used. */ 2052e68589bSMark F. Adams /* Not needed if -N switch used or number of point attributes is zero: */ 2062e68589bSMark F. Adams mid.pointattributelist = 0; 2072e68589bSMark F. Adams mid.pointmarkerlist = 0; /* Not needed if -N or -B switch used. */ 2082e68589bSMark F. Adams mid.trianglelist = 0; /* Not needed if -E switch used. */ 2092e68589bSMark F. Adams /* Not needed if -E switch used or number of triangle attributes is zero: */ 2102e68589bSMark F. Adams mid.triangleattributelist = 0; 2112e68589bSMark F. Adams mid.neighborlist = 0; /* Needed only if -n switch used. */ 2122e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P not used: */ 2132e68589bSMark F. Adams mid.segmentlist = 0; 2142e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P and -B not used: */ 2152e68589bSMark F. Adams mid.segmentmarkerlist = 0; 2162e68589bSMark F. Adams mid.edgelist = 0; /* Needed only if -e switch used. */ 2172e68589bSMark F. Adams mid.edgemarkerlist = 0; /* Needed if -e used and -B not used. */ 2182e68589bSMark F. Adams mid.numberoftriangles = 0; 2192e68589bSMark F. Adams 2202e68589bSMark F. Adams /* Triangulate the points. Switches are chosen to read and write a */ 2212e68589bSMark F. Adams /* PSLG (p), preserve the convex hull (c), number everything from */ 2222e68589bSMark F. Adams /* zero (z), assign a regional attribute to each element (A), and */ 2232e68589bSMark F. Adams /* produce an edge list (e), a Voronoi diagram (v), and a triangle */ 2242e68589bSMark F. Adams /* neighbor list (n). */ 2252e68589bSMark F. Adams if (nselected_2 != 0) { /* inactive processor */ 2262e68589bSMark F. Adams char args[] = "npczQ"; /* c is needed ? */ 2272e68589bSMark F. Adams triangulate(args, &in, &mid, (struct triangulateio*) NULL); 2282e68589bSMark F. Adams /* output .poly files for 'showme' */ 2292e68589bSMark F. Adams if (!PETSC_TRUE) { 2302e68589bSMark F. Adams static int level = 1; 2312e68589bSMark F. Adams FILE *file; char fname[32]; 2322e68589bSMark F. Adams 233c5df96a5SBarry Smith sprintf(fname,"C%d_%d.poly",level,rank); file = fopen(fname, "w"); 2342e68589bSMark F. Adams /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/ 2352e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); 2362e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2372e68589bSMark F. Adams for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid += 2) { 2382e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2392e68589bSMark F. Adams } 2402e68589bSMark F. Adams /*One line: <# of segments> <# of boundary markers (0 or 1)> */ 2412e68589bSMark F. Adams fprintf(file, "%d %d\n",0,0); 2422e68589bSMark F. Adams /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */ 2432e68589bSMark F. Adams /* One line: <# of holes> */ 2442e68589bSMark F. Adams fprintf(file, "%d\n",0); 2452e68589bSMark F. Adams /* Following lines: <hole #> <x> <y> */ 2462e68589bSMark F. Adams /* Optional line: <# of regional attributes and/or area constraints> */ 2472e68589bSMark F. Adams /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */ 2482e68589bSMark F. Adams fclose(file); 2492e68589bSMark F. Adams 2502e68589bSMark F. Adams /* elems */ 251c5df96a5SBarry Smith sprintf(fname,"C%d_%d.ele",level,rank); file = fopen(fname, "w"); 2522e68589bSMark F. Adams /* First line: <# of triangles> <nodes per triangle> <# of attributes> */ 2532e68589bSMark F. Adams fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0); 2542e68589bSMark F. Adams /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */ 2552e68589bSMark F. Adams for (kk=0,sid=0; kk<mid.numberoftriangles; kk++,sid += 3) { 2562e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]); 2572e68589bSMark F. Adams } 2582e68589bSMark F. Adams fclose(file); 2592e68589bSMark F. Adams 260c5df96a5SBarry Smith sprintf(fname,"C%d_%d.node",level,rank); file = fopen(fname, "w"); 2612e68589bSMark F. Adams /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */ 2622e68589bSMark F. Adams /* fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); */ 2632e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",nPlotPts,2,0,0); 2642e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2652e68589bSMark F. Adams for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid+=2) { 2662e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2672e68589bSMark F. Adams } 2682e68589bSMark F. Adams 2692e68589bSMark F. Adams sid /= 2; 2702e68589bSMark F. Adams for (jj=0; jj<nFineLoc; jj++) { 2712e68589bSMark F. Adams PetscBool sel = PETSC_TRUE; 2722e68589bSMark F. Adams for (kk=0; kk<nselected_2 && sel; kk++) { 2732e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 2742e68589bSMark F. Adams if (lid == jj) sel = PETSC_FALSE; 2752e68589bSMark F. Adams } 2762fa5cd67SKarl Rupp if (sel) fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]); 2772e68589bSMark F. Adams } 2782e68589bSMark F. Adams fclose(file); 27971959b99SBarry Smith if (sid != nPlotPts) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != nPlotPts %D",sid,nPlotPts); 2802e68589bSMark F. Adams level++; 2812e68589bSMark F. Adams } 2822e68589bSMark F. Adams } 2830cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 2840cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 2852e68589bSMark F. Adams #endif 2862e68589bSMark F. Adams { /* form P - setup some maps */ 2870cbbd2e1SMark F. Adams PetscInt clid,mm,*nTri,*node_tri; 2882e68589bSMark F. Adams 289e632b94dSBarry Smith ierr = PetscMalloc2(nselected_2, &node_tri,nselected_2, &nTri);CHKERRQ(ierr); 2902e68589bSMark F. Adams 2912e68589bSMark F. Adams /* need list of triangles on node */ 2922e68589bSMark F. Adams for (kk=0; kk<nselected_2; kk++) nTri[kk] = 0; 2932e68589bSMark F. Adams for (tid=0,kk=0; tid<mid.numberoftriangles; tid++) { 2942e68589bSMark F. Adams for (jj=0; jj<3; jj++) { 2952e68589bSMark F. Adams PetscInt cid = mid.trianglelist[kk++]; 2962e68589bSMark F. Adams if (nTri[cid] == 0) node_tri[cid] = tid; 2972e68589bSMark F. Adams nTri[cid]++; 2982e68589bSMark F. Adams } 2992e68589bSMark F. Adams } 3002e68589bSMark F. Adams #define EPS 1.e-12 3012e68589bSMark F. Adams /* find points and set prolongation */ 3020cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nFineLoc; mm++) { 303e78576d6SMark F. Adams PetscBool ise; 304e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise);CHKERRQ(ierr); 305e78576d6SMark F. Adams if (!ise) { 3060cbbd2e1SMark F. Adams const PetscInt lid = mm; 3072e68589bSMark F. Adams PetscScalar AA[3][3]; 3082e68589bSMark F. Adams PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO; 309539c167fSBarry Smith PetscCDIntNd *pos; 310539c167fSBarry Smith 311e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos);CHKERRQ(ierr); 312e78576d6SMark F. Adams while (pos) { 313e78576d6SMark F. Adams PetscInt flid; 314484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &flid);CHKERRQ(ierr); 315e78576d6SMark F. Adams ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos);CHKERRQ(ierr); 3160cbbd2e1SMark F. Adams 3172e68589bSMark F. Adams if (flid < nFineLoc) { /* could be a ghost */ 3182e68589bSMark F. Adams PetscInt bestTID = -1; PetscReal best_alpha = 1.e10; 3192e68589bSMark F. Adams const PetscInt fgid = flid + myFine0; 3202e68589bSMark F. Adams /* compute shape function for gid */ 321a2f3521dSMark F. Adams const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0}; 3222e68589bSMark F. Adams PetscBool haveit =PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3]; 3232fa5cd67SKarl Rupp 3242e68589bSMark F. Adams /* look for it */ 3250cbbd2e1SMark F. Adams for (tid = node_tri[clid], jj=0; 3262e68589bSMark F. Adams jj < 5 && !haveit && tid != -1; 3272e68589bSMark F. Adams jj++) { 3282e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3292e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3302e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 331a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3322e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3332e68589bSMark F. Adams } 3342e68589bSMark F. Adams 3352e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = (PetscScalar)fcoord[tt]; 3362e68589bSMark F. Adams 3372e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 3388b83055fSJed Brown PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3392e68589bSMark F. Adams { 3402e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal lowest=1.e10; 3412e68589bSMark F. Adams for (tt = 0, idx = 0; tt < 3; tt++) { 3422e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE; 3432e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) < lowest) { 3442e68589bSMark F. Adams lowest = PetscRealPart(alpha[tt]); 3452e68589bSMark F. Adams idx = tt; 3462e68589bSMark F. Adams } 3472e68589bSMark F. Adams } 3482e68589bSMark F. Adams haveit = have; 3492e68589bSMark F. Adams } 3502e68589bSMark F. Adams tid = mid.neighborlist[3*tid + idx]; 3512e68589bSMark F. Adams } 3522e68589bSMark F. Adams 3532e68589bSMark F. Adams if (!haveit) { 3542e68589bSMark F. Adams /* brute force */ 3552e68589bSMark F. Adams for (tid=0; tid<mid.numberoftriangles && !haveit; tid++) { 3562e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3572e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3582e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 359a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3602e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3612e68589bSMark F. Adams } 3622e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt]; 3632e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 3648b83055fSJed Brown PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3652e68589bSMark F. Adams { 3662e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal worst=0.0, v; 3672e68589bSMark F. Adams for (tt=0; tt<3 && have; tt++) { 3682e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE; 3692e68589bSMark F. Adams if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v; 3702e68589bSMark F. Adams } 3712e68589bSMark F. Adams if (worst < best_alpha) { 3722e68589bSMark F. Adams best_alpha = worst; bestTID = tid; 3732e68589bSMark F. Adams } 3742e68589bSMark F. Adams haveit = have; 3752e68589bSMark F. Adams } 3762e68589bSMark F. Adams } 3772e68589bSMark F. Adams } 3782e68589bSMark F. Adams if (!haveit) { 3792e68589bSMark F. Adams if (best_alpha > *a_worst_best) *a_worst_best = best_alpha; 3802e68589bSMark F. Adams /* use best one */ 3812e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3822e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*bestTID + tt]; 3832e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 384a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3852e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3862e68589bSMark F. Adams } 3872e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt]; 3882e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 3898b83055fSJed Brown PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3902e68589bSMark F. Adams } 3912e68589bSMark F. Adams 3922e68589bSMark F. Adams /* put in row of P */ 3932e68589bSMark F. Adams for (idx=0; idx<3; idx++) { 3942e68589bSMark F. Adams PetscScalar shp = alpha[idx]; 3952e68589bSMark F. Adams if (PetscAbs(PetscRealPart(shp)) > 1.e-6) { 3962e68589bSMark F. Adams PetscInt cgid = crsGID[clids[idx]]; 3972e68589bSMark F. Adams PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */ 3982e68589bSMark F. Adams for (tt=0; tt < bs; tt++, ii++, jj++) { 3992e68589bSMark F. Adams ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES);CHKERRQ(ierr); 4002e68589bSMark F. Adams } 4012e68589bSMark F. Adams } 4022e68589bSMark F. Adams } 4032e68589bSMark F. Adams } 4040cbbd2e1SMark F. Adams } /* aggregates iterations */ 4050cbbd2e1SMark F. Adams clid++; 4060cbbd2e1SMark F. Adams } /* a coarse agg */ 4070cbbd2e1SMark F. Adams } /* for all fine nodes */ 4080cbbd2e1SMark F. Adams 4092e68589bSMark F. Adams ierr = ISRestoreIndices(selected_2, &selected_idx_2);CHKERRQ(ierr); 4102e68589bSMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4112e68589bSMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4122e68589bSMark F. Adams 413e632b94dSBarry Smith ierr = PetscFree2(node_tri,nTri);CHKERRQ(ierr); 4142e68589bSMark F. Adams } 4150cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4160cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 4172e68589bSMark F. Adams #endif 4182e68589bSMark F. Adams free(mid.trianglelist); 4192e68589bSMark F. Adams free(mid.neighborlist); 4202e68589bSMark F. Adams ierr = PetscFree(in.pointlist);CHKERRQ(ierr); 4212e68589bSMark F. Adams PetscFunctionReturn(0); 4222e68589bSMark F. Adams #else 4233b4367a7SBarry Smith SETERRQ(PetscObjectComm((PetscObject)a_Prol),PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG"); 4242e68589bSMark F. Adams #endif 4252e68589bSMark F. Adams } 4262e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 4272e68589bSMark F. Adams /* 4282e68589bSMark F. Adams getGIDsOnSquareGraph - square graph, get 4292e68589bSMark F. Adams 4302e68589bSMark F. Adams Input Parameter: 4310cbbd2e1SMark F. Adams . nselected_1 - selected local indices (includes ghosts in input Gmat1) 432b43b03e9SMark F. Adams . clid_lid_1 - [nselected_1] lids of selected nodes 4332e68589bSMark F. Adams . Gmat1 - graph that goes with 'selected_1' 4342e68589bSMark F. Adams Output Parameter: 4352e68589bSMark F. Adams . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2) 4362e68589bSMark F. Adams . a_Gmat_2 - graph that is squared of 'Gmat_1' 4372e68589bSMark F. Adams . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes 4382e68589bSMark F. Adams */ 4392adfe9d3SBarry Smith static PetscErrorCode getGIDsOnSquareGraph(PetscInt nselected_1,const PetscInt clid_lid_1[],const Mat Gmat1,IS *a_selected_2,Mat *a_Gmat_2,PetscInt **a_crsGID) 4402e68589bSMark F. Adams { 4412e68589bSMark F. Adams PetscErrorCode ierr; 44273911c69SBarry Smith PetscMPIInt size; 443b43b03e9SMark F. Adams PetscInt *crsGID, kk,my0,Iend,nloc; 4443b4367a7SBarry Smith MPI_Comm comm; 4452e68589bSMark F. Adams 4462e68589bSMark F. Adams PetscFunctionBegin; 4473b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr); 4483b4367a7SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 4492e68589bSMark F. Adams ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend);CHKERRQ(ierr); /* AIJ */ 4502e68589bSMark F. Adams nloc = Iend - my0; /* this does not change */ 4512e68589bSMark F. Adams 452c5df96a5SBarry Smith if (size == 1) { /* not much to do in serial */ 453785e854fSJed Brown ierr = PetscMalloc1(nselected_1, &crsGID);CHKERRQ(ierr); 454b43b03e9SMark F. Adams for (kk=0; kk<nselected_1; kk++) crsGID[kk] = kk; 4552e68589bSMark F. Adams *a_Gmat_2 = 0; 456806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);CHKERRQ(ierr); 457806fa848SBarry Smith } else { 458b43b03e9SMark F. Adams PetscInt idx,num_fine_ghosts,num_crs_ghost,myCrs0; 4592e68589bSMark F. Adams Mat_MPIAIJ *mpimat2; 4602e68589bSMark F. Adams Mat Gmat2; 4612e68589bSMark F. Adams Vec locState; 4622e68589bSMark F. Adams PetscScalar *cpcol_state; 4632e68589bSMark F. Adams 4642e68589bSMark F. Adams /* scan my coarse zero gid, set 'lid_state' with coarse GID */ 465b43b03e9SMark F. Adams kk = nselected_1; 466dbd2ff41SBarry Smith ierr = MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr); 467b43b03e9SMark F. Adams myCrs0 -= nselected_1; 4682e68589bSMark F. Adams 469b43b03e9SMark F. Adams if (a_Gmat_2) { /* output */ 4702e68589bSMark F. Adams /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */ 4712e68589bSMark F. Adams ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr); 4722e68589bSMark F. Adams *a_Gmat_2 = Gmat2; /* output */ 473806fa848SBarry Smith } else Gmat2 = Gmat1; /* use local to get crsGIDs at least */ 4742e68589bSMark F. Adams /* get coarse grid GIDS for selected (locals and ghosts) */ 4752e68589bSMark F. Adams mpimat2 = (Mat_MPIAIJ*)Gmat2->data; 4762a7a6963SBarry Smith ierr = MatCreateVecs(Gmat2, &locState, 0);CHKERRQ(ierr); 4772e68589bSMark F. Adams ierr = VecSet(locState, (PetscScalar)(PetscReal)(-1));CHKERRQ(ierr); /* set with UNKNOWN state */ 478b43b03e9SMark F. Adams for (kk=0; kk<nselected_1; kk++) { 479b43b03e9SMark F. Adams PetscInt fgid = clid_lid_1[kk] + my0; 4802e68589bSMark F. Adams PetscScalar v = (PetscScalar)(kk+myCrs0); 4812e68589bSMark F. Adams ierr = VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES);CHKERRQ(ierr); /* set with PID */ 4822e68589bSMark F. Adams } 4832e68589bSMark F. Adams ierr = VecAssemblyBegin(locState);CHKERRQ(ierr); 4842e68589bSMark F. Adams ierr = VecAssemblyEnd(locState);CHKERRQ(ierr); 4852e68589bSMark F. Adams ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4862e68589bSMark F. Adams ierr = VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4872e68589bSMark F. Adams ierr = VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts);CHKERRQ(ierr); 4882e68589bSMark F. Adams ierr = VecGetArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr); 4892e68589bSMark F. Adams for (kk=0,num_crs_ghost=0; kk<num_fine_ghosts; kk++) { 4902e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++; 4912e68589bSMark F. Adams } 492854ce69bSBarry Smith ierr = PetscMalloc1(nselected_1+num_crs_ghost, &crsGID);CHKERRQ(ierr); /* output */ 4932e68589bSMark F. Adams { 4942e68589bSMark F. Adams PetscInt *selected_set; 495854ce69bSBarry Smith ierr = PetscMalloc1(nselected_1+num_crs_ghost, &selected_set);CHKERRQ(ierr); 4962e68589bSMark F. Adams /* do ghost of 'crsGID' */ 497b43b03e9SMark F. Adams for (kk=0,idx=nselected_1; kk<num_fine_ghosts; kk++) { 4982e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) { 4992e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5002e68589bSMark F. Adams selected_set[idx] = nloc + kk; 5012e68589bSMark F. Adams crsGID[idx++] = cgid; 5022e68589bSMark F. Adams } 5032e68589bSMark F. Adams } 50471959b99SBarry Smith if (idx != (nselected_1+num_crs_ghost)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != (nselected_1 %D + num_crs_ghost %D)",idx,nselected_1,num_crs_ghost); 5052e68589bSMark F. Adams ierr = VecRestoreArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr); 5062e68589bSMark F. Adams /* do locals in 'crsGID' */ 5072e68589bSMark F. Adams ierr = VecGetArray(locState, &cpcol_state);CHKERRQ(ierr); 5082e68589bSMark F. Adams for (kk=0,idx=0; kk<nloc; kk++) { 5092e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) { 5102e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5112e68589bSMark F. Adams selected_set[idx] = kk; 5122e68589bSMark F. Adams crsGID[idx++] = cgid; 5132e68589bSMark F. Adams } 5142e68589bSMark F. Adams } 51571959b99SBarry Smith if (idx != nselected_1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != nselected_1 %D",idx,nselected_1); 5162e68589bSMark F. Adams ierr = VecRestoreArray(locState, &cpcol_state);CHKERRQ(ierr); 5172e68589bSMark F. Adams 5182e68589bSMark F. Adams if (a_selected_2 != 0) { /* output */ 519806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);CHKERRQ(ierr); 520806fa848SBarry Smith } else { 5212e68589bSMark F. Adams ierr = PetscFree(selected_set);CHKERRQ(ierr); 5222e68589bSMark F. Adams } 5230cbbd2e1SMark F. Adams } 5242e68589bSMark F. Adams ierr = VecDestroy(&locState);CHKERRQ(ierr); 5252e68589bSMark F. Adams } 5262e68589bSMark F. Adams *a_crsGID = crsGID; /* output */ 5272e68589bSMark F. Adams PetscFunctionReturn(0); 5282e68589bSMark F. Adams } 5292e68589bSMark F. Adams 5302e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 5312e68589bSMark F. Adams /* 532fd1112cbSBarry Smith PCGAMGGraph_GEO 5332e68589bSMark F. Adams 5342e68589bSMark F. Adams Input Parameter: 5352e68589bSMark F. Adams . pc - this 5362e68589bSMark F. Adams . Amat - matrix on this fine level 5372e68589bSMark F. Adams Output Parameter: 538c8b0795cSMark F. Adams . a_Gmat 5392e68589bSMark F. Adams */ 5402adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_GEO(PC pc,Mat Amat,Mat *a_Gmat) 541c8b0795cSMark F. Adams { 542c8b0795cSMark F. Adams PetscErrorCode ierr; 543c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 544c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 545c1eae691SMark Adams const PetscReal vfilter = pc_gamg->threshold[0]; 5463b4367a7SBarry Smith MPI_Comm comm; 547c8b0795cSMark F. Adams Mat Gmat; 5480cbbd2e1SMark F. Adams PetscBool set,flg,symm; 5496e111a19SKarl Rupp 550c8b0795cSMark F. Adams PetscFunctionBegin; 5513b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 552fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGGraph_GEO,0,0,0,0);CHKERRQ(ierr); 553c8b0795cSMark F. Adams 5540cbbd2e1SMark F. Adams ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); 555263489e9SJed Brown symm = (PetscBool)!(set && flg); 5560cbbd2e1SMark F. Adams 5572d7fac45SMark F. Adams ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr); 558bf4339c2SBarry Smith ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr); 559c8b0795cSMark F. Adams 560c8b0795cSMark F. Adams *a_Gmat = Gmat; 561fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGGraph_GEO,0,0,0,0);CHKERRQ(ierr); 562c8b0795cSMark F. Adams PetscFunctionReturn(0); 563c8b0795cSMark F. Adams } 564c8b0795cSMark F. Adams 565c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 566c8b0795cSMark F. Adams /* 567fd1112cbSBarry Smith PCGAMGCoarsen_GEO 568c8b0795cSMark F. Adams 569c8b0795cSMark F. Adams Input Parameter: 570e0940f08SMark F. Adams . a_pc - this 571e0940f08SMark F. Adams . a_Gmat - graph 572c8b0795cSMark F. Adams Output Parameter: 573c8b0795cSMark F. Adams . a_llist_parent - linked list from selected indices for data locality only 574c8b0795cSMark F. Adams */ 575fd1112cbSBarry Smith PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc,Mat *a_Gmat,PetscCoarsenData **a_llist_parent) 576c8b0795cSMark F. Adams { 577c8b0795cSMark F. Adams PetscErrorCode ierr; 578c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,kk,Ii,ncols; 5790cbbd2e1SMark F. Adams IS perm; 580c8b0795cSMark F. Adams GAMGNode *gnodes; 581c8b0795cSMark F. Adams PetscInt *permute; 582e0940f08SMark F. Adams Mat Gmat = *a_Gmat; 5833b4367a7SBarry Smith MPI_Comm comm; 584b43b03e9SMark F. Adams MatCoarsen crs; 585c8b0795cSMark F. Adams 586c8b0795cSMark F. Adams PetscFunctionBegin; 5873b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)a_pc,&comm);CHKERRQ(ierr); 5880cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr); 589c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr); 590c8b0795cSMark F. Adams nloc = (Iend-Istart); 591c8b0795cSMark F. Adams 592c8b0795cSMark F. Adams /* create random permutation with sort for geo-mg */ 593785e854fSJed Brown ierr = PetscMalloc1(nloc, &gnodes);CHKERRQ(ierr); 594785e854fSJed Brown ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr); 595c8b0795cSMark F. Adams 596c8b0795cSMark F. Adams for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */ 597c8b0795cSMark F. Adams ierr = MatGetRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr); 598c8b0795cSMark F. Adams { 599c8b0795cSMark F. Adams PetscInt lid = Ii - Istart; 600c8b0795cSMark F. Adams gnodes[lid].lid = lid; 601c8b0795cSMark F. Adams gnodes[lid].degree = ncols; 602c8b0795cSMark F. Adams } 603c8b0795cSMark F. Adams ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr); 604c8b0795cSMark F. Adams } 605c8b0795cSMark F. Adams if (PETSC_TRUE) { 606e2c00dcbSBarry Smith PetscRandom rand; 607c8b0795cSMark F. Adams PetscBool *bIndexSet; 608e2c00dcbSBarry Smith PetscReal rr; 609e2c00dcbSBarry Smith PetscInt iSwapIndex; 610e2c00dcbSBarry Smith 611e2c00dcbSBarry Smith ierr = PetscRandomCreate(comm,&rand);CHKERRQ(ierr); 612e2c00dcbSBarry Smith ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr); 6132fa5cd67SKarl Rupp for (Ii = 0; Ii < nloc; Ii++) { 614e2c00dcbSBarry Smith ierr = PetscRandomGetValueReal(rand,&rr);CHKERRQ(ierr); 615e2c00dcbSBarry Smith iSwapIndex = (PetscInt) (rr*nloc); 6162fa5cd67SKarl Rupp if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 617c8b0795cSMark F. Adams GAMGNode iTemp = gnodes[iSwapIndex]; 618c8b0795cSMark F. Adams gnodes[iSwapIndex] = gnodes[Ii]; 619c8b0795cSMark F. Adams gnodes[Ii] = iTemp; 620c8b0795cSMark F. Adams bIndexSet[Ii] = PETSC_TRUE; 621c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 622c8b0795cSMark F. Adams } 623c8b0795cSMark F. Adams } 624e2c00dcbSBarry Smith ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 625c8b0795cSMark F. Adams ierr = PetscFree(bIndexSet);CHKERRQ(ierr); 626c8b0795cSMark F. Adams } 627c8b0795cSMark F. Adams /* only sort locals */ 6280cbbd2e1SMark F. Adams qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare); 629c8b0795cSMark F. Adams /* create IS of permutation */ 6302fa5cd67SKarl Rupp for (kk=0; kk<nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */ 631c8b0795cSMark F. Adams ierr = PetscFree(gnodes);CHKERRQ(ierr); 632b817416eSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr); 633c8b0795cSMark F. Adams 634c8b0795cSMark F. Adams /* get MIS aggs */ 635b43b03e9SMark F. Adams 6363b4367a7SBarry Smith ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr); 637b43b03e9SMark F. Adams ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); 638b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr); 639b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr); 640b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs(crs, PETSC_FALSE);CHKERRQ(ierr); 641b43b03e9SMark F. Adams ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 6420cbbd2e1SMark F. Adams ierr = MatCoarsenGetData(crs, a_llist_parent);CHKERRQ(ierr); 643b43b03e9SMark F. Adams ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 644c8b0795cSMark F. Adams 645c8b0795cSMark F. Adams ierr = ISDestroy(&perm);CHKERRQ(ierr); 6460cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr); 647c8b0795cSMark F. Adams PetscFunctionReturn(0); 648c8b0795cSMark F. Adams } 649c8b0795cSMark F. Adams 650c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 651c8b0795cSMark F. Adams /* 6520cbbd2e1SMark F. Adams PCGAMGProlongator_GEO 653c8b0795cSMark F. Adams 654c8b0795cSMark F. Adams Input Parameter: 655c8b0795cSMark F. Adams . pc - this 656c8b0795cSMark F. Adams . Amat - matrix on this fine level 657c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 6580cbbd2e1SMark F. Adams . selected_1 - [nselected] 6590cbbd2e1SMark F. Adams . agg_lists - [nselected] 660c8b0795cSMark F. Adams Output Parameter: 661c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 662c8b0795cSMark F. Adams */ 6632adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_GEO(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 6642e68589bSMark F. Adams { 6652e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 6662e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 667c8b0795cSMark F. Adams const PetscInt dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols; 6682e68589bSMark F. Adams PetscErrorCode ierr; 669b43b03e9SMark F. Adams PetscInt Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid; 670c8b0795cSMark F. Adams Mat Prol; 671c5df96a5SBarry Smith PetscMPIInt rank, size; 6723b4367a7SBarry Smith MPI_Comm comm; 6730cbbd2e1SMark F. Adams IS selected_2,selected_1; 6742e68589bSMark F. Adams const PetscInt *selected_idx; 675d9558ea9SBarry Smith MatType mtype; 6762e68589bSMark F. Adams 6772e68589bSMark F. Adams PetscFunctionBegin; 6783b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 6790cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr); 6803b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 6813b4367a7SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 6822e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 6832e68589bSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 68471959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 68571959b99SBarry Smith if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) % bs %D",Iend,Istart,bs); 6862e68589bSMark F. Adams 6872e68589bSMark F. Adams /* get 'nLocalSelected' */ 68841b27cdeSMark F. Adams ierr = PetscCDGetMIS(agg_lists, &selected_1);CHKERRQ(ierr); 689b43b03e9SMark F. Adams ierr = ISGetSize(selected_1, &jj);CHKERRQ(ierr); 690785e854fSJed Brown ierr = PetscMalloc1(jj, &clid_flid);CHKERRQ(ierr); 6912e68589bSMark F. Adams ierr = ISGetIndices(selected_1, &selected_idx);CHKERRQ(ierr); 692c8b0795cSMark F. Adams for (kk=0,nLocalSelected=0; kk<jj; kk++) { 6932e68589bSMark F. Adams PetscInt lid = selected_idx[kk]; 694b43b03e9SMark F. Adams if (lid<nloc) { 6950cbbd2e1SMark F. Adams ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr); 6962fa5cd67SKarl Rupp if (ncols>1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */ 6970cbbd2e1SMark F. Adams ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr); 698b43b03e9SMark F. Adams } 6992e68589bSMark F. Adams } 7002e68589bSMark F. Adams ierr = ISRestoreIndices(selected_1, &selected_idx);CHKERRQ(ierr); 701a2f3521dSMark F. Adams ierr = ISDestroy(&selected_1);CHKERRQ(ierr); /* this is selected_1 in serial */ 7022e68589bSMark F. Adams 703d9558ea9SBarry Smith /* create prolongator matrix */ 704d9558ea9SBarry Smith ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr); 7053b4367a7SBarry Smith ierr = MatCreate(comm, &Prol);CHKERRQ(ierr); 706806fa848SBarry Smith ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 707a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Prol, bs, bs);CHKERRQ(ierr); 708d9558ea9SBarry Smith ierr = MatSetType(Prol, mtype);CHKERRQ(ierr); 7090298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(Prol,3*data_cols,NULL);CHKERRQ(ierr); 7100298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(Prol,3*data_cols,NULL,3*data_cols,NULL);CHKERRQ(ierr); 7112e68589bSMark F. Adams 712c8b0795cSMark F. Adams /* can get all points "removed" - but not on geomg */ 7132e68589bSMark F. Adams ierr = MatGetSize(Prol, &kk, &jj);CHKERRQ(ierr); 7147f66b68fSMark Adams if (!jj) { 715bf4339c2SBarry Smith ierr = PetscInfo(pc,"ERROE: no selected points on coarse grid\n");CHKERRQ(ierr); 716b43b03e9SMark F. Adams ierr = PetscFree(clid_flid);CHKERRQ(ierr); 7172e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 7180298fd71SBarry Smith *a_P_out = NULL; /* out */ 7192e68589bSMark F. Adams PetscFunctionReturn(0); 7202e68589bSMark F. Adams } 7212e68589bSMark F. Adams 7222e68589bSMark F. Adams { 7232e68589bSMark F. Adams PetscReal *coords; 724a2f3521dSMark F. Adams PetscInt data_stride; 7250298fd71SBarry Smith PetscInt *crsGID = NULL; 7262e68589bSMark F. Adams Mat Gmat2; 7272e68589bSMark F. Adams 72871959b99SBarry Smith if (dim != data_cols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim %D != data_cols %D",dim,data_cols); 7292e68589bSMark F. Adams /* grow ghost data for better coarse grid cover of fine grid */ 7300cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7310cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 7322e68589bSMark F. Adams #endif 733a2f3521dSMark F. Adams /* messy method, squares graph and gets some data */ 734806fa848SBarry Smith ierr = getGIDsOnSquareGraph(nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID);CHKERRQ(ierr); 7352e68589bSMark F. Adams /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */ 7360cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7370cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 7382e68589bSMark F. Adams #endif 7392e68589bSMark F. Adams /* create global vector of coorindates in 'coords' */ 740c5df96a5SBarry Smith if (size > 1) { 741806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords);CHKERRQ(ierr); 742806fa848SBarry Smith } else { 743c8b0795cSMark F. Adams coords = (PetscReal*)pc_gamg->data; 744a2f3521dSMark F. Adams data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols; 7452e68589bSMark F. Adams } 7462e68589bSMark F. Adams ierr = MatDestroy(&Gmat2);CHKERRQ(ierr); 7472e68589bSMark F. Adams 7482e68589bSMark F. Adams /* triangulate */ 7492e68589bSMark F. Adams if (dim == 2) { 750c8b0795cSMark F. Adams PetscReal metric,tm; 7510cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7520cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr); 7532e68589bSMark F. Adams #endif 754806fa848SBarry Smith ierr = triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric);CHKERRQ(ierr); 7550cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7560cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr); 7572e68589bSMark F. Adams #endif 7582e68589bSMark F. Adams ierr = PetscFree(crsGID);CHKERRQ(ierr); 7592e68589bSMark F. Adams 7602e68589bSMark F. Adams /* clean up and create coordinates for coarse grid (output) */ 761c5df96a5SBarry Smith if (size > 1) ierr = PetscFree(coords);CHKERRQ(ierr); 7622e68589bSMark F. Adams 763b2566f29SBarry Smith ierr = MPIU_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm);CHKERRQ(ierr); 764c8b0795cSMark F. Adams if (tm > 1.) { /* needs to be globalized - should not happen */ 765bf4339c2SBarry Smith ierr = PetscInfo1(pc," failed metric for coarse grid %e\n",(double)tm);CHKERRQ(ierr); 7662e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 767806fa848SBarry Smith } else if (metric > .0) { 768bf4339c2SBarry Smith ierr = PetscInfo1(pc,"worst metric for coarse grid = %e\n",(double)metric);CHKERRQ(ierr); 7692e68589bSMark F. Adams } 7703b4367a7SBarry Smith } else SETERRQ(comm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG"); 7712e68589bSMark F. Adams { /* create next coords - output */ 7722e68589bSMark F. Adams PetscReal *crs_crds; 773785e854fSJed Brown ierr = PetscMalloc1(dim*nLocalSelected, &crs_crds);CHKERRQ(ierr); 7742e68589bSMark F. Adams for (kk=0; kk<nLocalSelected; kk++) { /* grab local select nodes to promote - output */ 775b43b03e9SMark F. Adams PetscInt lid = clid_flid[kk]; 776c8b0795cSMark F. Adams for (jj=0; jj<dim; jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid]; 7772e68589bSMark F. Adams } 778c8b0795cSMark F. Adams 779c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 780c8b0795cSMark F. Adams pc_gamg->data = crs_crds; /* out */ 781c8b0795cSMark F. Adams pc_gamg->data_sz = dim*nLocalSelected; 7822e68589bSMark F. Adams } 783a2f3521dSMark F. Adams ierr = ISDestroy(&selected_2);CHKERRQ(ierr); 7842e68589bSMark F. Adams } 785a2f3521dSMark F. Adams 7862e68589bSMark F. Adams *a_P_out = Prol; /* out */ 787b43b03e9SMark F. Adams ierr = PetscFree(clid_flid);CHKERRQ(ierr); 7880cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr); 789c8b0795cSMark F. Adams PetscFunctionReturn(0); 790c8b0795cSMark F. Adams } 791c8b0795cSMark F. Adams 7929b8ffb57SJed Brown static PetscErrorCode PCDestroy_GAMG_GEO(PC pc) 7939b8ffb57SJed Brown { 7949b8ffb57SJed Brown PetscErrorCode ierr; 7959b8ffb57SJed Brown 7969b8ffb57SJed Brown PetscFunctionBegin; 7979b8ffb57SJed Brown ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr); 7989b8ffb57SJed Brown PetscFunctionReturn(0); 7999b8ffb57SJed Brown } 8009b8ffb57SJed Brown 801c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 802c8b0795cSMark F. Adams /* 803c8b0795cSMark F. Adams PCCreateGAMG_GEO 804c8b0795cSMark F. Adams 805c8b0795cSMark F. Adams Input Parameter: 806c8b0795cSMark F. Adams . pc - 807c8b0795cSMark F. Adams */ 808c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_GEO(PC pc) 809c8b0795cSMark F. Adams { 810c8b0795cSMark F. Adams PetscErrorCode ierr; 811c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 812c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 813c8b0795cSMark F. Adams 814c8b0795cSMark F. Adams PetscFunctionBegin; 8151ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO; 8169b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_GEO; 817c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 818c8b0795cSMark F. Adams 819c8b0795cSMark F. Adams /* set internal function pointers */ 820fd1112cbSBarry Smith pc_gamg->ops->graph = PCGAMGGraph_GEO; 821fd1112cbSBarry Smith pc_gamg->ops->coarsen = PCGAMGCoarsen_GEO; 8221ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_GEO; 823fd1112cbSBarry Smith pc_gamg->ops->optprolongator = NULL; 8241ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_GEO; 825c8b0795cSMark F. Adams 826bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_GEO);CHKERRQ(ierr); 8272e68589bSMark F. Adams PetscFunctionReturn(0); 8282e68589bSMark F. Adams } 829