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) 82e68589bSMark F. Adams #define REAL PetscReal 92e68589bSMark F. Adams #include <triangle.h> 102e68589bSMark F. Adams #endif 112e68589bSMark F. Adams 122e68589bSMark F. Adams #include <petscblaslapack.h> 132e68589bSMark F. Adams 14c8b0795cSMark F. Adams /* Private context for the GAMG preconditioner */ 15c8b0795cSMark F. Adams typedef struct { 16c8b0795cSMark F. Adams PetscInt lid; /* local vertex index */ 17c8b0795cSMark F. Adams PetscInt degree; /* vertex degree */ 18c8b0795cSMark F. Adams } GAMGNode; 192fa5cd67SKarl Rupp 200cbbd2e1SMark F. Adams int petsc_geo_mg_compare(const void *a, const void *b) 21c8b0795cSMark F. Adams { 22c8b0795cSMark F. Adams return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree); 23c8b0795cSMark F. Adams } 24c8b0795cSMark F. Adams 252e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 262e68589bSMark F. Adams /* 272e68589bSMark F. Adams PCSetCoordinates_GEO 282e68589bSMark F. Adams 292e68589bSMark F. Adams Input Parameter: 302e68589bSMark F. Adams . pc - the preconditioner context 312e68589bSMark F. Adams */ 322e68589bSMark F. Adams #undef __FUNCT__ 332e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_GEO" 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; 392e68589bSMark F. Adams PetscInt arrsz,bs,my0,kk,ii,nloc,Iend; 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); 45a2f3521dSMark F. Adams 462e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr); 472e68589bSMark F. Adams nloc = (Iend-my0)/bs; 48a2f3521dSMark F. Adams 49ce94432eSBarry Smith if (nloc!=a_nloc) SETERRQ2(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Stokes not supported nloc = %d %d.",a_nloc,nloc); 50ce94432eSBarry Smith if ((Iend-my0)%bs!=0) SETERRQ1(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 512e68589bSMark F. Adams 52c8b0795cSMark F. Adams pc_gamg->data_cell_rows = 1; 537f66b68fSMark Adams if (!coords && nloc > 0) SETERRQ(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'."); 54c8b0795cSMark F. Adams pc_gamg->data_cell_cols = ndm; /* coordinates */ 552e68589bSMark F. Adams 56c8b0795cSMark F. Adams arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 572e68589bSMark F. Adams 582e68589bSMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 597f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 602e68589bSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 61854ce69bSBarry Smith ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr); 622e68589bSMark F. Adams } 632e68589bSMark F. Adams for (kk=0; kk<arrsz; kk++) pc_gamg->data[kk] = -999.; 642e68589bSMark F. Adams pc_gamg->data[arrsz] = -99.; 652e68589bSMark F. Adams /* copy data in - column oriented */ 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 } 7171959b99SBarry 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]); 722e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 732e68589bSMark F. Adams PetscFunctionReturn(0); 742e68589bSMark F. Adams } 752e68589bSMark F. Adams 762e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 772e68589bSMark F. Adams /* 782e68589bSMark F. Adams PCSetData_GEO 792e68589bSMark F. Adams 802e68589bSMark F. Adams Input Parameter: 812e68589bSMark F. Adams . pc - 822e68589bSMark F. Adams */ 832e68589bSMark F. Adams #undef __FUNCT__ 842e68589bSMark F. Adams #define __FUNCT__ "PCSetData_GEO" 85b8cd405aSMark F. Adams PetscErrorCode PCSetData_GEO(PC pc, Mat m) 862e68589bSMark F. Adams { 872e68589bSMark F. Adams PetscFunctionBegin; 88ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"GEO MG needs coordinates"); 892e68589bSMark F. Adams } 902e68589bSMark F. Adams 912e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 922e68589bSMark F. Adams /* 932e68589bSMark F. Adams PCSetFromOptions_GEO 942e68589bSMark F. Adams 952e68589bSMark F. Adams Input Parameter: 962e68589bSMark F. Adams . pc - 972e68589bSMark F. Adams */ 982e68589bSMark F. Adams #undef __FUNCT__ 992e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GEO" 1004416b707SBarry Smith PetscErrorCode PCSetFromOptions_GEO(PetscOptionItems *PetscOptionsObject,PC pc) 1012e68589bSMark F. Adams { 1022e68589bSMark F. Adams PetscErrorCode ierr; 1032e68589bSMark F. Adams 1042e68589bSMark F. Adams PetscFunctionBegin; 105e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-GEO options");CHKERRQ(ierr); 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 /* CHKERRQ(ierr); */ 1162e68589bSMark F. Adams } 1172e68589bSMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1182e68589bSMark F. Adams PetscFunctionReturn(0); 1192e68589bSMark F. Adams } 1202e68589bSMark F. Adams 1212e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1222e68589bSMark F. Adams /* 1232e68589bSMark F. Adams triangulateAndFormProl 1242e68589bSMark F. Adams 1252e68589bSMark F. Adams Input Parameter: 1262e68589bSMark F. Adams . selected_2 - list of selected local ID, includes selected ghosts 127a2f3521dSMark F. Adams . data_stride - 128a2f3521dSMark F. Adams . coords[2*data_stride] - column vector of local coordinates w/ ghosts 1292adfe9d3SBarry Smith . nselected_1 - selected IDs that go with base (1) graph includes selected ghosts 1300cbbd2e1SMark F. Adams . clid_lid_1[nselected_1] - lids of selected (c) nodes ??????????? 1312adfe9d3SBarry Smith . agg_lists_1 - list of aggregates selected_1 vertices of aggregate unselected vertices 1322e68589bSMark F. Adams . crsGID[selected.size()] - global index for prolongation operator 1332e68589bSMark F. Adams . bs - block size 1342e68589bSMark F. Adams Output Parameter: 1352e68589bSMark F. Adams . a_Prol - prolongation operator 1362e68589bSMark F. Adams . a_worst_best - measure of worst missed fine vertex, 0 is no misses 1372e68589bSMark F. Adams */ 1382e68589bSMark F. Adams #undef __FUNCT__ 1392e68589bSMark F. Adams #define __FUNCT__ "triangulateAndFormProl" 1402adfe9d3SBarry 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, 1412adfe9d3SBarry Smith const PetscInt crsGID[],PetscInt bs,Mat a_Prol,PetscReal *a_worst_best) 1422e68589bSMark F. Adams { 1432e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE) 1442e68589bSMark F. Adams PetscErrorCode ierr; 145b43b03e9SMark F. Adams PetscInt jj,tid,tt,idx,nselected_2; 1462e68589bSMark F. Adams struct triangulateio in,mid; 1470cbbd2e1SMark F. Adams const PetscInt *selected_idx_2; 14873911c69SBarry Smith PetscMPIInt rank; 1492e68589bSMark F. Adams PetscInt Istart,Iend,nFineLoc,myFine0; 1502e68589bSMark F. Adams int kk,nPlotPts,sid; 1513b4367a7SBarry Smith MPI_Comm comm; 152c8b0795cSMark F. Adams PetscReal tm; 153c8b0795cSMark F. Adams 1546e111a19SKarl Rupp PetscFunctionBegin; 1553b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr); 1563b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 157c8b0795cSMark F. Adams ierr = ISGetSize(selected_2, &nselected_2);CHKERRQ(ierr); 1582e68589bSMark F. Adams if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */ 1592e68589bSMark F. Adams *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */ 160806fa848SBarry Smith } else *a_worst_best = 0.0; 161b2566f29SBarry Smith ierr = MPIU_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm);CHKERRQ(ierr); 162c8b0795cSMark F. Adams if (tm > 0.0) { 163c8b0795cSMark F. Adams *a_worst_best = 100.0; 1642e68589bSMark F. Adams PetscFunctionReturn(0); 1652e68589bSMark F. Adams } 1662e68589bSMark F. Adams ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr); 1672e68589bSMark F. Adams nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs; 1682e68589bSMark F. Adams nPlotPts = nFineLoc; /* locals */ 1692e68589bSMark F. Adams /* traingle */ 1702e68589bSMark F. Adams /* Define input points - in*/ 1712e68589bSMark F. Adams in.numberofpoints = nselected_2; 1722e68589bSMark F. Adams in.numberofpointattributes = 0; 1732e68589bSMark F. Adams /* get nselected points */ 174e632b94dSBarry Smith ierr = PetscMalloc1(2*nselected_2, &in.pointlist);CHKERRQ(ierr); 1752e68589bSMark F. Adams ierr = ISGetIndices(selected_2, &selected_idx_2);CHKERRQ(ierr); 1762e68589bSMark F. Adams 1772e68589bSMark F. Adams for (kk=0,sid=0; kk<nselected_2; kk++,sid += 2) { 1782e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 1792e68589bSMark F. Adams in.pointlist[sid] = coords[lid]; 180a2f3521dSMark F. Adams in.pointlist[sid+1] = coords[data_stride + lid]; 1812e68589bSMark F. Adams if (lid>=nFineLoc) nPlotPts++; 1822e68589bSMark F. Adams } 18371959b99SBarry Smith if (sid != 2*nselected_2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != 2*nselected_2 %D",sid,nselected_2); 1842e68589bSMark F. Adams 1852e68589bSMark F. Adams in.numberofsegments = 0; 1862e68589bSMark F. Adams in.numberofedges = 0; 1872e68589bSMark F. Adams in.numberofholes = 0; 1882e68589bSMark F. Adams in.numberofregions = 0; 1892e68589bSMark F. Adams in.trianglelist = 0; 1902e68589bSMark F. Adams in.segmentmarkerlist = 0; 1912e68589bSMark F. Adams in.pointattributelist = 0; 1922e68589bSMark F. Adams in.pointmarkerlist = 0; 1932e68589bSMark F. Adams in.triangleattributelist = 0; 1942e68589bSMark F. Adams in.trianglearealist = 0; 1952e68589bSMark F. Adams in.segmentlist = 0; 1962e68589bSMark F. Adams in.holelist = 0; 1972e68589bSMark F. Adams in.regionlist = 0; 1982e68589bSMark F. Adams in.edgelist = 0; 1992e68589bSMark F. Adams in.edgemarkerlist = 0; 2002e68589bSMark F. Adams in.normlist = 0; 2012fa5cd67SKarl Rupp 2022e68589bSMark F. Adams /* triangulate */ 2032e68589bSMark F. Adams mid.pointlist = 0; /* Not needed if -N switch used. */ 2042e68589bSMark F. Adams /* Not needed if -N switch used or number of point attributes is zero: */ 2052e68589bSMark F. Adams mid.pointattributelist = 0; 2062e68589bSMark F. Adams mid.pointmarkerlist = 0; /* Not needed if -N or -B switch used. */ 2072e68589bSMark F. Adams mid.trianglelist = 0; /* Not needed if -E switch used. */ 2082e68589bSMark F. Adams /* Not needed if -E switch used or number of triangle attributes is zero: */ 2092e68589bSMark F. Adams mid.triangleattributelist = 0; 2102e68589bSMark F. Adams mid.neighborlist = 0; /* Needed only if -n switch used. */ 2112e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P not used: */ 2122e68589bSMark F. Adams mid.segmentlist = 0; 2132e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P and -B not used: */ 2142e68589bSMark F. Adams mid.segmentmarkerlist = 0; 2152e68589bSMark F. Adams mid.edgelist = 0; /* Needed only if -e switch used. */ 2162e68589bSMark F. Adams mid.edgemarkerlist = 0; /* Needed if -e used and -B not used. */ 2172e68589bSMark F. Adams mid.numberoftriangles = 0; 2182e68589bSMark F. Adams 2192e68589bSMark F. Adams /* Triangulate the points. Switches are chosen to read and write a */ 2202e68589bSMark F. Adams /* PSLG (p), preserve the convex hull (c), number everything from */ 2212e68589bSMark F. Adams /* zero (z), assign a regional attribute to each element (A), and */ 2222e68589bSMark F. Adams /* produce an edge list (e), a Voronoi diagram (v), and a triangle */ 2232e68589bSMark F. Adams /* neighbor list (n). */ 2242e68589bSMark F. Adams if (nselected_2 != 0) { /* inactive processor */ 2252e68589bSMark F. Adams char args[] = "npczQ"; /* c is needed ? */ 2262e68589bSMark F. Adams triangulate(args, &in, &mid, (struct triangulateio*) NULL); 2272e68589bSMark F. Adams /* output .poly files for 'showme' */ 2282e68589bSMark F. Adams if (!PETSC_TRUE) { 2292e68589bSMark F. Adams static int level = 1; 2302e68589bSMark F. Adams FILE *file; char fname[32]; 2312e68589bSMark F. Adams 232c5df96a5SBarry Smith sprintf(fname,"C%d_%d.poly",level,rank); file = fopen(fname, "w"); 2332e68589bSMark F. Adams /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/ 2342e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); 2352e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2362e68589bSMark F. Adams for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid += 2) { 2372e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2382e68589bSMark F. Adams } 2392e68589bSMark F. Adams /*One line: <# of segments> <# of boundary markers (0 or 1)> */ 2402e68589bSMark F. Adams fprintf(file, "%d %d\n",0,0); 2412e68589bSMark F. Adams /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */ 2422e68589bSMark F. Adams /* One line: <# of holes> */ 2432e68589bSMark F. Adams fprintf(file, "%d\n",0); 2442e68589bSMark F. Adams /* Following lines: <hole #> <x> <y> */ 2452e68589bSMark F. Adams /* Optional line: <# of regional attributes and/or area constraints> */ 2462e68589bSMark F. Adams /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */ 2472e68589bSMark F. Adams fclose(file); 2482e68589bSMark F. Adams 2492e68589bSMark F. Adams /* elems */ 250c5df96a5SBarry Smith sprintf(fname,"C%d_%d.ele",level,rank); file = fopen(fname, "w"); 2512e68589bSMark F. Adams /* First line: <# of triangles> <nodes per triangle> <# of attributes> */ 2522e68589bSMark F. Adams fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0); 2532e68589bSMark F. Adams /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */ 2542e68589bSMark F. Adams for (kk=0,sid=0; kk<mid.numberoftriangles; kk++,sid += 3) { 2552e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]); 2562e68589bSMark F. Adams } 2572e68589bSMark F. Adams fclose(file); 2582e68589bSMark F. Adams 259c5df96a5SBarry Smith sprintf(fname,"C%d_%d.node",level,rank); file = fopen(fname, "w"); 2602e68589bSMark F. Adams /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */ 2612e68589bSMark F. Adams /* fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); */ 2622e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",nPlotPts,2,0,0); 2632e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2642e68589bSMark F. Adams for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid+=2) { 2652e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2662e68589bSMark F. Adams } 2672e68589bSMark F. Adams 2682e68589bSMark F. Adams sid /= 2; 2692e68589bSMark F. Adams for (jj=0; jj<nFineLoc; jj++) { 2702e68589bSMark F. Adams PetscBool sel = PETSC_TRUE; 2712e68589bSMark F. Adams for (kk=0; kk<nselected_2 && sel; kk++) { 2722e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 2732e68589bSMark F. Adams if (lid == jj) sel = PETSC_FALSE; 2742e68589bSMark F. Adams } 2752fa5cd67SKarl Rupp if (sel) fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]); 2762e68589bSMark F. Adams } 2772e68589bSMark F. Adams fclose(file); 27871959b99SBarry Smith if (sid != nPlotPts) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != nPlotPts %D",sid,nPlotPts); 2792e68589bSMark F. Adams level++; 2802e68589bSMark F. Adams } 2812e68589bSMark F. Adams } 2820cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 2830cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 2842e68589bSMark F. Adams #endif 2852e68589bSMark F. Adams { /* form P - setup some maps */ 2860cbbd2e1SMark F. Adams PetscInt clid,mm,*nTri,*node_tri; 2872e68589bSMark F. Adams 288e632b94dSBarry Smith ierr = PetscMalloc2(nselected_2, &node_tri,nselected_2, &nTri);CHKERRQ(ierr); 2892e68589bSMark F. Adams 2902e68589bSMark F. Adams /* need list of triangles on node */ 2912e68589bSMark F. Adams for (kk=0; kk<nselected_2; kk++) nTri[kk] = 0; 2922e68589bSMark F. Adams for (tid=0,kk=0; tid<mid.numberoftriangles; tid++) { 2932e68589bSMark F. Adams for (jj=0; jj<3; jj++) { 2942e68589bSMark F. Adams PetscInt cid = mid.trianglelist[kk++]; 2952e68589bSMark F. Adams if (nTri[cid] == 0) node_tri[cid] = tid; 2962e68589bSMark F. Adams nTri[cid]++; 2972e68589bSMark F. Adams } 2982e68589bSMark F. Adams } 2992e68589bSMark F. Adams #define EPS 1.e-12 3002e68589bSMark F. Adams /* find points and set prolongation */ 3010cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nFineLoc; mm++) { 302e78576d6SMark F. Adams PetscBool ise; 303e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise);CHKERRQ(ierr); 304e78576d6SMark F. Adams if (!ise) { 3050cbbd2e1SMark F. Adams const PetscInt lid = mm; 3062e68589bSMark F. Adams PetscScalar AA[3][3]; 3072e68589bSMark F. Adams PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO; 308539c167fSBarry Smith PetscCDIntNd *pos; 309539c167fSBarry Smith 310e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos);CHKERRQ(ierr); 311e78576d6SMark F. Adams while (pos) { 312e78576d6SMark F. Adams PetscInt flid; 313484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &flid);CHKERRQ(ierr); 314e78576d6SMark F. Adams ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos);CHKERRQ(ierr); 3150cbbd2e1SMark F. Adams 3162e68589bSMark F. Adams if (flid < nFineLoc) { /* could be a ghost */ 3172e68589bSMark F. Adams PetscInt bestTID = -1; PetscReal best_alpha = 1.e10; 3182e68589bSMark F. Adams const PetscInt fgid = flid + myFine0; 3192e68589bSMark F. Adams /* compute shape function for gid */ 320a2f3521dSMark F. Adams const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0}; 3212e68589bSMark F. Adams PetscBool haveit =PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3]; 3222fa5cd67SKarl Rupp 3232e68589bSMark F. Adams /* look for it */ 3240cbbd2e1SMark F. Adams for (tid = node_tri[clid], jj=0; 3252e68589bSMark F. Adams jj < 5 && !haveit && tid != -1; 3262e68589bSMark F. Adams jj++) { 3272e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3282e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3292e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 330a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3312e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3322e68589bSMark F. Adams } 3332e68589bSMark F. Adams 3342e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = (PetscScalar)fcoord[tt]; 3352e68589bSMark F. Adams 3362e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 3378b83055fSJed Brown PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3382e68589bSMark F. Adams { 3392e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal lowest=1.e10; 3402e68589bSMark F. Adams for (tt = 0, idx = 0; tt < 3; tt++) { 3412e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE; 3422e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) < lowest) { 3432e68589bSMark F. Adams lowest = PetscRealPart(alpha[tt]); 3442e68589bSMark F. Adams idx = tt; 3452e68589bSMark F. Adams } 3462e68589bSMark F. Adams } 3472e68589bSMark F. Adams haveit = have; 3482e68589bSMark F. Adams } 3492e68589bSMark F. Adams tid = mid.neighborlist[3*tid + idx]; 3502e68589bSMark F. Adams } 3512e68589bSMark F. Adams 3522e68589bSMark F. Adams if (!haveit) { 3532e68589bSMark F. Adams /* brute force */ 3542e68589bSMark F. Adams for (tid=0; tid<mid.numberoftriangles && !haveit; tid++) { 3552e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3562e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3572e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 358a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3592e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3602e68589bSMark F. Adams } 3612e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt]; 3622e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 3638b83055fSJed Brown PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3642e68589bSMark F. Adams { 3652e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal worst=0.0, v; 3662e68589bSMark F. Adams for (tt=0; tt<3 && have; tt++) { 3672e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE; 3682e68589bSMark F. Adams if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v; 3692e68589bSMark F. Adams } 3702e68589bSMark F. Adams if (worst < best_alpha) { 3712e68589bSMark F. Adams best_alpha = worst; bestTID = tid; 3722e68589bSMark F. Adams } 3732e68589bSMark F. Adams haveit = have; 3742e68589bSMark F. Adams } 3752e68589bSMark F. Adams } 3762e68589bSMark F. Adams } 3772e68589bSMark F. Adams if (!haveit) { 3782e68589bSMark F. Adams if (best_alpha > *a_worst_best) *a_worst_best = best_alpha; 3792e68589bSMark F. Adams /* use best one */ 3802e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3812e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*bestTID + tt]; 3822e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 383a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3842e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3852e68589bSMark F. Adams } 3862e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt]; 3872e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 3888b83055fSJed Brown PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3892e68589bSMark F. Adams } 3902e68589bSMark F. Adams 3912e68589bSMark F. Adams /* put in row of P */ 3922e68589bSMark F. Adams for (idx=0; idx<3; idx++) { 3932e68589bSMark F. Adams PetscScalar shp = alpha[idx]; 3942e68589bSMark F. Adams if (PetscAbs(PetscRealPart(shp)) > 1.e-6) { 3952e68589bSMark F. Adams PetscInt cgid = crsGID[clids[idx]]; 3962e68589bSMark F. Adams PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */ 3972e68589bSMark F. Adams for (tt=0; tt < bs; tt++, ii++, jj++) { 3982e68589bSMark F. Adams ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES);CHKERRQ(ierr); 3992e68589bSMark F. Adams } 4002e68589bSMark F. Adams } 4012e68589bSMark F. Adams } 4022e68589bSMark F. Adams } 4030cbbd2e1SMark F. Adams } /* aggregates iterations */ 4040cbbd2e1SMark F. Adams clid++; 4050cbbd2e1SMark F. Adams } /* a coarse agg */ 4060cbbd2e1SMark F. Adams } /* for all fine nodes */ 4070cbbd2e1SMark F. Adams 4082e68589bSMark F. Adams ierr = ISRestoreIndices(selected_2, &selected_idx_2);CHKERRQ(ierr); 4092e68589bSMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4102e68589bSMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4112e68589bSMark F. Adams 412e632b94dSBarry Smith ierr = PetscFree2(node_tri,nTri);CHKERRQ(ierr); 4132e68589bSMark F. Adams } 4140cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4150cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 4162e68589bSMark F. Adams #endif 4172e68589bSMark F. Adams free(mid.trianglelist); 4182e68589bSMark F. Adams free(mid.neighborlist); 4192e68589bSMark F. Adams ierr = PetscFree(in.pointlist);CHKERRQ(ierr); 4202e68589bSMark F. Adams PetscFunctionReturn(0); 4212e68589bSMark F. Adams #else 4223b4367a7SBarry Smith SETERRQ(PetscObjectComm((PetscObject)a_Prol),PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG"); 4232e68589bSMark F. Adams #endif 4242e68589bSMark F. Adams } 4252e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 4262e68589bSMark F. Adams /* 4272e68589bSMark F. Adams getGIDsOnSquareGraph - square graph, get 4282e68589bSMark F. Adams 4292e68589bSMark F. Adams Input Parameter: 4300cbbd2e1SMark F. Adams . nselected_1 - selected local indices (includes ghosts in input Gmat1) 431b43b03e9SMark F. Adams . clid_lid_1 - [nselected_1] lids of selected nodes 4322e68589bSMark F. Adams . Gmat1 - graph that goes with 'selected_1' 4332e68589bSMark F. Adams Output Parameter: 4342e68589bSMark F. Adams . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2) 4352e68589bSMark F. Adams . a_Gmat_2 - graph that is squared of 'Gmat_1' 4362e68589bSMark F. Adams . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes 4372e68589bSMark F. Adams */ 4382e68589bSMark F. Adams #undef __FUNCT__ 4392e68589bSMark F. Adams #define __FUNCT__ "getGIDsOnSquareGraph" 4402adfe9d3SBarry 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) 4412e68589bSMark F. Adams { 4422e68589bSMark F. Adams PetscErrorCode ierr; 44373911c69SBarry Smith PetscMPIInt size; 444b43b03e9SMark F. Adams PetscInt *crsGID, kk,my0,Iend,nloc; 4453b4367a7SBarry Smith MPI_Comm comm; 4462e68589bSMark F. Adams 4472e68589bSMark F. Adams PetscFunctionBegin; 4483b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr); 4493b4367a7SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 4502e68589bSMark F. Adams ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend);CHKERRQ(ierr); /* AIJ */ 4512e68589bSMark F. Adams nloc = Iend - my0; /* this does not change */ 4522e68589bSMark F. Adams 453c5df96a5SBarry Smith if (size == 1) { /* not much to do in serial */ 454785e854fSJed Brown ierr = PetscMalloc1(nselected_1, &crsGID);CHKERRQ(ierr); 455b43b03e9SMark F. Adams for (kk=0; kk<nselected_1; kk++) crsGID[kk] = kk; 4562e68589bSMark F. Adams *a_Gmat_2 = 0; 457806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);CHKERRQ(ierr); 458806fa848SBarry Smith } else { 459b43b03e9SMark F. Adams PetscInt idx,num_fine_ghosts,num_crs_ghost,myCrs0; 4602e68589bSMark F. Adams Mat_MPIAIJ *mpimat2; 4612e68589bSMark F. Adams Mat Gmat2; 4622e68589bSMark F. Adams Vec locState; 4632e68589bSMark F. Adams PetscScalar *cpcol_state; 4642e68589bSMark F. Adams 4652e68589bSMark F. Adams /* scan my coarse zero gid, set 'lid_state' with coarse GID */ 466b43b03e9SMark F. Adams kk = nselected_1; 467dbd2ff41SBarry Smith ierr = MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr); 468b43b03e9SMark F. Adams myCrs0 -= nselected_1; 4692e68589bSMark F. Adams 470b43b03e9SMark F. Adams if (a_Gmat_2) { /* output */ 4712e68589bSMark F. Adams /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */ 4722e68589bSMark F. Adams ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr); 4732e68589bSMark F. Adams *a_Gmat_2 = Gmat2; /* output */ 474806fa848SBarry Smith } else Gmat2 = Gmat1; /* use local to get crsGIDs at least */ 4752e68589bSMark F. Adams /* get coarse grid GIDS for selected (locals and ghosts) */ 4762e68589bSMark F. Adams mpimat2 = (Mat_MPIAIJ*)Gmat2->data; 4772a7a6963SBarry Smith ierr = MatCreateVecs(Gmat2, &locState, 0);CHKERRQ(ierr); 4782e68589bSMark F. Adams ierr = VecSet(locState, (PetscScalar)(PetscReal)(-1));CHKERRQ(ierr); /* set with UNKNOWN state */ 479b43b03e9SMark F. Adams for (kk=0; kk<nselected_1; kk++) { 480b43b03e9SMark F. Adams PetscInt fgid = clid_lid_1[kk] + my0; 4812e68589bSMark F. Adams PetscScalar v = (PetscScalar)(kk+myCrs0); 4822e68589bSMark F. Adams ierr = VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES);CHKERRQ(ierr); /* set with PID */ 4832e68589bSMark F. Adams } 4842e68589bSMark F. Adams ierr = VecAssemblyBegin(locState);CHKERRQ(ierr); 4852e68589bSMark F. Adams ierr = VecAssemblyEnd(locState);CHKERRQ(ierr); 4862e68589bSMark F. Adams ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4872e68589bSMark F. Adams ierr = VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4882e68589bSMark F. Adams ierr = VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts);CHKERRQ(ierr); 4892e68589bSMark F. Adams ierr = VecGetArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr); 4902e68589bSMark F. Adams for (kk=0,num_crs_ghost=0; kk<num_fine_ghosts; kk++) { 4912e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++; 4922e68589bSMark F. Adams } 493854ce69bSBarry Smith ierr = PetscMalloc1(nselected_1+num_crs_ghost, &crsGID);CHKERRQ(ierr); /* output */ 4942e68589bSMark F. Adams { 4952e68589bSMark F. Adams PetscInt *selected_set; 496854ce69bSBarry Smith ierr = PetscMalloc1(nselected_1+num_crs_ghost, &selected_set);CHKERRQ(ierr); 4972e68589bSMark F. Adams /* do ghost of 'crsGID' */ 498b43b03e9SMark F. Adams for (kk=0,idx=nselected_1; kk<num_fine_ghosts; kk++) { 4992e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) { 5002e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5012e68589bSMark F. Adams selected_set[idx] = nloc + kk; 5022e68589bSMark F. Adams crsGID[idx++] = cgid; 5032e68589bSMark F. Adams } 5042e68589bSMark F. Adams } 50571959b99SBarry 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); 5062e68589bSMark F. Adams ierr = VecRestoreArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr); 5072e68589bSMark F. Adams /* do locals in 'crsGID' */ 5082e68589bSMark F. Adams ierr = VecGetArray(locState, &cpcol_state);CHKERRQ(ierr); 5092e68589bSMark F. Adams for (kk=0,idx=0; kk<nloc; kk++) { 5102e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) { 5112e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5122e68589bSMark F. Adams selected_set[idx] = kk; 5132e68589bSMark F. Adams crsGID[idx++] = cgid; 5142e68589bSMark F. Adams } 5152e68589bSMark F. Adams } 51671959b99SBarry Smith if (idx != nselected_1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != nselected_1 %D",idx,nselected_1); 5172e68589bSMark F. Adams ierr = VecRestoreArray(locState, &cpcol_state);CHKERRQ(ierr); 5182e68589bSMark F. Adams 5192e68589bSMark F. Adams if (a_selected_2 != 0) { /* output */ 520806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);CHKERRQ(ierr); 521806fa848SBarry Smith } else { 5222e68589bSMark F. Adams ierr = PetscFree(selected_set);CHKERRQ(ierr); 5232e68589bSMark F. Adams } 5240cbbd2e1SMark F. Adams } 5252e68589bSMark F. Adams ierr = VecDestroy(&locState);CHKERRQ(ierr); 5262e68589bSMark F. Adams } 5272e68589bSMark F. Adams *a_crsGID = crsGID; /* output */ 5282e68589bSMark F. Adams PetscFunctionReturn(0); 5292e68589bSMark F. Adams } 5302e68589bSMark F. Adams 5312e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 5322e68589bSMark F. Adams /* 533fd1112cbSBarry Smith PCGAMGGraph_GEO 5342e68589bSMark F. Adams 5352e68589bSMark F. Adams Input Parameter: 5362e68589bSMark F. Adams . pc - this 5372e68589bSMark F. Adams . Amat - matrix on this fine level 5382e68589bSMark F. Adams Output Parameter: 539c8b0795cSMark F. Adams . a_Gmat 5402e68589bSMark F. Adams */ 5412e68589bSMark F. Adams #undef __FUNCT__ 542fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGGraph_GEO" 5432adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_GEO(PC pc,Mat Amat,Mat *a_Gmat) 544c8b0795cSMark F. Adams { 545c8b0795cSMark F. Adams PetscErrorCode ierr; 546c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 547c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 548*c1eae691SMark Adams const PetscReal vfilter = pc_gamg->threshold[0]; 5493b4367a7SBarry Smith MPI_Comm comm; 550c8b0795cSMark F. Adams Mat Gmat; 5510cbbd2e1SMark F. Adams PetscBool set,flg,symm; 5526e111a19SKarl Rupp 553c8b0795cSMark F. Adams PetscFunctionBegin; 5543b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 555fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGGraph_GEO,0,0,0,0);CHKERRQ(ierr); 556c8b0795cSMark F. Adams 5570cbbd2e1SMark F. Adams ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); 558263489e9SJed Brown symm = (PetscBool)!(set && flg); 5590cbbd2e1SMark F. Adams 5602d7fac45SMark F. Adams ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr); 561bf4339c2SBarry Smith ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr); 562c8b0795cSMark F. Adams 563c8b0795cSMark F. Adams *a_Gmat = Gmat; 564fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGGraph_GEO,0,0,0,0);CHKERRQ(ierr); 565c8b0795cSMark F. Adams PetscFunctionReturn(0); 566c8b0795cSMark F. Adams } 567c8b0795cSMark F. Adams 568c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 569c8b0795cSMark F. Adams /* 570fd1112cbSBarry Smith PCGAMGCoarsen_GEO 571c8b0795cSMark F. Adams 572c8b0795cSMark F. Adams Input Parameter: 573e0940f08SMark F. Adams . a_pc - this 574e0940f08SMark F. Adams . a_Gmat - graph 575c8b0795cSMark F. Adams Output Parameter: 576c8b0795cSMark F. Adams . a_llist_parent - linked list from selected indices for data locality only 577c8b0795cSMark F. Adams */ 578c8b0795cSMark F. Adams #undef __FUNCT__ 579fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGCoarsen_GEO" 580fd1112cbSBarry Smith PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc,Mat *a_Gmat,PetscCoarsenData **a_llist_parent) 581c8b0795cSMark F. Adams { 582c8b0795cSMark F. Adams PetscErrorCode ierr; 583c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,kk,Ii,ncols; 5840cbbd2e1SMark F. Adams IS perm; 585c8b0795cSMark F. Adams GAMGNode *gnodes; 586c8b0795cSMark F. Adams PetscInt *permute; 587e0940f08SMark F. Adams Mat Gmat = *a_Gmat; 5883b4367a7SBarry Smith MPI_Comm comm; 589b43b03e9SMark F. Adams MatCoarsen crs; 590c8b0795cSMark F. Adams 591c8b0795cSMark F. Adams PetscFunctionBegin; 5923b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)a_pc,&comm);CHKERRQ(ierr); 5930cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr); 594c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr); 595c8b0795cSMark F. Adams nloc = (Iend-Istart); 596c8b0795cSMark F. Adams 597c8b0795cSMark F. Adams /* create random permutation with sort for geo-mg */ 598785e854fSJed Brown ierr = PetscMalloc1(nloc, &gnodes);CHKERRQ(ierr); 599785e854fSJed Brown ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr); 600c8b0795cSMark F. Adams 601c8b0795cSMark F. Adams for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */ 602c8b0795cSMark F. Adams ierr = MatGetRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr); 603c8b0795cSMark F. Adams { 604c8b0795cSMark F. Adams PetscInt lid = Ii - Istart; 605c8b0795cSMark F. Adams gnodes[lid].lid = lid; 606c8b0795cSMark F. Adams gnodes[lid].degree = ncols; 607c8b0795cSMark F. Adams } 608c8b0795cSMark F. Adams ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr); 609c8b0795cSMark F. Adams } 610c8b0795cSMark F. Adams if (PETSC_TRUE) { 611e2c00dcbSBarry Smith PetscRandom rand; 612c8b0795cSMark F. Adams PetscBool *bIndexSet; 613e2c00dcbSBarry Smith PetscReal rr; 614e2c00dcbSBarry Smith PetscInt iSwapIndex; 615e2c00dcbSBarry Smith 616e2c00dcbSBarry Smith ierr = PetscRandomCreate(comm,&rand);CHKERRQ(ierr); 617e2c00dcbSBarry Smith ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr); 6182fa5cd67SKarl Rupp for (Ii = 0; Ii < nloc; Ii++) { 619e2c00dcbSBarry Smith ierr = PetscRandomGetValueReal(rand,&rr);CHKERRQ(ierr); 620e2c00dcbSBarry Smith iSwapIndex = (PetscInt) (rr*nloc); 6212fa5cd67SKarl Rupp if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 622c8b0795cSMark F. Adams GAMGNode iTemp = gnodes[iSwapIndex]; 623c8b0795cSMark F. Adams gnodes[iSwapIndex] = gnodes[Ii]; 624c8b0795cSMark F. Adams gnodes[Ii] = iTemp; 625c8b0795cSMark F. Adams bIndexSet[Ii] = PETSC_TRUE; 626c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 627c8b0795cSMark F. Adams } 628c8b0795cSMark F. Adams } 629e2c00dcbSBarry Smith ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 630c8b0795cSMark F. Adams ierr = PetscFree(bIndexSet);CHKERRQ(ierr); 631c8b0795cSMark F. Adams } 632c8b0795cSMark F. Adams /* only sort locals */ 6330cbbd2e1SMark F. Adams qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare); 634c8b0795cSMark F. Adams /* create IS of permutation */ 6352fa5cd67SKarl Rupp for (kk=0; kk<nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */ 636806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr); 637c8b0795cSMark F. Adams 638c8b0795cSMark F. Adams ierr = PetscFree(gnodes);CHKERRQ(ierr); 639c8b0795cSMark F. Adams 640c8b0795cSMark F. Adams /* get MIS aggs */ 641b43b03e9SMark F. Adams 6423b4367a7SBarry Smith ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr); 643b43b03e9SMark F. Adams ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); 644b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr); 645b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr); 646b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs(crs, PETSC_FALSE);CHKERRQ(ierr); 647b43b03e9SMark F. Adams ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 6480cbbd2e1SMark F. Adams ierr = MatCoarsenGetData(crs, a_llist_parent);CHKERRQ(ierr); 649b43b03e9SMark F. Adams ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 650c8b0795cSMark F. Adams 651c8b0795cSMark F. Adams ierr = ISDestroy(&perm);CHKERRQ(ierr); 6520cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr); 653c8b0795cSMark F. Adams PetscFunctionReturn(0); 654c8b0795cSMark F. Adams } 655c8b0795cSMark F. Adams 656c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 657c8b0795cSMark F. Adams /* 6580cbbd2e1SMark F. Adams PCGAMGProlongator_GEO 659c8b0795cSMark F. Adams 660c8b0795cSMark F. Adams Input Parameter: 661c8b0795cSMark F. Adams . pc - this 662c8b0795cSMark F. Adams . Amat - matrix on this fine level 663c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 6640cbbd2e1SMark F. Adams . selected_1 - [nselected] 6650cbbd2e1SMark F. Adams . agg_lists - [nselected] 666c8b0795cSMark F. Adams Output Parameter: 667c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 668c8b0795cSMark F. Adams */ 669c8b0795cSMark F. Adams #undef __FUNCT__ 6700cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_GEO" 6712adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_GEO(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 6722e68589bSMark F. Adams { 6732e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 6742e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 675c8b0795cSMark F. Adams const PetscInt dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols; 6762e68589bSMark F. Adams PetscErrorCode ierr; 677b43b03e9SMark F. Adams PetscInt Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid; 678c8b0795cSMark F. Adams Mat Prol; 679c5df96a5SBarry Smith PetscMPIInt rank, size; 6803b4367a7SBarry Smith MPI_Comm comm; 6810cbbd2e1SMark F. Adams IS selected_2,selected_1; 6822e68589bSMark F. Adams const PetscInt *selected_idx; 683d9558ea9SBarry Smith MatType mtype; 6842e68589bSMark F. Adams 6852e68589bSMark F. Adams PetscFunctionBegin; 6863b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 6870cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr); 6883b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 6893b4367a7SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 6902e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 6912e68589bSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 69271959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 69371959b99SBarry Smith if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) % bs %D",Iend,Istart,bs); 6942e68589bSMark F. Adams 6952e68589bSMark F. Adams /* get 'nLocalSelected' */ 69641b27cdeSMark F. Adams ierr = PetscCDGetMIS(agg_lists, &selected_1);CHKERRQ(ierr); 697b43b03e9SMark F. Adams ierr = ISGetSize(selected_1, &jj);CHKERRQ(ierr); 698785e854fSJed Brown ierr = PetscMalloc1(jj, &clid_flid);CHKERRQ(ierr); 6992e68589bSMark F. Adams ierr = ISGetIndices(selected_1, &selected_idx);CHKERRQ(ierr); 700c8b0795cSMark F. Adams for (kk=0,nLocalSelected=0; kk<jj; kk++) { 7012e68589bSMark F. Adams PetscInt lid = selected_idx[kk]; 702b43b03e9SMark F. Adams if (lid<nloc) { 7030cbbd2e1SMark F. Adams ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr); 7042fa5cd67SKarl Rupp if (ncols>1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */ 7050cbbd2e1SMark F. Adams ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr); 706b43b03e9SMark F. Adams } 7072e68589bSMark F. Adams } 7082e68589bSMark F. Adams ierr = ISRestoreIndices(selected_1, &selected_idx);CHKERRQ(ierr); 709a2f3521dSMark F. Adams ierr = ISDestroy(&selected_1);CHKERRQ(ierr); /* this is selected_1 in serial */ 7102e68589bSMark F. Adams 711d9558ea9SBarry Smith /* create prolongator matrix */ 712d9558ea9SBarry Smith ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr); 7133b4367a7SBarry Smith ierr = MatCreate(comm, &Prol);CHKERRQ(ierr); 714806fa848SBarry Smith ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 715a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Prol, bs, bs);CHKERRQ(ierr); 716d9558ea9SBarry Smith ierr = MatSetType(Prol, mtype);CHKERRQ(ierr); 7170298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(Prol,3*data_cols,NULL);CHKERRQ(ierr); 7180298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(Prol,3*data_cols,NULL,3*data_cols,NULL);CHKERRQ(ierr); 7192e68589bSMark F. Adams 720c8b0795cSMark F. Adams /* can get all points "removed" - but not on geomg */ 7212e68589bSMark F. Adams ierr = MatGetSize(Prol, &kk, &jj);CHKERRQ(ierr); 7227f66b68fSMark Adams if (!jj) { 723bf4339c2SBarry Smith ierr = PetscInfo(pc,"ERROE: no selected points on coarse grid\n");CHKERRQ(ierr); 724b43b03e9SMark F. Adams ierr = PetscFree(clid_flid);CHKERRQ(ierr); 7252e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 7260298fd71SBarry Smith *a_P_out = NULL; /* out */ 7272e68589bSMark F. Adams PetscFunctionReturn(0); 7282e68589bSMark F. Adams } 7292e68589bSMark F. Adams 7302e68589bSMark F. Adams { 7312e68589bSMark F. Adams PetscReal *coords; 732a2f3521dSMark F. Adams PetscInt data_stride; 7330298fd71SBarry Smith PetscInt *crsGID = NULL; 7342e68589bSMark F. Adams Mat Gmat2; 7352e68589bSMark F. Adams 73671959b99SBarry Smith if (dim != data_cols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim %D != data_cols %D",dim,data_cols); 7372e68589bSMark F. Adams /* grow ghost data for better coarse grid cover of fine grid */ 7380cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7390cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 7402e68589bSMark F. Adams #endif 741a2f3521dSMark F. Adams /* messy method, squares graph and gets some data */ 742806fa848SBarry Smith ierr = getGIDsOnSquareGraph(nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID);CHKERRQ(ierr); 7432e68589bSMark F. Adams /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */ 7440cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7450cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 7462e68589bSMark F. Adams #endif 7472e68589bSMark F. Adams /* create global vector of coorindates in 'coords' */ 748c5df96a5SBarry Smith if (size > 1) { 749806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords);CHKERRQ(ierr); 750806fa848SBarry Smith } else { 751c8b0795cSMark F. Adams coords = (PetscReal*)pc_gamg->data; 752a2f3521dSMark F. Adams data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols; 7532e68589bSMark F. Adams } 7542e68589bSMark F. Adams ierr = MatDestroy(&Gmat2);CHKERRQ(ierr); 7552e68589bSMark F. Adams 7562e68589bSMark F. Adams /* triangulate */ 7572e68589bSMark F. Adams if (dim == 2) { 758c8b0795cSMark F. Adams PetscReal metric,tm; 7590cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7600cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr); 7612e68589bSMark F. Adams #endif 762806fa848SBarry Smith ierr = triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric);CHKERRQ(ierr); 7630cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7640cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr); 7652e68589bSMark F. Adams #endif 7662e68589bSMark F. Adams ierr = PetscFree(crsGID);CHKERRQ(ierr); 7672e68589bSMark F. Adams 7682e68589bSMark F. Adams /* clean up and create coordinates for coarse grid (output) */ 769c5df96a5SBarry Smith if (size > 1) ierr = PetscFree(coords);CHKERRQ(ierr); 7702e68589bSMark F. Adams 771b2566f29SBarry Smith ierr = MPIU_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm);CHKERRQ(ierr); 772c8b0795cSMark F. Adams if (tm > 1.) { /* needs to be globalized - should not happen */ 773bf4339c2SBarry Smith ierr = PetscInfo1(pc," failed metric for coarse grid %e\n",(double)tm);CHKERRQ(ierr); 7742e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 775806fa848SBarry Smith } else if (metric > .0) { 776bf4339c2SBarry Smith ierr = PetscInfo1(pc,"worst metric for coarse grid = %e\n",(double)metric);CHKERRQ(ierr); 7772e68589bSMark F. Adams } 7783b4367a7SBarry Smith } else SETERRQ(comm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG"); 7792e68589bSMark F. Adams { /* create next coords - output */ 7802e68589bSMark F. Adams PetscReal *crs_crds; 781785e854fSJed Brown ierr = PetscMalloc1(dim*nLocalSelected, &crs_crds);CHKERRQ(ierr); 7822e68589bSMark F. Adams for (kk=0; kk<nLocalSelected; kk++) { /* grab local select nodes to promote - output */ 783b43b03e9SMark F. Adams PetscInt lid = clid_flid[kk]; 784c8b0795cSMark F. Adams for (jj=0; jj<dim; jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid]; 7852e68589bSMark F. Adams } 786c8b0795cSMark F. Adams 787c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 788c8b0795cSMark F. Adams pc_gamg->data = crs_crds; /* out */ 789c8b0795cSMark F. Adams pc_gamg->data_sz = dim*nLocalSelected; 7902e68589bSMark F. Adams } 791a2f3521dSMark F. Adams ierr = ISDestroy(&selected_2);CHKERRQ(ierr); 7922e68589bSMark F. Adams } 793a2f3521dSMark F. Adams 7942e68589bSMark F. Adams *a_P_out = Prol; /* out */ 795b43b03e9SMark F. Adams ierr = PetscFree(clid_flid);CHKERRQ(ierr); 7960cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr); 797c8b0795cSMark F. Adams PetscFunctionReturn(0); 798c8b0795cSMark F. Adams } 799c8b0795cSMark F. Adams 8009b8ffb57SJed Brown #undef __FUNCT__ 8019b8ffb57SJed Brown #define __FUNCT__ "PCDestroy_GAMG_GEO" 8029b8ffb57SJed Brown static PetscErrorCode PCDestroy_GAMG_GEO(PC pc) 8039b8ffb57SJed Brown { 8049b8ffb57SJed Brown PetscErrorCode ierr; 8059b8ffb57SJed Brown 8069b8ffb57SJed Brown PetscFunctionBegin; 8079b8ffb57SJed Brown ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr); 8089b8ffb57SJed Brown PetscFunctionReturn(0); 8099b8ffb57SJed Brown } 8109b8ffb57SJed Brown 811c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 812c8b0795cSMark F. Adams /* 813c8b0795cSMark F. Adams PCCreateGAMG_GEO 814c8b0795cSMark F. Adams 815c8b0795cSMark F. Adams Input Parameter: 816c8b0795cSMark F. Adams . pc - 817c8b0795cSMark F. Adams */ 818c8b0795cSMark F. Adams #undef __FUNCT__ 819c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_GEO" 820c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_GEO(PC pc) 821c8b0795cSMark F. Adams { 822c8b0795cSMark F. Adams PetscErrorCode ierr; 823c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 824c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 825c8b0795cSMark F. Adams 826c8b0795cSMark F. Adams PetscFunctionBegin; 8271ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO; 8289b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_GEO; 829c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 830c8b0795cSMark F. Adams 831c8b0795cSMark F. Adams /* set internal function pointers */ 832fd1112cbSBarry Smith pc_gamg->ops->graph = PCGAMGGraph_GEO; 833fd1112cbSBarry Smith pc_gamg->ops->coarsen = PCGAMGCoarsen_GEO; 8341ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_GEO; 835fd1112cbSBarry Smith pc_gamg->ops->optprolongator = NULL; 8361ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_GEO; 837c8b0795cSMark F. Adams 838bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_GEO);CHKERRQ(ierr); 8392e68589bSMark F. Adams PetscFunctionReturn(0); 8402e68589bSMark F. Adams } 841