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