xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision d8c859f08f8da0a4ca0e6817f896afec0dc432fc)
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   // #define USE_R
17996434bc1SMark F. Adams #ifdef USE_R
18096434bc1SMark F. Adams   /* make R wih brute force for now */
18196434bc1SMark F. Adams   ierr = MatTranspose( Pold, Pnew );
18296434bc1SMark F. Adams   ierr = MatDestroy( &Pold );  CHKERRQ(ierr);
18396434bc1SMark F. Adams   ierr = MatRARt( a_Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
18496434bc1SMark F. Adams   Pold = Pnew;
18596434bc1SMark F. Adams #else
186038e3b61SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
18796434bc1SMark F. Adams #endif
188038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
189038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
190038e3b61SMark F. Adams   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
191038e3b61SMark F. Adams 
192676e1743SMark F. Adams   if( s_avoid_repart ) {
19322063be5SMark F. Adams     *a_Amat_crs = Cmat; /* output */
19422063be5SMark F. Adams   }
19522063be5SMark F. Adams   else {
19622063be5SMark F. Adams     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
1975ef31b24SMark F. Adams     Mat              adj;
19822063be5SMark F. Adams     const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols;
199b4fbaa2aSMark F. Adams     const PetscInt  stride0=ncrs0*a_ndata_rows;
20092a756f0SMark F. Adams     PetscInt        is_sz,*isnewproc_idx,ii,jj,kk,strideNew,*tidx;
201c9a0b8beSMark F. Adams     /* create sub communicator  */
202c9a0b8beSMark F. Adams     MPI_Comm        cm,new_comm;
203afc97cdcSMark F. Adams     MPI_Group       wg, g2;
204fd3c6afaSMark F. Adams     PetscInt       *counts,inpe;
205fd3c6afaSMark F. Adams     PetscMPIInt    *ranks;
20622063be5SMark F. Adams     IS              isscat;
20722063be5SMark F. Adams     PetscScalar    *array;
20822063be5SMark F. Adams     Vec             src_crd, dest_crd;
20922063be5SMark F. Adams     PetscReal      *data = *a_coarse_data;
21022063be5SMark F. Adams     VecScatter      vecscat;
211b4fbaa2aSMark F. Adams     IS  isnewproc;
212e33ef3b1SMark F. Adams 
21322063be5SMark F. Adams     /* get number of PEs to make active, reduce */
21422063be5SMark F. Adams     ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
215676e1743SMark F. Adams     new_npe = neq/s_min_eq_proc; /* hardwire min. number of eq/proc */
21622063be5SMark F. Adams     if( new_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1;
21722063be5SMark F. Adams     else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */
21822063be5SMark F. Adams 
21992a756f0SMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr);
220fd3c6afaSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
22192a756f0SMark F. Adams 
222fd3c6afaSMark F. Adams     ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr);
223afc97cdcSMark F. Adams     assert(counts[mype]==ncrs0);
224afc97cdcSMark F. Adams     /* count real active pes */
22522063be5SMark F. Adams     for( nactive = jj = 0 ; jj < npe ; jj++) {
226afc97cdcSMark F. Adams       if( counts[jj] != 0 ) {
22722063be5SMark F. Adams 	ranks[nactive++] = jj;
228afc97cdcSMark F. Adams         }
229afc97cdcSMark F. Adams     }
2303303bcf9Sadams 
2313303bcf9Sadams     if (nactive < new_npe) new_npe = nactive; /* this can happen with empty input procs */
23222063be5SMark F. Adams 
23358471d46SMark F. Adams #ifdef VERBOSE
23422063be5SMark 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);
23558471d46SMark F. Adams #endif
23658471d46SMark F. Adams 
23722063be5SMark F. Adams     *a_nactive_proc = new_npe; /* output */
2382f03bc48SMark F. Adams 
239afc97cdcSMark F. Adams     ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr);
24022063be5SMark F. Adams     ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr);
241afc97cdcSMark F. Adams     ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr);
242b4fbaa2aSMark F. Adams 
243afc97cdcSMark F. Adams     if( cm != MPI_COMM_NULL ) {
244b4fbaa2aSMark F. Adams       assert(ncrs0 != 0);
245c9a0b8beSMark F. Adams       ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr);
246c9a0b8beSMark F. Adams       ierr = MPI_Comm_free( &cm );                             CHKERRQ(ierr);
247afc97cdcSMark F. Adams     }
248b4fbaa2aSMark F. Adams     else assert(ncrs0 == 0);
249b4fbaa2aSMark F. Adams 
250afc97cdcSMark F. Adams     ierr = MPI_Group_free( &wg );                            CHKERRQ(ierr);
251afc97cdcSMark F. Adams     ierr = MPI_Group_free( &g2 );                            CHKERRQ(ierr);
252c9a0b8beSMark F. Adams 
2535ef31b24SMark F. Adams     /* MatPartitioningApply call MatConvert, which is collective */
254b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
255b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
256b4fbaa2aSMark F. Adams #endif
257038e3b61SMark F. Adams     if( a_cbs == 1) {
258038e3b61SMark F. Adams       ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
259eb07cef2SMark F. Adams     }
260eb07cef2SMark F. Adams     else{
261038e3b61SMark F. Adams       /* make a scalar matrix to partition */
262eb07cef2SMark F. Adams       Mat tMat;
26358471d46SMark F. Adams       PetscInt ncols,jj,Ii;
264b4fbaa2aSMark F. Adams       const PetscScalar *vals;
265b4fbaa2aSMark F. Adams       const PetscInt *idx;
26658471d46SMark F. Adams       PetscInt *d_nnz;
2675ee9c036SSatish Balay       static int llev = 0;
268b4fbaa2aSMark F. Adams 
26958471d46SMark F. Adams       ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
27058471d46SMark F. Adams       for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) {
27158471d46SMark F. Adams         ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
27258471d46SMark F. Adams         d_nnz[jj] = ncols/a_cbs;
27358471d46SMark F. Adams         if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
27458471d46SMark F. Adams         ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
27558471d46SMark F. Adams       }
2766876a03eSMark F. Adams 
277eb07cef2SMark F. Adams       ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
278eb07cef2SMark F. Adams                               PETSC_DETERMINE, PETSC_DETERMINE,
27958471d46SMark F. Adams                               0, d_nnz, 0, d_nnz,
280eb07cef2SMark F. Adams                               &tMat );
2816876a03eSMark F. Adams       CHKERRQ(ierr);
28258471d46SMark F. Adams       ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
283eb07cef2SMark F. Adams 
28422063be5SMark F. Adams       for ( ii = Istart0; ii < Iend0; ii++ ) {
28522063be5SMark F. Adams         PetscInt dest_row = ii/a_cbs;
28622063be5SMark F. Adams         ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
287eb07cef2SMark F. Adams         for( jj = 0 ; jj < ncols ; jj++ ){
288038e3b61SMark F. Adams           PetscInt dest_col = idx[jj]/a_cbs;
289eb07cef2SMark F. Adams           PetscScalar v = 1.0;
290eb07cef2SMark F. Adams           ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
291eb07cef2SMark F. Adams         }
29222063be5SMark F. Adams         ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
293eb07cef2SMark F. Adams       }
294eb07cef2SMark F. Adams       ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
295eb07cef2SMark F. Adams       ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
296eb07cef2SMark F. Adams 
297b4fbaa2aSMark F. Adams       if( llev++ == -1 ) {
298b4fbaa2aSMark F. Adams 	PetscViewer viewer; char fname[32];
299b4fbaa2aSMark F. Adams 	sprintf(fname,"part_mat_%d.mat",llev);
300b4fbaa2aSMark F. Adams 	PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
301b4fbaa2aSMark F. Adams 	ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
302b4fbaa2aSMark F. Adams 	ierr = PetscViewerDestroy( &viewer );
303b4fbaa2aSMark F. Adams       }
304b4fbaa2aSMark F. Adams 
305eb07cef2SMark F. Adams       ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
306eb07cef2SMark F. Adams 
307eb07cef2SMark F. Adams       ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
308eb07cef2SMark F. Adams     }
309afc97cdcSMark F. Adams     if( ncrs0 != 0 ){
310b4fbaa2aSMark F. Adams       const PetscInt *is_idx;
311b4fbaa2aSMark F. Adams       MatPartitioning  mpart;
3125ef31b24SMark F. Adams       /* hack to fix global data that pmetis.c uses in 'adj' */
31322063be5SMark F. Adams       for( nactive = jj = 0 ; jj < npe ; jj++ ) {
314afc97cdcSMark F. Adams 	if( counts[jj] != 0 ) {
31522063be5SMark F. Adams 	  adj->rmap->range[nactive++] = adj->rmap->range[jj];
316afc97cdcSMark F. Adams 	}
3175ef31b24SMark F. Adams       }
31822063be5SMark F. Adams       adj->rmap->range[nactive] = adj->rmap->range[npe];
3192f03bc48SMark F. Adams 
3205ef31b24SMark F. Adams       ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr);
3215ef31b24SMark F. Adams       ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
32211e60469SMark F. Adams       ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
3235ef31b24SMark F. Adams       ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
32411e60469SMark F. Adams       ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr);
32511e60469SMark F. Adams       ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
32622063be5SMark F. Adams 
3275ef31b24SMark F. Adams       /* collect IS info */
3285ef31b24SMark F. Adams       ierr = ISGetLocalSize( isnewproc, &is_sz );       CHKERRQ(ierr);
329038e3b61SMark F. Adams       ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr);
3305ef31b24SMark F. Adams       ierr = ISGetIndices( isnewproc, &is_idx );        CHKERRQ(ierr);
331b4fbaa2aSMark F. Adams       /* spread partitioning across machine - best way ??? */
332ecd57171SMark F. Adams       NN = 1; /*npe/new_npe;*/
333b4fbaa2aSMark F. Adams       for( kk = jj = 0 ; kk < is_sz ; kk++ ){
334afc97cdcSMark F. Adams         for( ii = 0 ; ii < a_cbs ; ii++, jj++ ) {
33522063be5SMark F. Adams           isnewproc_idx[jj] = is_idx[kk] * NN; /* distribution */
336eb07cef2SMark F. Adams         }
3375ef31b24SMark F. Adams       }
3385ef31b24SMark F. Adams       ierr = ISRestoreIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
3395ef31b24SMark F. Adams       ierr = ISDestroy( &isnewproc );                    CHKERRQ(ierr);
3404624a670SBarry Smith       ierr = PetscCommDestroy( &new_comm );              CHKERRQ(ierr);
341b4fbaa2aSMark F. Adams 
342b4fbaa2aSMark F. Adams       is_sz *= a_cbs;
3435ef31b24SMark F. Adams     }
3445ef31b24SMark F. Adams     else{
3455ef31b24SMark F. Adams       isnewproc_idx = 0;
3465ef31b24SMark F. Adams       is_sz = 0;
3475ef31b24SMark F. Adams     }
3485ef31b24SMark F. Adams     ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
3495ef31b24SMark F. Adams     ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc );
350afc97cdcSMark F. Adams     if( isnewproc_idx != 0 ) {
3515ef31b24SMark F. Adams       ierr = PetscFree( isnewproc_idx );  CHKERRQ(ierr);
3525ef31b24SMark F. Adams     }
353e33ef3b1SMark F. Adams 
35411e60469SMark F. Adams     /*
35511e60469SMark F. Adams      Create an index set from the isnewproc index set to indicate the mapping TO
35611e60469SMark F. Adams      */
35711e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
35811e60469SMark F. Adams     /*
35911e60469SMark F. Adams      Determine how many elements are assigned to each processor
36011e60469SMark F. Adams      */
361fd3c6afaSMark F. Adams     inpe = npe;
362fd3c6afaSMark F. Adams     ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr);
36311e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
364038e3b61SMark F. Adams     ncrs_new = counts[mype]/a_cbs;
36522063be5SMark F. Adams     strideNew = ncrs_new*a_ndata_rows;
366b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
367b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
368b4fbaa2aSMark F. Adams #endif
36922063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
37011e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
371d3d6bff4SMark F. Adams     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
372470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
37311e60469SMark F. Adams     /*
374d3d6bff4SMark F. Adams      There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
375d3d6bff4SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
37611e60469SMark F. Adams      */
37792a756f0SMark F. Adams     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
37811e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
379038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
380d3d6bff4SMark F. Adams       PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
381d3d6bff4SMark F. Adams       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
38211e60469SMark F. Adams     }
383038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
384d3d6bff4SMark F. Adams     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
38511e60469SMark F. Adams     CHKERRQ(ierr);
38692a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
38711e60469SMark F. Adams     /*
38811e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
38911e60469SMark F. Adams      */
390d3d6bff4SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
391d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
392d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs0 ; ii++) {
393d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
394d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
395676e1743SMark F. Adams           PetscScalar tt = (PetscScalar)data[ix];
396676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
397d3d6bff4SMark F. Adams 	}
398038e3b61SMark F. Adams       }
399eb07cef2SMark F. Adams     }
400eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
401eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
40211e60469SMark F. Adams     /*
40311e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
40411e60469SMark F. Adams       to the correct processor
40511e60469SMark F. Adams     */
40611e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
40711e60469SMark F. Adams     CHKERRQ(ierr);
40811e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
40911e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
41011e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
41111e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
41211e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
41311e60469SMark F. Adams     /*
41411e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
41511e60469SMark F. Adams     */
416eb07cef2SMark F. Adams     ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
417d3d6bff4SMark F. Adams     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
418eb07cef2SMark F. Adams 
41911e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
420eb07cef2SMark F. Adams     data = *a_coarse_data;
421d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
422d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs_new ; ii++) {
423d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
424d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
425d3d6bff4SMark F. Adams 	  data[ix] = PetscRealPart(array[jx]);
426d3d6bff4SMark F. Adams 	  array[jx] = 1.e300;
427d3d6bff4SMark F. Adams 	}
428038e3b61SMark F. Adams       }
429038e3b61SMark F. Adams     }
43011e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
43111e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
43211e60469SMark F. Adams     /*
43311e60469SMark F. Adams       Invert for MatGetSubMatrix
43411e60469SMark F. Adams     */
435038e3b61SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
43611e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
43711e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
43811e60469SMark F. Adams     /* A_crs output */
439038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
44011e60469SMark F. Adams     CHKERRQ(ierr);
441eb07cef2SMark F. Adams 
442038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
443e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
444038e3b61SMark F. Adams     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
445eb07cef2SMark F. Adams 
44611e60469SMark F. Adams     /* prolongator */
44711e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
44811e60469SMark F. Adams     {
44911e60469SMark F. Adams       IS findices;
45011e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
45196434bc1SMark F. Adams #ifdef USE_R
45296434bc1SMark F. Adams       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
45396434bc1SMark F. Adams #else
45411e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
45596434bc1SMark F. Adams #endif
45611e60469SMark F. Adams       CHKERRQ(ierr);
457d61a3a7dSMark F. Adams 
45811e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
45911e60469SMark F. Adams     }
460d61a3a7dSMark F. Adams 
4613530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
4623530afc2SMark F. Adams     *a_P_inout = Pnew; /* output */
46392a756f0SMark F. Adams 
46411e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
46592a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
46692a756f0SMark F. Adams     ierr = PetscFree( ranks );  CHKERRQ(ierr);
467e33ef3b1SMark F. Adams   }
4685b89ad90SMark F. Adams 
4695b89ad90SMark F. Adams   PetscFunctionReturn(0);
4705b89ad90SMark F. Adams }
4715b89ad90SMark F. Adams 
4725b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4735b89ad90SMark F. Adams /*
4745b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4755b89ad90SMark F. Adams                     by setting data structures and options.
4765b89ad90SMark F. Adams 
4775b89ad90SMark F. Adams    Input Parameter:
4785b89ad90SMark F. Adams .  pc - the preconditioner context
4795b89ad90SMark F. Adams 
4805b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4815b89ad90SMark F. Adams 
4825b89ad90SMark F. Adams    Notes:
4835b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4845b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4855b89ad90SMark F. Adams */
4865b89ad90SMark F. Adams #undef __FUNCT__
4875b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
488eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
4895b89ad90SMark F. Adams {
4905b89ad90SMark F. Adams   PetscErrorCode  ierr;
491eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
4925b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
49358471d46SMark F. Adams   PC_MG_Levels   **mglevels = mg->levels;
494eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
495d3d6bff4SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
496eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
4973530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
498038e3b61SMark F. Adams   PetscBool        isOK;
499587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
500587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
501737a81a9SMark F. Adams   MatInfo          info;
5025ef31b24SMark F. Adams 
5035b89ad90SMark F. Adams   PetscFunctionBegin;
50403a628feSMark F. Adams   pc_gamg->m_count++;
50558471d46SMark F. Adams   if( a_pc->setupcalled > 0 ) {
50603a628feSMark F. Adams     /* just do Galerkin grids */
50758471d46SMark F. Adams     Mat B,dA,dB;
50858471d46SMark F. Adams 
50958471d46SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
51058471d46SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
51158471d46SMark F. Adams 
51258471d46SMark F. Adams     /* currently only handle case where mat and pmat are the same on coarser levels */
51358471d46SMark F. Adams     ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
51458471d46SMark F. Adams     /* (re)set to get dirty flag */
51558471d46SMark F. Adams     ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
51658471d46SMark F. Adams     ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr);
51758471d46SMark F. Adams 
51858471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-2; level>-1; level--) {
51958471d46SMark F. Adams       ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
52003a628feSMark F. Adams       /* the first time through the matrix structure has changed from repartitioning */
52103a628feSMark F. Adams       if( pc_gamg->m_count == 2 ) {
52203a628feSMark F. Adams         ierr = MatDestroy( &B );  CHKERRQ(ierr);
52303a628feSMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
52403a628feSMark F. Adams         mglevels[level]->A = B;
52503a628feSMark F. Adams       }
52603a628feSMark F. Adams       else {
52758471d46SMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
52803a628feSMark F. Adams       }
52958471d46SMark F. Adams       ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
53058471d46SMark F. Adams       dB = B;
53158471d46SMark F. Adams       /* setup KSP/PC */
53258471d46SMark F. Adams       ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
53358471d46SMark F. Adams     }
53458471d46SMark F. Adams 
53558471d46SMark F. Adams #define PRINT_MATS !PETSC_TRUE
53658471d46SMark F. Adams     /* plot levels - A */
53758471d46SMark F. Adams     if( PRINT_MATS ) {
53858471d46SMark F. Adams       for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){
53958471d46SMark F. Adams         PetscViewer viewer;
54058471d46SMark F. Adams         char fname[32]; KSP smoother; Mat Tmat, TTm;
54158471d46SMark F. Adams         ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
54258471d46SMark F. Adams         ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr);
54303a628feSMark F. Adams         sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
54458471d46SMark F. Adams         ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
54558471d46SMark F. Adams         ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
54658471d46SMark F. Adams         ierr = MatView( Tmat, viewer ); CHKERRQ(ierr);
54758471d46SMark F. Adams         ierr = PetscViewerDestroy( &viewer );
54858471d46SMark F. Adams       }
54958471d46SMark F. Adams     }
55058471d46SMark F. Adams 
55103a628feSMark F. Adams     a_pc->setupcalled = 2;
55203a628feSMark F. Adams 
55358471d46SMark F. Adams     PetscFunctionReturn(0);
554eb07cef2SMark F. Adams   }
555baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
556baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
5575b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
5585b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
5593530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
560eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
561eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
562eb07cef2SMark F. Adams 
563038e3b61SMark F. Adams   /* get data of not around */
564038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
565038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
566eb07cef2SMark F. Adams   }
567eb07cef2SMark F. Adams   data = pc_gamg->m_data;
568038e3b61SMark F. Adams 
569eb07cef2SMark F. Adams   /* Get A_i and R_i */
570737a81a9SMark F. Adams   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
57158471d46SMark F. Adams #ifdef VERBOSE
5722c047e2dSMark 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",
5732c047e2dSMark F. Adams 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,
5742c047e2dSMark F. Adams 	      (int)(info.nz_used/(PetscReal)N),npe);
57558471d46SMark F. Adams #endif
5768f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
577785cba28SMark F. Adams         level < GAMG_MAXLEVELS-1 && (level==0 || M>TOP_GRID_LIM) && (npe==1 || nactivepe>1);
5780205a208SMark F. Adams         level++ ){
5795b89ad90SMark F. Adams     level1 = level + 1;
580b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
581b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
582b4fbaa2aSMark F. Adams #endif
583b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
584b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
585b4fbaa2aSMark F. Adams #endif
586470e26f8SMark F. Adams     ierr = createProlongation(Aarr[level], data, pc_gamg->m_dim, pc_gamg->m_data_cols, pc_gamg->m_method,
5873542efc5SMark F. Adams                               level, s_threshold, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level] );
5885b89ad90SMark F. Adams     CHKERRQ(ierr);
589d3d6bff4SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
590b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
591b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
592b4fbaa2aSMark F. Adams #endif
593baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
594baa4e9faSMark F. Adams     if( isOK ) {
595b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
596b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
597b4fbaa2aSMark F. Adams #endif
598470e26f8SMark F. Adams       ierr = partitionLevel( Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs,
599eb07cef2SMark F. Adams                              &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
6003530afc2SMark F. Adams       CHKERRQ(ierr);
601b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
602b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
603b4fbaa2aSMark F. Adams #endif
6043530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
605737a81a9SMark F. Adams       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
60658471d46SMark F. Adams #ifdef VERBOSE
60758471d46SMark 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",
6082c047e2dSMark F. Adams 		  mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols,
6092c047e2dSMark F. Adams 		  (int)(info.nz_used/(PetscReal)N),nactivepe);
61058471d46SMark F. Adams #endif
611e33ef3b1SMark F. Adams       /* coarse grids with SA can have zero row/cols from singleton aggregates */
61258471d46SMark F. Adams       /* aggregation method should gaurrentee this does not happen! */
613737a81a9SMark F. Adams 
61458471d46SMark F. Adams #ifdef VERBOSE
615e33ef3b1SMark F. Adams       if( PETSC_TRUE ){
616785cba28SMark F. Adams         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
617785cba28SMark F. Adams         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
618e33ef3b1SMark F. Adams         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
619785cba28SMark F. Adams         nloceq = Iend-Istart;
620e33ef3b1SMark F. Adams         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
621e33ef3b1SMark F. Adams         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
622e33ef3b1SMark F. Adams         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
623785cba28SMark F. Adams         for(kk=0;kk<nloceq;kk++){
624e33ef3b1SMark F. Adams           if(data_arr[kk]==0.0) {
625785cba28SMark F. Adams             id = kk + Istart;
626e33ef3b1SMark F. Adams             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
627e33ef3b1SMark F. Adams             CHKERRQ(ierr);
628ecd57171SMark F. Adams             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1);
629e33ef3b1SMark F. Adams           }
630e33ef3b1SMark F. Adams         }
631e33ef3b1SMark F. Adams         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
632e33ef3b1SMark F. Adams         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
633e33ef3b1SMark F. Adams         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
634e33ef3b1SMark F. Adams         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
635e33ef3b1SMark F. Adams       }
63658471d46SMark F. Adams #endif
637baa4e9faSMark F. Adams     }
638baa4e9faSMark F. Adams     else{
639be544d3cSMark F. Adams       coarse_data = 0;
640baa4e9faSMark F. Adams       break;
641baa4e9faSMark F. Adams     }
642eb07cef2SMark F. Adams     data = coarse_data;
643b4fbaa2aSMark F. Adams 
644b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
645b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
646b4fbaa2aSMark F. Adams #endif
6475b89ad90SMark F. Adams   }
648be544d3cSMark F. Adams   if( coarse_data ) {
649eb07cef2SMark F. Adams     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
650be544d3cSMark F. Adams   }
65158471d46SMark F. Adams #ifdef VERBOSE
652baa4e9faSMark F. Adams   PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
65358471d46SMark F. Adams #endif
6545b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
6555b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
6565b89ad90SMark F. Adams   fine_level = level;
657eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
6585b89ad90SMark F. Adams 
659fc4362bfSMark F. Adams   /* set default smoothers */
660587fa25dSMark F. Adams   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
661587fa25dSMark F. Adams         lidx <= fine_level;
662587fa25dSMark F. Adams         lidx++, level--) {
6635745e0f5SMark F. Adams     PetscReal emax, emin;
6645b89ad90SMark F. Adams     KSP smoother; PC subpc;
665587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
6665b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
667038e3b61SMark F. Adams     if( emaxs[level] > 0.0 ) emax=emaxs[level];
668587fa25dSMark F. Adams     else{ /* eigen estimate 'emax' */
669587fa25dSMark F. Adams       KSP eksp; Mat Lmat = Aarr[level];
670fc4362bfSMark F. Adams       Vec bb, xx; PC pc;
671038e3b61SMark F. Adams 
6725745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
6735745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
674fc4362bfSMark F. Adams       {
675fc4362bfSMark F. Adams 	PetscRandom    rctx;
676fc4362bfSMark F. Adams 	ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
677fc4362bfSMark F. Adams 	ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
678fc4362bfSMark F. Adams 	ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
679fc4362bfSMark F. Adams 	ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
6805b89ad90SMark F. Adams       }
681fc4362bfSMark F. Adams       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
682e33ef3b1SMark F. Adams       ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr);
683fc4362bfSMark F. Adams       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
68458471d46SMark F. Adams       ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
685fc4362bfSMark F. Adams       ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
68658471d46SMark F. Adams       ierr = PCSetType( pc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); /* should be same as eigen estimates op. */
687038e3b61SMark F. Adams       ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
688fc4362bfSMark F. Adams       CHKERRQ(ierr);
689fc4362bfSMark F. Adams       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
690fc4362bfSMark F. Adams 
691fc4362bfSMark F. Adams       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
692fc4362bfSMark F. Adams       ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
693fc4362bfSMark F. Adams       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
694fc4362bfSMark F. Adams       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
695fc4362bfSMark F. Adams       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
696fc4362bfSMark F. Adams       ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
69758471d46SMark F. Adams #ifdef VERBOSE
698*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);
69958471d46SMark F. Adams #endif
700fc4362bfSMark F. Adams     }
701038e3b61SMark F. Adams     {
702038e3b61SMark F. Adams       PetscInt N1, N0, tt;
703038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
704038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
705785cba28SMark F. Adams       emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */
706038e3b61SMark F. Adams       emax *= 1.05;
707038e3b61SMark F. Adams     }
708038e3b61SMark F. Adams 
70958471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );
710fc4362bfSMark F. Adams     ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
7110e1b4bd6SMark F. Adams     /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */
7125745e0f5SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
713e33ef3b1SMark F. Adams     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
7145745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
7155745e0f5SMark F. Adams   }
7165745e0f5SMark F. Adams   {
717ecd57171SMark F. Adams     /* coarse grid */
718ecd57171SMark F. Adams     KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
7195745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
720eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
72158471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
7225745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
723ecd57171SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
724ecd57171SMark F. Adams     ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
725ecd57171SMark F. Adams     ierr = PCSetUp( subpc ); CHKERRQ(ierr);
726ecd57171SMark F. Adams     ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
727ecd57171SMark F. Adams     assert(ii==1);
728ecd57171SMark F. Adams     ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
729ecd57171SMark F. Adams     ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
730fc4362bfSMark F. Adams   }
731737a81a9SMark F. Adams 
732fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
733eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
7345b89ad90SMark F. Adams   {
7355b89ad90SMark F. Adams     PetscBool galerkin;
736eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
7375b89ad90SMark F. Adams     if(galerkin){
7385b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
7395b89ad90SMark F. Adams     }
7405b89ad90SMark F. Adams   }
7415745e0f5SMark F. Adams 
74258471d46SMark F. Adams   /* plot levels - R/P */
74358471d46SMark F. Adams   if( PRINT_MATS ) {
74458471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-1;level>0;level--){
74558471d46SMark F. Adams       PetscViewer viewer;
74658471d46SMark F. Adams       char fname[32];
74703a628feSMark F. Adams       sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
74811e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
7495b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
750038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
7515b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
75203a628feSMark F. Adams       sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
753e33ef3b1SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
754e33ef3b1SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
755e33ef3b1SMark F. Adams       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
756e33ef3b1SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
7575b89ad90SMark F. Adams     }
75858471d46SMark F. Adams   }
75958471d46SMark F. Adams 
76058471d46SMark F. Adams   /* set interpolation between the levels, clean up */
76158471d46SMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
76258471d46SMark F. Adams        lidx<fine_level;
76358471d46SMark F. Adams        lidx++, level--){
76458471d46SMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
765587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
766587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
7675b89ad90SMark F. Adams   }
7685745e0f5SMark F. Adams 
7695b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
770eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
771eb07cef2SMark F. Adams   ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
77203a628feSMark F. Adams   a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
77358471d46SMark F. Adams 
77458471d46SMark F. Adams   {
77558471d46SMark F. Adams     KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
77658471d46SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
77758471d46SMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
77858471d46SMark F. Adams     ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
77958471d46SMark F. Adams   }
780e33ef3b1SMark F. Adams 
7815b89ad90SMark F. Adams   PetscFunctionReturn(0);
7825b89ad90SMark F. Adams }
7835b89ad90SMark F. Adams 
784eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7855b89ad90SMark F. Adams /*
7865b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7875b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7885b89ad90SMark F. Adams 
7895b89ad90SMark F. Adams    Input Parameter:
7905b89ad90SMark F. Adams .  pc - the preconditioner context
7915b89ad90SMark F. Adams 
7925b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7935b89ad90SMark F. Adams */
7945b89ad90SMark F. Adams #undef __FUNCT__
7955b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
7965b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
7975b89ad90SMark F. Adams {
7985b89ad90SMark F. Adams   PetscErrorCode  ierr;
7995b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
8005b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
8015b89ad90SMark F. Adams 
8025b89ad90SMark F. Adams   PetscFunctionBegin;
8035b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8045b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8055b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8065b89ad90SMark F. Adams   PetscFunctionReturn(0);
8075b89ad90SMark F. Adams }
8085b89ad90SMark F. Adams 
809676e1743SMark F. Adams 
810676e1743SMark F. Adams #undef __FUNCT__
811676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
812676e1743SMark F. Adams /*@
813676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
814676e1743SMark F. Adams    processor reduction.
815676e1743SMark F. Adams 
816676e1743SMark F. Adams    Not Collective on PC
817676e1743SMark F. Adams 
818676e1743SMark F. Adams    Input Parameters:
819676e1743SMark F. Adams .  pc - the preconditioner context
820676e1743SMark F. Adams 
821676e1743SMark F. Adams 
822676e1743SMark F. Adams    Options Database Key:
823676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
824676e1743SMark F. Adams 
825676e1743SMark F. Adams    Level: intermediate
826676e1743SMark F. Adams 
827676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
828676e1743SMark F. Adams 
829676e1743SMark F. Adams .seealso: ()
830676e1743SMark F. Adams @*/
831676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
832676e1743SMark F. Adams {
833676e1743SMark F. Adams   PetscErrorCode ierr;
834676e1743SMark F. Adams 
835676e1743SMark F. Adams   PetscFunctionBegin;
836676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
837676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
838676e1743SMark F. Adams   PetscFunctionReturn(0);
839676e1743SMark F. Adams }
840676e1743SMark F. Adams 
841676e1743SMark F. Adams EXTERN_C_BEGIN
842676e1743SMark F. Adams #undef __FUNCT__
843676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
844676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
845676e1743SMark F. Adams {
846676e1743SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
847676e1743SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
848676e1743SMark F. Adams 
849676e1743SMark F. Adams   PetscFunctionBegin;
850676e1743SMark F. Adams   if(n>0) s_min_eq_proc = n;
851676e1743SMark F. Adams   PetscFunctionReturn(0);
852676e1743SMark F. Adams }
853676e1743SMark F. Adams EXTERN_C_END
854676e1743SMark F. Adams 
855676e1743SMark F. Adams #undef __FUNCT__
856676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning"
857676e1743SMark F. Adams /*@
858676e1743SMark F. Adams    PCGAMGAvoidRepartitioning - Do not repartition the coarse grids
859676e1743SMark F. Adams 
860676e1743SMark F. Adams    Collective on PC
861676e1743SMark F. Adams 
862676e1743SMark F. Adams    Input Parameters:
863676e1743SMark F. Adams .  pc - the preconditioner context
864676e1743SMark F. Adams 
865676e1743SMark F. Adams 
866676e1743SMark F. Adams    Options Database Key:
867676e1743SMark F. Adams .  -pc_gamg_avoid_repartitioning
868676e1743SMark F. Adams 
869676e1743SMark F. Adams    Level: intermediate
870676e1743SMark F. Adams 
871676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
872676e1743SMark F. Adams 
873676e1743SMark F. Adams .seealso: ()
874676e1743SMark F. Adams @*/
875676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning(PC pc, PetscBool n)
876676e1743SMark F. Adams {
877676e1743SMark F. Adams   PetscErrorCode ierr;
878676e1743SMark F. Adams 
879676e1743SMark F. Adams   PetscFunctionBegin;
880676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
881676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGAvoidRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
882676e1743SMark F. Adams   PetscFunctionReturn(0);
883676e1743SMark F. Adams }
884676e1743SMark F. Adams 
885676e1743SMark F. Adams EXTERN_C_BEGIN
886676e1743SMark F. Adams #undef __FUNCT__
887676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning_GAMG"
888676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning_GAMG(PC pc, PetscBool n)
889676e1743SMark F. Adams {
890676e1743SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
891676e1743SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
892676e1743SMark F. Adams 
893676e1743SMark F. Adams   PetscFunctionBegin;
894676e1743SMark F. Adams   s_avoid_repart = n;
895676e1743SMark F. Adams   PetscFunctionReturn(0);
896676e1743SMark F. Adams }
897676e1743SMark F. Adams EXTERN_C_END
898676e1743SMark F. Adams 
899676e1743SMark F. Adams #undef __FUNCT__
9003542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
9013542efc5SMark F. Adams /*@
9023542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
9033542efc5SMark F. Adams 
9043542efc5SMark F. Adams    Not collective on PC
9053542efc5SMark F. Adams 
9063542efc5SMark F. Adams    Input Parameters:
9073542efc5SMark F. Adams .  pc - the preconditioner context
9083542efc5SMark F. Adams 
9093542efc5SMark F. Adams 
9103542efc5SMark F. Adams    Options Database Key:
9113542efc5SMark F. Adams .  -pc_gamg_threshold
9123542efc5SMark F. Adams 
9133542efc5SMark F. Adams    Level: intermediate
9143542efc5SMark F. Adams 
9153542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
9163542efc5SMark F. Adams 
9173542efc5SMark F. Adams .seealso: ()
9183542efc5SMark F. Adams @*/
9193542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
9203542efc5SMark F. Adams {
9213542efc5SMark F. Adams   PetscErrorCode ierr;
9223542efc5SMark F. Adams 
9233542efc5SMark F. Adams   PetscFunctionBegin;
9243542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9253542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
9263542efc5SMark F. Adams   PetscFunctionReturn(0);
9273542efc5SMark F. Adams }
9283542efc5SMark F. Adams 
9293542efc5SMark F. Adams EXTERN_C_BEGIN
9303542efc5SMark F. Adams #undef __FUNCT__
9313542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
9323542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
9333542efc5SMark F. Adams {
9343542efc5SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
9353542efc5SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
9363542efc5SMark F. Adams 
9373542efc5SMark F. Adams   PetscFunctionBegin;
9383542efc5SMark F. Adams   s_threshold = n;
9393542efc5SMark F. Adams   PetscFunctionReturn(0);
9403542efc5SMark F. Adams }
9413542efc5SMark F. Adams EXTERN_C_END
9423542efc5SMark F. Adams 
9433542efc5SMark F. Adams #undef __FUNCT__
944676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType"
945676e1743SMark F. Adams /*@
946676e1743SMark F. Adams    PCGAMGSetSolverType - Set solution method.
947676e1743SMark F. Adams 
948676e1743SMark F. Adams    Collective on PC
949676e1743SMark F. Adams 
950676e1743SMark F. Adams    Input Parameters:
951676e1743SMark F. Adams .  pc - the preconditioner context
952676e1743SMark F. Adams 
953676e1743SMark F. Adams 
954676e1743SMark F. Adams    Options Database Key:
9553542efc5SMark F. Adams .  -pc_gamg_type
956676e1743SMark F. Adams 
957676e1743SMark F. Adams    Level: intermediate
958676e1743SMark F. Adams 
959676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
960676e1743SMark F. Adams 
961676e1743SMark F. Adams .seealso: ()
962676e1743SMark F. Adams @*/
963676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz )
964676e1743SMark F. Adams {
965676e1743SMark F. Adams   PetscErrorCode ierr;
966676e1743SMark F. Adams 
967676e1743SMark F. Adams   PetscFunctionBegin;
968676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
969676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz));
970676e1743SMark F. Adams   CHKERRQ(ierr);
971676e1743SMark F. Adams   PetscFunctionReturn(0);
972676e1743SMark F. Adams }
973676e1743SMark F. Adams 
974676e1743SMark F. Adams EXTERN_C_BEGIN
975676e1743SMark F. Adams #undef __FUNCT__
976676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG"
977676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz )
978676e1743SMark F. Adams {
979676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
980676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
981676e1743SMark F. Adams 
982676e1743SMark F. Adams   PetscFunctionBegin;
983676e1743SMark F. Adams   if(sz < 64) strcpy(pc_gamg->m_type,str);
984676e1743SMark F. Adams   PetscFunctionReturn(0);
985676e1743SMark F. Adams }
986676e1743SMark F. Adams EXTERN_C_END
987676e1743SMark F. Adams 
9885b89ad90SMark F. Adams #undef __FUNCT__
9895b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
9905b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
9915b89ad90SMark F. Adams {
992676e1743SMark F. Adams   PetscErrorCode  ierr;
993676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
994676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
995676e1743SMark F. Adams   PetscBool flag;
9965b89ad90SMark F. Adams 
9975b89ad90SMark F. Adams   PetscFunctionBegin;
998676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
999676e1743SMark F. Adams   {
1000676e1743SMark F. Adams     ierr = PetscOptionsString("-pc_gamg_type",
1001470e26f8SMark F. Adams                               "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid (default)",
1002676e1743SMark F. Adams                               "PCGAMGSetSolverType",
1003676e1743SMark F. Adams                               pc_gamg->m_type,
1004676e1743SMark F. Adams                               pc_gamg->m_type,
1005676e1743SMark F. Adams                               64,
1006676e1743SMark F. Adams                               &flag );
1007676e1743SMark F. Adams     CHKERRQ(ierr);
1008*d8c859f0SMark F. Adams     if( flag && pc_gamg->m_data != 0 ) {
1009*d8c859f0SMark F. Adams       if( (strcmp(pc_gamg->m_type,"sa")==0 && pc_gamg->m_method != 2) ||
1010*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"pa")==0 && pc_gamg->m_method != 1) ||
1011*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"geo")==0 && pc_gamg->m_method != 0) ) {
1012*d8c859f0SMark F. Adams         SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "PCSetFromOptions called of PCSetCoordinates (with new method, after data was created)");
1013*d8c859f0SMark F. Adams       }
1014*d8c859f0SMark F. Adams     }
10152a44bfbaSMark F. Adams 
1016470e26f8SMark F. Adams     if (flag && strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2;
1017470e26f8SMark F. Adams     else if (flag && strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1;
1018470e26f8SMark F. Adams     else pc_gamg->m_method = 0;
10192a44bfbaSMark F. Adams 
10203542efc5SMark F. Adams     /* common (static) variable */
1021676e1743SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_avoid_repartitioning",
1022676e1743SMark F. Adams                             "Do not repartion coarse grids (false)",
1023676e1743SMark F. Adams                             "PCGAMGAvoidRepartitioning",
1024676e1743SMark F. Adams                             s_avoid_repart,
1025676e1743SMark F. Adams                             &s_avoid_repart,
1026676e1743SMark F. Adams                             &flag);
1027676e1743SMark F. Adams     CHKERRQ(ierr);
1028676e1743SMark F. Adams 
10293542efc5SMark F. Adams     /* common (static) variable */
1030676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1031676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1032676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
1033676e1743SMark F. Adams                            s_min_eq_proc,
1034676e1743SMark F. Adams                            &s_min_eq_proc,
1035676e1743SMark F. Adams                            &flag );
1036676e1743SMark F. Adams     CHKERRQ(ierr);
10373542efc5SMark F. Adams 
10383542efc5SMark F. Adams     /* common (static) variable */
10393542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
10403542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
10413542efc5SMark F. Adams                             "PCGAMGSetThreshold",
10423542efc5SMark F. Adams                             s_threshold,
10433542efc5SMark F. Adams                             &s_threshold,
10443542efc5SMark F. Adams                             &flag );
10453542efc5SMark F. Adams     CHKERRQ(ierr);
1046676e1743SMark F. Adams   }
1047676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1048676e1743SMark F. Adams 
10495b89ad90SMark F. Adams   PetscFunctionReturn(0);
10505b89ad90SMark F. Adams }
10515b89ad90SMark F. Adams 
10525b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
10535b89ad90SMark F. Adams /*
10545b89ad90SMark F. Adams  PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG
10555b89ad90SMark F. Adams 
10565b89ad90SMark F. Adams    Input Parameter:
10575b89ad90SMark F. Adams .  pc - the preconditioner context
10585b89ad90SMark F. Adams 
10595b89ad90SMark F. Adams    Application Interface Routine: PCCreate()
10605b89ad90SMark F. Adams 
10615b89ad90SMark F. Adams   */
10625b89ad90SMark F. Adams  /* MC
10635b89ad90SMark F. Adams      PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide
10645b89ad90SMark F. Adams        fine grid discretization matrix and coordinates on the fine grid.
10655b89ad90SMark F. Adams 
10665b89ad90SMark F. Adams    Options Database Key:
10675b89ad90SMark F. Adams    Multigrid options(inherited)
10685b89ad90SMark F. Adams +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
10695b89ad90SMark F. Adams .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
10705b89ad90SMark F. Adams .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
10715b89ad90SMark F. Adams    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
10725b89ad90SMark F. Adams    GAMG options:
10735b89ad90SMark F. Adams 
10745b89ad90SMark F. Adams    Level: intermediate
10755b89ad90SMark F. Adams   Concepts: multigrid
10765b89ad90SMark F. Adams 
10775b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
10785b89ad90SMark F. Adams            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
10795b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
10805b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
10815b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
10825b89ad90SMark F. Adams M */
10835b89ad90SMark F. Adams 
10845b89ad90SMark F. Adams EXTERN_C_BEGIN
10855b89ad90SMark F. Adams #undef __FUNCT__
10865b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
10875b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
10885b89ad90SMark F. Adams {
10895b89ad90SMark F. Adams   PetscErrorCode  ierr;
10905b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
10915b89ad90SMark F. Adams   PC_MG           *mg;
10925ef31b24SMark F. Adams   PetscClassId     cookie;
10935ee9c036SSatish Balay #if defined PETSC_USE_LOG
10945ee9c036SSatish Balay   static int count = 0;
10955ee9c036SSatish Balay #endif
10965b89ad90SMark F. Adams 
10975b89ad90SMark F. Adams   PetscFunctionBegin;
10985b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
10995b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
11005b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
11015b89ad90SMark F. Adams 
11025b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
11035b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
110403a628feSMark F. Adams   pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0;
11055b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
11065b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
11075b89ad90SMark F. Adams 
11085b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
11095b89ad90SMark F. Adams 
11105b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
11115b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
11125b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
11135b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
11145b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
11155b89ad90SMark F. Adams 
11165b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
11175b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
11185b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
1119c97e1df0SMark F. Adams 					    PCSetCoordinates_GAMG);
1120c97e1df0SMark F. Adams   CHKERRQ(ierr);
1121c97e1df0SMark F. Adams 
1122676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1123676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1124676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1125676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1126676e1743SMark F. Adams   CHKERRQ(ierr);
1127676e1743SMark F. Adams 
1128676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1129676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_C",
1130676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_GAMG",
1131676e1743SMark F. Adams 					    PCGAMGAvoidRepartitioning_GAMG);
1132676e1743SMark F. Adams   CHKERRQ(ierr);
1133676e1743SMark F. Adams 
1134676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
11353542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
11363542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
11373542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
11383542efc5SMark F. Adams   CHKERRQ(ierr);
11393542efc5SMark F. Adams 
11403542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1141676e1743SMark F. Adams 					    "PCGAMGSetSolverType_C",
1142676e1743SMark F. Adams 					    "PCGAMGSetSolverType_GAMG",
1143676e1743SMark F. Adams 					    PCGAMGSetSolverType_GAMG);
1144676e1743SMark F. Adams   CHKERRQ(ierr);
1145c97e1df0SMark F. Adams 
1146b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
1147785cba28SMark F. Adams   if( count++ == 0 ) {
11485ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
1149b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
11502a44bfbaSMark F. Adams     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
11512a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]);
11522a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]);
11532a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]);
1154b4fbaa2aSMark F. Adams     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1155b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1156b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1157b4fbaa2aSMark F. Adams     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1158b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_events[SET7]);
1159b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]); */
1160b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1161b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1162b4fbaa2aSMark F. Adams     PetscLogEventRegister("  PL repartition", cookie, &gamg_setup_events[SET12]);
1163b4fbaa2aSMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1164b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1165b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
11665ef31b24SMark F. Adams 
1167b4fbaa2aSMark F. Adams     /* create timer stages */
1168b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1169b4fbaa2aSMark F. Adams     {
1170b4fbaa2aSMark F. Adams       char str[32];
1171b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1172b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1173b4fbaa2aSMark F. Adams       PetscInt lidx;
1174b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1175b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1176b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1177b4fbaa2aSMark F. Adams       }
1178b4fbaa2aSMark F. Adams     }
1179b4fbaa2aSMark F. Adams #endif
1180b4fbaa2aSMark F. Adams   }
1181b4fbaa2aSMark F. Adams #endif
11825b89ad90SMark F. Adams   PetscFunctionReturn(0);
11835b89ad90SMark F. Adams }
11845b89ad90SMark F. Adams EXTERN_C_END
1185