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