xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision 302f38e8b74e96239826bae8d5bb9abb65139735)
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"
36*302f38e8SMark F. Adams PetscErrorCode PCSetCoordinates_GEO( PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords )
372e68589bSMark F. Adams {
382e68589bSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
392e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
402e68589bSMark F. Adams   PetscErrorCode ierr;
412e68589bSMark F. Adams   PetscInt       arrsz,bs,my0,kk,ii,nloc,Iend;
422e68589bSMark F. Adams   Mat            Amat = pc->pmat;
432e68589bSMark F. Adams 
442e68589bSMark F. Adams   PetscFunctionBegin;
452e68589bSMark F. Adams   PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 );
462e68589bSMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
472e68589bSMark F. Adams   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
482e68589bSMark F. Adams   nloc = (Iend-my0)/bs;
49*302f38e8SMark F. Adams   if(nloc!=a_nloc)SETERRQ2(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Stokes not supported nloc = %d %d.",a_nloc,nloc);
502e68589bSMark F. Adams   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
512e68589bSMark F. Adams 
52c8b0795cSMark F. Adams   pc_gamg->data_cell_rows = 1;
532e68589bSMark F. Adams   if( coords==0 && nloc > 0 ) {
542e68589bSMark F. Adams     SETERRQ(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
552e68589bSMark F. Adams   }
56c8b0795cSMark F. Adams   pc_gamg->data_cell_cols = ndm; /* coordinates */
572e68589bSMark F. Adams 
58c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
592e68589bSMark F. Adams 
602e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
612e68589bSMark F. Adams   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
622e68589bSMark F. Adams     ierr = PetscFree( pc_gamg->data );  CHKERRQ(ierr);
632e68589bSMark F. Adams     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr);
642e68589bSMark F. Adams   }
652e68589bSMark F. Adams   for(kk=0;kk<arrsz;kk++)pc_gamg->data[kk] = -999.;
662e68589bSMark F. Adams   pc_gamg->data[arrsz] = -99.;
672e68589bSMark F. Adams   /* copy data in - column oriented */
682e68589bSMark F. Adams   for( kk = 0 ; kk < nloc ; kk++ ){
692e68589bSMark F. Adams     for( ii = 0 ; ii < ndm ; ii++ ) {
702e68589bSMark F. Adams       pc_gamg->data[ii*nloc + kk] =  coords[kk*ndm + ii];
712e68589bSMark F. Adams     }
722e68589bSMark F. Adams   }
732e68589bSMark F. Adams   assert(pc_gamg->data[arrsz] == -99.);
742e68589bSMark F. Adams 
752e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
762e68589bSMark F. Adams 
772e68589bSMark F. Adams   PetscFunctionReturn(0);
782e68589bSMark F. Adams }
792e68589bSMark F. Adams EXTERN_C_END
802e68589bSMark F. Adams 
812e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
822e68589bSMark F. Adams /*
832e68589bSMark F. Adams    PCSetData_GEO
842e68589bSMark F. Adams 
852e68589bSMark F. Adams   Input Parameter:
862e68589bSMark F. Adams    . pc -
872e68589bSMark F. Adams */
882e68589bSMark F. Adams #undef __FUNCT__
892e68589bSMark F. Adams #define __FUNCT__ "PCSetData_GEO"
90b8cd405aSMark F. Adams PetscErrorCode PCSetData_GEO( PC pc, Mat m )
912e68589bSMark F. Adams {
922e68589bSMark F. Adams   PetscFunctionBegin;
932e68589bSMark F. Adams   SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"GEO MG needs coordinates");
942e68589bSMark F. Adams }
952e68589bSMark F. Adams 
962e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
972e68589bSMark F. Adams /*
982e68589bSMark F. Adams    PCSetFromOptions_GEO
992e68589bSMark F. Adams 
1002e68589bSMark F. Adams   Input Parameter:
1012e68589bSMark F. Adams    . pc -
1022e68589bSMark F. Adams */
1032e68589bSMark F. Adams #undef __FUNCT__
1042e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GEO"
1052e68589bSMark F. Adams PetscErrorCode PCSetFromOptions_GEO( PC pc )
1062e68589bSMark F. Adams {
1072e68589bSMark F. Adams   PetscErrorCode  ierr;
1082e68589bSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1092e68589bSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1102e68589bSMark F. Adams 
1112e68589bSMark F. Adams   PetscFunctionBegin;
1122e68589bSMark F. Adams   ierr = PetscOptionsHead("GAMG-GEO options"); CHKERRQ(ierr);
1132e68589bSMark F. Adams   {
1142e68589bSMark F. Adams     /* -pc_gamg_sa_nsmooths */
1152e68589bSMark F. Adams     /* pc_gamg_sa->smooths = 0; */
1162e68589bSMark F. Adams     /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */
1172e68589bSMark F. Adams     /*                        "smoothing steps for smoothed aggregation, usually 1 (0)", */
1182e68589bSMark F. Adams     /*                        "PCGAMGSetNSmooths_AGG", */
1192e68589bSMark F. Adams     /*                        pc_gamg_sa->smooths, */
1202e68589bSMark F. Adams     /*                        &pc_gamg_sa->smooths, */
1212e68589bSMark F. Adams     /*                        &flag);  */
1222e68589bSMark F. Adams     /* CHKERRQ(ierr); */
1232e68589bSMark F. Adams   }
1242e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1252e68589bSMark F. Adams 
1262e68589bSMark F. Adams   /* call base class */
1272e68589bSMark F. Adams   ierr = PCSetFromOptions_GAMG( pc ); CHKERRQ(ierr);
1282e68589bSMark F. Adams 
1292e68589bSMark F. Adams   if( pc_gamg->verbose ) {
13041b27cdeSMark F. Adams     MPI_Comm  wcomm = ((PetscObject)pc)->comm;
13141b27cdeSMark F. Adams     PetscPrintf(wcomm,"[%d]%s done\n",0,__FUNCT__);
1322e68589bSMark F. Adams   }
1332e68589bSMark F. Adams 
1342e68589bSMark F. Adams   PetscFunctionReturn(0);
1352e68589bSMark F. Adams }
1362e68589bSMark F. Adams 
1372e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1382e68589bSMark F. Adams /*
1392e68589bSMark F. Adams  triangulateAndFormProl
1402e68589bSMark F. Adams 
1412e68589bSMark F. Adams    Input Parameter:
1422e68589bSMark F. Adams    . selected_2 - list of selected local ID, includes selected ghosts
1432e68589bSMark F. Adams    . nnodes -
1442e68589bSMark F. Adams    . coords[2*nnodes] - column vector of local coordinates w/ ghosts
145b43b03e9SMark F. Adams    . nselected_1 - selected IDs that go with base (1) graph
1460cbbd2e1SMark F. Adams    . clid_lid_1[nselected_1] - lids of selected (c) nodes   ???????????
1470cbbd2e1SMark F. Adams    . agg_lists_1 - list of aggregates
1482e68589bSMark F. Adams    . crsGID[selected.size()] - global index for prolongation operator
1492e68589bSMark F. Adams    . bs - block size
1502e68589bSMark F. Adams   Output Parameter:
1512e68589bSMark F. Adams    . a_Prol - prolongation operator
1522e68589bSMark F. Adams    . a_worst_best - measure of worst missed fine vertex, 0 is no misses
1532e68589bSMark F. Adams */
1542e68589bSMark F. Adams #undef __FUNCT__
1552e68589bSMark F. Adams #define __FUNCT__ "triangulateAndFormProl"
1560cbbd2e1SMark F. Adams static PetscErrorCode triangulateAndFormProl( IS  selected_2, /* list of selected local ID, includes selected ghosts */
1572e68589bSMark F. Adams                                               const PetscInt nnodes,
1582e68589bSMark F. Adams                                               const PetscReal coords[], /* column vector of local coordinates w/ ghosts */
159b43b03e9SMark F. Adams                                               const PetscInt nselected_1, /* list of selected local ID, includes selected ghosts */
160b43b03e9SMark F. Adams                                               const PetscInt clid_lid_1[],
1610cbbd2e1SMark F. Adams                                               const PetscCoarsenData *agg_lists_1, /* selected_1 vertices of aggregate unselected vertices */
1622e68589bSMark F. Adams                                               const PetscInt crsGID[],
1632e68589bSMark F. Adams                                               const PetscInt bs,
1642e68589bSMark F. Adams                                               Mat a_Prol, /* prolongation operator (output) */
1652e68589bSMark F. Adams                                               PetscReal *a_worst_best /* measure of worst missed fine vertex, 0 is no misses */
1662e68589bSMark F. Adams                                               )
1672e68589bSMark F. Adams {
1682e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE)
1692e68589bSMark F. Adams   PetscErrorCode       ierr;
170b43b03e9SMark F. Adams   PetscInt             jj,tid,tt,idx,nselected_2;
1712e68589bSMark F. Adams   struct triangulateio in,mid;
1720cbbd2e1SMark F. Adams   const PetscInt      *selected_idx_2;
1732e68589bSMark F. Adams   PetscMPIInt          mype,npe;
1742e68589bSMark F. Adams   PetscInt             Istart,Iend,nFineLoc,myFine0;
1752e68589bSMark F. Adams   int                  kk,nPlotPts,sid;
176c8b0795cSMark F. Adams   MPI_Comm             wcomm = ((PetscObject)a_Prol)->comm;
177c8b0795cSMark F. Adams   PetscReal            tm;
1782e68589bSMark F. Adams   PetscFunctionBegin;
179c8b0795cSMark F. Adams 
180c8b0795cSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);    CHKERRQ(ierr);
181c8b0795cSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);     CHKERRQ(ierr);
182c8b0795cSMark F. Adams   ierr = ISGetSize( selected_2, &nselected_2 );        CHKERRQ(ierr);
1832e68589bSMark F. Adams   if(nselected_2 == 1 || nselected_2 == 2 ){ /* 0 happens on idle processors */
1842e68589bSMark F. Adams     *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */
185c8b0795cSMark F. Adams   }
186c8b0795cSMark F. Adams   else *a_worst_best = 0.0;
187c8b0795cSMark F. Adams   ierr = MPI_Allreduce( a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm );  CHKERRQ(ierr);
188c8b0795cSMark F. Adams   if( tm > 0.0 ) {
189c8b0795cSMark F. Adams     *a_worst_best = 100.0;
1902e68589bSMark F. Adams     PetscFunctionReturn(0);
1912e68589bSMark F. Adams   }
1922e68589bSMark F. Adams   ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend );  CHKERRQ(ierr);
1932e68589bSMark F. Adams   nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs;
1942e68589bSMark F. Adams   nPlotPts = nFineLoc; /* locals */
1952e68589bSMark F. Adams   /* traingle */
1962e68589bSMark F. Adams   /* Define input points - in*/
1972e68589bSMark F. Adams   in.numberofpoints = nselected_2;
1982e68589bSMark F. Adams   in.numberofpointattributes = 0;
1992e68589bSMark F. Adams   /* get nselected points */
2002e68589bSMark F. Adams   ierr = PetscMalloc( 2*(nselected_2)*sizeof(REAL), &in.pointlist ); CHKERRQ(ierr);
2012e68589bSMark F. Adams   ierr = ISGetIndices( selected_2, &selected_idx_2 );     CHKERRQ(ierr);
2022e68589bSMark F. Adams 
2032e68589bSMark F. Adams   for(kk=0,sid=0;kk<nselected_2;kk++,sid += 2){
2042e68589bSMark F. Adams     PetscInt lid = selected_idx_2[kk];
2052e68589bSMark F. Adams     in.pointlist[sid] = coords[lid];
2062e68589bSMark F. Adams     in.pointlist[sid+1] = coords[nnodes + lid];
2072e68589bSMark F. Adams     if(lid>=nFineLoc) nPlotPts++;
2082e68589bSMark F. Adams   }
2092e68589bSMark F. Adams   assert(sid==2*nselected_2);
2102e68589bSMark F. Adams 
2112e68589bSMark F. Adams   in.numberofsegments = 0;
2122e68589bSMark F. Adams   in.numberofedges = 0;
2132e68589bSMark F. Adams   in.numberofholes = 0;
2142e68589bSMark F. Adams   in.numberofregions = 0;
2152e68589bSMark F. Adams   in.trianglelist = 0;
2162e68589bSMark F. Adams   in.segmentmarkerlist = 0;
2172e68589bSMark F. Adams   in.pointattributelist = 0;
2182e68589bSMark F. Adams   in.pointmarkerlist = 0;
2192e68589bSMark F. Adams   in.triangleattributelist = 0;
2202e68589bSMark F. Adams   in.trianglearealist = 0;
2212e68589bSMark F. Adams   in.segmentlist = 0;
2222e68589bSMark F. Adams   in.holelist = 0;
2232e68589bSMark F. Adams   in.regionlist = 0;
2242e68589bSMark F. Adams   in.edgelist = 0;
2252e68589bSMark F. Adams   in.edgemarkerlist = 0;
2262e68589bSMark F. Adams   in.normlist = 0;
2272e68589bSMark F. Adams   /* triangulate */
2282e68589bSMark F. Adams   mid.pointlist = 0;            /* Not needed if -N switch used. */
2292e68589bSMark F. Adams   /* Not needed if -N switch used or number of point attributes is zero: */
2302e68589bSMark F. Adams   mid.pointattributelist = 0;
2312e68589bSMark F. Adams   mid.pointmarkerlist = 0; /* Not needed if -N or -B switch used. */
2322e68589bSMark F. Adams   mid.trianglelist = 0;          /* Not needed if -E switch used. */
2332e68589bSMark F. Adams   /* Not needed if -E switch used or number of triangle attributes is zero: */
2342e68589bSMark F. Adams   mid.triangleattributelist = 0;
2352e68589bSMark F. Adams   mid.neighborlist = 0;         /* Needed only if -n switch used. */
2362e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P not used: */
2372e68589bSMark F. Adams   mid.segmentlist = 0;
2382e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
2392e68589bSMark F. Adams   mid.segmentmarkerlist = 0;
2402e68589bSMark F. Adams   mid.edgelist = 0;             /* Needed only if -e switch used. */
2412e68589bSMark F. Adams   mid.edgemarkerlist = 0;   /* Needed if -e used and -B not used. */
2422e68589bSMark F. Adams   mid.numberoftriangles = 0;
2432e68589bSMark F. Adams 
2442e68589bSMark F. Adams   /* Triangulate the points.  Switches are chosen to read and write a  */
2452e68589bSMark F. Adams   /*   PSLG (p), preserve the convex hull (c), number everything from  */
2462e68589bSMark F. Adams   /*   zero (z), assign a regional attribute to each element (A), and  */
2472e68589bSMark F. Adams   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
2482e68589bSMark F. Adams   /*   neighbor list (n).                                            */
2492e68589bSMark F. Adams   if(nselected_2 != 0){ /* inactive processor */
2502e68589bSMark F. Adams     char args[] = "npczQ"; /* c is needed ? */
2512e68589bSMark F. Adams     triangulate(args, &in, &mid, (struct triangulateio *) NULL );
2522e68589bSMark F. Adams     /* output .poly files for 'showme' */
2532e68589bSMark F. Adams     if( !PETSC_TRUE ) {
2542e68589bSMark F. Adams       static int level = 1;
2552e68589bSMark F. Adams       FILE *file; char fname[32];
2562e68589bSMark F. Adams 
2572e68589bSMark F. Adams       sprintf(fname,"C%d_%d.poly",level,mype); file = fopen(fname, "w");
2582e68589bSMark F. Adams       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
2592e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0);
2602e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2612e68589bSMark F. Adams       for(kk=0,sid=0;kk<in.numberofpoints;kk++,sid += 2){
2622e68589bSMark F. Adams         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
2632e68589bSMark F. Adams       }
2642e68589bSMark F. Adams       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
2652e68589bSMark F. Adams       fprintf(file, "%d  %d\n",0,0);
2662e68589bSMark F. Adams       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
2672e68589bSMark F. Adams       /* One line: <# of holes> */
2682e68589bSMark F. Adams       fprintf(file, "%d\n",0);
2692e68589bSMark F. Adams       /* Following lines: <hole #> <x> <y> */
2702e68589bSMark F. Adams       /* Optional line: <# of regional attributes and/or area constraints> */
2712e68589bSMark F. Adams       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
2722e68589bSMark F. Adams       fclose(file);
2732e68589bSMark F. Adams 
2742e68589bSMark F. Adams       /* elems */
2752e68589bSMark F. Adams       sprintf(fname,"C%d_%d.ele",level,mype); file = fopen(fname, "w");
2762e68589bSMark F. Adams       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
2772e68589bSMark F. Adams       fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0);
2782e68589bSMark F. Adams       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
2792e68589bSMark F. Adams       for(kk=0,sid=0;kk<mid.numberoftriangles;kk++,sid += 3){
2802e68589bSMark F. Adams         fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]);
2812e68589bSMark F. Adams       }
2822e68589bSMark F. Adams       fclose(file);
2832e68589bSMark F. Adams 
2842e68589bSMark F. Adams       sprintf(fname,"C%d_%d.node",level,mype); file = fopen(fname, "w");
2852e68589bSMark F. Adams       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
2862e68589bSMark F. Adams       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
2872e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n",nPlotPts,2,0,0);
2882e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2892e68589bSMark F. Adams       for(kk=0,sid=0;kk<in.numberofpoints;kk++,sid+=2){
2902e68589bSMark F. Adams         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
2912e68589bSMark F. Adams       }
2922e68589bSMark F. Adams 
2932e68589bSMark F. Adams       sid /= 2;
2942e68589bSMark F. Adams       for(jj=0;jj<nFineLoc;jj++){
2952e68589bSMark F. Adams         PetscBool sel = PETSC_TRUE;
2962e68589bSMark F. Adams         for( kk=0 ; kk<nselected_2 && sel ; kk++ ){
2972e68589bSMark F. Adams           PetscInt lid = selected_idx_2[kk];
2982e68589bSMark F. Adams           if( lid == jj ) sel = PETSC_FALSE;
2992e68589bSMark F. Adams         }
3002e68589bSMark F. Adams         if( sel ) {
3012e68589bSMark F. Adams           fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[nnodes + jj]);
3022e68589bSMark F. Adams         }
3032e68589bSMark F. Adams       }
3042e68589bSMark F. Adams       fclose(file);
3052e68589bSMark F. Adams       assert(sid==nPlotPts);
3062e68589bSMark F. Adams       level++;
3072e68589bSMark F. Adams     }
3082e68589bSMark F. Adams   }
3090cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3100cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
3112e68589bSMark F. Adams #endif
3122e68589bSMark F. Adams   { /* form P - setup some maps */
3130cbbd2e1SMark F. Adams     PetscInt clid,mm,*nTri,*node_tri;
3142e68589bSMark F. Adams 
3152e68589bSMark F. Adams     ierr = PetscMalloc( nselected_2*sizeof(PetscInt), &node_tri ); CHKERRQ(ierr);
3162e68589bSMark F. Adams     ierr = PetscMalloc( nselected_2*sizeof(PetscInt), &nTri ); CHKERRQ(ierr);
3172e68589bSMark F. Adams 
3182e68589bSMark F. Adams     /* need list of triangles on node */
3192e68589bSMark F. Adams     for(kk=0;kk<nselected_2;kk++) nTri[kk] = 0;
3202e68589bSMark F. Adams     for(tid=0,kk=0;tid<mid.numberoftriangles;tid++){
3212e68589bSMark F. Adams       for(jj=0;jj<3;jj++) {
3222e68589bSMark F. Adams         PetscInt cid = mid.trianglelist[kk++];
3232e68589bSMark F. Adams         if( nTri[cid] == 0 ) node_tri[cid] = tid;
3242e68589bSMark F. Adams         nTri[cid]++;
3252e68589bSMark F. Adams       }
3262e68589bSMark F. Adams     }
3272e68589bSMark F. Adams #define EPS 1.e-12
3282e68589bSMark F. Adams     /* find points and set prolongation */
3290cbbd2e1SMark F. Adams     for( mm = clid = 0 ; mm < nFineLoc ; mm++ ){
330e78576d6SMark F. Adams       PetscBool ise;
331e78576d6SMark F. Adams       ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise); CHKERRQ(ierr);
332e78576d6SMark F. Adams       if( !ise ) {
3330cbbd2e1SMark F. Adams         const PetscInt lid = mm;
3340cbbd2e1SMark F. Adams         //for(clid_iterator=0;clid_iterator<nselected_1;clid_iterator++){
3350cbbd2e1SMark F. Adams         //PetscInt flid = clid_lid_1[clid_iterator]; assert(flid != -1);
3362e68589bSMark F. Adams         PetscScalar AA[3][3];
3372e68589bSMark F. Adams         PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO;
33841b27cdeSMark F. Adams         PetscCDPos         pos;
339e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos); CHKERRQ(ierr);
340e78576d6SMark F. Adams         while(pos){
341e78576d6SMark F. Adams           PetscInt flid;
342ffc955d6SMark F. Adams           ierr = PetscLLNGetID( pos, &flid ); CHKERRQ(ierr);
343e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos); CHKERRQ(ierr);
3440cbbd2e1SMark F. Adams 
3452e68589bSMark F. Adams           if( flid < nFineLoc ) {  /* could be a ghost */
3462e68589bSMark F. Adams             PetscInt bestTID = -1; PetscReal best_alpha = 1.e10;
3472e68589bSMark F. Adams             const PetscInt fgid = flid + myFine0;
3482e68589bSMark F. Adams             /* compute shape function for gid */
3492e68589bSMark F. Adams             const PetscReal fcoord[3] = {coords[flid],coords[nnodes+flid],1.0};
3502e68589bSMark F. Adams             PetscBool haveit=PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3];
3512e68589bSMark F. Adams             /* look for it */
3520cbbd2e1SMark F. Adams             for( tid = node_tri[clid], jj=0;
3532e68589bSMark F. Adams                  jj < 5 && !haveit && tid != -1;
3542e68589bSMark F. Adams                  jj++ ){
3552e68589bSMark F. Adams               for(tt=0;tt<3;tt++){
3562e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3*tid + tt];
3572e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
3582e68589bSMark F. Adams                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0;
3592e68589bSMark F. Adams                 clids[tt] = cid2; /* store for interp */
3602e68589bSMark F. Adams               }
3612e68589bSMark F. Adams 
3622e68589bSMark F. Adams               for(tt=0;tt<3;tt++) alpha[tt] = (PetscScalar)fcoord[tt];
3632e68589bSMark F. Adams 
3642e68589bSMark F. Adams               /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */
3652e68589bSMark F. Adams               LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
3662e68589bSMark F. Adams               {
3672e68589bSMark F. Adams                 PetscBool have=PETSC_TRUE;  PetscReal lowest=1.e10;
3682e68589bSMark F. Adams                 for( tt = 0, idx = 0 ; tt < 3 ; tt++ ) {
3692e68589bSMark F. Adams                   if( PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS ) have = PETSC_FALSE;
3702e68589bSMark F. Adams                   if( PetscRealPart(alpha[tt]) < lowest ){
3712e68589bSMark F. Adams                     lowest = PetscRealPart(alpha[tt]);
3722e68589bSMark F. Adams                     idx = tt;
3732e68589bSMark F. Adams                   }
3742e68589bSMark F. Adams                 }
3752e68589bSMark F. Adams                 haveit = have;
3762e68589bSMark F. Adams               }
3772e68589bSMark F. Adams               tid = mid.neighborlist[3*tid + idx];
3782e68589bSMark F. Adams             }
3792e68589bSMark F. Adams 
3802e68589bSMark F. Adams             if( !haveit ) {
3812e68589bSMark F. Adams               /* brute force */
3822e68589bSMark F. Adams               for(tid=0 ; tid<mid.numberoftriangles && !haveit ; tid++ ){
3832e68589bSMark F. Adams                 for(tt=0;tt<3;tt++){
3842e68589bSMark F. Adams                   PetscInt cid2 = mid.trianglelist[3*tid + tt];
3852e68589bSMark F. Adams                   PetscInt lid2 = selected_idx_2[cid2];
3862e68589bSMark F. Adams                   AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0;
3872e68589bSMark F. Adams                   clids[tt] = cid2; /* store for interp */
3882e68589bSMark F. Adams                 }
3892e68589bSMark F. Adams                 for(tt=0;tt<3;tt++) alpha[tt] = fcoord[tt];
3902e68589bSMark F. Adams                 /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */
3912e68589bSMark F. Adams                 LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
3922e68589bSMark F. Adams                 {
3932e68589bSMark F. Adams                   PetscBool have=PETSC_TRUE;  PetscReal worst=0.0, v;
3942e68589bSMark F. Adams                   for(tt=0; tt<3 && have ;tt++) {
3952e68589bSMark F. Adams                     if( PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS ) have=PETSC_FALSE;
3962e68589bSMark F. Adams                     if( (v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst ) worst = v;
3972e68589bSMark F. Adams                   }
3982e68589bSMark F. Adams                   if( worst < best_alpha ) {
3992e68589bSMark F. Adams                     best_alpha = worst; bestTID = tid;
4002e68589bSMark F. Adams                   }
4012e68589bSMark F. Adams                   haveit = have;
4022e68589bSMark F. Adams                 }
4032e68589bSMark F. Adams               }
4042e68589bSMark F. Adams             }
4052e68589bSMark F. Adams             if( !haveit ) {
4062e68589bSMark F. Adams               if( best_alpha > *a_worst_best ) *a_worst_best = best_alpha;
4072e68589bSMark F. Adams               /* use best one */
4082e68589bSMark F. Adams               for(tt=0;tt<3;tt++){
4092e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3*bestTID + tt];
4102e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
4112e68589bSMark F. Adams                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0;
4122e68589bSMark F. Adams                 clids[tt] = cid2; /* store for interp */
4132e68589bSMark F. Adams               }
4142e68589bSMark F. Adams               for(tt=0;tt<3;tt++) alpha[tt] = fcoord[tt];
4152e68589bSMark F. Adams               /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */
4162e68589bSMark F. Adams               LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
4172e68589bSMark F. Adams             }
4182e68589bSMark F. Adams 
4192e68589bSMark F. Adams             /* put in row of P */
4202e68589bSMark F. Adams             for(idx=0;idx<3;idx++){
4212e68589bSMark F. Adams               PetscScalar shp = alpha[idx];
4222e68589bSMark F. Adams               if( PetscAbs(PetscRealPart(shp)) > 1.e-6 ) {
4232e68589bSMark F. Adams                 PetscInt cgid = crsGID[clids[idx]];
4242e68589bSMark F. Adams                 PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */
4252e68589bSMark F. Adams                 for(tt=0 ; tt < bs ; tt++, ii++, jj++ ){
4262e68589bSMark F. Adams                   ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES); CHKERRQ(ierr);
4272e68589bSMark F. Adams                 }
4282e68589bSMark F. Adams               }
4292e68589bSMark F. Adams             }
4302e68589bSMark F. Adams           }
4310cbbd2e1SMark F. Adams         } /* aggregates iterations */
4320cbbd2e1SMark F. Adams         clid++;
4330cbbd2e1SMark F. Adams       } /* a coarse agg */
4340cbbd2e1SMark F. Adams     } /* for all fine nodes */
4350cbbd2e1SMark F. Adams 
4362e68589bSMark F. Adams     ierr = ISRestoreIndices( selected_2, &selected_idx_2 );     CHKERRQ(ierr);
4372e68589bSMark F. Adams     ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4382e68589bSMark F. Adams     ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4392e68589bSMark F. Adams 
4402e68589bSMark F. Adams     ierr = PetscFree( node_tri );  CHKERRQ(ierr);
4412e68589bSMark F. Adams     ierr = PetscFree( nTri );  CHKERRQ(ierr);
4422e68589bSMark F. Adams   }
4430cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4440cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
4452e68589bSMark F. Adams #endif
4462e68589bSMark F. Adams   free( mid.trianglelist );
4472e68589bSMark F. Adams   free( mid.neighborlist );
4482e68589bSMark F. Adams   ierr = PetscFree( in.pointlist );  CHKERRQ(ierr);
4492e68589bSMark F. Adams 
4502e68589bSMark F. Adams   PetscFunctionReturn(0);
4512e68589bSMark F. Adams #else
4522e68589bSMark F. Adams   SETERRQ(((PetscObject)a_Prol)->comm,PETSC_ERR_LIB,"configure with TRIANGLE to use geometric MG");
4532e68589bSMark F. Adams #endif
4542e68589bSMark F. Adams }
4552e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
4562e68589bSMark F. Adams /*
4572e68589bSMark F. Adams    getGIDsOnSquareGraph - square graph, get
4582e68589bSMark F. Adams 
4592e68589bSMark F. Adams    Input Parameter:
4600cbbd2e1SMark F. Adams    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
461b43b03e9SMark F. Adams    . clid_lid_1 - [nselected_1] lids of selected nodes
4622e68589bSMark F. Adams    . Gmat1 - graph that goes with 'selected_1'
4632e68589bSMark F. Adams    Output Parameter:
4642e68589bSMark F. Adams    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
4652e68589bSMark F. Adams    . a_Gmat_2 - graph that is squared of 'Gmat_1'
4662e68589bSMark F. Adams    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
4672e68589bSMark F. Adams */
4682e68589bSMark F. Adams #undef __FUNCT__
4692e68589bSMark F. Adams #define __FUNCT__ "getGIDsOnSquareGraph"
4700cbbd2e1SMark F. Adams static PetscErrorCode getGIDsOnSquareGraph( const PetscInt nselected_1,
471b43b03e9SMark F. Adams                                             const PetscInt clid_lid_1[],
4722e68589bSMark F. Adams                                             const Mat Gmat1,
4732e68589bSMark F. Adams                                             IS *a_selected_2,
4742e68589bSMark F. Adams                                             Mat *a_Gmat_2,
4752e68589bSMark F. Adams                                             PetscInt **a_crsGID
4762e68589bSMark F. Adams                                             )
4772e68589bSMark F. Adams {
4782e68589bSMark F. Adams   PetscErrorCode ierr;
4792e68589bSMark F. Adams   PetscMPIInt    mype,npe;
480b43b03e9SMark F. Adams   PetscInt       *crsGID, kk,my0,Iend,nloc;
4812e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Gmat1)->comm;
4822e68589bSMark F. Adams 
4832e68589bSMark F. Adams   PetscFunctionBegin;
4842e68589bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
4852e68589bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
4862e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend); CHKERRQ(ierr); /* AIJ */
4872e68589bSMark F. Adams   nloc = Iend - my0; /* this does not change */
4882e68589bSMark F. Adams 
4892e68589bSMark F. Adams   if (npe == 1) { /* not much to do in serial */
490b43b03e9SMark F. Adams     ierr = PetscMalloc( nselected_1*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr);
491b43b03e9SMark F. Adams     for(kk=0;kk<nselected_1;kk++) crsGID[kk] = kk;
4922e68589bSMark F. Adams     *a_Gmat_2 = 0;
493b43b03e9SMark F. Adams     ierr = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);
494b43b03e9SMark F. Adams     CHKERRQ(ierr);
4952e68589bSMark F. Adams   }
4962e68589bSMark F. Adams   else {
497b43b03e9SMark F. Adams     PetscInt      idx,num_fine_ghosts,num_crs_ghost,myCrs0;
4982e68589bSMark F. Adams     Mat_MPIAIJ   *mpimat2;
4992e68589bSMark F. Adams     Mat           Gmat2;
5002e68589bSMark F. Adams     Vec           locState;
5012e68589bSMark F. Adams     PetscScalar   *cpcol_state;
5022e68589bSMark F. Adams 
5032e68589bSMark F. Adams     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
504b43b03e9SMark F. Adams     kk = nselected_1;
505b43b03e9SMark F. Adams     MPI_Scan( &kk, &myCrs0, 1, MPIU_INT, MPIU_SUM, wcomm );
506b43b03e9SMark F. Adams     myCrs0 -= nselected_1;
5072e68589bSMark F. Adams 
508b43b03e9SMark F. Adams     if( a_Gmat_2 ) { /* output */
5092e68589bSMark F. Adams       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
5102e68589bSMark F. Adams       ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 );   CHKERRQ(ierr);
5112e68589bSMark F. Adams       *a_Gmat_2 = Gmat2; /* output */
5122e68589bSMark F. Adams     }
5132e68589bSMark F. Adams     else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
5142e68589bSMark F. Adams     /* get coarse grid GIDS for selected (locals and ghosts) */
5152e68589bSMark F. Adams     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
5162e68589bSMark F. Adams     ierr = MatGetVecs( Gmat2, &locState, 0 );         CHKERRQ(ierr);
5172e68589bSMark F. Adams     ierr = VecSet( locState, (PetscScalar)(PetscReal)(-1) );  CHKERRQ(ierr); /* set with UNKNOWN state */
518b43b03e9SMark F. Adams     for(kk=0;kk<nselected_1;kk++){
519b43b03e9SMark F. Adams       PetscInt fgid = clid_lid_1[kk] + my0;
5202e68589bSMark F. Adams       PetscScalar v = (PetscScalar)(kk+myCrs0);
5212e68589bSMark F. Adams       ierr = VecSetValues( locState, 1, &fgid, &v, INSERT_VALUES );  CHKERRQ(ierr); /* set with PID */
5222e68589bSMark F. Adams     }
5232e68589bSMark F. Adams     ierr = VecAssemblyBegin( locState ); CHKERRQ(ierr);
5242e68589bSMark F. Adams     ierr = VecAssemblyEnd( locState ); CHKERRQ(ierr);
5252e68589bSMark F. Adams     ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5262e68589bSMark F. Adams     ierr =   VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5272e68589bSMark F. Adams     ierr = VecGetLocalSize( mpimat2->lvec, &num_fine_ghosts ); CHKERRQ(ierr);
5282e68589bSMark F. Adams     ierr = VecGetArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr);
5292e68589bSMark F. Adams     for(kk=0,num_crs_ghost=0;kk<num_fine_ghosts;kk++){
5302e68589bSMark F. Adams       if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ) num_crs_ghost++;
5312e68589bSMark F. Adams     }
532b43b03e9SMark F. Adams     ierr = PetscMalloc( (nselected_1+num_crs_ghost)*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr); /* output */
5332e68589bSMark F. Adams     {
5342e68589bSMark F. Adams       PetscInt *selected_set;
535b43b03e9SMark F. Adams       ierr = PetscMalloc( (nselected_1+num_crs_ghost)*sizeof(PetscInt), &selected_set ); CHKERRQ(ierr);
5362e68589bSMark F. Adams       /* do ghost of 'crsGID' */
537b43b03e9SMark F. Adams       for(kk=0,idx=nselected_1;kk<num_fine_ghosts;kk++){
5382e68589bSMark F. Adams         if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){
5392e68589bSMark F. Adams           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
5402e68589bSMark F. Adams           selected_set[idx] = nloc + kk;
5412e68589bSMark F. Adams           crsGID[idx++] = cgid;
5422e68589bSMark F. Adams         }
5432e68589bSMark F. Adams       }
544b43b03e9SMark F. Adams       assert(idx==(nselected_1+num_crs_ghost));
5452e68589bSMark F. Adams       ierr = VecRestoreArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr);
5462e68589bSMark F. Adams       /* do locals in 'crsGID' */
5472e68589bSMark F. Adams       ierr = VecGetArray( locState, &cpcol_state ); CHKERRQ(ierr);
5482e68589bSMark F. Adams       for(kk=0,idx=0;kk<nloc;kk++){
5492e68589bSMark F. Adams         if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){
5502e68589bSMark F. Adams           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
5512e68589bSMark F. Adams           selected_set[idx] = kk;
5522e68589bSMark F. Adams           crsGID[idx++] = cgid;
5532e68589bSMark F. Adams         }
5542e68589bSMark F. Adams       }
555b43b03e9SMark F. Adams       assert(idx==nselected_1);
5562e68589bSMark F. Adams       ierr = VecRestoreArray( locState, &cpcol_state ); CHKERRQ(ierr);
5572e68589bSMark F. Adams 
5582e68589bSMark F. Adams       if( a_selected_2 != 0 ) { /* output */
5590cbbd2e1SMark F. Adams         ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);
5602e68589bSMark F. Adams         CHKERRQ(ierr);
5612e68589bSMark F. Adams       }
5620cbbd2e1SMark F. Adams       else {
5632e68589bSMark F. Adams         ierr = PetscFree( selected_set );  CHKERRQ(ierr);
5642e68589bSMark F. Adams       }
5650cbbd2e1SMark F. Adams     }
5662e68589bSMark F. Adams     ierr = VecDestroy( &locState );                    CHKERRQ(ierr);
5672e68589bSMark F. Adams   }
5682e68589bSMark F. Adams   *a_crsGID = crsGID; /* output */
5692e68589bSMark F. Adams 
5702e68589bSMark F. Adams   PetscFunctionReturn(0);
5712e68589bSMark F. Adams }
5722e68589bSMark F. Adams 
5732e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
5742e68589bSMark F. Adams /*
575c8b0795cSMark F. Adams    PCGAMGgraph_GEO
5762e68589bSMark F. Adams 
5772e68589bSMark F. Adams   Input Parameter:
5782e68589bSMark F. Adams    . pc - this
5792e68589bSMark F. Adams    . Amat - matrix on this fine level
5802e68589bSMark F. Adams   Output Parameter:
581c8b0795cSMark F. Adams    . a_Gmat
5822e68589bSMark F. Adams */
5832e68589bSMark F. Adams #undef __FUNCT__
584c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_GEO"
585c8b0795cSMark F. Adams PetscErrorCode PCGAMGgraph_GEO( PC pc,
5862e68589bSMark F. Adams                                 const Mat Amat,
587c8b0795cSMark F. Adams                                 Mat *a_Gmat
588c8b0795cSMark F. Adams                                 )
589c8b0795cSMark F. Adams {
590c8b0795cSMark F. Adams   PetscErrorCode ierr;
591c8b0795cSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
592c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
593c8b0795cSMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
594c8b0795cSMark F. Adams   const PetscReal vfilter = pc_gamg->threshold;
595c8b0795cSMark F. Adams   PetscMPIInt    mype,npe;
596c8b0795cSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
597c8b0795cSMark F. Adams   Mat            Gmat;
5980cbbd2e1SMark F. Adams   PetscBool  set,flg,symm;
599c8b0795cSMark F. Adams   PetscFunctionBegin;
6000cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
6010cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr);
6020cbbd2e1SMark F. Adams #endif
603c8b0795cSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
604c8b0795cSMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
605c8b0795cSMark F. Adams 
6060cbbd2e1SMark F. Adams   ierr = MatIsSymmetricKnown(Amat, &set, &flg);        CHKERRQ(ierr);
607263489e9SJed Brown   symm = (PetscBool)!(set && flg);
6080cbbd2e1SMark F. Adams 
6092d7fac45SMark F. Adams   ierr  = PCGAMGCreateGraph( Amat, &Gmat ); CHKERRQ( ierr );
6102d7fac45SMark F. Adams   ierr  = PCGAMGFilterGraph( &Gmat, vfilter, symm, verbose ); CHKERRQ( ierr );
611c8b0795cSMark F. Adams 
612c8b0795cSMark F. Adams   *a_Gmat = Gmat;
6130cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
6140cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr);
6150cbbd2e1SMark F. Adams #endif
616c8b0795cSMark F. Adams   PetscFunctionReturn(0);
617c8b0795cSMark F. Adams }
618c8b0795cSMark F. Adams 
619c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
620c8b0795cSMark F. Adams /*
621c8b0795cSMark F. Adams    PCGAMGcoarsen_GEO
622c8b0795cSMark F. Adams 
623c8b0795cSMark F. Adams   Input Parameter:
624e0940f08SMark F. Adams    . a_pc - this
625e0940f08SMark F. Adams    . a_Gmat - graph
626c8b0795cSMark F. Adams   Output Parameter:
627c8b0795cSMark F. Adams    . a_llist_parent - linked list from selected indices for data locality only
628c8b0795cSMark F. Adams */
629c8b0795cSMark F. Adams #undef __FUNCT__
630c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGcoarsen_GEO"
631e0940f08SMark F. Adams PetscErrorCode PCGAMGcoarsen_GEO( PC a_pc,
632e0940f08SMark F. Adams                                   Mat *a_Gmat,
6330cbbd2e1SMark F. Adams                                   PetscCoarsenData **a_llist_parent
634c8b0795cSMark F. Adams                                   )
635c8b0795cSMark F. Adams {
636c8b0795cSMark F. Adams   PetscErrorCode ierr;
637e0940f08SMark F. Adams   PC_MG          *mg = (PC_MG*)a_pc->data;
638c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
639c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,kk,Ii,ncols;
640c8b0795cSMark F. Adams   PetscMPIInt    mype,npe;
6410cbbd2e1SMark F. Adams   IS             perm;
642c8b0795cSMark F. Adams   GAMGNode *gnodes;
643c8b0795cSMark F. Adams   PetscInt *permute;
644e0940f08SMark F. Adams   Mat       Gmat = *a_Gmat;
645c8b0795cSMark F. Adams   MPI_Comm  wcomm = ((PetscObject)Gmat)->comm;
646b43b03e9SMark F. Adams   MatCoarsen crs;
647c8b0795cSMark F. Adams 
648c8b0795cSMark F. Adams   PetscFunctionBegin;
6490cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
6500cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
6510cbbd2e1SMark F. Adams #endif
652c8b0795cSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
653c8b0795cSMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
654c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange( Gmat, &Istart, &Iend ); CHKERRQ(ierr);
655c8b0795cSMark F. Adams   nloc = (Iend-Istart);
656c8b0795cSMark F. Adams 
657c8b0795cSMark F. Adams   /* create random permutation with sort for geo-mg */
658c8b0795cSMark F. Adams   ierr = PetscMalloc( nloc*sizeof(GAMGNode), &gnodes ); CHKERRQ(ierr);
659c8b0795cSMark F. Adams   ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr);
660c8b0795cSMark F. Adams 
661c8b0795cSMark F. Adams   for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */
662c8b0795cSMark F. Adams     ierr = MatGetRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr);
663c8b0795cSMark F. Adams     {
664c8b0795cSMark F. Adams       PetscInt lid = Ii - Istart;
665c8b0795cSMark F. Adams       gnodes[lid].lid = lid;
666c8b0795cSMark F. Adams       gnodes[lid].degree = ncols;
667c8b0795cSMark F. Adams     }
668c8b0795cSMark F. Adams     ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr);
669c8b0795cSMark F. Adams   }
670c8b0795cSMark F. Adams   /* randomize */
671c8b0795cSMark F. Adams   srand(1); /* make deterministic */
672c8b0795cSMark F. Adams   if( PETSC_TRUE ) {
673c8b0795cSMark F. Adams     PetscBool *bIndexSet;
674c8b0795cSMark F. Adams     ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr);
675c8b0795cSMark F. Adams     for ( Ii = 0; Ii < nloc ; Ii++) bIndexSet[Ii] = PETSC_FALSE;
676c8b0795cSMark F. Adams     for ( Ii = 0; Ii < nloc ; Ii++)
677c8b0795cSMark F. Adams     {
678c8b0795cSMark F. Adams       PetscInt iSwapIndex = rand()%nloc;
679c8b0795cSMark F. Adams       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii)
680c8b0795cSMark F. Adams       {
681c8b0795cSMark F. Adams         GAMGNode iTemp = gnodes[iSwapIndex];
682c8b0795cSMark F. Adams         gnodes[iSwapIndex] = gnodes[Ii];
683c8b0795cSMark F. Adams         gnodes[Ii] = iTemp;
684c8b0795cSMark F. Adams         bIndexSet[Ii] = PETSC_TRUE;
685c8b0795cSMark F. Adams         bIndexSet[iSwapIndex] = PETSC_TRUE;
686c8b0795cSMark F. Adams       }
687c8b0795cSMark F. Adams     }
688c8b0795cSMark F. Adams     ierr = PetscFree( bIndexSet );  CHKERRQ(ierr);
689c8b0795cSMark F. Adams   }
690c8b0795cSMark F. Adams   /* only sort locals */
6910cbbd2e1SMark F. Adams   qsort( gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare );
692c8b0795cSMark F. Adams   /* create IS of permutation */
693c8b0795cSMark F. Adams   for(kk=0;kk<nloc;kk++) { /* locals only */
694c8b0795cSMark F. Adams     permute[kk] = gnodes[kk].lid;
695c8b0795cSMark F. Adams   }
6960cbbd2e1SMark F. Adams   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);
697c8b0795cSMark F. Adams   CHKERRQ(ierr);
698c8b0795cSMark F. Adams 
699c8b0795cSMark F. Adams   ierr = PetscFree( gnodes );  CHKERRQ(ierr);
700c8b0795cSMark F. Adams 
701c8b0795cSMark F. Adams   /* get MIS aggs */
702b43b03e9SMark F. Adams 
703b43b03e9SMark F. Adams   ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr);
704b43b03e9SMark F. Adams   ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr);
705b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr);
706b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency( crs, Gmat ); CHKERRQ(ierr);
707b43b03e9SMark F. Adams   ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr);
708b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs( crs, PETSC_FALSE ); CHKERRQ(ierr);
709b43b03e9SMark F. Adams   ierr = MatCoarsenApply( crs ); CHKERRQ(ierr);
7100cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData( crs, a_llist_parent ); CHKERRQ(ierr);
711b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr);
712c8b0795cSMark F. Adams 
713c8b0795cSMark F. Adams   ierr = ISDestroy( &perm );                    CHKERRQ(ierr);
7140cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
7150cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
7160cbbd2e1SMark F. Adams #endif
717c8b0795cSMark F. Adams   PetscFunctionReturn(0);
718c8b0795cSMark F. Adams }
719c8b0795cSMark F. Adams 
720c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
721c8b0795cSMark F. Adams /*
7220cbbd2e1SMark F. Adams  PCGAMGProlongator_GEO
723c8b0795cSMark F. Adams 
724c8b0795cSMark F. Adams  Input Parameter:
725c8b0795cSMark F. Adams  . pc - this
726c8b0795cSMark F. Adams  . Amat - matrix on this fine level
727c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
7280cbbd2e1SMark F. Adams  . selected_1 - [nselected]
7290cbbd2e1SMark F. Adams  . agg_lists - [nselected]
730c8b0795cSMark F. Adams  Output Parameter:
731c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
732c8b0795cSMark F. Adams  */
733c8b0795cSMark F. Adams #undef __FUNCT__
7340cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_GEO"
7350cbbd2e1SMark F. Adams PetscErrorCode PCGAMGProlongator_GEO( PC pc,
736c8b0795cSMark F. Adams                                       const Mat Amat,
737c8b0795cSMark F. Adams                                       const Mat Gmat,
7380cbbd2e1SMark F. Adams                                       PetscCoarsenData *agg_lists,
739c8b0795cSMark F. Adams                                       Mat *a_P_out
7402e68589bSMark F. Adams                                       )
7412e68589bSMark F. Adams {
7422e68589bSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
7432e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
7442e68589bSMark F. Adams   const PetscInt  verbose = pc_gamg->verbose;
745c8b0795cSMark F. Adams   const PetscInt  dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
7462e68589bSMark F. Adams   PetscErrorCode ierr;
747b43b03e9SMark F. Adams   PetscInt       Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid;
748c8b0795cSMark F. Adams   Mat            Prol;
7492e68589bSMark F. Adams   PetscMPIInt    mype, npe;
7502e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
7510cbbd2e1SMark F. Adams   IS             selected_2,selected_1;
7522e68589bSMark F. Adams   const PetscInt *selected_idx;
7532e68589bSMark F. Adams 
7542e68589bSMark F. Adams   PetscFunctionBegin;
7550cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
7560cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
7570cbbd2e1SMark F. Adams #endif
7582e68589bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
7592e68589bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
7602e68589bSMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
7612e68589bSMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
762b43b03e9SMark F. Adams   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
7632e68589bSMark F. Adams 
7642e68589bSMark F. Adams   /* get 'nLocalSelected' */
76541b27cdeSMark F. Adams   ierr = PetscCDGetMIS( agg_lists, &selected_1 );        CHKERRQ(ierr);
766b43b03e9SMark F. Adams   ierr = ISGetSize( selected_1, &jj );               CHKERRQ(ierr);
767b43b03e9SMark F. Adams   ierr = PetscMalloc( jj*sizeof(PetscInt), &clid_flid ); CHKERRQ(ierr);
7682e68589bSMark F. Adams   ierr = ISGetIndices( selected_1, &selected_idx );     CHKERRQ(ierr);
769c8b0795cSMark F. Adams   for(kk=0,nLocalSelected=0;kk<jj;kk++) {
7702e68589bSMark F. Adams     PetscInt lid = selected_idx[kk];
771b43b03e9SMark F. Adams     if( lid<nloc ) {
7720cbbd2e1SMark F. Adams       ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0); CHKERRQ(ierr);
773b43b03e9SMark F. Adams       if( ncols>1 ) { /* fiter out singletons */
774b43b03e9SMark F. Adams         clid_flid[nLocalSelected++] = lid;
775b43b03e9SMark F. Adams       }
7760cbbd2e1SMark F. Adams       else assert(0); /* filtered in coarsening */
7770cbbd2e1SMark F. Adams       ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0); CHKERRQ(ierr);
778b43b03e9SMark F. Adams     }
7792e68589bSMark F. Adams   }
7802e68589bSMark F. Adams   ierr = ISRestoreIndices( selected_1, &selected_idx );     CHKERRQ(ierr);
7812e68589bSMark F. Adams 
7822e68589bSMark F. Adams   /* create prolongator, create P matrix */
78369b1f4b7SBarry Smith   ierr = MatCreateAIJ(wcomm,
7842e68589bSMark F. Adams                       nloc*bs, nLocalSelected*bs,
7852e68589bSMark F. Adams                       PETSC_DETERMINE, PETSC_DETERMINE,
7862e68589bSMark F. Adams                       3*data_cols, PETSC_NULL, /* don't have a good way to set this!!! */
7872e68589bSMark F. Adams                       3*data_cols, PETSC_NULL,
7882e68589bSMark F. Adams                       &Prol );
7892e68589bSMark F. Adams   CHKERRQ(ierr);
7902e68589bSMark F. Adams 
791c8b0795cSMark F. Adams   /* can get all points "removed" - but not on geomg */
7922e68589bSMark F. Adams   ierr =  MatGetSize( Prol, &kk, &jj ); CHKERRQ(ierr);
7932e68589bSMark F. Adams   if( jj==0 ) {
7942e68589bSMark F. Adams     if( verbose ) {
79541b27cdeSMark F. Adams       PetscPrintf(wcomm,"[%d]%s ERROE: no selected points on coarse grid\n",mype,__FUNCT__);
7962e68589bSMark F. Adams     }
797b43b03e9SMark F. Adams     ierr = PetscFree( clid_flid );  CHKERRQ(ierr);
7982e68589bSMark F. Adams     ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
7992e68589bSMark F. Adams     *a_P_out = PETSC_NULL;  /* out */
8002e68589bSMark F. Adams     PetscFunctionReturn(0);
8012e68589bSMark F. Adams   }
8022e68589bSMark F. Adams 
8032e68589bSMark F. Adams   {
8042e68589bSMark F. Adams     PetscReal *coords;
8052e68589bSMark F. Adams     PetscInt nnodes;
806ffc955d6SMark F. Adams     PetscInt  *crsGID = PETSC_NULL;
8072e68589bSMark F. Adams     Mat        Gmat2;
8082e68589bSMark F. Adams 
8092e68589bSMark F. Adams     assert(dim==data_cols);
8102e68589bSMark F. Adams     /* grow ghost data for better coarse grid cover of fine grid */
8110cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
8120cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
8132e68589bSMark F. Adams #endif
814b43b03e9SMark F. Adams     ierr = getGIDsOnSquareGraph( nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID );
8152e68589bSMark F. Adams     CHKERRQ(ierr);
8162e68589bSMark F. Adams     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
8170cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
8180cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
8192e68589bSMark F. Adams #endif
8202e68589bSMark F. Adams     /* create global vector of coorindates in 'coords' */
8212e68589bSMark F. Adams     if (npe > 1) {
8220cbbd2e1SMark F. Adams       ierr = PCGAMGGetDataWithGhosts( Gmat2, dim, pc_gamg->data, &nnodes, &coords );
8232e68589bSMark F. Adams       CHKERRQ(ierr);
8242e68589bSMark F. Adams     }
8252e68589bSMark F. Adams     else {
826c8b0795cSMark F. Adams       coords = (PetscReal*)pc_gamg->data;
8272e68589bSMark F. Adams       nnodes = nloc;
8282e68589bSMark F. Adams     }
8292e68589bSMark F. Adams     ierr = MatDestroy( &Gmat2 );  CHKERRQ(ierr);
8302e68589bSMark F. Adams 
8312e68589bSMark F. Adams     /* triangulate */
8322e68589bSMark F. Adams     if( dim == 2 ) {
833c8b0795cSMark F. Adams       PetscReal metric,tm;
8340cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
8350cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
8362e68589bSMark F. Adams #endif
8372e68589bSMark F. Adams       ierr = triangulateAndFormProl( selected_2, nnodes, coords,
8380cbbd2e1SMark F. Adams                                      nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric );
8392e68589bSMark F. Adams       CHKERRQ(ierr);
8400cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
8410cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0); CHKERRQ(ierr);
8422e68589bSMark F. Adams #endif
8432e68589bSMark F. Adams       ierr = PetscFree( crsGID );  CHKERRQ(ierr);
8442e68589bSMark F. Adams 
8452e68589bSMark F. Adams       /* clean up and create coordinates for coarse grid (output) */
8462e68589bSMark F. Adams       if (npe > 1) ierr = PetscFree( coords ); CHKERRQ(ierr);
8472e68589bSMark F. Adams 
848c8b0795cSMark F. Adams       ierr = MPI_Allreduce( &metric, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm );  CHKERRQ(ierr);
849c8b0795cSMark F. Adams       if( tm > 1. ) { /* needs to be globalized - should not happen */
8502e68589bSMark F. Adams         if( verbose ) {
85141b27cdeSMark F. Adams           PetscPrintf(wcomm,"[%d]%s failed metric for coarse grid %e\n",mype,__FUNCT__,tm);
8522e68589bSMark F. Adams         }
8532e68589bSMark F. Adams         ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
8542e68589bSMark F. Adams         Prol = PETSC_NULL;
8552e68589bSMark F. Adams       }
8562e68589bSMark F. Adams       else if( metric > .0 ) {
8572e68589bSMark F. Adams         if( verbose ) {
85841b27cdeSMark F. Adams           PetscPrintf(wcomm,"[%d]%s worst metric for coarse grid = %e\n",mype,__FUNCT__,metric);
8592e68589bSMark F. Adams         }
8602e68589bSMark F. Adams       }
8612e68589bSMark F. Adams     } else {
8622e68589bSMark F. Adams       SETERRQ(wcomm,PETSC_ERR_LIB,"3D not implemented for 'geo' AMG");
8632e68589bSMark F. Adams     }
8642e68589bSMark F. Adams     { /* create next coords - output */
8652e68589bSMark F. Adams       PetscReal *crs_crds;
8662e68589bSMark F. Adams       ierr = PetscMalloc( dim*nLocalSelected*sizeof(PetscReal), &crs_crds );
8672e68589bSMark F. Adams       CHKERRQ(ierr);
8682e68589bSMark F. Adams       for(kk=0;kk<nLocalSelected;kk++){/* grab local select nodes to promote - output */
869b43b03e9SMark F. Adams         PetscInt lid = clid_flid[kk];
870c8b0795cSMark F. Adams         for(jj=0;jj<dim;jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
8712e68589bSMark F. Adams       }
872c8b0795cSMark F. Adams 
873c8b0795cSMark F. Adams       ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
874c8b0795cSMark F. Adams       pc_gamg->data = crs_crds; /* out */
875c8b0795cSMark F. Adams       pc_gamg->data_sz = dim*nLocalSelected;
8762e68589bSMark F. Adams     }
8772e68589bSMark F. Adams     ierr = ISDestroy( &selected_2 ); CHKERRQ(ierr); /* this is selected_1 in serial */
8782e68589bSMark F. Adams   }
8792e68589bSMark F. Adams   *a_P_out = Prol;  /* out */
880b43b03e9SMark F. Adams   ierr = PetscFree( clid_flid );  CHKERRQ(ierr);
8810cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
8820cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
8830cbbd2e1SMark F. Adams #endif
884c8b0795cSMark F. Adams   PetscFunctionReturn(0);
885c8b0795cSMark F. Adams }
886c8b0795cSMark F. Adams 
887c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
888c8b0795cSMark F. Adams /*
889c8b0795cSMark F. Adams  PCCreateGAMG_GEO
890c8b0795cSMark F. Adams 
891c8b0795cSMark F. Adams   Input Parameter:
892c8b0795cSMark F. Adams    . pc -
893c8b0795cSMark F. Adams */
894c8b0795cSMark F. Adams #undef __FUNCT__
895c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_GEO"
896c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_GEO( PC pc )
897c8b0795cSMark F. Adams {
898c8b0795cSMark F. Adams   PetscErrorCode  ierr;
899c8b0795cSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
900c8b0795cSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
901c8b0795cSMark F. Adams 
902c8b0795cSMark F. Adams   PetscFunctionBegin;
903c8b0795cSMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GEO;
904c8b0795cSMark F. Adams   /* pc->ops->destroy        = PCDestroy_GEO; */
905c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
906c8b0795cSMark F. Adams 
907c8b0795cSMark F. Adams   /* set internal function pointers */
908c8b0795cSMark F. Adams   pc_gamg->graph = PCGAMGgraph_GEO;
909c8b0795cSMark F. Adams   pc_gamg->coarsen = PCGAMGcoarsen_GEO;
9100cbbd2e1SMark F. Adams   pc_gamg->prolongator = PCGAMGProlongator_GEO;
911c8b0795cSMark F. Adams   pc_gamg->optprol = 0;
912c8b0795cSMark F. Adams 
913c8b0795cSMark F. Adams   pc_gamg->createdefaultdata = PCSetData_GEO;
914c8b0795cSMark F. Adams 
915c8b0795cSMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
916c8b0795cSMark F. Adams                                             "PCSetCoordinates_C",
917c8b0795cSMark F. Adams                                             "PCSetCoordinates_GEO",
918c8b0795cSMark F. Adams                                             PCSetCoordinates_GEO);CHKERRQ(ierr);
9192e68589bSMark F. Adams 
9202e68589bSMark F. Adams   PetscFunctionReturn(0);
9212e68589bSMark F. Adams }
922