xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision fc4362bfe084354807e64d80cc2dddcb1c5ddd53)
15b89ad90SMark F. Adams /*
25b89ad90SMark F. Adams  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
45ef31b24SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>
5f96513f1SMatthew G Knepley 
6f96513f1SMatthew G Knepley PetscLogEvent gamg_setup_stages[NUM_SET];
7f96513f1SMatthew G Knepley 
85b89ad90SMark F. Adams /* Private context for the GAMG preconditioner */
95b89ad90SMark F. Adams typedef struct gamg_TAG {
105b89ad90SMark F. Adams   PetscInt       m_dim;
115b89ad90SMark F. Adams   PetscInt       m_Nlevels;
125b89ad90SMark F. Adams   PetscInt       m_data_sz;
135b89ad90SMark F. Adams   PetscReal     *m_data; /* blocked vector of vertex data on fine grid (coordinates) */
145b89ad90SMark F. Adams } PC_GAMG;
155b89ad90SMark F. Adams 
16*fc4362bfSMark F. Adams #define TOP_GRID_LIM 100
175ef31b24SMark F. Adams 
185b89ad90SMark F. Adams /* -----------------------------------------------------------------------------*/
195b89ad90SMark F. Adams #undef __FUNCT__
205b89ad90SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
215b89ad90SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
225b89ad90SMark F. Adams {
235b89ad90SMark F. Adams   PetscErrorCode  ierr;
245b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
255b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
265b89ad90SMark F. Adams 
275b89ad90SMark F. Adams   PetscFunctionBegin;
285b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr);
295b89ad90SMark F. Adams   PetscFunctionReturn(0);
305b89ad90SMark F. Adams }
315b89ad90SMark F. Adams 
325b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
335b89ad90SMark F. Adams /*
345b89ad90SMark F. Adams    PCSetCoordinates_GAMG
355b89ad90SMark F. Adams 
365b89ad90SMark F. Adams    Input Parameter:
375b89ad90SMark F. Adams    .  pc - the preconditioner context
385b89ad90SMark F. Adams */
39a92563c5SMark F. Adams EXTERN_C_BEGIN
405b89ad90SMark F. Adams #undef __FUNCT__
415b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG"
426c237d78SBarry Smith PetscErrorCode PCSetCoordinates_GAMG(PC pc, const int ndm,PetscReal *coords )
435b89ad90SMark F. Adams {
445b89ad90SMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
455b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
466c237d78SBarry Smith   PetscErrorCode ierr;
476c237d78SBarry Smith   PetscInt       bs, my0, tt;
486c237d78SBarry Smith   Mat            mat = pc->pmat;
496c237d78SBarry Smith   PetscInt       arrsz;
505b89ad90SMark F. Adams 
515b89ad90SMark F. Adams   PetscFunctionBegin;
525b89ad90SMark F. Adams   ierr  = MatGetBlockSize( mat, &bs );               CHKERRQ( ierr );
535b89ad90SMark F. Adams   ierr  = MatGetOwnershipRange( mat, &my0, &tt ); CHKERRQ(ierr);
546c237d78SBarry Smith   arrsz = (tt-my0)/bs*ndm;
555b89ad90SMark F. Adams 
565b89ad90SMark F. Adams   // put coordinates
576c237d78SBarry Smith   if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) {
585b89ad90SMark F. Adams     ierr = PetscFree( pc_gamg->m_data );  CHKERRQ(ierr);
596c237d78SBarry Smith     ierr = PetscMalloc(arrsz*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr);
605b89ad90SMark F. Adams   }
615b89ad90SMark F. Adams 
625b89ad90SMark F. Adams   /* copy data in */
635b89ad90SMark F. Adams   for(tt=0;tt<arrsz;tt++){
645b89ad90SMark F. Adams     pc_gamg->m_data[tt] = coords[tt];
655b89ad90SMark F. Adams   }
665b89ad90SMark F. Adams   pc_gamg->m_data_sz = arrsz;
675b89ad90SMark F. Adams   pc_gamg->m_dim = ndm;
685b89ad90SMark F. Adams 
695b89ad90SMark F. Adams   PetscFunctionReturn(0);
705b89ad90SMark F. Adams }
71a92563c5SMark F. Adams EXTERN_C_END
725b89ad90SMark F. Adams 
735b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
745b89ad90SMark F. Adams /*
7511e60469SMark F. Adams    partitionLevel
765b89ad90SMark F. Adams 
775b89ad90SMark F. Adams    Input Parameter:
783530afc2SMark F. Adams    . a_Amat_fine - matrix on this fine (k) level
7911e60469SMark F. Adams    . a_dime - 2 or 3
803530afc2SMark F. Adams    In/Output Parameter:
813530afc2SMark F. Adams    . a_P_inout - prolongation operator to the next level (k-1)
823530afc2SMark F. Adams    . a_coarse_crds - coordinates that need to be moved
833530afc2SMark F. Adams    . a_active_proc - number of active procs
8411e60469SMark F. Adams    Output Parameter:
853530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
865b89ad90SMark F. Adams */
875cb416c2SMark F. Adams 
885b89ad90SMark F. Adams #undef __FUNCT__
8911e60469SMark F. Adams #define __FUNCT__ "partitionLevel"
903530afc2SMark F. Adams PetscErrorCode partitionLevel( Mat a_Amat_fine,
9111e60469SMark F. Adams                                PetscInt a_dim,
923530afc2SMark F. Adams                                Mat *a_P_inout,
933530afc2SMark F. Adams                                PetscReal **a_coarse_crds,
943530afc2SMark F. Adams                                PetscMPIInt *a_active_proc,
953530afc2SMark F. Adams                                Mat *a_Amat_crs
9611e60469SMark F. Adams                             )
975b89ad90SMark F. Adams {
985b89ad90SMark F. Adams   PetscErrorCode   ierr;
993530afc2SMark F. Adams   Mat              Amat, Pnew, Pold = *a_P_inout;
10011e60469SMark F. Adams   IS               new_indices,isnum;
1013530afc2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
1025ef31b24SMark F. Adams   PetscMPIInt      nactive_procs,mype,npe;
10311e60469SMark F. Adams   PetscInt         Istart,Iend,Istart0,Iend0,ncrs0,ncrs_new,bs=1; /* bs ??? */
1045b89ad90SMark F. Adams 
1055b89ad90SMark F. Adams   PetscFunctionBegin;
10611e60469SMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
10711e60469SMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
10811e60469SMark F. Adams   /* RAP */
1093530afc2SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Amat ); CHKERRQ(ierr);
1103530afc2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart0, &Iend0 );    CHKERRQ(ierr); /* x2 size of 'a_coarse_crds' */
11111e60469SMark F. Adams   ncrs0 = (Iend0 - Istart0)/bs;
112acadaa72SMark F. Adams 
11311e60469SMark F. Adams   /* Repartition Amat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
11411e60469SMark F. Adams   {
1153530afc2SMark F. Adams     PetscInt        neq,N,counts[npe];
11611e60469SMark F. Adams     IS              isnewproc;
1175ef31b24SMark F. Adams     PetscMPIInt     new_npe,targ_npe;
1183530afc2SMark F. Adams 
1193530afc2SMark F. Adams     ierr = MatGetSize( Amat, &neq, &N );CHKERRQ(ierr);
1203530afc2SMark F. Adams #define MIN_EQ_PROC 100
1215ef31b24SMark F. Adams     nactive_procs = *a_active_proc;
1223530afc2SMark F. Adams     targ_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */
1235ef31b24SMark F. Adams     if( targ_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1; /* chop coarsest grid */
1245ef31b24SMark F. Adams     else if (targ_npe >= nactive_procs ) new_npe = nactive_procs; /* no change */
1253530afc2SMark F. Adams     else {
1265ef31b24SMark F. Adams       PetscMPIInt     factstart,fact;
1273530afc2SMark F. Adams       new_npe = -9999;
1285ef31b24SMark F. Adams       factstart = nactive_procs;
1293530afc2SMark F. Adams       for(fact=factstart;fact>0;fact--){ /* try to find a better number of procs */
1305ef31b24SMark F. Adams         if( nactive_procs%fact==0 && neq/(nactive_procs/fact) > MIN_EQ_PROC ) {
1315ef31b24SMark F. Adams           new_npe = nactive_procs/fact;
1323530afc2SMark F. Adams         }
1333530afc2SMark F. Adams       }
1343530afc2SMark F. Adams       assert(new_npe != -9999);
1353530afc2SMark F. Adams     }
1363530afc2SMark F. Adams     *a_active_proc = new_npe; /* output for next time */
1375ef31b24SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t\t%s num active procs = %d (N loc = %d)\n",__FUNCT__,new_npe,ncrs0);
13811e60469SMark F. Adams     { /* partition: get 'isnewproc' */
13911e60469SMark F. Adams       MatPartitioning  mpart;
1405ef31b24SMark F. Adams       Mat              adj;
1415ef31b24SMark F. Adams       const PetscInt  *is_idx;
1425ef31b24SMark F. Adams       PetscInt         is_sz,kk,jj,old_fact=(npe/nactive_procs),new_fact=(npe/new_npe),*isnewproc_idx;
143c9a0b8beSMark F. Adams       /* create sub communicator  */
144c9a0b8beSMark F. Adams       MPI_Comm cm,new_comm;
1455ef31b24SMark F. Adams       int membershipKey = mype % old_fact;
146c9a0b8beSMark F. Adams       ierr = MPI_Comm_split(wcomm, membershipKey, mype, &cm); CHKERRQ(ierr);
147c9a0b8beSMark F. Adams       ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr);
148c9a0b8beSMark F. Adams       ierr = MPI_Comm_free( &cm );                            CHKERRQ(ierr);
149c9a0b8beSMark F. Adams 
1505ef31b24SMark F. Adams       /* MatPartitioningApply call MatConvert, which is collective */
1515ef31b24SMark F. Adams       ierr = MatConvert(Amat,MATMPIADJ,MAT_INITIAL_MATRIX,&adj);CHKERRQ(ierr);
1525ef31b24SMark F. Adams       if( membershipKey == 0 ) {
1535ef31b24SMark F. Adams         /* hack to fix global data that pmetis.c uses in 'adj' */
1545ef31b24SMark F. Adams         for( kk=0 , jj=0 ; kk<=npe ; jj++, kk += old_fact ) {
1555ef31b24SMark F. Adams           adj->rmap->range[jj] = adj->rmap->range[kk];
1565ef31b24SMark F. Adams         }
1575ef31b24SMark F. Adams         ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr);
1585ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
15911e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr);
1605ef31b24SMark F. Adams         ierr = MatPartitioningSetNParts( mpart, new_npe ); CHKERRQ(ierr);
16111e60469SMark F. Adams         ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr);
16211e60469SMark F. Adams         ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr);
1635ef31b24SMark F. Adams         /* collect IS info */
1645ef31b24SMark F. Adams         ierr = ISGetLocalSize( isnewproc, &is_sz );        CHKERRQ(ierr);
1655ef31b24SMark F. Adams         ierr = PetscMalloc( is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr);
1665ef31b24SMark F. Adams         ierr = ISGetIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
1675ef31b24SMark F. Adams         for(kk=0;kk<is_sz;kk++){
168c9a0b8beSMark F. Adams           /* spread partitioning across machine - probably the right thing to do but machine specific */
1695ef31b24SMark F. Adams           isnewproc_idx[kk] = is_idx[kk] * new_fact;
1705ef31b24SMark F. Adams         }
1715ef31b24SMark F. Adams         ierr = ISRestoreIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
1725ef31b24SMark F. Adams         ierr = ISDestroy( &isnewproc );                    CHKERRQ(ierr);
1735ef31b24SMark F. Adams       }
1745ef31b24SMark F. Adams       else {
1755ef31b24SMark F. Adams         isnewproc_idx = 0;
1765ef31b24SMark F. Adams         is_sz = 0;
1775ef31b24SMark F. Adams       }
1785ef31b24SMark F. Adams       ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
1795ef31b24SMark F. Adams       ierr = MPI_Comm_free( &new_comm );    CHKERRQ(ierr);
1805ef31b24SMark F. Adams       ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc );
1815ef31b24SMark F. Adams       if( membershipKey == 0 ) {
1825ef31b24SMark F. Adams         ierr = PetscFree( isnewproc_idx );  CHKERRQ(ierr);
1835ef31b24SMark F. Adams       }
18411e60469SMark F. Adams     }
185c9a0b8beSMark F. Adams 
18611e60469SMark F. Adams     /*
18711e60469SMark F. Adams      Create an index set from the isnewproc index set to indicate the mapping TO
18811e60469SMark F. Adams      */
18911e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
19011e60469SMark F. Adams     /*
19111e60469SMark F. Adams      Determine how many elements are assigned to each processor
19211e60469SMark F. Adams      */
19311e60469SMark F. Adams     ierr = ISPartitioningCount( isnewproc, npe, counts ); CHKERRQ(ierr);
19411e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
19511e60469SMark F. Adams     ncrs_new = counts[mype];
19611e60469SMark F. Adams   }
1975ef31b24SMark F. Adams 
19811e60469SMark F. Adams   { /* Create a vector to contain the newly ordered element information */
19911e60469SMark F. Adams     const PetscInt *idx;
20011e60469SMark F. Adams     PetscInt        i,j,k;
20111e60469SMark F. Adams     IS              isscat;
20211e60469SMark F. Adams     PetscScalar    *array;
20311e60469SMark F. Adams     Vec             src_crd, dest_crd;
2043530afc2SMark F. Adams     PetscReal      *coords = *a_coarse_crds;
20511e60469SMark F. Adams     VecScatter      vecscat;
206acadaa72SMark F. Adams 
20711e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
20811e60469SMark F. Adams     ierr = VecSetSizes( dest_crd, a_dim*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
20911e60469SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/
21011e60469SMark F. Adams     /*
21111e60469SMark F. Adams      There are three data items per cell (element), the integer vertex numbers of its three
21211e60469SMark F. Adams      coordinates (we convert to double to use the scatter) (one can think
21311e60469SMark F. Adams      of the vectors of having a block size of 3, then there is one index in idx[] for each element)
21411e60469SMark F. Adams      */
21511e60469SMark F. Adams     {
21611e60469SMark F. Adams       PetscInt tidx[ncrs0*a_dim];
21711e60469SMark F. Adams       ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
21811e60469SMark F. Adams       for (i=0,j=0; i<ncrs0; i++) {
21911e60469SMark F. Adams         for (k=0; k<a_dim; k++) tidx[j++] = idx[i]*a_dim + k;
22011e60469SMark F. Adams       }
22111e60469SMark F. Adams       ierr = ISCreateGeneral( wcomm, a_dim*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
22211e60469SMark F. Adams       CHKERRQ(ierr);
22311e60469SMark F. Adams       ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
22411e60469SMark F. Adams     }
22511e60469SMark F. Adams     /*
22611e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
22711e60469SMark F. Adams      */
22811e60469SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, a_dim*ncrs0, &src_crd ); CHKERRQ(ierr);
22911e60469SMark F. Adams     ierr = VecGetArray( src_crd, &array );                        CHKERRQ(ierr);
23011e60469SMark F. Adams     for (i=0; i<a_dim*ncrs0; i++) array[i] = coords[i];
23111e60469SMark F. Adams     ierr = VecRestoreArray( src_crd, &array );                     CHKERRQ(ierr);
23211e60469SMark F. Adams     /*
23311e60469SMark F. Adams      Scatter the element vertex information (still in the original vertex ordering)
23411e60469SMark F. Adams      to the correct processor
23511e60469SMark F. Adams      */
23611e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
23711e60469SMark F. Adams     CHKERRQ(ierr);
23811e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
23911e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
24011e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
24111e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
24211e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
24311e60469SMark F. Adams     /*
24411e60469SMark F. Adams      Put the element vertex data into a new allocation of the gdata->ele
24511e60469SMark F. Adams      */
2463530afc2SMark F. Adams     ierr = PetscFree( *a_coarse_crds );    CHKERRQ(ierr);
2473530afc2SMark F. Adams     ierr = PetscMalloc( a_dim*ncrs_new*sizeof(PetscReal), a_coarse_crds );    CHKERRQ(ierr);
2483530afc2SMark F. Adams     coords = *a_coarse_crds; /* convient */
24911e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
25011e60469SMark F. Adams     for (i=0; i<a_dim*ncrs_new; i++) coords[i] = PetscRealPart(array[i]);
25111e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
25211e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
25311e60469SMark F. Adams   }
25411e60469SMark F. Adams   /*
25511e60469SMark F. Adams    Invert for MatGetSubMatrix
25611e60469SMark F. Adams    */
25711e60469SMark F. Adams   ierr = ISInvertPermutation( isnum, ncrs_new, &new_indices ); CHKERRQ(ierr);
25811e60469SMark F. Adams   ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
25911e60469SMark F. Adams   ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
26011e60469SMark F. Adams   /* A_crs output */
2613530afc2SMark F. Adams   ierr = MatGetSubMatrix( Amat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
26211e60469SMark F. Adams   CHKERRQ(ierr);
26311e60469SMark F. Adams   ierr = MatDestroy( &Amat ); CHKERRQ(ierr);
2643530afc2SMark F. Adams   Amat = *a_Amat_crs;
26511e60469SMark F. Adams   /* prolongator */
26611e60469SMark F. Adams   ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
26711e60469SMark F. Adams   {
26811e60469SMark F. Adams     IS findices;
26911e60469SMark F. Adams     ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
27011e60469SMark F. Adams     ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
27111e60469SMark F. Adams     CHKERRQ(ierr);
27211e60469SMark F. Adams     ierr = ISDestroy( &findices ); CHKERRQ(ierr);
27311e60469SMark F. Adams   }
2743530afc2SMark F. Adams   ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
2753530afc2SMark F. Adams   *a_P_inout = Pnew; /* output */
276acadaa72SMark F. Adams 
27711e60469SMark F. Adams   ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
2785b89ad90SMark F. Adams 
2795b89ad90SMark F. Adams   PetscFunctionReturn(0);
2805b89ad90SMark F. Adams }
2815b89ad90SMark F. Adams 
2825ef31b24SMark F. Adams #define GAMG_MAXLEVELS 20
2835ef31b24SMark F. Adams #if defined(PETSC_USE_LOG)
2845ef31b24SMark F. Adams PetscLogStage  gamg_stages[20];
2855ef31b24SMark F. Adams #endif
2865b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
2875b89ad90SMark F. Adams /*
2885b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
2895b89ad90SMark F. Adams                     by setting data structures and options.
2905b89ad90SMark F. Adams 
2915b89ad90SMark F. Adams    Input Parameter:
2925b89ad90SMark F. Adams .  pc - the preconditioner context
2935b89ad90SMark F. Adams 
2945b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
2955b89ad90SMark F. Adams 
2965b89ad90SMark F. Adams    Notes:
2975b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
2985b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
2995b89ad90SMark F. Adams */
3005b89ad90SMark F. Adams #undef __FUNCT__
3015b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
3025b89ad90SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC pc )
3035b89ad90SMark F. Adams {
3045b89ad90SMark F. Adams   PetscErrorCode  ierr;
3055b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
3065b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
3075b89ad90SMark F. Adams   Mat              Amat = pc->mat, Pmat = pc->pmat;
3085b89ad90SMark F. Adams   PetscBool        isSeq, isMPI;
3095b89ad90SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, lidx;
3105b89ad90SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
3113530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
3120205a208SMark F. Adams   PetscBool isOK;
3135b89ad90SMark F. Adams 
3145ef31b24SMark F. Adams   Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];  PetscReal *coarse_crds = 0, *crds = pc_gamg->m_data;
3155ef31b24SMark F. Adams 
3165b89ad90SMark F. Adams   PetscFunctionBegin;
317baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
318baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
3195b89ad90SMark F. Adams   if (pc->setupcalled){
3205b89ad90SMark F. Adams     /* no state data in GAMG to destroy (now) */
3215b89ad90SMark F. Adams     ierr = PCReset_MG(pc); CHKERRQ(ierr);
3225b89ad90SMark F. Adams   }
3236c237d78SBarry Smith   if (!pc_gamg->m_data) SETERRQ(wcomm,PETSC_ERR_SUP,"PCSetUp_GAMG called before PCSetCoordinates");
3245b89ad90SMark F. Adams   /* setup special features of PCGAMG */
3255b89ad90SMark F. Adams   ierr = PetscTypeCompare((PetscObject)Amat, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
3265b89ad90SMark F. Adams   ierr = PetscTypeCompare((PetscObject)Amat, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
3275b89ad90SMark F. Adams   if (isMPI) {
3285b89ad90SMark F. Adams   } else if (isSeq) {
3295b89ad90SMark F. Adams   } else SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with GAMG. GAMG can only handle AIJ matrices.",((PetscObject)Amat)->type_name);
3305b89ad90SMark F. Adams 
3315b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
3325b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
3335b89ad90SMark F. Adams   if(bs!=1) SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "GAMG only supports scalar prblems bs = '%d'.",bs);
3345b89ad90SMark F. Adams 
3355b89ad90SMark F. Adams   /* Get A_i and R_i */
3363530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
3373530afc2SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"[%d]%s level %d N=%d\n",0,__FUNCT__,0,N);
3380205a208SMark F. Adams   for (level=0, Aarr[0] = Pmat, nactivepe = npe;
3395ef31b24SMark F. Adams        level < GAMG_MAXLEVELS-1 && (level==0 || M>TOP_GRID_LIM); /* hard wired stopping logic */
3400205a208SMark F. Adams        level++ ) {
3415b89ad90SMark F. Adams     level1 = level + 1;
3425ef31b24SMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr);
3430f9369f8SMark F. Adams     ierr = createProlongation( Aarr[level], crds, pc_gamg->m_dim,
3445ef31b24SMark F. Adams                                &Parr[level1], &coarse_crds, &isOK );
3455b89ad90SMark F. Adams     CHKERRQ(ierr);
3465ef31b24SMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr);
3475ef31b24SMark F. Adams 
3485b89ad90SMark F. Adams     ierr = PetscFree( crds ); CHKERRQ( ierr );
349baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
350baa4e9faSMark F. Adams     if( isOK ) {
3515ef31b24SMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr);
3525ef31b24SMark F. Adams       ierr = partitionLevel( Aarr[level], pc_gamg->m_dim, &Parr[level1], &coarse_crds, &nactivepe, &Aarr[level1] );
3533530afc2SMark F. Adams       CHKERRQ(ierr);
3545ef31b24SMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr);
3553530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
3563530afc2SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"[%d]%s done level %d N=%d, active pe = %d\n",0,__FUNCT__,level1,N,nactivepe);
357baa4e9faSMark F. Adams     }
358baa4e9faSMark F. Adams     else{
359baa4e9faSMark F. Adams       break;
360baa4e9faSMark F. Adams     }
36111e60469SMark F. Adams     crds = coarse_crds;
3625b89ad90SMark F. Adams   }
3635b89ad90SMark F. Adams   ierr = PetscFree( coarse_crds ); CHKERRQ( ierr );
364baa4e9faSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
3655b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
3665b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
3675b89ad90SMark F. Adams   fine_level = level;
3685b89ad90SMark F. Adams   ierr = PCMGSetLevels(pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
3695b89ad90SMark F. Adams 
370*fc4362bfSMark F. Adams   /* set default smoothers */
371*fc4362bfSMark F. Adams   for (level=1,lidx=pc_gamg->m_Nlevels-2;
372*fc4362bfSMark F. Adams        level<=fine_level;
373*fc4362bfSMark F. Adams        level++,lidx--) {
3749e66a1f7SMark F. Adams     PetscReal emax = 3.0, emin;
3755b89ad90SMark F. Adams     KSP smoother; PC subpc;
3765b89ad90SMark F. Adams     ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
3775b89ad90SMark F. Adams     ierr = KSPSetType(smoother,KSPCHEBYCHEV);CHKERRQ(ierr);
3785b89ad90SMark F. Adams     ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
3795b89ad90SMark F. Adams     ierr = PCSetType(subpc,PCJACOBI);CHKERRQ(ierr);
380*fc4362bfSMark F. Adams     { /* eigen estimate 'emax' */
381*fc4362bfSMark F. Adams       KSP eksp; Mat Amat = Aarr[lidx];
382*fc4362bfSMark F. Adams       Vec bb, xx; PC pc;
383*fc4362bfSMark F. Adams       PetscInt N1, N0, tt;
384*fc4362bfSMark F. Adams       ierr = MatGetVecs( Amat, &bb, 0 );         CHKERRQ(ierr);
385*fc4362bfSMark F. Adams       ierr = MatGetVecs( Amat, &xx, 0 );         CHKERRQ(ierr);
386*fc4362bfSMark F. Adams       {
387*fc4362bfSMark F. Adams 	PetscRandom    rctx;
388*fc4362bfSMark F. Adams 	ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
389*fc4362bfSMark F. Adams 	ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
390*fc4362bfSMark F. Adams 	ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
391*fc4362bfSMark F. Adams 	ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
3925b89ad90SMark F. Adams       }
393*fc4362bfSMark F. Adams       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
394*fc4362bfSMark F. Adams       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
395*fc4362bfSMark F. Adams       ierr = KSPSetOperators( eksp, Amat, Amat, DIFFERENT_NONZERO_PATTERN ); CHKERRQ( ierr );
396*fc4362bfSMark F. Adams       ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
397*fc4362bfSMark F. Adams       ierr = PCSetType( pc, PCJACOBI ); CHKERRQ(ierr); /* should be same as above */
398*fc4362bfSMark F. Adams       ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 5 );
399*fc4362bfSMark F. Adams       CHKERRQ(ierr);
400*fc4362bfSMark F. Adams       ierr = KSPSetConvergenceTest( eksp, KSPSkipConverged, 0, 0 ); CHKERRQ(ierr);
401*fc4362bfSMark F. Adams       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
402*fc4362bfSMark F. Adams 
403*fc4362bfSMark F. Adams       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
404*fc4362bfSMark F. Adams       ierr = KSPSolve(eksp,bb,xx); CHKERRQ(ierr);
405*fc4362bfSMark F. Adams       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
406*fc4362bfSMark F. Adams 
407*fc4362bfSMark F. Adams       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
408*fc4362bfSMark F. Adams       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
409*fc4362bfSMark F. Adams       ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
410*fc4362bfSMark F. Adams 
411*fc4362bfSMark F. Adams       ierr = MatGetSize( Amat, &N1, &tt );         CHKERRQ(ierr);
412*fc4362bfSMark F. Adams       ierr = MatGetSize( Aarr[lidx+1], &N0, &tt );CHKERRQ(ierr);
413*fc4362bfSMark F. Adams       emin = emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */
414*fc4362bfSMark F. Adams     }
415*fc4362bfSMark F. Adams     ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
416*fc4362bfSMark F. Adams   }
417*fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
418*fc4362bfSMark F. Adams   ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr);
4195b89ad90SMark F. Adams   {
4205b89ad90SMark F. Adams     PetscBool galerkin;
4215b89ad90SMark F. Adams     ierr = PCMGGetGalerkin( pc,  &galerkin); CHKERRQ(ierr);
4225b89ad90SMark F. Adams     if(galerkin){
4235b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
4245b89ad90SMark F. Adams     }
4255b89ad90SMark F. Adams   }
4265ef31b24SMark F. Adams   {
4275ef31b24SMark F. Adams     char str[32];
4285ef31b24SMark F. Adams     sprintf(str,"MG Level %d (%d)",0,pc_gamg->m_Nlevels-1);
4295ef31b24SMark F. Adams     PetscLogStageRegister(str, &gamg_stages[fine_level]);
4305ef31b24SMark F. Adams   }
4315b89ad90SMark F. Adams   /* create coarse level and the interpolation between the levels */
432*fc4362bfSMark F. Adams   for (level=0,lidx=pc_gamg->m_Nlevels-1;
433*fc4362bfSMark F. Adams        level<fine_level;
434*fc4362bfSMark F. Adams        level++,lidx--){
4355b89ad90SMark F. Adams     level1 = level + 1;
4365ef31b24SMark F. Adams     ierr = PCMGSetInterpolation(pc,level1,Parr[lidx]);CHKERRQ(ierr);
4376c237d78SBarry Smith     if( !PETSC_TRUE ) {
43811e60469SMark F. Adams       PetscViewer viewer; char fname[32];
43911e60469SMark F. Adams       sprintf(fname,"Amat_%d.m",lidx);
44011e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
4415b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
44211e60469SMark F. Adams       ierr = MatView( Aarr[lidx], viewer ); CHKERRQ(ierr);
4435b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
4445b89ad90SMark F. Adams     }
4455ef31b24SMark F. Adams     ierr = MatDestroy( &Parr[lidx] );  CHKERRQ(ierr);
446a92563c5SMark F. Adams     {
4475b89ad90SMark F. Adams       KSP smoother;
4485b89ad90SMark F. Adams       ierr = PCMGGetSmoother(pc,level,&smoother); CHKERRQ(ierr);
4495b89ad90SMark F. Adams       ierr = KSPSetOperators( smoother, Aarr[lidx], Aarr[lidx], DIFFERENT_NONZERO_PATTERN );
4505b89ad90SMark F. Adams       CHKERRQ(ierr);
4515ef31b24SMark F. Adams       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
4525ef31b24SMark F. Adams     }
453a92563c5SMark F. Adams     ierr = MatDestroy( &Aarr[lidx] );  CHKERRQ(ierr);
4545ef31b24SMark F. Adams     {
4555ef31b24SMark F. Adams       char str[32];
4565ef31b24SMark F. Adams       sprintf(str,"MG Level %d (%d)",level+1,lidx-1);
4575ef31b24SMark F. Adams       PetscLogStageRegister(str, &gamg_stages[lidx-1]);
458a92563c5SMark F. Adams     }
4595b89ad90SMark F. Adams   }
4605b89ad90SMark F. Adams   { /* fine level (no P) */
4615b89ad90SMark F. Adams     KSP smoother;
4625b89ad90SMark F. Adams     ierr = PCMGGetSmoother(pc,fine_level,&smoother); CHKERRQ(ierr);
4635b89ad90SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], DIFFERENT_NONZERO_PATTERN );
4645b89ad90SMark F. Adams     CHKERRQ(ierr);
4655ef31b24SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
4665b89ad90SMark F. Adams   }
4675b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
4685b89ad90SMark F. Adams   pc->setupcalled = 0;
4695b89ad90SMark F. Adams   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
4705b89ad90SMark F. Adams   PetscFunctionReturn(0);
4715b89ad90SMark F. Adams }
4725b89ad90SMark F. Adams 
4735b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4745b89ad90SMark F. Adams /*
4755b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
4765b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
4775b89ad90SMark F. Adams 
4785b89ad90SMark F. Adams    Input Parameter:
4795b89ad90SMark F. Adams .  pc - the preconditioner context
4805b89ad90SMark F. Adams 
4815b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
4825b89ad90SMark F. Adams */
4835b89ad90SMark F. Adams #undef __FUNCT__
4845b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
4855b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
4865b89ad90SMark F. Adams {
4875b89ad90SMark F. Adams   PetscErrorCode  ierr;
4885b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
4895b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
4905b89ad90SMark F. Adams 
4915b89ad90SMark F. Adams   PetscFunctionBegin;
4925b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
4935b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
4945b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
4955b89ad90SMark F. Adams   PetscFunctionReturn(0);
4965b89ad90SMark F. Adams }
4975b89ad90SMark F. Adams 
4985b89ad90SMark F. Adams #undef __FUNCT__
4995b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
5005b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
5015b89ad90SMark F. Adams {
5025b89ad90SMark F. Adams   /* PetscErrorCode  ierr; */
5035b89ad90SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
5045b89ad90SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
5055b89ad90SMark F. Adams   /* MPI_Comm        comm = ((PetscObject)pc)->comm; */
5065b89ad90SMark F. Adams 
5075b89ad90SMark F. Adams   PetscFunctionBegin;
5085b89ad90SMark F. Adams   PetscFunctionReturn(0);
5095b89ad90SMark F. Adams }
5105b89ad90SMark F. Adams 
5115b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5125b89ad90SMark F. Adams /*
5135b89ad90SMark F. Adams  PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG
5145b89ad90SMark F. Adams 
5155b89ad90SMark F. Adams    Input Parameter:
5165b89ad90SMark F. Adams .  pc - the preconditioner context
5175b89ad90SMark F. Adams 
5185b89ad90SMark F. Adams    Application Interface Routine: PCCreate()
5195b89ad90SMark F. Adams 
5205b89ad90SMark F. Adams   */
5215b89ad90SMark F. Adams  /* MC
5225b89ad90SMark F. Adams      PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide
5235b89ad90SMark F. Adams        fine grid discretization matrix and coordinates on the fine grid.
5245b89ad90SMark F. Adams 
5255b89ad90SMark F. Adams    Options Database Key:
5265b89ad90SMark F. Adams    Multigrid options(inherited)
5275b89ad90SMark F. Adams +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
5285b89ad90SMark F. Adams .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
5295b89ad90SMark F. Adams .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
5305b89ad90SMark F. Adams    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
5315b89ad90SMark F. Adams    GAMG options:
5325b89ad90SMark F. Adams 
5335b89ad90SMark F. Adams    Level: intermediate
5345b89ad90SMark F. Adams   Concepts: multigrid
5355b89ad90SMark F. Adams 
5365b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
5375b89ad90SMark F. Adams            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
5385b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
5395b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
5405b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
5415b89ad90SMark F. Adams M */
5425b89ad90SMark F. Adams 
5435b89ad90SMark F. Adams EXTERN_C_BEGIN
5445b89ad90SMark F. Adams #undef __FUNCT__
5455b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
5465b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
5475b89ad90SMark F. Adams {
5485b89ad90SMark F. Adams   PetscErrorCode  ierr;
5495b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
5505b89ad90SMark F. Adams   PC_MG           *mg;
5515ef31b24SMark F. Adams   PetscClassId     cookie;
5525b89ad90SMark F. Adams 
5535b89ad90SMark F. Adams   PetscFunctionBegin;
5545b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
5555b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
5565b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
5575b89ad90SMark F. Adams 
5585b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
5595b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
5605b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
5615b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
5625b89ad90SMark F. Adams 
5635b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
5645b89ad90SMark F. Adams 
5655b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
5665b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
5675b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
5685b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
5695b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
5705b89ad90SMark F. Adams 
5715b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
5725b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
5735b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
5745b89ad90SMark F. Adams 					    PCSetCoordinates_GAMG);CHKERRQ(ierr);
5755ef31b24SMark F. Adams 
5765ef31b24SMark F. Adams   PetscClassIdRegister("GAMG Setup",&cookie);
5775ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-createProl", cookie, &gamg_setup_stages[SET1]);
5785ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-partLevel", cookie, &gamg_setup_stages[SET2]);
5795ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-MIS Graph", cookie, &gamg_setup_stages[SET3]);
5805ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-MIS-Agg", cookie, &gamg_setup_stages[SET4]);
5815ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-growSupp", cookie, &gamg_setup_stages[SET5]);
5825ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-tri-Prol", cookie, &gamg_setup_stages[SET6]);
5835ef31b24SMark F. Adams   PetscLogEventRegister("GAMG-find-prol", cookie, &gamg_setup_stages[FIND_V]);
5845ef31b24SMark F. Adams 
5855b89ad90SMark F. Adams   PetscFunctionReturn(0);
5865b89ad90SMark F. Adams }
5875b89ad90SMark F. Adams EXTERN_C_END
588