xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision b43b03e97385fe8ad94dfdfbf5044a665bcc78db)
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 int 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   PC_MG_Levels   **mglevels = mg->levels;
429   Mat              Amat = pc->mat, Pmat = pc->pmat;
430   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
431   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
432   PetscMPIInt      mype,npe,nactivepe;
433   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
434   PetscReal        emaxs[GAMG_MAXLEVELS];
435   MatInfo          info;
436 
437   PetscFunctionBegin;
438   pc_gamg->setup_count++;
439 
440   if( pc->setupcalled > 0 ) {
441     /* just do Galerkin grids */
442     Mat B,dA,dB;
443 
444     /* PCSetUp_MG seems to insists on setting this to GMRES */
445     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
446 
447     if( pc_gamg->Nlevels > 1 ) {
448       /* currently only handle case where mat and pmat are the same on coarser levels */
449       ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
450       /* (re)set to get dirty flag */
451       ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
452       ierr = KSPSetUp( mglevels[pc_gamg->Nlevels-1]->smoothd ); CHKERRQ(ierr);
453 
454       for (level=pc_gamg->Nlevels-2; level>-1; level--) {
455         ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
456         /* the first time through the matrix structure has changed from repartitioning */
457         if( pc_gamg->setup_count == 2 ) {
458           ierr = MatDestroy( &B );  CHKERRQ(ierr);
459           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
460           mglevels[level]->A = B;
461         }
462         else {
463           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
464         }
465         ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
466         dB = B;
467         /* setup KSP/PC */
468         ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
469       }
470     }
471 
472     pc->setupcalled = 2;
473 
474     PetscFunctionReturn(0);
475   }
476 
477   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
478   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
479   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
480   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
481   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
482   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
483   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
484 
485   if( pc_gamg->data == 0 && nloc > 0 ) {
486     if(!pc_gamg->createdefaultdata){
487       SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"'createdefaultdata' not set?!?! need to support NULL data!!!");
488     }
489     ierr = pc_gamg->createdefaultdata( pc ); CHKERRQ(ierr);
490   }
491 
492   /* Get A_i and R_i */
493   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
494   if (pc_gamg->verbose) {
495     PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
496                 mype,__FUNCT__,0,N,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
497                 (int)(info.nz_used/(PetscReal)N),npe);
498   }
499   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
500         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
501         level++ ){
502     level1 = level + 1;
503 #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
504     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
505 #endif
506 #if defined PETSC_USE_LOG
507     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr);
508 #endif
509     { /* construct prolongator */
510       Mat Gmat;   assert(pc_gamg->graph);
511       IS selected, llist;   assert(pc_gamg->coarsen);
512 
513       ierr = pc_gamg->graph( pc, Aarr[level], &Gmat ); CHKERRQ(ierr);
514       ierr = pc_gamg->coarsen( pc, Gmat, &selected, &llist ); CHKERRQ(ierr);
515       ierr = pc_gamg->prolongator( pc, Aarr[level], Gmat, selected, llist, &Parr[level1] );
516       CHKERRQ(ierr);
517 
518       if( Parr[level1] ){
519           /* get new block size of coarse matrices */
520         if( pc_gamg->col_bs_id != -1 ){
521           PetscBool flag;
522           ierr = PetscObjectComposedDataGetInt( (PetscObject)Parr[level1], pc_gamg->col_bs_id, bs, flag );
523           CHKERRQ( ierr ); assert(flag);
524         }
525       }
526 
527       if( pc_gamg->optprol && Parr[level1] ){
528         /* smooth */
529         ierr = pc_gamg->optprol( pc, Aarr[level], &Parr[level1] ); CHKERRQ(ierr);
530       }
531 
532       ierr = MatDestroy( &Gmat );  CHKERRQ(ierr);
533       ierr = ISDestroy( &selected );  CHKERRQ(ierr);
534       ierr = ISDestroy( &llist );  CHKERRQ(ierr);
535     }
536 #if defined PETSC_USE_LOG
537     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
538 #endif
539     /* cache eigen estimate */
540     if( pc_gamg->emax_id != -1 ){
541       PetscBool flag;
542       ierr = PetscObjectComposedDataGetReal( (PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag );
543       CHKERRQ( ierr );
544       if( !flag ) emaxs[level] = -1.;
545     }
546     else emaxs[level] = -1.;
547     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
548     if( !Parr[level1] ) {
549       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);
550       break;
551     }
552     else {
553 #if defined PETSC_USE_LOG
554       ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
555 #endif
556       ierr = createLevel( pc, Aarr[level], bs,
557                           (PetscBool)(level==pc_gamg->Nlevels-2),
558                           &Parr[level1], &nactivepe, &Aarr[level1] );
559       CHKERRQ(ierr);
560 #if defined PETSC_USE_LOG
561       ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
562 #endif
563       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
564       ierr = MatGetInfo( Aarr[level1], MAT_GLOBAL_SUM, &info ); CHKERRQ(ierr);
565       if (pc_gamg->verbose){
566         PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
567                     mype,__FUNCT__,(int)level1,N,pc_gamg->data_cell_cols,
568                     (int)(info.nz_used/(PetscReal)N),nactivepe);
569       }
570       /* stop if one node */
571       if( M/pc_gamg->data_cell_cols < 2 ) {
572         level++;
573         break;
574       }
575       /* Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; */
576       /* v = 1.e-10; /\* LU factor has hard wired numbers for small diags so this needs to match (yuk) *\/ */
577       /* ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); */
578       /* nloceq = Iend-Istart; */
579       /* ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr); */
580       /* ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr); */
581       /* ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr); */
582       /* for(kk=0;kk<nloceq;kk++){ */
583       /*   if(data_arr[kk]==0.0) { */
584       /*     id = kk + Istart;  */
585       /*     ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); */
586       /*     CHKERRQ(ierr); */
587       /*     PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1); */
588       /*   } */
589       /* } */
590       /* ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); */
591       /* ierr = VecDestroy( &diag );                CHKERRQ(ierr); */
592       /* ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
593       /* ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
594     }
595 
596 #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
597     ierr = PetscLogStagePop(); CHKERRQ( ierr );
598 #endif
599   } /* levels */
600 
601   if( pc_gamg->data ) {
602     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
603     pc_gamg->data = 0;
604   }
605 
606   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
607   pc_gamg->Nlevels = level + 1;
608   fine_level = level;
609   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
610 
611   if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if  */
612   /* set default smoothers */
613   for ( lidx = 1, level = pc_gamg->Nlevels-2;
614         lidx <= fine_level;
615         lidx++, level--) {
616     PetscReal emax, emin;
617     KSP smoother; PC subpc;
618     PetscBool isCheb;
619     /* set defaults */
620     ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
621     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
622     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
623     /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */
624     ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr);
625     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
626     /* overide defaults with input parameters */
627     ierr = KSPSetFromOptions( smoother ); CHKERRQ(ierr);
628 
629     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
630     /* do my own cheby */
631     ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr);
632 
633     if( isCheb ) {
634       ierr = PetscTypeCompare( (PetscObject)subpc, PCJACOBI, &isCheb ); CHKERRQ(ierr);
635       if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
636       else{ /* eigen estimate 'emax' */
637         KSP eksp; Mat Lmat = Aarr[level];
638         Vec bb, xx; PC pc;
639         const PCType type;
640 
641         ierr = PCGetType( subpc, &type );   CHKERRQ(ierr);
642         ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
643         ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
644         {
645           PetscRandom    rctx;
646           ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
647           ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
648           ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
649           ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
650         }
651         ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
652         ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
653         CHKERRQ(ierr);
654         ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
655 
656         ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
657         ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
658 
659         ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
660         ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
661         ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
662         ierr = PCSetType( pc, type ); CHKERRQ(ierr); /* should be same as eigen estimates op. */
663 
664         ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
665 
666 	/* solve - keep stuff out of logging */
667 	ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
668 	ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
669         ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
670 	ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
671 	ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
672 
673         ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
674 
675         ierr = VecDestroy( &xx );       CHKERRQ(ierr);
676         ierr = VecDestroy( &bb );       CHKERRQ(ierr);
677         ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
678 
679         if (pc_gamg->verbose) {
680           PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e\n",
681                       __FUNCT__,emax,emin);
682         }
683       }
684       {
685         PetscInt N1, N0, tt;
686         ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
687         ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
688         /* heuristic - is this crap? */
689         emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
690         emax *= 1.05;
691       }
692       ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
693     }
694   }
695   {
696     /* coarse grid */
697     KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
698     Mat Lmat = Aarr[pc_gamg->Nlevels-1];
699     ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
700     ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
701     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
702     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
703     ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
704     ierr = PCSetUp( subpc ); CHKERRQ(ierr);
705     ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
706     assert(ii==1);
707     ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
708     ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
709   }
710 
711   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
712   ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
713   ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr);
714   ierr = PetscOptionsEnd();CHKERRQ(ierr);
715   if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
716 
717   /* set interpolation between the levels, clean up */
718   for (lidx=0,level=pc_gamg->Nlevels-1;
719        lidx<fine_level;
720        lidx++, level--){
721     ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr);
722     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
723     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
724   }
725 
726   /* setupcalled is set to 0 so that MG is setup from scratch */
727   pc->setupcalled = 0;
728   ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
729   pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
730 
731   {
732     KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
733     ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
734     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
735     ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
736   }
737   }
738   else {
739     KSP smoother;
740     if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__);
741     ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
742     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
743     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
744     /* setupcalled is set to 0 so that MG is setup from scratch */
745     pc->setupcalled = 0;
746     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
747     pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
748   }
749   PetscFunctionReturn(0);
750 }
751 
752 /* ------------------------------------------------------------------------- */
753 /*
754    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
755    that was created with PCCreate_GAMG().
756 
757    Input Parameter:
758 .  pc - the preconditioner context
759 
760    Application Interface Routine: PCDestroy()
761 */
762 #undef __FUNCT__
763 #define __FUNCT__ "PCDestroy_GAMG"
764 PetscErrorCode PCDestroy_GAMG( PC pc )
765 {
766   PetscErrorCode  ierr;
767   PC_MG           *mg = (PC_MG*)pc->data;
768   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
769 
770   PetscFunctionBegin;
771   ierr = PCReset_GAMG( pc );CHKERRQ(ierr);
772   ierr = PetscFree( pc_gamg );CHKERRQ(ierr);
773   ierr = PCDestroy_MG( pc );CHKERRQ(ierr);
774   PetscFunctionReturn(0);
775 }
776 
777 
778 #undef __FUNCT__
779 #define __FUNCT__ "PCGAMGSetProcEqLim"
780 /*@
781    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
782    processor reduction.
783 
784    Not Collective on PC
785 
786    Input Parameters:
787 .  pc - the preconditioner context
788 
789 
790    Options Database Key:
791 .  -pc_gamg_process_eq_limit
792 
793    Level: intermediate
794 
795    Concepts: Unstructured multrigrid preconditioner
796 
797 .seealso: ()
798 @*/
799 PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
800 {
801   PetscErrorCode ierr;
802 
803   PetscFunctionBegin;
804   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
805   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
806   PetscFunctionReturn(0);
807 }
808 
809 EXTERN_C_BEGIN
810 #undef __FUNCT__
811 #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
812 PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
813 {
814   PC_MG           *mg = (PC_MG*)pc->data;
815   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
816 
817   PetscFunctionBegin;
818   if(n>0) pc_gamg->min_eq_proc = n;
819   PetscFunctionReturn(0);
820 }
821 EXTERN_C_END
822 
823 #undef __FUNCT__
824 #define __FUNCT__ "PCGAMGSetCoarseEqLim"
825 /*@
826    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
827 
828  Collective on PC
829 
830    Input Parameters:
831 .  pc - the preconditioner context
832 
833 
834    Options Database Key:
835 .  -pc_gamg_coarse_eq_limit
836 
837    Level: intermediate
838 
839    Concepts: Unstructured multrigrid preconditioner
840 
841 .seealso: ()
842  @*/
843 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
844 {
845   PetscErrorCode ierr;
846 
847   PetscFunctionBegin;
848   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
849   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
850   PetscFunctionReturn(0);
851 }
852 
853 EXTERN_C_BEGIN
854 #undef __FUNCT__
855 #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
856 PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
857 {
858   PC_MG           *mg = (PC_MG*)pc->data;
859   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
860 
861   PetscFunctionBegin;
862   if(n>0) pc_gamg->coarse_eq_limit = n;
863   PetscFunctionReturn(0);
864 }
865 EXTERN_C_END
866 
867 #undef __FUNCT__
868 #define __FUNCT__ "PCGAMGSetRepartitioning"
869 /*@
870    PCGAMGSetRepartitioning - Repartition the coarse grids
871 
872    Collective on PC
873 
874    Input Parameters:
875 .  pc - the preconditioner context
876 
877 
878    Options Database Key:
879 .  -pc_gamg_repartition
880 
881    Level: intermediate
882 
883    Concepts: Unstructured multrigrid preconditioner
884 
885 .seealso: ()
886 @*/
887 PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
888 {
889   PetscErrorCode ierr;
890 
891   PetscFunctionBegin;
892   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
893   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
894   PetscFunctionReturn(0);
895 }
896 
897 EXTERN_C_BEGIN
898 #undef __FUNCT__
899 #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
900 PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
901 {
902   PC_MG           *mg = (PC_MG*)pc->data;
903   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
904 
905   PetscFunctionBegin;
906   pc_gamg->repart = n;
907   PetscFunctionReturn(0);
908 }
909 EXTERN_C_END
910 
911 #undef __FUNCT__
912 #define __FUNCT__ "PCGAMGSetNlevels"
913 /*@
914    PCGAMGSetNlevels -
915 
916    Not collective on PC
917 
918    Input Parameters:
919 .  pc - the preconditioner context
920 
921 
922    Options Database Key:
923 .  -pc_mg_levels
924 
925    Level: intermediate
926 
927    Concepts: Unstructured multrigrid preconditioner
928 
929 .seealso: ()
930 @*/
931 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
932 {
933   PetscErrorCode ierr;
934 
935   PetscFunctionBegin;
936   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
937   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
938   PetscFunctionReturn(0);
939 }
940 
941 EXTERN_C_BEGIN
942 #undef __FUNCT__
943 #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
944 PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
945 {
946   PC_MG           *mg = (PC_MG*)pc->data;
947   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
948 
949   PetscFunctionBegin;
950   pc_gamg->Nlevels = n;
951   PetscFunctionReturn(0);
952 }
953 EXTERN_C_END
954 
955 #undef __FUNCT__
956 #define __FUNCT__ "PCGAMGSetThreshold"
957 /*@
958    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
959 
960    Not collective on PC
961 
962    Input Parameters:
963 .  pc - the preconditioner context
964 
965 
966    Options Database Key:
967 .  -pc_gamg_threshold
968 
969    Level: intermediate
970 
971    Concepts: Unstructured multrigrid preconditioner
972 
973 .seealso: ()
974 @*/
975 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
976 {
977   PetscErrorCode ierr;
978 
979   PetscFunctionBegin;
980   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
981   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
982   PetscFunctionReturn(0);
983 }
984 
985 EXTERN_C_BEGIN
986 #undef __FUNCT__
987 #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
988 PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
989 {
990   PC_MG           *mg = (PC_MG*)pc->data;
991   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
992 
993   PetscFunctionBegin;
994   pc_gamg->threshold = n;
995   PetscFunctionReturn(0);
996 }
997 EXTERN_C_END
998 
999 #undef __FUNCT__
1000 #define __FUNCT__ "PCGAMGSetType"
1001 /*@
1002    PCGAMGSetType - Set solution method - calls sub create method
1003 
1004    Collective on PC
1005 
1006    Input Parameters:
1007 .  pc - the preconditioner context
1008 
1009 
1010    Options Database Key:
1011 .  -pc_gamg_type
1012 
1013    Level: intermediate
1014 
1015    Concepts: Unstructured multrigrid preconditioner
1016 
1017 .seealso: ()
1018 @*/
1019 PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type )
1020 {
1021   PetscErrorCode ierr;
1022 
1023   PetscFunctionBegin;
1024   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1025   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type));
1026   CHKERRQ(ierr);
1027   PetscFunctionReturn(0);
1028 }
1029 
1030 EXTERN_C_BEGIN
1031 #undef __FUNCT__
1032 #define __FUNCT__ "PCGAMGSetType_GAMG"
1033 PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type )
1034 {
1035   PetscErrorCode ierr,(*r)(PC);
1036 
1037   PetscFunctionBegin;
1038   ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r);
1039   CHKERRQ(ierr);
1040 
1041   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1042 
1043   /* call sub create method */
1044   ierr = (*r)(pc); CHKERRQ(ierr);
1045 
1046   PetscFunctionReturn(0);
1047 }
1048 EXTERN_C_END
1049 
1050 #undef __FUNCT__
1051 #define __FUNCT__ "PCSetFromOptions_GAMG"
1052 PetscErrorCode PCSetFromOptions_GAMG( PC pc )
1053 {
1054   PetscErrorCode  ierr;
1055   PC_MG           *mg = (PC_MG*)pc->data;
1056   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1057   PetscBool        flag;
1058   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
1059 
1060   PetscFunctionBegin;
1061   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1062   {
1063     /* -pc_gamg_verbose */
1064     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
1065                            "none", pc_gamg->verbose,
1066                            &pc_gamg->verbose, PETSC_NULL );
1067     CHKERRQ(ierr);
1068 
1069     /* -pc_gamg_repartition */
1070     ierr = PetscOptionsBool("-pc_gamg_repartition",
1071                             "Repartion coarse grids (false)",
1072                             "PCGAMGRepartitioning",
1073                             pc_gamg->repart,
1074                             &pc_gamg->repart,
1075                             &flag);
1076     CHKERRQ(ierr);
1077 
1078     /* -pc_gamg_process_eq_limit */
1079     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1080                            "Limit (goal) on number of equations per process on coarse grids",
1081                            "PCGAMGSetProcEqLim",
1082                            pc_gamg->min_eq_proc,
1083                            &pc_gamg->min_eq_proc,
1084                            &flag );
1085     CHKERRQ(ierr);
1086 
1087     /* -pc_gamg_coarse_eq_limit */
1088     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1089                            "Limit on number of equations for the coarse grid",
1090                            "PCGAMGSetCoarseEqLim",
1091                            pc_gamg->coarse_eq_limit,
1092                            &pc_gamg->coarse_eq_limit,
1093                            &flag );
1094     CHKERRQ(ierr);
1095 
1096     /* -pc_gamg_threshold */
1097     ierr = PetscOptionsReal("-pc_gamg_threshold",
1098                             "Relative threshold to use for dropping edges in aggregation graph",
1099                             "PCGAMGSetThreshold",
1100                             pc_gamg->threshold,
1101                             &pc_gamg->threshold,
1102                             &flag );
1103     CHKERRQ(ierr);
1104     if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);
1105 
1106     ierr = PetscOptionsInt("-pc_mg_levels",
1107                            "Set number of MG levels",
1108                            "PCGAMGSetNlevels",
1109                            pc_gamg->Nlevels,
1110                            &pc_gamg->Nlevels,
1111                            &flag );
1112   }
1113   ierr = PetscOptionsTail();CHKERRQ(ierr);
1114 
1115   PetscFunctionReturn(0);
1116 }
1117 
1118 /* -------------------------------------------------------------------------- */
1119 /*MC
1120      PCCreate_GAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
1121        - This is the entry point to GAMG, registered in pcregis.c
1122 
1123    Options Database Keys:
1124    Multigrid options(inherited)
1125 +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1126 .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1127 .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1128 -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
1129 
1130   Level: intermediate
1131 
1132   Concepts: multigrid
1133 
1134 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1135            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1136            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1137            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1138            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1139 M*/
1140 EXTERN_C_BEGIN
1141 #undef __FUNCT__
1142 #define __FUNCT__ "PCCreate_GAMG"
1143 PetscErrorCode  PCCreate_GAMG( PC pc )
1144 {
1145   PetscErrorCode  ierr;
1146   PC_GAMG         *pc_gamg;
1147   PC_MG           *mg;
1148   PetscClassId     cookie;
1149 
1150 #if defined PETSC_USE_LOG
1151   static long count = 0;
1152 #endif
1153 
1154   PetscFunctionBegin;
1155   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1156   ierr = PCSetType( pc, PCMG );  CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
1157   ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr);
1158 
1159   /* create a supporting struct and attach it to pc */
1160   ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr);
1161   mg = (PC_MG*)pc->data;
1162   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
1163   mg->innerctx = pc_gamg;
1164 
1165   pc_gamg->setup_count = 0;
1166   /* these should be in subctx but repartitioning needs simple arrays */
1167   pc_gamg->data_sz = 0;
1168   pc_gamg->data = 0;
1169 
1170   /* register AMG type */
1171   if( !GAMGList ){
1172     ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
1173     ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
1174   }
1175 
1176   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1177   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1178   pc->ops->setup          = PCSetUp_GAMG;
1179   pc->ops->reset          = PCReset_GAMG;
1180   pc->ops->destroy        = PCDestroy_GAMG;
1181 
1182   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1183 					    "PCGAMGSetProcEqLim_C",
1184 					    "PCGAMGSetProcEqLim_GAMG",
1185 					    PCGAMGSetProcEqLim_GAMG);
1186   CHKERRQ(ierr);
1187 
1188   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1189 					    "PCGAMGSetCoarseEqLim_C",
1190 					    "PCGAMGSetCoarseEqLim_GAMG",
1191 					    PCGAMGSetCoarseEqLim_GAMG);
1192   CHKERRQ(ierr);
1193 
1194   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1195 					    "PCGAMGSetRepartitioning_C",
1196 					    "PCGAMGSetRepartitioning_GAMG",
1197 					    PCGAMGSetRepartitioning_GAMG);
1198   CHKERRQ(ierr);
1199 
1200   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1201 					    "PCGAMGSetThreshold_C",
1202 					    "PCGAMGSetThreshold_GAMG",
1203 					    PCGAMGSetThreshold_GAMG);
1204   CHKERRQ(ierr);
1205 
1206   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1207 					    "PCGAMGSetType_C",
1208 					    "PCGAMGSetType_GAMG",
1209 					    PCGAMGSetType_GAMG);
1210   CHKERRQ(ierr);
1211 
1212   pc_gamg->repart = PETSC_FALSE;
1213   pc_gamg->min_eq_proc = 100;
1214   pc_gamg->coarse_eq_limit = 800;
1215   pc_gamg->threshold = 0.05;
1216   pc_gamg->Nlevels = GAMG_MAXLEVELS;
1217   pc_gamg->verbose = 0;
1218   pc_gamg->emax_id = -1;
1219   pc_gamg->col_bs_id = -1;
1220 
1221   /* instantiate derived type */
1222   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1223   {
1224     char tname[256] = GAMGAGG;
1225     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",
1226                             GAMGList, tname, tname, sizeof(tname), PETSC_NULL );
1227     CHKERRQ(ierr);
1228     ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr);
1229   }
1230   ierr = PetscOptionsTail();CHKERRQ(ierr);
1231 
1232 #if defined PETSC_USE_LOG
1233   if( count++ == 0 ) {
1234     PetscClassIdRegister("GAMG Setup",&cookie);
1235     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
1236     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
1237     /* PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]); */
1238     /* PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]); */
1239     /* PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]); */
1240     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1241     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1242     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1243     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1244     PetscLogEventRegister("  SA: colect data", cookie, &gamg_setup_events[SET7]);
1245     PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]);
1246     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1247     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1248     PetscLogEventRegister("  repartition", cookie, &gamg_setup_events[SET12]);
1249     PetscLogEventRegister("  Invert-Sort", cookie, &gamg_setup_events[SET13]);
1250     PetscLogEventRegister("  Move A", cookie, &gamg_setup_events[SET14]);
1251     PetscLogEventRegister("  Move P", cookie, &gamg_setup_events[SET15]);
1252 
1253     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1254     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1255     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
1256     /* create timer stages */
1257 #if defined GAMG_STAGES
1258     {
1259       char str[32];
1260       sprintf(str,"MG Level %d (finest)",0);
1261       PetscLogStageRegister(str, &gamg_stages[0]);
1262       PetscInt lidx;
1263       for (lidx=1;lidx<9;lidx++){
1264 	sprintf(str,"MG Level %d",lidx);
1265 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1266       }
1267     }
1268 #endif
1269   }
1270 #endif
1271   PetscFunctionReturn(0);
1272 }
1273 EXTERN_C_END
1274