xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 3542efc54050ac0c2e879cbbaee517cefd02e6b1)
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;
19*3542efc5SMark F. Adams static PetscInt s_min_eq_proc = 600;
20*3542efc5SMark F. Adams static PetscReal s_threshold = 0.001;
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 */
178038e3b61SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
179737a81a9SMark F. Adams 
180038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
181038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
182038e3b61SMark F. Adams   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
183038e3b61SMark F. Adams 
184676e1743SMark F. Adams   if( s_avoid_repart ) {
18522063be5SMark F. Adams     *a_Amat_crs = Cmat; /* output */
18622063be5SMark F. Adams   }
18722063be5SMark F. Adams   else {
18822063be5SMark F. Adams     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
1895ef31b24SMark F. Adams     Mat              adj;
19022063be5SMark F. Adams     const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols;
191b4fbaa2aSMark F. Adams     const PetscInt  stride0=ncrs0*a_ndata_rows;
19292a756f0SMark F. Adams     PetscInt        is_sz,*isnewproc_idx,ii,jj,kk,strideNew,*tidx;
193c9a0b8beSMark F. Adams     /* create sub communicator  */
194c9a0b8beSMark F. Adams     MPI_Comm        cm,new_comm;
195afc97cdcSMark F. Adams     MPI_Group       wg, g2;
196fd3c6afaSMark F. Adams     PetscInt       *counts,inpe;
197fd3c6afaSMark F. Adams     PetscMPIInt    *ranks;
19822063be5SMark F. Adams     IS              isscat;
19922063be5SMark F. Adams     PetscScalar    *array;
20022063be5SMark F. Adams     Vec             src_crd, dest_crd;
20122063be5SMark F. Adams     PetscReal      *data = *a_coarse_data;
20222063be5SMark F. Adams     VecScatter      vecscat;
203b4fbaa2aSMark F. Adams     IS  isnewproc;
204e33ef3b1SMark F. Adams 
20522063be5SMark F. Adams     /* get number of PEs to make active, reduce */
20622063be5SMark F. Adams     ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
207676e1743SMark F. Adams     new_npe = neq/s_min_eq_proc; /* hardwire min. number of eq/proc */
20822063be5SMark F. Adams     if( new_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1;
20922063be5SMark F. Adams     else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */
21022063be5SMark F. Adams 
21192a756f0SMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr);
212fd3c6afaSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
21392a756f0SMark F. Adams 
214fd3c6afaSMark F. Adams     ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr);
215afc97cdcSMark F. Adams     assert(counts[mype]==ncrs0);
216afc97cdcSMark F. Adams     /* count real active pes */
21722063be5SMark F. Adams     for( nactive = jj = 0 ; jj < npe ; jj++) {
218afc97cdcSMark F. Adams       if( counts[jj] != 0 ) {
21922063be5SMark F. Adams 	ranks[nactive++] = jj;
220afc97cdcSMark F. Adams         }
221afc97cdcSMark F. Adams     }
2223303bcf9Sadams 
2233303bcf9Sadams     if (nactive < new_npe) new_npe = nactive; /* this can happen with empty input procs */
22422063be5SMark F. Adams 
22558471d46SMark F. Adams #ifdef VERBOSE
22622063be5SMark 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);
22758471d46SMark F. Adams #endif
22858471d46SMark F. Adams 
22922063be5SMark F. Adams     *a_nactive_proc = new_npe; /* output */
2302f03bc48SMark F. Adams 
231afc97cdcSMark F. Adams     ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr);
23222063be5SMark F. Adams     ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr);
233afc97cdcSMark F. Adams     ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr);
234b4fbaa2aSMark F. Adams 
235afc97cdcSMark F. Adams     if( cm != MPI_COMM_NULL ) {
236b4fbaa2aSMark F. Adams       assert(ncrs0 != 0);
237c9a0b8beSMark F. Adams       ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr);
238c9a0b8beSMark F. Adams       ierr = MPI_Comm_free( &cm );                             CHKERRQ(ierr);
239afc97cdcSMark F. Adams     }
240b4fbaa2aSMark F. Adams     else assert(ncrs0 == 0);
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;
259b4fbaa2aSMark F. Adams 
26058471d46SMark F. Adams       ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
26158471d46SMark F. Adams       for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) {
26258471d46SMark F. Adams         ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
26358471d46SMark F. Adams         d_nnz[jj] = ncols/a_cbs;
26458471d46SMark F. Adams         if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
26558471d46SMark F. Adams         ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
26658471d46SMark F. Adams       }
2676876a03eSMark F. Adams 
268eb07cef2SMark F. Adams       ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
269eb07cef2SMark F. Adams                               PETSC_DETERMINE, PETSC_DETERMINE,
27058471d46SMark F. Adams                               0, d_nnz, 0, d_nnz,
271eb07cef2SMark F. Adams                               &tMat );
2726876a03eSMark F. Adams       CHKERRQ(ierr);
27358471d46SMark F. Adams       ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
274eb07cef2SMark F. Adams 
27522063be5SMark F. Adams       for ( ii = Istart0; ii < Iend0; ii++ ) {
27622063be5SMark F. Adams         PetscInt dest_row = ii/a_cbs;
27722063be5SMark F. Adams         ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
278eb07cef2SMark F. Adams         for( jj = 0 ; jj < ncols ; jj++ ){
279038e3b61SMark F. Adams           PetscInt dest_col = idx[jj]/a_cbs;
280eb07cef2SMark F. Adams           PetscScalar v = 1.0;
281eb07cef2SMark F. Adams           ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
282eb07cef2SMark F. Adams         }
28322063be5SMark F. Adams         ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
284eb07cef2SMark F. Adams       }
285eb07cef2SMark F. Adams       ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
286eb07cef2SMark F. Adams       ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
287eb07cef2SMark F. Adams 
288b4fbaa2aSMark F. Adams       static int llev = 0;
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 
3125ef31b24SMark F. Adams       ierr = MatPartitioningCreate( new_comm, &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);
3324624a670SBarry Smith       ierr = PetscCommDestroy( &new_comm );              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);
443d61a3a7dSMark F. Adams 
44411e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
44511e60469SMark F. Adams       CHKERRQ(ierr);
446d61a3a7dSMark F. Adams 
44711e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
44811e60469SMark F. Adams     }
449d61a3a7dSMark F. Adams 
4503530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
4513530afc2SMark F. Adams     *a_P_inout = Pnew; /* output */
45292a756f0SMark F. Adams 
45311e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
45492a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
45592a756f0SMark F. Adams     ierr = PetscFree( ranks );  CHKERRQ(ierr);
456e33ef3b1SMark F. Adams   }
4575b89ad90SMark F. Adams 
4585b89ad90SMark F. Adams   PetscFunctionReturn(0);
4595b89ad90SMark F. Adams }
4605b89ad90SMark F. Adams 
4615b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4625b89ad90SMark F. Adams /*
4635b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4645b89ad90SMark F. Adams                     by setting data structures and options.
4655b89ad90SMark F. Adams 
4665b89ad90SMark F. Adams    Input Parameter:
4675b89ad90SMark F. Adams .  pc - the preconditioner context
4685b89ad90SMark F. Adams 
4695b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4705b89ad90SMark F. Adams 
4715b89ad90SMark F. Adams    Notes:
4725b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4735b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4745b89ad90SMark F. Adams */
4755b89ad90SMark F. Adams #undef __FUNCT__
4765b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
477eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
4785b89ad90SMark F. Adams {
4795b89ad90SMark F. Adams   PetscErrorCode  ierr;
480eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
4815b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
48258471d46SMark F. Adams   PC_MG_Levels   **mglevels = mg->levels;
483eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
484d3d6bff4SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
485eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
4863530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
487038e3b61SMark F. Adams   PetscBool        isOK;
488587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
489587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
490737a81a9SMark F. Adams   MatInfo          info;
4915ef31b24SMark F. Adams 
4925b89ad90SMark F. Adams   PetscFunctionBegin;
49303a628feSMark F. Adams   pc_gamg->m_count++;
49458471d46SMark F. Adams   if( a_pc->setupcalled > 0 ) {
49503a628feSMark F. Adams     /* just do Galerkin grids */
49658471d46SMark F. Adams     Mat B,dA,dB;
49758471d46SMark F. Adams 
49858471d46SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
49958471d46SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
50058471d46SMark F. Adams 
50158471d46SMark F. Adams     /* currently only handle case where mat and pmat are the same on coarser levels */
50258471d46SMark F. Adams     ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
50358471d46SMark F. Adams     /* (re)set to get dirty flag */
50458471d46SMark F. Adams     ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
50558471d46SMark F. Adams     ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr);
50658471d46SMark F. Adams 
50758471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-2; level>-1; level--) {
50858471d46SMark F. Adams       ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
50903a628feSMark F. Adams       /* the first time through the matrix structure has changed from repartitioning */
51003a628feSMark F. Adams       if( pc_gamg->m_count == 2 ) {
51103a628feSMark F. Adams         ierr = MatDestroy( &B );  CHKERRQ(ierr);
51203a628feSMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
51303a628feSMark F. Adams         mglevels[level]->A = B;
51403a628feSMark F. Adams       }
51503a628feSMark F. Adams       else {
51658471d46SMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
51703a628feSMark F. Adams       }
51858471d46SMark F. Adams       ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
51958471d46SMark F. Adams       dB = B;
52058471d46SMark F. Adams       /* setup KSP/PC */
52158471d46SMark F. Adams       ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
52258471d46SMark F. Adams     }
52358471d46SMark F. Adams 
52458471d46SMark F. Adams #define PRINT_MATS !PETSC_TRUE
52558471d46SMark F. Adams     /* plot levels - A */
52658471d46SMark F. Adams     if( PRINT_MATS ) {
52758471d46SMark F. Adams       for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){
52858471d46SMark F. Adams         PetscViewer viewer;
52958471d46SMark F. Adams         char fname[32]; KSP smoother; Mat Tmat, TTm;
53058471d46SMark F. Adams         ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
53158471d46SMark F. Adams         ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr);
53203a628feSMark F. Adams         sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
53358471d46SMark F. Adams         ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
53458471d46SMark F. Adams         ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
53558471d46SMark F. Adams         ierr = MatView( Tmat, viewer ); CHKERRQ(ierr);
53658471d46SMark F. Adams         ierr = PetscViewerDestroy( &viewer );
53758471d46SMark F. Adams       }
53858471d46SMark F. Adams     }
53958471d46SMark F. Adams 
54003a628feSMark F. Adams     a_pc->setupcalled = 2;
54103a628feSMark F. Adams 
54258471d46SMark F. Adams     PetscFunctionReturn(0);
543eb07cef2SMark F. Adams   }
544baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
545baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
5465b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
5475b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
5483530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
549eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
550eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
551eb07cef2SMark F. Adams 
552038e3b61SMark F. Adams   /* get data of not around */
553038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
554038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
555eb07cef2SMark F. Adams   }
556eb07cef2SMark F. Adams   data = pc_gamg->m_data;
557038e3b61SMark F. Adams 
558eb07cef2SMark F. Adams   /* Get A_i and R_i */
559737a81a9SMark F. Adams   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
56058471d46SMark F. Adams #ifdef VERBOSE
5612c047e2dSMark 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",
5622c047e2dSMark F. Adams 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,
5632c047e2dSMark F. Adams 	      (int)(info.nz_used/(PetscReal)N),npe);
56458471d46SMark F. Adams #endif
5658f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
566785cba28SMark F. Adams         level < GAMG_MAXLEVELS-1 && (level==0 || M>TOP_GRID_LIM) && (npe==1 || nactivepe>1);
5670205a208SMark F. Adams         level++ ){
5685b89ad90SMark F. Adams     level1 = level + 1;
569b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
570b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
571b4fbaa2aSMark F. Adams #endif
572b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
573b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
574b4fbaa2aSMark F. Adams #endif
575470e26f8SMark F. Adams     ierr = createProlongation(Aarr[level], data, pc_gamg->m_dim, pc_gamg->m_data_cols, pc_gamg->m_method,
576*3542efc5SMark F. Adams                               level, s_threshold, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level] );
5775b89ad90SMark F. Adams     CHKERRQ(ierr);
578d3d6bff4SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
579b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
580b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
581b4fbaa2aSMark F. Adams #endif
582baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
583baa4e9faSMark F. Adams     if( isOK ) {
584b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
585b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
586b4fbaa2aSMark F. Adams #endif
587470e26f8SMark F. Adams       ierr = partitionLevel( Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs,
588eb07cef2SMark F. Adams                              &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
5893530afc2SMark F. Adams       CHKERRQ(ierr);
590b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
591b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
592b4fbaa2aSMark F. Adams #endif
5933530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
594737a81a9SMark F. Adams       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
59558471d46SMark F. Adams #ifdef VERBOSE
59658471d46SMark 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",
5972c047e2dSMark F. Adams 		  mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols,
5982c047e2dSMark F. Adams 		  (int)(info.nz_used/(PetscReal)N),nactivepe);
59958471d46SMark F. Adams #endif
600e33ef3b1SMark F. Adams       /* coarse grids with SA can have zero row/cols from singleton aggregates */
60158471d46SMark F. Adams       /* aggregation method should gaurrentee this does not happen! */
602737a81a9SMark F. Adams 
60358471d46SMark F. Adams #ifdef VERBOSE
604e33ef3b1SMark F. Adams       if( PETSC_TRUE ){
605785cba28SMark F. Adams         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
606785cba28SMark F. Adams         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
607e33ef3b1SMark F. Adams         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
608785cba28SMark F. Adams         nloceq = Iend-Istart;
609e33ef3b1SMark F. Adams         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
610e33ef3b1SMark F. Adams         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
611e33ef3b1SMark F. Adams         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
612785cba28SMark F. Adams         for(kk=0;kk<nloceq;kk++){
613e33ef3b1SMark F. Adams           if(data_arr[kk]==0.0) {
614785cba28SMark F. Adams             id = kk + Istart;
615e33ef3b1SMark F. Adams             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
616e33ef3b1SMark F. Adams             CHKERRQ(ierr);
617ecd57171SMark F. Adams             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1);
618e33ef3b1SMark F. Adams           }
619e33ef3b1SMark F. Adams         }
620e33ef3b1SMark F. Adams         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
621e33ef3b1SMark F. Adams         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
622e33ef3b1SMark F. Adams         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
623e33ef3b1SMark F. Adams         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
624e33ef3b1SMark F. Adams       }
62558471d46SMark F. Adams #endif
626baa4e9faSMark F. Adams     }
627baa4e9faSMark F. Adams     else{
628be544d3cSMark F. Adams       coarse_data = 0;
629baa4e9faSMark F. Adams       break;
630baa4e9faSMark F. Adams     }
631eb07cef2SMark F. Adams     data = coarse_data;
632b4fbaa2aSMark F. Adams 
633b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
634b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
635b4fbaa2aSMark F. Adams #endif
6365b89ad90SMark F. Adams   }
637be544d3cSMark F. Adams   if( coarse_data ) {
638eb07cef2SMark F. Adams     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
639be544d3cSMark F. Adams   }
64058471d46SMark F. Adams #ifdef VERBOSE
641baa4e9faSMark F. Adams   PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
64258471d46SMark F. Adams #endif
6435b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
6445b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
6455b89ad90SMark F. Adams   fine_level = level;
646eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
6475b89ad90SMark F. Adams 
648fc4362bfSMark F. Adams   /* set default smoothers */
649587fa25dSMark F. Adams   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
650587fa25dSMark F. Adams         lidx <= fine_level;
651587fa25dSMark F. Adams         lidx++, level--) {
6525745e0f5SMark F. Adams     PetscReal emax, emin;
6535b89ad90SMark F. Adams     KSP smoother; PC subpc;
654587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
6555b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
656038e3b61SMark F. Adams     if( emaxs[level] > 0.0 ) emax=emaxs[level];
657587fa25dSMark F. Adams     else{ /* eigen estimate 'emax' */
658587fa25dSMark F. Adams       KSP eksp; Mat Lmat = Aarr[level];
659fc4362bfSMark F. Adams       Vec bb, xx; PC pc;
660038e3b61SMark F. Adams 
6615745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
6625745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
663fc4362bfSMark F. Adams       {
664fc4362bfSMark F. Adams 	PetscRandom    rctx;
665fc4362bfSMark F. Adams 	ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
666fc4362bfSMark F. Adams 	ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
667fc4362bfSMark F. Adams 	ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
668fc4362bfSMark F. Adams 	ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
6695b89ad90SMark F. Adams       }
670fc4362bfSMark F. Adams       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
671e33ef3b1SMark F. Adams       ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr);
672fc4362bfSMark F. Adams       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
67358471d46SMark F. Adams       ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
674fc4362bfSMark F. Adams       ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
67558471d46SMark F. Adams       ierr = PCSetType( pc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); /* should be same as eigen estimates op. */
676038e3b61SMark F. Adams       ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
677fc4362bfSMark F. Adams       CHKERRQ(ierr);
678fc4362bfSMark F. Adams       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
679fc4362bfSMark F. Adams 
680fc4362bfSMark F. Adams       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
681fc4362bfSMark F. Adams       ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
682fc4362bfSMark F. Adams       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
683fc4362bfSMark F. Adams       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
684fc4362bfSMark F. Adams       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
685fc4362bfSMark F. Adams       ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
68658471d46SMark F. Adams #ifdef VERBOSE
6871fddbf69SMark F. Adams       PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s max eigen=%e min=%e PC=%s\n",__FUNCT__,emax,emin,PETSC_GAMG_SMOOTHER);
68858471d46SMark F. Adams #endif
689fc4362bfSMark F. Adams     }
690038e3b61SMark F. Adams     {
691038e3b61SMark F. Adams       PetscInt N1, N0, tt;
692038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
693038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
694785cba28SMark F. Adams       emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */
695038e3b61SMark F. Adams       emax *= 1.05;
696038e3b61SMark F. Adams     }
697038e3b61SMark F. Adams 
69858471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );
699fc4362bfSMark F. Adams     ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
7000e1b4bd6SMark F. Adams     /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */
7015745e0f5SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
702e33ef3b1SMark F. Adams     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
7035745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
7045745e0f5SMark F. Adams   }
7055745e0f5SMark F. Adams   {
706ecd57171SMark F. Adams     /* coarse grid */
707ecd57171SMark F. Adams     KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
7085745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
709eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
71058471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
7115745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
712ecd57171SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
713ecd57171SMark F. Adams     ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
714ecd57171SMark F. Adams     ierr = PCSetUp( subpc ); CHKERRQ(ierr);
715ecd57171SMark F. Adams     ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
716ecd57171SMark F. Adams     assert(ii==1);
717ecd57171SMark F. Adams     ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
718ecd57171SMark F. Adams     ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
719fc4362bfSMark F. Adams   }
720737a81a9SMark F. Adams 
721fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
722eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
7235b89ad90SMark F. Adams   {
7245b89ad90SMark F. Adams     PetscBool galerkin;
725eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
7265b89ad90SMark F. Adams     if(galerkin){
7275b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
7285b89ad90SMark F. Adams     }
7295b89ad90SMark F. Adams   }
7305745e0f5SMark F. Adams 
73158471d46SMark F. Adams   /* plot levels - R/P */
73258471d46SMark F. Adams   if( PRINT_MATS ) {
73358471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-1;level>0;level--){
73458471d46SMark F. Adams       PetscViewer viewer;
73558471d46SMark F. Adams       char fname[32];
73603a628feSMark F. Adams       sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
73711e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
7385b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
739038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
7405b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
74103a628feSMark F. Adams       sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
742e33ef3b1SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
743e33ef3b1SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
744e33ef3b1SMark F. Adams       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
745e33ef3b1SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
7465b89ad90SMark F. Adams     }
74758471d46SMark F. Adams   }
74858471d46SMark F. Adams 
74958471d46SMark F. Adams   /* set interpolation between the levels, clean up */
75058471d46SMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
75158471d46SMark F. Adams        lidx<fine_level;
75258471d46SMark F. Adams        lidx++, level--){
75358471d46SMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
754587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
755587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
7565b89ad90SMark F. Adams   }
7575745e0f5SMark F. Adams 
7585b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
759eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
760eb07cef2SMark F. Adams   ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
76103a628feSMark F. Adams   a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
76258471d46SMark F. Adams 
76358471d46SMark F. Adams   {
76458471d46SMark F. Adams     KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
76558471d46SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
76658471d46SMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
76758471d46SMark F. Adams     ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
76858471d46SMark F. Adams   }
769e33ef3b1SMark F. Adams 
7705b89ad90SMark F. Adams   PetscFunctionReturn(0);
7715b89ad90SMark F. Adams }
7725b89ad90SMark F. Adams 
773eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7745b89ad90SMark F. Adams /*
7755b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7765b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7775b89ad90SMark F. Adams 
7785b89ad90SMark F. Adams    Input Parameter:
7795b89ad90SMark F. Adams .  pc - the preconditioner context
7805b89ad90SMark F. Adams 
7815b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7825b89ad90SMark F. Adams */
7835b89ad90SMark F. Adams #undef __FUNCT__
7845b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
7855b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
7865b89ad90SMark F. Adams {
7875b89ad90SMark F. Adams   PetscErrorCode  ierr;
7885b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
7895b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
7905b89ad90SMark F. Adams 
7915b89ad90SMark F. Adams   PetscFunctionBegin;
7925b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7935b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7945b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7955b89ad90SMark F. Adams   PetscFunctionReturn(0);
7965b89ad90SMark F. Adams }
7975b89ad90SMark F. Adams 
798676e1743SMark F. Adams 
799676e1743SMark F. Adams #undef __FUNCT__
800676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
801676e1743SMark F. Adams /*@
802676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
803676e1743SMark F. Adams    processor reduction.
804676e1743SMark F. Adams 
805676e1743SMark F. Adams    Not Collective on PC
806676e1743SMark F. Adams 
807676e1743SMark F. Adams    Input Parameters:
808676e1743SMark F. Adams .  pc - the preconditioner context
809676e1743SMark F. Adams 
810676e1743SMark F. Adams 
811676e1743SMark F. Adams    Options Database Key:
812676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
813676e1743SMark F. Adams 
814676e1743SMark F. Adams    Level: intermediate
815676e1743SMark F. Adams 
816676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
817676e1743SMark F. Adams 
818676e1743SMark F. Adams .seealso: ()
819676e1743SMark F. Adams @*/
820676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
821676e1743SMark F. Adams {
822676e1743SMark F. Adams   PetscErrorCode ierr;
823676e1743SMark F. Adams 
824676e1743SMark F. Adams   PetscFunctionBegin;
825676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
826676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
827676e1743SMark F. Adams   PetscFunctionReturn(0);
828676e1743SMark F. Adams }
829676e1743SMark F. Adams 
830676e1743SMark F. Adams EXTERN_C_BEGIN
831676e1743SMark F. Adams #undef __FUNCT__
832676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
833676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
834676e1743SMark F. Adams {
835676e1743SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
836676e1743SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
837676e1743SMark F. Adams 
838676e1743SMark F. Adams   PetscFunctionBegin;
839676e1743SMark F. Adams   if(n>0) s_min_eq_proc = n;
840676e1743SMark F. Adams   PetscFunctionReturn(0);
841676e1743SMark F. Adams }
842676e1743SMark F. Adams EXTERN_C_END
843676e1743SMark F. Adams 
844676e1743SMark F. Adams #undef __FUNCT__
845676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning"
846676e1743SMark F. Adams /*@
847676e1743SMark F. Adams    PCGAMGAvoidRepartitioning - Do not repartition the coarse grids
848676e1743SMark F. Adams 
849676e1743SMark F. Adams    Collective on PC
850676e1743SMark F. Adams 
851676e1743SMark F. Adams    Input Parameters:
852676e1743SMark F. Adams .  pc - the preconditioner context
853676e1743SMark F. Adams 
854676e1743SMark F. Adams 
855676e1743SMark F. Adams    Options Database Key:
856676e1743SMark F. Adams .  -pc_gamg_avoid_repartitioning
857676e1743SMark F. Adams 
858676e1743SMark F. Adams    Level: intermediate
859676e1743SMark F. Adams 
860676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
861676e1743SMark F. Adams 
862676e1743SMark F. Adams .seealso: ()
863676e1743SMark F. Adams @*/
864676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning(PC pc, PetscBool n)
865676e1743SMark F. Adams {
866676e1743SMark F. Adams   PetscErrorCode ierr;
867676e1743SMark F. Adams 
868676e1743SMark F. Adams   PetscFunctionBegin;
869676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
870676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGAvoidRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
871676e1743SMark F. Adams   PetscFunctionReturn(0);
872676e1743SMark F. Adams }
873676e1743SMark F. Adams 
874676e1743SMark F. Adams EXTERN_C_BEGIN
875676e1743SMark F. Adams #undef __FUNCT__
876676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning_GAMG"
877676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning_GAMG(PC pc, PetscBool n)
878676e1743SMark F. Adams {
879676e1743SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
880676e1743SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
881676e1743SMark F. Adams 
882676e1743SMark F. Adams   PetscFunctionBegin;
883676e1743SMark F. Adams   s_avoid_repart = n;
884676e1743SMark F. Adams   PetscFunctionReturn(0);
885676e1743SMark F. Adams }
886676e1743SMark F. Adams EXTERN_C_END
887676e1743SMark F. Adams 
888676e1743SMark F. Adams #undef __FUNCT__
889*3542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
890*3542efc5SMark F. Adams /*@
891*3542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
892*3542efc5SMark F. Adams 
893*3542efc5SMark F. Adams    Not collective on PC
894*3542efc5SMark F. Adams 
895*3542efc5SMark F. Adams    Input Parameters:
896*3542efc5SMark F. Adams .  pc - the preconditioner context
897*3542efc5SMark F. Adams 
898*3542efc5SMark F. Adams 
899*3542efc5SMark F. Adams    Options Database Key:
900*3542efc5SMark F. Adams .  -pc_gamg_threshold
901*3542efc5SMark F. Adams 
902*3542efc5SMark F. Adams    Level: intermediate
903*3542efc5SMark F. Adams 
904*3542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
905*3542efc5SMark F. Adams 
906*3542efc5SMark F. Adams .seealso: ()
907*3542efc5SMark F. Adams @*/
908*3542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
909*3542efc5SMark F. Adams {
910*3542efc5SMark F. Adams   PetscErrorCode ierr;
911*3542efc5SMark F. Adams 
912*3542efc5SMark F. Adams   PetscFunctionBegin;
913*3542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
914*3542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
915*3542efc5SMark F. Adams   PetscFunctionReturn(0);
916*3542efc5SMark F. Adams }
917*3542efc5SMark F. Adams 
918*3542efc5SMark F. Adams EXTERN_C_BEGIN
919*3542efc5SMark F. Adams #undef __FUNCT__
920*3542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
921*3542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
922*3542efc5SMark F. Adams {
923*3542efc5SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
924*3542efc5SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
925*3542efc5SMark F. Adams 
926*3542efc5SMark F. Adams   PetscFunctionBegin;
927*3542efc5SMark F. Adams   s_threshold = n;
928*3542efc5SMark F. Adams   PetscFunctionReturn(0);
929*3542efc5SMark F. Adams }
930*3542efc5SMark F. Adams EXTERN_C_END
931*3542efc5SMark F. Adams 
932*3542efc5SMark F. Adams #undef __FUNCT__
933676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType"
934676e1743SMark F. Adams /*@
935676e1743SMark F. Adams    PCGAMGSetSolverType - Set solution method.
936676e1743SMark F. Adams 
937676e1743SMark F. Adams    Collective on PC
938676e1743SMark F. Adams 
939676e1743SMark F. Adams    Input Parameters:
940676e1743SMark F. Adams .  pc - the preconditioner context
941676e1743SMark F. Adams 
942676e1743SMark F. Adams 
943676e1743SMark F. Adams    Options Database Key:
944*3542efc5SMark F. Adams .  -pc_gamg_type
945676e1743SMark F. Adams 
946676e1743SMark F. Adams    Level: intermediate
947676e1743SMark F. Adams 
948676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
949676e1743SMark F. Adams 
950676e1743SMark F. Adams .seealso: ()
951676e1743SMark F. Adams @*/
952676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz )
953676e1743SMark F. Adams {
954676e1743SMark F. Adams   PetscErrorCode ierr;
955676e1743SMark F. Adams 
956676e1743SMark F. Adams   PetscFunctionBegin;
957676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
958676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz));
959676e1743SMark F. Adams   CHKERRQ(ierr);
960676e1743SMark F. Adams   PetscFunctionReturn(0);
961676e1743SMark F. Adams }
962676e1743SMark F. Adams 
963676e1743SMark F. Adams EXTERN_C_BEGIN
964676e1743SMark F. Adams #undef __FUNCT__
965676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG"
966676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz )
967676e1743SMark F. Adams {
968676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
969676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
970676e1743SMark F. Adams 
971676e1743SMark F. Adams   PetscFunctionBegin;
972676e1743SMark F. Adams   if(sz < 64) strcpy(pc_gamg->m_type,str);
973676e1743SMark F. Adams   PetscFunctionReturn(0);
974676e1743SMark F. Adams }
975676e1743SMark F. Adams EXTERN_C_END
976676e1743SMark F. Adams 
9775b89ad90SMark F. Adams #undef __FUNCT__
9785b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
9795b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
9805b89ad90SMark F. Adams {
981676e1743SMark F. Adams   PetscErrorCode  ierr;
982676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
983676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
984676e1743SMark F. Adams   PetscBool flag;
9855b89ad90SMark F. Adams 
9865b89ad90SMark F. Adams   PetscFunctionBegin;
987676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
988676e1743SMark F. Adams   {
989676e1743SMark F. Adams     ierr = PetscOptionsString("-pc_gamg_type",
990470e26f8SMark F. Adams                               "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid (default)",
991676e1743SMark F. Adams                               "PCGAMGSetSolverType",
992676e1743SMark F. Adams                               pc_gamg->m_type,
993676e1743SMark F. Adams                               pc_gamg->m_type,
994676e1743SMark F. Adams                               64,
995676e1743SMark F. Adams                               &flag );
996676e1743SMark F. Adams     CHKERRQ(ierr);
9972a44bfbaSMark F. Adams 
998470e26f8SMark F. Adams     if (flag && strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2;
999470e26f8SMark F. Adams     else if (flag && strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1;
1000470e26f8SMark F. Adams     else pc_gamg->m_method = 0;
10012a44bfbaSMark F. Adams 
1002*3542efc5SMark F. Adams     /* common (static) variable */
1003676e1743SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_avoid_repartitioning",
1004676e1743SMark F. Adams                             "Do not repartion coarse grids (false)",
1005676e1743SMark F. Adams                             "PCGAMGAvoidRepartitioning",
1006676e1743SMark F. Adams                             s_avoid_repart,
1007676e1743SMark F. Adams                             &s_avoid_repart,
1008676e1743SMark F. Adams                             &flag);
1009676e1743SMark F. Adams     CHKERRQ(ierr);
1010676e1743SMark F. Adams 
1011*3542efc5SMark F. Adams     /* common (static) variable */
1012676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1013676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1014676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
1015676e1743SMark F. Adams                            s_min_eq_proc,
1016676e1743SMark F. Adams                            &s_min_eq_proc,
1017676e1743SMark F. Adams                            &flag );
1018676e1743SMark F. Adams     CHKERRQ(ierr);
1019*3542efc5SMark F. Adams 
1020*3542efc5SMark F. Adams     /* common (static) variable */
1021*3542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
1022*3542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
1023*3542efc5SMark F. Adams                             "PCGAMGSetThreshold",
1024*3542efc5SMark F. Adams                             s_threshold,
1025*3542efc5SMark F. Adams                             &s_threshold,
1026*3542efc5SMark F. Adams                             &flag );
1027*3542efc5SMark F. Adams     CHKERRQ(ierr);
1028676e1743SMark F. Adams   }
1029676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1030676e1743SMark F. Adams 
10315b89ad90SMark F. Adams   PetscFunctionReturn(0);
10325b89ad90SMark F. Adams }
10335b89ad90SMark F. Adams 
10345b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
10355b89ad90SMark F. Adams /*
10365b89ad90SMark F. Adams  PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG
10375b89ad90SMark F. Adams 
10385b89ad90SMark F. Adams    Input Parameter:
10395b89ad90SMark F. Adams .  pc - the preconditioner context
10405b89ad90SMark F. Adams 
10415b89ad90SMark F. Adams    Application Interface Routine: PCCreate()
10425b89ad90SMark F. Adams 
10435b89ad90SMark F. Adams   */
10445b89ad90SMark F. Adams  /* MC
10455b89ad90SMark F. Adams      PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide
10465b89ad90SMark F. Adams        fine grid discretization matrix and coordinates on the fine grid.
10475b89ad90SMark F. Adams 
10485b89ad90SMark F. Adams    Options Database Key:
10495b89ad90SMark F. Adams    Multigrid options(inherited)
10505b89ad90SMark F. Adams +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
10515b89ad90SMark F. Adams .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
10525b89ad90SMark F. Adams .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
10535b89ad90SMark F. Adams    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
10545b89ad90SMark F. Adams    GAMG options:
10555b89ad90SMark F. Adams 
10565b89ad90SMark F. Adams    Level: intermediate
10575b89ad90SMark F. Adams   Concepts: multigrid
10585b89ad90SMark F. Adams 
10595b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
10605b89ad90SMark F. Adams            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
10615b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
10625b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
10635b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
10645b89ad90SMark F. Adams M */
10655b89ad90SMark F. Adams 
10665b89ad90SMark F. Adams EXTERN_C_BEGIN
10675b89ad90SMark F. Adams #undef __FUNCT__
10685b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
10695b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
10705b89ad90SMark F. Adams {
10715b89ad90SMark F. Adams   PetscErrorCode  ierr;
10725b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
10735b89ad90SMark F. Adams   PC_MG           *mg;
10745ef31b24SMark F. Adams   PetscClassId     cookie;
10755b89ad90SMark F. Adams 
10765b89ad90SMark F. Adams   PetscFunctionBegin;
10775b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
10785b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
10795b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
10805b89ad90SMark F. Adams 
10815b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
10825b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
108303a628feSMark F. Adams   pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0;
10845b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
10855b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
10865b89ad90SMark F. Adams 
10875b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
10885b89ad90SMark F. Adams 
10895b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
10905b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
10915b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
10925b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
10935b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
10945b89ad90SMark F. Adams 
10955b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
10965b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
10975b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
1098c97e1df0SMark F. Adams 					    PCSetCoordinates_GAMG);
1099c97e1df0SMark F. Adams   CHKERRQ(ierr);
1100c97e1df0SMark F. Adams 
1101676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1102676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1103676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1104676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1105676e1743SMark F. Adams   CHKERRQ(ierr);
1106676e1743SMark F. Adams 
1107676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1108676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_C",
1109676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_GAMG",
1110676e1743SMark F. Adams 					    PCGAMGAvoidRepartitioning_GAMG);
1111676e1743SMark F. Adams   CHKERRQ(ierr);
1112676e1743SMark F. Adams 
1113676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1114*3542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
1115*3542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
1116*3542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
1117*3542efc5SMark F. Adams   CHKERRQ(ierr);
1118*3542efc5SMark F. Adams 
1119*3542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1120676e1743SMark F. Adams 					    "PCGAMGSetSolverType_C",
1121676e1743SMark F. Adams 					    "PCGAMGSetSolverType_GAMG",
1122676e1743SMark F. Adams 					    PCGAMGSetSolverType_GAMG);
1123676e1743SMark F. Adams   CHKERRQ(ierr);
1124c97e1df0SMark F. Adams 
1125b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
1126785cba28SMark F. Adams   static int count = 0;
1127785cba28SMark F. Adams   if( count++ == 0 ) {
11285ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
1129b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
11302a44bfbaSMark F. Adams     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
11312a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]);
11322a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]);
11332a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]);
1134b4fbaa2aSMark F. Adams     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1135b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1136b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1137b4fbaa2aSMark F. Adams     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1138b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_events[SET7]);
1139b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]); */
1140b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1141b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1142b4fbaa2aSMark F. Adams     PetscLogEventRegister("  PL repartition", cookie, &gamg_setup_events[SET12]);
1143b4fbaa2aSMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1144b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1145b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
11465ef31b24SMark F. Adams 
1147b4fbaa2aSMark F. Adams     /* create timer stages */
1148b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1149b4fbaa2aSMark F. Adams     {
1150b4fbaa2aSMark F. Adams       char str[32];
1151b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1152b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1153b4fbaa2aSMark F. Adams       PetscInt lidx;
1154b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1155b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1156b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1157b4fbaa2aSMark F. Adams       }
1158b4fbaa2aSMark F. Adams     }
1159b4fbaa2aSMark F. Adams #endif
1160b4fbaa2aSMark F. Adams   }
1161b4fbaa2aSMark F. Adams #endif
11625b89ad90SMark F. Adams   PetscFunctionReturn(0);
11635b89ad90SMark F. Adams }
11645b89ad90SMark F. Adams EXTERN_C_END
1165