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*/ 6b45d2f2cSJed Brown #include <petsc-private/kspimpl.h> 72e68589bSMark F. Adams 82e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE) 92e68589bSMark F. Adams #define REAL PetscReal 102e68589bSMark F. Adams #include <triangle.h> 112e68589bSMark F. Adams #endif 122e68589bSMark F. Adams 132e68589bSMark F. Adams #include <petscblaslapack.h> 142e68589bSMark F. Adams 15c8b0795cSMark F. Adams /* Private context for the GAMG preconditioner */ 16c8b0795cSMark F. Adams typedef struct { 17c8b0795cSMark F. Adams PetscInt lid; /* local vertex index */ 18c8b0795cSMark F. Adams PetscInt degree; /* vertex degree */ 19c8b0795cSMark F. Adams } GAMGNode; 202fa5cd67SKarl Rupp 210cbbd2e1SMark F. Adams int petsc_geo_mg_compare(const void *a, const void *b) 22c8b0795cSMark F. Adams { 23c8b0795cSMark F. Adams return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree); 24c8b0795cSMark F. Adams } 25c8b0795cSMark F. Adams 262e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 272e68589bSMark F. Adams /* 282e68589bSMark F. Adams PCSetCoordinates_GEO 292e68589bSMark F. Adams 302e68589bSMark F. Adams Input Parameter: 312e68589bSMark F. Adams . pc - the preconditioner context 322e68589bSMark F. Adams */ 332e68589bSMark F. Adams EXTERN_C_BEGIN 342e68589bSMark F. Adams #undef __FUNCT__ 352e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_GEO" 36302f38e8SMark F. Adams PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 372e68589bSMark F. Adams { 382e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 392e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 402e68589bSMark F. Adams PetscErrorCode ierr; 412e68589bSMark F. Adams PetscInt arrsz,bs,my0,kk,ii,nloc,Iend; 422e68589bSMark F. Adams Mat Amat = pc->pmat; 432e68589bSMark F. Adams 442e68589bSMark F. Adams PetscFunctionBegin; 452e68589bSMark F. Adams PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1); 462e68589bSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 47a2f3521dSMark F. Adams 482e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr); 492e68589bSMark F. Adams nloc = (Iend-my0)/bs; 50a2f3521dSMark F. Adams 51*ce94432eSBarry Smith if (nloc!=a_nloc) SETERRQ2(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Stokes not supported nloc = %d %d.",a_nloc,nloc); 52*ce94432eSBarry Smith if ((Iend-my0)%bs!=0) SETERRQ1(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 532e68589bSMark F. Adams 54c8b0795cSMark F. Adams pc_gamg->data_cell_rows = 1; 55*ce94432eSBarry Smith if (coords==0 && nloc > 0) SETERRQ(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'."); 56c8b0795cSMark F. Adams pc_gamg->data_cell_cols = ndm; /* coordinates */ 572e68589bSMark F. Adams 58c8b0795cSMark F. Adams arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 592e68589bSMark F. Adams 602e68589bSMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 612e68589bSMark F. Adams if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) { 622e68589bSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 632e68589bSMark F. Adams ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr); 642e68589bSMark F. Adams } 652e68589bSMark F. Adams for (kk=0; kk<arrsz; kk++) pc_gamg->data[kk] = -999.; 662e68589bSMark F. Adams pc_gamg->data[arrsz] = -99.; 672e68589bSMark F. Adams /* copy data in - column oriented */ 682e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) { 692e68589bSMark F. Adams for (ii = 0; ii < ndm; ii++) { 702e68589bSMark F. Adams pc_gamg->data[ii*nloc + kk] = coords[kk*ndm + ii]; 712e68589bSMark F. Adams } 722e68589bSMark F. Adams } 7371959b99SBarry 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]); 742e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 752e68589bSMark F. Adams PetscFunctionReturn(0); 762e68589bSMark F. Adams } 772e68589bSMark F. Adams EXTERN_C_END 782e68589bSMark F. Adams 792e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 802e68589bSMark F. Adams /* 812e68589bSMark F. Adams PCSetData_GEO 822e68589bSMark F. Adams 832e68589bSMark F. Adams Input Parameter: 842e68589bSMark F. Adams . pc - 852e68589bSMark F. Adams */ 862e68589bSMark F. Adams #undef __FUNCT__ 872e68589bSMark F. Adams #define __FUNCT__ "PCSetData_GEO" 88b8cd405aSMark F. Adams PetscErrorCode PCSetData_GEO(PC pc, Mat m) 892e68589bSMark F. Adams { 902e68589bSMark F. Adams PetscFunctionBegin; 91*ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"GEO MG needs coordinates"); 922e68589bSMark F. Adams } 932e68589bSMark F. Adams 942e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 952e68589bSMark F. Adams /* 962e68589bSMark F. Adams PCSetFromOptions_GEO 972e68589bSMark F. Adams 982e68589bSMark F. Adams Input Parameter: 992e68589bSMark F. Adams . pc - 1002e68589bSMark F. Adams */ 1012e68589bSMark F. Adams #undef __FUNCT__ 1022e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GEO" 1032e68589bSMark F. Adams PetscErrorCode PCSetFromOptions_GEO(PC pc) 1042e68589bSMark F. Adams { 1052e68589bSMark F. Adams PetscErrorCode ierr; 1062e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1072e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1082e68589bSMark F. Adams 1092e68589bSMark F. Adams PetscFunctionBegin; 1102e68589bSMark F. Adams ierr = PetscOptionsHead("GAMG-GEO options");CHKERRQ(ierr); 1112e68589bSMark F. Adams { 1122e68589bSMark F. Adams /* -pc_gamg_sa_nsmooths */ 1132e68589bSMark F. Adams /* pc_gamg_sa->smooths = 0; */ 1142e68589bSMark F. Adams /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */ 1152e68589bSMark F. Adams /* "smoothing steps for smoothed aggregation, usually 1 (0)", */ 1162e68589bSMark F. Adams /* "PCGAMGSetNSmooths_AGG", */ 1172e68589bSMark F. Adams /* pc_gamg_sa->smooths, */ 1182e68589bSMark F. Adams /* &pc_gamg_sa->smooths, */ 1192e68589bSMark F. Adams /* &flag); */ 1202e68589bSMark F. Adams /* CHKERRQ(ierr); */ 1212e68589bSMark F. Adams } 1222e68589bSMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1232e68589bSMark F. Adams 1242e68589bSMark F. Adams /* call base class */ 1252e68589bSMark F. Adams ierr = PCSetFromOptions_GAMG(pc);CHKERRQ(ierr); 1262e68589bSMark F. Adams 1272e68589bSMark F. Adams if (pc_gamg->verbose) { 128*ce94432eSBarry Smith PetscPrintf(PetscObjectComm((PetscObject)pc),"[%d]%s done\n",0,__FUNCT__); 1292e68589bSMark F. Adams } 1302e68589bSMark F. Adams PetscFunctionReturn(0); 1312e68589bSMark F. Adams } 1322e68589bSMark F. Adams 1332e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1342e68589bSMark F. Adams /* 1352e68589bSMark F. Adams triangulateAndFormProl 1362e68589bSMark F. Adams 1372e68589bSMark F. Adams Input Parameter: 1382e68589bSMark F. Adams . selected_2 - list of selected local ID, includes selected ghosts 139a2f3521dSMark F. Adams . data_stride - 140a2f3521dSMark F. Adams . coords[2*data_stride] - column vector of local coordinates w/ ghosts 141b43b03e9SMark F. Adams . nselected_1 - selected IDs that go with base (1) graph 1420cbbd2e1SMark F. Adams . clid_lid_1[nselected_1] - lids of selected (c) nodes ??????????? 1430cbbd2e1SMark F. Adams . agg_lists_1 - list of aggregates 1442e68589bSMark F. Adams . crsGID[selected.size()] - global index for prolongation operator 1452e68589bSMark F. Adams . bs - block size 1462e68589bSMark F. Adams Output Parameter: 1472e68589bSMark F. Adams . a_Prol - prolongation operator 1482e68589bSMark F. Adams . a_worst_best - measure of worst missed fine vertex, 0 is no misses 1492e68589bSMark F. Adams */ 1502e68589bSMark F. Adams #undef __FUNCT__ 1512e68589bSMark F. Adams #define __FUNCT__ "triangulateAndFormProl" 1520cbbd2e1SMark F. Adams static PetscErrorCode triangulateAndFormProl(IS selected_2, /* list of selected local ID, includes selected ghosts */ 153a2f3521dSMark F. Adams const PetscInt data_stride, 1542e68589bSMark F. Adams const PetscReal coords[], /* column vector of local coordinates w/ ghosts */ 155b43b03e9SMark F. Adams const PetscInt nselected_1, /* list of selected local ID, includes selected ghosts */ 156b43b03e9SMark F. Adams const PetscInt clid_lid_1[], 1570cbbd2e1SMark F. Adams const PetscCoarsenData *agg_lists_1, /* selected_1 vertices of aggregate unselected vertices */ 1582e68589bSMark F. Adams const PetscInt crsGID[], 1592e68589bSMark F. Adams const PetscInt bs, 1602e68589bSMark F. Adams Mat a_Prol, /* prolongation operator (output) */ 1611147fc2aSKarl Rupp PetscReal *a_worst_best) /* measure of worst missed fine vertex, 0 is no misses */ 1622e68589bSMark F. Adams { 1632e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE) 1642e68589bSMark F. Adams PetscErrorCode ierr; 165b43b03e9SMark F. Adams PetscInt jj,tid,tt,idx,nselected_2; 1662e68589bSMark F. Adams struct triangulateio in,mid; 1670cbbd2e1SMark F. Adams const PetscInt *selected_idx_2; 168c5df96a5SBarry Smith PetscMPIInt rank,size; 1692e68589bSMark F. Adams PetscInt Istart,Iend,nFineLoc,myFine0; 1702e68589bSMark F. Adams int kk,nPlotPts,sid; 171c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Prol)->comm; 172c8b0795cSMark F. Adams PetscReal tm; 173c8b0795cSMark F. Adams 1746e111a19SKarl Rupp PetscFunctionBegin; 175c5df96a5SBarry Smith ierr = MPI_Comm_rank(wcomm,&rank);CHKERRQ(ierr); 176c5df96a5SBarry Smith ierr = MPI_Comm_size(wcomm,&size);CHKERRQ(ierr); 177c8b0795cSMark F. Adams ierr = ISGetSize(selected_2, &nselected_2);CHKERRQ(ierr); 1782e68589bSMark F. Adams if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */ 1792e68589bSMark F. Adams *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */ 180806fa848SBarry Smith } else *a_worst_best = 0.0; 181c8b0795cSMark F. Adams ierr = MPI_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm);CHKERRQ(ierr); 182c8b0795cSMark F. Adams if (tm > 0.0) { 183c8b0795cSMark F. Adams *a_worst_best = 100.0; 1842e68589bSMark F. Adams PetscFunctionReturn(0); 1852e68589bSMark F. Adams } 1862e68589bSMark F. Adams ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr); 1872e68589bSMark F. Adams nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs; 1882e68589bSMark F. Adams nPlotPts = nFineLoc; /* locals */ 1892e68589bSMark F. Adams /* traingle */ 1902e68589bSMark F. Adams /* Define input points - in*/ 1912e68589bSMark F. Adams in.numberofpoints = nselected_2; 1922e68589bSMark F. Adams in.numberofpointattributes = 0; 1932e68589bSMark F. Adams /* get nselected points */ 1942e68589bSMark F. Adams ierr = PetscMalloc(2*(nselected_2)*sizeof(REAL), &in.pointlist);CHKERRQ(ierr); 1952e68589bSMark F. Adams ierr = ISGetIndices(selected_2, &selected_idx_2);CHKERRQ(ierr); 1962e68589bSMark F. Adams 1972e68589bSMark F. Adams for (kk=0,sid=0; kk<nselected_2; kk++,sid += 2) { 1982e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 1992e68589bSMark F. Adams in.pointlist[sid] = coords[lid]; 200a2f3521dSMark F. Adams in.pointlist[sid+1] = coords[data_stride + lid]; 2012e68589bSMark F. Adams if (lid>=nFineLoc) nPlotPts++; 2022e68589bSMark F. Adams } 20371959b99SBarry Smith if (sid != 2*nselected_2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != 2*nselected_2 %D",sid,nselected_2); 2042e68589bSMark F. Adams 2052e68589bSMark F. Adams in.numberofsegments = 0; 2062e68589bSMark F. Adams in.numberofedges = 0; 2072e68589bSMark F. Adams in.numberofholes = 0; 2082e68589bSMark F. Adams in.numberofregions = 0; 2092e68589bSMark F. Adams in.trianglelist = 0; 2102e68589bSMark F. Adams in.segmentmarkerlist = 0; 2112e68589bSMark F. Adams in.pointattributelist = 0; 2122e68589bSMark F. Adams in.pointmarkerlist = 0; 2132e68589bSMark F. Adams in.triangleattributelist = 0; 2142e68589bSMark F. Adams in.trianglearealist = 0; 2152e68589bSMark F. Adams in.segmentlist = 0; 2162e68589bSMark F. Adams in.holelist = 0; 2172e68589bSMark F. Adams in.regionlist = 0; 2182e68589bSMark F. Adams in.edgelist = 0; 2192e68589bSMark F. Adams in.edgemarkerlist = 0; 2202e68589bSMark F. Adams in.normlist = 0; 2212fa5cd67SKarl Rupp 2222e68589bSMark F. Adams /* triangulate */ 2232e68589bSMark F. Adams mid.pointlist = 0; /* Not needed if -N switch used. */ 2242e68589bSMark F. Adams /* Not needed if -N switch used or number of point attributes is zero: */ 2252e68589bSMark F. Adams mid.pointattributelist = 0; 2262e68589bSMark F. Adams mid.pointmarkerlist = 0; /* Not needed if -N or -B switch used. */ 2272e68589bSMark F. Adams mid.trianglelist = 0; /* Not needed if -E switch used. */ 2282e68589bSMark F. Adams /* Not needed if -E switch used or number of triangle attributes is zero: */ 2292e68589bSMark F. Adams mid.triangleattributelist = 0; 2302e68589bSMark F. Adams mid.neighborlist = 0; /* Needed only if -n switch used. */ 2312e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P not used: */ 2322e68589bSMark F. Adams mid.segmentlist = 0; 2332e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P and -B not used: */ 2342e68589bSMark F. Adams mid.segmentmarkerlist = 0; 2352e68589bSMark F. Adams mid.edgelist = 0; /* Needed only if -e switch used. */ 2362e68589bSMark F. Adams mid.edgemarkerlist = 0; /* Needed if -e used and -B not used. */ 2372e68589bSMark F. Adams mid.numberoftriangles = 0; 2382e68589bSMark F. Adams 2392e68589bSMark F. Adams /* Triangulate the points. Switches are chosen to read and write a */ 2402e68589bSMark F. Adams /* PSLG (p), preserve the convex hull (c), number everything from */ 2412e68589bSMark F. Adams /* zero (z), assign a regional attribute to each element (A), and */ 2422e68589bSMark F. Adams /* produce an edge list (e), a Voronoi diagram (v), and a triangle */ 2432e68589bSMark F. Adams /* neighbor list (n). */ 2442e68589bSMark F. Adams if (nselected_2 != 0) { /* inactive processor */ 2452e68589bSMark F. Adams char args[] = "npczQ"; /* c is needed ? */ 2462e68589bSMark F. Adams triangulate(args, &in, &mid, (struct triangulateio*) NULL); 2472e68589bSMark F. Adams /* output .poly files for 'showme' */ 2482e68589bSMark F. Adams if (!PETSC_TRUE) { 2492e68589bSMark F. Adams static int level = 1; 2502e68589bSMark F. Adams FILE *file; char fname[32]; 2512e68589bSMark F. Adams 252c5df96a5SBarry Smith sprintf(fname,"C%d_%d.poly",level,rank); file = fopen(fname, "w"); 2532e68589bSMark F. Adams /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/ 2542e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); 2552e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2562e68589bSMark F. Adams for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid += 2) { 2572e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2582e68589bSMark F. Adams } 2592e68589bSMark F. Adams /*One line: <# of segments> <# of boundary markers (0 or 1)> */ 2602e68589bSMark F. Adams fprintf(file, "%d %d\n",0,0); 2612e68589bSMark F. Adams /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */ 2622e68589bSMark F. Adams /* One line: <# of holes> */ 2632e68589bSMark F. Adams fprintf(file, "%d\n",0); 2642e68589bSMark F. Adams /* Following lines: <hole #> <x> <y> */ 2652e68589bSMark F. Adams /* Optional line: <# of regional attributes and/or area constraints> */ 2662e68589bSMark F. Adams /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */ 2672e68589bSMark F. Adams fclose(file); 2682e68589bSMark F. Adams 2692e68589bSMark F. Adams /* elems */ 270c5df96a5SBarry Smith sprintf(fname,"C%d_%d.ele",level,rank); file = fopen(fname, "w"); 2712e68589bSMark F. Adams /* First line: <# of triangles> <nodes per triangle> <# of attributes> */ 2722e68589bSMark F. Adams fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0); 2732e68589bSMark F. Adams /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */ 2742e68589bSMark F. Adams for (kk=0,sid=0; kk<mid.numberoftriangles; kk++,sid += 3) { 2752e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]); 2762e68589bSMark F. Adams } 2772e68589bSMark F. Adams fclose(file); 2782e68589bSMark F. Adams 279c5df96a5SBarry Smith sprintf(fname,"C%d_%d.node",level,rank); file = fopen(fname, "w"); 2802e68589bSMark F. Adams /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */ 2812e68589bSMark F. Adams /* fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); */ 2822e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",nPlotPts,2,0,0); 2832e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2842e68589bSMark F. Adams for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid+=2) { 2852e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2862e68589bSMark F. Adams } 2872e68589bSMark F. Adams 2882e68589bSMark F. Adams sid /= 2; 2892e68589bSMark F. Adams for (jj=0; jj<nFineLoc; jj++) { 2902e68589bSMark F. Adams PetscBool sel = PETSC_TRUE; 2912e68589bSMark F. Adams for (kk=0; kk<nselected_2 && sel; kk++) { 2922e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 2932e68589bSMark F. Adams if (lid == jj) sel = PETSC_FALSE; 2942e68589bSMark F. Adams } 2952fa5cd67SKarl Rupp if (sel) fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]); 2962e68589bSMark F. Adams } 2972e68589bSMark F. Adams fclose(file); 29871959b99SBarry Smith if (sid != nPlotPts) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != nPlotPts %D",sid,nPlotPts); 2992e68589bSMark F. Adams level++; 3002e68589bSMark F. Adams } 3012e68589bSMark F. Adams } 3020cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3030cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 3042e68589bSMark F. Adams #endif 3052e68589bSMark F. Adams { /* form P - setup some maps */ 3060cbbd2e1SMark F. Adams PetscInt clid,mm,*nTri,*node_tri; 3072e68589bSMark F. Adams 3082e68589bSMark F. Adams ierr = PetscMalloc(nselected_2*sizeof(PetscInt), &node_tri);CHKERRQ(ierr); 3092e68589bSMark F. Adams ierr = PetscMalloc(nselected_2*sizeof(PetscInt), &nTri);CHKERRQ(ierr); 3102e68589bSMark F. Adams 3112e68589bSMark F. Adams /* need list of triangles on node */ 3122e68589bSMark F. Adams for (kk=0; kk<nselected_2; kk++) nTri[kk] = 0; 3132e68589bSMark F. Adams for (tid=0,kk=0; tid<mid.numberoftriangles; tid++) { 3142e68589bSMark F. Adams for (jj=0; jj<3; jj++) { 3152e68589bSMark F. Adams PetscInt cid = mid.trianglelist[kk++]; 3162e68589bSMark F. Adams if (nTri[cid] == 0) node_tri[cid] = tid; 3172e68589bSMark F. Adams nTri[cid]++; 3182e68589bSMark F. Adams } 3192e68589bSMark F. Adams } 3202e68589bSMark F. Adams #define EPS 1.e-12 3212e68589bSMark F. Adams /* find points and set prolongation */ 3220cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nFineLoc; mm++) { 323e78576d6SMark F. Adams PetscBool ise; 324e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise);CHKERRQ(ierr); 325e78576d6SMark F. Adams if (!ise) { 3260cbbd2e1SMark F. Adams const PetscInt lid = mm; 327c41eeb4cSKarl Rupp /* for (clid_iterator=0;clid_iterator<nselected_1;clid_iterator++) { */ 3282e68589bSMark F. Adams PetscScalar AA[3][3]; 3292e68589bSMark F. Adams PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO; 33041b27cdeSMark F. Adams PetscCDPos pos; 331e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos);CHKERRQ(ierr); 332e78576d6SMark F. Adams while (pos) { 333e78576d6SMark F. Adams PetscInt flid; 334ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &flid);CHKERRQ(ierr); 335e78576d6SMark F. Adams ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos);CHKERRQ(ierr); 3360cbbd2e1SMark F. Adams 3372e68589bSMark F. Adams if (flid < nFineLoc) { /* could be a ghost */ 3382e68589bSMark F. Adams PetscInt bestTID = -1; PetscReal best_alpha = 1.e10; 3392e68589bSMark F. Adams const PetscInt fgid = flid + myFine0; 3402e68589bSMark F. Adams /* compute shape function for gid */ 341a2f3521dSMark F. Adams const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0}; 3422e68589bSMark F. Adams PetscBool haveit =PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3]; 3432fa5cd67SKarl Rupp 3442e68589bSMark F. Adams /* look for it */ 3450cbbd2e1SMark F. Adams for (tid = node_tri[clid], jj=0; 3462e68589bSMark F. Adams jj < 5 && !haveit && tid != -1; 3472e68589bSMark F. Adams jj++) { 3482e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3492e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3502e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 351a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3522e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3532e68589bSMark F. Adams } 3542e68589bSMark F. Adams 3552e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = (PetscScalar)fcoord[tt]; 3562e68589bSMark F. Adams 3572e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 358f75e95b9SBarry Smith PetscStackCall("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3592e68589bSMark F. Adams { 3602e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal lowest=1.e10; 3612e68589bSMark F. Adams for (tt = 0, idx = 0; tt < 3; tt++) { 3622e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE; 3632e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) < lowest) { 3642e68589bSMark F. Adams lowest = PetscRealPart(alpha[tt]); 3652e68589bSMark F. Adams idx = tt; 3662e68589bSMark F. Adams } 3672e68589bSMark F. Adams } 3682e68589bSMark F. Adams haveit = have; 3692e68589bSMark F. Adams } 3702e68589bSMark F. Adams tid = mid.neighborlist[3*tid + idx]; 3712e68589bSMark F. Adams } 3722e68589bSMark F. Adams 3732e68589bSMark F. Adams if (!haveit) { 3742e68589bSMark F. Adams /* brute force */ 3752e68589bSMark F. Adams for (tid=0; tid<mid.numberoftriangles && !haveit; tid++) { 3762e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 3772e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3782e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 379a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 3802e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3812e68589bSMark F. Adams } 3822e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt]; 3832e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 384f75e95b9SBarry Smith PetscStackCall("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 3852e68589bSMark F. Adams { 3862e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal worst=0.0, v; 3872e68589bSMark F. Adams for (tt=0; tt<3 && have; tt++) { 3882e68589bSMark F. Adams if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE; 3892e68589bSMark F. Adams if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v; 3902e68589bSMark F. Adams } 3912e68589bSMark F. Adams if (worst < best_alpha) { 3922e68589bSMark F. Adams best_alpha = worst; bestTID = tid; 3932e68589bSMark F. Adams } 3942e68589bSMark F. Adams haveit = have; 3952e68589bSMark F. Adams } 3962e68589bSMark F. Adams } 3972e68589bSMark F. Adams } 3982e68589bSMark F. Adams if (!haveit) { 3992e68589bSMark F. Adams if (best_alpha > *a_worst_best) *a_worst_best = best_alpha; 4002e68589bSMark F. Adams /* use best one */ 4012e68589bSMark F. Adams for (tt=0; tt<3; tt++) { 4022e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*bestTID + tt]; 4032e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 404a2f3521dSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0; 4052e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 4062e68589bSMark F. Adams } 4072e68589bSMark F. Adams for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt]; 4082e68589bSMark F. Adams /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */ 409f75e95b9SBarry Smith PetscStackCall("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO)); 4102e68589bSMark F. Adams } 4112e68589bSMark F. Adams 4122e68589bSMark F. Adams /* put in row of P */ 4132e68589bSMark F. Adams for (idx=0; idx<3; idx++) { 4142e68589bSMark F. Adams PetscScalar shp = alpha[idx]; 4152e68589bSMark F. Adams if (PetscAbs(PetscRealPart(shp)) > 1.e-6) { 4162e68589bSMark F. Adams PetscInt cgid = crsGID[clids[idx]]; 4172e68589bSMark F. Adams PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */ 4182e68589bSMark F. Adams for (tt=0; tt < bs; tt++, ii++, jj++) { 4192e68589bSMark F. Adams ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES);CHKERRQ(ierr); 4202e68589bSMark F. Adams } 4212e68589bSMark F. Adams } 4222e68589bSMark F. Adams } 4232e68589bSMark F. Adams } 4240cbbd2e1SMark F. Adams } /* aggregates iterations */ 4250cbbd2e1SMark F. Adams clid++; 4260cbbd2e1SMark F. Adams } /* a coarse agg */ 4270cbbd2e1SMark F. Adams } /* for all fine nodes */ 4280cbbd2e1SMark F. Adams 4292e68589bSMark F. Adams ierr = ISRestoreIndices(selected_2, &selected_idx_2);CHKERRQ(ierr); 4302e68589bSMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4312e68589bSMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4322e68589bSMark F. Adams 4332e68589bSMark F. Adams ierr = PetscFree(node_tri);CHKERRQ(ierr); 4342e68589bSMark F. Adams ierr = PetscFree(nTri);CHKERRQ(ierr); 4352e68589bSMark F. Adams } 4360cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4370cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 4382e68589bSMark F. Adams #endif 4392e68589bSMark F. Adams free(mid.trianglelist); 4402e68589bSMark F. Adams free(mid.neighborlist); 4412e68589bSMark F. Adams ierr = PetscFree(in.pointlist);CHKERRQ(ierr); 4422e68589bSMark F. Adams PetscFunctionReturn(0); 4432e68589bSMark F. Adams #else 444c666adbfSMark F. Adams SETERRQ(((PetscObject)a_Prol)->comm,PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG"); 4452e68589bSMark F. Adams #endif 4462e68589bSMark F. Adams } 4472e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 4482e68589bSMark F. Adams /* 4492e68589bSMark F. Adams getGIDsOnSquareGraph - square graph, get 4502e68589bSMark F. Adams 4512e68589bSMark F. Adams Input Parameter: 4520cbbd2e1SMark F. Adams . nselected_1 - selected local indices (includes ghosts in input Gmat1) 453b43b03e9SMark F. Adams . clid_lid_1 - [nselected_1] lids of selected nodes 4542e68589bSMark F. Adams . Gmat1 - graph that goes with 'selected_1' 4552e68589bSMark F. Adams Output Parameter: 4562e68589bSMark F. Adams . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2) 4572e68589bSMark F. Adams . a_Gmat_2 - graph that is squared of 'Gmat_1' 4582e68589bSMark F. Adams . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes 4592e68589bSMark F. Adams */ 4602e68589bSMark F. Adams #undef __FUNCT__ 4612e68589bSMark F. Adams #define __FUNCT__ "getGIDsOnSquareGraph" 4622fa5cd67SKarl Rupp static PetscErrorCode getGIDsOnSquareGraph(const PetscInt nselected_1,const PetscInt clid_lid_1[],const Mat Gmat1,IS *a_selected_2,Mat *a_Gmat_2,PetscInt **a_crsGID) 4632e68589bSMark F. Adams { 4642e68589bSMark F. Adams PetscErrorCode ierr; 465c5df96a5SBarry Smith PetscMPIInt rank,size; 466b43b03e9SMark F. Adams PetscInt *crsGID, kk,my0,Iend,nloc; 4672e68589bSMark F. Adams MPI_Comm wcomm = ((PetscObject)Gmat1)->comm; 4682e68589bSMark F. Adams 4692e68589bSMark F. Adams PetscFunctionBegin; 470c5df96a5SBarry Smith ierr = MPI_Comm_rank(wcomm,&rank);CHKERRQ(ierr); 471c5df96a5SBarry Smith ierr = MPI_Comm_size(wcomm,&size);CHKERRQ(ierr); 4722e68589bSMark F. Adams ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend);CHKERRQ(ierr); /* AIJ */ 4732e68589bSMark F. Adams nloc = Iend - my0; /* this does not change */ 4742e68589bSMark F. Adams 475c5df96a5SBarry Smith if (size == 1) { /* not much to do in serial */ 476b43b03e9SMark F. Adams ierr = PetscMalloc(nselected_1*sizeof(PetscInt), &crsGID);CHKERRQ(ierr); 477b43b03e9SMark F. Adams for (kk=0; kk<nselected_1; kk++) crsGID[kk] = kk; 4782e68589bSMark F. Adams *a_Gmat_2 = 0; 479806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);CHKERRQ(ierr); 480806fa848SBarry Smith } else { 481b43b03e9SMark F. Adams PetscInt idx,num_fine_ghosts,num_crs_ghost,myCrs0; 4822e68589bSMark F. Adams Mat_MPIAIJ *mpimat2; 4832e68589bSMark F. Adams Mat Gmat2; 4842e68589bSMark F. Adams Vec locState; 4852e68589bSMark F. Adams PetscScalar *cpcol_state; 4862e68589bSMark F. Adams 4872e68589bSMark F. Adams /* scan my coarse zero gid, set 'lid_state' with coarse GID */ 488b43b03e9SMark F. Adams kk = nselected_1; 489b43b03e9SMark F. Adams MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPIU_SUM, wcomm); 490b43b03e9SMark F. Adams myCrs0 -= nselected_1; 4912e68589bSMark F. Adams 492b43b03e9SMark F. Adams if (a_Gmat_2) { /* output */ 4932e68589bSMark F. Adams /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */ 4942e68589bSMark F. Adams ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr); 4952e68589bSMark F. Adams *a_Gmat_2 = Gmat2; /* output */ 496806fa848SBarry Smith } else Gmat2 = Gmat1; /* use local to get crsGIDs at least */ 4972e68589bSMark F. Adams /* get coarse grid GIDS for selected (locals and ghosts) */ 4982e68589bSMark F. Adams mpimat2 = (Mat_MPIAIJ*)Gmat2->data; 4992e68589bSMark F. Adams ierr = MatGetVecs(Gmat2, &locState, 0);CHKERRQ(ierr); 5002e68589bSMark F. Adams ierr = VecSet(locState, (PetscScalar)(PetscReal)(-1));CHKERRQ(ierr); /* set with UNKNOWN state */ 501b43b03e9SMark F. Adams for (kk=0; kk<nselected_1; kk++) { 502b43b03e9SMark F. Adams PetscInt fgid = clid_lid_1[kk] + my0; 5032e68589bSMark F. Adams PetscScalar v = (PetscScalar)(kk+myCrs0); 5042e68589bSMark F. Adams ierr = VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES);CHKERRQ(ierr); /* set with PID */ 5052e68589bSMark F. Adams } 5062e68589bSMark F. Adams ierr = VecAssemblyBegin(locState);CHKERRQ(ierr); 5072e68589bSMark F. Adams ierr = VecAssemblyEnd(locState);CHKERRQ(ierr); 5082e68589bSMark F. Adams ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5092e68589bSMark F. Adams ierr = VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5102e68589bSMark F. Adams ierr = VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts);CHKERRQ(ierr); 5112e68589bSMark F. Adams ierr = VecGetArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr); 5122e68589bSMark F. Adams for (kk=0,num_crs_ghost=0; kk<num_fine_ghosts; kk++) { 5132e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++; 5142e68589bSMark F. Adams } 515b43b03e9SMark F. Adams ierr = PetscMalloc((nselected_1+num_crs_ghost)*sizeof(PetscInt), &crsGID);CHKERRQ(ierr); /* output */ 5162e68589bSMark F. Adams { 5172e68589bSMark F. Adams PetscInt *selected_set; 518b43b03e9SMark F. Adams ierr = PetscMalloc((nselected_1+num_crs_ghost)*sizeof(PetscInt), &selected_set);CHKERRQ(ierr); 5192e68589bSMark F. Adams /* do ghost of 'crsGID' */ 520b43b03e9SMark F. Adams for (kk=0,idx=nselected_1; kk<num_fine_ghosts; kk++) { 5212e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) { 5222e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5232e68589bSMark F. Adams selected_set[idx] = nloc + kk; 5242e68589bSMark F. Adams crsGID[idx++] = cgid; 5252e68589bSMark F. Adams } 5262e68589bSMark F. Adams } 52771959b99SBarry 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); 5282e68589bSMark F. Adams ierr = VecRestoreArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr); 5292e68589bSMark F. Adams /* do locals in 'crsGID' */ 5302e68589bSMark F. Adams ierr = VecGetArray(locState, &cpcol_state);CHKERRQ(ierr); 5312e68589bSMark F. Adams for (kk=0,idx=0; kk<nloc; kk++) { 5322e68589bSMark F. Adams if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) { 5332e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5342e68589bSMark F. Adams selected_set[idx] = kk; 5352e68589bSMark F. Adams crsGID[idx++] = cgid; 5362e68589bSMark F. Adams } 5372e68589bSMark F. Adams } 53871959b99SBarry Smith if (idx != nselected_1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != nselected_1 %D",idx,nselected_1); 5392e68589bSMark F. Adams ierr = VecRestoreArray(locState, &cpcol_state);CHKERRQ(ierr); 5402e68589bSMark F. Adams 5412e68589bSMark F. Adams if (a_selected_2 != 0) { /* output */ 542806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);CHKERRQ(ierr); 543806fa848SBarry Smith } else { 5442e68589bSMark F. Adams ierr = PetscFree(selected_set);CHKERRQ(ierr); 5452e68589bSMark F. Adams } 5460cbbd2e1SMark F. Adams } 5472e68589bSMark F. Adams ierr = VecDestroy(&locState);CHKERRQ(ierr); 5482e68589bSMark F. Adams } 5492e68589bSMark F. Adams *a_crsGID = crsGID; /* output */ 5502e68589bSMark F. Adams PetscFunctionReturn(0); 5512e68589bSMark F. Adams } 5522e68589bSMark F. Adams 5532e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 5542e68589bSMark F. Adams /* 555c8b0795cSMark F. Adams PCGAMGgraph_GEO 5562e68589bSMark F. Adams 5572e68589bSMark F. Adams Input Parameter: 5582e68589bSMark F. Adams . pc - this 5592e68589bSMark F. Adams . Amat - matrix on this fine level 5602e68589bSMark F. Adams Output Parameter: 561c8b0795cSMark F. Adams . a_Gmat 5622e68589bSMark F. Adams */ 5632e68589bSMark F. Adams #undef __FUNCT__ 564c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_GEO" 5652fa5cd67SKarl Rupp PetscErrorCode PCGAMGgraph_GEO(PC pc,const Mat Amat,Mat *a_Gmat) 566c8b0795cSMark F. Adams { 567c8b0795cSMark F. Adams PetscErrorCode ierr; 568c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 569c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 570c8b0795cSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 571c8b0795cSMark F. Adams const PetscReal vfilter = pc_gamg->threshold; 572c5df96a5SBarry Smith PetscMPIInt rank,size; 573*ce94432eSBarry Smith MPI_Comm wcomm; 574c8b0795cSMark F. Adams Mat Gmat; 5750cbbd2e1SMark F. Adams PetscBool set,flg,symm; 5766e111a19SKarl Rupp 577c8b0795cSMark F. Adams PetscFunctionBegin; 578*ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&wcomm);CHKERRQ(ierr); 5790cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 5800cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr); 5810cbbd2e1SMark F. Adams #endif 582c5df96a5SBarry Smith ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr); 583c5df96a5SBarry Smith ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr); 584c8b0795cSMark F. Adams 5850cbbd2e1SMark F. Adams ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); 586263489e9SJed Brown symm = (PetscBool)!(set && flg); 5870cbbd2e1SMark F. Adams 5882d7fac45SMark F. Adams ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr); 5892d7fac45SMark F. Adams ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm, verbose);CHKERRQ(ierr); 590c8b0795cSMark F. Adams 591c8b0795cSMark F. Adams *a_Gmat = Gmat; 5920cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 5930cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr); 5940cbbd2e1SMark F. Adams #endif 595c8b0795cSMark F. Adams PetscFunctionReturn(0); 596c8b0795cSMark F. Adams } 597c8b0795cSMark F. Adams 598c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 599c8b0795cSMark F. Adams /* 600c8b0795cSMark F. Adams PCGAMGcoarsen_GEO 601c8b0795cSMark F. Adams 602c8b0795cSMark F. Adams Input Parameter: 603e0940f08SMark F. Adams . a_pc - this 604e0940f08SMark F. Adams . a_Gmat - graph 605c8b0795cSMark F. Adams Output Parameter: 606c8b0795cSMark F. Adams . a_llist_parent - linked list from selected indices for data locality only 607c8b0795cSMark F. Adams */ 608c8b0795cSMark F. Adams #undef __FUNCT__ 609c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGcoarsen_GEO" 6102fa5cd67SKarl Rupp PetscErrorCode PCGAMGcoarsen_GEO(PC a_pc,Mat *a_Gmat,PetscCoarsenData **a_llist_parent) 611c8b0795cSMark F. Adams { 612c8b0795cSMark F. Adams PetscErrorCode ierr; 613e0940f08SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 614c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 615c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,kk,Ii,ncols; 616c5df96a5SBarry Smith PetscMPIInt rank,size; 6170cbbd2e1SMark F. Adams IS perm; 618c8b0795cSMark F. Adams GAMGNode *gnodes; 619c8b0795cSMark F. Adams PetscInt *permute; 620e0940f08SMark F. Adams Mat Gmat = *a_Gmat; 621*ce94432eSBarry Smith MPI_Comm wcomm; 622b43b03e9SMark F. Adams MatCoarsen crs; 623c8b0795cSMark F. Adams 624c8b0795cSMark F. Adams PetscFunctionBegin; 625*ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)a_pc,&wcomm);CHKERRQ(ierr); 6260cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 6270cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr); 6280cbbd2e1SMark F. Adams #endif 629c5df96a5SBarry Smith ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr); 630c5df96a5SBarry Smith ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr); 631c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr); 632c8b0795cSMark F. Adams nloc = (Iend-Istart); 633c8b0795cSMark F. Adams 634c8b0795cSMark F. Adams /* create random permutation with sort for geo-mg */ 635c8b0795cSMark F. Adams ierr = PetscMalloc(nloc*sizeof(GAMGNode), &gnodes);CHKERRQ(ierr); 636c8b0795cSMark F. Adams ierr = PetscMalloc(nloc*sizeof(PetscInt), &permute);CHKERRQ(ierr); 637c8b0795cSMark F. Adams 638c8b0795cSMark F. Adams for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */ 639c8b0795cSMark F. Adams ierr = MatGetRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr); 640c8b0795cSMark F. Adams { 641c8b0795cSMark F. Adams PetscInt lid = Ii - Istart; 642c8b0795cSMark F. Adams gnodes[lid].lid = lid; 643c8b0795cSMark F. Adams gnodes[lid].degree = ncols; 644c8b0795cSMark F. Adams } 645c8b0795cSMark F. Adams ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr); 646c8b0795cSMark F. Adams } 647c8b0795cSMark F. Adams /* randomize */ 648c8b0795cSMark F. Adams srand(1); /* make deterministic */ 649c8b0795cSMark F. Adams if (PETSC_TRUE) { 650c8b0795cSMark F. Adams PetscBool *bIndexSet; 651c8b0795cSMark F. Adams ierr = PetscMalloc(nloc*sizeof(PetscBool), &bIndexSet);CHKERRQ(ierr); 652c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) bIndexSet[Ii] = PETSC_FALSE; 6532fa5cd67SKarl Rupp for (Ii = 0; Ii < nloc; Ii++) { 654c8b0795cSMark F. Adams PetscInt iSwapIndex = rand()%nloc; 6552fa5cd67SKarl Rupp if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 656c8b0795cSMark F. Adams GAMGNode iTemp = gnodes[iSwapIndex]; 657c8b0795cSMark F. Adams gnodes[iSwapIndex] = gnodes[Ii]; 658c8b0795cSMark F. Adams gnodes[Ii] = iTemp; 659c8b0795cSMark F. Adams bIndexSet[Ii] = PETSC_TRUE; 660c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 661c8b0795cSMark F. Adams } 662c8b0795cSMark F. Adams } 663c8b0795cSMark F. Adams ierr = PetscFree(bIndexSet);CHKERRQ(ierr); 664c8b0795cSMark F. Adams } 665c8b0795cSMark F. Adams /* only sort locals */ 6660cbbd2e1SMark F. Adams qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare); 667c8b0795cSMark F. Adams /* create IS of permutation */ 6682fa5cd67SKarl Rupp for (kk=0; kk<nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */ 669806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr); 670c8b0795cSMark F. Adams 671c8b0795cSMark F. Adams ierr = PetscFree(gnodes);CHKERRQ(ierr); 672c8b0795cSMark F. Adams 673c8b0795cSMark F. Adams /* get MIS aggs */ 674b43b03e9SMark F. Adams 675b43b03e9SMark F. Adams ierr = MatCoarsenCreate(wcomm, &crs);CHKERRQ(ierr); 676b43b03e9SMark F. Adams ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); 677b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr); 678b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr); 679b43b03e9SMark F. Adams ierr = MatCoarsenSetVerbose(crs, pc_gamg->verbose);CHKERRQ(ierr); 680b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs(crs, PETSC_FALSE);CHKERRQ(ierr); 681b43b03e9SMark F. Adams ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 6820cbbd2e1SMark F. Adams ierr = MatCoarsenGetData(crs, a_llist_parent);CHKERRQ(ierr); 683b43b03e9SMark F. Adams ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 684c8b0795cSMark F. Adams 685c8b0795cSMark F. Adams ierr = ISDestroy(&perm);CHKERRQ(ierr); 6860cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 6870cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr); 6880cbbd2e1SMark F. Adams #endif 689c8b0795cSMark F. Adams PetscFunctionReturn(0); 690c8b0795cSMark F. Adams } 691c8b0795cSMark F. Adams 692c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 693c8b0795cSMark F. Adams /* 6940cbbd2e1SMark F. Adams PCGAMGProlongator_GEO 695c8b0795cSMark F. Adams 696c8b0795cSMark F. Adams Input Parameter: 697c8b0795cSMark F. Adams . pc - this 698c8b0795cSMark F. Adams . Amat - matrix on this fine level 699c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 7000cbbd2e1SMark F. Adams . selected_1 - [nselected] 7010cbbd2e1SMark F. Adams . agg_lists - [nselected] 702c8b0795cSMark F. Adams Output Parameter: 703c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 704c8b0795cSMark F. Adams */ 705c8b0795cSMark F. Adams #undef __FUNCT__ 7060cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_GEO" 7072fa5cd67SKarl Rupp PetscErrorCode PCGAMGProlongator_GEO(PC pc,const Mat Amat,const Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 7082e68589bSMark F. Adams { 7092e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 7102e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 7112e68589bSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 712c8b0795cSMark F. Adams const PetscInt dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols; 7132e68589bSMark F. Adams PetscErrorCode ierr; 714b43b03e9SMark F. Adams PetscInt Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid; 715c8b0795cSMark F. Adams Mat Prol; 716c5df96a5SBarry Smith PetscMPIInt rank, size; 717*ce94432eSBarry Smith MPI_Comm wcomm; 7180cbbd2e1SMark F. Adams IS selected_2,selected_1; 7192e68589bSMark F. Adams const PetscInt *selected_idx; 7202e68589bSMark F. Adams 7212e68589bSMark F. Adams PetscFunctionBegin; 722*ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&wcomm);CHKERRQ(ierr); 7230cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 7240cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr); 7250cbbd2e1SMark F. Adams #endif 726c5df96a5SBarry Smith ierr = MPI_Comm_rank(wcomm,&rank);CHKERRQ(ierr); 727c5df96a5SBarry Smith ierr = MPI_Comm_size(wcomm,&size);CHKERRQ(ierr); 7282e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 7292e68589bSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 73071959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 73171959b99SBarry Smith if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) % bs %D",Iend,Istart,bs); 7322e68589bSMark F. Adams 7332e68589bSMark F. Adams /* get 'nLocalSelected' */ 73441b27cdeSMark F. Adams ierr = PetscCDGetMIS(agg_lists, &selected_1);CHKERRQ(ierr); 735b43b03e9SMark F. Adams ierr = ISGetSize(selected_1, &jj);CHKERRQ(ierr); 736b43b03e9SMark F. Adams ierr = PetscMalloc(jj*sizeof(PetscInt), &clid_flid);CHKERRQ(ierr); 7372e68589bSMark F. Adams ierr = ISGetIndices(selected_1, &selected_idx);CHKERRQ(ierr); 738c8b0795cSMark F. Adams for (kk=0,nLocalSelected=0; kk<jj; kk++) { 7392e68589bSMark F. Adams PetscInt lid = selected_idx[kk]; 740b43b03e9SMark F. Adams if (lid<nloc) { 7410cbbd2e1SMark F. Adams ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr); 7422fa5cd67SKarl Rupp if (ncols>1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */ 7430cbbd2e1SMark F. Adams ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr); 744b43b03e9SMark F. Adams } 7452e68589bSMark F. Adams } 7462e68589bSMark F. Adams ierr = ISRestoreIndices(selected_1, &selected_idx);CHKERRQ(ierr); 747a2f3521dSMark F. Adams ierr = ISDestroy(&selected_1);CHKERRQ(ierr); /* this is selected_1 in serial */ 7482e68589bSMark F. Adams 7492e68589bSMark F. Adams /* create prolongator, create P matrix */ 750a2f3521dSMark F. Adams ierr = MatCreate(wcomm, &Prol);CHKERRQ(ierr); 751806fa848SBarry Smith ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 752a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Prol, bs, bs);CHKERRQ(ierr); 753a2f3521dSMark F. Adams ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr); 7540298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(Prol,3*data_cols,NULL);CHKERRQ(ierr); 7550298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(Prol,3*data_cols,NULL,3*data_cols,NULL);CHKERRQ(ierr); 756a2f3521dSMark F. Adams /* ierr = MatCreateAIJ(wcomm, */ 757a2f3521dSMark F. Adams /* nloc*bs, nLocalSelected*bs, */ 758a2f3521dSMark F. Adams /* PETSC_DETERMINE, PETSC_DETERMINE, */ 7590298fd71SBarry Smith /* 3*data_cols, NULL, */ 7600298fd71SBarry Smith /* 3*data_cols, NULL, */ 761a2f3521dSMark F. Adams /* &Prol); */ 762a2f3521dSMark F. Adams /* CHKERRQ(ierr); */ 7632e68589bSMark F. Adams 764c8b0795cSMark F. Adams /* can get all points "removed" - but not on geomg */ 7652e68589bSMark F. Adams ierr = MatGetSize(Prol, &kk, &jj);CHKERRQ(ierr); 7662e68589bSMark F. Adams if (jj==0) { 7672fa5cd67SKarl Rupp if (verbose) PetscPrintf(wcomm,"[%d]%s ERROE: no selected points on coarse grid\n",rank,__FUNCT__); 768b43b03e9SMark F. Adams ierr = PetscFree(clid_flid);CHKERRQ(ierr); 7692e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 7700298fd71SBarry Smith *a_P_out = NULL; /* out */ 7712e68589bSMark F. Adams PetscFunctionReturn(0); 7722e68589bSMark F. Adams } 7732e68589bSMark F. Adams 7742e68589bSMark F. Adams { 7752e68589bSMark F. Adams PetscReal *coords; 776a2f3521dSMark F. Adams PetscInt data_stride; 7770298fd71SBarry Smith PetscInt *crsGID = NULL; 7782e68589bSMark F. Adams Mat Gmat2; 7792e68589bSMark F. Adams 78071959b99SBarry Smith if (dim != data_cols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim %D != data_cols %D",dim,data_cols); 7812e68589bSMark F. Adams /* grow ghost data for better coarse grid cover of fine grid */ 7820cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7830cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 7842e68589bSMark F. Adams #endif 785a2f3521dSMark F. Adams /* messy method, squares graph and gets some data */ 786806fa848SBarry Smith ierr = getGIDsOnSquareGraph(nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID);CHKERRQ(ierr); 7872e68589bSMark F. Adams /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */ 7880cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 7890cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 7902e68589bSMark F. Adams #endif 7912e68589bSMark F. Adams /* create global vector of coorindates in 'coords' */ 792c5df96a5SBarry Smith if (size > 1) { 793806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords);CHKERRQ(ierr); 794806fa848SBarry Smith } else { 795c8b0795cSMark F. Adams coords = (PetscReal*)pc_gamg->data; 796a2f3521dSMark F. Adams data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols; 7972e68589bSMark F. Adams } 7982e68589bSMark F. Adams ierr = MatDestroy(&Gmat2);CHKERRQ(ierr); 7992e68589bSMark F. Adams 8002e68589bSMark F. Adams /* triangulate */ 8012e68589bSMark F. Adams if (dim == 2) { 802c8b0795cSMark F. Adams PetscReal metric,tm; 8030cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 8040cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr); 8052e68589bSMark F. Adams #endif 806806fa848SBarry Smith ierr = triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric);CHKERRQ(ierr); 8070cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 8080cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr); 8092e68589bSMark F. Adams #endif 8102e68589bSMark F. Adams ierr = PetscFree(crsGID);CHKERRQ(ierr); 8112e68589bSMark F. Adams 8122e68589bSMark F. Adams /* clean up and create coordinates for coarse grid (output) */ 813c5df96a5SBarry Smith if (size > 1) ierr = PetscFree(coords);CHKERRQ(ierr); 8142e68589bSMark F. Adams 815c8b0795cSMark F. Adams ierr = MPI_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm);CHKERRQ(ierr); 816c8b0795cSMark F. Adams if (tm > 1.) { /* needs to be globalized - should not happen */ 8172fa5cd67SKarl Rupp if (verbose) PetscPrintf(wcomm,"[%d]%s failed metric for coarse grid %e\n",rank,__FUNCT__,tm); 8182e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 8190298fd71SBarry Smith Prol = NULL; 820806fa848SBarry Smith } else if (metric > .0) { 8212fa5cd67SKarl Rupp if (verbose) PetscPrintf(wcomm,"[%d]%s worst metric for coarse grid = %e\n",rank,__FUNCT__,metric); 8222e68589bSMark F. Adams } 823f23aa3ddSBarry Smith } else SETERRQ(wcomm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG"); 8242e68589bSMark F. Adams { /* create next coords - output */ 8252e68589bSMark F. Adams PetscReal *crs_crds; 826806fa848SBarry Smith ierr = PetscMalloc(dim*nLocalSelected*sizeof(PetscReal), &crs_crds);CHKERRQ(ierr); 8272e68589bSMark F. Adams for (kk=0; kk<nLocalSelected; kk++) { /* grab local select nodes to promote - output */ 828b43b03e9SMark F. Adams PetscInt lid = clid_flid[kk]; 829c8b0795cSMark F. Adams for (jj=0; jj<dim; jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid]; 8302e68589bSMark F. Adams } 831c8b0795cSMark F. Adams 832c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 833c8b0795cSMark F. Adams pc_gamg->data = crs_crds; /* out */ 834c8b0795cSMark F. Adams pc_gamg->data_sz = dim*nLocalSelected; 8352e68589bSMark F. Adams } 836a2f3521dSMark F. Adams ierr = ISDestroy(&selected_2);CHKERRQ(ierr); 8372e68589bSMark F. Adams } 838a2f3521dSMark F. Adams 8392e68589bSMark F. Adams *a_P_out = Prol; /* out */ 840b43b03e9SMark F. Adams ierr = PetscFree(clid_flid);CHKERRQ(ierr); 8410cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 8420cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr); 8430cbbd2e1SMark F. Adams #endif 844c8b0795cSMark F. Adams PetscFunctionReturn(0); 845c8b0795cSMark F. Adams } 846c8b0795cSMark F. Adams 847c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 848c8b0795cSMark F. Adams /* 849c8b0795cSMark F. Adams PCCreateGAMG_GEO 850c8b0795cSMark F. Adams 851c8b0795cSMark F. Adams Input Parameter: 852c8b0795cSMark F. Adams . pc - 853c8b0795cSMark F. Adams */ 854c8b0795cSMark F. Adams #undef __FUNCT__ 855c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_GEO" 856c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_GEO(PC pc) 857c8b0795cSMark F. Adams { 858c8b0795cSMark F. Adams PetscErrorCode ierr; 859c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 860c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 861c8b0795cSMark F. Adams 862c8b0795cSMark F. Adams PetscFunctionBegin; 863c8b0795cSMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GEO; 864c8b0795cSMark F. Adams /* pc->ops->destroy = PCDestroy_GEO; */ 865c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 866c8b0795cSMark F. Adams 867c8b0795cSMark F. Adams /* set internal function pointers */ 868c8b0795cSMark F. Adams pc_gamg->graph = PCGAMGgraph_GEO; 869c8b0795cSMark F. Adams pc_gamg->coarsen = PCGAMGcoarsen_GEO; 8700cbbd2e1SMark F. Adams pc_gamg->prolongator = PCGAMGProlongator_GEO; 871c8b0795cSMark F. Adams pc_gamg->optprol = 0; 872a2f3521dSMark F. Adams pc_gamg->formkktprol = 0; 873c8b0795cSMark F. Adams 874c8b0795cSMark F. Adams pc_gamg->createdefaultdata = PCSetData_GEO; 875c8b0795cSMark F. Adams 876c8b0795cSMark F. Adams ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 877c8b0795cSMark F. Adams "PCSetCoordinates_C", 878c8b0795cSMark F. Adams "PCSetCoordinates_GEO", 879c8b0795cSMark F. Adams PCSetCoordinates_GEO);CHKERRQ(ierr); 8802e68589bSMark F. Adams PetscFunctionReturn(0); 8812e68589bSMark F. Adams } 882