xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision c20e42284aa1d345fa3c09863aa2c890811fb5bd)
15b89ad90SMark F. Adams /*
25b89ad90SMark F. Adams  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4313a3e24SSatish Balay #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
5313a3e24SSatish Balay #include "private/matimpl.h"
6f96513f1SMatthew G Knepley 
7b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
8b4fbaa2aSMark F. Adams PetscLogEvent gamg_setup_events[NUM_SET];
9b4fbaa2aSMark F. Adams #endif
10b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
11b4fbaa2aSMark F. Adams 
12b8fd24d8SMark F. Adams /*#define GAMG_STAGES*/
13b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
14b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
15b4fbaa2aSMark F. Adams #endif
16f96513f1SMatthew G Knepley 
175b89ad90SMark F. Adams /* Private context for the GAMG preconditioner */
185b89ad90SMark F. Adams typedef struct gamg_TAG {
195b89ad90SMark F. Adams   PetscInt       m_dim;
205b89ad90SMark F. Adams   PetscInt       m_Nlevels;
215b89ad90SMark F. Adams   PetscInt       m_data_sz;
22d3d6bff4SMark F. Adams   PetscInt       m_data_rows;
23d3d6bff4SMark F. Adams   PetscInt       m_data_cols;
2403a628feSMark F. Adams   PetscInt       m_count;
252c047e2dSMark F. Adams   PetscInt       m_method; /* 0: geomg; 1: plain agg 'pa'; 2: smoothed agg 'sa' */
265b89ad90SMark F. Adams   PetscReal     *m_data; /* blocked vector of vertex data on fine grid (coordinates) */
27676e1743SMark F. Adams   char           m_type[64];
28c20e4228SMark F. Adams   PetscBool      m_avoid_repart;
29c20e4228SMark F. Adams   PetscInt       m_min_eq_proc;
30c20e4228SMark F. Adams   PetscReal      m_threshold;
315b89ad90SMark F. Adams } PC_GAMG;
325b89ad90SMark F. Adams 
335b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
345b89ad90SMark F. Adams /*
355b89ad90SMark F. Adams    PCSetCoordinates_GAMG
365b89ad90SMark F. Adams 
375b89ad90SMark F. Adams    Input Parameter:
385b89ad90SMark F. Adams    .  pc - the preconditioner context
395b89ad90SMark F. Adams */
40a92563c5SMark F. Adams EXTERN_C_BEGIN
415b89ad90SMark F. Adams #undef __FUNCT__
425b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG"
43eb07cef2SMark F. Adams PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords )
445b89ad90SMark F. Adams {
45eb07cef2SMark F. Adams   PC_MG          *mg = (PC_MG*)a_pc->data;
465b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
476c237d78SBarry Smith   PetscErrorCode ierr;
48d3d6bff4SMark F. Adams   PetscInt       arrsz,bs,my0,kk,ii,jj,nloc,Iend;
49038e3b61SMark F. Adams   Mat            Amat = a_pc->pmat;
505b89ad90SMark F. Adams 
515b89ad90SMark F. Adams   PetscFunctionBegin;
5258471d46SMark F. Adams   PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 );
53038e3b61SMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
54d3d6bff4SMark F. Adams   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
55d3d6bff4SMark F. Adams   nloc = (Iend-my0)/bs;
56d3d6bff4SMark F. Adams   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
57038e3b61SMark F. Adams 
58d3d6bff4SMark F. Adams   pc_gamg->m_data_rows = 1;
592a44bfbaSMark F. Adams   if(a_coords==0 && pc_gamg->m_method==0) pc_gamg->m_method = 2; /* use SA if no coords */
60470e26f8SMark F. Adams   if( pc_gamg->m_method==0 ) pc_gamg->m_data_cols = a_ndm; /* coordinates */
61038e3b61SMark F. Adams   else{ /* SA: null space vectors */
62d3d6bff4SMark F. Adams     if(a_coords != 0 && bs==1 ) pc_gamg->m_data_cols = 1; /* scalar w/ coords and SA (not needed) */
63d3d6bff4SMark F. Adams     else if(a_coords != 0 ) pc_gamg->m_data_cols = (a_ndm==2 ? 3 : 6); /* elasticity */
64d3d6bff4SMark F. Adams     else pc_gamg->m_data_cols = bs; /* no data, force SA with constant null space vectors */
65d3d6bff4SMark F. Adams     pc_gamg->m_data_rows = bs;
66038e3b61SMark F. Adams   }
67d3d6bff4SMark F. Adams   arrsz = nloc*pc_gamg->m_data_rows*pc_gamg->m_data_cols;
685b89ad90SMark F. Adams 
69038e3b61SMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
706e3e101aSMark F. Adams   if (pc_gamg->m_data==0 || (pc_gamg->m_data_sz != arrsz)) {
715b89ad90SMark F. Adams     ierr = PetscFree( pc_gamg->m_data );  CHKERRQ(ierr);
7233a2b366SJed Brown     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->m_data ); CHKERRQ(ierr);
735b89ad90SMark F. Adams   }
74038e3b61SMark F. Adams   for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.;
75eb07cef2SMark F. Adams   pc_gamg->m_data[arrsz] = -99.;
76038e3b61SMark F. Adams   /* copy data in - column oriented */
77470e26f8SMark F. Adams   if( pc_gamg->m_method != 0 ) {
78d3d6bff4SMark F. Adams     const PetscInt M = Iend - my0;
79038e3b61SMark F. Adams     for(kk=0;kk<nloc;kk++){
80038e3b61SMark F. Adams       PetscReal *data = &pc_gamg->m_data[kk*bs];
81d3d6bff4SMark F. Adams       if( pc_gamg->m_data_cols==1 ) *data = 1.0;
82038e3b61SMark F. Adams       else {
83038e3b61SMark F. Adams         for(ii=0;ii<bs;ii++)
84038e3b61SMark F. Adams 	  for(jj=0;jj<bs;jj++)
85038e3b61SMark F. Adams 	    if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */
86038e3b61SMark F. Adams 	    else data[ii*M + jj] = 0.0;
87eb07cef2SMark F. Adams         if( a_coords != 0 ) {
88038e3b61SMark F. Adams           if( a_ndm == 2 ){ /* rotational modes */
89038e3b61SMark F. Adams             data += 2*M;
90038e3b61SMark F. Adams             data[0] = -a_coords[2*kk+1];
91038e3b61SMark F. Adams             data[1] =  a_coords[2*kk];
925b89ad90SMark F. Adams           }
93eb07cef2SMark F. Adams           else {
94038e3b61SMark F. Adams             data += 3*M;
95038e3b61SMark F. Adams             data[0] = 0.0;               data[M+0] =  a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1];
96038e3b61SMark F. Adams             data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0;               data[2*M+1] =  a_coords[3*kk];
97038e3b61SMark F. Adams             data[2] =  a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk];   data[2*M+2] = 0.0;
98038e3b61SMark F. Adams           }
99eb07cef2SMark F. Adams         }
100eb07cef2SMark F. Adams       }
101eb07cef2SMark F. Adams     }
102eb07cef2SMark F. Adams   }
103eb07cef2SMark F. Adams   else {
104038e3b61SMark F. Adams     for( kk = 0 ; kk < nloc ; kk++ ){
105038e3b61SMark F. Adams       for( ii = 0 ; ii < a_ndm ; ii++ ) {
106038e3b61SMark F. Adams         pc_gamg->m_data[ii*nloc + kk] =  a_coords[kk*a_ndm + ii];
107eb07cef2SMark F. Adams       }
108eb07cef2SMark F. Adams     }
109038e3b61SMark F. Adams   }
110eb07cef2SMark F. Adams   assert(pc_gamg->m_data[arrsz] == -99.);
111038e3b61SMark F. Adams 
1125b89ad90SMark F. Adams   pc_gamg->m_data_sz = arrsz;
113eb07cef2SMark F. Adams   pc_gamg->m_dim = a_ndm;
1145b89ad90SMark F. Adams 
1155b89ad90SMark F. Adams   PetscFunctionReturn(0);
1165b89ad90SMark F. Adams }
117a92563c5SMark F. Adams EXTERN_C_END
1185b89ad90SMark F. Adams 
119d3d6bff4SMark F. Adams 
120d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
121d3d6bff4SMark F. Adams #undef __FUNCT__
122d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
123d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
124d3d6bff4SMark F. Adams {
125d3d6bff4SMark F. Adams   PetscErrorCode  ierr;
126d3d6bff4SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
127d3d6bff4SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
128d3d6bff4SMark F. Adams 
129d3d6bff4SMark F. Adams   PetscFunctionBegin;
13058471d46SMark F. Adams   if( pc_gamg->m_data != 0 ) { /* this should not happen, cleaned up in SetUp */
131d3d6bff4SMark F. Adams     ierr = PetscFree(pc_gamg->m_data); CHKERRQ(ierr);
13258471d46SMark F. Adams   }
133d3d6bff4SMark F. Adams   pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0;
134d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
135d3d6bff4SMark F. Adams }
136d3d6bff4SMark F. Adams 
1375b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1385b89ad90SMark F. Adams /*
13911e60469SMark F. Adams    partitionLevel
1405b89ad90SMark F. Adams 
1415b89ad90SMark F. Adams    Input Parameter:
1423530afc2SMark F. Adams    . a_Amat_fine - matrix on this fine (k) level
143d3d6bff4SMark F. Adams    . a_ndata_rows - size of data to move (coarse grid)
144d3d6bff4SMark F. Adams    . a_ndata_cols - size of data to move (coarse grid)
1453530afc2SMark F. Adams    In/Output Parameter:
1463530afc2SMark F. Adams    . a_P_inout - prolongation operator to the next level (k-1)
147eb07cef2SMark F. Adams    . a_coarse_data - data that need to be moved
148afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
14911e60469SMark F. Adams    Output Parameter:
1503530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
151c20e4228SMark F. Adams    . a_avoid_repart -
152c20e4228SMark F. Adams    . a_min_eq_proc - goal for number of eqs. to have per proc. after proc reduction
1535b89ad90SMark F. Adams */
1545cb416c2SMark F. Adams 
1555b89ad90SMark F. Adams #undef __FUNCT__
15611e60469SMark F. Adams #define __FUNCT__ "partitionLevel"
1573530afc2SMark F. Adams PetscErrorCode partitionLevel( Mat a_Amat_fine,
158d3d6bff4SMark F. Adams                                PetscInt a_ndata_rows,
159d3d6bff4SMark F. Adams                                PetscInt a_ndata_cols,
160038e3b61SMark F. Adams 			       PetscInt a_cbs,
1613530afc2SMark F. Adams                                Mat *a_P_inout,
162eb07cef2SMark F. Adams                                PetscReal **a_coarse_data,
163afc97cdcSMark F. Adams                                PetscMPIInt *a_nactive_proc,
164c20e4228SMark F. Adams                                Mat *a_Amat_crs,
165c20e4228SMark F. Adams                                PetscBool a_avoid_repart,
166c20e4228SMark F. Adams                                PetscInt a_min_eq_proc
16711e60469SMark F. Adams                                )
1685b89ad90SMark F. Adams {
1695b89ad90SMark F. Adams   PetscErrorCode   ierr;
170038e3b61SMark F. Adams   Mat              Cmat,Pnew,Pold=*a_P_inout;
17111e60469SMark F. Adams   IS               new_indices,isnum;
1723530afc2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
173fa31c87dSMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive;
174fa31c87dSMark F. Adams   PetscInt         neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0;
1755b89ad90SMark F. Adams 
1765b89ad90SMark F. Adams   PetscFunctionBegin;
17711e60469SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
17811e60469SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
17911e60469SMark F. Adams   /* RAP */
18096434bc1SMark F. Adams #ifdef USE_R
18196434bc1SMark F. Adams   /* make R wih brute force for now */
18296434bc1SMark F. Adams   ierr = MatTranspose( Pold, Pnew );
18396434bc1SMark F. Adams   ierr = MatDestroy( &Pold );  CHKERRQ(ierr);
18496434bc1SMark F. Adams   ierr = MatRARt( a_Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
18596434bc1SMark F. Adams   Pold = Pnew;
18696434bc1SMark F. Adams #else
187038e3b61SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
18896434bc1SMark F. Adams #endif
189038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
190038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
191038e3b61SMark F. Adams   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
192038e3b61SMark F. Adams 
193c20e4228SMark F. Adams   if( a_avoid_repart) {
19422063be5SMark F. Adams     *a_Amat_crs = Cmat; /* output */
19522063be5SMark F. Adams   }
19622063be5SMark F. Adams   else {
19722063be5SMark F. Adams     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
1985ef31b24SMark F. Adams     Mat              adj;
19922063be5SMark F. Adams     const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols;
200b4fbaa2aSMark F. Adams     const PetscInt  stride0=ncrs0*a_ndata_rows;
20192a756f0SMark F. Adams     PetscInt        is_sz,*isnewproc_idx,ii,jj,kk,strideNew,*tidx;
202c9a0b8beSMark F. Adams     /* create sub communicator  */
20371d60426SBarry Smith     MPI_Comm        cm;
204afc97cdcSMark F. Adams     MPI_Group       wg, g2;
205fd3c6afaSMark F. Adams     PetscInt       *counts,inpe;
206fd3c6afaSMark F. Adams     PetscMPIInt    *ranks;
20722063be5SMark F. Adams     IS              isscat;
20822063be5SMark F. Adams     PetscScalar    *array;
20922063be5SMark F. Adams     Vec             src_crd, dest_crd;
21022063be5SMark F. Adams     PetscReal      *data = *a_coarse_data;
21122063be5SMark F. Adams     VecScatter      vecscat;
212b4fbaa2aSMark F. Adams     IS  isnewproc;
213e33ef3b1SMark F. Adams 
21422063be5SMark F. Adams     /* get number of PEs to make active, reduce */
21522063be5SMark F. Adams     ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
216c20e4228SMark F. Adams     new_npe = neq/a_min_eq_proc; /* hardwire min. number of eq/proc */
217c20e4228SMark F. Adams     if( new_npe == 0 || neq < 2*a_min_eq_proc ) new_npe = 1;
21822063be5SMark F. Adams     else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */
21922063be5SMark F. Adams 
22092a756f0SMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr);
221fd3c6afaSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
22292a756f0SMark F. Adams 
223fd3c6afaSMark F. Adams     ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr);
224afc97cdcSMark F. Adams     assert(counts[mype]==ncrs0);
225afc97cdcSMark F. Adams     /* count real active pes */
22622063be5SMark F. Adams     for( nactive = jj = 0 ; jj < npe ; jj++) {
227afc97cdcSMark F. Adams       if( counts[jj] != 0 ) {
22822063be5SMark F. Adams 	ranks[nactive++] = jj;
229afc97cdcSMark F. Adams         }
230afc97cdcSMark F. Adams     }
2313303bcf9Sadams 
2323303bcf9Sadams     if (nactive < new_npe) new_npe = nactive; /* this can happen with empty input procs */
23322063be5SMark F. Adams 
23458471d46SMark F. Adams #ifdef VERBOSE
23522063be5SMark F. Adams     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);
23658471d46SMark F. Adams #endif
23758471d46SMark F. Adams 
23822063be5SMark F. Adams     *a_nactive_proc = new_npe; /* output */
2392f03bc48SMark F. Adams 
240afc97cdcSMark F. Adams     ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr);
24122063be5SMark F. Adams     ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr);
242afc97cdcSMark F. Adams     ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr);
243b4fbaa2aSMark F. Adams 
244afc97cdcSMark F. Adams     ierr = MPI_Group_free( &wg );                            CHKERRQ(ierr);
245afc97cdcSMark F. Adams     ierr = MPI_Group_free( &g2 );                            CHKERRQ(ierr);
246c9a0b8beSMark F. Adams 
2475ef31b24SMark F. Adams     /* MatPartitioningApply call MatConvert, which is collective */
248b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
249b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
250b4fbaa2aSMark F. Adams #endif
251038e3b61SMark F. Adams     if( a_cbs == 1) {
252038e3b61SMark F. Adams       ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
253eb07cef2SMark F. Adams     }
254eb07cef2SMark F. Adams     else{
255038e3b61SMark F. Adams       /* make a scalar matrix to partition */
256eb07cef2SMark F. Adams       Mat tMat;
25758471d46SMark F. Adams       PetscInt ncols,jj,Ii;
258b4fbaa2aSMark F. Adams       const PetscScalar *vals;
259b4fbaa2aSMark F. Adams       const PetscInt *idx;
26058471d46SMark F. Adams       PetscInt *d_nnz;
2615ee9c036SSatish Balay       static int llev = 0;
262b4fbaa2aSMark F. Adams 
26358471d46SMark F. Adams       ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
26458471d46SMark F. Adams       for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) {
26558471d46SMark F. Adams         ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
26658471d46SMark F. Adams         d_nnz[jj] = ncols/a_cbs;
26758471d46SMark F. Adams         if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
26858471d46SMark F. Adams         ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
26958471d46SMark F. Adams       }
2706876a03eSMark F. Adams 
271eb07cef2SMark F. Adams       ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
272eb07cef2SMark F. Adams                               PETSC_DETERMINE, PETSC_DETERMINE,
27358471d46SMark F. Adams                               0, d_nnz, 0, d_nnz,
274eb07cef2SMark F. Adams                               &tMat );
2756876a03eSMark F. Adams       CHKERRQ(ierr);
27658471d46SMark F. Adams       ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
277eb07cef2SMark F. Adams 
27822063be5SMark F. Adams       for ( ii = Istart0; ii < Iend0; ii++ ) {
27922063be5SMark F. Adams         PetscInt dest_row = ii/a_cbs;
28022063be5SMark F. Adams         ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
281eb07cef2SMark F. Adams         for( jj = 0 ; jj < ncols ; jj++ ){
282038e3b61SMark F. Adams           PetscInt dest_col = idx[jj]/a_cbs;
283eb07cef2SMark F. Adams           PetscScalar v = 1.0;
284eb07cef2SMark F. Adams           ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
285eb07cef2SMark F. Adams         }
28622063be5SMark F. Adams         ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
287eb07cef2SMark F. Adams       }
288eb07cef2SMark F. Adams       ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
289eb07cef2SMark F. Adams       ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
290eb07cef2SMark F. Adams 
291b4fbaa2aSMark F. Adams       if( llev++ == -1 ) {
292b4fbaa2aSMark F. Adams 	PetscViewer viewer; char fname[32];
293b4fbaa2aSMark F. Adams 	sprintf(fname,"part_mat_%d.mat",llev);
294b4fbaa2aSMark F. Adams 	PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
295b4fbaa2aSMark F. Adams 	ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
296b4fbaa2aSMark F. Adams 	ierr = PetscViewerDestroy( &viewer );
297b4fbaa2aSMark F. Adams       }
298b4fbaa2aSMark F. Adams 
299eb07cef2SMark F. Adams       ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
300eb07cef2SMark F. Adams 
301eb07cef2SMark F. Adams       ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
302eb07cef2SMark F. Adams     }
303afc97cdcSMark F. Adams     if( ncrs0 != 0 ){
304b4fbaa2aSMark F. Adams       const PetscInt *is_idx;
305b4fbaa2aSMark F. Adams       MatPartitioning  mpart;
3065ef31b24SMark F. Adams       /* hack to fix global data that pmetis.c uses in 'adj' */
30722063be5SMark F. Adams       for( nactive = jj = 0 ; jj < npe ; jj++ ) {
308afc97cdcSMark F. Adams 	if( counts[jj] != 0 ) {
30922063be5SMark F. Adams 	  adj->rmap->range[nactive++] = adj->rmap->range[jj];
310afc97cdcSMark F. Adams 	}
3115ef31b24SMark F. Adams       }
31222063be5SMark F. Adams       adj->rmap->range[nactive] = adj->rmap->range[npe];
3132f03bc48SMark F. Adams 
31471d60426SBarry Smith       ierr = MatPartitioningCreate( cm, &mpart ); CHKERRQ(ierr);
3155ef31b24SMark F. Adams       ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
31611e60469SMark F. Adams       ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
3175ef31b24SMark F. Adams       ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
31811e60469SMark F. Adams       ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr);
31911e60469SMark F. Adams       ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
32022063be5SMark F. Adams 
3215ef31b24SMark F. Adams       /* collect IS info */
3225ef31b24SMark F. Adams       ierr = ISGetLocalSize( isnewproc, &is_sz );       CHKERRQ(ierr);
323038e3b61SMark F. Adams       ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr);
3245ef31b24SMark F. Adams       ierr = ISGetIndices( isnewproc, &is_idx );        CHKERRQ(ierr);
325b4fbaa2aSMark F. Adams       /* spread partitioning across machine - best way ??? */
326ecd57171SMark F. Adams       NN = 1; /*npe/new_npe;*/
327b4fbaa2aSMark F. Adams       for( kk = jj = 0 ; kk < is_sz ; kk++ ){
328afc97cdcSMark F. Adams         for( ii = 0 ; ii < a_cbs ; ii++, jj++ ) {
32922063be5SMark F. Adams           isnewproc_idx[jj] = is_idx[kk] * NN; /* distribution */
330eb07cef2SMark F. Adams         }
3315ef31b24SMark F. Adams       }
3325ef31b24SMark F. Adams       ierr = ISRestoreIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
3335ef31b24SMark F. Adams       ierr = ISDestroy( &isnewproc );                    CHKERRQ(ierr);
33471d60426SBarry Smith       ierr = MPI_Comm_free( &cm );              CHKERRQ(ierr);
335b4fbaa2aSMark F. Adams 
336b4fbaa2aSMark F. Adams       is_sz *= a_cbs;
3375ef31b24SMark F. Adams     }
3385ef31b24SMark F. Adams     else{
3395ef31b24SMark F. Adams       isnewproc_idx = 0;
3405ef31b24SMark F. Adams       is_sz = 0;
3415ef31b24SMark F. Adams     }
3425ef31b24SMark F. Adams     ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
3435ef31b24SMark F. Adams     ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc );
344afc97cdcSMark F. Adams     if( isnewproc_idx != 0 ) {
3455ef31b24SMark F. Adams       ierr = PetscFree( isnewproc_idx );  CHKERRQ(ierr);
3465ef31b24SMark F. Adams     }
347e33ef3b1SMark F. Adams 
34811e60469SMark F. Adams     /*
34911e60469SMark F. Adams      Create an index set from the isnewproc index set to indicate the mapping TO
35011e60469SMark F. Adams      */
35111e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
35211e60469SMark F. Adams     /*
35311e60469SMark F. Adams      Determine how many elements are assigned to each processor
35411e60469SMark F. Adams      */
355fd3c6afaSMark F. Adams     inpe = npe;
356fd3c6afaSMark F. Adams     ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr);
35711e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
358038e3b61SMark F. Adams     ncrs_new = counts[mype]/a_cbs;
35922063be5SMark F. Adams     strideNew = ncrs_new*a_ndata_rows;
360b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
361b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
362b4fbaa2aSMark F. Adams #endif
36322063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
36411e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
365d3d6bff4SMark F. Adams     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
366470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
36711e60469SMark F. Adams     /*
368d3d6bff4SMark F. Adams      There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
369d3d6bff4SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
37011e60469SMark F. Adams      */
37192a756f0SMark F. Adams     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
37211e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
373038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
374d3d6bff4SMark F. Adams       PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
375d3d6bff4SMark F. Adams       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
37611e60469SMark F. Adams     }
377038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
378d3d6bff4SMark F. Adams     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
37911e60469SMark F. Adams     CHKERRQ(ierr);
38092a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
38111e60469SMark F. Adams     /*
38211e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
38311e60469SMark F. Adams      */
384d3d6bff4SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
385d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
386d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs0 ; ii++) {
387d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
388d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
389676e1743SMark F. Adams           PetscScalar tt = (PetscScalar)data[ix];
390676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
391d3d6bff4SMark F. Adams 	}
392038e3b61SMark F. Adams       }
393eb07cef2SMark F. Adams     }
394eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
395eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
39611e60469SMark F. Adams     /*
39711e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
39811e60469SMark F. Adams       to the correct processor
39911e60469SMark F. Adams     */
40011e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
40111e60469SMark F. Adams     CHKERRQ(ierr);
40211e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
40311e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40411e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40511e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
40611e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
40711e60469SMark F. Adams     /*
40811e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
40911e60469SMark F. Adams     */
410eb07cef2SMark F. Adams     ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
411d3d6bff4SMark F. Adams     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
412eb07cef2SMark F. Adams 
41311e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
414eb07cef2SMark F. Adams     data = *a_coarse_data;
415d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
416d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs_new ; ii++) {
417d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
418d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
419d3d6bff4SMark F. Adams 	  data[ix] = PetscRealPart(array[jx]);
420d3d6bff4SMark F. Adams 	  array[jx] = 1.e300;
421d3d6bff4SMark F. Adams 	}
422038e3b61SMark F. Adams       }
423038e3b61SMark F. Adams     }
42411e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
42511e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
42611e60469SMark F. Adams     /*
42711e60469SMark F. Adams       Invert for MatGetSubMatrix
42811e60469SMark F. Adams     */
429038e3b61SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
43011e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
43111e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
43211e60469SMark F. Adams     /* A_crs output */
433038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
43411e60469SMark F. Adams     CHKERRQ(ierr);
435eb07cef2SMark F. Adams 
436038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
437e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
438038e3b61SMark F. Adams     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
439eb07cef2SMark F. Adams 
44011e60469SMark F. Adams     /* prolongator */
44111e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
44211e60469SMark F. Adams     {
44311e60469SMark F. Adams       IS findices;
44411e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
44596434bc1SMark F. Adams #ifdef USE_R
44696434bc1SMark F. Adams       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
44796434bc1SMark F. Adams #else
44811e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
44996434bc1SMark F. Adams #endif
45011e60469SMark F. Adams       CHKERRQ(ierr);
451d61a3a7dSMark F. Adams 
45211e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
45311e60469SMark F. Adams     }
454d61a3a7dSMark F. Adams 
4553530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
4563530afc2SMark F. Adams     *a_P_inout = Pnew; /* output */
45792a756f0SMark F. Adams 
45811e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
45992a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
46092a756f0SMark F. Adams     ierr = PetscFree( ranks );  CHKERRQ(ierr);
461e33ef3b1SMark F. Adams   }
4625b89ad90SMark F. Adams 
4635b89ad90SMark F. Adams   PetscFunctionReturn(0);
4645b89ad90SMark F. Adams }
4655b89ad90SMark F. Adams 
4665b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4675b89ad90SMark F. Adams /*
4685b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4695b89ad90SMark F. Adams                     by setting data structures and options.
4705b89ad90SMark F. Adams 
4715b89ad90SMark F. Adams    Input Parameter:
4725b89ad90SMark F. Adams .  pc - the preconditioner context
4735b89ad90SMark F. Adams 
4745b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4755b89ad90SMark F. Adams 
4765b89ad90SMark F. Adams    Notes:
4775b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4785b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4795b89ad90SMark F. Adams */
4805b89ad90SMark F. Adams #undef __FUNCT__
4815b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
482eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
4835b89ad90SMark F. Adams {
4845b89ad90SMark F. Adams   PetscErrorCode  ierr;
485eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
4865b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
48758471d46SMark F. Adams   PC_MG_Levels   **mglevels = mg->levels;
488eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
489d3d6bff4SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
490eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
4913530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
492038e3b61SMark F. Adams   PetscBool        isOK;
493587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
494587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
495737a81a9SMark F. Adams   MatInfo          info;
4965ef31b24SMark F. Adams 
4975b89ad90SMark F. Adams   PetscFunctionBegin;
49803a628feSMark F. Adams   pc_gamg->m_count++;
49958471d46SMark F. Adams   if( a_pc->setupcalled > 0 ) {
50003a628feSMark F. Adams     /* just do Galerkin grids */
50158471d46SMark F. Adams     Mat B,dA,dB;
50258471d46SMark F. Adams 
50358471d46SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
50458471d46SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
50558471d46SMark F. Adams 
50658471d46SMark F. Adams     /* currently only handle case where mat and pmat are the same on coarser levels */
50758471d46SMark F. Adams     ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
50858471d46SMark F. Adams     /* (re)set to get dirty flag */
50958471d46SMark F. Adams     ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
51058471d46SMark F. Adams     ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr);
51158471d46SMark F. Adams 
51258471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-2; level>-1; level--) {
51358471d46SMark F. Adams       ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
51403a628feSMark F. Adams       /* the first time through the matrix structure has changed from repartitioning */
51503a628feSMark F. Adams       if( pc_gamg->m_count == 2 ) {
51603a628feSMark F. Adams         ierr = MatDestroy( &B );  CHKERRQ(ierr);
51703a628feSMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
51803a628feSMark F. Adams         mglevels[level]->A = B;
51903a628feSMark F. Adams       }
52003a628feSMark F. Adams       else {
52158471d46SMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
52203a628feSMark F. Adams       }
52358471d46SMark F. Adams       ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
52458471d46SMark F. Adams       dB = B;
52558471d46SMark F. Adams       /* setup KSP/PC */
52658471d46SMark F. Adams       ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
52758471d46SMark F. Adams     }
52858471d46SMark F. Adams 
52958471d46SMark F. Adams #define PRINT_MATS !PETSC_TRUE
53058471d46SMark F. Adams     /* plot levels - A */
53158471d46SMark F. Adams     if( PRINT_MATS ) {
53258471d46SMark F. Adams       for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){
53358471d46SMark F. Adams         PetscViewer viewer;
53458471d46SMark F. Adams         char fname[32]; KSP smoother; Mat Tmat, TTm;
53558471d46SMark F. Adams         ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
53658471d46SMark F. Adams         ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr);
53703a628feSMark F. Adams         sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
53858471d46SMark F. Adams         ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
53958471d46SMark F. Adams         ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
54058471d46SMark F. Adams         ierr = MatView( Tmat, viewer ); CHKERRQ(ierr);
54158471d46SMark F. Adams         ierr = PetscViewerDestroy( &viewer );
54258471d46SMark F. Adams       }
54358471d46SMark F. Adams     }
54458471d46SMark F. Adams 
54503a628feSMark F. Adams     a_pc->setupcalled = 2;
54603a628feSMark F. Adams 
54758471d46SMark F. Adams     PetscFunctionReturn(0);
548eb07cef2SMark F. Adams   }
549baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
550baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
5515b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
5525b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
5533530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
554eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
555eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
556eb07cef2SMark F. Adams 
557038e3b61SMark F. Adams   /* get data of not around */
558038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
559038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
560eb07cef2SMark F. Adams   }
561eb07cef2SMark F. Adams   data = pc_gamg->m_data;
562038e3b61SMark F. Adams 
563eb07cef2SMark F. Adams   /* Get A_i and R_i */
564737a81a9SMark F. Adams   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
56558471d46SMark F. Adams #ifdef VERBOSE
5662c047e2dSMark F. Adams   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",
5672c047e2dSMark F. Adams 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,
5682c047e2dSMark F. Adams 	      (int)(info.nz_used/(PetscReal)N),npe);
56958471d46SMark F. Adams #endif
5708f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
571c20e4228SMark F. Adams         level < (GAMG_MAXLEVELS-1) && (level==0 || M>2*pc_gamg->m_min_eq_proc); /* && (npe==1 || nactivepe>1); */
5720205a208SMark F. Adams         level++ ){
5735b89ad90SMark F. Adams     level1 = level + 1;
574b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
575b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
576b4fbaa2aSMark F. Adams #endif
577b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
578b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
579b4fbaa2aSMark F. Adams #endif
580470e26f8SMark F. Adams     ierr = createProlongation(Aarr[level], data, pc_gamg->m_dim, pc_gamg->m_data_cols, pc_gamg->m_method,
581c20e4228SMark F. Adams                               level, pc_gamg->m_threshold, &bs, &Parr[level1], &coarse_data, &isOK,
582c20e4228SMark F. Adams                               &emaxs[level] );
5835b89ad90SMark F. Adams     CHKERRQ(ierr);
584d3d6bff4SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
585b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
586b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
587b4fbaa2aSMark F. Adams #endif
588baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
589baa4e9faSMark F. Adams     if( isOK ) {
590b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
591b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
592b4fbaa2aSMark F. Adams #endif
593470e26f8SMark F. Adams       ierr = partitionLevel( Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs,
594c20e4228SMark F. Adams                              &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1],
595c20e4228SMark F. Adams                              pc_gamg->m_avoid_repart, pc_gamg->m_min_eq_proc );
5963530afc2SMark F. Adams       CHKERRQ(ierr);
597b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
598b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
599b4fbaa2aSMark F. Adams #endif
6003530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
601737a81a9SMark F. Adams       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
60258471d46SMark F. Adams #ifdef VERBOSE
60358471d46SMark F. Adams       PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
6042c047e2dSMark F. Adams 		  mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols,
6052c047e2dSMark F. Adams 		  (int)(info.nz_used/(PetscReal)N),nactivepe);
60658471d46SMark F. Adams #endif
607e33ef3b1SMark F. Adams       /* coarse grids with SA can have zero row/cols from singleton aggregates */
60858471d46SMark F. Adams       /* aggregation method should gaurrentee this does not happen! */
609737a81a9SMark F. Adams 
61058471d46SMark F. Adams #ifdef VERBOSE
611e33ef3b1SMark F. Adams       if( PETSC_TRUE ){
612785cba28SMark F. Adams         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
613785cba28SMark F. Adams         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
614e33ef3b1SMark F. Adams         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
615785cba28SMark F. Adams         nloceq = Iend-Istart;
616e33ef3b1SMark F. Adams         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
617e33ef3b1SMark F. Adams         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
618e33ef3b1SMark F. Adams         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
619785cba28SMark F. Adams         for(kk=0;kk<nloceq;kk++){
620e33ef3b1SMark F. Adams           if(data_arr[kk]==0.0) {
621785cba28SMark F. Adams             id = kk + Istart;
622e33ef3b1SMark F. Adams             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
623e33ef3b1SMark F. Adams             CHKERRQ(ierr);
624ecd57171SMark F. Adams             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1);
625e33ef3b1SMark F. Adams           }
626e33ef3b1SMark F. Adams         }
627e33ef3b1SMark F. Adams         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
628e33ef3b1SMark F. Adams         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
629e33ef3b1SMark F. Adams         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
630e33ef3b1SMark F. Adams         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
631e33ef3b1SMark F. Adams       }
63258471d46SMark F. Adams #endif
633baa4e9faSMark F. Adams     }
634baa4e9faSMark F. Adams     else{
635be544d3cSMark F. Adams       coarse_data = 0;
636baa4e9faSMark F. Adams       break;
637baa4e9faSMark F. Adams     }
638eb07cef2SMark F. Adams     data = coarse_data;
639b4fbaa2aSMark F. Adams 
640b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
641b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
642b4fbaa2aSMark F. Adams #endif
6435b89ad90SMark F. Adams   }
644be544d3cSMark F. Adams   if( coarse_data ) {
645eb07cef2SMark F. Adams     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
646be544d3cSMark F. Adams   }
64758471d46SMark F. Adams #ifdef VERBOSE
648baa4e9faSMark F. Adams   PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
64958471d46SMark F. Adams #endif
6505b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
6515b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
6525b89ad90SMark F. Adams   fine_level = level;
653eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
6545b89ad90SMark F. Adams 
655fc4362bfSMark F. Adams   /* set default smoothers */
656587fa25dSMark F. Adams   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
657587fa25dSMark F. Adams         lidx <= fine_level;
658587fa25dSMark F. Adams         lidx++, level--) {
6595745e0f5SMark F. Adams     PetscReal emax, emin;
6605b89ad90SMark F. Adams     KSP smoother; PC subpc;
661587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
6625b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
663038e3b61SMark F. Adams     if( emaxs[level] > 0.0 ) emax=emaxs[level];
664587fa25dSMark F. Adams     else{ /* eigen estimate 'emax' */
665587fa25dSMark F. Adams       KSP eksp; Mat Lmat = Aarr[level];
666fc4362bfSMark F. Adams       Vec bb, xx; PC pc;
667038e3b61SMark F. Adams 
6685745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
6695745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
670fc4362bfSMark F. Adams       {
671fc4362bfSMark F. Adams 	PetscRandom    rctx;
672fc4362bfSMark F. Adams 	ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
673fc4362bfSMark F. Adams 	ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
674fc4362bfSMark F. Adams 	ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
675fc4362bfSMark F. Adams 	ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
6765b89ad90SMark F. Adams       }
677fc4362bfSMark F. Adams       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
678e33ef3b1SMark F. Adams       ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr);
679fc4362bfSMark F. Adams       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
68058471d46SMark F. Adams       ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
681fc4362bfSMark F. Adams       ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
68258471d46SMark F. Adams       ierr = PCSetType( pc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); /* should be same as eigen estimates op. */
683038e3b61SMark F. Adams       ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
684fc4362bfSMark F. Adams       CHKERRQ(ierr);
685fc4362bfSMark F. Adams       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
686fc4362bfSMark F. Adams 
687fc4362bfSMark F. Adams       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
688fc4362bfSMark F. Adams       ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
689fc4362bfSMark F. Adams       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
690fc4362bfSMark F. Adams       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
691fc4362bfSMark F. Adams       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
692fc4362bfSMark F. Adams       ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
69358471d46SMark F. Adams #ifdef VERBOSE
694*d8c859f0SMark F. Adams       PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s PC setup max eigen=%e min=%e PC=%s\n",__FUNCT__,emax,emin,PETSC_GAMG_SMOOTHER);
69558471d46SMark F. Adams #endif
696fc4362bfSMark F. Adams     }
697038e3b61SMark F. Adams     {
698038e3b61SMark F. Adams       PetscInt N1, N0, tt;
699038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
700038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
701785cba28SMark F. Adams       emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */
702038e3b61SMark F. Adams       emax *= 1.05;
703038e3b61SMark F. Adams     }
704038e3b61SMark F. Adams 
70558471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );
706fc4362bfSMark F. Adams     ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
7070e1b4bd6SMark F. Adams     /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */
7085745e0f5SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
709e33ef3b1SMark F. Adams     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
7105745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
7115745e0f5SMark F. Adams   }
7125745e0f5SMark F. Adams   {
713ecd57171SMark F. Adams     /* coarse grid */
714ecd57171SMark F. Adams     KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
7155745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
716eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
71758471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
7185745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
719ecd57171SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
720ecd57171SMark F. Adams     ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
721ecd57171SMark F. Adams     ierr = PCSetUp( subpc ); CHKERRQ(ierr);
722ecd57171SMark F. Adams     ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
723ecd57171SMark F. Adams     assert(ii==1);
724ecd57171SMark F. Adams     ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
725ecd57171SMark F. Adams     ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
726fc4362bfSMark F. Adams   }
727737a81a9SMark F. Adams 
728fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
729eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
7305b89ad90SMark F. Adams   {
7315b89ad90SMark F. Adams     PetscBool galerkin;
732eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
7335b89ad90SMark F. Adams     if(galerkin){
7345b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
7355b89ad90SMark F. Adams     }
7365b89ad90SMark F. Adams   }
7375745e0f5SMark F. Adams 
73858471d46SMark F. Adams   /* plot levels - R/P */
73958471d46SMark F. Adams   if( PRINT_MATS ) {
74058471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-1;level>0;level--){
74158471d46SMark F. Adams       PetscViewer viewer;
74258471d46SMark F. Adams       char fname[32];
74303a628feSMark F. Adams       sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
74411e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
7455b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
746038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
7475b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
74803a628feSMark F. Adams       sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
749e33ef3b1SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
750e33ef3b1SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
751e33ef3b1SMark F. Adams       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
752e33ef3b1SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
7535b89ad90SMark F. Adams     }
75458471d46SMark F. Adams   }
75558471d46SMark F. Adams 
75658471d46SMark F. Adams   /* set interpolation between the levels, clean up */
75758471d46SMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
75858471d46SMark F. Adams        lidx<fine_level;
75958471d46SMark F. Adams        lidx++, level--){
76058471d46SMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
761587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
762587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
7635b89ad90SMark F. Adams   }
7645745e0f5SMark F. Adams 
7655b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
766eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
767eb07cef2SMark F. Adams   ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
76803a628feSMark F. Adams   a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
76958471d46SMark F. Adams 
77058471d46SMark F. Adams   {
77158471d46SMark F. Adams     KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
77258471d46SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
77358471d46SMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
77458471d46SMark F. Adams     ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
77558471d46SMark F. Adams   }
776e33ef3b1SMark F. Adams 
7775b89ad90SMark F. Adams   PetscFunctionReturn(0);
7785b89ad90SMark F. Adams }
7795b89ad90SMark F. Adams 
780eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7815b89ad90SMark F. Adams /*
7825b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7835b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7845b89ad90SMark F. Adams 
7855b89ad90SMark F. Adams    Input Parameter:
7865b89ad90SMark F. Adams .  pc - the preconditioner context
7875b89ad90SMark F. Adams 
7885b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7895b89ad90SMark F. Adams */
7905b89ad90SMark F. Adams #undef __FUNCT__
7915b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
7925b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
7935b89ad90SMark F. Adams {
7945b89ad90SMark F. Adams   PetscErrorCode  ierr;
7955b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
7965b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
7975b89ad90SMark F. Adams 
7985b89ad90SMark F. Adams   PetscFunctionBegin;
7995b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8005b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8015b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8025b89ad90SMark F. Adams   PetscFunctionReturn(0);
8035b89ad90SMark F. Adams }
8045b89ad90SMark F. Adams 
805676e1743SMark F. Adams 
806676e1743SMark F. Adams #undef __FUNCT__
807676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
808676e1743SMark F. Adams /*@
809676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
810676e1743SMark F. Adams    processor reduction.
811676e1743SMark F. Adams 
812676e1743SMark F. Adams    Not Collective on PC
813676e1743SMark F. Adams 
814676e1743SMark F. Adams    Input Parameters:
815676e1743SMark F. Adams .  pc - the preconditioner context
816676e1743SMark F. Adams 
817676e1743SMark F. Adams 
818676e1743SMark F. Adams    Options Database Key:
819676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
820676e1743SMark F. Adams 
821676e1743SMark F. Adams    Level: intermediate
822676e1743SMark F. Adams 
823676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
824676e1743SMark F. Adams 
825676e1743SMark F. Adams .seealso: ()
826676e1743SMark F. Adams @*/
827676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
828676e1743SMark F. Adams {
829676e1743SMark F. Adams   PetscErrorCode ierr;
830676e1743SMark F. Adams 
831676e1743SMark F. Adams   PetscFunctionBegin;
832676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
833676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
834676e1743SMark F. Adams   PetscFunctionReturn(0);
835676e1743SMark F. Adams }
836676e1743SMark F. Adams 
837676e1743SMark F. Adams EXTERN_C_BEGIN
838676e1743SMark F. Adams #undef __FUNCT__
839676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
840676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
841676e1743SMark F. Adams {
842c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
843c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
844676e1743SMark F. Adams 
845676e1743SMark F. Adams   PetscFunctionBegin;
846c20e4228SMark F. Adams   if(n>0) pc_gamg->m_min_eq_proc = n;
847676e1743SMark F. Adams   PetscFunctionReturn(0);
848676e1743SMark F. Adams }
849676e1743SMark F. Adams EXTERN_C_END
850676e1743SMark F. Adams 
851676e1743SMark F. Adams #undef __FUNCT__
852676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning"
853676e1743SMark F. Adams /*@
854676e1743SMark F. Adams    PCGAMGAvoidRepartitioning - Do not repartition the coarse grids
855676e1743SMark F. Adams 
856676e1743SMark F. Adams    Collective on PC
857676e1743SMark F. Adams 
858676e1743SMark F. Adams    Input Parameters:
859676e1743SMark F. Adams .  pc - the preconditioner context
860676e1743SMark F. Adams 
861676e1743SMark F. Adams 
862676e1743SMark F. Adams    Options Database Key:
863676e1743SMark F. Adams .  -pc_gamg_avoid_repartitioning
864676e1743SMark F. Adams 
865676e1743SMark F. Adams    Level: intermediate
866676e1743SMark F. Adams 
867676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
868676e1743SMark F. Adams 
869676e1743SMark F. Adams .seealso: ()
870676e1743SMark F. Adams @*/
871676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning(PC pc, PetscBool n)
872676e1743SMark F. Adams {
873676e1743SMark F. Adams   PetscErrorCode ierr;
874676e1743SMark F. Adams 
875676e1743SMark F. Adams   PetscFunctionBegin;
876676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
877676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGAvoidRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
878676e1743SMark F. Adams   PetscFunctionReturn(0);
879676e1743SMark F. Adams }
880676e1743SMark F. Adams 
881676e1743SMark F. Adams EXTERN_C_BEGIN
882676e1743SMark F. Adams #undef __FUNCT__
883676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning_GAMG"
884676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning_GAMG(PC pc, PetscBool n)
885676e1743SMark F. Adams {
886c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
887c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
888676e1743SMark F. Adams 
889676e1743SMark F. Adams   PetscFunctionBegin;
890c20e4228SMark F. Adams   pc_gamg->m_avoid_repart = n;
891676e1743SMark F. Adams   PetscFunctionReturn(0);
892676e1743SMark F. Adams }
893676e1743SMark F. Adams EXTERN_C_END
894676e1743SMark F. Adams 
895676e1743SMark F. Adams #undef __FUNCT__
8963542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
8973542efc5SMark F. Adams /*@
8983542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
8993542efc5SMark F. Adams 
9003542efc5SMark F. Adams    Not collective on PC
9013542efc5SMark F. Adams 
9023542efc5SMark F. Adams    Input Parameters:
9033542efc5SMark F. Adams .  pc - the preconditioner context
9043542efc5SMark F. Adams 
9053542efc5SMark F. Adams 
9063542efc5SMark F. Adams    Options Database Key:
9073542efc5SMark F. Adams .  -pc_gamg_threshold
9083542efc5SMark F. Adams 
9093542efc5SMark F. Adams    Level: intermediate
9103542efc5SMark F. Adams 
9113542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
9123542efc5SMark F. Adams 
9133542efc5SMark F. Adams .seealso: ()
9143542efc5SMark F. Adams @*/
9153542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
9163542efc5SMark F. Adams {
9173542efc5SMark F. Adams   PetscErrorCode ierr;
9183542efc5SMark F. Adams 
9193542efc5SMark F. Adams   PetscFunctionBegin;
9203542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9213542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
9223542efc5SMark F. Adams   PetscFunctionReturn(0);
9233542efc5SMark F. Adams }
9243542efc5SMark F. Adams 
9253542efc5SMark F. Adams EXTERN_C_BEGIN
9263542efc5SMark F. Adams #undef __FUNCT__
9273542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
9283542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
9293542efc5SMark F. Adams {
930c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
931c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
9323542efc5SMark F. Adams 
9333542efc5SMark F. Adams   PetscFunctionBegin;
934c20e4228SMark F. Adams   pc_gamg->m_threshold = n;
9353542efc5SMark F. Adams   PetscFunctionReturn(0);
9363542efc5SMark F. Adams }
9373542efc5SMark F. Adams EXTERN_C_END
9383542efc5SMark F. Adams 
9393542efc5SMark F. Adams #undef __FUNCT__
940676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType"
941676e1743SMark F. Adams /*@
942676e1743SMark F. Adams    PCGAMGSetSolverType - Set solution method.
943676e1743SMark F. Adams 
944676e1743SMark F. Adams    Collective on PC
945676e1743SMark F. Adams 
946676e1743SMark F. Adams    Input Parameters:
947676e1743SMark F. Adams .  pc - the preconditioner context
948676e1743SMark F. Adams 
949676e1743SMark F. Adams 
950676e1743SMark F. Adams    Options Database Key:
9513542efc5SMark F. Adams .  -pc_gamg_type
952676e1743SMark F. Adams 
953676e1743SMark F. Adams    Level: intermediate
954676e1743SMark F. Adams 
955676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
956676e1743SMark F. Adams 
957676e1743SMark F. Adams .seealso: ()
958676e1743SMark F. Adams @*/
959676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz )
960676e1743SMark F. Adams {
961676e1743SMark F. Adams   PetscErrorCode ierr;
962676e1743SMark F. Adams 
963676e1743SMark F. Adams   PetscFunctionBegin;
964676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
965676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz));
966676e1743SMark F. Adams   CHKERRQ(ierr);
967676e1743SMark F. Adams   PetscFunctionReturn(0);
968676e1743SMark F. Adams }
969676e1743SMark F. Adams 
970676e1743SMark F. Adams EXTERN_C_BEGIN
971676e1743SMark F. Adams #undef __FUNCT__
972676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG"
973676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz )
974676e1743SMark F. Adams {
975676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
976676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
977676e1743SMark F. Adams 
978676e1743SMark F. Adams   PetscFunctionBegin;
979676e1743SMark F. Adams   if(sz < 64) strcpy(pc_gamg->m_type,str);
980676e1743SMark F. Adams   PetscFunctionReturn(0);
981676e1743SMark F. Adams }
982676e1743SMark F. Adams EXTERN_C_END
983676e1743SMark F. Adams 
9845b89ad90SMark F. Adams #undef __FUNCT__
9855b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
9865b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
9875b89ad90SMark F. Adams {
988676e1743SMark F. Adams   PetscErrorCode  ierr;
989676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
990676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
991676e1743SMark F. Adams   PetscBool flag;
9925b89ad90SMark F. Adams 
9935b89ad90SMark F. Adams   PetscFunctionBegin;
994676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
995676e1743SMark F. Adams   {
996676e1743SMark F. Adams     ierr = PetscOptionsString("-pc_gamg_type",
997470e26f8SMark F. Adams                               "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid (default)",
998676e1743SMark F. Adams                               "PCGAMGSetSolverType",
999676e1743SMark F. Adams                               pc_gamg->m_type,
1000676e1743SMark F. Adams                               pc_gamg->m_type,
1001676e1743SMark F. Adams                               64,
1002676e1743SMark F. Adams                               &flag );
1003676e1743SMark F. Adams     CHKERRQ(ierr);
1004*d8c859f0SMark F. Adams     if( flag && pc_gamg->m_data != 0 ) {
1005*d8c859f0SMark F. Adams       if( (strcmp(pc_gamg->m_type,"sa")==0 && pc_gamg->m_method != 2) ||
1006*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"pa")==0 && pc_gamg->m_method != 1) ||
1007*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"geo")==0 && pc_gamg->m_method != 0) ) {
1008*d8c859f0SMark F. Adams         SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "PCSetFromOptions called of PCSetCoordinates (with new method, after data was created)");
1009*d8c859f0SMark F. Adams       }
1010*d8c859f0SMark F. Adams     }
10112a44bfbaSMark F. Adams 
1012470e26f8SMark F. Adams     if (flag && strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2;
1013470e26f8SMark F. Adams     else if (flag && strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1;
1014470e26f8SMark F. Adams     else pc_gamg->m_method = 0;
10152a44bfbaSMark F. Adams 
1016c20e4228SMark F. Adams     /* -pc_gamg_avoid_repartitioning */
1017c20e4228SMark F. Adams     pc_gamg->m_avoid_repart = PETSC_FALSE;
1018676e1743SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_avoid_repartitioning",
1019676e1743SMark F. Adams                             "Do not repartion coarse grids (false)",
1020676e1743SMark F. Adams                             "PCGAMGAvoidRepartitioning",
1021c20e4228SMark F. Adams                             pc_gamg->m_avoid_repart,
1022c20e4228SMark F. Adams                             &pc_gamg->m_avoid_repart,
1023676e1743SMark F. Adams                             &flag);
1024676e1743SMark F. Adams     CHKERRQ(ierr);
1025676e1743SMark F. Adams 
1026c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1027c20e4228SMark F. Adams     pc_gamg->m_min_eq_proc = 600;
1028676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1029676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1030676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
1031c20e4228SMark F. Adams                            pc_gamg->m_min_eq_proc,
1032c20e4228SMark F. Adams                            &pc_gamg->m_min_eq_proc,
1033676e1743SMark F. Adams                            &flag );
1034676e1743SMark F. Adams     CHKERRQ(ierr);
10353542efc5SMark F. Adams 
1036c20e4228SMark F. Adams     /* -pc_gamg_threshold */
1037c20e4228SMark F. Adams     pc_gamg->m_threshold = 0.05;
10383542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
10393542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
10403542efc5SMark F. Adams                             "PCGAMGSetThreshold",
1041c20e4228SMark F. Adams                             pc_gamg->m_threshold,
1042c20e4228SMark F. Adams                             &pc_gamg->m_threshold,
10433542efc5SMark F. Adams                             &flag );
10443542efc5SMark F. Adams     CHKERRQ(ierr);
1045676e1743SMark F. Adams   }
1046676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1047676e1743SMark F. Adams 
10485b89ad90SMark F. Adams   PetscFunctionReturn(0);
10495b89ad90SMark F. Adams }
10505b89ad90SMark F. Adams 
10515b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
10525b89ad90SMark F. Adams /*MC
1053dc76db92SMark F. Adams      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning. This preconditioner currently has two
1054dc76db92SMark F. Adams            AMG methods: 1) an unstructured geometric method, which requires that you provide coordinates for each
1055dc76db92SMark F. Adams            vertex, and 2) smoothed aggregation.  Smoothed aggregation (SA) can work without coordinates but it
1056dc76db92SMark F. Adams            will generate some common non-trivial null spaces if coordinates are provided.  The input fine grid matrix
1057dc76db92SMark F. Adams            must have the block size set for 'system' problems (with multiple dof per vertex/cell) to work properly.
1058dc76db92SMark F. Adams            SA will generate rotational rigid body mode null space vectors, in addition to the trivial translational
1059dc76db92SMark F. Adams            modes, when coordinates are provide in 2D and 3D.
10605b89ad90SMark F. Adams 
1061280d9858SJed Brown    Options Database Keys:
10625b89ad90SMark F. Adams    Multigrid options(inherited)
1063280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1064280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1065280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1066280d9858SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
10675b89ad90SMark F. Adams 
10685b89ad90SMark F. Adams   Level: intermediate
1069280d9858SJed Brown 
10705b89ad90SMark F. Adams   Concepts: multigrid
10715b89ad90SMark F. Adams 
10725b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1073280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
10745b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
10755b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
10765b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
10775b89ad90SMark F. Adams M*/
10785b89ad90SMark F. Adams 
10795b89ad90SMark F. Adams EXTERN_C_BEGIN
10805b89ad90SMark F. Adams #undef __FUNCT__
10815b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
10825b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
10835b89ad90SMark F. Adams {
10845b89ad90SMark F. Adams   PetscErrorCode  ierr;
10855b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
10865b89ad90SMark F. Adams   PC_MG           *mg;
10875ef31b24SMark F. Adams   PetscClassId     cookie;
10885ee9c036SSatish Balay #if defined PETSC_USE_LOG
10895ee9c036SSatish Balay   static int count = 0;
10905ee9c036SSatish Balay #endif
10915b89ad90SMark F. Adams 
10925b89ad90SMark F. Adams   PetscFunctionBegin;
10935b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
10945b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
10955b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
10965b89ad90SMark F. Adams 
10975b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
10985b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
109903a628feSMark F. Adams   pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0;
11005b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
11015b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
11025b89ad90SMark F. Adams 
11035b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
11045b89ad90SMark F. Adams 
11055b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
11065b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
11075b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
11085b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
11095b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
11105b89ad90SMark F. Adams 
11115b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
11125b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
11135b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
1114c97e1df0SMark F. Adams 					    PCSetCoordinates_GAMG);
1115c97e1df0SMark F. Adams   CHKERRQ(ierr);
1116c97e1df0SMark F. Adams 
1117676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1118676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1119676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1120676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1121676e1743SMark F. Adams   CHKERRQ(ierr);
1122676e1743SMark F. Adams 
1123676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1124676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_C",
1125676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_GAMG",
1126676e1743SMark F. Adams 					    PCGAMGAvoidRepartitioning_GAMG);
1127676e1743SMark F. Adams   CHKERRQ(ierr);
1128676e1743SMark F. Adams 
1129676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
11303542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
11313542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
11323542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
11333542efc5SMark F. Adams   CHKERRQ(ierr);
11343542efc5SMark F. Adams 
11353542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1136676e1743SMark F. Adams 					    "PCGAMGSetSolverType_C",
1137676e1743SMark F. Adams 					    "PCGAMGSetSolverType_GAMG",
1138676e1743SMark F. Adams 					    PCGAMGSetSolverType_GAMG);
1139676e1743SMark F. Adams   CHKERRQ(ierr);
1140c97e1df0SMark F. Adams 
1141b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
1142785cba28SMark F. Adams   if( count++ == 0 ) {
11435ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
1144b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
11452a44bfbaSMark F. Adams     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
11462a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]);
11472a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]);
11482a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]);
1149b4fbaa2aSMark F. Adams     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1150b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1151b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1152b4fbaa2aSMark F. Adams     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1153b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_events[SET7]);
1154b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]); */
1155b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1156b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1157b4fbaa2aSMark F. Adams     PetscLogEventRegister("  PL repartition", cookie, &gamg_setup_events[SET12]);
1158b4fbaa2aSMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1159b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1160b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
11615ef31b24SMark F. Adams 
1162b4fbaa2aSMark F. Adams     /* create timer stages */
1163b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1164b4fbaa2aSMark F. Adams     {
1165b4fbaa2aSMark F. Adams       char str[32];
1166b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1167b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1168b4fbaa2aSMark F. Adams       PetscInt lidx;
1169b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1170b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1171b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1172b4fbaa2aSMark F. Adams       }
1173b4fbaa2aSMark F. Adams     }
1174b4fbaa2aSMark F. Adams #endif
1175b4fbaa2aSMark F. Adams   }
1176b4fbaa2aSMark F. Adams #endif
11775b89ad90SMark F. Adams   PetscFunctionReturn(0);
11785b89ad90SMark F. Adams }
11795b89ad90SMark F. Adams EXTERN_C_END
1180