xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision bf4339c22951fe2ed4cd378ee6bed5e3925f33c6)
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 #include <petscblaslapack.h>
92e68589bSMark F. Adams 
102e68589bSMark F. Adams typedef struct {
11c8b0795cSMark F. Adams   PetscInt  nsmooths;
12c8b0795cSMark F. Adams   PetscBool sym_graph;
13ef4ad70eSMark F. Adams   PetscBool square_graph;
142e68589bSMark F. Adams } PC_GAMG_AGG;
152e68589bSMark F. Adams 
162e68589bSMark F. Adams #undef __FUNCT__
172e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths"
182e68589bSMark F. Adams /*@
192e68589bSMark F. Adams    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical)
202e68589bSMark F. Adams 
212e68589bSMark F. Adams    Not Collective on PC
222e68589bSMark F. Adams 
232e68589bSMark F. Adams    Input Parameters:
242e68589bSMark F. Adams .  pc - the preconditioner context
252e68589bSMark F. Adams 
262e68589bSMark F. Adams    Options Database Key:
272e68589bSMark F. Adams .  -pc_gamg_agg_nsmooths
282e68589bSMark F. Adams 
292e68589bSMark F. Adams    Level: intermediate
302e68589bSMark F. Adams 
312e68589bSMark F. Adams    Concepts: Aggregation AMG preconditioner
322e68589bSMark F. Adams 
332e68589bSMark F. Adams .seealso: ()
342e68589bSMark F. Adams @*/
352e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
362e68589bSMark F. Adams {
372e68589bSMark F. Adams   PetscErrorCode ierr;
382e68589bSMark F. Adams 
392e68589bSMark F. Adams   PetscFunctionBegin;
402e68589bSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
412e68589bSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
422e68589bSMark F. Adams   PetscFunctionReturn(0);
432e68589bSMark F. Adams }
442e68589bSMark F. Adams 
452e68589bSMark F. Adams #undef __FUNCT__
462e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths_GAMG"
472e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths_GAMG(PC pc, PetscInt n)
482e68589bSMark F. Adams {
492e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
502e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
51c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
522e68589bSMark F. Adams 
532e68589bSMark F. Adams   PetscFunctionBegin;
54c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
55c8b0795cSMark F. Adams   PetscFunctionReturn(0);
56c8b0795cSMark F. Adams }
57c8b0795cSMark F. Adams 
58c8b0795cSMark F. Adams #undef __FUNCT__
59c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph"
60c8b0795cSMark F. Adams /*@
61c8b0795cSMark F. Adams    PCGAMGSetSymGraph -
62c8b0795cSMark F. Adams 
63c8b0795cSMark F. Adams    Not Collective on PC
64c8b0795cSMark F. Adams 
65c8b0795cSMark F. Adams    Input Parameters:
66c8b0795cSMark F. Adams .  pc - the preconditioner context
67c8b0795cSMark F. Adams 
68c8b0795cSMark F. Adams    Options Database Key:
69c8b0795cSMark F. Adams .  -pc_gamg_sym_graph
70c8b0795cSMark F. Adams 
71c8b0795cSMark F. Adams    Level: intermediate
72c8b0795cSMark F. Adams 
73c8b0795cSMark F. Adams    Concepts: Aggregation AMG preconditioner
74c8b0795cSMark F. Adams 
75c8b0795cSMark F. Adams .seealso: ()
76c8b0795cSMark F. Adams @*/
77c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n)
78c8b0795cSMark F. Adams {
79c8b0795cSMark F. Adams   PetscErrorCode ierr;
80c8b0795cSMark F. Adams 
81c8b0795cSMark F. Adams   PetscFunctionBegin;
82c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
83c8b0795cSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
84c8b0795cSMark F. Adams   PetscFunctionReturn(0);
85c8b0795cSMark F. Adams }
86c8b0795cSMark F. Adams 
87c8b0795cSMark F. Adams #undef __FUNCT__
88c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph_GAMG"
89c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph_GAMG(PC pc, PetscBool n)
90c8b0795cSMark F. Adams {
91c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
92c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
93c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
94c8b0795cSMark F. Adams 
95c8b0795cSMark F. Adams   PetscFunctionBegin;
96c8b0795cSMark F. Adams   pc_gamg_agg->sym_graph = n;
972e68589bSMark F. Adams   PetscFunctionReturn(0);
982e68589bSMark F. Adams }
992e68589bSMark F. Adams 
100ef4ad70eSMark F. Adams #undef __FUNCT__
101ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph"
102ef4ad70eSMark F. Adams /*@
103ef4ad70eSMark F. Adams    PCGAMGSetSquareGraph -
104ef4ad70eSMark F. Adams 
105ef4ad70eSMark F. Adams    Not Collective on PC
106ef4ad70eSMark F. Adams 
107ef4ad70eSMark F. Adams    Input Parameters:
108ef4ad70eSMark F. Adams .  pc - the preconditioner context
109ef4ad70eSMark F. Adams 
110ef4ad70eSMark F. Adams    Options Database Key:
111ef4ad70eSMark F. Adams .  -pc_gamg_square_graph
112ef4ad70eSMark F. Adams 
113ef4ad70eSMark F. Adams    Level: intermediate
114ef4ad70eSMark F. Adams 
115ef4ad70eSMark F. Adams    Concepts: Aggregation AMG preconditioner
116ef4ad70eSMark F. Adams 
117ef4ad70eSMark F. Adams .seealso: ()
118ef4ad70eSMark F. Adams @*/
119ef4ad70eSMark F. Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscBool n)
120ef4ad70eSMark F. Adams {
121ef4ad70eSMark F. Adams   PetscErrorCode ierr;
122ef4ad70eSMark F. Adams 
123ef4ad70eSMark F. Adams   PetscFunctionBegin;
124ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
125ef4ad70eSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
126ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
127ef4ad70eSMark F. Adams }
128ef4ad70eSMark F. Adams 
129ef4ad70eSMark F. Adams #undef __FUNCT__
130ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph_GAMG"
131ef4ad70eSMark F. Adams PetscErrorCode PCGAMGSetSquareGraph_GAMG(PC pc, PetscBool n)
132ef4ad70eSMark F. Adams {
133ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
134ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
135ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
136ef4ad70eSMark F. Adams 
137ef4ad70eSMark F. Adams   PetscFunctionBegin;
138ef4ad70eSMark F. Adams   pc_gamg_agg->square_graph = n;
139ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
140ef4ad70eSMark F. Adams }
141ef4ad70eSMark F. Adams 
1422e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1432e68589bSMark F. Adams /*
1442e68589bSMark F. Adams    PCSetFromOptions_GAMG_AGG
1452e68589bSMark F. Adams 
1462e68589bSMark F. Adams   Input Parameter:
1472e68589bSMark F. Adams    . pc -
1482e68589bSMark F. Adams */
1492e68589bSMark F. Adams #undef __FUNCT__
1502e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG_AGG"
1518c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptions *PetscOptionsObject,PC pc)
1522e68589bSMark F. Adams {
1532e68589bSMark F. Adams   PetscErrorCode ierr;
1542e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1552e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
156c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1572e68589bSMark F. Adams 
1582e68589bSMark F. Adams   PetscFunctionBegin;
159e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr);
1602e68589bSMark F. Adams   {
1612e68589bSMark F. Adams     /* -pc_gamg_agg_nsmooths */
162d3042614SMark Adams     pc_gamg_agg->nsmooths = 1;
1632fa5cd67SKarl Rupp 
1648afaa268SBarry Smith     ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths","smoothing steps for smoothed aggregation, usually 1","PCGAMGSetNSmooths",pc_gamg_agg->nsmooths,&pc_gamg_agg->nsmooths,NULL);CHKERRQ(ierr);
165c8b0795cSMark F. Adams 
166c8b0795cSMark F. Adams     /* -pc_gamg_sym_graph */
167c8b0795cSMark F. Adams     pc_gamg_agg->sym_graph = PETSC_FALSE;
1682fa5cd67SKarl Rupp 
1698afaa268SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_sym_graph","Set for asymmetric matrices","PCGAMGSetSymGraph",pc_gamg_agg->sym_graph,&pc_gamg_agg->sym_graph,NULL);CHKERRQ(ierr);
170ef4ad70eSMark F. Adams 
171ef4ad70eSMark F. Adams     /* -pc_gamg_square_graph */
172ef4ad70eSMark F. Adams     pc_gamg_agg->square_graph = PETSC_TRUE;
1732fa5cd67SKarl Rupp 
1748afaa268SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_square_graph","For faster coarsening and lower coarse grid complexity","PCGAMGSetSquareGraph",pc_gamg_agg->square_graph,&pc_gamg_agg->square_graph,NULL);CHKERRQ(ierr);
1752e68589bSMark F. Adams   }
1762e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1772e68589bSMark F. Adams   PetscFunctionReturn(0);
1782e68589bSMark F. Adams }
1792e68589bSMark F. Adams 
1802e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1812e68589bSMark F. Adams /*
1822e68589bSMark F. Adams    PCDestroy_AGG
1832e68589bSMark F. Adams 
1842e68589bSMark F. Adams   Input Parameter:
1852e68589bSMark F. Adams    . pc -
1862e68589bSMark F. Adams */
1872e68589bSMark F. Adams #undef __FUNCT__
1889b8ffb57SJed Brown #define __FUNCT__ "PCDestroy_GAMG_AGG"
1899b8ffb57SJed Brown PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
1902e68589bSMark F. Adams {
1912e68589bSMark F. Adams   PetscErrorCode ierr;
1922e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1932e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1942e68589bSMark F. Adams 
1952e68589bSMark F. Adams   PetscFunctionBegin;
1969b8ffb57SJed Brown   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
1973ae0bb68SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
1982e68589bSMark F. Adams   PetscFunctionReturn(0);
1992e68589bSMark F. Adams }
2002e68589bSMark F. Adams 
2012e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2022e68589bSMark F. Adams /*
2032e68589bSMark F. Adams    PCSetCoordinates_AGG
204302f38e8SMark F. Adams      - collective
2052e68589bSMark F. Adams 
2062e68589bSMark F. Adams    Input Parameter:
2072e68589bSMark F. Adams    . pc - the preconditioner context
208a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
209302f38e8SMark F. Adams    . a_nloc - number of vertices local
210302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
2112e68589bSMark F. Adams */
2121e6b0712SBarry Smith 
2132e68589bSMark F. Adams #undef __FUNCT__
2142e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG"
2151e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
2162e68589bSMark F. Adams {
2172e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
2182e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
2192e68589bSMark F. Adams   PetscErrorCode ierr;
22069344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
221a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
2222e68589bSMark F. Adams 
2232e68589bSMark F. Adams   PetscFunctionBegin;
224a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
225a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
226302f38e8SMark F. Adams   nloc = a_nloc;
2272e68589bSMark F. Adams 
2282e68589bSMark F. Adams   /* SA: null space vectors */
22969344418SMark F. Adams   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
23069344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
231a2f3521dSMark F. Adams   else if (coords) {
232c666adbfSMark F. Adams     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf);
23369344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
23469344418SMark F. Adams     if (ndm != ndf) {
235c666adbfSMark F. Adams       if (pc_gamg->data_cell_cols != ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Don't know how to create null space for ndm=%d, ndf=%d.  Use MatSetNearNullSpace.",ndm,ndf);
236a2f3521dSMark F. Adams     }
237806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
23871959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
23971959b99SBarry Smith   if (pc_gamg->data_cell_cols <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data_cell_cols %D <= 0",pc_gamg->data_cell_cols);
240c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2412e68589bSMark F. Adams 
2422e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
2432e68589bSMark F. Adams   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
2442e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
245854ce69bSBarry Smith     ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr);
2462e68589bSMark F. Adams   }
2472e68589bSMark F. Adams   /* copy data in - column oriented */
2482e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
24969344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
25069344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
251c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2522e68589bSMark F. Adams     else {
25369344418SMark F. Adams       /* translational modes */
2542fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
2552fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
25669344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2572e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2582fa5cd67SKarl Rupp         }
2592fa5cd67SKarl Rupp       }
26069344418SMark F. Adams 
26169344418SMark F. Adams       /* rotational modes */
2622e68589bSMark F. Adams       if (coords) {
26369344418SMark F. Adams         if (ndm == 2) {
2642e68589bSMark F. Adams           data   += 2*M;
2652e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2662e68589bSMark F. Adams           data[1] =  coords[2*kk];
267806fa848SBarry Smith         } else {
2682e68589bSMark F. Adams           data   += 3*M;
2692e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2702e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
2712e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2722e68589bSMark F. Adams         }
2732e68589bSMark F. Adams       }
2742e68589bSMark F. Adams     }
2752e68589bSMark F. Adams   }
2762e68589bSMark F. Adams 
2772e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2782e68589bSMark F. Adams   PetscFunctionReturn(0);
2792e68589bSMark F. Adams }
2802e68589bSMark F. Adams 
281b43b03e9SMark F. Adams typedef PetscInt NState;
282b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
283b43b03e9SMark F. Adams static const NState DELETED =-1;
284b43b03e9SMark F. Adams static const NState REMOVED =-3;
285b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
286b43b03e9SMark F. Adams 
287c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
288c8b0795cSMark F. Adams /*
289b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
290b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
291c8b0795cSMark F. Adams 
292c8b0795cSMark F. Adams    Input Parameter:
293c8b0795cSMark F. Adams    . Gmat_2 - glabal matrix of graph (data not defined)
294c8b0795cSMark F. Adams    . Gmat_1 - base graph to grab with
295c8b0795cSMark F. Adams    Input/Output Parameter:
2960cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids)
297c8b0795cSMark F. Adams */
298c8b0795cSMark F. Adams #undef __FUNCT__
299c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs"
3000cbbd2e1SMark F. Adams static PetscErrorCode smoothAggs(const Mat Gmat_2, /* base (squared) graph */
3010cbbd2e1SMark F. Adams                                  const Mat Gmat_1,  /* base graph */
3020cbbd2e1SMark F. Adams                                  /* const IS selected_2, [nselected local] selected vertices */
3031147fc2aSKarl Rupp                                  PetscCoarsenData *aggs_2)  /* [nselected local] global ID of aggregate */
304c8b0795cSMark F. Adams {
305c8b0795cSMark F. Adams   PetscErrorCode ierr;
306c8b0795cSMark F. Adams   PetscBool      isMPI;
307c8b0795cSMark F. Adams   Mat_SeqAIJ     *matA_1, *matB_1=0, *matA_2, *matB_2=0;
3083b4367a7SBarry Smith   MPI_Comm       comm;
3090cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
310c8b0795cSMark F. Adams   Mat_MPIAIJ     *mpimat_2 = 0, *mpimat_1=0;
311c8b0795cSMark F. Adams   const PetscInt nloc      = Gmat_2->rmap->n;
3120cbbd2e1SMark F. Adams   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
3130cbbd2e1SMark F. Adams   PetscInt       *lid_cprowID_1;
314c8b0795cSMark F. Adams   NState         *lid_state;
3150cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
316c8b0795cSMark F. Adams 
317c8b0795cSMark F. Adams   PetscFunctionBegin;
3183b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr);
319c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
320c8b0795cSMark F. Adams 
3210cbbd2e1SMark F. Adams   if (PETSC_FALSE) {
322c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0;
323c8b0795cSMark F. Adams     sprintf(fname,"Gmat2_%d.m",llev++);
3243b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
325c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
326c8b0795cSMark F. Adams     ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr);
327c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
328c8b0795cSMark F. Adams   }
329c8b0795cSMark F. Adams 
330c8b0795cSMark F. Adams   /* get submatrices */
331251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
332c8b0795cSMark F. Adams   if (isMPI) {
333c8b0795cSMark F. Adams     /* grab matrix objects */
334c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
335c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
336c8b0795cSMark F. Adams     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
337c8b0795cSMark F. Adams     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
338c8b0795cSMark F. Adams     matA_2   = (Mat_SeqAIJ*)mpimat_2->A->data;
339c8b0795cSMark F. Adams     matB_2   = (Mat_SeqAIJ*)mpimat_2->B->data;
340c8b0795cSMark F. Adams 
341c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
34211e456e1SBarry Smith     ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
343c8b0795cSMark F. Adams 
344785e854fSJed Brown     ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr);
3450cbbd2e1SMark F. Adams     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
346c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
347c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
348c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
349c8b0795cSMark F. Adams     }
350806fa848SBarry Smith   } else {
351c8b0795cSMark F. Adams     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
352c8b0795cSMark F. Adams     matA_2        = (Mat_SeqAIJ*)Gmat_2->data;
3530298fd71SBarry Smith     lid_cprowID_1 = NULL;
354c8b0795cSMark F. Adams   }
35578a438d6SMark Adams   if (nloc>0) {
35671959b99SBarry Smith     if (!(matA_1 && !matA_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_1 && !matA_1->compressedrow.use)");
35771959b99SBarry Smith     if (!(matB_1==0 || matB_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_1==0 || matB_1->compressedrow.use)");
35871959b99SBarry Smith     if (!(matA_2 && !matA_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_2 && !matA_2->compressedrow.use)");
35971959b99SBarry Smith     if (!(matB_2==0 || matB_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_2==0 || matB_2->compressedrow.use)");
36078a438d6SMark Adams   }
361c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
362785e854fSJed Brown   ierr = PetscMalloc1(nloc, &lid_state);CHKERRQ(ierr);
363785e854fSJed Brown   ierr = PetscMalloc1(nloc, &lid_parent_gid);CHKERRQ(ierr);
364c8b0795cSMark F. Adams   for (lid = 0; lid < nloc; lid++) {
3650cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
366c8b0795cSMark F. Adams     lid_state[lid]      = DELETED;
367c8b0795cSMark F. Adams   }
3680cbbd2e1SMark F. Adams 
3690cbbd2e1SMark F. Adams   /* set lid_state */
3700cbbd2e1SMark F. Adams   for (lid = 0; lid < nloc; lid++) {
37141b27cdeSMark F. Adams     PetscCDPos pos;
372e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
373e78576d6SMark F. Adams     if (pos) {
374e78576d6SMark F. Adams       PetscInt gid1;
3752fa5cd67SKarl Rupp 
37671959b99SBarry Smith       ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
37771959b99SBarry Smith       if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0);
3780cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
379b43b03e9SMark F. Adams     }
380b43b03e9SMark F. Adams   }
3810cbbd2e1SMark F. Adams 
3820cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
383c8b0795cSMark F. Adams   for (lid=kk=0; lid<nloc; lid++) {
384c8b0795cSMark F. Adams     NState state = lid_state[lid];
385c8b0795cSMark F. Adams     if (IS_SELECTED(state)) {
38641b27cdeSMark F. Adams       PetscCDPos pos;
387e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
388e78576d6SMark F. Adams       while (pos) {
389e78576d6SMark F. Adams         PetscInt gid1;
390ffc955d6SMark F. Adams         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
391e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
392e78576d6SMark F. Adams 
3932fa5cd67SKarl Rupp         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
394c8b0795cSMark F. Adams       }
3950cbbd2e1SMark F. Adams     }
3960cbbd2e1SMark F. Adams   }
3970cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
398c8b0795cSMark F. Adams   if (isMPI) {
399c8b0795cSMark F. Adams     Vec tempVec;
400c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
4012a7a6963SBarry Smith     ierr = MatCreateVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr);
402c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
403c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
404c8b0795cSMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
405c8b0795cSMark F. Adams     }
406c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
407c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
408806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
409806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
410c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
411c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
412806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
413806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
414c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
4150cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
4160cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4170cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
4180cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
4190cbbd2e1SMark F. Adams     }
4200cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
4210cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
4220cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
423806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
424806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4250cbbd2e1SMark F. Adams     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
4260cbbd2e1SMark F. Adams 
427c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
428c8b0795cSMark F. Adams   } /* ismpi */
429c8b0795cSMark F. Adams 
430c8b0795cSMark F. Adams   /* doit */
431c8b0795cSMark F. Adams   for (lid=0; lid<nloc; lid++) {
432c8b0795cSMark F. Adams     NState state = lid_state[lid];
4330cbbd2e1SMark F. Adams     if (IS_SELECTED(state)) {
4340cbbd2e1SMark F. Adams       /* steal locals */
435c8b0795cSMark F. Adams       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
436c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
437c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
4380cbbd2e1SMark F. Adams         PetscInt lidj   = idx[j], sgid;
439c8b0795cSMark F. Adams         NState   statej = lid_state[lidj];
4400cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
4410cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
4420cbbd2e1SMark F. Adams           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
4430cbbd2e1SMark F. Adams             PetscInt   hav=0,slid=sgid-my0,gidj=lidj+my0;
4440298fd71SBarry Smith             PetscCDPos pos,last=NULL;
445c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
446e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
447e78576d6SMark F. Adams             while (pos) {
448e78576d6SMark F. Adams               PetscInt gid;
449ffc955d6SMark F. Adams               ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
4500cbbd2e1SMark F. Adams               if (gid == gidj) {
45171959b99SBarry Smith                 if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
45241b27cdeSMark F. Adams                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
45341b27cdeSMark F. Adams                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
4540cbbd2e1SMark F. Adams                 hav  = 1;
455c8b0795cSMark F. Adams                 break;
456806fa848SBarry Smith               } else last = pos;
457e78576d6SMark F. Adams 
458e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
459c8b0795cSMark F. Adams             }
460c8b0795cSMark F. Adams             if (hav!=1) {
46171959b99SBarry Smith               if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
462c666adbfSMark F. Adams               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
463c8b0795cSMark F. Adams             }
464806fa848SBarry Smith           } else {            /* I'm stealing this local, owned by a ghost */
4652fa5cd67SKarl Rupp             if (sgid != -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Have un-symmetric graph (apparently). Use '-pc_gamg_sym_graph true' to symetrize the graph or '-pc_gamg_threshold 0.0' if the matrix is structurally symmetric.");
46641b27cdeSMark F. Adams             ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr);
467c8b0795cSMark F. Adams           }
468c8b0795cSMark F. Adams         }
4690cbbd2e1SMark F. Adams       } /* local neighbors */
470806fa848SBarry Smith     } else if (state == DELETED && lid_cprowID_1) {
4710cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
472c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
473c8b0795cSMark F. Adams       if ((ix=lid_cprowID_1[lid]) != -1) {
474c8b0795cSMark F. Adams         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
475c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
476c8b0795cSMark F. Adams         for (j=0; j<n; j++) {
477c8b0795cSMark F. Adams           PetscInt cpid   = idx[j];
478c8b0795cSMark F. Adams           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
479c8b0795cSMark F. Adams           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
4800cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
4810cbbd2e1SMark F. Adams             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
4820cbbd2e1SMark F. Adams               PetscInt   hav=0,oldslidj=sgidold-my0;
4830298fd71SBarry Smith               PetscCDPos pos,last=NULL;
4840cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
485e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
486e78576d6SMark F. Adams               while (pos) {
487e78576d6SMark F. Adams                 PetscInt gid;
488ffc955d6SMark F. Adams                 ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
4890cbbd2e1SMark F. Adams                 if (lid+my0 == gid) {
4900cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
49171959b99SBarry Smith                   if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
49241b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr);
4930cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
4940cbbd2e1SMark F. Adams                   hav = 1;
495c8b0795cSMark F. Adams                   break;
496806fa848SBarry Smith                 } else last = pos;
497e78576d6SMark F. Adams 
498e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
499c8b0795cSMark F. Adams               }
500c8b0795cSMark F. Adams               if (hav!=1) {
501c666adbfSMark F. Adams                 if (hav==0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
502c666adbfSMark F. Adams                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
503c8b0795cSMark F. Adams               }
504806fa848SBarry Smith             } else {
5050cbbd2e1SMark F. Adams               /* ghosts remove this later */
5060cbbd2e1SMark F. Adams             }
507c8b0795cSMark F. Adams           }
508c8b0795cSMark F. Adams         }
509c8b0795cSMark F. Adams       }
510c8b0795cSMark F. Adams     } /* selected/deleted */
511c8b0795cSMark F. Adams   } /* node loop */
512c8b0795cSMark F. Adams 
513c8b0795cSMark F. Adams   if (isMPI) {
5140cbbd2e1SMark F. Adams     PetscScalar   *cpcol_2_parent,*cpcol_2_gid;
5150cbbd2e1SMark F. Adams     Vec           tempVec,ghostgids2,ghostparents2;
5160cbbd2e1SMark F. Adams     PetscInt      cpid,nghost_2;
5170cbbd2e1SMark F. Adams     GAMGHashTable gid_cpid;
518c8b0795cSMark F. Adams 
5190cbbd2e1SMark F. Adams     ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr);
5202a7a6963SBarry Smith     ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr);
5210cbbd2e1SMark F. Adams 
5220cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
523c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5240cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr);
525c8b0795cSMark F. Adams     }
526c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
527c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5280cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr);
529806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
530806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5310cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
5320cbbd2e1SMark F. Adams 
5330cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
5340cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5350cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
5360cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
5370cbbd2e1SMark F. Adams     }
5380cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
5390cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5400cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr);
541806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
542806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5430cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
5440cbbd2e1SMark F. Adams 
545c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
546c8b0795cSMark F. Adams 
5470cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
548f10ff945SMark F. Adams     ierr = GAMGTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr);
5490cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5500cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
5510cbbd2e1SMark F. Adams       if (state==DELETED) {
5520cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5530cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
5540cbbd2e1SMark F. Adams         if (sgid_old == -1 && sgid_new != -1) {
5550cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5560cbbd2e1SMark F. Adams           ierr = GAMGTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr);
5570cbbd2e1SMark F. Adams         }
5580cbbd2e1SMark F. Adams       }
5590cbbd2e1SMark F. Adams     }
560c8b0795cSMark F. Adams 
5610cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
562c8b0795cSMark F. Adams     for (lid=0; lid<nloc; lid++) {
563c8b0795cSMark F. Adams       NState state = lid_state[lid];
564c8b0795cSMark F. Adams       if (IS_SELECTED(state)) {
5650298fd71SBarry Smith         PetscCDPos pos,last=NULL;
566c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
567e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
568e78576d6SMark F. Adams         while (pos) {
569e78576d6SMark F. Adams           PetscInt gid;
570ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
571e78576d6SMark F. Adams 
5720cbbd2e1SMark F. Adams           if (gid < my0 || gid >= Iend) {
5730cbbd2e1SMark F. Adams             ierr = GAMGTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr);
5740cbbd2e1SMark F. Adams             if (cpid != -1) {
5750cbbd2e1SMark F. Adams               /* a moved ghost - */
5760cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
57741b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr);
578806fa848SBarry Smith             } else last = pos;
579806fa848SBarry Smith           } else last = pos;
580e78576d6SMark F. Adams 
581e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
582c8b0795cSMark F. Adams         } /* loop over list of deleted */
583c8b0795cSMark F. Adams       } /* selected */
584c8b0795cSMark F. Adams     }
5850cbbd2e1SMark F. Adams     ierr = GAMGTableDestroy(&gid_cpid);CHKERRQ(ierr);
586c8b0795cSMark F. Adams 
5870cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
5880cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5890cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5900cbbd2e1SMark F. Adams       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
5910cbbd2e1SMark F. Adams         PetscInt   gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5920cbbd2e1SMark F. Adams         PetscInt   slid_new=sgid_new-my0,hav=0;
59341b27cdeSMark F. Adams         PetscCDPos pos;
5940cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
595e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
596e78576d6SMark F. Adams         while (pos) {
597e78576d6SMark F. Adams           PetscInt gidj;
598ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &gidj);CHKERRQ(ierr);
599e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
600e78576d6SMark F. Adams 
6010cbbd2e1SMark F. Adams           if (gidj == gid) { hav = 1; break; }
602c8b0795cSMark F. Adams         }
603c8b0795cSMark F. Adams         if (hav != 1) {
604ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
60541b27cdeSMark F. Adams           ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr);
606c8b0795cSMark F. Adams         }
607c8b0795cSMark F. Adams       }
608c8b0795cSMark F. Adams     }
609c8b0795cSMark F. Adams 
6100cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
6110cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
6120cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
6130cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
614c8b0795cSMark F. Adams     ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr);
6150cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr);
6160cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr);
6170cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr);
618c8b0795cSMark F. Adams   }
619c8b0795cSMark F. Adams 
6200cbbd2e1SMark F. Adams   ierr = PetscFree(lid_parent_gid);CHKERRQ(ierr);
621c8b0795cSMark F. Adams   ierr = PetscFree(lid_state);CHKERRQ(ierr);
622c8b0795cSMark F. Adams   PetscFunctionReturn(0);
623c8b0795cSMark F. Adams }
6242e68589bSMark F. Adams 
6252e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6262e68589bSMark F. Adams /*
627a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
628a2f3521dSMark F. Adams       Looks in Mat for near null space.
629a2f3521dSMark F. Adams       Does not work for Stokes
6302e68589bSMark F. Adams 
6312e68589bSMark F. Adams   Input Parameter:
6322e68589bSMark F. Adams    . pc -
633a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
6342e68589bSMark F. Adams */
6352e68589bSMark F. Adams #undef __FUNCT__
6362e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG"
637b8cd405aSMark F. Adams PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
6382e68589bSMark F. Adams {
6392e68589bSMark F. Adams   PetscErrorCode ierr;
640b8cd405aSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
641b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
642b8cd405aSMark F. Adams   MatNullSpace   mnull;
643b8cd405aSMark F. Adams 
6442e68589bSMark F. Adams   PetscFunctionBegin;
645b8cd405aSMark F. Adams   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
646b8cd405aSMark F. Adams   if (!mnull) {
647a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
6489d1b15c3SMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
64971959b99SBarry Smith     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);
65071959b99SBarry Smith     if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs);
6510298fd71SBarry Smith     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr);
652806fa848SBarry Smith   } else {
653b8cd405aSMark F. Adams     PetscReal *nullvec;
654b8cd405aSMark F. Adams     PetscBool has_const;
655b8cd405aSMark F. Adams     PetscInt  i,j,mlocal,nvec,bs;
656b8cd405aSMark F. Adams     const Vec *vecs; const PetscScalar *v;
6570298fd71SBarry Smith     ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr);
658b8cd405aSMark F. Adams     ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr);
659a0dea87dSMark Adams     pc_gamg->data_sz = (nvec+!!has_const)*mlocal;
660785e854fSJed Brown     ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr);
661b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
662b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
663b8cd405aSMark F. Adams       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
664b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
665b8cd405aSMark F. Adams       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
666b8cd405aSMark F. Adams     }
667b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
668b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
6692fa5cd67SKarl Rupp 
67006e133e7SMark Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
6712fa5cd67SKarl Rupp 
672b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
673b8cd405aSMark F. Adams   }
6742e68589bSMark F. Adams   PetscFunctionReturn(0);
6752e68589bSMark F. Adams }
6762e68589bSMark F. Adams 
6772e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6782e68589bSMark F. Adams /*
6792e68589bSMark F. Adams  formProl0
6802e68589bSMark F. Adams 
6812e68589bSMark F. Adams    Input Parameter:
6820cbbd2e1SMark F. Adams    . agg_llists - list of arrays with aggregates
6832e68589bSMark F. Adams    . bs - block size
6840cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
6850cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
6862e68589bSMark F. Adams    . data_stride - bs*(nloc nodes + ghost nodes)
6872e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
6882e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
6892e68589bSMark F. Adams   Output Parameter:
6902e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
6912e68589bSMark F. Adams    . a_Prol - prolongation operator
6922e68589bSMark F. Adams */
6932e68589bSMark F. Adams #undef __FUNCT__
6942e68589bSMark F. Adams #define __FUNCT__ "formProl0"
6950cbbd2e1SMark F. Adams static PetscErrorCode formProl0(const PetscCoarsenData *agg_llists, /* list from selected vertices of aggregate unselected vertices */
6960cbbd2e1SMark F. Adams                                 const PetscInt bs,          /* (row) block size */
6970cbbd2e1SMark F. Adams                                 const PetscInt nSAvec,      /* column bs */
6980cbbd2e1SMark F. Adams                                 const PetscInt my0crs,      /* global index of start of locals */
6990cbbd2e1SMark F. Adams                                 const PetscInt data_stride, /* (nloc+nghost)*bs */
7000cbbd2e1SMark F. Adams                                 PetscReal      data_in[],   /* [data_stride][nSAvec] */
7010cbbd2e1SMark F. Adams                                 const PetscInt flid_fgid[], /* [data_stride/bs] */
7022e68589bSMark F. Adams                                 PetscReal **a_data_out,
7031147fc2aSKarl Rupp                                 Mat a_Prol) /* prolongation operator (output)*/
7042e68589bSMark F. Adams {
7052e68589bSMark F. Adams   PetscErrorCode ierr;
7060cbbd2e1SMark F. Adams   PetscInt       Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
7073b4367a7SBarry Smith   MPI_Comm       comm;
70873911c69SBarry Smith   PetscMPIInt    rank;
7092e68589bSMark F. Adams   PetscReal      *out_data;
71041b27cdeSMark F. Adams   PetscCDPos     pos;
7110cbbd2e1SMark F. Adams   GAMGHashTable  fgid_flid;
7120cbbd2e1SMark F. Adams 
713797e13b7SMark F. Adams /* #define OUT_AGGS */
714519f805aSKarl Rupp #if defined(OUT_AGGS)
7150298fd71SBarry Smith   static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM;
7169057884aSMark F. Adams #endif
7172e68589bSMark F. Adams 
7182e68589bSMark F. Adams   PetscFunctionBegin;
7193b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
7203b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
7212e68589bSMark F. Adams   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
72271959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
72371959b99SBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Iend %D - Istart %d must be divisible by bs %D",Iend,Istart,bs);
7240cbbd2e1SMark F. Adams   Iend   /= bs;
7250cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
7262e68589bSMark F. Adams 
727f10ff945SMark F. Adams   ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
7280cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
7290cbbd2e1SMark F. Adams     ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
7302e68589bSMark F. Adams   }
7312e68589bSMark F. Adams 
732519f805aSKarl Rupp #if defined(OUT_AGGS)
733c5df96a5SBarry Smith   sprintf(fname,"aggs_%d_%d.m",llev++,rank);
7342fa5cd67SKarl Rupp   if (llev==1) file = fopen(fname,"w");
7350cbbd2e1SMark F. Adams   MatGetSize(a_Prol, &pM, &jj);
7360cbbd2e1SMark F. Adams #endif
7370cbbd2e1SMark F. Adams 
7380cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
7390cbbd2e1SMark F. Adams   for (nSelected=mm=0; mm<nloc; mm++) {
740e78576d6SMark F. Adams     PetscBool ise;
741e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
742e78576d6SMark F. Adams     if (!ise) nSelected++;
7430cbbd2e1SMark F. Adams   }
7440cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
74571959b99SBarry Smith   if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D  != my0crs %D",ii,nSAvec,my0crs);
74671959b99SBarry Smith   if (nSelected != (jj-ii)/nSAvec) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nSelected %D != (jj %D - ii %D)/nSAvec %D",nSelected,jj,ii,nSAvec);
7470cbbd2e1SMark F. Adams 
7482e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
7490cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
7502fa5cd67SKarl Rupp 
751785e854fSJed Brown   ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr);
75233812677SSatish Balay   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL;
7530cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
7542e68589bSMark F. Adams 
7552e68589bSMark F. Adams   /* find points and set prolongation */
756c8b0795cSMark F. Adams   minsz = 100;
7572e68589bSMark F. Adams   ndone = 0;
7580cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
759e78576d6SMark F. Adams     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
760e78576d6SMark F. Adams     if (jj > 0) {
7610cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
7620cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
7630cbbd2e1SMark F. Adams       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
7642e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
7652e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
7662e68589bSMark F. Adams       PetscInt       *fids;
76765d7b583SSatish Balay       PetscReal      *data;
7680cbbd2e1SMark F. Adams       /* count agg */
7690cbbd2e1SMark F. Adams       if (asz<minsz) minsz = asz;
7700cbbd2e1SMark F. Adams 
7710cbbd2e1SMark F. Adams       /* get block */
772854ce69bSBarry Smith       ierr = PetscMalloc1(Mdata*N, &qqc);CHKERRQ(ierr);
773854ce69bSBarry Smith       ierr = PetscMalloc1(M*N, &qqr);CHKERRQ(ierr);
774785e854fSJed Brown       ierr = PetscMalloc1(N, &TAU);CHKERRQ(ierr);
775785e854fSJed Brown       ierr = PetscMalloc1(LWORK, &WORK);CHKERRQ(ierr);
776785e854fSJed Brown       ierr = PetscMalloc1(M, &fids);CHKERRQ(ierr);
7772e68589bSMark F. Adams 
7782e68589bSMark F. Adams       aggID = 0;
779e78576d6SMark F. Adams       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
780e78576d6SMark F. Adams       while (pos) {
781e78576d6SMark F. Adams         PetscInt gid1;
782ffc955d6SMark F. Adams         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
783e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
784e78576d6SMark F. Adams 
7850cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
7860cbbd2e1SMark F. Adams         else {
7870cbbd2e1SMark F. Adams           ierr = GAMGTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
78871959b99SBarry Smith           if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
7890cbbd2e1SMark F. Adams         }
7902e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
79165d7b583SSatish Balay         data = &data_in[flid*bs];
7922e68589bSMark F. Adams         for (kk = ii = 0; ii < bs; ii++) {
7932e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
7940cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
7950cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
7962e68589bSMark F. Adams           }
7972e68589bSMark F. Adams         }
798519f805aSKarl Rupp #if defined(OUT_AGGS)
799b2a4f308SMark F. Adams         if (llev==1) {
8009057884aSMark F. Adams           char     str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
8010cbbd2e1SMark F. Adams           PetscInt MM,pi,pj;
802c5df96a5SBarry Smith           str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11];
803f7620de1SMatthew G Knepley           MM      = (PetscInt)(PetscSqrtReal((PetscReal)pM));
8040cbbd2e1SMark F. Adams           pj      = gid1/MM; pi = gid1%MM;
805b2a4f308SMark F. Adams           fprintf(file,str,(double)pi,(double)pj);
806b2a4f308SMark F. Adams           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
8079057884aSMark F. Adams         }
8089057884aSMark F. Adams #endif
8092e68589bSMark F. Adams         /* set fine IDs */
8102e68589bSMark F. Adams         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
8112e68589bSMark F. Adams 
8122e68589bSMark F. Adams         aggID++;
8130cbbd2e1SMark F. Adams       }
8142e68589bSMark F. Adams 
8152e68589bSMark F. Adams       /* pad with zeros */
8162e68589bSMark F. Adams       for (ii = asz*bs; ii < Mdata; ii++) {
8172e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
8182e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
8192e68589bSMark F. Adams         }
8202e68589bSMark F. Adams       }
8212e68589bSMark F. Adams 
8222e68589bSMark F. Adams       ndone += aggID;
8232e68589bSMark F. Adams       /* QR */
82484df3f90SSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
8258b83055fSJed Brown       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
82684df3f90SSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
827d23427c4SJed Brown       if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
8282e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
8292e68589bSMark F. Adams       {
8302e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
8312e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
8322e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
83333812677SSatish Balay             if (data[jj*out_data_stride + ii] != PETSC_MAX_REAL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"data[jj*out_data_stride + ii] != %e",PETSC_MAX_REAL);
8340cbbd2e1SMark F. Adams            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
8350cbbd2e1SMark F. Adams            else data[jj*out_data_stride + ii] = 0.;
8362e68589bSMark F. Adams           }
8372e68589bSMark F. Adams         }
8382e68589bSMark F. Adams       }
8392e68589bSMark F. Adams 
8402e68589bSMark F. Adams       /* get Q - row oriented */
8418b83055fSJed Brown       PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
842c666adbfSMark F. Adams       if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
8432e68589bSMark F. Adams 
8442e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
8452e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
8462e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
8472e68589bSMark F. Adams         }
8482e68589bSMark F. Adams       }
8492e68589bSMark F. Adams 
8502e68589bSMark F. Adams       /* add diagonal block of P0 */
851c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
852c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
853c8b0795cSMark F. Adams       }
8542e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
8552e68589bSMark F. Adams 
8562e68589bSMark F. Adams       ierr = PetscFree(qqc);CHKERRQ(ierr);
8572e68589bSMark F. Adams       ierr = PetscFree(qqr);CHKERRQ(ierr);
8582e68589bSMark F. Adams       ierr = PetscFree(TAU);CHKERRQ(ierr);
8592e68589bSMark F. Adams       ierr = PetscFree(WORK);CHKERRQ(ierr);
8602e68589bSMark F. Adams       ierr = PetscFree(fids);CHKERRQ(ierr);
861b43b03e9SMark F. Adams       clid++;
8620cbbd2e1SMark F. Adams     } /* coarse agg */
8630cbbd2e1SMark F. Adams   } /* for all fine nodes */
8640cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8650cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8662e68589bSMark F. Adams 
8673b4367a7SBarry Smith /* ierr = MPI_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */
8682e68589bSMark F. Adams /* MatGetSize(a_Prol, &kk, &jj); */
8693b4367a7SBarry Smith /* ierr = MPI_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */
8703b4367a7SBarry Smith /* PetscPrintf(comm," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",rank,__FUNCT__,ii,kk/bs,ndone,jj); */
8712e68589bSMark F. Adams 
872519f805aSKarl Rupp #if defined(OUT_AGGS)
873b2a4f308SMark F. Adams   if (llev==1) fclose(file);
8749057884aSMark F. Adams #endif
8750cbbd2e1SMark F. Adams   ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr);
8762e68589bSMark F. Adams   PetscFunctionReturn(0);
8772e68589bSMark F. Adams }
8782e68589bSMark F. Adams 
8795adeb434SBarry Smith #undef __FUNCT__
8805adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG_AGG"
8815adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer)
8825adeb434SBarry Smith {
8835adeb434SBarry Smith   PetscErrorCode ierr;
8845adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
8855adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
8865adeb434SBarry Smith   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
8875adeb434SBarry Smith 
8885adeb434SBarry Smith   PetscFunctionBegin;
8895adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"      AGG specific options\n");CHKERRQ(ierr);
8905adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr);
8915adeb434SBarry Smith   PetscFunctionReturn(0);
8925adeb434SBarry Smith }
8935adeb434SBarry Smith 
8942e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
8952e68589bSMark F. Adams /*
896fd1112cbSBarry Smith    PCGAMGGraph_AGG
8972e68589bSMark F. Adams 
8982e68589bSMark F. Adams   Input Parameter:
8992e68589bSMark F. Adams    . pc - this
9002e68589bSMark F. Adams    . Amat - matrix on this fine level
9012e68589bSMark F. Adams   Output Parameter:
902c8b0795cSMark F. Adams    . a_Gmat -
9032e68589bSMark F. Adams */
9042e68589bSMark F. Adams #undef __FUNCT__
905fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGGraph_AGG"
906fd1112cbSBarry Smith PetscErrorCode PCGAMGGraph_AGG(PC pc,const Mat Amat,Mat *a_Gmat)
907c8b0795cSMark F. Adams {
908c8b0795cSMark F. Adams   PetscErrorCode            ierr;
909c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
910c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
911c8b0795cSMark F. Adams   const PetscReal           vfilter      = pc_gamg->threshold;
912c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
913e0940f08SMark F. Adams   Mat                       Gmat;
9143b4367a7SBarry Smith   MPI_Comm                  comm;
91567744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
916c8b0795cSMark F. Adams 
917c8b0795cSMark F. Adams   PetscFunctionBegin;
9183b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
919fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
920c8b0795cSMark F. Adams 
92167744fedSMark F. Adams   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
922c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
9230cbbd2e1SMark F. Adams 
9242d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
925*bf4339c2SBarry Smith   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr);
926c8b0795cSMark F. Adams 
927e0940f08SMark F. Adams   *a_Gmat = Gmat;
928c8b0795cSMark F. Adams 
929fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
930c8b0795cSMark F. Adams   PetscFunctionReturn(0);
931c8b0795cSMark F. Adams }
932c8b0795cSMark F. Adams 
933c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
934c8b0795cSMark F. Adams /*
935b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
936c8b0795cSMark F. Adams 
937c8b0795cSMark F. Adams   Input Parameter:
938e0940f08SMark F. Adams    . a_pc - this
939e0940f08SMark F. Adams   Input/Output Parameter:
9400cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
941c8b0795cSMark F. Adams   Output Parameter:
9420cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
943c8b0795cSMark F. Adams */
944c8b0795cSMark F. Adams #undef __FUNCT__
945b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG"
9462fa5cd67SKarl Rupp PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
947c8b0795cSMark F. Adams {
948c8b0795cSMark F. Adams   PetscErrorCode ierr;
949e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
950c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
951c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
9520cbbd2e1SMark F. Adams   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
9530cbbd2e1SMark F. Adams   IS             perm;
954c8b0795cSMark F. Adams   PetscInt       Ii,nloc,bs,n,m;
955c8b0795cSMark F. Adams   PetscInt       *permute;
956c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
957b43b03e9SMark F. Adams   MatCoarsen     crs;
9583b4367a7SBarry Smith   MPI_Comm       comm;
95973911c69SBarry Smith   PetscMPIInt    rank;
960c8b0795cSMark F. Adams 
961c8b0795cSMark F. Adams   PetscFunctionBegin;
9620cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
9633b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
9643b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
965e0940f08SMark F. Adams   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
96671959b99SBarry Smith   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr);
96771959b99SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs);
968c8b0795cSMark F. Adams   nloc = n/bs;
969c8b0795cSMark F. Adams 
97057d29afaSToby Isaac   if (pc_gamg->firstCoarsen && pc_gamg_agg->square_graph) { /* we don't want to square coarse grids */
971878e152fSMark F. Adams     /* ierr = MatMatTransposeMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2); */
972806fa848SBarry Smith     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
973806fa848SBarry Smith   } else Gmat2 = Gmat1;
974c8b0795cSMark F. Adams 
975c8b0795cSMark F. Adams   /* get MIS aggs */
976c8b0795cSMark F. Adams   /* randomize */
977785e854fSJed Brown   ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr);
978785e854fSJed Brown   ierr = PetscMalloc1(nloc, &bIndexSet);CHKERRQ(ierr);
979c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
980c8b0795cSMark F. Adams     bIndexSet[Ii] = PETSC_FALSE;
981c8b0795cSMark F. Adams     permute[Ii]   = Ii;
982c8b0795cSMark F. Adams   }
983c8b0795cSMark F. Adams   srand(1); /* make deterministic */
984c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
985c8b0795cSMark F. Adams     PetscInt iSwapIndex = rand()%nloc;
986c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
987c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
988c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
989c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
990c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
991c8b0795cSMark F. Adams     }
992c8b0795cSMark F. Adams   }
993c8b0795cSMark F. Adams   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
994c8b0795cSMark F. Adams 
995806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
9960cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
9970cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
998b43b03e9SMark F. Adams #endif
9993b4367a7SBarry Smith   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
10009057884aSMark F. Adams   /* ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); */
10019057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1002b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
1003b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
1004b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1005b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
10060cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
1007b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1008b43b03e9SMark F. Adams 
1009c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
1010c8b0795cSMark F. Adams   ierr = PetscFree(permute);CHKERRQ(ierr);
10110cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10120cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1013b43b03e9SMark F. Adams #endif
10149f3f12c8SMark F. Adams 
1015c8b0795cSMark F. Adams   /* smooth aggs */
1016e0940f08SMark F. Adams   if (Gmat2 != Gmat1) {
10170cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10180cbbd2e1SMark F. Adams     ierr     = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
1019c8b0795cSMark F. Adams     ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
1020e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
102141b27cdeSMark F. Adams     ierr     = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10223b4367a7SBarry Smith     if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1023806fa848SBarry Smith   } else {
10240cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
1025f10ff945SMark F. Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenAppply) */
102641b27cdeSMark F. Adams     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10270cbbd2e1SMark F. Adams     if (mat) {
10280cbbd2e1SMark F. Adams       ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
10290cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
10300cbbd2e1SMark F. Adams     }
10310cbbd2e1SMark F. Adams   }
10320cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
1033c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1034c8b0795cSMark F. Adams }
1035c8b0795cSMark F. Adams 
1036c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1037c8b0795cSMark F. Adams /*
10380cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
1039c8b0795cSMark F. Adams 
1040c8b0795cSMark F. Adams  Input Parameter:
1041c8b0795cSMark F. Adams  . pc - this
1042c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1043c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
10440cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1045c8b0795cSMark F. Adams  Output Parameter:
1046c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1047c8b0795cSMark F. Adams  */
1048c8b0795cSMark F. Adams #undef __FUNCT__
10490cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG"
10502fa5cd67SKarl Rupp PetscErrorCode PCGAMGProlongator_AGG(PC pc,const Mat Amat,const Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
10512e68589bSMark F. Adams {
10522e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
10532e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
1054c8b0795cSMark F. Adams   const PetscInt data_cols = pc_gamg->data_cell_cols;
10552e68589bSMark F. Adams   PetscErrorCode ierr;
1056c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1057c8b0795cSMark F. Adams   Mat            Prol;
1058c5df96a5SBarry Smith   PetscMPIInt    rank, size;
10593b4367a7SBarry Smith   MPI_Comm       comm;
10600cbbd2e1SMark F. Adams   const PetscInt col_bs = data_cols;
1061c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
1062c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
10632e68589bSMark F. Adams 
10642e68589bSMark F. Adams   PetscFunctionBegin;
10653b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
10660cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
10673b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10683b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
10692e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1070c8b0795cSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
107171959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
107271959b99SBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs);
10732e68589bSMark F. Adams 
10742e68589bSMark F. Adams   /* get 'nLocalSelected' */
10750cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1076e78576d6SMark F. Adams     PetscBool ise;
10770cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
1078e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1079e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
10802e68589bSMark F. Adams   }
10812e68589bSMark F. Adams 
10822e68589bSMark F. Adams   /* create prolongator, create P matrix */
10833b4367a7SBarry Smith   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
1084806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1085a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
1086a2f3521dSMark F. Adams   ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr);
10870298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(Prol, data_cols, NULL);CHKERRQ(ierr);
10880298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(Prol,data_cols, NULL,data_cols, NULL);CHKERRQ(ierr);
1089a2f3521dSMark F. Adams   /* nloc*bs, nLocalSelected*col_bs, */
1090a2f3521dSMark F. Adams   /* PETSC_DETERMINE, PETSC_DETERMINE, */
10910298fd71SBarry Smith   /* data_cols, NULL, data_cols, NULL, */
1092a2f3521dSMark F. Adams   /* &Prol); */
10932e68589bSMark F. Adams 
10942e68589bSMark F. Adams   /* can get all points "removed" */
1095c8b0795cSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
1096c8b0795cSMark F. Adams   if (ii==0) {
1097*bf4339c2SBarry Smith     ierr = PetscInfo(pc,"no selected points on coarse grid\n");CHKERRQ(ierr);
10982e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
10990298fd71SBarry Smith     *a_P_out = NULL;  /* out */
11002e68589bSMark F. Adams     PetscFunctionReturn(0);
11012e68589bSMark F. Adams   }
1102*bf4339c2SBarry Smith   ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr);
1103c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
11040cbbd2e1SMark F. Adams 
110571959b99SBarry Smith   if ((kk-myCrs0) % col_bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D -myCrs0 %D) not divisible by col_bs %D",kk,myCrs0,col_bs);
1106c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
110771959b99SBarry Smith   if ((kk/col_bs-myCrs0) != nLocalSelected) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D/col_bs %D - myCrs0 %D) != nLocalSelected %D)",kk,col_bs,myCrs0,nLocalSelected);
11082e68589bSMark F. Adams 
11092e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
11100cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11110cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11122e68589bSMark F. Adams #endif
1113c5df96a5SBarry Smith   if (size > 1) { /*  */
11142e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
1115785e854fSJed Brown     ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr);
11162e68589bSMark F. Adams     for (jj = 0; jj < data_cols; jj++) {
1117c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1118a2f3521dSMark F. Adams         PetscInt        ii,stride;
1119c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
11202fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
11212fa5cd67SKarl Rupp 
1122806fa848SBarry Smith         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1123a2f3521dSMark F. Adams 
11242e68589bSMark F. Adams         if (jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
1125785e854fSJed Brown           ierr    = PetscMalloc1(stride*bs*data_cols, &data_w_ghost);CHKERRQ(ierr);
1126a2f3521dSMark F. Adams           nbnodes = bs*stride;
11272e68589bSMark F. Adams         }
1128a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
11292fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
11302e68589bSMark F. Adams         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
11312e68589bSMark F. Adams       }
11322e68589bSMark F. Adams     }
11332e68589bSMark F. Adams     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1134806fa848SBarry Smith   } else {
1135c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
1136c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
11372e68589bSMark F. Adams   }
11382e68589bSMark F. Adams 
11392e68589bSMark F. Adams   /* get P0 */
1140c5df96a5SBarry Smith   if (size > 1) {
11412e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1142a2f3521dSMark F. Adams     PetscInt  stride;
11432e68589bSMark F. Adams 
1144785e854fSJed Brown     ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr);
11452e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1146806fa848SBarry Smith     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1147785e854fSJed Brown     ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr);
1148a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
11492e68589bSMark F. Adams     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1150a2f3521dSMark F. Adams 
115171959b99SBarry Smith     if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs);
11522e68589bSMark F. Adams     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1153806fa848SBarry Smith   } else {
1154785e854fSJed Brown     ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr);
11552e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
11562e68589bSMark F. Adams   }
11570cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11580cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11590cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
11602e68589bSMark F. Adams #endif
1161c8b0795cSMark F. Adams   {
11620298fd71SBarry Smith     PetscReal *data_out = NULL;
1163806fa848SBarry Smith     ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1164c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
11652fa5cd67SKarl Rupp 
1166c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
1167c8b0795cSMark F. Adams     pc_gamg->data_cell_rows = data_cols;
1168c8b0795cSMark F. Adams     pc_gamg->data_sz        = data_cols*data_cols*nLocalSelected;
1169c8b0795cSMark F. Adams   }
11700cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11710cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1172c8b0795cSMark F. Adams #endif
1173c5df96a5SBarry Smith   if (size > 1) ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);
11742e68589bSMark F. Adams   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
11752e68589bSMark F. Adams 
1176c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
1177ed4e82eaSPeter Brune 
11780cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1179c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1180c8b0795cSMark F. Adams }
1181c8b0795cSMark F. Adams 
1182c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1183c8b0795cSMark F. Adams /*
1184fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
1185c8b0795cSMark F. Adams 
1186c8b0795cSMark F. Adams   Input Parameter:
1187c8b0795cSMark F. Adams    . pc - this
1188c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1189c8b0795cSMark F. Adams  In/Output Parameter:
1190c8b0795cSMark F. Adams    . a_P_out - prolongation operator to the next level
1191c8b0795cSMark F. Adams */
1192c8b0795cSMark F. Adams #undef __FUNCT__
1193fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGOptProlongator_AGG"
1194fd1112cbSBarry Smith PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,const Mat Amat,Mat *a_P)
1195c8b0795cSMark F. Adams {
1196c8b0795cSMark F. Adams   PetscErrorCode ierr;
1197c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1198c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1199c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1200c8b0795cSMark F. Adams   PetscInt       jj;
1201c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
12023b4367a7SBarry Smith   MPI_Comm       comm;
1203c8b0795cSMark F. Adams 
1204c8b0795cSMark F. Adams   PetscFunctionBegin;
12053b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
1206fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1207c8b0795cSMark F. Adams 
12082e68589bSMark F. Adams   /* smooth P0 */
1209c8b0795cSMark F. Adams   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
12102e68589bSMark F. Adams     Mat       tMat;
12112e68589bSMark F. Adams     Vec       diag;
12122e68589bSMark F. Adams     PetscReal alpha, emax, emin;
12130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12140cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12152e68589bSMark F. Adams #endif
12162e68589bSMark F. Adams     if (jj == 0) {
12172e68589bSMark F. Adams       KSP eksp;
12182e68589bSMark F. Adams       Vec bb, xx;
1219da509fc8SJed Brown       PC  epc;
12202a7a6963SBarry Smith       ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr);
12212a7a6963SBarry Smith       ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr);
12222e68589bSMark F. Adams       {
12232e68589bSMark F. Adams         PetscRandom rctx;
12243b4367a7SBarry Smith         ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
12252e68589bSMark F. Adams         ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
12262e68589bSMark F. Adams         ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
12272e68589bSMark F. Adams         ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
12282e68589bSMark F. Adams       }
1229e696ed0bSMark F. Adams 
1230e696ed0bSMark F. Adams       /* zeroing out BC rows -- needed for crazy matrices */
1231e696ed0bSMark F. Adams       {
1232e696ed0bSMark F. Adams         PetscInt    Istart,Iend,ncols,jj,Ii;
1233e696ed0bSMark F. Adams         PetscScalar zero = 0.0;
1234e696ed0bSMark F. Adams         ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1235e696ed0bSMark F. Adams         for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
1236e696ed0bSMark F. Adams           ierr = MatGetRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1237e696ed0bSMark F. Adams           if (ncols <= 1) {
1238e696ed0bSMark F. Adams             ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
1239e696ed0bSMark F. Adams           }
1240e696ed0bSMark F. Adams           ierr = MatRestoreRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1241e696ed0bSMark F. Adams         }
1242e696ed0bSMark F. Adams         ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
1243e696ed0bSMark F. Adams         ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
1244e696ed0bSMark F. Adams       }
1245e696ed0bSMark F. Adams 
12463b4367a7SBarry Smith       ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr);
1247806fa848SBarry Smith       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
12482e68589bSMark F. Adams       ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1249da509fc8SJed Brown       ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
1250c2436cacSMark F. Adams       ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
1251c2436cacSMark F. Adams       ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
1252c2436cacSMark F. Adams 
1253c2436cacSMark F. Adams       ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
125423ee1639SBarry Smith       ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr);
12552e68589bSMark F. Adams       ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
12562e68589bSMark F. Adams 
1257da509fc8SJed Brown       ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr);
1258da509fc8SJed Brown       ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1259c2436cacSMark F. Adams 
12602e68589bSMark F. Adams       /* solve - keep stuff out of logging */
12612e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
12622e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
12632e68589bSMark F. Adams       ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
12642e68589bSMark F. Adams       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
12652e68589bSMark F. Adams       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
12662e68589bSMark F. Adams 
12672e68589bSMark F. Adams       ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
1268*bf4339c2SBarry Smith       ierr = PetscInfo3(pc,"\t\t\t smooth P0: max eigen=%e min=%e PC=%s\n",emax,emin,PCJACOBI);CHKERRQ(ierr);
12692e68589bSMark F. Adams       ierr = VecDestroy(&xx);CHKERRQ(ierr);
12702e68589bSMark F. Adams       ierr = VecDestroy(&bb);CHKERRQ(ierr);
12712e68589bSMark F. Adams       ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
12722e68589bSMark F. Adams 
12732e68589bSMark F. Adams       if (pc_gamg->emax_id == -1) {
12743bf036e2SBarry Smith         ierr = PetscObjectComposedDataRegister(&pc_gamg->emax_id);CHKERRQ(ierr);
127571959b99SBarry Smith         if (pc_gamg->emax_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->emax_id == -1");
12762e68589bSMark F. Adams       }
12772e68589bSMark F. Adams       ierr = PetscObjectComposedDataSetScalar((PetscObject)Amat, pc_gamg->emax_id, emax);CHKERRQ(ierr);
12782e68589bSMark F. Adams     }
12792e68589bSMark F. Adams 
12802e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
12812e68589bSMark F. Adams     ierr  = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
12822a7a6963SBarry Smith     ierr  = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr);
12832e68589bSMark F. Adams     ierr  = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
12842e68589bSMark F. Adams     ierr  = VecReciprocal(diag);CHKERRQ(ierr);
12852e68589bSMark F. Adams     ierr  = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
12862e68589bSMark F. Adams     ierr  = VecDestroy(&diag);CHKERRQ(ierr);
1287b4ec6429SMark F. Adams     alpha = -1.4/emax;
12882e68589bSMark F. Adams     ierr  = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
12892e68589bSMark F. Adams     ierr  = MatDestroy(&Prol);CHKERRQ(ierr);
12902e68589bSMark F. Adams     Prol  = tMat;
12910cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12920cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12932e68589bSMark F. Adams #endif
12942e68589bSMark F. Adams   }
1295fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1296c8b0795cSMark F. Adams   *a_P = Prol;
1297c8b0795cSMark F. Adams   PetscFunctionReturn(0);
12982e68589bSMark F. Adams }
12992e68589bSMark F. Adams 
1300c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1301c8b0795cSMark F. Adams /*
1302c8b0795cSMark F. Adams    PCCreateGAMG_AGG
13032e68589bSMark F. Adams 
1304c8b0795cSMark F. Adams   Input Parameter:
1305c8b0795cSMark F. Adams    . pc -
1306c8b0795cSMark F. Adams */
1307c8b0795cSMark F. Adams #undef __FUNCT__
1308c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG"
1309c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1310c8b0795cSMark F. Adams {
1311c8b0795cSMark F. Adams   PetscErrorCode ierr;
1312c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1313c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1314c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
13152e68589bSMark F. Adams 
1316c8b0795cSMark F. Adams   PetscFunctionBegin;
1317c8b0795cSMark F. Adams   /* create sub context for SA */
1318b00a9115SJed Brown   ierr            = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr);
1319c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1320c8b0795cSMark F. Adams 
13211ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
13229b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1323c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1324c8b0795cSMark F. Adams 
1325c8b0795cSMark F. Adams   /* set internal function pointers */
1326fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
13271ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
13281ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
1329fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
13301ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
13315adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
1332c8b0795cSMark F. Adams 
1333bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr);
13342e68589bSMark F. Adams   PetscFunctionReturn(0);
13352e68589bSMark F. Adams }
1336