xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision dc76db9228eb457bb7d236429860079c152c6c94)
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 */
18676e1743SMark F. Adams static PetscBool s_avoid_repart = PETSC_FALSE;
193542efc5SMark F. Adams static PetscInt s_min_eq_proc = 600;
20*d8c859f0SMark F. Adams static PetscReal s_threshold = 0.05;
215b89ad90SMark F. Adams typedef struct gamg_TAG {
225b89ad90SMark F. Adams   PetscInt       m_dim;
235b89ad90SMark F. Adams   PetscInt       m_Nlevels;
245b89ad90SMark F. Adams   PetscInt       m_data_sz;
25d3d6bff4SMark F. Adams   PetscInt       m_data_rows;
26d3d6bff4SMark F. Adams   PetscInt       m_data_cols;
2703a628feSMark F. Adams   PetscInt       m_count;
282c047e2dSMark F. Adams   PetscInt       m_method; /* 0: geomg; 1: plain agg 'pa'; 2: smoothed agg 'sa' */
295b89ad90SMark F. Adams   PetscReal     *m_data; /* blocked vector of vertex data on fine grid (coordinates) */
30676e1743SMark F. Adams   char           m_type[64];
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);
72eb07cef2SMark F. Adams     ierr = PetscMalloc((arrsz+1)*sizeof(double), &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)
1515b89ad90SMark F. Adams */
1525cb416c2SMark F. Adams 
153676e1743SMark F. Adams #define TOP_GRID_LIM 2*s_min_eq_proc /* this will happen anyway */
15422063be5SMark 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,
1643530afc2SMark F. Adams                                Mat *a_Amat_crs
16511e60469SMark F. Adams                                )
1665b89ad90SMark F. Adams {
1675b89ad90SMark F. Adams   PetscErrorCode   ierr;
168038e3b61SMark F. Adams   Mat              Cmat,Pnew,Pold=*a_P_inout;
16911e60469SMark F. Adams   IS               new_indices,isnum;
1703530afc2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
171fa31c87dSMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive;
172fa31c87dSMark F. Adams   PetscInt         neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0;
1735b89ad90SMark F. Adams 
1745b89ad90SMark F. Adams   PetscFunctionBegin;
17511e60469SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
17611e60469SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
17711e60469SMark F. Adams   /* RAP */
17896434bc1SMark F. Adams #ifdef USE_R
17996434bc1SMark F. Adams   /* make R wih brute force for now */
18096434bc1SMark F. Adams   ierr = MatTranspose( Pold, Pnew );
18196434bc1SMark F. Adams   ierr = MatDestroy( &Pold );  CHKERRQ(ierr);
18296434bc1SMark F. Adams   ierr = MatRARt( a_Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
18396434bc1SMark F. Adams   Pold = Pnew;
18496434bc1SMark F. Adams #else
185038e3b61SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
18696434bc1SMark F. Adams #endif
187038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
188038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
189038e3b61SMark F. Adams   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
190038e3b61SMark F. Adams 
19171d60426SBarry Smith   if( s_avoid_repart) {
19222063be5SMark F. Adams     *a_Amat_crs = Cmat; /* output */
19322063be5SMark F. Adams   }
19422063be5SMark F. Adams   else {
19522063be5SMark F. Adams     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
1965ef31b24SMark F. Adams     Mat              adj;
19722063be5SMark F. Adams     const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols;
198b4fbaa2aSMark F. Adams     const PetscInt  stride0=ncrs0*a_ndata_rows;
19992a756f0SMark F. Adams     PetscInt        is_sz,*isnewproc_idx,ii,jj,kk,strideNew,*tidx;
200c9a0b8beSMark F. Adams     /* create sub communicator  */
20171d60426SBarry Smith     MPI_Comm        cm;
202afc97cdcSMark F. Adams     MPI_Group       wg, g2;
203fd3c6afaSMark F. Adams     PetscInt       *counts,inpe;
204fd3c6afaSMark F. Adams     PetscMPIInt    *ranks;
20522063be5SMark F. Adams     IS              isscat;
20622063be5SMark F. Adams     PetscScalar    *array;
20722063be5SMark F. Adams     Vec             src_crd, dest_crd;
20822063be5SMark F. Adams     PetscReal      *data = *a_coarse_data;
20922063be5SMark F. Adams     VecScatter      vecscat;
210b4fbaa2aSMark F. Adams     IS  isnewproc;
211e33ef3b1SMark F. Adams 
21222063be5SMark F. Adams     /* get number of PEs to make active, reduce */
21322063be5SMark F. Adams     ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
214676e1743SMark F. Adams     new_npe = neq/s_min_eq_proc; /* hardwire min. number of eq/proc */
21522063be5SMark F. Adams     if( new_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1;
21622063be5SMark F. Adams     else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */
21722063be5SMark F. Adams 
21892a756f0SMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr);
219fd3c6afaSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
22092a756f0SMark F. Adams 
221fd3c6afaSMark F. Adams     ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr);
222afc97cdcSMark F. Adams     assert(counts[mype]==ncrs0);
223afc97cdcSMark F. Adams     /* count real active pes */
22422063be5SMark F. Adams     for( nactive = jj = 0 ; jj < npe ; jj++) {
225afc97cdcSMark F. Adams       if( counts[jj] != 0 ) {
22622063be5SMark F. Adams 	ranks[nactive++] = jj;
227afc97cdcSMark F. Adams         }
228afc97cdcSMark F. Adams     }
2293303bcf9Sadams 
2303303bcf9Sadams     if (nactive < new_npe) new_npe = nactive; /* this can happen with empty input procs */
23122063be5SMark F. Adams 
23258471d46SMark F. Adams #ifdef VERBOSE
23322063be5SMark 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);
23458471d46SMark F. Adams #endif
23558471d46SMark F. Adams 
23622063be5SMark F. Adams     *a_nactive_proc = new_npe; /* output */
2372f03bc48SMark F. Adams 
238afc97cdcSMark F. Adams     ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr);
23922063be5SMark F. Adams     ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr);
240afc97cdcSMark F. Adams     ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr);
241b4fbaa2aSMark F. Adams 
242afc97cdcSMark F. Adams     ierr = MPI_Group_free( &wg );                            CHKERRQ(ierr);
243afc97cdcSMark F. Adams     ierr = MPI_Group_free( &g2 );                            CHKERRQ(ierr);
244c9a0b8beSMark F. Adams 
2455ef31b24SMark F. Adams     /* MatPartitioningApply call MatConvert, which is collective */
246b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
247b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
248b4fbaa2aSMark F. Adams #endif
249038e3b61SMark F. Adams     if( a_cbs == 1) {
250038e3b61SMark F. Adams       ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
251eb07cef2SMark F. Adams     }
252eb07cef2SMark F. Adams     else{
253038e3b61SMark F. Adams       /* make a scalar matrix to partition */
254eb07cef2SMark F. Adams       Mat tMat;
25558471d46SMark F. Adams       PetscInt ncols,jj,Ii;
256b4fbaa2aSMark F. Adams       const PetscScalar *vals;
257b4fbaa2aSMark F. Adams       const PetscInt *idx;
25858471d46SMark F. Adams       PetscInt *d_nnz;
2595ee9c036SSatish Balay       static int llev = 0;
260b4fbaa2aSMark F. Adams 
26158471d46SMark F. Adams       ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
26258471d46SMark F. Adams       for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) {
26358471d46SMark F. Adams         ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
26458471d46SMark F. Adams         d_nnz[jj] = ncols/a_cbs;
26558471d46SMark F. Adams         if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
26658471d46SMark F. Adams         ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
26758471d46SMark F. Adams       }
2686876a03eSMark F. Adams 
269eb07cef2SMark F. Adams       ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
270eb07cef2SMark F. Adams                               PETSC_DETERMINE, PETSC_DETERMINE,
27158471d46SMark F. Adams                               0, d_nnz, 0, d_nnz,
272eb07cef2SMark F. Adams                               &tMat );
2736876a03eSMark F. Adams       CHKERRQ(ierr);
27458471d46SMark F. Adams       ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
275eb07cef2SMark F. Adams 
27622063be5SMark F. Adams       for ( ii = Istart0; ii < Iend0; ii++ ) {
27722063be5SMark F. Adams         PetscInt dest_row = ii/a_cbs;
27822063be5SMark F. Adams         ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
279eb07cef2SMark F. Adams         for( jj = 0 ; jj < ncols ; jj++ ){
280038e3b61SMark F. Adams           PetscInt dest_col = idx[jj]/a_cbs;
281eb07cef2SMark F. Adams           PetscScalar v = 1.0;
282eb07cef2SMark F. Adams           ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
283eb07cef2SMark F. Adams         }
28422063be5SMark F. Adams         ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
285eb07cef2SMark F. Adams       }
286eb07cef2SMark F. Adams       ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
287eb07cef2SMark F. Adams       ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
288eb07cef2SMark F. Adams 
289b4fbaa2aSMark F. Adams       if( llev++ == -1 ) {
290b4fbaa2aSMark F. Adams 	PetscViewer viewer; char fname[32];
291b4fbaa2aSMark F. Adams 	sprintf(fname,"part_mat_%d.mat",llev);
292b4fbaa2aSMark F. Adams 	PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
293b4fbaa2aSMark F. Adams 	ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
294b4fbaa2aSMark F. Adams 	ierr = PetscViewerDestroy( &viewer );
295b4fbaa2aSMark F. Adams       }
296b4fbaa2aSMark F. Adams 
297eb07cef2SMark F. Adams       ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
298eb07cef2SMark F. Adams 
299eb07cef2SMark F. Adams       ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
300eb07cef2SMark F. Adams     }
301afc97cdcSMark F. Adams     if( ncrs0 != 0 ){
302b4fbaa2aSMark F. Adams       const PetscInt *is_idx;
303b4fbaa2aSMark F. Adams       MatPartitioning  mpart;
3045ef31b24SMark F. Adams       /* hack to fix global data that pmetis.c uses in 'adj' */
30522063be5SMark F. Adams       for( nactive = jj = 0 ; jj < npe ; jj++ ) {
306afc97cdcSMark F. Adams 	if( counts[jj] != 0 ) {
30722063be5SMark F. Adams 	  adj->rmap->range[nactive++] = adj->rmap->range[jj];
308afc97cdcSMark F. Adams 	}
3095ef31b24SMark F. Adams       }
31022063be5SMark F. Adams       adj->rmap->range[nactive] = adj->rmap->range[npe];
3112f03bc48SMark F. Adams 
31271d60426SBarry Smith       ierr = MatPartitioningCreate( cm, &mpart ); CHKERRQ(ierr);
3135ef31b24SMark F. Adams       ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
31411e60469SMark F. Adams       ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
3155ef31b24SMark F. Adams       ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
31611e60469SMark F. Adams       ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr);
31711e60469SMark F. Adams       ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
31822063be5SMark F. Adams 
3195ef31b24SMark F. Adams       /* collect IS info */
3205ef31b24SMark F. Adams       ierr = ISGetLocalSize( isnewproc, &is_sz );       CHKERRQ(ierr);
321038e3b61SMark F. Adams       ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr);
3225ef31b24SMark F. Adams       ierr = ISGetIndices( isnewproc, &is_idx );        CHKERRQ(ierr);
323b4fbaa2aSMark F. Adams       /* spread partitioning across machine - best way ??? */
324ecd57171SMark F. Adams       NN = 1; /*npe/new_npe;*/
325b4fbaa2aSMark F. Adams       for( kk = jj = 0 ; kk < is_sz ; kk++ ){
326afc97cdcSMark F. Adams         for( ii = 0 ; ii < a_cbs ; ii++, jj++ ) {
32722063be5SMark F. Adams           isnewproc_idx[jj] = is_idx[kk] * NN; /* distribution */
328eb07cef2SMark F. Adams         }
3295ef31b24SMark F. Adams       }
3305ef31b24SMark F. Adams       ierr = ISRestoreIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
3315ef31b24SMark F. Adams       ierr = ISDestroy( &isnewproc );                    CHKERRQ(ierr);
33271d60426SBarry Smith       ierr = MPI_Comm_free( &cm );              CHKERRQ(ierr);
333b4fbaa2aSMark F. Adams 
334b4fbaa2aSMark F. Adams       is_sz *= a_cbs;
3355ef31b24SMark F. Adams     }
3365ef31b24SMark F. Adams     else{
3375ef31b24SMark F. Adams       isnewproc_idx = 0;
3385ef31b24SMark F. Adams       is_sz = 0;
3395ef31b24SMark F. Adams     }
3405ef31b24SMark F. Adams     ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
3415ef31b24SMark F. Adams     ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc );
342afc97cdcSMark F. Adams     if( isnewproc_idx != 0 ) {
3435ef31b24SMark F. Adams       ierr = PetscFree( isnewproc_idx );  CHKERRQ(ierr);
3445ef31b24SMark F. Adams     }
345e33ef3b1SMark F. Adams 
34611e60469SMark F. Adams     /*
34711e60469SMark F. Adams      Create an index set from the isnewproc index set to indicate the mapping TO
34811e60469SMark F. Adams      */
34911e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
35011e60469SMark F. Adams     /*
35111e60469SMark F. Adams      Determine how many elements are assigned to each processor
35211e60469SMark F. Adams      */
353fd3c6afaSMark F. Adams     inpe = npe;
354fd3c6afaSMark F. Adams     ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr);
35511e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
356038e3b61SMark F. Adams     ncrs_new = counts[mype]/a_cbs;
35722063be5SMark F. Adams     strideNew = ncrs_new*a_ndata_rows;
358b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
359b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
360b4fbaa2aSMark F. Adams #endif
36122063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
36211e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
363d3d6bff4SMark F. Adams     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
364470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
36511e60469SMark F. Adams     /*
366d3d6bff4SMark F. Adams      There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
367d3d6bff4SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
36811e60469SMark F. Adams      */
36992a756f0SMark F. Adams     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
37011e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
371038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
372d3d6bff4SMark F. Adams       PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
373d3d6bff4SMark F. Adams       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
37411e60469SMark F. Adams     }
375038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
376d3d6bff4SMark F. Adams     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
37711e60469SMark F. Adams     CHKERRQ(ierr);
37892a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
37911e60469SMark F. Adams     /*
38011e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
38111e60469SMark F. Adams      */
382d3d6bff4SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
383d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
384d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs0 ; ii++) {
385d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
386d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
387676e1743SMark F. Adams           PetscScalar tt = (PetscScalar)data[ix];
388676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
389d3d6bff4SMark F. Adams 	}
390038e3b61SMark F. Adams       }
391eb07cef2SMark F. Adams     }
392eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
393eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
39411e60469SMark F. Adams     /*
39511e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
39611e60469SMark F. Adams       to the correct processor
39711e60469SMark F. Adams     */
39811e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
39911e60469SMark F. Adams     CHKERRQ(ierr);
40011e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
40111e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40211e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40311e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
40411e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
40511e60469SMark F. Adams     /*
40611e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
40711e60469SMark F. Adams     */
408eb07cef2SMark F. Adams     ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
409d3d6bff4SMark F. Adams     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
410eb07cef2SMark F. Adams 
41111e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
412eb07cef2SMark F. Adams     data = *a_coarse_data;
413d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
414d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs_new ; ii++) {
415d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
416d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
417d3d6bff4SMark F. Adams 	  data[ix] = PetscRealPart(array[jx]);
418d3d6bff4SMark F. Adams 	  array[jx] = 1.e300;
419d3d6bff4SMark F. Adams 	}
420038e3b61SMark F. Adams       }
421038e3b61SMark F. Adams     }
42211e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
42311e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
42411e60469SMark F. Adams     /*
42511e60469SMark F. Adams       Invert for MatGetSubMatrix
42611e60469SMark F. Adams     */
427038e3b61SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
42811e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
42911e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
43011e60469SMark F. Adams     /* A_crs output */
431038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
43211e60469SMark F. Adams     CHKERRQ(ierr);
433eb07cef2SMark F. Adams 
434038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
435e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
436038e3b61SMark F. Adams     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
437eb07cef2SMark F. Adams 
43811e60469SMark F. Adams     /* prolongator */
43911e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
44011e60469SMark F. Adams     {
44111e60469SMark F. Adams       IS findices;
44211e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
44396434bc1SMark F. Adams #ifdef USE_R
44496434bc1SMark F. Adams       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
44596434bc1SMark F. Adams #else
44611e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
44796434bc1SMark F. Adams #endif
44811e60469SMark F. Adams       CHKERRQ(ierr);
449d61a3a7dSMark F. Adams 
45011e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
45111e60469SMark F. Adams     }
452d61a3a7dSMark F. Adams 
4533530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
4543530afc2SMark F. Adams     *a_P_inout = Pnew; /* output */
45592a756f0SMark F. Adams 
45611e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
45792a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
45892a756f0SMark F. Adams     ierr = PetscFree( ranks );  CHKERRQ(ierr);
459e33ef3b1SMark F. Adams   }
4605b89ad90SMark F. Adams 
4615b89ad90SMark F. Adams   PetscFunctionReturn(0);
4625b89ad90SMark F. Adams }
4635b89ad90SMark F. Adams 
4645b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4655b89ad90SMark F. Adams /*
4665b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4675b89ad90SMark F. Adams                     by setting data structures and options.
4685b89ad90SMark F. Adams 
4695b89ad90SMark F. Adams    Input Parameter:
4705b89ad90SMark F. Adams .  pc - the preconditioner context
4715b89ad90SMark F. Adams 
4725b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4735b89ad90SMark F. Adams 
4745b89ad90SMark F. Adams    Notes:
4755b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4765b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4775b89ad90SMark F. Adams */
4785b89ad90SMark F. Adams #undef __FUNCT__
4795b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
480eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
4815b89ad90SMark F. Adams {
4825b89ad90SMark F. Adams   PetscErrorCode  ierr;
483eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
4845b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
48558471d46SMark F. Adams   PC_MG_Levels   **mglevels = mg->levels;
486eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
487d3d6bff4SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
488eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
4893530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
490038e3b61SMark F. Adams   PetscBool        isOK;
491587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
492587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
493737a81a9SMark F. Adams   MatInfo          info;
4945ef31b24SMark F. Adams 
4955b89ad90SMark F. Adams   PetscFunctionBegin;
49603a628feSMark F. Adams   pc_gamg->m_count++;
49758471d46SMark F. Adams   if( a_pc->setupcalled > 0 ) {
49803a628feSMark F. Adams     /* just do Galerkin grids */
49958471d46SMark F. Adams     Mat B,dA,dB;
50058471d46SMark F. Adams 
50158471d46SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
50258471d46SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
50358471d46SMark F. Adams 
50458471d46SMark F. Adams     /* currently only handle case where mat and pmat are the same on coarser levels */
50558471d46SMark F. Adams     ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
50658471d46SMark F. Adams     /* (re)set to get dirty flag */
50758471d46SMark F. Adams     ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
50858471d46SMark F. Adams     ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr);
50958471d46SMark F. Adams 
51058471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-2; level>-1; level--) {
51158471d46SMark F. Adams       ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
51203a628feSMark F. Adams       /* the first time through the matrix structure has changed from repartitioning */
51303a628feSMark F. Adams       if( pc_gamg->m_count == 2 ) {
51403a628feSMark F. Adams         ierr = MatDestroy( &B );  CHKERRQ(ierr);
51503a628feSMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
51603a628feSMark F. Adams         mglevels[level]->A = B;
51703a628feSMark F. Adams       }
51803a628feSMark F. Adams       else {
51958471d46SMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
52003a628feSMark F. Adams       }
52158471d46SMark F. Adams       ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
52258471d46SMark F. Adams       dB = B;
52358471d46SMark F. Adams       /* setup KSP/PC */
52458471d46SMark F. Adams       ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
52558471d46SMark F. Adams     }
52658471d46SMark F. Adams 
52758471d46SMark F. Adams #define PRINT_MATS !PETSC_TRUE
52858471d46SMark F. Adams     /* plot levels - A */
52958471d46SMark F. Adams     if( PRINT_MATS ) {
53058471d46SMark F. Adams       for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){
53158471d46SMark F. Adams         PetscViewer viewer;
53258471d46SMark F. Adams         char fname[32]; KSP smoother; Mat Tmat, TTm;
53358471d46SMark F. Adams         ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
53458471d46SMark F. Adams         ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr);
53503a628feSMark F. Adams         sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
53658471d46SMark F. Adams         ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
53758471d46SMark F. Adams         ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
53858471d46SMark F. Adams         ierr = MatView( Tmat, viewer ); CHKERRQ(ierr);
53958471d46SMark F. Adams         ierr = PetscViewerDestroy( &viewer );
54058471d46SMark F. Adams       }
54158471d46SMark F. Adams     }
54258471d46SMark F. Adams 
54303a628feSMark F. Adams     a_pc->setupcalled = 2;
54403a628feSMark F. Adams 
54558471d46SMark F. Adams     PetscFunctionReturn(0);
546eb07cef2SMark F. Adams   }
547baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
548baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
5495b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
5505b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
5513530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
552eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
553eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
554eb07cef2SMark F. Adams 
555038e3b61SMark F. Adams   /* get data of not around */
556038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
557038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
558eb07cef2SMark F. Adams   }
559eb07cef2SMark F. Adams   data = pc_gamg->m_data;
560038e3b61SMark F. Adams 
561eb07cef2SMark F. Adams   /* Get A_i and R_i */
562737a81a9SMark F. Adams   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
56358471d46SMark F. Adams #ifdef VERBOSE
5642c047e2dSMark 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",
5652c047e2dSMark F. Adams 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,
5662c047e2dSMark F. Adams 	      (int)(info.nz_used/(PetscReal)N),npe);
56758471d46SMark F. Adams #endif
5688f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
56969f22f16SMark F. Adams         level < (GAMG_MAXLEVELS-1) && (level==0 || M>TOP_GRID_LIM); /* && (npe==1 || nactivepe>1); */
5700205a208SMark F. Adams         level++ ){
5715b89ad90SMark F. Adams     level1 = level + 1;
572b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
573b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
574b4fbaa2aSMark F. Adams #endif
575b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
576b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
577b4fbaa2aSMark F. Adams #endif
578470e26f8SMark F. Adams     ierr = createProlongation(Aarr[level], data, pc_gamg->m_dim, pc_gamg->m_data_cols, pc_gamg->m_method,
5793542efc5SMark F. Adams                               level, s_threshold, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level] );
5805b89ad90SMark F. Adams     CHKERRQ(ierr);
581d3d6bff4SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
582b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
583b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
584b4fbaa2aSMark F. Adams #endif
585baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
586baa4e9faSMark F. Adams     if( isOK ) {
587b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
588b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
589b4fbaa2aSMark F. Adams #endif
590470e26f8SMark F. Adams       ierr = partitionLevel( Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs,
591eb07cef2SMark F. Adams                              &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
5923530afc2SMark F. Adams       CHKERRQ(ierr);
593b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
594b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
595b4fbaa2aSMark F. Adams #endif
5963530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
597737a81a9SMark F. Adams       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
59858471d46SMark F. Adams #ifdef VERBOSE
59958471d46SMark 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",
6002c047e2dSMark F. Adams 		  mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols,
6012c047e2dSMark F. Adams 		  (int)(info.nz_used/(PetscReal)N),nactivepe);
60258471d46SMark F. Adams #endif
603e33ef3b1SMark F. Adams       /* coarse grids with SA can have zero row/cols from singleton aggregates */
60458471d46SMark F. Adams       /* aggregation method should gaurrentee this does not happen! */
605737a81a9SMark F. Adams 
60658471d46SMark F. Adams #ifdef VERBOSE
607e33ef3b1SMark F. Adams       if( PETSC_TRUE ){
608785cba28SMark F. Adams         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
609785cba28SMark F. Adams         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
610e33ef3b1SMark F. Adams         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
611785cba28SMark F. Adams         nloceq = Iend-Istart;
612e33ef3b1SMark F. Adams         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
613e33ef3b1SMark F. Adams         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
614e33ef3b1SMark F. Adams         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
615785cba28SMark F. Adams         for(kk=0;kk<nloceq;kk++){
616e33ef3b1SMark F. Adams           if(data_arr[kk]==0.0) {
617785cba28SMark F. Adams             id = kk + Istart;
618e33ef3b1SMark F. Adams             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
619e33ef3b1SMark F. Adams             CHKERRQ(ierr);
620ecd57171SMark F. Adams             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1);
621e33ef3b1SMark F. Adams           }
622e33ef3b1SMark F. Adams         }
623e33ef3b1SMark F. Adams         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
624e33ef3b1SMark F. Adams         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
625e33ef3b1SMark F. Adams         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
626e33ef3b1SMark F. Adams         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
627e33ef3b1SMark F. Adams       }
62858471d46SMark F. Adams #endif
629baa4e9faSMark F. Adams     }
630baa4e9faSMark F. Adams     else{
631be544d3cSMark F. Adams       coarse_data = 0;
632baa4e9faSMark F. Adams       break;
633baa4e9faSMark F. Adams     }
634eb07cef2SMark F. Adams     data = coarse_data;
635b4fbaa2aSMark F. Adams 
636b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
637b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
638b4fbaa2aSMark F. Adams #endif
6395b89ad90SMark F. Adams   }
640be544d3cSMark F. Adams   if( coarse_data ) {
641eb07cef2SMark F. Adams     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
642be544d3cSMark F. Adams   }
64358471d46SMark F. Adams #ifdef VERBOSE
644baa4e9faSMark F. Adams   PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
64558471d46SMark F. Adams #endif
6465b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
6475b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
6485b89ad90SMark F. Adams   fine_level = level;
649eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
6505b89ad90SMark F. Adams 
651fc4362bfSMark F. Adams   /* set default smoothers */
652587fa25dSMark F. Adams   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
653587fa25dSMark F. Adams         lidx <= fine_level;
654587fa25dSMark F. Adams         lidx++, level--) {
6555745e0f5SMark F. Adams     PetscReal emax, emin;
6565b89ad90SMark F. Adams     KSP smoother; PC subpc;
657587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
6585b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
659038e3b61SMark F. Adams     if( emaxs[level] > 0.0 ) emax=emaxs[level];
660587fa25dSMark F. Adams     else{ /* eigen estimate 'emax' */
661587fa25dSMark F. Adams       KSP eksp; Mat Lmat = Aarr[level];
662fc4362bfSMark F. Adams       Vec bb, xx; PC pc;
663038e3b61SMark F. Adams 
6645745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
6655745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
666fc4362bfSMark F. Adams       {
667fc4362bfSMark F. Adams 	PetscRandom    rctx;
668fc4362bfSMark F. Adams 	ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
669fc4362bfSMark F. Adams 	ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
670fc4362bfSMark F. Adams 	ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
671fc4362bfSMark F. Adams 	ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
6725b89ad90SMark F. Adams       }
673fc4362bfSMark F. Adams       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
674e33ef3b1SMark F. Adams       ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr);
675fc4362bfSMark F. Adams       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
67658471d46SMark F. Adams       ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
677fc4362bfSMark F. Adams       ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
67858471d46SMark F. Adams       ierr = PCSetType( pc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); /* should be same as eigen estimates op. */
679038e3b61SMark F. Adams       ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
680fc4362bfSMark F. Adams       CHKERRQ(ierr);
681fc4362bfSMark F. Adams       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
682fc4362bfSMark F. Adams 
683fc4362bfSMark F. Adams       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
684fc4362bfSMark F. Adams       ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
685fc4362bfSMark F. Adams       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
686fc4362bfSMark F. Adams       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
687fc4362bfSMark F. Adams       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
688fc4362bfSMark F. Adams       ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
68958471d46SMark F. Adams #ifdef VERBOSE
690*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);
69158471d46SMark F. Adams #endif
692fc4362bfSMark F. Adams     }
693038e3b61SMark F. Adams     {
694038e3b61SMark F. Adams       PetscInt N1, N0, tt;
695038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
696038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
697785cba28SMark F. Adams       emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */
698038e3b61SMark F. Adams       emax *= 1.05;
699038e3b61SMark F. Adams     }
700038e3b61SMark F. Adams 
70158471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );
702fc4362bfSMark F. Adams     ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
7030e1b4bd6SMark F. Adams     /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */
7045745e0f5SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
705e33ef3b1SMark F. Adams     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
7065745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
7075745e0f5SMark F. Adams   }
7085745e0f5SMark F. Adams   {
709ecd57171SMark F. Adams     /* coarse grid */
710ecd57171SMark F. Adams     KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
7115745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
712eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
71358471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
7145745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
715ecd57171SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
716ecd57171SMark F. Adams     ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
717ecd57171SMark F. Adams     ierr = PCSetUp( subpc ); CHKERRQ(ierr);
718ecd57171SMark F. Adams     ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
719ecd57171SMark F. Adams     assert(ii==1);
720ecd57171SMark F. Adams     ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
721ecd57171SMark F. Adams     ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
722fc4362bfSMark F. Adams   }
723737a81a9SMark F. Adams 
724fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
725eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
7265b89ad90SMark F. Adams   {
7275b89ad90SMark F. Adams     PetscBool galerkin;
728eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
7295b89ad90SMark F. Adams     if(galerkin){
7305b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
7315b89ad90SMark F. Adams     }
7325b89ad90SMark F. Adams   }
7335745e0f5SMark F. Adams 
73458471d46SMark F. Adams   /* plot levels - R/P */
73558471d46SMark F. Adams   if( PRINT_MATS ) {
73658471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-1;level>0;level--){
73758471d46SMark F. Adams       PetscViewer viewer;
73858471d46SMark F. Adams       char fname[32];
73903a628feSMark F. Adams       sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
74011e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
7415b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
742038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
7435b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
74403a628feSMark F. Adams       sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
745e33ef3b1SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
746e33ef3b1SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
747e33ef3b1SMark F. Adams       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
748e33ef3b1SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
7495b89ad90SMark F. Adams     }
75058471d46SMark F. Adams   }
75158471d46SMark F. Adams 
75258471d46SMark F. Adams   /* set interpolation between the levels, clean up */
75358471d46SMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
75458471d46SMark F. Adams        lidx<fine_level;
75558471d46SMark F. Adams        lidx++, level--){
75658471d46SMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
757587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
758587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
7595b89ad90SMark F. Adams   }
7605745e0f5SMark F. Adams 
7615b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
762eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
763eb07cef2SMark F. Adams   ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
76403a628feSMark F. Adams   a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
76558471d46SMark F. Adams 
76658471d46SMark F. Adams   {
76758471d46SMark F. Adams     KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
76858471d46SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
76958471d46SMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
77058471d46SMark F. Adams     ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
77158471d46SMark F. Adams   }
772e33ef3b1SMark F. Adams 
7735b89ad90SMark F. Adams   PetscFunctionReturn(0);
7745b89ad90SMark F. Adams }
7755b89ad90SMark F. Adams 
776eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7775b89ad90SMark F. Adams /*
7785b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7795b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7805b89ad90SMark F. Adams 
7815b89ad90SMark F. Adams    Input Parameter:
7825b89ad90SMark F. Adams .  pc - the preconditioner context
7835b89ad90SMark F. Adams 
7845b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7855b89ad90SMark F. Adams */
7865b89ad90SMark F. Adams #undef __FUNCT__
7875b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
7885b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
7895b89ad90SMark F. Adams {
7905b89ad90SMark F. Adams   PetscErrorCode  ierr;
7915b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
7925b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
7935b89ad90SMark F. Adams 
7945b89ad90SMark F. Adams   PetscFunctionBegin;
7955b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7965b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7975b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7985b89ad90SMark F. Adams   PetscFunctionReturn(0);
7995b89ad90SMark F. Adams }
8005b89ad90SMark F. Adams 
801676e1743SMark F. Adams 
802676e1743SMark F. Adams #undef __FUNCT__
803676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
804676e1743SMark F. Adams /*@
805676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
806676e1743SMark F. Adams    processor reduction.
807676e1743SMark F. Adams 
808676e1743SMark F. Adams    Not Collective on PC
809676e1743SMark F. Adams 
810676e1743SMark F. Adams    Input Parameters:
811676e1743SMark F. Adams .  pc - the preconditioner context
812676e1743SMark F. Adams 
813676e1743SMark F. Adams 
814676e1743SMark F. Adams    Options Database Key:
815676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
816676e1743SMark F. Adams 
817676e1743SMark F. Adams    Level: intermediate
818676e1743SMark F. Adams 
819676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
820676e1743SMark F. Adams 
821676e1743SMark F. Adams .seealso: ()
822676e1743SMark F. Adams @*/
823676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
824676e1743SMark F. Adams {
825676e1743SMark F. Adams   PetscErrorCode ierr;
826676e1743SMark F. Adams 
827676e1743SMark F. Adams   PetscFunctionBegin;
828676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
829676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
830676e1743SMark F. Adams   PetscFunctionReturn(0);
831676e1743SMark F. Adams }
832676e1743SMark F. Adams 
833676e1743SMark F. Adams EXTERN_C_BEGIN
834676e1743SMark F. Adams #undef __FUNCT__
835676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
836676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
837676e1743SMark F. Adams {
838676e1743SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
839676e1743SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
840676e1743SMark F. Adams 
841676e1743SMark F. Adams   PetscFunctionBegin;
842676e1743SMark F. Adams   if(n>0) s_min_eq_proc = n;
843676e1743SMark F. Adams   PetscFunctionReturn(0);
844676e1743SMark F. Adams }
845676e1743SMark F. Adams EXTERN_C_END
846676e1743SMark F. Adams 
847676e1743SMark F. Adams #undef __FUNCT__
848676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning"
849676e1743SMark F. Adams /*@
850676e1743SMark F. Adams    PCGAMGAvoidRepartitioning - Do not repartition the coarse grids
851676e1743SMark F. Adams 
852676e1743SMark F. Adams    Collective on PC
853676e1743SMark F. Adams 
854676e1743SMark F. Adams    Input Parameters:
855676e1743SMark F. Adams .  pc - the preconditioner context
856676e1743SMark F. Adams 
857676e1743SMark F. Adams 
858676e1743SMark F. Adams    Options Database Key:
859676e1743SMark F. Adams .  -pc_gamg_avoid_repartitioning
860676e1743SMark F. Adams 
861676e1743SMark F. Adams    Level: intermediate
862676e1743SMark F. Adams 
863676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
864676e1743SMark F. Adams 
865676e1743SMark F. Adams .seealso: ()
866676e1743SMark F. Adams @*/
867676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning(PC pc, PetscBool n)
868676e1743SMark F. Adams {
869676e1743SMark F. Adams   PetscErrorCode ierr;
870676e1743SMark F. Adams 
871676e1743SMark F. Adams   PetscFunctionBegin;
872676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
873676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGAvoidRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
874676e1743SMark F. Adams   PetscFunctionReturn(0);
875676e1743SMark F. Adams }
876676e1743SMark F. Adams 
877676e1743SMark F. Adams EXTERN_C_BEGIN
878676e1743SMark F. Adams #undef __FUNCT__
879676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning_GAMG"
880676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning_GAMG(PC pc, PetscBool n)
881676e1743SMark F. Adams {
882676e1743SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
883676e1743SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
884676e1743SMark F. Adams 
885676e1743SMark F. Adams   PetscFunctionBegin;
886676e1743SMark F. Adams   s_avoid_repart = n;
887676e1743SMark F. Adams   PetscFunctionReturn(0);
888676e1743SMark F. Adams }
889676e1743SMark F. Adams EXTERN_C_END
890676e1743SMark F. Adams 
891676e1743SMark F. Adams #undef __FUNCT__
8923542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
8933542efc5SMark F. Adams /*@
8943542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
8953542efc5SMark F. Adams 
8963542efc5SMark F. Adams    Not collective on PC
8973542efc5SMark F. Adams 
8983542efc5SMark F. Adams    Input Parameters:
8993542efc5SMark F. Adams .  pc - the preconditioner context
9003542efc5SMark F. Adams 
9013542efc5SMark F. Adams 
9023542efc5SMark F. Adams    Options Database Key:
9033542efc5SMark F. Adams .  -pc_gamg_threshold
9043542efc5SMark F. Adams 
9053542efc5SMark F. Adams    Level: intermediate
9063542efc5SMark F. Adams 
9073542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
9083542efc5SMark F. Adams 
9093542efc5SMark F. Adams .seealso: ()
9103542efc5SMark F. Adams @*/
9113542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
9123542efc5SMark F. Adams {
9133542efc5SMark F. Adams   PetscErrorCode ierr;
9143542efc5SMark F. Adams 
9153542efc5SMark F. Adams   PetscFunctionBegin;
9163542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9173542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
9183542efc5SMark F. Adams   PetscFunctionReturn(0);
9193542efc5SMark F. Adams }
9203542efc5SMark F. Adams 
9213542efc5SMark F. Adams EXTERN_C_BEGIN
9223542efc5SMark F. Adams #undef __FUNCT__
9233542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
9243542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
9253542efc5SMark F. Adams {
9263542efc5SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
9273542efc5SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
9283542efc5SMark F. Adams 
9293542efc5SMark F. Adams   PetscFunctionBegin;
9303542efc5SMark F. Adams   s_threshold = n;
9313542efc5SMark F. Adams   PetscFunctionReturn(0);
9323542efc5SMark F. Adams }
9333542efc5SMark F. Adams EXTERN_C_END
9343542efc5SMark F. Adams 
9353542efc5SMark F. Adams #undef __FUNCT__
936676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType"
937676e1743SMark F. Adams /*@
938676e1743SMark F. Adams    PCGAMGSetSolverType - Set solution method.
939676e1743SMark F. Adams 
940676e1743SMark F. Adams    Collective on PC
941676e1743SMark F. Adams 
942676e1743SMark F. Adams    Input Parameters:
943676e1743SMark F. Adams .  pc - the preconditioner context
944676e1743SMark F. Adams 
945676e1743SMark F. Adams 
946676e1743SMark F. Adams    Options Database Key:
9473542efc5SMark F. Adams .  -pc_gamg_type
948676e1743SMark F. Adams 
949676e1743SMark F. Adams    Level: intermediate
950676e1743SMark F. Adams 
951676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
952676e1743SMark F. Adams 
953676e1743SMark F. Adams .seealso: ()
954676e1743SMark F. Adams @*/
955676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz )
956676e1743SMark F. Adams {
957676e1743SMark F. Adams   PetscErrorCode ierr;
958676e1743SMark F. Adams 
959676e1743SMark F. Adams   PetscFunctionBegin;
960676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
961676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz));
962676e1743SMark F. Adams   CHKERRQ(ierr);
963676e1743SMark F. Adams   PetscFunctionReturn(0);
964676e1743SMark F. Adams }
965676e1743SMark F. Adams 
966676e1743SMark F. Adams EXTERN_C_BEGIN
967676e1743SMark F. Adams #undef __FUNCT__
968676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG"
969676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz )
970676e1743SMark F. Adams {
971676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
972676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
973676e1743SMark F. Adams 
974676e1743SMark F. Adams   PetscFunctionBegin;
975676e1743SMark F. Adams   if(sz < 64) strcpy(pc_gamg->m_type,str);
976676e1743SMark F. Adams   PetscFunctionReturn(0);
977676e1743SMark F. Adams }
978676e1743SMark F. Adams EXTERN_C_END
979676e1743SMark F. Adams 
9805b89ad90SMark F. Adams #undef __FUNCT__
9815b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
9825b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
9835b89ad90SMark F. Adams {
984676e1743SMark F. Adams   PetscErrorCode  ierr;
985676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
986676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
987676e1743SMark F. Adams   PetscBool flag;
9885b89ad90SMark F. Adams 
9895b89ad90SMark F. Adams   PetscFunctionBegin;
990676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
991676e1743SMark F. Adams   {
992676e1743SMark F. Adams     ierr = PetscOptionsString("-pc_gamg_type",
993470e26f8SMark F. Adams                               "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid (default)",
994676e1743SMark F. Adams                               "PCGAMGSetSolverType",
995676e1743SMark F. Adams                               pc_gamg->m_type,
996676e1743SMark F. Adams                               pc_gamg->m_type,
997676e1743SMark F. Adams                               64,
998676e1743SMark F. Adams                               &flag );
999676e1743SMark F. Adams     CHKERRQ(ierr);
1000*d8c859f0SMark F. Adams     if( flag && pc_gamg->m_data != 0 ) {
1001*d8c859f0SMark F. Adams       if( (strcmp(pc_gamg->m_type,"sa")==0 && pc_gamg->m_method != 2) ||
1002*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"pa")==0 && pc_gamg->m_method != 1) ||
1003*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"geo")==0 && pc_gamg->m_method != 0) ) {
1004*d8c859f0SMark F. Adams         SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "PCSetFromOptions called of PCSetCoordinates (with new method, after data was created)");
1005*d8c859f0SMark F. Adams       }
1006*d8c859f0SMark F. Adams     }
10072a44bfbaSMark F. Adams 
1008470e26f8SMark F. Adams     if (flag && strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2;
1009470e26f8SMark F. Adams     else if (flag && strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1;
1010470e26f8SMark F. Adams     else pc_gamg->m_method = 0;
10112a44bfbaSMark F. Adams 
10123542efc5SMark F. Adams     /* common (static) variable */
1013676e1743SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_avoid_repartitioning",
1014676e1743SMark F. Adams                             "Do not repartion coarse grids (false)",
1015676e1743SMark F. Adams                             "PCGAMGAvoidRepartitioning",
1016676e1743SMark F. Adams                             s_avoid_repart,
1017676e1743SMark F. Adams                             &s_avoid_repart,
1018676e1743SMark F. Adams                             &flag);
1019676e1743SMark F. Adams     CHKERRQ(ierr);
1020676e1743SMark F. Adams 
10213542efc5SMark F. Adams     /* common (static) variable */
1022676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1023676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1024676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
1025676e1743SMark F. Adams                            s_min_eq_proc,
1026676e1743SMark F. Adams                            &s_min_eq_proc,
1027676e1743SMark F. Adams                            &flag );
1028676e1743SMark F. Adams     CHKERRQ(ierr);
10293542efc5SMark F. Adams 
10303542efc5SMark F. Adams     /* common (static) variable */
10313542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
10323542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
10333542efc5SMark F. Adams                             "PCGAMGSetThreshold",
10343542efc5SMark F. Adams                             s_threshold,
10353542efc5SMark F. Adams                             &s_threshold,
10363542efc5SMark F. Adams                             &flag );
10373542efc5SMark F. Adams     CHKERRQ(ierr);
1038676e1743SMark F. Adams   }
1039676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1040676e1743SMark F. Adams 
10415b89ad90SMark F. Adams   PetscFunctionReturn(0);
10425b89ad90SMark F. Adams }
10435b89ad90SMark F. Adams 
10445b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
10455b89ad90SMark F. Adams /*
10465b89ad90SMark F. Adams  PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG
10475b89ad90SMark F. Adams 
10485b89ad90SMark F. Adams    Input Parameter:
10495b89ad90SMark F. Adams .  pc - the preconditioner context
10505b89ad90SMark F. Adams 
10515b89ad90SMark F. Adams    Application Interface Routine: PCCreate()
10525b89ad90SMark F. Adams 
10535b89ad90SMark F. Adams   */
10545b89ad90SMark F. Adams  /* MC
1055dc76db92SMark F. Adams      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning. This preconditioner currently has two
1056dc76db92SMark F. Adams            AMG methods: 1) an unstructured geometric method, which requires that you provide coordinates for each
1057dc76db92SMark F. Adams            vertex, and 2) smoothed aggregation.  Smoothed aggregation (SA) can work without coordinates but it
1058dc76db92SMark F. Adams            will generate some common non-trivial null spaces if coordinates are provided.  The input fine grid matrix
1059dc76db92SMark F. Adams            must have the block size set for 'system' problems (with multiple dof per vertex/cell) to work properly.
1060dc76db92SMark F. Adams            SA will generate rotational rigid body mode null space vectors, in addition to the trivial translational
1061dc76db92SMark F. Adams            modes, when coordinates are provide in 2D and 3D.
10625b89ad90SMark F. Adams 
10635b89ad90SMark F. Adams    Options Database Key:
10645b89ad90SMark F. Adams    Multigrid options(inherited)
10655b89ad90SMark F. Adams +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
10665b89ad90SMark F. Adams .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
10675b89ad90SMark F. Adams .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
10685b89ad90SMark F. Adams    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
10695b89ad90SMark F. Adams    GAMG options:
10705b89ad90SMark F. Adams 
10715b89ad90SMark F. Adams    Level: intermediate
10725b89ad90SMark F. Adams   Concepts: multigrid
10735b89ad90SMark F. Adams 
10745b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
10755b89ad90SMark F. Adams            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
10765b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
10775b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
10785b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
10795b89ad90SMark F. Adams M */
10805b89ad90SMark F. Adams 
10815b89ad90SMark F. Adams EXTERN_C_BEGIN
10825b89ad90SMark F. Adams #undef __FUNCT__
10835b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
10845b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
10855b89ad90SMark F. Adams {
10865b89ad90SMark F. Adams   PetscErrorCode  ierr;
10875b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
10885b89ad90SMark F. Adams   PC_MG           *mg;
10895ef31b24SMark F. Adams   PetscClassId     cookie;
10905ee9c036SSatish Balay #if defined PETSC_USE_LOG
10915ee9c036SSatish Balay   static int count = 0;
10925ee9c036SSatish Balay #endif
10935b89ad90SMark F. Adams 
10945b89ad90SMark F. Adams   PetscFunctionBegin;
10955b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
10965b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
10975b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
10985b89ad90SMark F. Adams 
10995b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
11005b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
110103a628feSMark F. Adams   pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0;
11025b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
11035b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
11045b89ad90SMark F. Adams 
11055b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
11065b89ad90SMark F. Adams 
11075b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
11085b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
11095b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
11105b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
11115b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
11125b89ad90SMark F. Adams 
11135b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
11145b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
11155b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
1116c97e1df0SMark F. Adams 					    PCSetCoordinates_GAMG);
1117c97e1df0SMark F. Adams   CHKERRQ(ierr);
1118c97e1df0SMark F. Adams 
1119676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1120676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1121676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1122676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1123676e1743SMark F. Adams   CHKERRQ(ierr);
1124676e1743SMark F. Adams 
1125676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1126676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_C",
1127676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_GAMG",
1128676e1743SMark F. Adams 					    PCGAMGAvoidRepartitioning_GAMG);
1129676e1743SMark F. Adams   CHKERRQ(ierr);
1130676e1743SMark F. Adams 
1131676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
11323542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
11333542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
11343542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
11353542efc5SMark F. Adams   CHKERRQ(ierr);
11363542efc5SMark F. Adams 
11373542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1138676e1743SMark F. Adams 					    "PCGAMGSetSolverType_C",
1139676e1743SMark F. Adams 					    "PCGAMGSetSolverType_GAMG",
1140676e1743SMark F. Adams 					    PCGAMGSetSolverType_GAMG);
1141676e1743SMark F. Adams   CHKERRQ(ierr);
1142c97e1df0SMark F. Adams 
1143b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
1144785cba28SMark F. Adams   if( count++ == 0 ) {
11455ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
1146b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
11472a44bfbaSMark F. Adams     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
11482a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]);
11492a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]);
11502a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]);
1151b4fbaa2aSMark F. Adams     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1152b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1153b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1154b4fbaa2aSMark F. Adams     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1155b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_events[SET7]);
1156b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]); */
1157b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1158b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1159b4fbaa2aSMark F. Adams     PetscLogEventRegister("  PL repartition", cookie, &gamg_setup_events[SET12]);
1160b4fbaa2aSMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1161b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1162b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
11635ef31b24SMark F. Adams 
1164b4fbaa2aSMark F. Adams     /* create timer stages */
1165b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1166b4fbaa2aSMark F. Adams     {
1167b4fbaa2aSMark F. Adams       char str[32];
1168b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1169b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1170b4fbaa2aSMark F. Adams       PetscInt lidx;
1171b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1172b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1173b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1174b4fbaa2aSMark F. Adams       }
1175b4fbaa2aSMark F. Adams     }
1176b4fbaa2aSMark F. Adams #endif
1177b4fbaa2aSMark F. Adams   }
1178b4fbaa2aSMark F. Adams #endif
11795b89ad90SMark F. Adams   PetscFunctionReturn(0);
11805b89ad90SMark F. Adams }
11815b89ad90SMark F. Adams EXTERN_C_END
1182