xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 26f2ff8fd9bc4ea607dadf72d7897bc7a997deee)
1 /*
2  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
3  */
4 #include <../src/ksp/pc/impls/gamg/gamg.h>
5 
6 PetscLogEvent gamg_setup_stages[NUM_SET];
7 
8 /* Private context for the GAMG preconditioner */
9 typedef struct gamg_TAG {
10   PetscInt       m_dim;
11   PetscInt       m_Nlevels;
12   PetscInt       m_data_sz;
13   PetscInt       m_data_rows;
14   PetscInt       m_data_cols;
15   PetscBool      m_useSA;
16   PetscReal     *m_data; /* blocked vector of vertex data on fine grid (coordinates) */
17 } PC_GAMG;
18 
19 /* -------------------------------------------------------------------------- */
20 /*
21    PCSetCoordinates_GAMG
22 
23    Input Parameter:
24    .  pc - the preconditioner context
25 */
26 EXTERN_C_BEGIN
27 #undef __FUNCT__
28 #define __FUNCT__ "PCSetCoordinates_GAMG"
29 PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords )
30 {
31   PC_MG          *mg = (PC_MG*)a_pc->data;
32   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
33   PetscErrorCode ierr;
34   PetscInt       arrsz,bs,my0,kk,ii,jj,nloc,Iend;
35   Mat            Amat = a_pc->pmat;
36   PetscBool      flag;
37   char           str[64];
38 
39   PetscFunctionBegin;
40   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
41   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
42   nloc = (Iend-my0)/bs;
43   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
44 
45   ierr  = PetscOptionsGetString(PETSC_NULL,"-pc_gamg_type",str,64,&flag);    CHKERRQ( ierr );
46   pc_gamg->m_useSA = (PetscBool)(flag && strcmp(str,"sa") == 0);
47 
48   pc_gamg->m_data_rows = 1;
49   if(a_coords == 0) pc_gamg->m_useSA = PETSC_TRUE; /* use SA if no data */
50   if( !pc_gamg->m_useSA ) pc_gamg->m_data_cols = a_ndm; /* coordinates */
51   else{ /* SA: null space vectors */
52     if(a_coords != 0 && bs==1 ) pc_gamg->m_data_cols = 1; /* scalar w/ coords and SA (not needed) */
53     else if(a_coords != 0) pc_gamg->m_data_cols = (a_ndm==2 ? 3 : 6); /* elasticity */
54     else pc_gamg->m_data_cols = bs; /* no data, force SA with constant null space vectors */
55     pc_gamg->m_data_rows = bs;
56   }
57   arrsz = nloc*pc_gamg->m_data_rows*pc_gamg->m_data_cols;
58 
59   /* create data - syntactic sugar that should be refactored at some point */
60   if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) {
61     ierr = PetscFree( pc_gamg->m_data );  CHKERRQ(ierr);
62     ierr = PetscMalloc((arrsz+1)*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr);
63   }
64   for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.;
65   pc_gamg->m_data[arrsz] = -99.;
66   /* copy data in - column oriented */
67   if( pc_gamg->m_useSA ) {
68     const PetscInt M = Iend - my0;
69     for(kk=0;kk<nloc;kk++){
70       PetscReal *data = &pc_gamg->m_data[kk*bs];
71       if( pc_gamg->m_data_cols==1 ) *data = 1.0;
72       else {
73         for(ii=0;ii<bs;ii++)
74 	  for(jj=0;jj<bs;jj++)
75 	    if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */
76 	    else data[ii*M + jj] = 0.0;
77         if( a_coords != 0 ) {
78           if( a_ndm == 2 ){ /* rotational modes */
79             data += 2*M;
80             data[0] = -a_coords[2*kk+1];
81             data[1] =  a_coords[2*kk];
82           }
83           else {
84             data += 3*M;
85             data[0] = 0.0;               data[M+0] =  a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1];
86             data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0;               data[2*M+1] =  a_coords[3*kk];
87             data[2] =  a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk];   data[2*M+2] = 0.0;
88           }
89         }
90       }
91     }
92   }
93   else {
94     for( kk = 0 ; kk < nloc ; kk++ ){
95       for( ii = 0 ; ii < a_ndm ; ii++ ) {
96         pc_gamg->m_data[ii*nloc + kk] =  a_coords[kk*a_ndm + ii];
97       }
98     }
99   }
100   assert(pc_gamg->m_data[arrsz] == -99.);
101 
102   pc_gamg->m_data_sz = arrsz;
103   pc_gamg->m_dim = a_ndm;
104 
105   PetscFunctionReturn(0);
106 }
107 EXTERN_C_END
108 
109 
110 /* -----------------------------------------------------------------------------*/
111 #undef __FUNCT__
112 #define __FUNCT__ "PCReset_GAMG"
113 PetscErrorCode PCReset_GAMG(PC pc)
114 {
115   PetscErrorCode  ierr;
116   PC_MG           *mg = (PC_MG*)pc->data;
117   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
118 
119   PetscFunctionBegin;
120   ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr);
121   pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0;
122   PetscFunctionReturn(0);
123 }
124 
125 /* -------------------------------------------------------------------------- */
126 /*
127    partitionLevel
128 
129    Input Parameter:
130    . a_Amat_fine - matrix on this fine (k) level
131    . a_ndata_rows - size of data to move (coarse grid)
132    . a_ndata_cols - size of data to move (coarse grid)
133    In/Output Parameter:
134    . a_P_inout - prolongation operator to the next level (k-1)
135    . a_coarse_data - data that need to be moved
136    . a_active_proc - number of active procs
137    Output Parameter:
138    . a_Amat_crs - coarse matrix that is created (k-1)
139 */
140 
141 #undef __FUNCT__
142 #define __FUNCT__ "partitionLevel"
143 PetscErrorCode partitionLevel( Mat a_Amat_fine,
144                                PetscInt a_ndata_rows,
145                                PetscInt a_ndata_cols,
146 			       PetscInt a_cbs,
147                                Mat *a_P_inout,
148                                PetscReal **a_coarse_data,
149                                PetscMPIInt *a_active_proc,
150                                Mat *a_Amat_crs
151                                )
152 {
153   PetscErrorCode   ierr;
154   Mat              Cmat,Pnew,Pold=*a_P_inout;
155   IS               new_indices,isnum;
156   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
157   PetscMPIInt      nactive_procs,mype,npe;
158   PetscInt         Istart,Iend,Istart0,Iend0,ncrs0,ncrs_new,fbs;
159   PetscInt         neq,NN;
160   PetscMPIInt      new_npe,targ_npe;
161   PetscBool        flag = PETSC_FALSE;
162 
163   PetscFunctionBegin;
164   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
165   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
166   ierr = MatGetBlockSize( a_Amat_fine, &fbs ); CHKERRQ(ierr);
167   /* RAP */
168   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
169 
170   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
171 
172   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
173   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
174 
175   /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
176   ierr = MatGetSize( Cmat, &neq, &NN );CHKERRQ(ierr);
177 #define MIN_EQ_PROC 2000
178   nactive_procs = *a_active_proc;
179   targ_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */
180 #define TOP_GRID_LIM 2000
181   if( targ_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1; /* chop coarsest grid */
182   else if (targ_npe >= nactive_procs ) new_npe = nactive_procs; /* no change */
183   else {
184     PetscMPIInt factstart,fact;
185     new_npe = -9999;
186     factstart = nactive_procs;
187     for(fact=factstart;fact>0;fact--){ /* try to find a better number of procs */
188       if( nactive_procs%fact==0 && neq/(nactive_procs/fact) > MIN_EQ_PROC ) {
189         new_npe = nactive_procs/fact;
190       }
191     }
192     if(new_npe == -9999) new_npe = nactive_procs;
193   }
194 
195   ierr  = PetscOptionsHasName(PETSC_NULL,"-pc_gamg_avoid_repartitioning",&flag);
196   CHKERRQ( ierr );
197   if( !flag ) { /* re-partition */
198     MatPartitioning  mpart;
199     Mat              adj;
200     const PetscInt  *is_idx;
201     PetscInt         is_sz,kk,jj,ii,old_fact=(npe/nactive_procs), *isnewproc_idx;
202     /* create sub communicator  */
203     MPI_Comm         cm,new_comm;
204     PetscInt         membershipKey = mype%old_fact, new_fact=(npe/new_npe),counts[npe];
205     IS               isnewproc;
206 
207     *a_active_proc = new_npe; /* output for next level */
208     ierr = MPI_Comm_split(wcomm, membershipKey, mype, &cm); CHKERRQ(ierr);
209     ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr);
210     ierr = MPI_Comm_free( &cm );                            CHKERRQ(ierr);
211 
212     /* MatPartitioningApply call MatConvert, which is collective */
213     ierr = PetscLogEventBegin(gamg_setup_stages[SET12],0,0,0,0);CHKERRQ(ierr);
214     if( a_cbs==1) {
215       ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
216     }
217     else {
218       /* make a scalar matrix to partition */
219       Mat tMat;
220       PetscInt Ii,ncols; const PetscScalar *vals; const PetscInt *idx;
221       MatInfo info;
222       ierr = MatGetInfo(Cmat,MAT_LOCAL,&info); CHKERRQ(ierr);
223       ncols = (PetscInt)info.nz_used/((ncrs0+1)*a_cbs*a_cbs)+1;
224 
225       ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
226                               PETSC_DETERMINE, PETSC_DETERMINE,
227                               2*ncols, PETSC_NULL, ncols, PETSC_NULL,
228                               &tMat );
229       CHKERRQ(ierr);
230 
231       for ( Ii = Istart0; Ii < Iend0; Ii++ ) {
232         PetscInt dest_row = Ii/a_cbs;
233         ierr = MatGetRow(Cmat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr);
234         for( jj = 0 ; jj < ncols ; jj++ ){
235           PetscInt dest_col = idx[jj]/a_cbs;
236           PetscScalar v = 1.0;
237           ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
238         }
239         ierr = MatRestoreRow(Cmat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr);
240       }
241       ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
242       ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
243 
244       ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
245 
246       ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
247     }
248     if( membershipKey == 0 ) {
249       if(ncrs0==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"zero local nodes -- increase min");
250       /* hack to fix global data that pmetis.c uses in 'adj' */
251       for( kk=0 , jj=0 ; kk<=npe ; jj++, kk += old_fact ) {
252         adj->rmap->range[jj] = adj->rmap->range[kk];
253       }
254       ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr);
255       ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
256       ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr);
257       ierr = MatPartitioningSetNParts( mpart, new_npe ); CHKERRQ(ierr);
258       ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr);
259       ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr);
260       /* collect IS info */
261       ierr = ISGetLocalSize( isnewproc, &is_sz );        CHKERRQ(ierr);
262       ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr);
263       ierr = ISGetIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
264       /* spread partitioning across machine - probably the right thing to do but machine spec. */
265       for(kk=0,jj=0;kk<is_sz;kk++){
266         for(ii=0 ; ii<a_cbs ; ii++, jj++ ) { /* expand for equation level by 'bs' */
267           isnewproc_idx[jj] = is_idx[kk] * new_fact;
268         }
269       }
270       ierr = ISRestoreIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
271       ierr = ISDestroy( &isnewproc );                    CHKERRQ(ierr);
272       is_sz *= a_cbs;
273     }
274     else {
275       isnewproc_idx = 0;
276       is_sz = 0;
277     }
278     ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
279     ierr = MPI_Comm_free( &new_comm );    CHKERRQ(ierr);
280     ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc );
281     if( membershipKey == 0 ) {
282       ierr = PetscFree( isnewproc_idx );  CHKERRQ(ierr);
283     }
284 
285     /*
286      Create an index set from the isnewproc index set to indicate the mapping TO
287      */
288     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
289     /*
290      Determine how many elements are assigned to each processor
291      */
292     ierr = ISPartitioningCount( isnewproc, npe, counts ); CHKERRQ(ierr);
293     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
294     ncrs_new = counts[mype]/a_cbs;
295     ierr = PetscLogEventEnd(gamg_setup_stages[SET12],0,0,0,0);   CHKERRQ(ierr);
296 
297     { /* Create a vector to contain the newly ordered element information */
298       const PetscInt *idx, data_sz=a_ndata_rows*a_ndata_cols;
299       const PetscInt  stride0=ncrs0*a_ndata_rows,strideNew=ncrs_new*a_ndata_rows;
300       PetscInt        ii,jj,kk;
301       IS              isscat;
302       PetscScalar    *array;
303       Vec             src_crd, dest_crd;
304       PetscReal      *data = *a_coarse_data;
305       VecScatter      vecscat;
306       PetscInt        tidx[ncrs0*data_sz];
307 
308       ierr = VecCreate( wcomm, &dest_crd );
309       ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
310       ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/
311       /*
312        There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
313        a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
314        */
315       ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
316       for(ii=0,jj=0; ii<ncrs0 ; ii++) {
317         PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
318         for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
319       }
320       ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
321       ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
322       CHKERRQ(ierr);
323       /*
324        Create a vector to contain the original vertex information for each element
325        */
326       ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
327       for( jj=0; jj<a_ndata_cols ; jj++ ) {
328         for( ii=0 ; ii<ncrs0 ; ii++) {
329           for( kk=0; kk<a_ndata_rows ; kk++ ) {
330             PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
331             ierr = VecSetValues( src_crd, 1, &jx, &data[ix], INSERT_VALUES );  CHKERRQ(ierr);
332           }
333         }
334       }
335       ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
336       ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
337       /*
338        Scatter the element vertex information (still in the original vertex ordering)
339        to the correct processor
340        */
341       ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
342       CHKERRQ(ierr);
343       ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
344       ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
345       ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
346       ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
347       ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
348       /*
349        Put the element vertex data into a new allocation of the gdata->ele
350        */
351       ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
352       ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
353 
354       ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
355       data = *a_coarse_data;
356       for( jj=0; jj<a_ndata_cols ; jj++ ) {
357         for( ii=0 ; ii<ncrs_new ; ii++) {
358           for( kk=0; kk<a_ndata_rows ; kk++ ) {
359             PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
360             data[ix] = PetscRealPart(array[jx]);
361             array[jx] = 1.e300;
362           }
363         }
364       }
365       ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
366       ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
367     }
368 
369     /*
370      Invert for MatGetSubMatrix
371      */
372     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
373     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
374     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
375     /* A_crs output */
376     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
377     CHKERRQ(ierr);
378 
379     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
380     Cmat = *a_Amat_crs; /* output */
381     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
382 
383     /* prolongator */
384     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
385     {
386       IS findices;
387       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
388       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
389       CHKERRQ(ierr);
390       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
391     }
392     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
393     *a_P_inout = Pnew; /* output */
394     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
395   }
396   else {
397     *a_Amat_crs = Cmat; /* output */
398   }
399 
400   PetscFunctionReturn(0);
401 }
402 
403 #define GAMG_MAXLEVELS 30
404 #if defined(PETSC_USE_LOG)
405 PetscLogStage  gamg_stages[20];
406 #endif
407 /* -------------------------------------------------------------------------- */
408 /*
409    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
410                     by setting data structures and options.
411 
412    Input Parameter:
413 .  pc - the preconditioner context
414 
415    Application Interface Routine: PCSetUp()
416 
417    Notes:
418    The interface routine PCSetUp() is not usually called directly by
419    the user, but instead is called by PCApply() if necessary.
420 */
421 #undef __FUNCT__
422 #define __FUNCT__ "PCSetUp_GAMG"
423 PetscErrorCode PCSetUp_GAMG( PC a_pc )
424 {
425   PetscErrorCode  ierr;
426   PC_MG           *mg = (PC_MG*)a_pc->data;
427   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
428   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
429   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
430   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
431   PetscMPIInt      mype,npe,nactivepe;
432   PetscBool        isOK;
433   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
434   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
435   MatInfo          info;
436 
437   PetscFunctionBegin;
438   if( a_pc->setupcalled ) {
439     /* no state data in GAMG to destroy */
440     ierr = PCReset_MG( a_pc ); CHKERRQ(ierr);
441   }
442   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
443   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
444   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
445   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
446   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
447   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
448   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
449 
450   /* get data of not around */
451   if( pc_gamg->m_data == 0 && nloc > 0 ) {
452     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
453   }
454   data = pc_gamg->m_data;
455 
456   /* Get A_i and R_i */
457   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
458   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",
459 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,(PetscInt)(info.nz_used/(PetscReal)N),npe);
460   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
461         level < GAMG_MAXLEVELS-1 && (level==0 || M>TOP_GRID_LIM) && (npe==1 || nactivepe>1);
462         level++ ){
463     level1 = level + 1;
464     ierr = PetscLogEventBegin(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr);
465     ierr = createProlongation(Aarr[level], data, pc_gamg->m_dim, pc_gamg->m_data_cols, pc_gamg->m_useSA,
466                               level, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level] );
467     CHKERRQ(ierr);
468     ierr = PetscFree( data ); CHKERRQ( ierr );
469     ierr = PetscLogEventEnd(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr);
470 
471     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
472 
473     if( isOK ) {
474       ierr = PetscLogEventBegin(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr);
475       ierr = partitionLevel( Aarr[level], pc_gamg->m_useSA ? bs : 1, pc_gamg->m_data_cols, bs,
476                              &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
477       CHKERRQ(ierr);
478       ierr = PetscLogEventEnd(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr);
479       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
480       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
481       PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, bs=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
482 		  mype,__FUNCT__,level1,N,bs,pc_gamg->m_data_cols,(PetscInt)(info.nz_used/(PetscReal)N),nactivepe);
483       /* coarse grids with SA can have zero row/cols from singleton aggregates */
484       /* aggregation method can probably gaurrentee this does not happen! - be safe for now */
485 
486       if( PETSC_TRUE ){
487         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
488         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
489         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
490         nloceq = Iend-Istart;
491         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
492         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
493         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
494         for(kk=0;kk<nloceq;kk++){
495           if(data_arr[kk]==0.0) {
496             id = kk + Istart;
497             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
498             CHKERRQ(ierr);
499             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added diag to zero (%d) on level %d \n",mype,__FUNCT__,id,level);
500           }
501         }
502         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
503         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
504         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
505         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
506       }
507     }
508     else{
509       coarse_data = 0;
510       break;
511     }
512     data = coarse_data;
513   }
514   if( coarse_data ) {
515     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
516   }
517   PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
518   pc_gamg->m_data = 0; /* destroyed coordinate data */
519   pc_gamg->m_Nlevels = level + 1;
520   fine_level = level;
521   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
522 
523   /* set default smoothers */
524   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
525         lidx <= fine_level;
526         lidx++, level--) {
527     PetscReal emax, emin;
528     KSP smoother; PC subpc;
529     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
530     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
531     if( emaxs[level] > 0.0 ) emax=emaxs[level];
532     else{ /* eigen estimate 'emax' */
533       KSP eksp; Mat Lmat = Aarr[level];
534       Vec bb, xx; PC pc;
535 
536       ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
537       ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
538       {
539 	PetscRandom    rctx;
540 	ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
541 	ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
542 	ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
543 	ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
544       }
545       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
546       ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr);
547       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
548       ierr = KSPSetOperators( eksp, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN ); CHKERRQ( ierr );
549       ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
550       ierr = PCSetType( pc, PCPBJACOBI ); CHKERRQ(ierr); /* should be same as above */
551       ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
552       CHKERRQ(ierr);
553       //ierr = KSPSetConvergenceTest( eksp, KSPSkipConverged, 0, 0 ); CHKERRQ(ierr);
554       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
555 
556       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
557       ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
558       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
559       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
560       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
561       ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
562       PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s max eigen=%e min=%e PC=%s\n",__FUNCT__,emax,emin,PETSC_GAMG_SMOOTHER);
563     }
564     {
565       PetscInt N1, N0, tt;
566       ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
567       ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
568       emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */
569       emax *= 1.05;
570 
571     }
572 
573     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], DIFFERENT_NONZERO_PATTERN );
574     ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
575     /*ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr);*/
576     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
577     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
578     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
579   }
580   {
581     KSP smoother; /* coarse grid */
582     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
583     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
584     ierr = KSPSetOperators( smoother, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN );
585     CHKERRQ(ierr);
586     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
587   }
588 
589   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
590   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
591   {
592     PetscBool galerkin;
593     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
594     if(galerkin){
595       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
596     }
597   }
598 
599   /* set interpolation between the levels, create timer stages, clean up */
600   if( PETSC_FALSE ) {
601     char str[32];
602     sprintf(str,"MG Level %d (%d)",0,pc_gamg->m_Nlevels-1);
603     PetscLogStageRegister(str, &gamg_stages[fine_level]);
604   }
605   for (lidx=0,level=pc_gamg->m_Nlevels-1;
606        lidx<fine_level;
607        lidx++, level--){
608     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
609     if( !PETSC_TRUE ) {
610       PetscViewer viewer; char fname[32];
611       sprintf(fname,"Pmat_%d.m",level);
612       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
613       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
614       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
615       ierr = PetscViewerDestroy( &viewer );
616       sprintf(fname,"Amat_%d.m",level);
617       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
618       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
619       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
620       ierr = PetscViewerDestroy( &viewer );
621     }
622     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
623     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
624     if( PETSC_FALSE ) {
625       char str[32];
626       sprintf(str,"MG Level %d (%d)",lidx+1,level-1);
627       PetscLogStageRegister(str, &gamg_stages[level-1]);
628     }
629   }
630 
631   /* setupcalled is set to 0 so that MG is setup from scratch */
632   a_pc->setupcalled = 0;
633   ierr = PCSetUp_MG(a_pc);CHKERRQ(ierr);
634 
635   PetscFunctionReturn(0);
636 }
637 
638 /* ------------------------------------------------------------------------- */
639 /*
640    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
641    that was created with PCCreate_GAMG().
642 
643    Input Parameter:
644 .  pc - the preconditioner context
645 
646    Application Interface Routine: PCDestroy()
647 */
648 #undef __FUNCT__
649 #define __FUNCT__ "PCDestroy_GAMG"
650 PetscErrorCode PCDestroy_GAMG(PC pc)
651 {
652   PetscErrorCode  ierr;
653   PC_MG           *mg = (PC_MG*)pc->data;
654   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
655 
656   PetscFunctionBegin;
657   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
658   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
659   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
660   PetscFunctionReturn(0);
661 }
662 
663 #undef __FUNCT__
664 #define __FUNCT__ "PCSetFromOptions_GAMG"
665 PetscErrorCode PCSetFromOptions_GAMG(PC pc)
666 {
667   /* PetscErrorCode  ierr; */
668   /* PC_MG           *mg = (PC_MG*)pc->data; */
669   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
670   /* MPI_Comm        comm = ((PetscObject)pc)->comm; */
671 
672   PetscFunctionBegin;
673   PetscFunctionReturn(0);
674 }
675 
676 /* -------------------------------------------------------------------------- */
677 /*
678  PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG
679 
680    Input Parameter:
681 .  pc - the preconditioner context
682 
683    Application Interface Routine: PCCreate()
684 
685   */
686  /* MC
687      PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide
688        fine grid discretization matrix and coordinates on the fine grid.
689 
690    Options Database Key:
691    Multigrid options(inherited)
692 +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
693 .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
694 .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
695    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
696    GAMG options:
697 
698    Level: intermediate
699   Concepts: multigrid
700 
701 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
702            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
703            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
704            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
705            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
706 M */
707 
708 EXTERN_C_BEGIN
709 #undef __FUNCT__
710 #define __FUNCT__ "PCCreate_GAMG"
711 PetscErrorCode  PCCreate_GAMG(PC pc)
712 {
713   PetscErrorCode  ierr;
714   PC_GAMG         *pc_gamg;
715   PC_MG           *mg;
716   PetscClassId     cookie;
717 
718   PetscFunctionBegin;
719   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
720   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
721   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
722 
723   /* create a supporting struct and attach it to pc */
724   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
725   mg = (PC_MG*)pc->data;
726   mg->innerctx = pc_gamg;
727 
728   pc_gamg->m_Nlevels    = -1;
729 
730   /* overwrite the pointers of PCMG by the functions of PCGAMG */
731   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
732   pc->ops->setup          = PCSetUp_GAMG;
733   pc->ops->reset          = PCReset_GAMG;
734   pc->ops->destroy        = PCDestroy_GAMG;
735 
736   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
737 					    "PCSetCoordinates_C",
738 					    "PCSetCoordinates_GAMG",
739 					    PCSetCoordinates_GAMG);CHKERRQ(ierr);
740   static int count = 0;
741   if( count++ == 0 ) {
742     PetscClassIdRegister("GAMG Setup",&cookie);
743     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_stages[SET1]);
744     PetscLogEventRegister(" make graph", cookie, &gamg_setup_stages[SET3]);
745     PetscLogEventRegister(" MIS/Agg", cookie, &gamg_setup_stages[SET4]);
746     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_stages[SET5]);
747     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_stages[SET6]);
748     PetscLogEventRegister("   search & set", cookie, &gamg_setup_stages[FIND_V]);
749     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_stages[SET7]);
750     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_stages[SET8]); */
751     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_stages[SET9]);
752     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_stages[SET2]);
753     PetscLogEventRegister(" PL repartition", cookie, &gamg_setup_stages[SET12]);
754     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_stages[SET13]); */
755     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_stages[SET10]); */
756     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_stages[SET11]); */
757   }
758 
759   PetscFunctionReturn(0);
760 }
761 EXTERN_C_END
762