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