xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision e213d8f197ea160ee2f05fc0bd1d6231459cf4e1)
1 /*
2  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
3  */
4 #include <../src/ksp/pc/impls/gamg/gamg.h>
5 
6 PetscLogEvent gamg_setup_stages[NUM_SET];
7 
8 /* Private context for the GAMG preconditioner */
9 typedef struct gamg_TAG {
10   PetscInt       m_dim;
11   PetscInt       m_Nlevels;
12   PetscInt       m_data_sz;
13   PetscInt       m_data_rows;
14   PetscInt       m_data_cols;
15   PetscBool      m_useSA;
16   PetscReal     *m_data; /* blocked vector of vertex data on fine grid (coordinates) */
17 } PC_GAMG;
18 
19 /* -------------------------------------------------------------------------- */
20 /*
21    PCSetCoordinates_GAMG
22 
23    Input Parameter:
24    .  pc - the preconditioner context
25 */
26 EXTERN_C_BEGIN
27 #undef __FUNCT__
28 #define __FUNCT__ "PCSetCoordinates_GAMG"
29 PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords )
30 {
31   PC_MG          *mg = (PC_MG*)a_pc->data;
32   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
33   PetscErrorCode ierr;
34   PetscInt       arrsz,bs,my0,kk,ii,jj,nloc,Iend;
35   Mat            Amat = a_pc->pmat;
36   PetscBool      flag;
37   char           str[16];
38 
39   PetscFunctionBegin;
40   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
41   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
42   nloc = (Iend-my0)/bs;
43   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
44 
45   ierr  = PetscOptionsGetString(PETSC_NULL,"-pc_gamg_type",str,16,&flag);    CHKERRQ( ierr );
46   pc_gamg->m_useSA = (PetscBool)(flag && strcmp(str,"sa") == 0);
47 
48   pc_gamg->m_data_rows = 1;
49   if(a_coords == 0) pc_gamg->m_useSA = PETSC_TRUE; /* use SA if no data */
50   if( !pc_gamg->m_useSA ) pc_gamg->m_data_cols = a_ndm; /* coordinates */
51   else{ /* SA: null space vectors */
52     if(a_coords != 0 && bs==1 ) pc_gamg->m_data_cols = 1; /* scalar w/ coords and SA (not needed) */
53     else if(a_coords != 0) pc_gamg->m_data_cols = (a_ndm==2 ? 3 : 6); /* elasticity */
54     else pc_gamg->m_data_cols = bs; /* no data, force SA with constant null space vectors */
55     pc_gamg->m_data_rows = bs;
56   }
57   arrsz = nloc*pc_gamg->m_data_rows*pc_gamg->m_data_cols;
58 
59   /* create data - syntactic sugar that should be refactored at some point */
60   if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) {
61     ierr = PetscFree( pc_gamg->m_data );  CHKERRQ(ierr);
62     ierr = PetscMalloc((arrsz+1)*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr);
63   }
64   for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.;
65   pc_gamg->m_data[arrsz] = -99.;
66   /* copy data in - column oriented */
67   if( pc_gamg->m_useSA ) {
68     const PetscInt M = Iend - my0;
69     for(kk=0;kk<nloc;kk++){
70       PetscReal *data = &pc_gamg->m_data[kk*bs];
71       if( pc_gamg->m_data_cols==1 ) *data = 1.0;
72       else {
73         for(ii=0;ii<bs;ii++)
74 	  for(jj=0;jj<bs;jj++)
75 	    if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */
76 	    else data[ii*M + jj] = 0.0;
77         if( a_coords != 0 ) {
78           if( a_ndm == 2 ){ /* rotational modes */
79             data += 2*M;
80             data[0] = -a_coords[2*kk+1];
81             data[1] =  a_coords[2*kk];
82           }
83           else {
84             data += 3*M;
85             data[0] = 0.0;               data[M+0] =  a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1];
86             data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0;               data[2*M+1] =  a_coords[3*kk];
87             data[2] =  a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk];   data[2*M+2] = 0.0;
88           }
89         }
90       }
91     }
92   }
93   else {
94     for( kk = 0 ; kk < nloc ; kk++ ){
95       for( ii = 0 ; ii < a_ndm ; ii++ ) {
96         pc_gamg->m_data[ii*nloc + kk] =  a_coords[kk*a_ndm + ii];
97       }
98     }
99   }
100 
101   assert(pc_gamg->m_data[arrsz] == -99.);
102   for(kk=0;kk<arrsz;kk++) assert(pc_gamg->m_data[kk] != -999.); // debug
103 
104   pc_gamg->m_data_sz = arrsz;
105   pc_gamg->m_dim = a_ndm;
106 
107   PetscFunctionReturn(0);
108 }
109 EXTERN_C_END
110 
111 
112 /* -----------------------------------------------------------------------------*/
113 #undef __FUNCT__
114 #define __FUNCT__ "PCReset_GAMG"
115 PetscErrorCode PCReset_GAMG(PC pc)
116 {
117   PetscErrorCode  ierr;
118   PC_MG           *mg = (PC_MG*)pc->data;
119   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
120 
121   PetscFunctionBegin;
122   ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr);
123   pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0;
124   PetscFunctionReturn(0);
125 }
126 
127 /* -------------------------------------------------------------------------- */
128 /*
129    partitionLevel
130 
131    Input Parameter:
132    . a_Amat_fine - matrix on this fine (k) level
133    . a_ndata_rows - size of data to move (coarse grid)
134    . a_ndata_cols - size of data to move (coarse grid)
135    In/Output Parameter:
136    . a_P_inout - prolongation operator to the next level (k-1)
137    . a_coarse_data - data that need to be moved
138    . a_active_proc - number of active procs
139    Output Parameter:
140    . a_Amat_crs - coarse matrix that is created (k-1)
141 */
142 
143 #undef __FUNCT__
144 #define __FUNCT__ "partitionLevel"
145 PetscErrorCode partitionLevel( Mat a_Amat_fine,
146                                PetscInt a_ndata_rows,
147                                PetscInt a_ndata_cols,
148 			       PetscInt a_cbs,
149                                Mat *a_P_inout,
150                                PetscReal **a_coarse_data,
151                                PetscMPIInt *a_active_proc,
152                                Mat *a_Amat_crs
153                                )
154 {
155   PetscErrorCode   ierr;
156   Mat              Cmat,Pnew,Pold=*a_P_inout;
157   IS               new_indices,isnum;
158   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
159   PetscMPIInt      nactive_procs,mype,npe;
160   PetscInt         Istart,Iend,Istart0,Iend0,ncrs0,ncrs_new,fbs;
161   PetscInt         neq,NN;
162   PetscMPIInt      new_npe,targ_npe;
163   PetscBool        flag = PETSC_FALSE;
164 
165   PetscFunctionBegin;
166   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
167   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
168   ierr = MatGetBlockSize( a_Amat_fine, &fbs ); CHKERRQ(ierr);
169   /* RAP */
170   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
171 
172   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
173 
174   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
175   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
176 
177   /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
178   ierr = MatGetSize( Cmat, &neq, &NN );CHKERRQ(ierr);
179 #define MIN_EQ_PROC 500
180   nactive_procs = *a_active_proc;
181   targ_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */
182 #define TOP_GRID_LIM 500
183   if( targ_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1; /* chop coarsest grid */
184   else if (targ_npe >= nactive_procs ) new_npe = nactive_procs; /* no change */
185   else {
186     PetscMPIInt factstart,fact;
187     new_npe = -9999;
188     factstart = nactive_procs;
189     for(fact=factstart;fact>0;fact--){ /* try to find a better number of procs */
190       if( nactive_procs%fact==0 && neq/(nactive_procs/fact) > MIN_EQ_PROC ) {
191         new_npe = nactive_procs/fact;
192       }
193     }
194     if(new_npe == -9999) new_npe = nactive_procs;
195   }
196 
197   ierr  = PetscOptionsHasName(PETSC_NULL,"-pc_gamg_avoid_repartitioning",&flag);
198   CHKERRQ( ierr );
199   if( !flag ) { /* re-partition */
200     MatPartitioning  mpart;
201     Mat              adj;
202     const PetscInt  *is_idx;
203     PetscInt         is_sz,kk,jj,ii,old_fact=(npe/nactive_procs), *isnewproc_idx;
204     /* create sub communicator  */
205     MPI_Comm         cm,new_comm;
206     PetscInt         membershipKey = mype%old_fact, new_fact=(npe/new_npe),counts[npe];
207     IS               isnewproc;
208 
209     *a_active_proc = new_npe; /* output for next level */
210     ierr = MPI_Comm_split(wcomm, membershipKey, mype, &cm); CHKERRQ(ierr);
211     ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr);
212     ierr = MPI_Comm_free( &cm );                            CHKERRQ(ierr);
213 
214     /* MatPartitioningApply call MatConvert, which is collective */
215     ierr = PetscLogEventBegin(gamg_setup_stages[SET12],0,0,0,0);CHKERRQ(ierr);
216     if( a_cbs==1) {
217       ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
218     }
219     else {
220       /* make a scalar matrix to partition */
221       Mat tMat;
222       PetscInt Ii,ncols; const PetscScalar *vals; const PetscInt *idx;
223       MatInfo info;
224       ierr = MatGetInfo(Cmat,MAT_LOCAL,&info); CHKERRQ(ierr);
225       ncols = (PetscInt)info.nz_used/((ncrs0+1)*a_cbs*a_cbs)+1;
226 
227       ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
228                               PETSC_DETERMINE, PETSC_DETERMINE,
229                               2*ncols, PETSC_NULL, ncols, PETSC_NULL,
230                               &tMat );
231       CHKERRQ(ierr);
232 
233       for ( Ii = Istart0; Ii < Iend0; Ii++ ) {
234         PetscInt dest_row = Ii/a_cbs;
235         ierr = MatGetRow(Cmat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr);
236         for( jj = 0 ; jj < ncols ; jj++ ){
237           PetscInt dest_col = idx[jj]/a_cbs;
238           PetscScalar v = 1.0;
239           ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
240         }
241         ierr = MatRestoreRow(Cmat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr);
242       }
243       ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
244       ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
245 
246       ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
247 
248       ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
249     }
250     if( membershipKey == 0 ) {
251       if(ncrs0==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"zero local nodes -- increase min");
252       /* hack to fix global data that pmetis.c uses in 'adj' */
253       for( kk=0 , jj=0 ; kk<=npe ; jj++, kk += old_fact ) {
254         adj->rmap->range[jj] = adj->rmap->range[kk];
255       }
256       ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr);
257       ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
258       ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr);
259       ierr = MatPartitioningSetNParts( mpart, new_npe ); CHKERRQ(ierr);
260       ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr);
261       ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr);
262       /* collect IS info */
263       ierr = ISGetLocalSize( isnewproc, &is_sz );        CHKERRQ(ierr);
264       ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr);
265       ierr = ISGetIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
266       /* spread partitioning across machine - probably the right thing to do but machine spec. */
267       for(kk=0,jj=0;kk<is_sz;kk++){
268         for(ii=0 ; ii<a_cbs ; ii++, jj++ ) { /* expand for equation level by 'bs' */
269           isnewproc_idx[jj] = is_idx[kk] * new_fact;
270         }
271       }
272       ierr = ISRestoreIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
273       ierr = ISDestroy( &isnewproc );                    CHKERRQ(ierr);
274       is_sz *= a_cbs;
275     }
276     else {
277       isnewproc_idx = 0;
278       is_sz = 0;
279     }
280     ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
281     ierr = MPI_Comm_free( &new_comm );    CHKERRQ(ierr);
282     ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc );
283     if( membershipKey == 0 ) {
284       ierr = PetscFree( isnewproc_idx );  CHKERRQ(ierr);
285     }
286 
287     /*
288      Create an index set from the isnewproc index set to indicate the mapping TO
289      */
290     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
291     /*
292      Determine how many elements are assigned to each processor
293      */
294     ierr = ISPartitioningCount( isnewproc, npe, counts ); CHKERRQ(ierr);
295     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
296     ncrs_new = counts[mype]/a_cbs;
297     ierr = PetscLogEventEnd(gamg_setup_stages[SET12],0,0,0,0);   CHKERRQ(ierr);
298 
299     { /* Create a vector to contain the newly ordered element information */
300       const PetscInt *idx, data_sz=a_ndata_rows*a_ndata_cols;
301       const PetscInt  stride0=ncrs0*a_ndata_rows,strideNew=ncrs_new*a_ndata_rows;
302       PetscInt        ii,jj,kk;
303       IS              isscat;
304       PetscScalar    *array;
305       Vec             src_crd, dest_crd;
306       PetscReal      *data = *a_coarse_data;
307       VecScatter      vecscat;
308       PetscInt        tidx[ncrs0*data_sz];
309 
310       ierr = VecCreate( wcomm, &dest_crd );
311       ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
312       ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/
313       /*
314        There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
315        a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
316        */
317       ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
318       for(ii=0,jj=0; ii<ncrs0 ; ii++) {
319         PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
320         for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
321       }
322       ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
323       ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
324       CHKERRQ(ierr);
325       /*
326        Create a vector to contain the original vertex information for each element
327        */
328       ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
329       for( jj=0; jj<a_ndata_cols ; jj++ ) {
330         for( ii=0 ; ii<ncrs0 ; ii++) {
331           for( kk=0; kk<a_ndata_rows ; kk++ ) {
332             PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
333             ierr = VecSetValues( src_crd, 1, &jx, &data[ix], INSERT_VALUES );  CHKERRQ(ierr);
334           }
335         }
336       }
337       ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
338       ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
339       /*
340        Scatter the element vertex information (still in the original vertex ordering)
341        to the correct processor
342        */
343       ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
344       CHKERRQ(ierr);
345       ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
346       ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
347       ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
348       ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
349       ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
350       /*
351        Put the element vertex data into a new allocation of the gdata->ele
352        */
353       ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
354       ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
355 
356       ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
357       data = *a_coarse_data;
358       for( jj=0; jj<a_ndata_cols ; jj++ ) {
359         for( ii=0 ; ii<ncrs_new ; ii++) {
360           for( kk=0; kk<a_ndata_rows ; kk++ ) {
361             PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
362             data[ix] = PetscRealPart(array[jx]);
363             array[jx] = 1.e300;
364           }
365         }
366       }
367       ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
368       ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
369     }
370 
371     /*
372      Invert for MatGetSubMatrix
373      */
374     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
375     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
376     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
377     /* A_crs output */
378     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
379     CHKERRQ(ierr);
380 
381     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
382     Cmat = *a_Amat_crs; /* output */
383     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
384 
385     /* prolongator */
386     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
387     {
388       IS findices;
389       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
390       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
391       CHKERRQ(ierr);
392       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
393     }
394     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
395     *a_P_inout = Pnew; /* output */
396     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
397   }
398   else {
399     *a_Amat_crs = Cmat; /* output */
400   }
401 
402   PetscFunctionReturn(0);
403 }
404 
405 #define GAMG_MAXLEVELS 30
406 #if defined(PETSC_USE_LOG)
407 PetscLogStage  gamg_stages[20];
408 #endif
409 /* -------------------------------------------------------------------------- */
410 /*
411    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
412                     by setting data structures and options.
413 
414    Input Parameter:
415 .  pc - the preconditioner context
416 
417    Application Interface Routine: PCSetUp()
418 
419    Notes:
420    The interface routine PCSetUp() is not usually called directly by
421    the user, but instead is called by PCApply() if necessary.
422 */
423 #undef __FUNCT__
424 #define __FUNCT__ "PCSetUp_GAMG"
425 PetscErrorCode PCSetUp_GAMG( PC a_pc )
426 {
427   PetscErrorCode  ierr;
428   PC_MG           *mg = (PC_MG*)a_pc->data;
429   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
430   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
431   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
432   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
433   PetscMPIInt      mype,npe,nactivepe;
434   PetscBool        isOK;
435   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
436   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
437   MatInfo          info;
438 
439   PetscFunctionBegin;
440   if( a_pc->setupcalled ) {
441     /* no state data in GAMG to destroy */
442     ierr = PCReset_MG( a_pc ); CHKERRQ(ierr);
443   }
444   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
445   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
446   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
447   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
448   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
449   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
450   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
451 
452   /* get data of not around */
453   if( pc_gamg->m_data == 0 && nloc > 0 ) {
454     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
455   }
456   data = pc_gamg->m_data;
457 
458   /* Get A_i and R_i */
459   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
460   PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
461 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,(PetscInt)info.nz_used/N,npe);
462   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
463         level < GAMG_MAXLEVELS-1 && (level==0 || M/bs>TOP_GRID_LIM) && (npe==1 || nactivepe>1);
464         level++ ){
465     level1 = level + 1;
466     ierr = PetscLogEventBegin(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr);
467     ierr = createProlongation(Aarr[level], data, pc_gamg->m_dim, pc_gamg->m_data_cols, pc_gamg->m_useSA,
468                               &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level] );
469     CHKERRQ(ierr);
470     ierr = PetscFree( data ); CHKERRQ( ierr );
471     ierr = PetscLogEventEnd(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr);
472 
473     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
474 
475     if( isOK ) {
476       ierr = PetscLogEventBegin(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr);
477       ierr = partitionLevel( Aarr[level], pc_gamg->m_useSA ? bs : 1, pc_gamg->m_data_cols, bs,
478                              &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
479       CHKERRQ(ierr);
480       ierr = PetscLogEventEnd(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr);
481       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
482       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
483       PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, bs=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
484 		  mype,__FUNCT__,level1,N,bs,pc_gamg->m_data_cols,(PetscInt)info.nz_used/N,nactivepe);
485       /* coarse grids with SA can have zero row/cols from singleton aggregates */
486       /* aggregation method can probably gaurrentee this does not happen! - be safe for now */
487 
488       if( PETSC_TRUE ){
489         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloc,id;
490         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
491         nloc = Iend-Istart;
492         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
493         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
494         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
495         for(kk=0;kk<nloc;kk++){
496           if(data_arr[kk]==0.0) {
497             id = kk + Istart; v = 1.e-1;
498             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
499             CHKERRQ(ierr);
500             /* PetscPrintf(PETSC_COMM_SELF,"[%d]%s warning: added diag to zero (%d) on level %d \n",mype,__FUNCT__,id,level); */
501           }
502         }
503         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
504         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
505         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
506         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
507       }
508     }
509     else{
510       break;
511     }
512     data = coarse_data;
513   }
514   ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
515   PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
516   pc_gamg->m_data = 0; /* destroyed coordinate data */
517   pc_gamg->m_Nlevels = level + 1;
518   fine_level = level;
519   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
520 
521   /* set default smoothers */
522   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
523         lidx <= fine_level;
524         lidx++, level--) {
525     PetscReal emax, emin;
526     KSP smoother; PC subpc;
527     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
528     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
529     if( emaxs[level] > 0.0 ) emax = emaxs[level];
530     else{ /* eigen estimate 'emax' */
531       KSP eksp; Mat Lmat = Aarr[level];
532       Vec bb, xx; PC pc;
533 
534       ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
535       ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
536       {
537 	PetscRandom    rctx;
538 	ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
539 	ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
540 	ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
541 	ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
542       }
543       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
544       ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr);
545       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
546       ierr = KSPSetOperators( eksp, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN ); CHKERRQ( ierr );
547       ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
548       ierr = PCSetType( pc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); /* should be same as above */
549       ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
550       CHKERRQ(ierr);
551       //ierr = KSPSetConvergenceTest( eksp, KSPSkipConverged, 0, 0 ); CHKERRQ(ierr);
552       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
553 
554       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
555       ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
556       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
557       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
558       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
559       ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
560       PetscPrintf(PETSC_COMM_WORLD,"%s max eigen = %e min = %e\n",__FUNCT__,emax,emin);
561     }
562     {
563       PetscInt N1, N0, tt;
564       ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
565       ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
566       emin = .5*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */
567       emax *= 1.05;
568     }
569 
570     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], DIFFERENT_NONZERO_PATTERN );
571     ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
572     ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr);
573     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
574     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
575     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
576   }
577   {
578     KSP smoother; /* coarse grid */
579     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
580     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
581     ierr = KSPSetOperators( smoother, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN );
582     CHKERRQ(ierr);
583     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
584   }
585 
586   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
587   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
588   {
589     PetscBool galerkin;
590     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
591     if(galerkin){
592       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
593     }
594   }
595 
596   /* set interpolation between the levels, create timer stages, clean up */
597   if( PETSC_FALSE ) {
598     char str[32];
599     sprintf(str,"MG Level %d (%d)",0,pc_gamg->m_Nlevels-1);
600     PetscLogStageRegister(str, &gamg_stages[fine_level]);
601   }
602   for (lidx=0,level=pc_gamg->m_Nlevels-1;
603        lidx<fine_level;
604        lidx++, level--){
605     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
606     if( !PETSC_TRUE ) {
607       PetscViewer viewer; char fname[32];
608       sprintf(fname,"Pmat_%d.m",level);
609       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
610       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
611       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
612       ierr = PetscViewerDestroy( &viewer );
613       sprintf(fname,"Amat_%d.m",level);
614       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
615       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
616       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
617       ierr = PetscViewerDestroy( &viewer );
618     }
619     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
620     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
621     if( PETSC_FALSE ) {
622       char str[32];
623       sprintf(str,"MG Level %d (%d)",lidx+1,level-1);
624       PetscLogStageRegister(str, &gamg_stages[level-1]);
625     }
626   }
627 
628   /* setupcalled is set to 0 so that MG is setup from scratch */
629   a_pc->setupcalled = 0;
630   ierr = PCSetUp_MG(a_pc);CHKERRQ(ierr);
631 
632   PetscFunctionReturn(0);
633 }
634 
635 /* ------------------------------------------------------------------------- */
636 /*
637    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
638    that was created with PCCreate_GAMG().
639 
640    Input Parameter:
641 .  pc - the preconditioner context
642 
643    Application Interface Routine: PCDestroy()
644 */
645 #undef __FUNCT__
646 #define __FUNCT__ "PCDestroy_GAMG"
647 PetscErrorCode PCDestroy_GAMG(PC pc)
648 {
649   PetscErrorCode  ierr;
650   PC_MG           *mg = (PC_MG*)pc->data;
651   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
652 
653   PetscFunctionBegin;
654   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
655   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
656   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
657   PetscFunctionReturn(0);
658 }
659 
660 #undef __FUNCT__
661 #define __FUNCT__ "PCSetFromOptions_GAMG"
662 PetscErrorCode PCSetFromOptions_GAMG(PC pc)
663 {
664   /* PetscErrorCode  ierr; */
665   /* PC_MG           *mg = (PC_MG*)pc->data; */
666   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
667   /* MPI_Comm        comm = ((PetscObject)pc)->comm; */
668 
669   PetscFunctionBegin;
670   PetscFunctionReturn(0);
671 }
672 
673 /* -------------------------------------------------------------------------- */
674 /*
675  PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG
676 
677    Input Parameter:
678 .  pc - the preconditioner context
679 
680    Application Interface Routine: PCCreate()
681 
682   */
683  /* MC
684      PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide
685        fine grid discretization matrix and coordinates on the fine grid.
686 
687    Options Database Key:
688    Multigrid options(inherited)
689 +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
690 .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
691 .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
692    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
693    GAMG options:
694 
695    Level: intermediate
696   Concepts: multigrid
697 
698 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
699            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
700            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
701            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
702            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
703 M */
704 
705 EXTERN_C_BEGIN
706 #undef __FUNCT__
707 #define __FUNCT__ "PCCreate_GAMG"
708 PetscErrorCode  PCCreate_GAMG(PC pc)
709 {
710   PetscErrorCode  ierr;
711   PC_GAMG         *pc_gamg;
712   PC_MG           *mg;
713   PetscClassId     cookie;
714 
715   PetscFunctionBegin;
716   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
717   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
718   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
719 
720   /* create a supporting struct and attach it to pc */
721   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
722   mg = (PC_MG*)pc->data;
723   mg->innerctx = pc_gamg;
724 
725   pc_gamg->m_Nlevels    = -1;
726 
727   /* overwrite the pointers of PCMG by the functions of PCGAMG */
728   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
729   pc->ops->setup          = PCSetUp_GAMG;
730   pc->ops->reset          = PCReset_GAMG;
731   pc->ops->destroy        = PCDestroy_GAMG;
732 
733   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
734 					    "PCSetCoordinates_C",
735 					    "PCSetCoordinates_GAMG",
736 					    PCSetCoordinates_GAMG);CHKERRQ(ierr);
737 
738   PetscClassIdRegister("GAMG Setup",&cookie);
739   PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_stages[SET1]);
740   PetscLogEventRegister(" make graph", cookie, &gamg_setup_stages[SET3]);
741   PetscLogEventRegister(" MIS/Agg", cookie, &gamg_setup_stages[SET4]);
742   PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_stages[SET5]);
743   PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_stages[SET6]);
744   PetscLogEventRegister("   search & set", cookie, &gamg_setup_stages[FIND_V]);
745   PetscLogEventRegister("  SA: init", cookie, &gamg_setup_stages[SET7]);
746   /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_stages[SET8]); */
747   PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_stages[SET9]);
748   PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_stages[SET2]);
749   PetscLogEventRegister(" PL repartition", cookie, &gamg_setup_stages[SET12]);
750   /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_stages[SET13]); */
751   /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_stages[SET10]); */
752   /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_stages[SET11]); */
753 
754   PetscFunctionReturn(0);
755 }
756 EXTERN_C_END
757