12e68589bSMark F. Adams /* 22e68589bSMark F. Adams GAMG geometric-algebric multiogrid PC - Mark Adams 2011 32e68589bSMark F. Adams */ 42e68589bSMark F. Adams 52e68589bSMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 62e68589bSMark F. Adams #include <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 16*c8b0795cSMark F. Adams /* Private context for the GAMG preconditioner */ 17*c8b0795cSMark F. Adams typedef struct{ 18*c8b0795cSMark F. Adams PetscInt lid; /* local vertex index */ 19*c8b0795cSMark F. Adams PetscInt degree; /* vertex degree */ 20*c8b0795cSMark F. Adams } GAMGNode; 21*c8b0795cSMark F. Adams int compare (const void *a, const void *b) 22*c8b0795cSMark F. Adams { 23*c8b0795cSMark F. Adams return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree); 24*c8b0795cSMark F. Adams } 25*c8b0795cSMark 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 51*c8b0795cSMark 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 } 55*c8b0795cSMark F. Adams pc_gamg->data_cell_cols = ndm; /* coordinates */ 562e68589bSMark F. Adams 57*c8b0795cSMark 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 ) { 1292e68589bSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"[%d]%s done\n",0,__FUNCT__); 1302e68589bSMark F. Adams } 1312e68589bSMark F. Adams 1322e68589bSMark F. Adams PetscFunctionReturn(0); 1332e68589bSMark F. Adams } 1342e68589bSMark F. Adams 1352e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1362e68589bSMark F. Adams /* 1372e68589bSMark F. Adams triangulateAndFormProl 1382e68589bSMark F. Adams 1392e68589bSMark F. Adams Input Parameter: 1402e68589bSMark F. Adams . selected_2 - list of selected local ID, includes selected ghosts 1412e68589bSMark F. Adams . nnodes - 1422e68589bSMark F. Adams . coords[2*nnodes] - column vector of local coordinates w/ ghosts 1432e68589bSMark F. Adams . selected_1 - selected IDs that go with base (1) graph 1442e68589bSMark F. Adams . locals_llist - linked list with (some) locality info of base graph 1452e68589bSMark F. Adams . crsGID[selected.size()] - global index for prolongation operator 1462e68589bSMark F. Adams . bs - block size 1472e68589bSMark F. Adams Output Parameter: 1482e68589bSMark F. Adams . a_Prol - prolongation operator 1492e68589bSMark F. Adams . a_worst_best - measure of worst missed fine vertex, 0 is no misses 1502e68589bSMark F. Adams */ 1512e68589bSMark F. Adams #undef __FUNCT__ 1522e68589bSMark F. Adams #define __FUNCT__ "triangulateAndFormProl" 1532e68589bSMark F. Adams PetscErrorCode triangulateAndFormProl( IS selected_2, /* list of selected local ID, includes selected ghosts */ 1542e68589bSMark F. Adams const PetscInt nnodes, 1552e68589bSMark F. Adams const PetscReal coords[], /* column vector of local coordinates w/ ghosts */ 1562e68589bSMark F. Adams IS selected_1, /* list of selected local ID, includes selected ghosts */ 1572e68589bSMark F. Adams IS locals_llist, /* linked list from selected 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) */ 1612e68589bSMark F. Adams PetscReal *a_worst_best /* measure of worst missed fine vertex, 0 is no misses */ 1622e68589bSMark F. Adams ) 1632e68589bSMark F. Adams { 1642e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE) 1652e68589bSMark F. Adams PetscErrorCode ierr; 1662e68589bSMark F. Adams PetscInt jj,tid,tt,idx,nselected_1,nselected_2; 1672e68589bSMark F. Adams struct triangulateio in,mid; 1682e68589bSMark F. Adams const PetscInt *selected_idx_1,*selected_idx_2,*llist_idx; 1692e68589bSMark F. Adams PetscMPIInt mype,npe; 1702e68589bSMark F. Adams PetscInt Istart,Iend,nFineLoc,myFine0; 1712e68589bSMark F. Adams int kk,nPlotPts,sid; 172*c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Prol)->comm; 173*c8b0795cSMark F. Adams PetscReal tm; 1742e68589bSMark F. Adams PetscFunctionBegin; 175*c8b0795cSMark F. Adams 176*c8b0795cSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype); CHKERRQ(ierr); 177*c8b0795cSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe); CHKERRQ(ierr); 178*c8b0795cSMark F. Adams ierr = ISGetSize( selected_1, &nselected_1 ); CHKERRQ(ierr); 179*c8b0795cSMark F. Adams ierr = ISGetSize( selected_2, &nselected_2 ); CHKERRQ(ierr); 1802e68589bSMark F. Adams if(nselected_2 == 1 || nselected_2 == 2 ){ /* 0 happens on idle processors */ 1812e68589bSMark F. Adams *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */ 182*c8b0795cSMark F. Adams } 183*c8b0795cSMark F. Adams else *a_worst_best = 0.0; 184*c8b0795cSMark F. Adams ierr = MPI_Allreduce( a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm ); CHKERRQ(ierr); 185*c8b0795cSMark F. Adams if( tm > 0.0 ) { 186*c8b0795cSMark F. Adams *a_worst_best = 100.0; 1872e68589bSMark F. Adams PetscFunctionReturn(0); 1882e68589bSMark F. Adams } 1892e68589bSMark F. Adams ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend ); CHKERRQ(ierr); 1902e68589bSMark F. Adams nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs; 1912e68589bSMark F. Adams nPlotPts = nFineLoc; /* locals */ 1922e68589bSMark F. Adams /* traingle */ 1932e68589bSMark F. Adams /* Define input points - in*/ 1942e68589bSMark F. Adams in.numberofpoints = nselected_2; 1952e68589bSMark F. Adams in.numberofpointattributes = 0; 1962e68589bSMark F. Adams /* get nselected points */ 1972e68589bSMark F. Adams ierr = PetscMalloc( 2*(nselected_2)*sizeof(REAL), &in.pointlist ); CHKERRQ(ierr); 1982e68589bSMark F. Adams ierr = ISGetIndices( selected_2, &selected_idx_2 ); CHKERRQ(ierr); 1992e68589bSMark F. Adams 2002e68589bSMark F. Adams for(kk=0,sid=0;kk<nselected_2;kk++,sid += 2){ 2012e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 2022e68589bSMark F. Adams in.pointlist[sid] = coords[lid]; 2032e68589bSMark F. Adams in.pointlist[sid+1] = coords[nnodes + lid]; 2042e68589bSMark F. Adams if(lid>=nFineLoc) nPlotPts++; 2052e68589bSMark F. Adams } 2062e68589bSMark F. Adams assert(sid==2*nselected_2); 2072e68589bSMark F. Adams 2082e68589bSMark F. Adams in.numberofsegments = 0; 2092e68589bSMark F. Adams in.numberofedges = 0; 2102e68589bSMark F. Adams in.numberofholes = 0; 2112e68589bSMark F. Adams in.numberofregions = 0; 2122e68589bSMark F. Adams in.trianglelist = 0; 2132e68589bSMark F. Adams in.segmentmarkerlist = 0; 2142e68589bSMark F. Adams in.pointattributelist = 0; 2152e68589bSMark F. Adams in.pointmarkerlist = 0; 2162e68589bSMark F. Adams in.triangleattributelist = 0; 2172e68589bSMark F. Adams in.trianglearealist = 0; 2182e68589bSMark F. Adams in.segmentlist = 0; 2192e68589bSMark F. Adams in.holelist = 0; 2202e68589bSMark F. Adams in.regionlist = 0; 2212e68589bSMark F. Adams in.edgelist = 0; 2222e68589bSMark F. Adams in.edgemarkerlist = 0; 2232e68589bSMark F. Adams in.normlist = 0; 2242e68589bSMark F. Adams /* triangulate */ 2252e68589bSMark F. Adams mid.pointlist = 0; /* Not needed if -N switch used. */ 2262e68589bSMark F. Adams /* Not needed if -N switch used or number of point attributes is zero: */ 2272e68589bSMark F. Adams mid.pointattributelist = 0; 2282e68589bSMark F. Adams mid.pointmarkerlist = 0; /* Not needed if -N or -B switch used. */ 2292e68589bSMark F. Adams mid.trianglelist = 0; /* Not needed if -E switch used. */ 2302e68589bSMark F. Adams /* Not needed if -E switch used or number of triangle attributes is zero: */ 2312e68589bSMark F. Adams mid.triangleattributelist = 0; 2322e68589bSMark F. Adams mid.neighborlist = 0; /* Needed only if -n switch used. */ 2332e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P not used: */ 2342e68589bSMark F. Adams mid.segmentlist = 0; 2352e68589bSMark F. Adams /* Needed only if segments are output (-p or -c) and -P and -B not used: */ 2362e68589bSMark F. Adams mid.segmentmarkerlist = 0; 2372e68589bSMark F. Adams mid.edgelist = 0; /* Needed only if -e switch used. */ 2382e68589bSMark F. Adams mid.edgemarkerlist = 0; /* Needed if -e used and -B not used. */ 2392e68589bSMark F. Adams mid.numberoftriangles = 0; 2402e68589bSMark F. Adams 2412e68589bSMark F. Adams /* Triangulate the points. Switches are chosen to read and write a */ 2422e68589bSMark F. Adams /* PSLG (p), preserve the convex hull (c), number everything from */ 2432e68589bSMark F. Adams /* zero (z), assign a regional attribute to each element (A), and */ 2442e68589bSMark F. Adams /* produce an edge list (e), a Voronoi diagram (v), and a triangle */ 2452e68589bSMark F. Adams /* neighbor list (n). */ 2462e68589bSMark F. Adams if(nselected_2 != 0){ /* inactive processor */ 2472e68589bSMark F. Adams char args[] = "npczQ"; /* c is needed ? */ 2482e68589bSMark F. Adams triangulate(args, &in, &mid, (struct triangulateio *) NULL ); 2492e68589bSMark F. Adams /* output .poly files for 'showme' */ 2502e68589bSMark F. Adams if( !PETSC_TRUE ) { 2512e68589bSMark F. Adams static int level = 1; 2522e68589bSMark F. Adams FILE *file; char fname[32]; 2532e68589bSMark F. Adams 2542e68589bSMark F. Adams sprintf(fname,"C%d_%d.poly",level,mype); file = fopen(fname, "w"); 2552e68589bSMark F. Adams /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/ 2562e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); 2572e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2582e68589bSMark F. Adams for(kk=0,sid=0;kk<in.numberofpoints;kk++,sid += 2){ 2592e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2602e68589bSMark F. Adams } 2612e68589bSMark F. Adams /*One line: <# of segments> <# of boundary markers (0 or 1)> */ 2622e68589bSMark F. Adams fprintf(file, "%d %d\n",0,0); 2632e68589bSMark F. Adams /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */ 2642e68589bSMark F. Adams /* One line: <# of holes> */ 2652e68589bSMark F. Adams fprintf(file, "%d\n",0); 2662e68589bSMark F. Adams /* Following lines: <hole #> <x> <y> */ 2672e68589bSMark F. Adams /* Optional line: <# of regional attributes and/or area constraints> */ 2682e68589bSMark F. Adams /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */ 2692e68589bSMark F. Adams fclose(file); 2702e68589bSMark F. Adams 2712e68589bSMark F. Adams /* elems */ 2722e68589bSMark F. Adams sprintf(fname,"C%d_%d.ele",level,mype); file = fopen(fname, "w"); 2732e68589bSMark F. Adams /* First line: <# of triangles> <nodes per triangle> <# of attributes> */ 2742e68589bSMark F. Adams fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0); 2752e68589bSMark F. Adams /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */ 2762e68589bSMark F. Adams for(kk=0,sid=0;kk<mid.numberoftriangles;kk++,sid += 3){ 2772e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]); 2782e68589bSMark F. Adams } 2792e68589bSMark F. Adams fclose(file); 2802e68589bSMark F. Adams 2812e68589bSMark F. Adams sprintf(fname,"C%d_%d.node",level,mype); file = fopen(fname, "w"); 2822e68589bSMark F. Adams /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */ 2832e68589bSMark F. Adams /* fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); */ 2842e68589bSMark F. Adams fprintf(file, "%d %d %d %d\n",nPlotPts,2,0,0); 2852e68589bSMark F. Adams /*Following lines: <vertex #> <x> <y> */ 2862e68589bSMark F. Adams for(kk=0,sid=0;kk<in.numberofpoints;kk++,sid+=2){ 2872e68589bSMark F. Adams fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 2882e68589bSMark F. Adams } 2892e68589bSMark F. Adams 2902e68589bSMark F. Adams sid /= 2; 2912e68589bSMark F. Adams for(jj=0;jj<nFineLoc;jj++){ 2922e68589bSMark F. Adams PetscBool sel = PETSC_TRUE; 2932e68589bSMark F. Adams for( kk=0 ; kk<nselected_2 && sel ; kk++ ){ 2942e68589bSMark F. Adams PetscInt lid = selected_idx_2[kk]; 2952e68589bSMark F. Adams if( lid == jj ) sel = PETSC_FALSE; 2962e68589bSMark F. Adams } 2972e68589bSMark F. Adams if( sel ) { 2982e68589bSMark F. Adams fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[nnodes + jj]); 2992e68589bSMark F. Adams } 3002e68589bSMark F. Adams } 3012e68589bSMark F. Adams fclose(file); 3022e68589bSMark F. Adams assert(sid==nPlotPts); 3032e68589bSMark F. Adams level++; 3042e68589bSMark F. Adams } 3052e68589bSMark F. Adams } 3062e68589bSMark F. Adams #if defined PETSC_USE_LOG 3072e68589bSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 3082e68589bSMark F. Adams #endif 3092e68589bSMark F. Adams { /* form P - setup some maps */ 3102e68589bSMark F. Adams PetscInt clid_iterator; 3112e68589bSMark F. Adams PetscInt *nTri, *node_tri; 3122e68589bSMark F. Adams 3132e68589bSMark F. Adams ierr = PetscMalloc( nselected_2*sizeof(PetscInt), &node_tri ); CHKERRQ(ierr); 3142e68589bSMark F. Adams ierr = PetscMalloc( nselected_2*sizeof(PetscInt), &nTri ); CHKERRQ(ierr); 3152e68589bSMark F. Adams 3162e68589bSMark F. Adams /* need list of triangles on node */ 3172e68589bSMark F. Adams for(kk=0;kk<nselected_2;kk++) nTri[kk] = 0; 3182e68589bSMark F. Adams for(tid=0,kk=0;tid<mid.numberoftriangles;tid++){ 3192e68589bSMark F. Adams for(jj=0;jj<3;jj++) { 3202e68589bSMark F. Adams PetscInt cid = mid.trianglelist[kk++]; 3212e68589bSMark F. Adams if( nTri[cid] == 0 ) node_tri[cid] = tid; 3222e68589bSMark F. Adams nTri[cid]++; 3232e68589bSMark F. Adams } 3242e68589bSMark F. Adams } 3252e68589bSMark F. Adams #define EPS 1.e-12 3262e68589bSMark F. Adams /* find points and set prolongation */ 3272e68589bSMark F. Adams ierr = ISGetIndices( selected_1, &selected_idx_1 ); CHKERRQ(ierr); 3282e68589bSMark F. Adams ierr = ISGetIndices( locals_llist, &llist_idx ); CHKERRQ(ierr); 3292e68589bSMark F. Adams for( clid_iterator = 0 ; clid_iterator < nselected_1 ; clid_iterator++ ){ 3302e68589bSMark F. Adams PetscInt flid = selected_idx_1[clid_iterator]; assert(flid != -1); 3312e68589bSMark F. Adams PetscScalar AA[3][3]; 3322e68589bSMark F. Adams PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO; 3332e68589bSMark F. Adams do{ 3342e68589bSMark F. Adams if( flid < nFineLoc ) { /* could be a ghost */ 3352e68589bSMark F. Adams PetscInt bestTID = -1; PetscReal best_alpha = 1.e10; 3362e68589bSMark F. Adams const PetscInt fgid = flid + myFine0; 3372e68589bSMark F. Adams /* compute shape function for gid */ 3382e68589bSMark F. Adams const PetscReal fcoord[3] = { coords[flid], coords[nnodes + flid], 1.0 }; 3392e68589bSMark F. Adams PetscBool haveit = PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3]; 3402e68589bSMark F. Adams /* look for it */ 3412e68589bSMark F. Adams for( tid = node_tri[clid_iterator], jj=0; 3422e68589bSMark F. Adams jj < 5 && !haveit && tid != -1; 3432e68589bSMark F. Adams jj++ ){ 3442e68589bSMark F. Adams for(tt=0;tt<3;tt++){ 3452e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3462e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 3472e68589bSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0; 3482e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3492e68589bSMark F. Adams } 3502e68589bSMark F. Adams 3512e68589bSMark F. Adams for(tt=0;tt<3;tt++) alpha[tt] = (PetscScalar)fcoord[tt]; 3522e68589bSMark F. Adams 3532e68589bSMark F. Adams /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */ 3542e68589bSMark F. Adams LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO); 3552e68589bSMark F. Adams { 3562e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal lowest=1.e10; 3572e68589bSMark F. Adams for( tt = 0, idx = 0 ; tt < 3 ; tt++ ) { 3582e68589bSMark F. Adams if( PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS ) have = PETSC_FALSE; 3592e68589bSMark F. Adams if( PetscRealPart(alpha[tt]) < lowest ){ 3602e68589bSMark F. Adams lowest = PetscRealPart(alpha[tt]); 3612e68589bSMark F. Adams idx = tt; 3622e68589bSMark F. Adams } 3632e68589bSMark F. Adams } 3642e68589bSMark F. Adams haveit = have; 3652e68589bSMark F. Adams } 3662e68589bSMark F. Adams tid = mid.neighborlist[3*tid + idx]; 3672e68589bSMark F. Adams } 3682e68589bSMark F. Adams 3692e68589bSMark F. Adams if( !haveit ) { 3702e68589bSMark F. Adams /* brute force */ 3712e68589bSMark F. Adams for(tid=0 ; tid<mid.numberoftriangles && !haveit ; tid++ ){ 3722e68589bSMark F. Adams for(tt=0;tt<3;tt++){ 3732e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*tid + tt]; 3742e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 3752e68589bSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0; 3762e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 3772e68589bSMark F. Adams } 3782e68589bSMark F. Adams for(tt=0;tt<3;tt++) alpha[tt] = fcoord[tt]; 3792e68589bSMark F. Adams /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */ 3802e68589bSMark F. Adams LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO); 3812e68589bSMark F. Adams { 3822e68589bSMark F. Adams PetscBool have=PETSC_TRUE; PetscReal worst=0.0, v; 3832e68589bSMark F. Adams for(tt=0; tt<3 && have ;tt++) { 3842e68589bSMark F. Adams if( PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS ) have=PETSC_FALSE; 3852e68589bSMark F. Adams if( (v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst ) worst = v; 3862e68589bSMark F. Adams } 3872e68589bSMark F. Adams if( worst < best_alpha ) { 3882e68589bSMark F. Adams best_alpha = worst; bestTID = tid; 3892e68589bSMark F. Adams } 3902e68589bSMark F. Adams haveit = have; 3912e68589bSMark F. Adams } 3922e68589bSMark F. Adams } 3932e68589bSMark F. Adams } 3942e68589bSMark F. Adams if( !haveit ) { 3952e68589bSMark F. Adams if( best_alpha > *a_worst_best ) *a_worst_best = best_alpha; 3962e68589bSMark F. Adams /* use best one */ 3972e68589bSMark F. Adams for(tt=0;tt<3;tt++){ 3982e68589bSMark F. Adams PetscInt cid2 = mid.trianglelist[3*bestTID + tt]; 3992e68589bSMark F. Adams PetscInt lid2 = selected_idx_2[cid2]; 4002e68589bSMark F. Adams AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0; 4012e68589bSMark F. Adams clids[tt] = cid2; /* store for interp */ 4022e68589bSMark F. Adams } 4032e68589bSMark F. Adams for(tt=0;tt<3;tt++) alpha[tt] = fcoord[tt]; 4042e68589bSMark F. Adams /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */ 4052e68589bSMark F. Adams LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO); 4062e68589bSMark F. Adams } 4072e68589bSMark F. Adams 4082e68589bSMark F. Adams /* put in row of P */ 4092e68589bSMark F. Adams for(idx=0;idx<3;idx++){ 4102e68589bSMark F. Adams PetscScalar shp = alpha[idx]; 4112e68589bSMark F. Adams if( PetscAbs(PetscRealPart(shp)) > 1.e-6 ) { 4122e68589bSMark F. Adams PetscInt cgid = crsGID[clids[idx]]; 4132e68589bSMark F. Adams PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */ 4142e68589bSMark F. Adams for(tt=0 ; tt < bs ; tt++, ii++, jj++ ){ 4152e68589bSMark F. Adams ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES); CHKERRQ(ierr); 4162e68589bSMark F. Adams } 4172e68589bSMark F. Adams } 4182e68589bSMark F. Adams } 4192e68589bSMark F. Adams } /* local vertex test */ 4202e68589bSMark F. Adams } while( (flid=llist_idx[flid]) != -1 ); 4212e68589bSMark F. Adams } 4222e68589bSMark F. Adams ierr = ISRestoreIndices( selected_2, &selected_idx_2 ); CHKERRQ(ierr); 4232e68589bSMark F. Adams ierr = ISRestoreIndices( selected_1, &selected_idx_1 ); CHKERRQ(ierr); 4242e68589bSMark F. Adams ierr = ISRestoreIndices( locals_llist, &llist_idx ); CHKERRQ(ierr); 4252e68589bSMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4262e68589bSMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4272e68589bSMark F. Adams 4282e68589bSMark F. Adams ierr = PetscFree( node_tri ); CHKERRQ(ierr); 4292e68589bSMark F. Adams ierr = PetscFree( nTri ); CHKERRQ(ierr); 4302e68589bSMark F. Adams } 4312e68589bSMark F. Adams #if defined PETSC_USE_LOG 4322e68589bSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 4332e68589bSMark F. Adams #endif 4342e68589bSMark F. Adams free( mid.trianglelist ); 4352e68589bSMark F. Adams free( mid.neighborlist ); 4362e68589bSMark F. Adams ierr = PetscFree( in.pointlist ); CHKERRQ(ierr); 4372e68589bSMark F. Adams 4382e68589bSMark F. Adams PetscFunctionReturn(0); 4392e68589bSMark F. Adams #else 4402e68589bSMark F. Adams SETERRQ(((PetscObject)a_Prol)->comm,PETSC_ERR_LIB,"configure with TRIANGLE to use geometric MG"); 4412e68589bSMark F. Adams #endif 4422e68589bSMark F. Adams } 4432e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 4442e68589bSMark F. Adams /* 4452e68589bSMark F. Adams getGIDsOnSquareGraph - square graph, get 4462e68589bSMark F. Adams 4472e68589bSMark F. Adams Input Parameter: 4482e68589bSMark F. Adams . selected_1 - selected local indices (includes ghosts in input Gmat_1) 4492e68589bSMark F. Adams . Gmat1 - graph that goes with 'selected_1' 4502e68589bSMark F. Adams Output Parameter: 4512e68589bSMark F. Adams . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2) 4522e68589bSMark F. Adams . a_Gmat_2 - graph that is squared of 'Gmat_1' 4532e68589bSMark F. Adams . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes 4542e68589bSMark F. Adams */ 4552e68589bSMark F. Adams #undef __FUNCT__ 4562e68589bSMark F. Adams #define __FUNCT__ "getGIDsOnSquareGraph" 4572e68589bSMark F. Adams PetscErrorCode getGIDsOnSquareGraph( const IS selected_1, 4582e68589bSMark F. Adams const Mat Gmat1, 4592e68589bSMark F. Adams IS *a_selected_2, 4602e68589bSMark F. Adams Mat *a_Gmat_2, 4612e68589bSMark F. Adams PetscInt **a_crsGID 4622e68589bSMark F. Adams ) 4632e68589bSMark F. Adams { 4642e68589bSMark F. Adams PetscErrorCode ierr; 4652e68589bSMark F. Adams PetscMPIInt mype,npe; 4662e68589bSMark F. Adams PetscInt *crsGID, kk,my0,Iend,nloc,nSelected_1; 4672e68589bSMark F. Adams const PetscInt *selected_idx; 4682e68589bSMark F. Adams MPI_Comm wcomm = ((PetscObject)Gmat1)->comm; 4692e68589bSMark F. Adams 4702e68589bSMark F. Adams PetscFunctionBegin; 4712e68589bSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 4722e68589bSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 4732e68589bSMark F. Adams ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend); CHKERRQ(ierr); /* AIJ */ 4742e68589bSMark F. Adams nloc = Iend - my0; /* this does not change */ 4752e68589bSMark F. Adams ierr = ISGetLocalSize( selected_1, &nSelected_1 ); CHKERRQ(ierr); 4762e68589bSMark F. Adams 4772e68589bSMark F. Adams if (npe == 1) { /* not much to do in serial */ 4782e68589bSMark F. Adams ierr = PetscMalloc( nSelected_1*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr); 4792e68589bSMark F. Adams for(kk=0;kk<nSelected_1;kk++) crsGID[kk] = kk; 4802e68589bSMark F. Adams *a_Gmat_2 = 0; 4812e68589bSMark F. Adams *a_selected_2 = selected_1; /* needed? */ 4822e68589bSMark F. Adams } 4832e68589bSMark F. Adams else { 4842e68589bSMark F. Adams PetscInt idx,num_fine_ghosts,num_crs_ghost,nLocalSelected,myCrs0; 4852e68589bSMark F. Adams Mat_MPIAIJ *mpimat2; 4862e68589bSMark F. Adams Mat Gmat2; 4872e68589bSMark F. Adams Vec locState; 4882e68589bSMark F. Adams PetscScalar *cpcol_state; 4892e68589bSMark F. Adams 4902e68589bSMark F. Adams /* get 'nLocalSelected' */ 4912e68589bSMark F. Adams ierr = ISGetIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 4922e68589bSMark F. Adams for(kk=0,nLocalSelected=0;kk<nSelected_1;kk++){ 4932e68589bSMark F. Adams PetscInt lid = selected_idx[kk]; 4942e68589bSMark F. Adams if(lid<nloc) nLocalSelected++; 4952e68589bSMark F. Adams } 4962e68589bSMark F. Adams ierr = ISRestoreIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 4972e68589bSMark F. Adams /* scan my coarse zero gid, set 'lid_state' with coarse GID */ 4982e68589bSMark F. Adams MPI_Scan( &nLocalSelected, &myCrs0, 1, MPIU_INT, MPIU_SUM, wcomm ); 4992e68589bSMark F. Adams myCrs0 -= nLocalSelected; 5002e68589bSMark F. Adams 5012e68589bSMark F. Adams if( a_Gmat_2 != 0 ) { /* output */ 5022e68589bSMark F. Adams /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */ 5032e68589bSMark F. Adams ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 ); CHKERRQ(ierr); 5042e68589bSMark F. Adams *a_Gmat_2 = Gmat2; /* output */ 5052e68589bSMark F. Adams } 5062e68589bSMark F. Adams else Gmat2 = Gmat1; /* use local to get crsGIDs at least */ 5072e68589bSMark F. Adams /* get coarse grid GIDS for selected (locals and ghosts) */ 5082e68589bSMark F. Adams mpimat2 = (Mat_MPIAIJ*)Gmat2->data; 5092e68589bSMark F. Adams ierr = MatGetVecs( Gmat2, &locState, 0 ); CHKERRQ(ierr); 5102e68589bSMark F. Adams ierr = VecSet( locState, (PetscScalar)(PetscReal)(-1) ); CHKERRQ(ierr); /* set with UNKNOWN state */ 5112e68589bSMark F. Adams ierr = ISGetIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 5122e68589bSMark F. Adams for(kk=0;kk<nLocalSelected;kk++){ 5132e68589bSMark F. Adams PetscInt fgid = selected_idx[kk] + my0; 5142e68589bSMark F. Adams PetscScalar v = (PetscScalar)(kk+myCrs0); 5152e68589bSMark F. Adams ierr = VecSetValues( locState, 1, &fgid, &v, INSERT_VALUES ); CHKERRQ(ierr); /* set with PID */ 5162e68589bSMark F. Adams } 5172e68589bSMark F. Adams ierr = ISRestoreIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 5182e68589bSMark F. Adams ierr = VecAssemblyBegin( locState ); CHKERRQ(ierr); 5192e68589bSMark F. Adams ierr = VecAssemblyEnd( locState ); CHKERRQ(ierr); 5202e68589bSMark F. Adams ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5212e68589bSMark F. Adams ierr = VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5222e68589bSMark F. Adams ierr = VecGetLocalSize( mpimat2->lvec, &num_fine_ghosts ); CHKERRQ(ierr); 5232e68589bSMark F. Adams ierr = VecGetArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr); 5242e68589bSMark F. Adams for(kk=0,num_crs_ghost=0;kk<num_fine_ghosts;kk++){ 5252e68589bSMark F. Adams if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ) num_crs_ghost++; 5262e68589bSMark F. Adams } 5272e68589bSMark F. Adams ierr = PetscMalloc( (nLocalSelected+num_crs_ghost)*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr); /* output */ 5282e68589bSMark F. Adams { 5292e68589bSMark F. Adams PetscInt *selected_set; 5302e68589bSMark F. Adams ierr = PetscMalloc( (nLocalSelected+num_crs_ghost)*sizeof(PetscInt), &selected_set ); CHKERRQ(ierr); 5312e68589bSMark F. Adams /* do ghost of 'crsGID' */ 5322e68589bSMark F. Adams for(kk=0,idx=nLocalSelected;kk<num_fine_ghosts;kk++){ 5332e68589bSMark F. Adams if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){ 5342e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5352e68589bSMark F. Adams selected_set[idx] = nloc + kk; 5362e68589bSMark F. Adams crsGID[idx++] = cgid; 5372e68589bSMark F. Adams } 5382e68589bSMark F. Adams } 5392e68589bSMark F. Adams assert(idx==(nLocalSelected+num_crs_ghost)); 5402e68589bSMark F. Adams ierr = VecRestoreArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr); 5412e68589bSMark F. Adams /* do locals in 'crsGID' */ 5422e68589bSMark F. Adams ierr = VecGetArray( locState, &cpcol_state ); CHKERRQ(ierr); 5432e68589bSMark F. Adams for(kk=0,idx=0;kk<nloc;kk++){ 5442e68589bSMark F. Adams if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){ 5452e68589bSMark F. Adams PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 5462e68589bSMark F. Adams selected_set[idx] = kk; 5472e68589bSMark F. Adams crsGID[idx++] = cgid; 5482e68589bSMark F. Adams } 5492e68589bSMark F. Adams } 5502e68589bSMark F. Adams assert(idx==nLocalSelected); 5512e68589bSMark F. Adams ierr = VecRestoreArray( locState, &cpcol_state ); CHKERRQ(ierr); 5522e68589bSMark F. Adams 5532e68589bSMark F. Adams if( a_selected_2 != 0 ) { /* output */ 5542e68589bSMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF,(nLocalSelected+num_crs_ghost),selected_set,PETSC_COPY_VALUES,a_selected_2); 5552e68589bSMark F. Adams CHKERRQ(ierr); 5562e68589bSMark F. Adams } 5572e68589bSMark F. Adams ierr = PetscFree( selected_set ); CHKERRQ(ierr); 5582e68589bSMark F. Adams } 5592e68589bSMark F. Adams ierr = VecDestroy( &locState ); CHKERRQ(ierr); 5602e68589bSMark F. Adams } 5612e68589bSMark F. Adams *a_crsGID = crsGID; /* output */ 5622e68589bSMark F. Adams 5632e68589bSMark F. Adams PetscFunctionReturn(0); 5642e68589bSMark F. Adams } 5652e68589bSMark F. Adams 5662e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 5672e68589bSMark F. Adams /* 568*c8b0795cSMark F. Adams PCGAMGgraph_GEO 5692e68589bSMark F. Adams 5702e68589bSMark F. Adams Input Parameter: 5712e68589bSMark F. Adams . pc - this 5722e68589bSMark F. Adams . Amat - matrix on this fine level 5732e68589bSMark F. Adams Output Parameter: 574*c8b0795cSMark F. Adams . a_Gmat 5752e68589bSMark F. Adams */ 5762e68589bSMark F. Adams #undef __FUNCT__ 577*c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_GEO" 578*c8b0795cSMark F. Adams PetscErrorCode PCGAMGgraph_GEO( PC pc, 5792e68589bSMark F. Adams const Mat Amat, 580*c8b0795cSMark F. Adams Mat *a_Gmat 581*c8b0795cSMark F. Adams ) 582*c8b0795cSMark F. Adams { 583*c8b0795cSMark F. Adams PetscErrorCode ierr; 584*c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 585*c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 586*c8b0795cSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 587*c8b0795cSMark F. Adams const PetscReal vfilter = pc_gamg->threshold; 588*c8b0795cSMark F. Adams PetscMPIInt mype,npe; 589*c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat)->comm; 590*c8b0795cSMark F. Adams Mat Gmat; 591*c8b0795cSMark F. Adams 592*c8b0795cSMark F. Adams PetscFunctionBegin; 593*c8b0795cSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 594*c8b0795cSMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 595*c8b0795cSMark F. Adams 596*c8b0795cSMark F. Adams ierr = createSimpleGraph( Amat, &Gmat ); CHKERRQ( ierr ); 597*c8b0795cSMark F. Adams ierr = scaleFilterGraph( &Gmat, vfilter, PETSC_TRUE, verbose ); CHKERRQ( ierr ); 598*c8b0795cSMark F. Adams 599*c8b0795cSMark F. Adams *a_Gmat = Gmat; 600*c8b0795cSMark F. Adams 601*c8b0795cSMark F. Adams PetscFunctionReturn(0); 602*c8b0795cSMark F. Adams } 603*c8b0795cSMark F. Adams 604*c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 605*c8b0795cSMark F. Adams /* 606*c8b0795cSMark F. Adams PCGAMGcoarsen_GEO 607*c8b0795cSMark F. Adams 608*c8b0795cSMark F. Adams Input Parameter: 609*c8b0795cSMark F. Adams . pc - this 610*c8b0795cSMark F. Adams . Gmat - graph 611*c8b0795cSMark F. Adams Output Parameter: 612*c8b0795cSMark F. Adams . a_selected - selected indices (local) 613*c8b0795cSMark F. Adams . a_llist_parent - linked list from selected indices for data locality only 614*c8b0795cSMark F. Adams */ 615*c8b0795cSMark F. Adams #undef __FUNCT__ 616*c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGcoarsen_GEO" 617*c8b0795cSMark F. Adams PetscErrorCode PCGAMGcoarsen_GEO( PC pc, 618*c8b0795cSMark F. Adams const Mat Gmat, 619*c8b0795cSMark F. Adams IS *a_selected, 620*c8b0795cSMark F. Adams IS *a_llist_parent 621*c8b0795cSMark F. Adams ) 622*c8b0795cSMark F. Adams { 623*c8b0795cSMark F. Adams PetscErrorCode ierr; 624*c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 625*c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 626*c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,kk,Ii,ncols; 627*c8b0795cSMark F. Adams PetscMPIInt mype,npe; 628*c8b0795cSMark F. Adams IS perm,selected,llist_parent; 629*c8b0795cSMark F. Adams GAMGNode *gnodes; 630*c8b0795cSMark F. Adams PetscInt *permute; 631*c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)Gmat)->comm; 632*c8b0795cSMark F. Adams 633*c8b0795cSMark F. Adams PetscFunctionBegin; 634*c8b0795cSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 635*c8b0795cSMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 636*c8b0795cSMark F. Adams ierr = MatGetOwnershipRange( Gmat, &Istart, &Iend ); CHKERRQ(ierr); 637*c8b0795cSMark F. Adams nloc = (Iend-Istart); 638*c8b0795cSMark F. Adams 639*c8b0795cSMark F. Adams /* create random permutation with sort for geo-mg */ 640*c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(GAMGNode), &gnodes ); CHKERRQ(ierr); 641*c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr); 642*c8b0795cSMark F. Adams 643*c8b0795cSMark F. Adams for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */ 644*c8b0795cSMark F. Adams ierr = MatGetRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr); 645*c8b0795cSMark F. Adams { 646*c8b0795cSMark F. Adams PetscInt lid = Ii - Istart; 647*c8b0795cSMark F. Adams gnodes[lid].lid = lid; 648*c8b0795cSMark F. Adams gnodes[lid].degree = ncols; 649*c8b0795cSMark F. Adams } 650*c8b0795cSMark F. Adams ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr); 651*c8b0795cSMark F. Adams } 652*c8b0795cSMark F. Adams /* randomize */ 653*c8b0795cSMark F. Adams srand(1); /* make deterministic */ 654*c8b0795cSMark F. Adams if( PETSC_TRUE ) { 655*c8b0795cSMark F. Adams PetscBool *bIndexSet; 656*c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr); 657*c8b0795cSMark F. Adams for ( Ii = 0; Ii < nloc ; Ii++) bIndexSet[Ii] = PETSC_FALSE; 658*c8b0795cSMark F. Adams for ( Ii = 0; Ii < nloc ; Ii++) 659*c8b0795cSMark F. Adams { 660*c8b0795cSMark F. Adams PetscInt iSwapIndex = rand()%nloc; 661*c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) 662*c8b0795cSMark F. Adams { 663*c8b0795cSMark F. Adams GAMGNode iTemp = gnodes[iSwapIndex]; 664*c8b0795cSMark F. Adams gnodes[iSwapIndex] = gnodes[Ii]; 665*c8b0795cSMark F. Adams gnodes[Ii] = iTemp; 666*c8b0795cSMark F. Adams bIndexSet[Ii] = PETSC_TRUE; 667*c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 668*c8b0795cSMark F. Adams } 669*c8b0795cSMark F. Adams } 670*c8b0795cSMark F. Adams ierr = PetscFree( bIndexSet ); CHKERRQ(ierr); 671*c8b0795cSMark F. Adams } 672*c8b0795cSMark F. Adams /* only sort locals */ 673*c8b0795cSMark F. Adams qsort( gnodes, nloc, sizeof(GAMGNode), compare ); 674*c8b0795cSMark F. Adams /* create IS of permutation */ 675*c8b0795cSMark F. Adams for(kk=0;kk<nloc;kk++) { /* locals only */ 676*c8b0795cSMark F. Adams permute[kk] = gnodes[kk].lid; 677*c8b0795cSMark F. Adams } 678*c8b0795cSMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_COPY_VALUES, &perm); 679*c8b0795cSMark F. Adams CHKERRQ(ierr); 680*c8b0795cSMark F. Adams 681*c8b0795cSMark F. Adams ierr = PetscFree( gnodes ); CHKERRQ(ierr); 682*c8b0795cSMark F. Adams ierr = PetscFree( permute ); CHKERRQ(ierr); 683*c8b0795cSMark F. Adams 684*c8b0795cSMark F. Adams /* get MIS aggs */ 685*c8b0795cSMark F. Adams ierr = maxIndSetAgg( perm, Gmat, PETSC_FALSE, pc_gamg->verbose, &selected, &llist_parent ); CHKERRQ(ierr); 686*c8b0795cSMark F. Adams 687*c8b0795cSMark F. Adams ierr = ISDestroy( &perm ); CHKERRQ(ierr); 688*c8b0795cSMark F. Adams 689*c8b0795cSMark F. Adams *a_selected = selected; 690*c8b0795cSMark F. Adams *a_llist_parent = llist_parent; 691*c8b0795cSMark F. Adams 692*c8b0795cSMark F. Adams PetscFunctionReturn(0); 693*c8b0795cSMark F. Adams } 694*c8b0795cSMark F. Adams 695*c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 696*c8b0795cSMark F. Adams /* 697*c8b0795cSMark F. Adams PCGAMGprolongator_GEO 698*c8b0795cSMark F. Adams 699*c8b0795cSMark F. Adams Input Parameter: 700*c8b0795cSMark F. Adams . pc - this 701*c8b0795cSMark F. Adams . Amat - matrix on this fine level 702*c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 703*c8b0795cSMark F. Adams . selected - [nselected inc. chosts] 704*c8b0795cSMark F. Adams . llist_1 - [nloc + Gmat.nghost] linked list 705*c8b0795cSMark F. Adams Output Parameter: 706*c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 707*c8b0795cSMark F. Adams */ 708*c8b0795cSMark F. Adams #undef __FUNCT__ 709*c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGprolongator_GEO" 710*c8b0795cSMark F. Adams PetscErrorCode PCGAMGprolongator_GEO( PC pc, 711*c8b0795cSMark F. Adams const Mat Amat, 712*c8b0795cSMark F. Adams const Mat Gmat, 713*c8b0795cSMark F. Adams IS selected_1, 714*c8b0795cSMark F. Adams IS llist_1, 715*c8b0795cSMark F. Adams Mat *a_P_out 7162e68589bSMark F. Adams ) 7172e68589bSMark F. Adams { 7182e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 7192e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 7202e68589bSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 721*c8b0795cSMark F. Adams const PetscInt dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols; 7222e68589bSMark F. Adams PetscErrorCode ierr; 723*c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,jj,kk,nLocalSelected,bs; 724*c8b0795cSMark F. Adams Mat Prol; 7252e68589bSMark F. Adams PetscMPIInt mype, npe; 7262e68589bSMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat)->comm; 727*c8b0795cSMark F. Adams IS selected_2; 7282e68589bSMark F. Adams const PetscInt *selected_idx; 7292e68589bSMark F. Adams 7302e68589bSMark F. Adams PetscFunctionBegin; 7312e68589bSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 7322e68589bSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 7332e68589bSMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 7342e68589bSMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 735*c8b0795cSMark F. Adams nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs==0); 7362e68589bSMark F. Adams 7372e68589bSMark F. Adams /* get 'nLocalSelected' */ 738*c8b0795cSMark F. Adams ierr = ISGetLocalSize( selected_1, &jj ); CHKERRQ(ierr); 7392e68589bSMark F. Adams ierr = ISGetIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 740*c8b0795cSMark F. Adams for(kk=0,nLocalSelected=0;kk<jj;kk++) { 7412e68589bSMark F. Adams PetscInt lid = selected_idx[kk]; 7422e68589bSMark F. Adams if(lid<nloc) nLocalSelected++; 7432e68589bSMark F. Adams } 7442e68589bSMark F. Adams ierr = ISRestoreIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 7452e68589bSMark F. Adams 7462e68589bSMark F. Adams /* create prolongator, create P matrix */ 7472e68589bSMark F. Adams ierr = MatCreateMPIAIJ(wcomm, 7482e68589bSMark F. Adams nloc*bs, nLocalSelected*bs, 7492e68589bSMark F. Adams PETSC_DETERMINE, PETSC_DETERMINE, 7502e68589bSMark F. Adams 3*data_cols, PETSC_NULL, /* don't have a good way to set this!!! */ 7512e68589bSMark F. Adams 3*data_cols, PETSC_NULL, 7522e68589bSMark F. Adams &Prol ); 7532e68589bSMark F. Adams CHKERRQ(ierr); 7542e68589bSMark F. Adams 755*c8b0795cSMark F. Adams /* can get all points "removed" - but not on geomg */ 7562e68589bSMark F. Adams ierr = MatGetSize( Prol, &kk, &jj ); CHKERRQ(ierr); 7572e68589bSMark F. Adams if( jj==0 ) { 7582e68589bSMark F. Adams if( verbose ) { 759*c8b0795cSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"[%d]%s ERROE: no selected points on coarse grid\n",mype,__FUNCT__); 7602e68589bSMark F. Adams } 7612e68589bSMark F. Adams ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 7622e68589bSMark F. Adams *a_P_out = PETSC_NULL; /* out */ 7632e68589bSMark F. Adams PetscFunctionReturn(0); 7642e68589bSMark F. Adams } 7652e68589bSMark F. Adams 7662e68589bSMark F. Adams { 7672e68589bSMark F. Adams PetscReal *coords; 7682e68589bSMark F. Adams PetscInt nnodes; 7692e68589bSMark F. Adams PetscInt *crsGID; 7702e68589bSMark F. Adams Mat Gmat2; 7712e68589bSMark F. Adams 7722e68589bSMark F. Adams assert(dim==data_cols); 7732e68589bSMark F. Adams /* grow ghost data for better coarse grid cover of fine grid */ 7742e68589bSMark F. Adams #if defined PETSC_USE_LOG 7752e68589bSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 7762e68589bSMark F. Adams #endif 7772e68589bSMark F. Adams ierr = getGIDsOnSquareGraph( selected_1, Gmat, &selected_2, &Gmat2, &crsGID ); 7782e68589bSMark F. Adams CHKERRQ(ierr); 7792e68589bSMark F. Adams /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */ 7802e68589bSMark F. Adams #if defined PETSC_USE_LOG 7812e68589bSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 7822e68589bSMark F. Adams #endif 7832e68589bSMark F. Adams /* create global vector of coorindates in 'coords' */ 7842e68589bSMark F. Adams if (npe > 1) { 785*c8b0795cSMark F. Adams ierr = getDataWithGhosts( Gmat2, dim, pc_gamg->data, &nnodes, &coords ); 7862e68589bSMark F. Adams CHKERRQ(ierr); 7872e68589bSMark F. Adams } 7882e68589bSMark F. Adams else { 789*c8b0795cSMark F. Adams coords = (PetscReal*)pc_gamg->data; 7902e68589bSMark F. Adams nnodes = nloc; 7912e68589bSMark F. Adams } 7922e68589bSMark F. Adams ierr = MatDestroy( &Gmat2 ); CHKERRQ(ierr); 7932e68589bSMark F. Adams 7942e68589bSMark F. Adams /* triangulate */ 7952e68589bSMark F. Adams if( dim == 2 ) { 796*c8b0795cSMark F. Adams PetscReal metric,tm; 7972e68589bSMark F. Adams #if defined PETSC_USE_LOG 7982e68589bSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr); 7992e68589bSMark F. Adams #endif 8002e68589bSMark F. Adams ierr = triangulateAndFormProl( selected_2, nnodes, coords, 8012e68589bSMark F. Adams selected_1, llist_1, crsGID, bs, Prol, &metric ); 8022e68589bSMark F. Adams CHKERRQ(ierr); 8032e68589bSMark F. Adams #if defined PETSC_USE_LOG 8042e68589bSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET6],0,0,0,0); CHKERRQ(ierr); 8052e68589bSMark F. Adams #endif 8062e68589bSMark F. Adams ierr = PetscFree( crsGID ); CHKERRQ(ierr); 8072e68589bSMark F. Adams 8082e68589bSMark F. Adams /* clean up and create coordinates for coarse grid (output) */ 8092e68589bSMark F. Adams if (npe > 1) ierr = PetscFree( coords ); CHKERRQ(ierr); 8102e68589bSMark F. Adams 811*c8b0795cSMark F. Adams ierr = MPI_Allreduce( &metric, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm ); CHKERRQ(ierr); 812*c8b0795cSMark F. Adams if( tm > 1. ) { /* needs to be globalized - should not happen */ 8132e68589bSMark F. Adams if( verbose ) { 814*c8b0795cSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"[%d]%s failed metric for coarse grid %e\n",mype,__FUNCT__,tm); 8152e68589bSMark F. Adams } 8162e68589bSMark F. Adams ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 8172e68589bSMark F. Adams Prol = PETSC_NULL; 8182e68589bSMark F. Adams } 8192e68589bSMark F. Adams else if( metric > .0 ) { 8202e68589bSMark F. Adams if( verbose ) { 821*c8b0795cSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"[%d]%s worst metric for coarse grid = %e\n",mype,__FUNCT__,metric); 8222e68589bSMark F. Adams } 8232e68589bSMark F. Adams } 8242e68589bSMark F. Adams } else { 8252e68589bSMark F. Adams SETERRQ(wcomm,PETSC_ERR_LIB,"3D not implemented for 'geo' AMG"); 8262e68589bSMark F. Adams } 8272e68589bSMark F. Adams { /* create next coords - output */ 8282e68589bSMark F. Adams PetscReal *crs_crds; 8292e68589bSMark F. Adams ierr = PetscMalloc( dim*nLocalSelected*sizeof(PetscReal), &crs_crds ); 8302e68589bSMark F. Adams CHKERRQ(ierr); 8312e68589bSMark F. Adams ierr = ISGetIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 8322e68589bSMark F. Adams for(kk=0;kk<nLocalSelected;kk++){/* grab local select nodes to promote - output */ 8332e68589bSMark F. Adams PetscInt lid = selected_idx[kk]; 834*c8b0795cSMark F. Adams for(jj=0;jj<dim;jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid]; 8352e68589bSMark F. Adams } 8362e68589bSMark F. Adams ierr = ISRestoreIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 837*c8b0795cSMark F. Adams 838*c8b0795cSMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 839*c8b0795cSMark F. Adams pc_gamg->data = crs_crds; /* out */ 840*c8b0795cSMark F. Adams pc_gamg->data_sz = dim*nLocalSelected; 8412e68589bSMark F. Adams } 8422e68589bSMark F. Adams if (npe > 1) { 8432e68589bSMark F. Adams ierr = ISDestroy( &selected_2 ); CHKERRQ(ierr); /* this is selected_1 in serial */ 8442e68589bSMark F. Adams } 8452e68589bSMark F. Adams } 8462e68589bSMark F. Adams *a_P_out = Prol; /* out */ 8472e68589bSMark F. Adams 848*c8b0795cSMark F. Adams PetscFunctionReturn(0); 849*c8b0795cSMark F. Adams } 850*c8b0795cSMark F. Adams 851*c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 852*c8b0795cSMark F. Adams /* 853*c8b0795cSMark F. Adams PCCreateGAMG_GEO 854*c8b0795cSMark F. Adams 855*c8b0795cSMark F. Adams Input Parameter: 856*c8b0795cSMark F. Adams . pc - 857*c8b0795cSMark F. Adams */ 858*c8b0795cSMark F. Adams #undef __FUNCT__ 859*c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_GEO" 860*c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_GEO( PC pc ) 861*c8b0795cSMark F. Adams { 862*c8b0795cSMark F. Adams PetscErrorCode ierr; 863*c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 864*c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 865*c8b0795cSMark F. Adams 866*c8b0795cSMark F. Adams PetscFunctionBegin; 867*c8b0795cSMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GEO; 868*c8b0795cSMark F. Adams /* pc->ops->destroy = PCDestroy_GEO; */ 869*c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 870*c8b0795cSMark F. Adams 871*c8b0795cSMark F. Adams /* set internal function pointers */ 872*c8b0795cSMark F. Adams pc_gamg->graph = PCGAMGgraph_GEO; 873*c8b0795cSMark F. Adams pc_gamg->coarsen = PCGAMGcoarsen_GEO; 874*c8b0795cSMark F. Adams pc_gamg->prolongator = PCGAMGprolongator_GEO; 875*c8b0795cSMark F. Adams pc_gamg->optprol = 0; 876*c8b0795cSMark F. Adams 877*c8b0795cSMark F. Adams pc_gamg->createdefaultdata = PCSetData_GEO; 878*c8b0795cSMark F. Adams 879*c8b0795cSMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 880*c8b0795cSMark F. Adams "PCSetCoordinates_C", 881*c8b0795cSMark F. Adams "PCSetCoordinates_GEO", 882*c8b0795cSMark F. Adams PCSetCoordinates_GEO);CHKERRQ(ierr); 8832e68589bSMark F. Adams 8842e68589bSMark F. Adams PetscFunctionReturn(0); 8852e68589bSMark F. Adams } 886