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