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 <assert.h> 142e68589bSMark F. Adams #include <petscblaslapack.h> 152e68589bSMark F. Adams 16c8b0795cSMark F. Adams /* Private context for the GAMG preconditioner */ 17c8b0795cSMark F. Adams typedef struct{ 18c8b0795cSMark F. Adams PetscInt lid; /* local vertex index */ 19c8b0795cSMark F. Adams PetscInt degree; /* vertex degree */ 20c8b0795cSMark F. Adams } GAMGNode; 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" 362e68589bSMark F. Adams PetscErrorCode PCSetCoordinates_GEO( PC pc, PetscInt ndm, 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 ); 472e68589bSMark F. Adams ierr = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr); 482e68589bSMark F. Adams nloc = (Iend-my0)/bs; 492e68589bSMark F. Adams if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 502e68589bSMark F. Adams 51c8b0795cSMark F. Adams pc_gamg->data_cell_rows = 1; 522e68589bSMark F. Adams if( coords==0 && nloc > 0 ) { 532e68589bSMark F. Adams SETERRQ(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'."); 542e68589bSMark F. Adams } 55c8b0795cSMark F. Adams pc_gamg->data_cell_cols = ndm; /* coordinates */ 562e68589bSMark F. Adams 57c8b0795cSMark F. Adams arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 582e68589bSMark F. Adams 592e68589bSMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 602e68589bSMark F. Adams if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) { 612e68589bSMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr); 622e68589bSMark F. Adams ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr); 632e68589bSMark F. Adams } 642e68589bSMark F. Adams for(kk=0;kk<arrsz;kk++)pc_gamg->data[kk] = -999.; 652e68589bSMark F. Adams pc_gamg->data[arrsz] = -99.; 662e68589bSMark F. Adams /* copy data in - column oriented */ 672e68589bSMark F. Adams for( kk = 0 ; kk < nloc ; kk++ ){ 682e68589bSMark F. Adams for( ii = 0 ; ii < ndm ; ii++ ) { 692e68589bSMark F. Adams pc_gamg->data[ii*nloc + kk] = coords[kk*ndm + ii]; 702e68589bSMark F. Adams } 712e68589bSMark F. Adams } 722e68589bSMark F. Adams assert(pc_gamg->data[arrsz] == -99.); 732e68589bSMark F. Adams 742e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 752e68589bSMark F. Adams 762e68589bSMark F. Adams PetscFunctionReturn(0); 772e68589bSMark F. Adams } 782e68589bSMark F. Adams EXTERN_C_END 792e68589bSMark F. Adams 802e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 812e68589bSMark F. Adams /* 822e68589bSMark F. Adams PCSetData_GEO 832e68589bSMark F. Adams 842e68589bSMark F. Adams Input Parameter: 852e68589bSMark F. Adams . pc - 862e68589bSMark F. Adams */ 872e68589bSMark F. Adams #undef __FUNCT__ 882e68589bSMark F. Adams #define __FUNCT__ "PCSetData_GEO" 892e68589bSMark F. Adams PetscErrorCode PCSetData_GEO( PC pc ) 902e68589bSMark F. Adams { 912e68589bSMark F. Adams PetscFunctionBegin; 922e68589bSMark F. Adams SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"GEO MG needs coordinates"); 932e68589bSMark F. Adams } 942e68589bSMark F. Adams 952e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 962e68589bSMark F. Adams /* 972e68589bSMark F. Adams PCSetFromOptions_GEO 982e68589bSMark F. Adams 992e68589bSMark F. Adams Input Parameter: 1002e68589bSMark F. Adams . pc - 1012e68589bSMark F. Adams */ 1022e68589bSMark F. Adams #undef __FUNCT__ 1032e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GEO" 1042e68589bSMark F. Adams PetscErrorCode PCSetFromOptions_GEO( PC pc ) 1052e68589bSMark F. Adams { 1062e68589bSMark F. Adams PetscErrorCode ierr; 1072e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1082e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1092e68589bSMark F. Adams 1102e68589bSMark F. Adams PetscFunctionBegin; 1112e68589bSMark F. Adams ierr = PetscOptionsHead("GAMG-GEO options"); CHKERRQ(ierr); 1122e68589bSMark F. Adams { 1132e68589bSMark F. Adams /* -pc_gamg_sa_nsmooths */ 1142e68589bSMark F. Adams /* pc_gamg_sa->smooths = 0; */ 1152e68589bSMark F. Adams /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */ 1162e68589bSMark F. Adams /* "smoothing steps for smoothed aggregation, usually 1 (0)", */ 1172e68589bSMark F. Adams /* "PCGAMGSetNSmooths_AGG", */ 1182e68589bSMark F. Adams /* pc_gamg_sa->smooths, */ 1192e68589bSMark F. Adams /* &pc_gamg_sa->smooths, */ 1202e68589bSMark F. Adams /* &flag); */ 1212e68589bSMark F. Adams /* CHKERRQ(ierr); */ 1222e68589bSMark F. Adams } 1232e68589bSMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1242e68589bSMark F. Adams 1252e68589bSMark F. Adams /* call base class */ 1262e68589bSMark F. Adams ierr = PCSetFromOptions_GAMG( pc ); CHKERRQ(ierr); 1272e68589bSMark F. Adams 1282e68589bSMark F. Adams if( pc_gamg->verbose ) { 129*41b27cdeSMark F. Adams MPI_Comm wcomm = ((PetscObject)pc)->comm; 130*41b27cdeSMark F. Adams PetscPrintf(wcomm,"[%d]%s done\n",0,__FUNCT__); 1312e68589bSMark F. Adams } 1322e68589bSMark F. Adams 1332e68589bSMark F. Adams PetscFunctionReturn(0); 1342e68589bSMark F. Adams } 1352e68589bSMark F. Adams 1362e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1372e68589bSMark F. Adams /* 1382e68589bSMark F. Adams triangulateAndFormProl 1392e68589bSMark F. Adams 1402e68589bSMark F. Adams Input Parameter: 1412e68589bSMark F. Adams . selected_2 - list of selected local ID, includes selected ghosts 1422e68589bSMark F. Adams . nnodes - 1432e68589bSMark F. Adams . coords[2*nnodes] - column vector of local coordinates w/ ghosts 144b43b03e9SMark F. Adams . nselected_1 - selected IDs that go with base (1) graph 1450cbbd2e1SMark F. Adams . clid_lid_1[nselected_1] - lids of selected (c) nodes ??????????? 1460cbbd2e1SMark F. Adams . agg_lists_1 - list of aggregates 1472e68589bSMark F. Adams . crsGID[selected.size()] - global index for prolongation operator 1482e68589bSMark F. Adams . bs - block size 1492e68589bSMark F. Adams Output Parameter: 1502e68589bSMark F. Adams . a_Prol - prolongation operator 1512e68589bSMark F. Adams . a_worst_best - measure of worst missed fine vertex, 0 is no misses 1522e68589bSMark F. Adams */ 1532e68589bSMark F. Adams #undef __FUNCT__ 1542e68589bSMark F. Adams #define __FUNCT__ "triangulateAndFormProl" 1550cbbd2e1SMark F. Adams static PetscErrorCode triangulateAndFormProl( IS selected_2, /* list of selected local ID, includes selected ghosts */ 1562e68589bSMark F. Adams const PetscInt nnodes, 1572e68589bSMark F. Adams const PetscReal coords[], /* column vector of local coordinates w/ ghosts */ 158b43b03e9SMark F. Adams const PetscInt nselected_1, /* list of selected local ID, includes selected ghosts */ 159b43b03e9SMark F. Adams const PetscInt clid_lid_1[], 1600cbbd2e1SMark F. Adams const PetscCoarsenData *agg_lists_1, /* selected_1 vertices of aggregate unselected vertices */ 1612e68589bSMark F. Adams const PetscInt crsGID[], 1622e68589bSMark F. Adams const PetscInt bs, 1632e68589bSMark F. Adams Mat a_Prol, /* prolongation operator (output) */ 1642e68589bSMark F. Adams PetscReal *a_worst_best /* measure of worst missed fine vertex, 0 is no misses */ 1652e68589bSMark F. Adams ) 1662e68589bSMark F. Adams { 1672e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE) 1682e68589bSMark F. Adams PetscErrorCode ierr; 169b43b03e9SMark F. Adams PetscInt jj,tid,tt,idx,nselected_2; 1702e68589bSMark F. Adams struct triangulateio in,mid; 1710cbbd2e1SMark F. Adams const PetscInt *selected_idx_2; 1722e68589bSMark F. Adams PetscMPIInt mype,npe; 1732e68589bSMark F. Adams PetscInt Istart,Iend,nFineLoc,myFine0; 1742e68589bSMark F. Adams int kk,nPlotPts,sid; 175c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Prol)->comm; 176c8b0795cSMark F. Adams PetscReal tm; 1772e68589bSMark F. Adams PetscFunctionBegin; 178c8b0795cSMark F. Adams 179c8b0795cSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype); CHKERRQ(ierr); 180c8b0795cSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe); CHKERRQ(ierr); 181c8b0795cSMark F. Adams ierr = ISGetSize( selected_2, &nselected_2 ); CHKERRQ(ierr); 1822e68589bSMark F. Adams if(nselected_2 == 1 || nselected_2 == 2 ){ /* 0 happens on idle processors */ 1832e68589bSMark F. Adams *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */ 184c8b0795cSMark F. Adams } 185c8b0795cSMark F. Adams else *a_worst_best = 0.0; 186c8b0795cSMark F. Adams ierr = MPI_Allreduce( a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm ); CHKERRQ(ierr); 187c8b0795cSMark F. Adams if( tm > 0.0 ) { 188c8b0795cSMark F. Adams *a_worst_best = 100.0; 1892e68589bSMark F. Adams PetscFunctionReturn(0); 1902e68589bSMark F. Adams } 1912e68589bSMark F. Adams ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend ); CHKERRQ(ierr); 1922e68589bSMark F. Adams nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs; 1932e68589bSMark F. Adams nPlotPts = nFineLoc; /* locals */ 1942e68589bSMark F. Adams /* traingle */ 1952e68589bSMark F. Adams /* Define input points - in*/ 1962e68589bSMark F. Adams in.numberofpoints = nselected_2; 1972e68589bSMark F. Adams in.numberofpointattributes = 0; 1982e68589bSMark F. Adams /* get nselected points */ 1992e68589bSMark F. Adams ierr = PetscMalloc( 2*(nselected_2)*sizeof(REAL), &in.pointlist ); CHKERRQ(ierr); 2002e68589bSMark F. Adams ierr = ISGetIndices( selected_2, &selected_idx_2 ); CHKERRQ(ierr); 2012e68589bSMark F. Adams 2022e68589bSMark F. Adams for(kk=0,sid=0;kk<nselected_2;kk++,sid += 2){ 2032e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 2042e68589bSMark F. Adams in.pointlist[sid] = coords[lid]; 2052e68589bSMark F. Adams in.pointlist[sid+1] = coords[nnodes + lid]; 2062e68589bSMark F. Adams if(lid>=nFineLoc) nPlotPts++; 2072e68589bSMark F. Adams } 2082e68589bSMark F. Adams assert(sid==2*nselected_2); 2092e68589bSMark F. Adams 2102e68589bSMark F. Adams in.numberofsegments = 0; 2112e68589bSMark F. Adams in.numberofedges = 0; 2122e68589bSMark F. Adams in.numberofholes = 0; 2132e68589bSMark F. Adams in.numberofregions = 0; 2142e68589bSMark F. Adams in.trianglelist = 0; 2152e68589bSMark F. Adams in.segmentmarkerlist = 0; 2162e68589bSMark F. Adams in.pointattributelist = 0; 2172e68589bSMark F. Adams in.pointmarkerlist = 0; 2182e68589bSMark F. Adams in.triangleattributelist = 0; 2192e68589bSMark F. Adams in.trianglearealist = 0; 2202e68589bSMark F. Adams in.segmentlist = 0; 2212e68589bSMark F. Adams in.holelist = 0; 2222e68589bSMark F. Adams in.regionlist = 0; 2232e68589bSMark F. Adams in.edgelist = 0; 2242e68589bSMark F. Adams in.edgemarkerlist = 0; 2252e68589bSMark F. Adams in.normlist = 0; 2262e68589bSMark F. Adams /* triangulate */ 2272e68589bSMark F. Adams mid.pointlist = 0; /* Not needed if -N switch used. */ 2282e68589bSMark F. Adams /* Not needed if -N switch used or number of point attributes is zero: */ 2292e68589bSMark F. Adams mid.pointattributelist = 0; 2302e68589bSMark F. Adams mid.pointmarkerlist = 0; /* Not needed if -N or -B switch used. */ 2312e68589bSMark F. Adams mid.trianglelist = 0; /* Not needed if -E switch used. */ 2322e68589bSMark F. Adams /* Not needed if -E switch used or number of triangle attributes is zero: */ 2332e68589bSMark F. Adams mid.triangleattributelist = 0; 2342e68589bSMark F. Adams mid.neighborlist = 0; /* Needed only if -n switch used. */ 2352e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P not used: */ 2362e68589bSMark F. Adams mid.segmentlist = 0; 2372e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P and -B not used: */ 2382e68589bSMark F. Adams mid.segmentmarkerlist = 0; 2392e68589bSMark F. Adams mid.edgelist = 0; /* Needed only if -e switch used. */ 2402e68589bSMark F. Adams mid.edgemarkerlist = 0; /* Needed if -e used and -B not used. */ 2412e68589bSMark F. Adams mid.numberoftriangles = 0; 2422e68589bSMark F. Adams 2432e68589bSMark F. Adams /* Triangulate the points. Switches are chosen to read and write a */ 2442e68589bSMark F. Adams /* PSLG (p), preserve the convex hull (c), number everything from */ 2452e68589bSMark F. Adams /* zero (z), assign a regional attribute to each element (A), and */ 2462e68589bSMark F. Adams /* produce an edge list (e), a Voronoi diagram (v), and a triangle */ 2472e68589bSMark F. Adams /* neighbor list (n). */ 2482e68589bSMark F. Adams if(nselected_2 != 0){ /* inactive processor */ 2492e68589bSMark F. Adams char args[] = "npczQ"; /* c is needed ? */ 2502e68589bSMark F. Adams triangulate(args, &in, &mid, (struct triangulateio *) NULL ); 2512e68589bSMark F. Adams /* output .poly files for 'showme' */ 2522e68589bSMark F. Adams if( !PETSC_TRUE ) { 2532e68589bSMark F. Adams static int level = 1; 2542e68589bSMark F. Adams FILE *file; char fname[32]; 2552e68589bSMark F. Adams 2562e68589bSMark F. Adams sprintf(fname,"C%d_%d.poly",level,mype); file = fopen(fname, "w"); 2572e68589bSMark F. Adams /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/ 2582e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); 2592e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2602e68589bSMark F. Adams for(kk=0,sid=0;kk<in.numberofpoints;kk++,sid += 2){ 2612e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2622e68589bSMark F. Adams } 2632e68589bSMark F. Adams /*One line: <# of segments> <# of boundary markers (0 or 1)> */ 2642e68589bSMark F. Adams fprintf(file, "%d %d\n",0,0); 2652e68589bSMark F. Adams /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */ 2662e68589bSMark F. Adams /* One line: <# of holes> */ 2672e68589bSMark F. Adams fprintf(file, "%d\n",0); 2682e68589bSMark F. Adams /* Following lines: <hole #> <x> <y> */ 2692e68589bSMark F. Adams /* Optional line: <# of regional attributes and/or area constraints> */ 2702e68589bSMark F. Adams /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */ 2712e68589bSMark F. Adams fclose(file); 2722e68589bSMark F. Adams 2732e68589bSMark F. Adams /* elems */ 2742e68589bSMark F. Adams sprintf(fname,"C%d_%d.ele",level,mype); file = fopen(fname, "w"); 2752e68589bSMark F. Adams /* First line: <# of triangles> <nodes per triangle> <# of attributes> */ 2762e68589bSMark F. Adams fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0); 2772e68589bSMark F. Adams /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */ 2782e68589bSMark F. Adams for(kk=0,sid=0;kk<mid.numberoftriangles;kk++,sid += 3){ 2792e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]); 2802e68589bSMark F. Adams } 2812e68589bSMark F. Adams fclose(file); 2822e68589bSMark F. Adams 2832e68589bSMark F. Adams sprintf(fname,"C%d_%d.node",level,mype); file = fopen(fname, "w"); 2842e68589bSMark F. Adams /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */ 2852e68589bSMark F. Adams /* fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); */ 2862e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",nPlotPts,2,0,0); 2872e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2882e68589bSMark F. Adams for(kk=0,sid=0;kk<in.numberofpoints;kk++,sid+=2){ 2892e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2902e68589bSMark F. Adams } 2912e68589bSMark F. Adams 2922e68589bSMark F. Adams sid /= 2; 2932e68589bSMark F. Adams for(jj=0;jj<nFineLoc;jj++){ 2942e68589bSMark F. Adams PetscBool sel = PETSC_TRUE; 2952e68589bSMark F. Adams for( kk=0 ; kk<nselected_2 && sel ; kk++ ){ 2962e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 2972e68589bSMark F. Adams if( lid == jj ) sel = PETSC_FALSE; 2982e68589bSMark F. Adams } 2992e68589bSMark F. Adams if( sel ) { 3002e68589bSMark F. Adams fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[nnodes + jj]); 3012e68589bSMark F. Adams } 3022e68589bSMark F. Adams } 3032e68589bSMark F. Adams fclose(file); 3042e68589bSMark F. Adams assert(sid==nPlotPts); 3052e68589bSMark F. Adams level++; 3062e68589bSMark F. Adams } 3072e68589bSMark F. Adams } 3080cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3090cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 3102e68589bSMark F. Adams #endif 3112e68589bSMark F. Adams { /* form P - setup some maps */ 3120cbbd2e1SMark F. Adams PetscInt clid,mm,*nTri,*node_tri; 3132e68589bSMark F. Adams 3142e68589bSMark F. Adams ierr = PetscMalloc( nselected_2*sizeof(PetscInt), &node_tri ); CHKERRQ(ierr); 3152e68589bSMark F. Adams ierr = PetscMalloc( nselected_2*sizeof(PetscInt), &nTri ); CHKERRQ(ierr); 3162e68589bSMark F. Adams 3172e68589bSMark F. Adams /* need list of triangles on node */ 3182e68589bSMark F. Adams for(kk=0;kk<nselected_2;kk++) nTri[kk] = 0; 3192e68589bSMark F. Adams for(tid=0,kk=0;tid<mid.numberoftriangles;tid++){ 3202e68589bSMark F. Adams for(jj=0;jj<3;jj++) { 3212e68589bSMark F. Adams PetscInt cid = mid.trianglelist[kk++]; 3222e68589bSMark F. Adams if( nTri[cid] == 0 ) node_tri[cid] = tid; 3232e68589bSMark F. Adams nTri[cid]++; 3242e68589bSMark F. Adams } 3252e68589bSMark F. Adams } 3262e68589bSMark F. Adams #define EPS 1.e-12 3272e68589bSMark F. Adams /* find points and set prolongation */ 3280cbbd2e1SMark F. Adams for( mm = clid = 0 ; mm < nFineLoc ; mm++ ){ 329*41b27cdeSMark F. Adams if( (jj=PetscCDSizeAt(agg_lists_1,mm)) > 0 ) { 3300cbbd2e1SMark F. Adams const PetscInt lid = mm; 3310cbbd2e1SMark F. Adams //for(clid_iterator=0;clid_iterator<nselected_1;clid_iterator++){ 3320cbbd2e1SMark F. Adams //PetscInt flid = clid_lid_1[clid_iterator]; assert(flid != -1); 3332e68589bSMark F. Adams PetscScalar AA[3][3]; 3342e68589bSMark F. Adams PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO; 335*41b27cdeSMark F. Adams PetscCDPos pos; 3360cbbd2e1SMark F. Adams 337*41b27cdeSMark F. Adams for( pos=PetscCDGetHeadPos(agg_lists_1,lid) ; 3380cbbd2e1SMark F. Adams pos ; 339*41b27cdeSMark F. Adams pos=PetscCDGetNextPos(agg_lists_1,lid,pos) ){ 3400cbbd2e1SMark F. Adams PetscInt flid = LLNGetID(pos); 3412e68589bSMark F. Adams if( flid < nFineLoc ) { /* could be a ghost */ 3422e68589bSMark F. Adams PetscInt bestTID = -1; PetscReal best_alpha = 1.e10; 3432e68589bSMark F. Adams const PetscInt fgid = flid + myFine0; 3442e68589bSMark F. Adams /* compute shape function for gid */ 3452e68589bSMark F. Adams const PetscReal fcoord[3] = {coords[flid],coords[nnodes+flid],1.0}; 3462e68589bSMark F. Adams PetscBool haveit=PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3]; 3472e68589bSMark F. Adams /* look for it */ 3480cbbd2e1SMark F. Adams for( tid = node_tri[clid], jj=0; 3492e68589bSMark F. Adams jj < 5 && !haveit && tid != -1; 3502e68589bSMark F. Adams jj++ ){ 3512e68589bSMark F. Adams for(tt=0;tt<3;tt++){ 3522e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3532e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 3542e68589bSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0; 3552e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3562e68589bSMark F. Adams } 3572e68589bSMark F. Adams 3582e68589bSMark F. Adams for(tt=0;tt<3;tt++) alpha[tt] = (PetscScalar)fcoord[tt]; 3592e68589bSMark F. Adams 3602e68589bSMark F. Adams /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */ 3612e68589bSMark F. Adams LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO); 3622e68589bSMark F. Adams { 3632e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal lowest=1.e10; 3642e68589bSMark F. Adams for( tt = 0, idx = 0 ; tt < 3 ; tt++ ) { 3652e68589bSMark F. Adams if( PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS ) have = PETSC_FALSE; 3662e68589bSMark F. Adams if( PetscRealPart(alpha[tt]) < lowest ){ 3672e68589bSMark F. Adams lowest = PetscRealPart(alpha[tt]); 3682e68589bSMark F. Adams idx = tt; 3692e68589bSMark F. Adams } 3702e68589bSMark F. Adams } 3712e68589bSMark F. Adams haveit = have; 3722e68589bSMark F. Adams } 3732e68589bSMark F. Adams tid = mid.neighborlist[3*tid + idx]; 3742e68589bSMark F. Adams } 3752e68589bSMark F. Adams 3762e68589bSMark F. Adams if( !haveit ) { 3772e68589bSMark F. Adams /* brute force */ 3782e68589bSMark F. Adams for(tid=0 ; tid<mid.numberoftriangles && !haveit ; tid++ ){ 3792e68589bSMark F. Adams for(tt=0;tt<3;tt++){ 3802e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3812e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 3822e68589bSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0; 3832e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3842e68589bSMark F. Adams } 3852e68589bSMark F. Adams for(tt=0;tt<3;tt++) alpha[tt] = fcoord[tt]; 3862e68589bSMark F. Adams /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */ 3872e68589bSMark F. Adams LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO); 3882e68589bSMark F. Adams { 3892e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal worst=0.0, v; 3902e68589bSMark F. Adams for(tt=0; tt<3 && have ;tt++) { 3912e68589bSMark F. Adams if( PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS ) have=PETSC_FALSE; 3922e68589bSMark F. Adams if( (v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst ) worst = v; 3932e68589bSMark F. Adams } 3942e68589bSMark F. Adams if( worst < best_alpha ) { 3952e68589bSMark F. Adams best_alpha = worst; bestTID = tid; 3962e68589bSMark F. Adams } 3972e68589bSMark F. Adams haveit = have; 3982e68589bSMark F. Adams } 3992e68589bSMark F. Adams } 4002e68589bSMark F. Adams } 4012e68589bSMark F. Adams if( !haveit ) { 4022e68589bSMark F. Adams if( best_alpha > *a_worst_best ) *a_worst_best = best_alpha; 4032e68589bSMark F. Adams /* use best one */ 4042e68589bSMark F. Adams for(tt=0;tt<3;tt++){ 4052e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*bestTID + tt]; 4062e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 4072e68589bSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0; 4082e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 4092e68589bSMark F. Adams } 4102e68589bSMark F. Adams for(tt=0;tt<3;tt++) alpha[tt] = fcoord[tt]; 4112e68589bSMark F. Adams /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */ 4122e68589bSMark F. Adams LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO); 4132e68589bSMark F. Adams } 4142e68589bSMark F. Adams 4152e68589bSMark F. Adams /* put in row of P */ 4162e68589bSMark F. Adams for(idx=0;idx<3;idx++){ 4172e68589bSMark F. Adams PetscScalar shp = alpha[idx]; 4182e68589bSMark F. Adams if( PetscAbs(PetscRealPart(shp)) > 1.e-6 ) { 4192e68589bSMark F. Adams PetscInt cgid = crsGID[clids[idx]]; 4202e68589bSMark F. Adams PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */ 4212e68589bSMark F. Adams for(tt=0 ; tt < bs ; tt++, ii++, jj++ ){ 4222e68589bSMark F. Adams ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES); CHKERRQ(ierr); 4232e68589bSMark F. Adams } 4242e68589bSMark F. Adams } 4252e68589bSMark F. Adams } 4262e68589bSMark F. Adams } 4270cbbd2e1SMark F. Adams } /* aggregates iterations */ 4280cbbd2e1SMark F. Adams clid++; 4290cbbd2e1SMark F. Adams } /* a coarse agg */ 4300cbbd2e1SMark F. Adams } /* for all fine nodes */ 4310cbbd2e1SMark F. Adams 4322e68589bSMark F. Adams ierr = ISRestoreIndices( selected_2, &selected_idx_2 ); CHKERRQ(ierr); 4332e68589bSMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4342e68589bSMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4352e68589bSMark F. Adams 4362e68589bSMark F. Adams ierr = PetscFree( node_tri ); CHKERRQ(ierr); 4372e68589bSMark F. Adams ierr = PetscFree( nTri ); CHKERRQ(ierr); 4382e68589bSMark F. Adams } 4390cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4400cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 4412e68589bSMark F. Adams #endif 4422e68589bSMark F. Adams free( mid.trianglelist ); 4432e68589bSMark F. Adams free( mid.neighborlist ); 4442e68589bSMark F. Adams ierr = PetscFree( in.pointlist ); CHKERRQ(ierr); 4452e68589bSMark F. Adams 4462e68589bSMark F. Adams PetscFunctionReturn(0); 4472e68589bSMark F. Adams #else 4482e68589bSMark F. Adams SETERRQ(((PetscObject)a_Prol)->comm,PETSC_ERR_LIB,"configure with TRIANGLE to use geometric MG"); 4492e68589bSMark F. Adams #endif 4502e68589bSMark F. Adams } 4512e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 4522e68589bSMark F. Adams /* 4532e68589bSMark F. Adams getGIDsOnSquareGraph - square graph, get 4542e68589bSMark F. Adams 4552e68589bSMark F. Adams Input Parameter: 4560cbbd2e1SMark F. Adams . nselected_1 - selected local indices (includes ghosts in input Gmat1) 457b43b03e9SMark F. Adams . clid_lid_1 - [nselected_1] lids of selected nodes 4582e68589bSMark F. Adams . Gmat1 - graph that goes with 'selected_1' 4592e68589bSMark F. Adams Output Parameter: 4602e68589bSMark F. Adams . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2) 4612e68589bSMark F. Adams . a_Gmat_2 - graph that is squared of 'Gmat_1' 4622e68589bSMark F. Adams . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes 4632e68589bSMark F. Adams */ 4642e68589bSMark F. Adams #undef __FUNCT__ 4652e68589bSMark F. Adams #define __FUNCT__ "getGIDsOnSquareGraph" 4660cbbd2e1SMark F. Adams static PetscErrorCode getGIDsOnSquareGraph( const PetscInt nselected_1, 467b43b03e9SMark F. Adams const PetscInt clid_lid_1[], 4682e68589bSMark F. Adams const Mat Gmat1, 4692e68589bSMark F. Adams IS *a_selected_2, 4702e68589bSMark F. Adams Mat *a_Gmat_2, 4712e68589bSMark F. Adams PetscInt **a_crsGID 4722e68589bSMark F. Adams ) 4732e68589bSMark F. Adams { 4742e68589bSMark F. Adams PetscErrorCode ierr; 4752e68589bSMark F. Adams PetscMPIInt mype,npe; 476b43b03e9SMark F. Adams PetscInt *crsGID, kk,my0,Iend,nloc; 4772e68589bSMark F. Adams MPI_Comm wcomm = ((PetscObject)Gmat1)->comm; 4782e68589bSMark F. Adams 4792e68589bSMark F. Adams PetscFunctionBegin; 4802e68589bSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 4812e68589bSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 4822e68589bSMark F. Adams ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend); CHKERRQ(ierr); /* AIJ */ 4832e68589bSMark F. Adams nloc = Iend - my0; /* this does not change */ 4842e68589bSMark F. Adams 4852e68589bSMark F. Adams if (npe == 1) { /* not much to do in serial */ 486b43b03e9SMark F. Adams ierr = PetscMalloc( nselected_1*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr); 487b43b03e9SMark F. Adams for(kk=0;kk<nselected_1;kk++) crsGID[kk] = kk; 4882e68589bSMark F. Adams *a_Gmat_2 = 0; 489b43b03e9SMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2); 490b43b03e9SMark F. Adams CHKERRQ(ierr); 4912e68589bSMark F. Adams } 4922e68589bSMark F. Adams else { 493b43b03e9SMark F. Adams PetscInt idx,num_fine_ghosts,num_crs_ghost,myCrs0; 4942e68589bSMark F. Adams Mat_MPIAIJ *mpimat2; 4952e68589bSMark F. Adams Mat Gmat2; 4962e68589bSMark F. Adams Vec locState; 4972e68589bSMark F. Adams PetscScalar *cpcol_state; 4982e68589bSMark F. Adams 4992e68589bSMark F. Adams /* scan my coarse zero gid, set 'lid_state' with coarse GID */ 500b43b03e9SMark F. Adams kk = nselected_1; 501b43b03e9SMark F. Adams MPI_Scan( &kk, &myCrs0, 1, MPIU_INT, MPIU_SUM, wcomm ); 502b43b03e9SMark F. Adams myCrs0 -= nselected_1; 5032e68589bSMark F. Adams 504b43b03e9SMark F. Adams if( a_Gmat_2 ) { /* output */ 5052e68589bSMark F. Adams /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */ 5062e68589bSMark F. Adams ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 ); CHKERRQ(ierr); 5072e68589bSMark F. Adams *a_Gmat_2 = Gmat2; /* output */ 5082e68589bSMark F. Adams } 5092e68589bSMark F. Adams else Gmat2 = Gmat1; /* use local to get crsGIDs at least */ 5102e68589bSMark F. Adams /* get coarse grid GIDS for selected (locals and ghosts) */ 5112e68589bSMark F. Adams mpimat2 = (Mat_MPIAIJ*)Gmat2->data; 5122e68589bSMark F. Adams ierr = MatGetVecs( Gmat2, &locState, 0 ); CHKERRQ(ierr); 5132e68589bSMark F. Adams ierr = VecSet( locState, (PetscScalar)(PetscReal)(-1) ); CHKERRQ(ierr); /* set with UNKNOWN state */ 514b43b03e9SMark F. Adams for(kk=0;kk<nselected_1;kk++){ 515b43b03e9SMark F. Adams PetscInt fgid = clid_lid_1[kk] + my0; 5162e68589bSMark F. Adams PetscScalar v = (PetscScalar)(kk+myCrs0); 5172e68589bSMark F. Adams ierr = VecSetValues( locState, 1, &fgid, &v, INSERT_VALUES ); CHKERRQ(ierr); /* set with PID */ 5182e68589bSMark F. Adams } 5192e68589bSMark F. Adams ierr = VecAssemblyBegin( locState ); CHKERRQ(ierr); 5202e68589bSMark F. Adams ierr = VecAssemblyEnd( locState ); CHKERRQ(ierr); 5212e68589bSMark F. Adams ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5222e68589bSMark F. Adams ierr = VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5232e68589bSMark F. Adams ierr = VecGetLocalSize( mpimat2->lvec, &num_fine_ghosts ); CHKERRQ(ierr); 5242e68589bSMark F. Adams ierr = VecGetArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr); 5252e68589bSMark F. Adams for(kk=0,num_crs_ghost=0;kk<num_fine_ghosts;kk++){ 5262e68589bSMark F. Adams if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ) num_crs_ghost++; 5272e68589bSMark F. Adams } 528b43b03e9SMark F. Adams ierr = PetscMalloc( (nselected_1+num_crs_ghost)*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr); /* output */ 5292e68589bSMark F. Adams { 5302e68589bSMark F. Adams PetscInt *selected_set; 531b43b03e9SMark F. Adams ierr = PetscMalloc( (nselected_1+num_crs_ghost)*sizeof(PetscInt), &selected_set ); CHKERRQ(ierr); 5322e68589bSMark F. Adams /* do ghost of 'crsGID' */ 533b43b03e9SMark F. Adams for(kk=0,idx=nselected_1;kk<num_fine_ghosts;kk++){ 5342e68589bSMark F. Adams if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){ 5352e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5362e68589bSMark F. Adams selected_set[idx] = nloc + kk; 5372e68589bSMark F. Adams crsGID[idx++] = cgid; 5382e68589bSMark F. Adams } 5392e68589bSMark F. Adams } 540b43b03e9SMark F. Adams assert(idx==(nselected_1+num_crs_ghost)); 5412e68589bSMark F. Adams ierr = VecRestoreArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr); 5422e68589bSMark F. Adams /* do locals in 'crsGID' */ 5432e68589bSMark F. Adams ierr = VecGetArray( locState, &cpcol_state ); CHKERRQ(ierr); 5442e68589bSMark F. Adams for(kk=0,idx=0;kk<nloc;kk++){ 5452e68589bSMark F. Adams if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){ 5462e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5472e68589bSMark F. Adams selected_set[idx] = kk; 5482e68589bSMark F. Adams crsGID[idx++] = cgid; 5492e68589bSMark F. Adams } 5502e68589bSMark F. Adams } 551b43b03e9SMark F. Adams assert(idx==nselected_1); 5522e68589bSMark F. Adams ierr = VecRestoreArray( locState, &cpcol_state ); CHKERRQ(ierr); 5532e68589bSMark F. Adams 5542e68589bSMark F. Adams if( a_selected_2 != 0 ) { /* output */ 5550cbbd2e1SMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2); 5562e68589bSMark F. Adams CHKERRQ(ierr); 5572e68589bSMark F. Adams } 5580cbbd2e1SMark F. Adams else { 5592e68589bSMark F. Adams ierr = PetscFree( selected_set ); CHKERRQ(ierr); 5602e68589bSMark F. Adams } 5610cbbd2e1SMark F. Adams } 5622e68589bSMark F. Adams ierr = VecDestroy( &locState ); CHKERRQ(ierr); 5632e68589bSMark F. Adams } 5642e68589bSMark F. Adams *a_crsGID = crsGID; /* output */ 5652e68589bSMark F. Adams 5662e68589bSMark F. Adams PetscFunctionReturn(0); 5672e68589bSMark F. Adams } 5682e68589bSMark F. Adams 5692e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 5702e68589bSMark F. Adams /* 571c8b0795cSMark F. Adams PCGAMGgraph_GEO 5722e68589bSMark F. Adams 5732e68589bSMark F. Adams Input Parameter: 5742e68589bSMark F. Adams . pc - this 5752e68589bSMark F. Adams . Amat - matrix on this fine level 5762e68589bSMark F. Adams Output Parameter: 577c8b0795cSMark F. Adams . a_Gmat 5782e68589bSMark F. Adams */ 5792e68589bSMark F. Adams #undef __FUNCT__ 580c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_GEO" 581c8b0795cSMark F. Adams PetscErrorCode PCGAMGgraph_GEO( PC pc, 5822e68589bSMark F. Adams const Mat Amat, 583c8b0795cSMark F. Adams Mat *a_Gmat 584c8b0795cSMark F. Adams ) 585c8b0795cSMark F. Adams { 586c8b0795cSMark F. Adams PetscErrorCode ierr; 587c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 588c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 589c8b0795cSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 590c8b0795cSMark F. Adams const PetscReal vfilter = pc_gamg->threshold; 591c8b0795cSMark F. Adams PetscMPIInt mype,npe; 592c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat)->comm; 593c8b0795cSMark F. Adams Mat Gmat; 5940cbbd2e1SMark F. Adams PetscBool set,flg,symm; 595c8b0795cSMark F. Adams PetscFunctionBegin; 5960cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 5970cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr); 5980cbbd2e1SMark F. Adams #endif 599c8b0795cSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 600c8b0795cSMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 601c8b0795cSMark F. Adams 6020cbbd2e1SMark F. Adams ierr = MatIsSymmetricKnown(Amat, &set, &flg); CHKERRQ(ierr); 603263489e9SJed Brown symm = (PetscBool)!(set && flg); 6040cbbd2e1SMark F. Adams 6050cbbd2e1SMark F. Adams ierr = PCGAMGCreateSimpleGraph( Amat, &Gmat ); CHKERRQ( ierr ); 6060cbbd2e1SMark F. Adams ierr = PCGAMGScaleFilterGraph( &Gmat, vfilter, symm, verbose ); CHKERRQ( ierr ); 607c8b0795cSMark F. Adams 608c8b0795cSMark F. Adams *a_Gmat = Gmat; 6090cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 6100cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr); 6110cbbd2e1SMark F. Adams #endif 612c8b0795cSMark F. Adams PetscFunctionReturn(0); 613c8b0795cSMark F. Adams } 614c8b0795cSMark F. Adams 615c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 616c8b0795cSMark F. Adams /* 617c8b0795cSMark F. Adams PCGAMGcoarsen_GEO 618c8b0795cSMark F. Adams 619c8b0795cSMark F. Adams Input Parameter: 620e0940f08SMark F. Adams . a_pc - this 621e0940f08SMark F. Adams . a_Gmat - graph 622c8b0795cSMark F. Adams Output Parameter: 623c8b0795cSMark F. Adams . a_llist_parent - linked list from selected indices for data locality only 624c8b0795cSMark F. Adams */ 625c8b0795cSMark F. Adams #undef __FUNCT__ 626c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGcoarsen_GEO" 627e0940f08SMark F. Adams PetscErrorCode PCGAMGcoarsen_GEO( PC a_pc, 628e0940f08SMark F. Adams Mat *a_Gmat, 6290cbbd2e1SMark F. Adams PetscCoarsenData **a_llist_parent 630c8b0795cSMark F. Adams ) 631c8b0795cSMark F. Adams { 632c8b0795cSMark F. Adams PetscErrorCode ierr; 633e0940f08SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 634c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 635c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,kk,Ii,ncols; 636c8b0795cSMark F. Adams PetscMPIInt mype,npe; 6370cbbd2e1SMark F. Adams IS perm; 638c8b0795cSMark F. Adams GAMGNode *gnodes; 639c8b0795cSMark F. Adams PetscInt *permute; 640e0940f08SMark F. Adams Mat Gmat = *a_Gmat; 641c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)Gmat)->comm; 642b43b03e9SMark F. Adams MatCoarsen crs; 643c8b0795cSMark F. Adams 644c8b0795cSMark F. Adams PetscFunctionBegin; 6450cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 6460cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr); 6470cbbd2e1SMark F. Adams #endif 648c8b0795cSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 649c8b0795cSMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 650c8b0795cSMark F. Adams ierr = MatGetOwnershipRange( Gmat, &Istart, &Iend ); CHKERRQ(ierr); 651c8b0795cSMark F. Adams nloc = (Iend-Istart); 652c8b0795cSMark F. Adams 653c8b0795cSMark F. Adams /* create random permutation with sort for geo-mg */ 654c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(GAMGNode), &gnodes ); CHKERRQ(ierr); 655c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr); 656c8b0795cSMark F. Adams 657c8b0795cSMark F. Adams for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */ 658c8b0795cSMark F. Adams ierr = MatGetRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr); 659c8b0795cSMark F. Adams { 660c8b0795cSMark F. Adams PetscInt lid = Ii - Istart; 661c8b0795cSMark F. Adams gnodes[lid].lid = lid; 662c8b0795cSMark F. Adams gnodes[lid].degree = ncols; 663c8b0795cSMark F. Adams } 664c8b0795cSMark F. Adams ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr); 665c8b0795cSMark F. Adams } 666c8b0795cSMark F. Adams /* randomize */ 667c8b0795cSMark F. Adams srand(1); /* make deterministic */ 668c8b0795cSMark F. Adams if( PETSC_TRUE ) { 669c8b0795cSMark F. Adams PetscBool *bIndexSet; 670c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr); 671c8b0795cSMark F. Adams for ( Ii = 0; Ii < nloc ; Ii++) bIndexSet[Ii] = PETSC_FALSE; 672c8b0795cSMark F. Adams for ( Ii = 0; Ii < nloc ; Ii++) 673c8b0795cSMark F. Adams { 674c8b0795cSMark F. Adams PetscInt iSwapIndex = rand()%nloc; 675c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) 676c8b0795cSMark F. Adams { 677c8b0795cSMark F. Adams GAMGNode iTemp = gnodes[iSwapIndex]; 678c8b0795cSMark F. Adams gnodes[iSwapIndex] = gnodes[Ii]; 679c8b0795cSMark F. Adams gnodes[Ii] = iTemp; 680c8b0795cSMark F. Adams bIndexSet[Ii] = PETSC_TRUE; 681c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 682c8b0795cSMark F. Adams } 683c8b0795cSMark F. Adams } 684c8b0795cSMark F. Adams ierr = PetscFree( bIndexSet ); CHKERRQ(ierr); 685c8b0795cSMark F. Adams } 686c8b0795cSMark F. Adams /* only sort locals */ 6870cbbd2e1SMark F. Adams qsort( gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare ); 688c8b0795cSMark F. Adams /* create IS of permutation */ 689c8b0795cSMark F. Adams for(kk=0;kk<nloc;kk++) { /* locals only */ 690c8b0795cSMark F. Adams permute[kk] = gnodes[kk].lid; 691c8b0795cSMark F. Adams } 6920cbbd2e1SMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm); 693c8b0795cSMark F. Adams CHKERRQ(ierr); 694c8b0795cSMark F. Adams 695c8b0795cSMark F. Adams ierr = PetscFree( gnodes ); CHKERRQ(ierr); 696c8b0795cSMark F. Adams 697c8b0795cSMark F. Adams /* get MIS aggs */ 698b43b03e9SMark F. Adams 699b43b03e9SMark F. Adams ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr); 700b43b03e9SMark F. Adams ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr); 701b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr); 702b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency( crs, Gmat ); CHKERRQ(ierr); 703b43b03e9SMark F. Adams ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr); 704b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs( crs, PETSC_FALSE ); CHKERRQ(ierr); 705b43b03e9SMark F. Adams ierr = MatCoarsenApply( crs ); CHKERRQ(ierr); 7060cbbd2e1SMark F. Adams ierr = MatCoarsenGetData( crs, a_llist_parent ); CHKERRQ(ierr); 707b43b03e9SMark F. Adams ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr); 708c8b0795cSMark F. Adams 709c8b0795cSMark F. Adams ierr = ISDestroy( &perm ); CHKERRQ(ierr); 7100cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 7110cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr); 7120cbbd2e1SMark F. Adams #endif 713c8b0795cSMark F. Adams PetscFunctionReturn(0); 714c8b0795cSMark F. Adams } 715c8b0795cSMark F. Adams 716c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 717c8b0795cSMark F. Adams /* 7180cbbd2e1SMark F. Adams PCGAMGProlongator_GEO 719c8b0795cSMark F. Adams 720c8b0795cSMark F. Adams Input Parameter: 721c8b0795cSMark F. Adams . pc - this 722c8b0795cSMark F. Adams . Amat - matrix on this fine level 723c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 7240cbbd2e1SMark F. Adams . selected_1 - [nselected] 7250cbbd2e1SMark F. Adams . agg_lists - [nselected] 726c8b0795cSMark F. Adams Output Parameter: 727c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 728c8b0795cSMark F. Adams */ 729c8b0795cSMark F. Adams #undef __FUNCT__ 7300cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_GEO" 7310cbbd2e1SMark F. Adams PetscErrorCode PCGAMGProlongator_GEO( PC pc, 732c8b0795cSMark F. Adams const Mat Amat, 733c8b0795cSMark F. Adams const Mat Gmat, 7340cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists, 735c8b0795cSMark F. Adams Mat *a_P_out 7362e68589bSMark F. Adams ) 7372e68589bSMark F. Adams { 7382e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 7392e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 7402e68589bSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 741c8b0795cSMark F. Adams const PetscInt dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols; 7422e68589bSMark F. Adams PetscErrorCode ierr; 743b43b03e9SMark F. Adams PetscInt Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid; 744c8b0795cSMark F. Adams Mat Prol; 7452e68589bSMark F. Adams PetscMPIInt mype, npe; 7462e68589bSMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat)->comm; 7470cbbd2e1SMark F. Adams IS selected_2,selected_1; 7482e68589bSMark F. Adams const PetscInt *selected_idx; 7492e68589bSMark F. Adams 7502e68589bSMark F. Adams PetscFunctionBegin; 7510cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 7520cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr); 7530cbbd2e1SMark F. Adams #endif 7542e68589bSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 7552e68589bSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 7562e68589bSMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 7572e68589bSMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 758b43b03e9SMark F. Adams nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0); 7592e68589bSMark F. Adams 7602e68589bSMark F. Adams /* get 'nLocalSelected' */ 761*41b27cdeSMark F. Adams ierr = PetscCDGetMIS( agg_lists, &selected_1 ); CHKERRQ(ierr); 762b43b03e9SMark F. Adams ierr = ISGetSize( selected_1, &jj ); CHKERRQ(ierr); 763b43b03e9SMark F. Adams ierr = PetscMalloc( jj*sizeof(PetscInt), &clid_flid ); CHKERRQ(ierr); 7642e68589bSMark F. Adams ierr = ISGetIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 765c8b0795cSMark F. Adams for(kk=0,nLocalSelected=0;kk<jj;kk++) { 7662e68589bSMark F. Adams PetscInt lid = selected_idx[kk]; 767b43b03e9SMark F. Adams if( lid<nloc ) { 7680cbbd2e1SMark F. Adams ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0); CHKERRQ(ierr); 769b43b03e9SMark F. Adams if( ncols>1 ) { /* fiter out singletons */ 770b43b03e9SMark F. Adams clid_flid[nLocalSelected++] = lid; 771b43b03e9SMark F. Adams } 7720cbbd2e1SMark F. Adams else assert(0); /* filtered in coarsening */ 7730cbbd2e1SMark F. Adams ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0); CHKERRQ(ierr); 774b43b03e9SMark F. Adams } 7752e68589bSMark F. Adams } 7762e68589bSMark F. Adams ierr = ISRestoreIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 7772e68589bSMark F. Adams 7782e68589bSMark F. Adams /* create prolongator, create P matrix */ 77969b1f4b7SBarry Smith ierr = MatCreateAIJ(wcomm, 7802e68589bSMark F. Adams nloc*bs, nLocalSelected*bs, 7812e68589bSMark F. Adams PETSC_DETERMINE, PETSC_DETERMINE, 7822e68589bSMark F. Adams 3*data_cols, PETSC_NULL, /* don't have a good way to set this!!! */ 7832e68589bSMark F. Adams 3*data_cols, PETSC_NULL, 7842e68589bSMark F. Adams &Prol ); 7852e68589bSMark F. Adams CHKERRQ(ierr); 7862e68589bSMark F. Adams 787c8b0795cSMark F. Adams /* can get all points "removed" - but not on geomg */ 7882e68589bSMark F. Adams ierr = MatGetSize( Prol, &kk, &jj ); CHKERRQ(ierr); 7892e68589bSMark F. Adams if( jj==0 ) { 7902e68589bSMark F. Adams if( verbose ) { 791*41b27cdeSMark F. Adams PetscPrintf(wcomm,"[%d]%s ERROE: no selected points on coarse grid\n",mype,__FUNCT__); 7922e68589bSMark F. Adams } 793b43b03e9SMark F. Adams ierr = PetscFree( clid_flid ); CHKERRQ(ierr); 7942e68589bSMark F. Adams ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 7952e68589bSMark F. Adams *a_P_out = PETSC_NULL; /* out */ 7962e68589bSMark F. Adams PetscFunctionReturn(0); 7972e68589bSMark F. Adams } 7982e68589bSMark F. Adams 7992e68589bSMark F. Adams { 8002e68589bSMark F. Adams PetscReal *coords; 8012e68589bSMark F. Adams PetscInt nnodes; 8022e68589bSMark F. Adams PetscInt *crsGID; 8032e68589bSMark F. Adams Mat Gmat2; 8042e68589bSMark F. Adams 8052e68589bSMark F. Adams assert(dim==data_cols); 8062e68589bSMark F. Adams /* grow ghost data for better coarse grid cover of fine grid */ 8070cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 8080cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 8092e68589bSMark F. Adams #endif 810b43b03e9SMark F. Adams ierr = getGIDsOnSquareGraph( nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID ); 8112e68589bSMark F. Adams CHKERRQ(ierr); 8122e68589bSMark F. Adams /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */ 8130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 8140cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 8152e68589bSMark F. Adams #endif 8162e68589bSMark F. Adams /* create global vector of coorindates in 'coords' */ 8172e68589bSMark F. Adams if (npe > 1) { 8180cbbd2e1SMark F. Adams ierr = PCGAMGGetDataWithGhosts( Gmat2, dim, pc_gamg->data, &nnodes, &coords ); 8192e68589bSMark F. Adams CHKERRQ(ierr); 8202e68589bSMark F. Adams } 8212e68589bSMark F. Adams else { 822c8b0795cSMark F. Adams coords = (PetscReal*)pc_gamg->data; 8232e68589bSMark F. Adams nnodes = nloc; 8242e68589bSMark F. Adams } 8252e68589bSMark F. Adams ierr = MatDestroy( &Gmat2 ); CHKERRQ(ierr); 8262e68589bSMark F. Adams 8272e68589bSMark F. Adams /* triangulate */ 8282e68589bSMark F. Adams if( dim == 2 ) { 829c8b0795cSMark F. Adams PetscReal metric,tm; 8300cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 8310cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr); 8322e68589bSMark F. Adams #endif 8332e68589bSMark F. Adams ierr = triangulateAndFormProl( selected_2, nnodes, coords, 8340cbbd2e1SMark F. Adams nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric ); 8352e68589bSMark F. Adams CHKERRQ(ierr); 8360cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 8370cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0); CHKERRQ(ierr); 8382e68589bSMark F. Adams #endif 8392e68589bSMark F. Adams ierr = PetscFree( crsGID ); CHKERRQ(ierr); 8402e68589bSMark F. Adams 8412e68589bSMark F. Adams /* clean up and create coordinates for coarse grid (output) */ 8422e68589bSMark F. Adams if (npe > 1) ierr = PetscFree( coords ); CHKERRQ(ierr); 8432e68589bSMark F. Adams 844c8b0795cSMark F. Adams ierr = MPI_Allreduce( &metric, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm ); CHKERRQ(ierr); 845c8b0795cSMark F. Adams if( tm > 1. ) { /* needs to be globalized - should not happen */ 8462e68589bSMark F. Adams if( verbose ) { 847*41b27cdeSMark F. Adams PetscPrintf(wcomm,"[%d]%s failed metric for coarse grid %e\n",mype,__FUNCT__,tm); 8482e68589bSMark F. Adams } 8492e68589bSMark F. Adams ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 8502e68589bSMark F. Adams Prol = PETSC_NULL; 8512e68589bSMark F. Adams } 8522e68589bSMark F. Adams else if( metric > .0 ) { 8532e68589bSMark F. Adams if( verbose ) { 854*41b27cdeSMark F. Adams PetscPrintf(wcomm,"[%d]%s worst metric for coarse grid = %e\n",mype,__FUNCT__,metric); 8552e68589bSMark F. Adams } 8562e68589bSMark F. Adams } 8572e68589bSMark F. Adams } else { 8582e68589bSMark F. Adams SETERRQ(wcomm,PETSC_ERR_LIB,"3D not implemented for 'geo' AMG"); 8592e68589bSMark F. Adams } 8602e68589bSMark F. Adams { /* create next coords - output */ 8612e68589bSMark F. Adams PetscReal *crs_crds; 8622e68589bSMark F. Adams ierr = PetscMalloc( dim*nLocalSelected*sizeof(PetscReal), &crs_crds ); 8632e68589bSMark F. Adams CHKERRQ(ierr); 8642e68589bSMark F. Adams for(kk=0;kk<nLocalSelected;kk++){/* grab local select nodes to promote - output */ 865b43b03e9SMark F. Adams PetscInt lid = clid_flid[kk]; 866c8b0795cSMark F. Adams for(jj=0;jj<dim;jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid]; 8672e68589bSMark F. Adams } 868c8b0795cSMark F. Adams 869c8b0795cSMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 870c8b0795cSMark F. Adams pc_gamg->data = crs_crds; /* out */ 871c8b0795cSMark F. Adams pc_gamg->data_sz = dim*nLocalSelected; 8722e68589bSMark F. Adams } 8732e68589bSMark F. Adams ierr = ISDestroy( &selected_2 ); CHKERRQ(ierr); /* this is selected_1 in serial */ 8742e68589bSMark F. Adams } 8752e68589bSMark F. Adams *a_P_out = Prol; /* out */ 876b43b03e9SMark F. Adams ierr = PetscFree( clid_flid ); CHKERRQ(ierr); 8770cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 8780cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr); 8790cbbd2e1SMark F. Adams #endif 880c8b0795cSMark F. Adams PetscFunctionReturn(0); 881c8b0795cSMark F. Adams } 882c8b0795cSMark F. Adams 883c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 884c8b0795cSMark F. Adams /* 885c8b0795cSMark F. Adams PCCreateGAMG_GEO 886c8b0795cSMark F. Adams 887c8b0795cSMark F. Adams Input Parameter: 888c8b0795cSMark F. Adams . pc - 889c8b0795cSMark F. Adams */ 890c8b0795cSMark F. Adams #undef __FUNCT__ 891c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_GEO" 892c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_GEO( PC pc ) 893c8b0795cSMark F. Adams { 894c8b0795cSMark F. Adams PetscErrorCode ierr; 895c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 896c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 897c8b0795cSMark F. Adams 898c8b0795cSMark F. Adams PetscFunctionBegin; 899c8b0795cSMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GEO; 900c8b0795cSMark F. Adams /* pc->ops->destroy = PCDestroy_GEO; */ 901c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 902c8b0795cSMark F. Adams 903c8b0795cSMark F. Adams /* set internal function pointers */ 904c8b0795cSMark F. Adams pc_gamg->graph = PCGAMGgraph_GEO; 905c8b0795cSMark F. Adams pc_gamg->coarsen = PCGAMGcoarsen_GEO; 9060cbbd2e1SMark F. Adams pc_gamg->prolongator = PCGAMGProlongator_GEO; 907c8b0795cSMark F. Adams pc_gamg->optprol = 0; 908c8b0795cSMark F. Adams 909c8b0795cSMark F. Adams pc_gamg->createdefaultdata = PCSetData_GEO; 910c8b0795cSMark F. Adams 911c8b0795cSMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 912c8b0795cSMark F. Adams "PCSetCoordinates_C", 913c8b0795cSMark F. Adams "PCSetCoordinates_GEO", 914c8b0795cSMark F. Adams PCSetCoordinates_GEO);CHKERRQ(ierr); 9152e68589bSMark F. Adams 9162e68589bSMark F. Adams PetscFunctionReturn(0); 9172e68589bSMark F. Adams } 918