xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision e4bbc6d007f0412834508dfd030a46bcf377ed51)
1 /*
2  GAMG geometric-algebric multigrid PC - Mark Adams 2011
3  */
4 #include "petsc-private/matimpl.h"
5 #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
6 #include <petsc-private/kspimpl.h>
7 
8 #if defined PETSC_GAMG_USE_LOG
9 PetscLogEvent petsc_gamg_setup_events[NUM_SET];
10 #endif
11 
12 #if defined PETSC_USE_LOG
13 PetscLogEvent PC_GAMGGgraph_AGG;
14 PetscLogEvent PC_GAMGGgraph_GEO;
15 PetscLogEvent PC_GAMGCoarsen_AGG;
16 PetscLogEvent PC_GAMGCoarsen_GEO;
17 PetscLogEvent PC_GAMGProlongator_AGG;
18 PetscLogEvent PC_GAMGProlongator_GEO;
19 PetscLogEvent PC_GAMGOptprol_AGG;
20 #endif
21 
22 #define GAMG_MAXLEVELS 30
23 
24 /* #define GAMG_STAGES */
25 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
26 static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
27 #endif
28 
29 static PetscFList GAMGList = 0;
30 
31 /* ----------------------------------------------------------------------------- */
32 #undef __FUNCT__
33 #define __FUNCT__ "PCReset_GAMG"
34 PetscErrorCode PCReset_GAMG(PC pc)
35 {
36   PetscErrorCode  ierr;
37   PC_MG           *mg = (PC_MG*)pc->data;
38   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
39 
40   PetscFunctionBegin;
41   if( pc_gamg->data != 0 ) { /* this should not happen, cleaned up in SetUp */
42     ierr = PetscFree(pc_gamg->data); CHKERRQ(ierr);
43   }
44   pc_gamg->data = 0; pc_gamg->data_sz = 0;
45   PetscFunctionReturn(0);
46 }
47 
48 /* -------------------------------------------------------------------------- */
49 /*
50    createLevel: create coarse op with RAP.  repartition and/or reduce number
51      of active processors.
52 
53    Input Parameter:
54    . pc - parameters
55    . Amat_fine - matrix on this fine (k) level
56    . cbs - coarse block size
57    . isLast -
58    In/Output Parameter:
59    . a_P_inout - prolongation operator to the next level (k-1)
60    . a_nactive_proc - number of active procs
61    Output Parameter:
62    . a_Amat_crs - coarse matrix that is created (k-1)
63 */
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "createLevel"
67 static PetscErrorCode createLevel( const PC pc,
68                             const Mat Amat_fine,
69                             const PetscInt cbs,
70                             const PetscBool isLast,
71                             Mat *a_P_inout,
72                             PetscMPIInt *a_nactive_proc,
73                             Mat *a_Amat_crs
74                             )
75 {
76   PC_MG           *mg = (PC_MG*)pc->data;
77   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
78   const PetscBool  repart = pc_gamg->repart;
79   const PetscInt   min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
80   PetscErrorCode   ierr;
81   Mat              Cmat,Pnew,Pold=*a_P_inout;
82   IS               new_indices,isnum;
83   MPI_Comm         wcomm = ((PetscObject)Amat_fine)->comm;
84   PetscMPIInt      mype,npe,new_npe,nactive = *a_nactive_proc;
85   PetscInt         neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0,rfactor;
86 
87   PetscFunctionBegin;
88   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
89   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
90 
91   /* RAP */
92 #ifdef USE_R
93   /* make R wih brute force for now */
94   ierr = MatTranspose( Pold, Pnew );
95   ierr = MatDestroy( &Pold );  CHKERRQ(ierr);
96   ierr = MatRARt( Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
97   Pold = Pnew;
98 #else
99   ierr = MatPtAP( Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
100 #endif
101   ierr = MatSetBlockSize( Cmat, cbs );      CHKERRQ(ierr);
102   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
103   ncrs0 = (Iend0-Istart0)/cbs; assert((Iend0-Istart0)%cbs == 0);
104 
105   /* get number of PEs to make active, reduce */
106   ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
107   new_npe = (PetscMPIInt)((float)neq/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
108   if( new_npe == 0 || neq < coarse_max ) new_npe = 1;
109   else if (new_npe >= nactive ) new_npe = nactive; /* no change, rare */
110   if( isLast ) new_npe = 1;
111 
112   if( !repart && new_npe==nactive ) {
113     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction */
114   }
115   else {
116     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
117     const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,data_sz=ndata_rows*ndata_cols;
118     const PetscInt  stride0=ncrs0*pc_gamg->data_cell_rows;
119     PetscInt        *counts,is_sz,*newproc_idx,ii,jj,kk,strideNew,*tidx,inpe,targetPE;
120     IS              isnewproc;
121     VecScatter      vecscat;
122     PetscScalar    *array;
123     Vec             src_crd, dest_crd;
124     IS              isscat;
125 
126     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
127 #if defined PETSC_GAMG_USE_LOG
128       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
129 #endif
130     if( repart ) {
131       /* create sub communicator  */
132       Mat             adj;
133 
134       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,neq);
135 
136       /* MatPartitioningApply call MatConvert, which is collective */
137       if( cbs == 1 ) {
138 	ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
139       }
140       else{
141 	/* make a scalar matrix to partition */
142 	Mat tMat;
143 	PetscInt ncols,jj,Ii;
144 	const PetscScalar *vals;
145 	const PetscInt *idx;
146 	PetscInt *d_nnz, *o_nnz;
147 	static PetscInt llev = 0;
148 
149 	ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
150 	ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr);
151 	for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += cbs, jj++ ) {
152 	  ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
153 	  d_nnz[jj] = ncols/cbs;
154 	  o_nnz[jj] = ncols/cbs;
155 	  ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
156 	  if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
157 	  if( o_nnz[jj] > (neq/cbs-ncrs0) ) o_nnz[jj] = neq/cbs-ncrs0;
158 	}
159 
160 	ierr = MatCreateAIJ( wcomm, ncrs0, ncrs0,
161                              PETSC_DETERMINE, PETSC_DETERMINE,
162                              0, d_nnz, 0, o_nnz,
163                              &tMat );
164 	CHKERRQ(ierr);
165 	ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
166 	ierr = PetscFree( o_nnz ); CHKERRQ(ierr);
167 
168 	for ( ii = Istart0; ii < Iend0; ii++ ) {
169 	  PetscInt dest_row = ii/cbs;
170 	  ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
171 	  for( jj = 0 ; jj < ncols ; jj++ ){
172 	    PetscInt dest_col = idx[jj]/cbs;
173 	    PetscScalar v = 1.0;
174 	    ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
175 	  }
176 	  ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
177 	}
178 	ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
179 	ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
180 
181 	if( llev++ == -1 ) {
182 	  PetscViewer viewer; char fname[32];
183 	  ierr = PetscSNPrintf(fname,sizeof fname,"part_mat_%D.mat",llev);CHKERRQ(ierr);
184 	  PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
185 	  ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
186 	  ierr = PetscViewerDestroy( &viewer );
187 	}
188 
189 	ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
190 
191 	ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
192       }
193 
194       { /* partition: get newproc_idx, set is_sz */
195 	char prefix[256];
196 	const char *pcpre;
197 	const PetscInt *is_idx;
198 	MatPartitioning  mpart;
199 	IS proc_is;
200 
201 	ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr);
202 	ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
203 	ierr = PCGetOptionsPrefix( pc, &pcpre );CHKERRQ(ierr);
204 	ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
205 	ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
206 	ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
207 	ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
208 	ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr);
209 	ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
210 
211 	/* collect IS info */
212 	ierr = ISGetLocalSize( proc_is, &is_sz );       CHKERRQ(ierr);
213 	ierr = PetscMalloc( cbs*is_sz*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr);
214 	ierr = ISGetIndices( proc_is, &is_idx );        CHKERRQ(ierr);
215 	NN = 1; /* bring to "front" of machine */
216 	/*NN = npe/new_npe;*/ /* spread partitioning across machine */
217 	for( kk = jj = 0 ; kk < is_sz ; kk++ ){
218 	  for( ii = 0 ; ii < cbs ; ii++, jj++ ){
219 	    newproc_idx[jj] = is_idx[kk] * NN; /* distribution */
220 	  }
221 	}
222 	ierr = ISRestoreIndices( proc_is, &is_idx );     CHKERRQ(ierr);
223 	ierr = ISDestroy( &proc_is );                    CHKERRQ(ierr);
224 
225 	is_sz *= cbs;
226       }
227       ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
228 
229       ierr = ISCreateGeneral( wcomm, is_sz, newproc_idx, PETSC_COPY_VALUES, &isnewproc );
230       CHKERRQ(ierr);
231       if( newproc_idx != 0 ) {
232 	ierr = PetscFree( newproc_idx );  CHKERRQ(ierr);
233       }
234     }
235     else { /* simple aggreagtion of parts */
236       /* find factor */
237       if( new_npe == 1 ) rfactor = npe; /* easy */
238       else {
239 	PetscReal best_fact = 0.;
240 	jj = -1;
241 	for( kk = 1 ; kk <= npe ; kk++ ){
242 	  if( npe%kk==0 ) { /* a candidate */
243 	    PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe;
244 	    if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */
245 	    if( fact > best_fact ) {
246 	      best_fact = fact; jj = kk;
247 	    }
248 	  }
249 	}
250 	if(jj!=-1) rfactor = jj;
251 	else rfactor = 1; /* prime? */
252       }
253       new_npe = npe/rfactor;
254 
255       if( new_npe==nactive ) {
256 	*a_Amat_crs = Cmat; /* output - no repartitioning or reduction */
257 	ierr = PetscFree( counts );  CHKERRQ(ierr);
258 	*a_nactive_proc = new_npe; /* output */
259 	if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq=%d\n",mype,__FUNCT__,new_npe,neq);
260 	PetscFunctionReturn(0);
261       }
262 
263       /* reduce - isnewproc */
264       targetPE = mype/rfactor;
265 
266       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s aggregate processors: npe: %d --> %d, neq=%d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,neq);
267       is_sz = ncrs0*cbs;
268       ierr = ISCreateStride( wcomm, is_sz, targetPE, 0, &isnewproc );
269       CHKERRQ(ierr);
270     }
271 
272     /*
273       Create an index set from the isnewproc index set to indicate the mapping TO
274     */
275     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
276     /*
277       Determine how many elements are assigned to each processor
278     */
279     inpe = npe;
280     ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr);
281     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
282     ncrs_new = counts[mype]/cbs;
283     strideNew = ncrs_new*ndata_rows;
284 #if defined PETSC_GAMG_USE_LOG
285       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
286 #endif
287     /* Create a vector to contain the newly ordered element information */
288     ierr = VecCreate( wcomm, &dest_crd );
289     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
290     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
291     /*
292      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
293      a block size of ...).  Note, ISs are expanded into equation space by 'cbs'.
294      */
295     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
296     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
297     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
298       PetscInt id = idx[ii*cbs]/cbs; /* get node back */
299       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
300     }
301     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
302     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
303     CHKERRQ(ierr);
304     ierr = PetscFree( tidx );  CHKERRQ(ierr);
305     /*
306      Create a vector to contain the original vertex information for each element
307      */
308     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
309     for( jj=0; jj<ndata_cols ; jj++ ) {
310       for( ii=0 ; ii<ncrs0 ; ii++) {
311 	for( kk=0; kk<ndata_rows ; kk++ ) {
312 	  PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*ndata_cols + jj;
313           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
314 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
315 	}
316       }
317     }
318     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
319     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
320     /*
321       Scatter the element vertex information (still in the original vertex ordering)
322       to the correct processor
323     */
324     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
325     CHKERRQ(ierr);
326     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
327     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
328     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
329     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
330     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
331     /*
332       Put the element vertex data into a new allocation of the gdata->ele
333     */
334     ierr = PetscFree( pc_gamg->data );    CHKERRQ(ierr);
335     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), &pc_gamg->data );    CHKERRQ(ierr);
336     pc_gamg->data_sz = data_sz*ncrs_new;
337     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
338     for( jj=0; jj<ndata_cols ; jj++ ) {
339       for( ii=0 ; ii<ncrs_new ; ii++) {
340 	for( kk=0; kk<ndata_rows ; kk++ ) {
341 	  PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*ndata_cols + jj;
342 	  pc_gamg->data[ix] = PetscRealPart(array[jx]);
343 	  array[jx] = 1.e300;
344 	}
345       }
346     }
347     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
348     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
349 #if defined PETSC_GAMG_USE_LOG
350     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
351 #endif
352     /*
353       Invert for MatGetSubMatrix
354     */
355     ierr = ISInvertPermutation( isnum, ncrs_new*cbs, &new_indices ); CHKERRQ(ierr);
356     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
357     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
358 #if defined PETSC_GAMG_USE_LOG
359     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
360     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
361 #endif
362     /* A_crs output */
363     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
364     CHKERRQ(ierr);
365 
366     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
367     Cmat = *a_Amat_crs; /* output */
368     ierr = MatSetBlockSize( Cmat, cbs );      CHKERRQ(ierr);
369 #if defined PETSC_GAMG_USE_LOG
370     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
371 #endif
372     /* prolongator */
373     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
374     {
375       IS findices;
376 #if defined PETSC_GAMG_USE_LOG
377       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
378 #endif
379       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
380 #ifdef USE_R
381       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
382 #else
383       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
384 #endif
385       CHKERRQ(ierr);
386       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
387 #if defined PETSC_GAMG_USE_LOG
388       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
389 #endif
390     }
391     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
392     *a_P_inout = Pnew; /* output - repartitioned */
393     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
394     ierr = PetscFree( counts );  CHKERRQ(ierr);
395   }
396 
397   *a_nactive_proc = new_npe; /* output */
398 
399   if( !PETSC_TRUE ) {
400     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
401     if(llev==0) {
402       sprintf(fname,"Cmat_%d.m",llev++);
403       PetscViewerASCIIOpen(wcomm,fname,&viewer);
404       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
405       ierr = MatView(Amat_fine, viewer ); CHKERRQ(ierr);
406       ierr = PetscViewerDestroy( &viewer );
407     }
408     sprintf(fname,"Cmat_%d.m",llev++);
409     PetscViewerASCIIOpen(wcomm,fname,&viewer);
410     ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
411     ierr = MatView(Cmat, viewer ); CHKERRQ(ierr);
412     ierr = PetscViewerDestroy( &viewer );
413   }
414 
415   PetscFunctionReturn(0);
416 }
417 
418 /* -------------------------------------------------------------------------- */
419 /*
420    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
421                     by setting data structures and options.
422 
423    Input Parameter:
424 .  pc - the preconditioner context
425 
426    Application Interface Routine: PCSetUp()
427 
428    Notes:
429    The interface routine PCSetUp() is not usually called directly by
430    the user, but instead is called by PCApply() if necessary.
431 */
432 #undef __FUNCT__
433 #define __FUNCT__ "PCSetUp_GAMG"
434 PetscErrorCode PCSetUp_GAMG( PC pc )
435 {
436   PetscErrorCode  ierr;
437   PC_MG           *mg = (PC_MG*)pc->data;
438   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
439   Mat              Pmat = pc->pmat;
440   PetscInt         fine_level,level,level1,M,N,bs,nloc,lidx,Istart,Iend,nASMBlocksArr[GAMG_MAXLEVELS];
441   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
442   PetscMPIInt      mype,npe,nactivepe;
443   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
444   PetscReal        emaxs[GAMG_MAXLEVELS];
445   IS              *ASMLocalIDsArr[GAMG_MAXLEVELS];
446   PetscLogDouble   nnz0=0,nnztot=0;
447   MatInfo          info;
448 
449   PetscFunctionBegin;
450   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
451   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
452   if( pc_gamg->setup_count++ > 0 ) {
453     PC_MG_Levels   **mglevels = mg->levels;
454     /* just do Galerkin grids */
455     Mat B,dA,dB;
456     assert(pc->setupcalled);
457 
458     if( pc_gamg->Nlevels > 1 ) {
459       /* currently only handle case where mat and pmat are the same on coarser levels */
460       ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
461       /* (re)set to get dirty flag */
462       ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
463 
464       for (level=pc_gamg->Nlevels-2; level>-1; level--) {
465         /* the first time through the matrix structure has changed from repartitioning */
466         if( pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
467           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
468           ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
469           mglevels[level]->A = B;
470         }
471         else {
472           ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
473           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
474         }
475         ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
476         dB = B;
477         }
478     }
479 
480     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
481 
482     /* PCSetUp_MG seems to insists on setting this to GMRES */
483     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
484 
485     PetscFunctionReturn(0);
486   }
487   assert(pc->setupcalled == 0);
488 
489   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
490   ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr);
491   ierr = MatGetSize( Pmat, &M, &N );CHKERRQ(ierr);
492   ierr = MatGetOwnershipRange( Pmat, &Istart, &Iend ); CHKERRQ(ierr);
493   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
494 
495   if( pc_gamg->data == 0 && nloc > 0 ) {
496     if(!pc_gamg->createdefaultdata){
497       SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"'createdefaultdata' not set?!?! need to support NULL data!!!");
498     }
499     ierr = pc_gamg->createdefaultdata( pc ); CHKERRQ(ierr);
500   }
501 
502   /* Get A_i and R_i */
503   if (pc_gamg->verbose) {
504     if(pc_gamg->verbose==1) ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);
505     else ierr =  MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);
506     CHKERRQ(ierr);
507     nnz0 = info.nz_used;
508     nnztot = info.nz_used;
509     PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
510                 mype,__FUNCT__,0,N,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
511                 (int)(nnz0/(PetscReal)N),npe);
512   }
513 
514   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
515         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
516         level++ ){
517     level1 = level + 1;
518 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
519     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
520 #endif
521 #if defined PETSC_GAMG_USE_LOG
522     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr);
523 #endif
524     { /* construct prolongator */
525       Mat Gmat;
526       PetscCoarsenData *agg_lists;
527       PetscInt fine_bs = bs;
528 
529       ierr = pc_gamg->graph( pc, Aarr[level], &Gmat ); CHKERRQ(ierr);
530       ierr = pc_gamg->coarsen( pc, &Gmat, &agg_lists ); CHKERRQ(ierr);
531       ierr = pc_gamg->prolongator( pc, Aarr[level], Gmat, agg_lists, &Parr[level1] );
532       CHKERRQ(ierr);
533 
534       if( Parr[level1] ){
535         /* get new block size of coarse matrices */
536         if( pc_gamg->col_bs_id != -1 ){
537           PetscBool flag;
538           ierr = PetscObjectComposedDataGetInt( (PetscObject)Parr[level1], pc_gamg->col_bs_id, bs, flag );
539           CHKERRQ( ierr ); assert(flag);
540         }
541       }
542 
543       if( pc_gamg->optprol && Parr[level1] ){
544         /* smooth */
545         ierr = pc_gamg->optprol( pc, Aarr[level], &Parr[level1] ); CHKERRQ(ierr);
546       }
547 
548       ierr = MatDestroy( &Gmat );      CHKERRQ(ierr);
549 
550       if ( pc_gamg->use_aggs_in_gasm ) {
551         ierr = PetscCDGetASMBlocks(agg_lists, fine_bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);  CHKERRQ(ierr);
552       }
553 
554       ierr = PetscCDDestroy( agg_lists );  CHKERRQ(ierr);
555     }
556 #if defined PETSC_GAMG_USE_LOG
557     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
558 #endif
559     /* cache eigen estimate */
560     if( pc_gamg->emax_id != -1 ){
561       PetscBool flag;
562       ierr = PetscObjectComposedDataGetReal( (PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag );
563       CHKERRQ( ierr );
564       if( !flag ) emaxs[level] = -1.;
565     }
566     else emaxs[level] = -1.;
567     if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
568     if( !Parr[level1] ) {
569       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);
570       break;
571     }
572 #if defined PETSC_GAMG_USE_LOG
573     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
574 #endif
575     ierr = createLevel( pc, Aarr[level], bs,
576                         (PetscBool)(level==pc_gamg->Nlevels-2),
577                         &Parr[level1], &nactivepe, &Aarr[level1] );
578     CHKERRQ(ierr);
579 #if defined PETSC_GAMG_USE_LOG
580     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
581 #endif
582     ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
583     if (pc_gamg->verbose){
584       PetscInt NN = M;
585       if(pc_gamg->verbose==1) {
586         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);
587         ierr = MatGetLocalSize( Aarr[level1], &NN, &N );CHKERRQ(ierr);
588       }
589       else ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info);
590       CHKERRQ(ierr);
591       nnztot += info.nz_used;
592       PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
593                   mype,__FUNCT__,(int)level1,N,pc_gamg->data_cell_cols,
594                   (int)(info.nz_used/(PetscReal)NN), nactivepe );
595       CHKERRQ(ierr);
596     }
597     /* stop if one node */
598     if( M/pc_gamg->data_cell_cols < 2 ) {
599       level++;
600       break;
601     }
602     /* Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; */
603     /* v = 1.e-10; /\* LU factor has hard wired numbers for small diags so this needs to match (yuk) *\/ */
604     /* ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); */
605     /* nloceq = Iend-Istart; */
606     /* ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr); */
607     /* ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr); */
608     /* ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr); */
609     /* for(kk=0;kk<nloceq;kk++){ */
610     /*   if(data_arr[kk]==0.0) { */
611     /*     id = kk + Istart; */
612     /*     ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); */
613     /*     CHKERRQ(ierr); */
614     /*     PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1); */
615     /*   } */
616     /* } */
617     /* ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); */
618     /* ierr = VecDestroy( &diag );                CHKERRQ(ierr); */
619     /* ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
620     /* ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
621 
622 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
623     ierr = PetscLogStagePop(); CHKERRQ( ierr );
624 #endif
625   } /* levels */
626 
627   if( pc_gamg->data ) {
628     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
629     pc_gamg->data = 0;
630   }
631 
632   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
633   pc_gamg->Nlevels = level + 1;
634   fine_level = level;
635   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
636 
637   /* simple setup */
638   if( !PETSC_TRUE ){
639     PC_MG_Levels **mglevels = mg->levels;
640     for (lidx=0,level=pc_gamg->Nlevels-1;
641          lidx<fine_level;
642          lidx++, level--){
643       ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr);
644       ierr = KSPSetOperators( mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );CHKERRQ(ierr);
645       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
646       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
647     }
648     ierr = KSPSetOperators( mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
649 
650     ierr = PCSetUp_MG( pc );  CHKERRQ( ierr );
651   }
652   else if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if one level */
653     /* set default smoothers & set operators */
654     for ( lidx = 1, level = pc_gamg->Nlevels-2;
655           lidx <= fine_level;
656           lidx++, level--) {
657       KSP smoother;
658       PC subpc;
659       /* set defaults */
660       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
661       ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
662       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
663 
664       /* override defaults and command line args (!) */
665       if ( pc_gamg->use_aggs_in_gasm ) {
666         if( PETSC_FALSE ){
667           PetscInt ii;
668           PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s Nblocks=%d\n",mype,__FUNCT__,nASMBlocksArr[level]);
669           for(ii=0;ii<nASMBlocksArr[level];ii++){
670             ISView(ASMLocalIDsArr[level][ii],PETSC_VIEWER_STDOUT_SELF);
671           }
672         }
673         PetscInt sz = nASMBlocksArr[level];
674         IS *is = ASMLocalIDsArr[level];
675         ierr = PCSetType( subpc, PCGASM ); CHKERRQ(ierr);
676         if(sz==0){
677           IS is;
678           PetscInt my0,kk;
679           ierr = MatGetOwnershipRange( Aarr[level], &my0, &kk ); CHKERRQ(ierr);
680           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is ); CHKERRQ(ierr);
681           ierr = PCGASMSetLocalSubdomains( subpc, 1, PETSC_NULL, &is ); CHKERRQ(ierr);
682         }
683         else {
684           ierr = PCGASMSetLocalSubdomains( subpc, sz, PETSC_NULL, is ); CHKERRQ(ierr);
685           ierr = PetscFree( is ); CHKERRQ(ierr);
686         }
687         ierr = PCGASMSetOverlap( subpc, 0 ); CHKERRQ(ierr);
688 
689         ASMLocalIDsArr[level] = PETSC_NULL;
690         nASMBlocksArr[level] = 0;
691         ierr = PCGASMSetType( subpc, PC_GASM_BASIC ); CHKERRQ(ierr);
692       }
693       else {
694         ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr);
695       }
696       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
697       /* set ops */
698       ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
699       ierr = PCMGSetInterpolation( pc, lidx, Parr[level+1] );CHKERRQ(ierr);
700     }
701     {
702       /* coarse grid */
703       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
704       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
705       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
706       ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
707       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
708       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
709       ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
710       ierr = PCSetUp( subpc ); CHKERRQ(ierr);
711       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
712       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
713       ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
714     }
715 
716     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
717     ierr = PetscObjectOptionsBegin( (PetscObject)pc );CHKERRQ(ierr);
718     ierr = PCSetFromOptions_MG( pc ); CHKERRQ(ierr);
719     ierr = PetscOptionsEnd();  CHKERRQ(ierr);
720     if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
721 
722     /* create cheby smoothers */
723     for ( lidx = 1, level = pc_gamg->Nlevels-2;
724           lidx <= fine_level;
725           lidx++, level--) {
726       KSP smoother;
727       PetscBool flag;
728       PC subpc;
729 
730       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
731 
732       /* do my own cheby */
733       ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &flag ); CHKERRQ(ierr);
734       if( flag ) {
735         PetscReal emax, emin;
736         ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
737         ierr = PetscTypeCompare( (PetscObject)subpc, PCJACOBI, &flag ); CHKERRQ(ierr);
738         if( flag && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
739         else{ /* eigen estimate 'emax' */
740           KSP eksp; Mat Lmat = Aarr[level];
741           Vec bb, xx;
742 
743           ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
744           ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
745           {
746             PetscRandom    rctx;
747             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
748             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
749             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
750             ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
751           }
752           ierr = KSPCreate( wcomm, &eksp );CHKERRQ(ierr);
753           ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
754           CHKERRQ(ierr);
755           ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
756 
757           ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
758           ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
759 
760           ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
761           ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
762           ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
763 
764           /* set PC type to be same as smoother, inc ref count, may not work */
765           ierr = PetscObjectReference( (PetscObject)subpc );   CHKERRQ(ierr);
766           ierr = KSPSetPC( eksp, subpc ); CHKERRQ( ierr );
767 
768           /* solve - keep stuff out of logging */
769           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
770           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
771           ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
772           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
773           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
774 
775           ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
776 
777           ierr = VecDestroy( &xx );       CHKERRQ(ierr);
778           ierr = VecDestroy( &bb );       CHKERRQ(ierr);
779 
780           ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
781 
782           if( pc_gamg->verbose > 0 ) {
783             PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d\n",__FUNCT__,emax,emin,level);
784           }
785         }
786         {
787           PetscInt N1, N0, tt;
788           ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
789           ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
790           /* heuristic - is this crap? */
791           emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
792           emax *= 1.05;
793         }
794         ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
795       } /* setup checby flag */
796     } /* non-coarse levels */
797 
798     /* clean up */
799     for(level=1;level<pc_gamg->Nlevels;level++){
800       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
801       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
802     }
803 
804     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
805 
806     {
807       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
808       ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
809       ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
810     }
811   }
812   else {
813     KSP smoother;
814     if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__);
815     ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
816     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
817     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
818     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
819   }
820 
821   PetscFunctionReturn(0);
822 }
823 
824 /* ------------------------------------------------------------------------- */
825 /*
826  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
827    that was created with PCCreate_GAMG().
828 
829    Input Parameter:
830 .  pc - the preconditioner context
831 
832    Application Interface Routine: PCDestroy()
833 */
834 #undef __FUNCT__
835 #define __FUNCT__ "PCDestroy_GAMG"
836 PetscErrorCode PCDestroy_GAMG( PC pc )
837 {
838   PetscErrorCode  ierr;
839   PC_MG           *mg = (PC_MG*)pc->data;
840   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
841 
842   PetscFunctionBegin;
843   ierr = PCReset_GAMG( pc );CHKERRQ(ierr);
844   ierr = PetscFree( pc_gamg );CHKERRQ(ierr);
845   ierr = PCDestroy_MG( pc );CHKERRQ(ierr);
846   PetscFunctionReturn(0);
847 }
848 
849 
850 #undef __FUNCT__
851 #define __FUNCT__ "PCGAMGSetProcEqLim"
852 /*@
853    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
854    processor reduction.
855 
856    Not Collective on PC
857 
858    Input Parameters:
859 .  pc - the preconditioner context
860 
861 
862    Options Database Key:
863 .  -pc_gamg_process_eq_limit
864 
865    Level: intermediate
866 
867    Concepts: Unstructured multrigrid preconditioner
868 
869 .seealso: ()
870 @*/
871 PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
872 {
873   PetscErrorCode ierr;
874 
875   PetscFunctionBegin;
876   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
877   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
878   PetscFunctionReturn(0);
879 }
880 
881 EXTERN_C_BEGIN
882 #undef __FUNCT__
883 #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
884 PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
885 {
886   PC_MG           *mg = (PC_MG*)pc->data;
887   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
888 
889   PetscFunctionBegin;
890   if(n>0) pc_gamg->min_eq_proc = n;
891   PetscFunctionReturn(0);
892 }
893 EXTERN_C_END
894 
895 #undef __FUNCT__
896 #define __FUNCT__ "PCGAMGSetCoarseEqLim"
897 /*@
898    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
899 
900  Collective on PC
901 
902    Input Parameters:
903 .  pc - the preconditioner context
904 
905 
906    Options Database Key:
907 .  -pc_gamg_coarse_eq_limit
908 
909    Level: intermediate
910 
911    Concepts: Unstructured multrigrid preconditioner
912 
913 .seealso: ()
914  @*/
915 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
916 {
917   PetscErrorCode ierr;
918 
919   PetscFunctionBegin;
920   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
921   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
922   PetscFunctionReturn(0);
923 }
924 
925 EXTERN_C_BEGIN
926 #undef __FUNCT__
927 #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
928 PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
929 {
930   PC_MG           *mg = (PC_MG*)pc->data;
931   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
932 
933   PetscFunctionBegin;
934   if(n>0) pc_gamg->coarse_eq_limit = n;
935   PetscFunctionReturn(0);
936 }
937 EXTERN_C_END
938 
939 #undef __FUNCT__
940 #define __FUNCT__ "PCGAMGSetRepartitioning"
941 /*@
942    PCGAMGSetRepartitioning - Repartition the coarse grids
943 
944    Collective on PC
945 
946    Input Parameters:
947 .  pc - the preconditioner context
948 
949 
950    Options Database Key:
951 .  -pc_gamg_repartition
952 
953    Level: intermediate
954 
955    Concepts: Unstructured multrigrid preconditioner
956 
957 .seealso: ()
958 @*/
959 PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
960 {
961   PetscErrorCode ierr;
962 
963   PetscFunctionBegin;
964   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
965   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
966   PetscFunctionReturn(0);
967 }
968 
969 EXTERN_C_BEGIN
970 #undef __FUNCT__
971 #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
972 PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
973 {
974   PC_MG           *mg = (PC_MG*)pc->data;
975   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
976 
977   PetscFunctionBegin;
978   pc_gamg->repart = n;
979   PetscFunctionReturn(0);
980 }
981 EXTERN_C_END
982 
983 #undef __FUNCT__
984 #define __FUNCT__ "PCGAMGSetUseASMAggs"
985 /*@
986    PCGAMGSetUseASMAggs -
987 
988    Collective on PC
989 
990    Input Parameters:
991 .  pc - the preconditioner context
992 
993 
994    Options Database Key:
995 .  -pc_gamg_use_agg_gasm
996 
997    Level: intermediate
998 
999    Concepts: Unstructured multrigrid preconditioner
1000 
1001 .seealso: ()
1002 @*/
1003 PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1004 {
1005   PetscErrorCode ierr;
1006 
1007   PetscFunctionBegin;
1008   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1009   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1010   PetscFunctionReturn(0);
1011 }
1012 
1013 EXTERN_C_BEGIN
1014 #undef __FUNCT__
1015 #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
1016 PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1017 {
1018   PC_MG           *mg = (PC_MG*)pc->data;
1019   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1020 
1021   PetscFunctionBegin;
1022   pc_gamg->use_aggs_in_gasm = n;
1023   PetscFunctionReturn(0);
1024 }
1025 EXTERN_C_END
1026 
1027 #undef __FUNCT__
1028 #define __FUNCT__ "PCGAMGSetNlevels"
1029 /*@
1030    PCGAMGSetNlevels -
1031 
1032    Not collective on PC
1033 
1034    Input Parameters:
1035 .  pc - the preconditioner context
1036 
1037 
1038    Options Database Key:
1039 .  -pc_mg_levels
1040 
1041    Level: intermediate
1042 
1043    Concepts: Unstructured multrigrid preconditioner
1044 
1045 .seealso: ()
1046 @*/
1047 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1048 {
1049   PetscErrorCode ierr;
1050 
1051   PetscFunctionBegin;
1052   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1053   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1054   PetscFunctionReturn(0);
1055 }
1056 
1057 EXTERN_C_BEGIN
1058 #undef __FUNCT__
1059 #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
1060 PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1061 {
1062   PC_MG           *mg = (PC_MG*)pc->data;
1063   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1064 
1065   PetscFunctionBegin;
1066   pc_gamg->Nlevels = n;
1067   PetscFunctionReturn(0);
1068 }
1069 EXTERN_C_END
1070 
1071 #undef __FUNCT__
1072 #define __FUNCT__ "PCGAMGSetThreshold"
1073 /*@
1074    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
1075 
1076    Not collective on PC
1077 
1078    Input Parameters:
1079 .  pc - the preconditioner context
1080 
1081 
1082    Options Database Key:
1083 .  -pc_gamg_threshold
1084 
1085    Level: intermediate
1086 
1087    Concepts: Unstructured multrigrid preconditioner
1088 
1089 .seealso: ()
1090 @*/
1091 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
1092 {
1093   PetscErrorCode ierr;
1094 
1095   PetscFunctionBegin;
1096   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1097   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
1098   PetscFunctionReturn(0);
1099 }
1100 
1101 EXTERN_C_BEGIN
1102 #undef __FUNCT__
1103 #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
1104 PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
1105 {
1106   PC_MG           *mg = (PC_MG*)pc->data;
1107   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1108 
1109   PetscFunctionBegin;
1110   pc_gamg->threshold = n;
1111   PetscFunctionReturn(0);
1112 }
1113 EXTERN_C_END
1114 
1115 #undef __FUNCT__
1116 #define __FUNCT__ "PCGAMGSetType"
1117 /*@
1118    PCGAMGSetType - Set solution method - calls sub create method
1119 
1120    Collective on PC
1121 
1122    Input Parameters:
1123 .  pc - the preconditioner context
1124 
1125 
1126    Options Database Key:
1127 .  -pc_gamg_type
1128 
1129    Level: intermediate
1130 
1131    Concepts: Unstructured multrigrid preconditioner
1132 
1133 .seealso: ()
1134 @*/
1135 PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type )
1136 {
1137   PetscErrorCode ierr;
1138 
1139   PetscFunctionBegin;
1140   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1141   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type));
1142   CHKERRQ(ierr);
1143   PetscFunctionReturn(0);
1144 }
1145 
1146 EXTERN_C_BEGIN
1147 #undef __FUNCT__
1148 #define __FUNCT__ "PCGAMGSetType_GAMG"
1149 PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type )
1150 {
1151   PetscErrorCode ierr,(*r)(PC);
1152 
1153   PetscFunctionBegin;
1154   ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r);
1155   CHKERRQ(ierr);
1156 
1157   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1158 
1159   /* call sub create method */
1160   ierr = (*r)(pc); CHKERRQ(ierr);
1161 
1162   PetscFunctionReturn(0);
1163 }
1164 EXTERN_C_END
1165 
1166 #undef __FUNCT__
1167 #define __FUNCT__ "PCSetFromOptions_GAMG"
1168 PetscErrorCode PCSetFromOptions_GAMG( PC pc )
1169 {
1170   PetscErrorCode  ierr;
1171   PC_MG           *mg = (PC_MG*)pc->data;
1172   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1173   PetscBool        flag;
1174   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
1175 
1176   PetscFunctionBegin;
1177   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1178   {
1179     /* -pc_gamg_verbose */
1180     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
1181                            "none", pc_gamg->verbose,
1182                            &pc_gamg->verbose, PETSC_NULL );
1183     CHKERRQ(ierr);
1184 
1185     /* -pc_gamg_repartition */
1186     ierr = PetscOptionsBool("-pc_gamg_repartition",
1187                             "Repartion coarse grids (false)",
1188                             "PCGAMGRepartitioning",
1189                             pc_gamg->repart,
1190                             &pc_gamg->repart,
1191                             &flag);
1192     CHKERRQ(ierr);
1193 
1194     /* -pc_gamg_use_agg_gasm */
1195     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1196                             "Use aggregation agragates for GASM smoother (false)",
1197                             "PCGAMGUseASMAggs",
1198                             pc_gamg->use_aggs_in_gasm,
1199                             &pc_gamg->use_aggs_in_gasm,
1200                             &flag);
1201     CHKERRQ(ierr);
1202 
1203     /* -pc_gamg_process_eq_limit */
1204     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1205                            "Limit (goal) on number of equations per process on coarse grids",
1206                            "PCGAMGSetProcEqLim",
1207                            pc_gamg->min_eq_proc,
1208                            &pc_gamg->min_eq_proc,
1209                            &flag );
1210     CHKERRQ(ierr);
1211 
1212     /* -pc_gamg_coarse_eq_limit */
1213     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1214                            "Limit on number of equations for the coarse grid",
1215                            "PCGAMGSetCoarseEqLim",
1216                            pc_gamg->coarse_eq_limit,
1217                            &pc_gamg->coarse_eq_limit,
1218                            &flag );
1219     CHKERRQ(ierr);
1220 
1221     /* -pc_gamg_threshold */
1222     ierr = PetscOptionsReal("-pc_gamg_threshold",
1223                             "Relative threshold to use for dropping edges in aggregation graph",
1224                             "PCGAMGSetThreshold",
1225                             pc_gamg->threshold,
1226                             &pc_gamg->threshold,
1227                             &flag );
1228     CHKERRQ(ierr);
1229     if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);
1230 
1231     ierr = PetscOptionsInt("-pc_mg_levels",
1232                            "Set number of MG levels",
1233                            "PCGAMGSetNlevels",
1234                            pc_gamg->Nlevels,
1235                            &pc_gamg->Nlevels,
1236                            &flag );
1237   }
1238   ierr = PetscOptionsTail();CHKERRQ(ierr);
1239 
1240   PetscFunctionReturn(0);
1241 }
1242 
1243 /* -------------------------------------------------------------------------- */
1244 /*MC
1245      PCCreate_GAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
1246        - This is the entry point to GAMG, registered in pcregis.c
1247 
1248    Options Database Keys:
1249    Multigrid options(inherited)
1250 +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1251 .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1252 .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1253 -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
1254 
1255   Level: intermediate
1256 
1257   Concepts: multigrid
1258 
1259 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1260            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1261            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1262            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1263            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1264 M*/
1265 EXTERN_C_BEGIN
1266 #undef __FUNCT__
1267 #define __FUNCT__ "PCCreate_GAMG"
1268 PetscErrorCode  PCCreate_GAMG( PC pc )
1269 {
1270   PetscErrorCode  ierr;
1271   PC_GAMG         *pc_gamg;
1272   PC_MG           *mg;
1273 
1274 #if defined PETSC_GAMG_USE_LOG
1275   static long count = 0;
1276 #endif
1277 
1278   PetscFunctionBegin;
1279   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1280   ierr = PCSetType( pc, PCMG );  CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
1281   ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr);
1282 
1283   /* create a supporting struct and attach it to pc */
1284   ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr);
1285   mg = (PC_MG*)pc->data;
1286   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
1287   mg->innerctx = pc_gamg;
1288 
1289   pc_gamg->setup_count = 0;
1290   /* these should be in subctx but repartitioning needs simple arrays */
1291   pc_gamg->data_sz = 0;
1292   pc_gamg->data = 0;
1293 
1294   /* register AMG type */
1295   if( !GAMGList ){
1296     ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
1297     ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
1298   }
1299 
1300   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1301   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1302   pc->ops->setup          = PCSetUp_GAMG;
1303   pc->ops->reset          = PCReset_GAMG;
1304   pc->ops->destroy        = PCDestroy_GAMG;
1305 
1306   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1307 					    "PCGAMGSetProcEqLim_C",
1308 					    "PCGAMGSetProcEqLim_GAMG",
1309 					    PCGAMGSetProcEqLim_GAMG);
1310   CHKERRQ(ierr);
1311 
1312   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1313 					    "PCGAMGSetCoarseEqLim_C",
1314 					    "PCGAMGSetCoarseEqLim_GAMG",
1315 					    PCGAMGSetCoarseEqLim_GAMG);
1316   CHKERRQ(ierr);
1317 
1318   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1319 					    "PCGAMGSetRepartitioning_C",
1320 					    "PCGAMGSetRepartitioning_GAMG",
1321 					    PCGAMGSetRepartitioning_GAMG);
1322   CHKERRQ(ierr);
1323 
1324   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1325 					    "PCGAMGSetUseASMAggs_C",
1326 					    "PCGAMGSetUseASMAggs_GAMG",
1327 					    PCGAMGSetUseASMAggs_GAMG);
1328   CHKERRQ(ierr);
1329 
1330   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1331 					    "PCGAMGSetThreshold_C",
1332 					    "PCGAMGSetThreshold_GAMG",
1333 					    PCGAMGSetThreshold_GAMG);
1334   CHKERRQ(ierr);
1335 
1336   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1337 					    "PCGAMGSetType_C",
1338 					    "PCGAMGSetType_GAMG",
1339 					    PCGAMGSetType_GAMG);
1340   CHKERRQ(ierr);
1341 
1342   pc_gamg->repart = PETSC_FALSE;
1343   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1344   pc_gamg->min_eq_proc = 100;
1345   pc_gamg->coarse_eq_limit = 800;
1346   pc_gamg->threshold = 0.05;
1347   pc_gamg->Nlevels = GAMG_MAXLEVELS;
1348   pc_gamg->verbose = 0;
1349   pc_gamg->emax_id = -1;
1350   pc_gamg->col_bs_id = -1;
1351 
1352   /* instantiate derived type */
1353   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1354   {
1355     char tname[256] = GAMGAGG;
1356     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",
1357                             GAMGList, tname, tname, sizeof(tname), PETSC_NULL );
1358     CHKERRQ(ierr);
1359     ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr);
1360   }
1361   ierr = PetscOptionsTail();CHKERRQ(ierr);
1362 
1363   /* private events */
1364 #if defined PETSC_GAMG_USE_LOG
1365   if( count++ == 0 ) {
1366     PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1367     PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1368     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1369     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1370     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1371     PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1372     PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1373     PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1374     PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1375     PetscLogEventRegister("  SA: colect data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1376     PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1377     PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1378     PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1379     PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1380     PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1381     PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1382     PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);
1383 
1384     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1385     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1386     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1387     /* create timer stages */
1388 #if defined GAMG_STAGES
1389     {
1390       char str[32];
1391       sprintf(str,"MG Level %d (finest)",0);
1392       PetscLogStageRegister(str, &gamg_stages[0]);
1393       PetscInt lidx;
1394       for (lidx=1;lidx<9;lidx++){
1395 	sprintf(str,"MG Level %d",lidx);
1396 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1397       }
1398     }
1399 #endif
1400   }
1401 #endif
1402   /* general events */
1403 #if defined PETSC_USE_LOG
1404   PetscLogEventRegister("PCGAMGgraph_AGG", PC_CLASSID, &PC_GAMGGgraph_AGG);
1405   PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);
1406   PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1407   PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1408   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1409   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1410   PetscLogEventRegister("PCGAMGProlOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);
1411 #endif
1412 
1413   PetscFunctionReturn(0);
1414 }
1415 EXTERN_C_END
1416