xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 302f38e8b74e96239826bae8d5bb9abb65139735)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4b45d2f2cSJed Brown #include "petsc-private/matimpl.h"
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
6b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
7f96513f1SMatthew G Knepley 
80cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
90cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
10b4fbaa2aSMark F. Adams #endif
110cbbd2e1SMark F. Adams 
120cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
130cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_AGG;
140cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_GEO;
150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGOptprol_AGG;
200cbbd2e1SMark F. Adams #endif
210cbbd2e1SMark F. Adams 
22b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
23b4fbaa2aSMark F. Adams 
24b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
250cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
26b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
27b4fbaa2aSMark F. Adams #endif
28f96513f1SMatthew G Knepley 
299d5b6da9SMark F. Adams static PetscFList GAMGList = 0;
309d5b6da9SMark F. Adams 
31d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
32d3d6bff4SMark F. Adams #undef __FUNCT__
33d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
34d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
35d3d6bff4SMark F. Adams {
36d3d6bff4SMark F. Adams   PetscErrorCode  ierr;
37d3d6bff4SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
38d3d6bff4SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
39d3d6bff4SMark F. Adams 
40d3d6bff4SMark F. Adams   PetscFunctionBegin;
419d5b6da9SMark F. Adams   if( pc_gamg->data != 0 ) { /* this should not happen, cleaned up in SetUp */
429d5b6da9SMark F. Adams     ierr = PetscFree(pc_gamg->data); CHKERRQ(ierr);
4358471d46SMark F. Adams   }
449d5b6da9SMark F. Adams   pc_gamg->data = 0; pc_gamg->data_sz = 0;
45d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
46d3d6bff4SMark F. Adams }
47d3d6bff4SMark F. Adams 
485b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
495b89ad90SMark F. Adams /*
50a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
51a147abb0SMark F. Adams      of active processors.
525b89ad90SMark F. Adams 
535b89ad90SMark F. Adams    Input Parameter:
549d5b6da9SMark F. Adams    . pc - parameters
559d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
569d5b6da9SMark F. Adams    . cbs - coarse block size
579d5b6da9SMark F. Adams    . isLast -
583530afc2SMark F. Adams    In/Output Parameter:
593530afc2SMark F. Adams    . a_P_inout - prolongation operator to the next level (k-1)
60afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
6111e60469SMark F. Adams    Output Parameter:
623530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
635b89ad90SMark F. Adams */
645cb416c2SMark F. Adams 
655b89ad90SMark F. Adams #undef __FUNCT__
668263b398SMark F. Adams #define __FUNCT__ "createLevel"
670cbbd2e1SMark F. Adams static PetscErrorCode createLevel( const PC pc,
689d5b6da9SMark F. Adams                             const Mat Amat_fine,
699d5b6da9SMark F. Adams                             const PetscInt cbs,
709d5b6da9SMark F. Adams                             const PetscBool isLast,
713530afc2SMark F. Adams                             Mat *a_P_inout,
72afc97cdcSMark F. Adams                             PetscMPIInt *a_nactive_proc,
73389730f3SMark F. Adams                             Mat *a_Amat_crs
7411e60469SMark F. Adams                             )
755b89ad90SMark F. Adams {
769d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
77486a8d0bSJed Brown   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
789d5b6da9SMark F. Adams   const PetscBool  repart = pc_gamg->repart;
799d5b6da9SMark F. Adams   const PetscInt   min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
805b89ad90SMark F. Adams   PetscErrorCode   ierr;
81038e3b61SMark F. Adams   Mat              Cmat,Pnew,Pold=*a_P_inout;
8211e60469SMark F. Adams   IS               new_indices,isnum;
839d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)Amat_fine)->comm;
845a9b9e01SMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive = *a_nactive_proc;
855a9b9e01SMark F. Adams   PetscInt         neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0,rfactor;
865b89ad90SMark F. Adams 
875b89ad90SMark F. Adams   PetscFunctionBegin;
8811e60469SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
8911e60469SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
90f6536408SMark F. Adams 
9111e60469SMark F. Adams   /* RAP */
9296434bc1SMark F. Adams #ifdef USE_R
9396434bc1SMark F. Adams   /* make R wih brute force for now */
9496434bc1SMark F. Adams   ierr = MatTranspose( Pold, Pnew );
9596434bc1SMark F. Adams   ierr = MatDestroy( &Pold );  CHKERRQ(ierr);
969d5b6da9SMark F. Adams   ierr = MatRARt( Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
9796434bc1SMark F. Adams   Pold = Pnew;
9896434bc1SMark F. Adams #else
999d5b6da9SMark F. Adams   ierr = MatPtAP( Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
10096434bc1SMark F. Adams #endif
1019d5b6da9SMark F. Adams   ierr = MatSetBlockSize( Cmat, cbs );      CHKERRQ(ierr);
102038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
1039d5b6da9SMark F. Adams   ncrs0 = (Iend0-Istart0)/cbs; assert((Iend0-Istart0)%cbs == 0);
104038e3b61SMark F. Adams 
105f852f58cSMark F. Adams   /* get number of PEs to make active, reduce */
106f852f58cSMark F. Adams   ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
10781708dccSMark F. Adams   new_npe = (PetscMPIInt)((float)neq/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
108f852f58cSMark F. Adams   if( new_npe == 0 || neq < coarse_max ) new_npe = 1;
1095a9b9e01SMark F. Adams   else if (new_npe >= nactive ) new_npe = nactive; /* no change, rare */
1109d5b6da9SMark F. Adams   if( isLast ) new_npe = 1;
111f852f58cSMark F. Adams 
1125a9b9e01SMark F. Adams   if( !repart && new_npe==nactive ) {
113a147abb0SMark F. Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction */
11422063be5SMark F. Adams   }
11522063be5SMark F. Adams   else {
11622063be5SMark F. Adams     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
117c8b0795cSMark F. Adams     const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,data_sz=ndata_rows*ndata_cols;
118c8b0795cSMark F. Adams     const PetscInt  stride0=ncrs0*pc_gamg->data_cell_rows;
1195a9b9e01SMark F. Adams     PetscInt        *counts,is_sz,*newproc_idx,ii,jj,kk,strideNew,*tidx,inpe,targetPE;
1205a9b9e01SMark F. Adams     IS              isnewproc;
1215a9b9e01SMark F. Adams     VecScatter      vecscat;
12222063be5SMark F. Adams     PetscScalar    *array;
12322063be5SMark F. Adams     Vec             src_crd, dest_crd;
1245a9b9e01SMark F. Adams     IS              isscat;
125e33ef3b1SMark F. Adams 
126fd3c6afaSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
1270cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1280cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
129b4fbaa2aSMark F. Adams #endif
1305a9b9e01SMark F. Adams     if( repart ) {
1315a9b9e01SMark F. Adams       /* create sub communicator  */
1325a9b9e01SMark F. Adams       Mat             adj;
1335a9b9e01SMark F. Adams 
134b43b03e9SMark F. Adams       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,neq);
1355a9b9e01SMark F. Adams 
1365a9b9e01SMark F. Adams       /* MatPartitioningApply call MatConvert, which is collective */
1379d5b6da9SMark F. Adams       if( cbs == 1 ) {
138038e3b61SMark F. Adams 	ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
139eb07cef2SMark F. Adams       }
140eb07cef2SMark F. Adams       else{
141038e3b61SMark F. Adams 	/* make a scalar matrix to partition */
142eb07cef2SMark F. Adams 	Mat tMat;
14358471d46SMark F. Adams 	PetscInt ncols,jj,Ii;
144b4fbaa2aSMark F. Adams 	const PetscScalar *vals;
145b4fbaa2aSMark F. Adams 	const PetscInt *idx;
1465f8cf99dSMark F. Adams 	PetscInt *d_nnz, *o_nnz;
1479057884aSMark F. Adams 	static PetscInt llev = 0;
148b4fbaa2aSMark F. Adams 
14958471d46SMark F. Adams 	ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
1505f8cf99dSMark F. Adams 	ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr);
1519d5b6da9SMark F. Adams 	for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += cbs, jj++ ) {
15258471d46SMark F. Adams 	  ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
1539d5b6da9SMark F. Adams 	  d_nnz[jj] = ncols/cbs;
1549d5b6da9SMark F. Adams 	  o_nnz[jj] = ncols/cbs;
15558471d46SMark F. Adams 	  ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
1565f8cf99dSMark F. Adams 	  if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
1579d5b6da9SMark F. Adams 	  if( o_nnz[jj] > (neq/cbs-ncrs0) ) o_nnz[jj] = neq/cbs-ncrs0;
15858471d46SMark F. Adams 	}
1596876a03eSMark F. Adams 
16069b1f4b7SBarry Smith 	ierr = MatCreateAIJ( wcomm, ncrs0, ncrs0,
161eb07cef2SMark F. Adams                              PETSC_DETERMINE, PETSC_DETERMINE,
1625f8cf99dSMark F. Adams                              0, d_nnz, 0, o_nnz,
163eb07cef2SMark F. Adams                              &tMat );
1646876a03eSMark F. Adams 	CHKERRQ(ierr);
16558471d46SMark F. Adams 	ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
1665f8cf99dSMark F. Adams 	ierr = PetscFree( o_nnz ); CHKERRQ(ierr);
167eb07cef2SMark F. Adams 
16822063be5SMark F. Adams 	for ( ii = Istart0; ii < Iend0; ii++ ) {
1699d5b6da9SMark F. Adams 	  PetscInt dest_row = ii/cbs;
17022063be5SMark F. Adams 	  ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
171eb07cef2SMark F. Adams 	  for( jj = 0 ; jj < ncols ; jj++ ){
1729d5b6da9SMark F. Adams 	    PetscInt dest_col = idx[jj]/cbs;
173eb07cef2SMark F. Adams 	    PetscScalar v = 1.0;
174eb07cef2SMark F. Adams 	    ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
175eb07cef2SMark F. Adams 	  }
17622063be5SMark F. Adams 	  ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
177eb07cef2SMark F. Adams 	}
178eb07cef2SMark F. Adams 	ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
179eb07cef2SMark F. Adams 	ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
180eb07cef2SMark F. Adams 
181b4fbaa2aSMark F. Adams 	if( llev++ == -1 ) {
182b4fbaa2aSMark F. Adams 	  PetscViewer viewer; char fname[32];
1834704b1adSJed Brown 	  ierr = PetscSNPrintf(fname,sizeof fname,"part_mat_%D.mat",llev);CHKERRQ(ierr);
184b4fbaa2aSMark F. Adams 	  PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
185b4fbaa2aSMark F. Adams 	  ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
186b4fbaa2aSMark F. Adams 	  ierr = PetscViewerDestroy( &viewer );
187b4fbaa2aSMark F. Adams 	}
188b4fbaa2aSMark F. Adams 
189eb07cef2SMark F. Adams 	ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
190eb07cef2SMark F. Adams 
191eb07cef2SMark F. Adams 	ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
192eb07cef2SMark F. Adams       }
193f150b916SMark F. Adams 
1945a9b9e01SMark F. Adams       { /* partition: get newproc_idx, set is_sz */
1955a9b9e01SMark F. Adams 	char prefix[256];
1965a9b9e01SMark F. Adams 	const char *pcpre;
197b4fbaa2aSMark F. Adams 	const PetscInt *is_idx;
198b4fbaa2aSMark F. Adams 	MatPartitioning  mpart;
199a4b7d37bSMark F. Adams 	IS proc_is;
2002f03bc48SMark F. Adams 
2015a9b9e01SMark F. Adams 	ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr);
2025ef31b24SMark F. Adams 	ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
2039d5b6da9SMark F. Adams 	ierr = PCGetOptionsPrefix( pc, &pcpre );CHKERRQ(ierr);
20459a0be82SJed Brown 	ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
20559a0be82SJed Brown 	ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
20611e60469SMark F. Adams 	ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
2075ef31b24SMark F. Adams 	ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
208a4b7d37bSMark F. Adams 	ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr);
20911e60469SMark F. Adams 	ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
2105a9b9e01SMark F. Adams 
2115ef31b24SMark F. Adams 	/* collect IS info */
212a4b7d37bSMark F. Adams 	ierr = ISGetLocalSize( proc_is, &is_sz );       CHKERRQ(ierr);
2139d5b6da9SMark F. Adams 	ierr = PetscMalloc( cbs*is_sz*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr);
214a4b7d37bSMark F. Adams 	ierr = ISGetIndices( proc_is, &is_idx );        CHKERRQ(ierr);
215a147abb0SMark F. Adams 	NN = 1; /* bring to "front" of machine */
216a147abb0SMark F. Adams 	/*NN = npe/new_npe;*/ /* spread partitioning across machine */
217b4fbaa2aSMark F. Adams 	for( kk = jj = 0 ; kk < is_sz ; kk++ ){
2189d5b6da9SMark F. Adams 	  for( ii = 0 ; ii < cbs ; ii++, jj++ ){
2198263b398SMark F. Adams 	    newproc_idx[jj] = is_idx[kk] * NN; /* distribution */
220eb07cef2SMark F. Adams 	  }
2215ef31b24SMark F. Adams 	}
222a4b7d37bSMark F. Adams 	ierr = ISRestoreIndices( proc_is, &is_idx );     CHKERRQ(ierr);
223a4b7d37bSMark F. Adams 	ierr = ISDestroy( &proc_is );                    CHKERRQ(ierr);
224b4fbaa2aSMark F. Adams 
2259d5b6da9SMark F. Adams 	is_sz *= cbs;
2265ef31b24SMark F. Adams       }
2275ef31b24SMark F. Adams       ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
2285a9b9e01SMark F. Adams 
2298263b398SMark F. Adams       ierr = ISCreateGeneral( wcomm, is_sz, newproc_idx, PETSC_COPY_VALUES, &isnewproc );
2305a9b9e01SMark F. Adams       CHKERRQ(ierr);
2318263b398SMark F. Adams       if( newproc_idx != 0 ) {
2328263b398SMark F. Adams 	ierr = PetscFree( newproc_idx );  CHKERRQ(ierr);
2335ef31b24SMark F. Adams       }
2345a9b9e01SMark F. Adams     }
2355a9b9e01SMark F. Adams     else { /* simple aggreagtion of parts */
2365a9b9e01SMark F. Adams       /* find factor */
2375a9b9e01SMark F. Adams       if( new_npe == 1 ) rfactor = npe; /* easy */
2385a9b9e01SMark F. Adams       else {
2395a9b9e01SMark F. Adams 	PetscReal best_fact = 0.;
2405a9b9e01SMark F. Adams 	jj = -1;
2415a9b9e01SMark F. Adams 	for( kk = 1 ; kk <= npe ; kk++ ){
2425a9b9e01SMark F. Adams 	  if( npe%kk==0 ) { /* a candidate */
2435a9b9e01SMark F. Adams 	    PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe;
2445a9b9e01SMark F. Adams 	    if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2455a9b9e01SMark F. Adams 	    if( fact > best_fact ) {
2465a9b9e01SMark F. Adams 	      best_fact = fact; jj = kk;
2475a9b9e01SMark F. Adams 	    }
2485a9b9e01SMark F. Adams 	  }
2495a9b9e01SMark F. Adams 	}
2505a9b9e01SMark F. Adams 	if(jj!=-1) rfactor = jj;
2515a9b9e01SMark F. Adams 	else rfactor = 1; /* prime? */
2525a9b9e01SMark F. Adams       }
2535a9b9e01SMark F. Adams       new_npe = npe/rfactor;
2545a9b9e01SMark F. Adams 
2555a9b9e01SMark F. Adams       if( new_npe==nactive ) {
2565a9b9e01SMark F. Adams 	*a_Amat_crs = Cmat; /* output - no repartitioning or reduction */
2575a9b9e01SMark F. Adams 	ierr = PetscFree( counts );  CHKERRQ(ierr);
2585a9b9e01SMark F. Adams 	*a_nactive_proc = new_npe; /* output */
259b43b03e9SMark F. Adams 	if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq=%d\n",mype,__FUNCT__,new_npe,neq);
2605a9b9e01SMark F. Adams 	PetscFunctionReturn(0);
2615a9b9e01SMark F. Adams       }
2625a9b9e01SMark F. Adams 
2635a9b9e01SMark F. Adams       /* reduce - isnewproc */
2645a9b9e01SMark F. Adams       targetPE = mype/rfactor;
2655a9b9e01SMark F. Adams 
266b43b03e9SMark F. Adams       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s aggregate processors: npe: %d --> %d, neq=%d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,neq);
2679d5b6da9SMark F. Adams       is_sz = ncrs0*cbs;
2685a9b9e01SMark F. Adams       ierr = ISCreateStride( wcomm, is_sz, targetPE, 0, &isnewproc );
2695a9b9e01SMark F. Adams       CHKERRQ(ierr);
2705a9b9e01SMark F. Adams     }
271e33ef3b1SMark F. Adams 
27211e60469SMark F. Adams     /*
27311e60469SMark F. Adams       Create an index set from the isnewproc index set to indicate the mapping TO
27411e60469SMark F. Adams     */
27511e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
27611e60469SMark F. Adams     /*
27711e60469SMark F. Adams       Determine how many elements are assigned to each processor
27811e60469SMark F. Adams     */
279fd3c6afaSMark F. Adams     inpe = npe;
280fd3c6afaSMark F. Adams     ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr);
28111e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
2829d5b6da9SMark F. Adams     ncrs_new = counts[mype]/cbs;
2839d5b6da9SMark F. Adams     strideNew = ncrs_new*ndata_rows;
2840cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2850cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
286b4fbaa2aSMark F. Adams #endif
28722063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
28811e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
289d3d6bff4SMark F. Adams     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
290470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
29111e60469SMark F. Adams     /*
2929d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
2939d5b6da9SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cbs'.
29411e60469SMark F. Adams      */
29592a756f0SMark F. Adams     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
29611e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
297038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
2989d5b6da9SMark F. Adams       PetscInt id = idx[ii*cbs]/cbs; /* get node back */
299d3d6bff4SMark F. Adams       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
30011e60469SMark F. Adams     }
301038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
302d3d6bff4SMark F. Adams     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
30311e60469SMark F. Adams     CHKERRQ(ierr);
30492a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
30511e60469SMark F. Adams     /*
30611e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
30711e60469SMark F. Adams      */
308d3d6bff4SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
3099d5b6da9SMark F. Adams     for( jj=0; jj<ndata_cols ; jj++ ) {
310d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs0 ; ii++) {
3119d5b6da9SMark F. Adams 	for( kk=0; kk<ndata_rows ; kk++ ) {
3129d5b6da9SMark F. Adams 	  PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*ndata_cols + jj;
313c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
314676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
315d3d6bff4SMark F. Adams 	}
316038e3b61SMark F. Adams       }
317eb07cef2SMark F. Adams     }
318eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
319eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
32011e60469SMark F. Adams     /*
32111e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
32211e60469SMark F. Adams       to the correct processor
32311e60469SMark F. Adams     */
32411e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
32511e60469SMark F. Adams     CHKERRQ(ierr);
32611e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
32711e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32811e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32911e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
33011e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
33111e60469SMark F. Adams     /*
33211e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
33311e60469SMark F. Adams     */
334c8b0795cSMark F. Adams     ierr = PetscFree( pc_gamg->data );    CHKERRQ(ierr);
335c8b0795cSMark F. Adams     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), &pc_gamg->data );    CHKERRQ(ierr);
336c8b0795cSMark F. Adams     pc_gamg->data_sz = data_sz*ncrs_new;
33711e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
3389d5b6da9SMark F. Adams     for( jj=0; jj<ndata_cols ; jj++ ) {
339d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs_new ; ii++) {
3409d5b6da9SMark F. Adams 	for( kk=0; kk<ndata_rows ; kk++ ) {
3419d5b6da9SMark F. Adams 	  PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*ndata_cols + jj;
342c8b0795cSMark F. Adams 	  pc_gamg->data[ix] = PetscRealPart(array[jx]);
343d3d6bff4SMark F. Adams 	  array[jx] = 1.e300;
344d3d6bff4SMark F. Adams 	}
345038e3b61SMark F. Adams       }
346038e3b61SMark F. Adams     }
34711e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
34811e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
3490cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3500cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
351ed3f9983SMark F. Adams #endif
35211e60469SMark F. Adams     /*
35311e60469SMark F. Adams       Invert for MatGetSubMatrix
35411e60469SMark F. Adams     */
3559d5b6da9SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*cbs, &new_indices ); CHKERRQ(ierr);
35611e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
35711e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
3580cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3590cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3600cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
361ed3f9983SMark F. Adams #endif
36211e60469SMark F. Adams     /* A_crs output */
363038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
36411e60469SMark F. Adams     CHKERRQ(ierr);
365eb07cef2SMark F. Adams 
366038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
367e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
3689d5b6da9SMark F. Adams     ierr = MatSetBlockSize( Cmat, cbs );      CHKERRQ(ierr);
3690cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3700cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
371ed3f9983SMark F. Adams #endif
37211e60469SMark F. Adams     /* prolongator */
37311e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
37411e60469SMark F. Adams     {
37511e60469SMark F. Adams       IS findices;
3760cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3770cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
378ed3f9983SMark F. Adams #endif
37911e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
38096434bc1SMark F. Adams #ifdef USE_R
38196434bc1SMark F. Adams       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
38296434bc1SMark F. Adams #else
38311e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
38496434bc1SMark F. Adams #endif
38511e60469SMark F. Adams       CHKERRQ(ierr);
38611e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
3870cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3880cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
389ed3f9983SMark F. Adams #endif
39011e60469SMark F. Adams     }
3913530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
392a147abb0SMark F. Adams     *a_P_inout = Pnew; /* output - repartitioned */
39311e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
39492a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
395e33ef3b1SMark F. Adams   }
3965b89ad90SMark F. Adams 
3975a9b9e01SMark F. Adams   *a_nactive_proc = new_npe; /* output */
3985a9b9e01SMark F. Adams 
399c8b0795cSMark F. Adams   if( !PETSC_TRUE ) {
400c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
401c8b0795cSMark F. Adams     if(llev==0) {
402c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
403c8b0795cSMark F. Adams       PetscViewerASCIIOpen(wcomm,fname,&viewer);
404c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
405c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer ); CHKERRQ(ierr);
406c8b0795cSMark F. Adams       ierr = PetscViewerDestroy( &viewer );
407c8b0795cSMark F. Adams     }
408c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
409c8b0795cSMark F. Adams     PetscViewerASCIIOpen(wcomm,fname,&viewer);
410c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
411c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer ); CHKERRQ(ierr);
412c8b0795cSMark F. Adams     ierr = PetscViewerDestroy( &viewer );
413c8b0795cSMark F. Adams   }
414c8b0795cSMark F. Adams 
4155b89ad90SMark F. Adams   PetscFunctionReturn(0);
4165b89ad90SMark F. Adams }
4175b89ad90SMark F. Adams 
4185b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4195b89ad90SMark F. Adams /*
4205b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4215b89ad90SMark F. Adams                     by setting data structures and options.
4225b89ad90SMark F. Adams 
4235b89ad90SMark F. Adams    Input Parameter:
4245b89ad90SMark F. Adams .  pc - the preconditioner context
4255b89ad90SMark F. Adams 
4265b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4275b89ad90SMark F. Adams 
4285b89ad90SMark F. Adams    Notes:
4295b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4305b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4315b89ad90SMark F. Adams */
4325b89ad90SMark F. Adams #undef __FUNCT__
4335b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4349d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC pc )
4355b89ad90SMark F. Adams {
4365b89ad90SMark F. Adams   PetscErrorCode  ierr;
4379d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
4385b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
4392adcac29SMark F. Adams   Mat              Pmat = pc->pmat;
440ffc955d6SMark F. Adams   PetscInt         fine_level,level,level1,M,N,bs,nloc,lidx,Istart,Iend,nASMBlocksArr[GAMG_MAXLEVELS];
4419d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
4423530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
443587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
444c8b0795cSMark F. Adams   PetscReal        emaxs[GAMG_MAXLEVELS];
445a94c3b12SMark F. Adams   IS              *ASMLocalIDsArr[GAMG_MAXLEVELS],removedEqs[GAMG_MAXLEVELS];
446a94c3b12SMark F. Adams   PetscInt         level_bs[GAMG_MAXLEVELS];
4470534fa05SHong Zhang   PetscLogDouble   nnz0=0,nnztot=0;
448737a81a9SMark F. Adams   MatInfo          info;
449*302f38e8SMark F. Adams   PetscBool        stokes = PETSC_FALSE;
4505ef31b24SMark F. Adams 
4515b89ad90SMark F. Adams   PetscFunctionBegin;
45284d3f75bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
45384d3f75bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
45484d3f75bSMark F. Adams   if( pc_gamg->setup_count++ > 0 ) {
45584d3f75bSMark F. Adams     PC_MG_Levels   **mglevels = mg->levels;
45603a628feSMark F. Adams     /* just do Galerkin grids */
45758471d46SMark F. Adams     Mat B,dA,dB;
458d5280255SMark F. Adams     assert(pc->setupcalled);
45958471d46SMark F. Adams 
4609d5b6da9SMark F. Adams     if( pc_gamg->Nlevels > 1 ) {
46158471d46SMark F. Adams       /* currently only handle case where mat and pmat are the same on coarser levels */
4629d5b6da9SMark F. Adams       ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
46358471d46SMark F. Adams       /* (re)set to get dirty flag */
4649d5b6da9SMark F. Adams       ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
46558471d46SMark F. Adams 
4669d5b6da9SMark F. Adams       for (level=pc_gamg->Nlevels-2; level>-1; level--) {
46703a628feSMark F. Adams         /* the first time through the matrix structure has changed from repartitioning */
46884d3f75bSMark F. Adams         if( pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
46903a628feSMark F. Adams           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
470084a8fe3SJed Brown           ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
47103a628feSMark F. Adams           mglevels[level]->A = B;
47203a628feSMark F. Adams         }
47303a628feSMark F. Adams         else {
474084a8fe3SJed Brown           ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
47558471d46SMark F. Adams           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
47603a628feSMark F. Adams         }
47758471d46SMark F. Adams         ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
47858471d46SMark F. Adams         dB = B;
47958471d46SMark F. Adams         }
4805f8cf99dSMark F. Adams     }
481d5280255SMark F. Adams 
482d5280255SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
483d5280255SMark F. Adams 
484d5280255SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
485d5280255SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
486d5280255SMark F. Adams 
48758471d46SMark F. Adams     PetscFunctionReturn(0);
488eb07cef2SMark F. Adams   }
48984d3f75bSMark F. Adams   assert(pc->setupcalled == 0);
490f6536408SMark F. Adams 
491*302f38e8SMark F. Adams   /* deal with Stokes */
492*302f38e8SMark F. Adams   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
493eb07cef2SMark F. Adams 
494*302f38e8SMark F. Adams   if( pc_gamg->data == 0 ) {
4959d5b6da9SMark F. Adams     if( !pc_gamg->createdefaultdata ){
496*302f38e8SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_LIB,"'createdefaultdata' not set?!?! need to support NULL data!!!");
497*302f38e8SMark F. Adams     }
498*302f38e8SMark F. Adams     if( stokes ) {
499*302f38e8SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_LIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
500eb07cef2SMark F. Adams     }
501b8cd405aSMark F. Adams     ierr = pc_gamg->createdefaultdata( pc, Pmat ); CHKERRQ(ierr);
5029d5b6da9SMark F. Adams   }
503038e3b61SMark F. Adams 
504*302f38e8SMark F. Adams   /* get basic dims */
505*302f38e8SMark F. Adams   if( stokes ) {
506*302f38e8SMark F. Adams     bs = pc_gamg->data_cell_rows;
507*302f38e8SMark F. Adams     nloc = pc_gamg->data_sz/pc_gamg->data_cell_cols/bs;
508*302f38e8SMark F. Adams     /* M,N,Iend,Istart,Pmat -- TODO */
509*302f38e8SMark F. Adams     ierr = MatGetBlockSize( Pmat, &M ); CHKERRQ(ierr); assert(M==bs);
510*302f38e8SMark F. Adams     ierr = MatGetOwnershipRange( Pmat, &Istart, &Iend ); CHKERRQ(ierr); assert((Iend-Istart)/bs == nloc);
511*302f38e8SMark F. Adams     ierr = MatGetSize( Pmat, &M, &N );CHKERRQ(ierr);
512*302f38e8SMark F. Adams   }
513*302f38e8SMark F. Adams   else {
514*302f38e8SMark F. Adams     ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr);
515*302f38e8SMark F. Adams     ierr = MatGetSize( Pmat, &M, &N );CHKERRQ(ierr);
516*302f38e8SMark F. Adams     ierr = MatGetOwnershipRange( Pmat, &Istart, &Iend ); CHKERRQ(ierr);
517*302f38e8SMark F. Adams     nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
518*302f38e8SMark F. Adams   }
519eb07cef2SMark F. Adams   /* Get A_i and R_i */
520c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
5212adcac29SMark F. Adams     if(pc_gamg->verbose==1) ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);
5222adcac29SMark F. Adams     else ierr =  MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);
523b2a4f308SMark F. Adams     CHKERRQ(ierr);
524b2a4f308SMark F. Adams     nnz0 = info.nz_used;
525b2a4f308SMark F. Adams     nnztot = info.nz_used;
526b43b03e9SMark F. Adams     PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
527c8b0795cSMark F. Adams                 mype,__FUNCT__,0,N,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
528b2a4f308SMark F. Adams                 (int)(nnz0/(PetscReal)N),npe);
529c8b0795cSMark F. Adams   }
53084d3f75bSMark F. Adams 
5318f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
5329d5b6da9SMark F. Adams         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
5330205a208SMark F. Adams         level++ ){
5345b89ad90SMark F. Adams     level1 = level + 1;
5350cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
536b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
537b4fbaa2aSMark F. Adams #endif
5380cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5390cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr);
540b4fbaa2aSMark F. Adams #endif
541c8b0795cSMark F. Adams     { /* construct prolongator */
542725b86d8SJed Brown       Mat Gmat;
5430cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
544c8b0795cSMark F. Adams 
545a94c3b12SMark F. Adams       level_bs[level] = bs;
546b43b03e9SMark F. Adams       ierr = pc_gamg->graph( pc, Aarr[level], &Gmat ); CHKERRQ(ierr);
5470cbbd2e1SMark F. Adams       ierr = pc_gamg->coarsen( pc, &Gmat, &agg_lists ); CHKERRQ(ierr);
5480cbbd2e1SMark F. Adams       ierr = pc_gamg->prolongator( pc, Aarr[level], Gmat, agg_lists, &Parr[level1] );
5495b89ad90SMark F. Adams       CHKERRQ(ierr);
550c8b0795cSMark F. Adams 
551c8b0795cSMark F. Adams       if( Parr[level1] ){
5529d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
553c8b0795cSMark F. Adams         if( pc_gamg->col_bs_id != -1 ){
5549d5b6da9SMark F. Adams           PetscBool flag;
5559d5b6da9SMark F. Adams           ierr = PetscObjectComposedDataGetInt( (PetscObject)Parr[level1], pc_gamg->col_bs_id, bs, flag );
556c8b0795cSMark F. Adams           CHKERRQ( ierr ); assert(flag);
5579d5b6da9SMark F. Adams         }
558c8b0795cSMark F. Adams       }
559c8b0795cSMark F. Adams 
560c8b0795cSMark F. Adams       if( pc_gamg->optprol && Parr[level1] ){
561c8b0795cSMark F. Adams         /* smooth */
562c8b0795cSMark F. Adams         ierr = pc_gamg->optprol( pc, Aarr[level], &Parr[level1] ); CHKERRQ(ierr);
563c8b0795cSMark F. Adams       }
564c8b0795cSMark F. Adams 
565c8b0795cSMark F. Adams       ierr = MatDestroy( &Gmat );      CHKERRQ(ierr);
566ffc955d6SMark F. Adams 
567ffc955d6SMark F. Adams       if ( pc_gamg->use_aggs_in_gasm ) {
568a94c3b12SMark F. Adams         ierr = PetscCDGetASMBlocks(agg_lists, level_bs[level], &nASMBlocksArr[level], &ASMLocalIDsArr[level]);  CHKERRQ(ierr);
569ffc955d6SMark F. Adams       }
570ffc955d6SMark F. Adams 
571a94c3b12SMark F. Adams       ierr = PetscCDGetRemovedIS( agg_lists, &removedEqs[level] );  CHKERRQ(ierr);
572a94c3b12SMark F. Adams 
57341b27cdeSMark F. Adams       ierr = PetscCDDestroy( agg_lists );  CHKERRQ(ierr);
574c8b0795cSMark F. Adams     }
5750cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5760cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
577c8b0795cSMark F. Adams #endif
5789d5b6da9SMark F. Adams     /* cache eigen estimate */
5799d5b6da9SMark F. Adams     if( pc_gamg->emax_id != -1 ){
5809d5b6da9SMark F. Adams       PetscBool flag;
581f73f8d2cSSatish Balay       ierr = PetscObjectComposedDataGetReal( (PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag );
5829d5b6da9SMark F. Adams       CHKERRQ( ierr );
5839d5b6da9SMark F. Adams       if( !flag ) emaxs[level] = -1.;
5849d5b6da9SMark F. Adams     }
5859d5b6da9SMark F. Adams     else emaxs[level] = -1.;
5862adcac29SMark F. Adams     if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
587c8b0795cSMark F. Adams     if( !Parr[level1] ) {
588b43b03e9SMark F. Adams       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);
589c8b0795cSMark F. Adams       break;
590c8b0795cSMark F. Adams     }
5910cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5920cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
593b4fbaa2aSMark F. Adams #endif
594c8b0795cSMark F. Adams     ierr = createLevel( pc, Aarr[level], bs,
5959d5b6da9SMark F. Adams                         (PetscBool)(level==pc_gamg->Nlevels-2),
596c8b0795cSMark F. Adams                         &Parr[level1], &nactivepe, &Aarr[level1] );
5973530afc2SMark F. Adams     CHKERRQ(ierr);
5980cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5990cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
600b4fbaa2aSMark F. Adams #endif
6013530afc2SMark F. Adams     ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
602c8b0795cSMark F. Adams     if (pc_gamg->verbose){
6030cbbd2e1SMark F. Adams       PetscInt NN = M;
6040cbbd2e1SMark F. Adams       if(pc_gamg->verbose==1) {
6050cbbd2e1SMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);
6060cbbd2e1SMark F. Adams         ierr = MatGetLocalSize( Aarr[level1], &NN, &N );CHKERRQ(ierr);
6070cbbd2e1SMark F. Adams       }
6080cbbd2e1SMark F. Adams       else ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info);
6090cbbd2e1SMark F. Adams       CHKERRQ(ierr);
6100cbbd2e1SMark F. Adams       nnztot += info.nz_used;
611b43b03e9SMark F. Adams       PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
612c8b0795cSMark F. Adams                   mype,__FUNCT__,(int)level1,N,pc_gamg->data_cell_cols,
6130cbbd2e1SMark F. Adams                   (int)(info.nz_used/(PetscReal)NN), nactivepe );
6140cbbd2e1SMark F. Adams       CHKERRQ(ierr);
615c8b0795cSMark F. Adams     }
61681708dccSMark F. Adams     /* stop if one node */
617c8b0795cSMark F. Adams     if( M/pc_gamg->data_cell_cols < 2 ) {
61881708dccSMark F. Adams       level++;
61981708dccSMark F. Adams       break;
62081708dccSMark F. Adams     }
621c8b0795cSMark F. Adams     /* Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; */
622c8b0795cSMark F. Adams     /* v = 1.e-10; /\* LU factor has hard wired numbers for small diags so this needs to match (yuk) *\/ */
623c8b0795cSMark F. Adams     /* ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); */
624c8b0795cSMark F. Adams     /* nloceq = Iend-Istart; */
625c8b0795cSMark F. Adams     /* ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr); */
626c8b0795cSMark F. Adams     /* ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr); */
627c8b0795cSMark F. Adams     /* ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr); */
628c8b0795cSMark F. Adams     /* for(kk=0;kk<nloceq;kk++){ */
629c8b0795cSMark F. Adams     /*   if(data_arr[kk]==0.0) { */
630c8b0795cSMark F. Adams     /*     id = kk + Istart; */
631c8b0795cSMark F. Adams     /*     ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); */
632c8b0795cSMark F. Adams     /*     CHKERRQ(ierr); */
633c8b0795cSMark F. Adams     /*     PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1); */
634c8b0795cSMark F. Adams     /*   } */
635c8b0795cSMark F. Adams     /* } */
636c8b0795cSMark F. Adams     /* ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); */
637c8b0795cSMark F. Adams     /* ierr = VecDestroy( &diag );                CHKERRQ(ierr); */
63884d3f75bSMark F. Adams     /* ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
63984d3f75bSMark F. Adams     /* ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
640b4fbaa2aSMark F. Adams 
6410cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
642b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
643b4fbaa2aSMark F. Adams #endif
644c8b0795cSMark F. Adams   } /* levels */
645c8b0795cSMark F. Adams 
646c8b0795cSMark F. Adams   if( pc_gamg->data ) {
647c8b0795cSMark F. Adams     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
648c8b0795cSMark F. Adams     pc_gamg->data = 0;
6495b89ad90SMark F. Adams   }
650c8b0795cSMark F. Adams 
6510cbbd2e1SMark F. Adams   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
6529d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6535b89ad90SMark F. Adams   fine_level = level;
6549d5b6da9SMark F. Adams   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
6555b89ad90SMark F. Adams 
65684d3f75bSMark F. Adams   /* simple setup */
65784d3f75bSMark F. Adams   if( !PETSC_TRUE ){
65884d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
65984d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
66084d3f75bSMark F. Adams          lidx<fine_level;
66184d3f75bSMark F. Adams          lidx++, level--){
66284d3f75bSMark F. Adams       ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr);
66384d3f75bSMark F. Adams       ierr = KSPSetOperators( mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );CHKERRQ(ierr);
66484d3f75bSMark F. Adams       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
66584d3f75bSMark F. Adams       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
66684d3f75bSMark F. Adams     }
66784d3f75bSMark F. Adams     ierr = KSPSetOperators( mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
66884d3f75bSMark F. Adams 
66984d3f75bSMark F. Adams     ierr = PCSetUp_MG( pc );  CHKERRQ( ierr );
67084d3f75bSMark F. Adams   }
67184d3f75bSMark F. Adams   else if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if one level */
672d5280255SMark F. Adams     /* set default smoothers & set operators */
6739d5b6da9SMark F. Adams     for ( lidx = 1, level = pc_gamg->Nlevels-2;
674587fa25dSMark F. Adams           lidx <= fine_level;
675587fa25dSMark F. Adams           lidx++, level--) {
676ffc955d6SMark F. Adams       KSP smoother;
677ffc955d6SMark F. Adams       PC subpc;
678f6536408SMark F. Adams       /* set defaults */
6799d5b6da9SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
6805b89ad90SMark F. Adams       ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
681f6536408SMark F. Adams       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
682ffc955d6SMark F. Adams 
683ffc955d6SMark F. Adams       /* override defaults and command line args (!) */
684ffc955d6SMark F. Adams       if ( pc_gamg->use_aggs_in_gasm ) {
6852d3561bbSSatish Balay         PetscInt sz;
6862d3561bbSSatish Balay         IS *is;
687ffc955d6SMark F. Adams         if( PETSC_FALSE ){
688ffc955d6SMark F. Adams           PetscInt ii;
689ffc955d6SMark F. Adams           PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s Nblocks=%d\n",mype,__FUNCT__,nASMBlocksArr[level]);
690ffc955d6SMark F. Adams           for(ii=0;ii<nASMBlocksArr[level];ii++){
691ffc955d6SMark F. Adams             ISView(ASMLocalIDsArr[level][ii],PETSC_VIEWER_STDOUT_SELF);
692ffc955d6SMark F. Adams           }
693ffc955d6SMark F. Adams         }
6942d3561bbSSatish Balay         sz = nASMBlocksArr[level];
6952d3561bbSSatish Balay         is = ASMLocalIDsArr[level];
696ffc955d6SMark F. Adams         ierr = PCSetType( subpc, PCGASM ); CHKERRQ(ierr);
697ffc955d6SMark F. Adams         if(sz==0){
698ffc955d6SMark F. Adams           IS is;
699ffc955d6SMark F. Adams           PetscInt my0,kk;
700ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange( Aarr[level], &my0, &kk ); CHKERRQ(ierr);
701ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is ); CHKERRQ(ierr);
702a94c3b12SMark F. Adams           ierr = PCGASMSetLocalSubdomains( subpc, 1, &is, PETSC_NULL ); CHKERRQ(ierr);
703a94c3b12SMark F. Adams           ierr = ISDestroy( &is ); CHKERRQ(ierr);
704ffc955d6SMark F. Adams         }
705ffc955d6SMark F. Adams         else {
706a94c3b12SMark F. Adams           PetscInt kk;
707a94c3b12SMark F. Adams           ierr = PCGASMSetLocalSubdomains( subpc, sz, is, PETSC_NULL ); CHKERRQ(ierr);
708a94c3b12SMark F. Adams           for(kk=0;kk<sz;kk++){
709a94c3b12SMark F. Adams             ierr = ISDestroy( &is[kk] ); CHKERRQ(ierr);
710a94c3b12SMark F. Adams           }
711ffc955d6SMark F. Adams           ierr = PetscFree( is ); CHKERRQ(ierr);
712ffc955d6SMark F. Adams         }
713ffc955d6SMark F. Adams         ierr = PCGASMSetOverlap( subpc, 0 ); CHKERRQ(ierr);
714ffc955d6SMark F. Adams 
715ffc955d6SMark F. Adams         ASMLocalIDsArr[level] = PETSC_NULL;
716ffc955d6SMark F. Adams         nASMBlocksArr[level] = 0;
717ffc955d6SMark F. Adams         ierr = PCGASMSetType( subpc, PC_GASM_BASIC ); CHKERRQ(ierr);
718ffc955d6SMark F. Adams       }
719ffc955d6SMark F. Adams       else {
7209d5b6da9SMark F. Adams         ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr);
721ffc955d6SMark F. Adams       }
722f6536408SMark F. Adams       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
723d5280255SMark F. Adams       /* set ops */
724f6536408SMark F. Adams       ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
725d5280255SMark F. Adams       ierr = PCMGSetInterpolation( pc, lidx, Parr[level+1] );CHKERRQ(ierr);
726d5280255SMark F. Adams     }
727d5280255SMark F. Adams     {
728d5280255SMark F. Adams       /* coarse grid */
729d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
730d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
731d5280255SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
732d5280255SMark F. Adams       ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
733d5280255SMark F. Adams       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
734d5280255SMark F. Adams       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
735d5280255SMark F. Adams       ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
736d5280255SMark F. Adams       ierr = PCSetUp( subpc ); CHKERRQ(ierr);
737d5280255SMark F. Adams       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
738d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
739d5280255SMark F. Adams       ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
740d5280255SMark F. Adams     }
741d5280255SMark F. Adams 
742d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
743d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin( (PetscObject)pc );CHKERRQ(ierr);
744d5280255SMark F. Adams     ierr = PCSetFromOptions_MG( pc ); CHKERRQ(ierr);
745d5280255SMark F. Adams     ierr = PetscOptionsEnd();  CHKERRQ(ierr);
746d5280255SMark F. Adams     if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
747d5280255SMark F. Adams 
748d5280255SMark F. Adams     /* create cheby smoothers */
749d5280255SMark F. Adams     for ( lidx = 1, level = pc_gamg->Nlevels-2;
750d5280255SMark F. Adams           lidx <= fine_level;
751d5280255SMark F. Adams           lidx++, level--) {
752d5280255SMark F. Adams       KSP smoother;
753ffc955d6SMark F. Adams       PetscBool flag;
754d5280255SMark F. Adams       PC subpc;
755d5280255SMark F. Adams 
756ffc955d6SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
757ffc955d6SMark F. Adams 
758ffc955d6SMark F. Adams       /* do my own cheby */
759ffc955d6SMark F. Adams       ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &flag ); CHKERRQ(ierr);
760ffc955d6SMark F. Adams       if( flag ) {
761ffc955d6SMark F. Adams         PetscReal emax, emin;
762d5280255SMark F. Adams         ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
763ffc955d6SMark F. Adams         ierr = PetscTypeCompare( (PetscObject)subpc, PCJACOBI, &flag ); CHKERRQ(ierr);
764ffc955d6SMark F. Adams         if( flag && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
765587fa25dSMark F. Adams         else{ /* eigen estimate 'emax' */
766587fa25dSMark F. Adams           KSP eksp; Mat Lmat = Aarr[level];
767b2a4f308SMark F. Adams           Vec bb, xx;
768038e3b61SMark F. Adams 
7695745e0f5SMark F. Adams           ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
7705745e0f5SMark F. Adams           ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
771fc4362bfSMark F. Adams           {
772fc4362bfSMark F. Adams             PetscRandom    rctx;
773fc4362bfSMark F. Adams             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
774fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
775fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
776fc4362bfSMark F. Adams             ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
7775b89ad90SMark F. Adams           }
778a94c3b12SMark F. Adams           if( removedEqs[level] ) {
779*302f38e8SMark F. Adams             /* being very careful - zeroing out BC rows (this is not done in agg.c estimates) */
780a94c3b12SMark F. Adams             PetscScalar *zeros;
781bea1fdf0SMark F. Adams             PetscInt ii,jj, *idx_bs, sz, bs=level_bs[level];
782bea1fdf0SMark F. Adams             const PetscInt *idx;
783a94c3b12SMark F. Adams             ierr = ISGetLocalSize( removedEqs[level], &sz ); CHKERRQ(ierr);
784a94c3b12SMark F. Adams             ierr = PetscMalloc( bs*sz*sizeof(PetscScalar), &zeros ); CHKERRQ(ierr);
785a94c3b12SMark F. Adams             for(ii=0;ii<bs*sz;ii++) zeros[ii] = 0.;
786a94c3b12SMark F. Adams             ierr = PetscMalloc( bs*sz*sizeof(PetscInt), &idx_bs ); CHKERRQ(ierr);
787a94c3b12SMark F. Adams             ierr = ISGetIndices( removedEqs[level], &idx); CHKERRQ(ierr);
788a94c3b12SMark F. Adams             for(ii=0;ii<sz;ii++) {
789a94c3b12SMark F. Adams               for(jj=0;jj<bs;jj++) {
790a94c3b12SMark F. Adams                 idx_bs[ii] = bs*idx[ii]+jj;
791a94c3b12SMark F. Adams               }
792a94c3b12SMark F. Adams             }
793a94c3b12SMark F. Adams             ierr = ISRestoreIndices( removedEqs[level], &idx ); CHKERRQ(ierr);
794a94c3b12SMark F. Adams             if( sz > 0 ) {
795a94c3b12SMark F. Adams               ierr = VecSetValues( bb, sz, idx_bs, zeros, INSERT_VALUES );  CHKERRQ(ierr);
796a94c3b12SMark F. Adams             }
797a94c3b12SMark F. Adams             ierr = PetscFree( idx_bs );  CHKERRQ(ierr);
798a94c3b12SMark F. Adams             ierr = PetscFree( zeros );  CHKERRQ(ierr);
799a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb); CHKERRQ(ierr);
800a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb); CHKERRQ(ierr);
801a94c3b12SMark F. Adams             ierr = ISDestroy( &removedEqs[level] );                    CHKERRQ(ierr);
802a94c3b12SMark F. Adams           }
803fc4362bfSMark F. Adams           ierr = KSPCreate( wcomm, &eksp );CHKERRQ(ierr);
804038e3b61SMark F. Adams           ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
805fc4362bfSMark F. Adams           CHKERRQ(ierr);
806fc4362bfSMark F. Adams           ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
8071a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
8081a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix( eksp, "gamg_est_");         CHKERRQ(ierr);
809f6536408SMark F. Adams           ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
810f6536408SMark F. Adams 
811f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
812f6536408SMark F. Adams           ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
813fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
8145a9b9e01SMark F. Adams 
815d3d0db20SJed Brown           /* set PC type to be same as smoother */
816ffc955d6SMark F. Adams           ierr = KSPSetPC( eksp, subpc ); CHKERRQ( ierr );
817b2a4f308SMark F. Adams 
8185a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
8195a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
8205a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
821fc4362bfSMark F. Adams           ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
8225a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
8235a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
8245a9b9e01SMark F. Adams 
825fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
8265a9b9e01SMark F. Adams 
827fc4362bfSMark F. Adams           ierr = VecDestroy( &xx );       CHKERRQ(ierr);
828fc4362bfSMark F. Adams           ierr = VecDestroy( &bb );       CHKERRQ(ierr);
829fc4362bfSMark F. Adams           ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
830f6536408SMark F. Adams 
831ffc955d6SMark F. Adams           if( pc_gamg->verbose > 0 ) {
832a94c3b12SMark F. Adams             PetscInt N1, tt;
833a94c3b12SMark F. Adams             ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
834a94c3b12SMark F. Adams             PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1);
835f6536408SMark F. Adams           }
836fc4362bfSMark F. Adams         }
837038e3b61SMark F. Adams         {
838038e3b61SMark F. Adams           PetscInt N1, N0, tt;
839038e3b61SMark F. Adams           ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
840038e3b61SMark F. Adams           ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
841f6536408SMark F. Adams           /* heuristic - is this crap? */
842f6536408SMark F. Adams           emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
843038e3b61SMark F. Adams           emax *= 1.05;
844038e3b61SMark F. Adams         }
845fc4362bfSMark F. Adams         ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
846ffc955d6SMark F. Adams       } /* setup checby flag */
847ffc955d6SMark F. Adams     } /* non-coarse levels */
848737a81a9SMark F. Adams 
849d5280255SMark F. Adams     /* clean up */
850d5280255SMark F. Adams     for(level=1;level<pc_gamg->Nlevels;level++){
851587fa25dSMark F. Adams       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
852587fa25dSMark F. Adams       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
8535b89ad90SMark F. Adams     }
8540cbbd2e1SMark F. Adams 
8550cbbd2e1SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
8560cbbd2e1SMark F. Adams 
857a94c3b12SMark F. Adams     if( PETSC_FALSE ){
85858471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
8599d5b6da9SMark F. Adams       ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
86058471d46SMark F. Adams       ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
86158471d46SMark F. Adams     }
8625f8cf99dSMark F. Adams   }
8635f8cf99dSMark F. Adams   else {
8645f8cf99dSMark F. Adams     KSP smoother;
865b43b03e9SMark F. Adams     if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__);
8669d5b6da9SMark F. Adams     ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
8675f8cf99dSMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
8685f8cf99dSMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
8699d5b6da9SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
8705f8cf99dSMark F. Adams   }
87184d3f75bSMark F. Adams 
8725b89ad90SMark F. Adams   PetscFunctionReturn(0);
8735b89ad90SMark F. Adams }
8745b89ad90SMark F. Adams 
875eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8765b89ad90SMark F. Adams /*
8775b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8785b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8795b89ad90SMark F. Adams 
8805b89ad90SMark F. Adams    Input Parameter:
8815b89ad90SMark F. Adams .  pc - the preconditioner context
8825b89ad90SMark F. Adams 
8835b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8845b89ad90SMark F. Adams */
8855b89ad90SMark F. Adams #undef __FUNCT__
8865b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
8875b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG( PC pc )
8885b89ad90SMark F. Adams {
8895b89ad90SMark F. Adams   PetscErrorCode  ierr;
8905b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
8915b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
8925b89ad90SMark F. Adams 
8935b89ad90SMark F. Adams   PetscFunctionBegin;
8945b89ad90SMark F. Adams   ierr = PCReset_GAMG( pc );CHKERRQ(ierr);
8955b89ad90SMark F. Adams   ierr = PetscFree( pc_gamg );CHKERRQ(ierr);
8965b89ad90SMark F. Adams   ierr = PCDestroy_MG( pc );CHKERRQ(ierr);
8975b89ad90SMark F. Adams   PetscFunctionReturn(0);
8985b89ad90SMark F. Adams }
8995b89ad90SMark F. Adams 
900676e1743SMark F. Adams 
901676e1743SMark F. Adams #undef __FUNCT__
902676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
903676e1743SMark F. Adams /*@
904676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
905676e1743SMark F. Adams    processor reduction.
906676e1743SMark F. Adams 
907676e1743SMark F. Adams    Not Collective on PC
908676e1743SMark F. Adams 
909676e1743SMark F. Adams    Input Parameters:
910676e1743SMark F. Adams .  pc - the preconditioner context
911676e1743SMark F. Adams 
912676e1743SMark F. Adams 
913676e1743SMark F. Adams    Options Database Key:
914676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
915676e1743SMark F. Adams 
916676e1743SMark F. Adams    Level: intermediate
917676e1743SMark F. Adams 
918676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
919676e1743SMark F. Adams 
920676e1743SMark F. Adams .seealso: ()
921676e1743SMark F. Adams @*/
922676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
923676e1743SMark F. Adams {
924676e1743SMark F. Adams   PetscErrorCode ierr;
925676e1743SMark F. Adams 
926676e1743SMark F. Adams   PetscFunctionBegin;
927676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
928676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
929676e1743SMark F. Adams   PetscFunctionReturn(0);
930676e1743SMark F. Adams }
931676e1743SMark F. Adams 
932676e1743SMark F. Adams EXTERN_C_BEGIN
933676e1743SMark F. Adams #undef __FUNCT__
934676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
935676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
936676e1743SMark F. Adams {
937c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
938c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
939676e1743SMark F. Adams 
940676e1743SMark F. Adams   PetscFunctionBegin;
9419d5b6da9SMark F. Adams   if(n>0) pc_gamg->min_eq_proc = n;
942676e1743SMark F. Adams   PetscFunctionReturn(0);
943676e1743SMark F. Adams }
944676e1743SMark F. Adams EXTERN_C_END
945676e1743SMark F. Adams 
946676e1743SMark F. Adams #undef __FUNCT__
947389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
948389730f3SMark F. Adams /*@
949389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
950389730f3SMark F. Adams 
951389730f3SMark F. Adams  Collective on PC
952389730f3SMark F. Adams 
953389730f3SMark F. Adams    Input Parameters:
954389730f3SMark F. Adams .  pc - the preconditioner context
955389730f3SMark F. Adams 
956389730f3SMark F. Adams 
957389730f3SMark F. Adams    Options Database Key:
958389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
959389730f3SMark F. Adams 
960389730f3SMark F. Adams    Level: intermediate
961389730f3SMark F. Adams 
962389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
963389730f3SMark F. Adams 
964389730f3SMark F. Adams .seealso: ()
965389730f3SMark F. Adams  @*/
966389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
967389730f3SMark F. Adams {
968389730f3SMark F. Adams   PetscErrorCode ierr;
969389730f3SMark F. Adams 
970389730f3SMark F. Adams   PetscFunctionBegin;
971389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
972389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
973389730f3SMark F. Adams   PetscFunctionReturn(0);
974389730f3SMark F. Adams }
975389730f3SMark F. Adams 
976389730f3SMark F. Adams EXTERN_C_BEGIN
977389730f3SMark F. Adams #undef __FUNCT__
978389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
979389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
980389730f3SMark F. Adams {
981389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
982389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
983389730f3SMark F. Adams 
984389730f3SMark F. Adams   PetscFunctionBegin;
9859d5b6da9SMark F. Adams   if(n>0) pc_gamg->coarse_eq_limit = n;
986389730f3SMark F. Adams   PetscFunctionReturn(0);
987389730f3SMark F. Adams }
988389730f3SMark F. Adams EXTERN_C_END
989389730f3SMark F. Adams 
990389730f3SMark F. Adams #undef __FUNCT__
9918263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
992676e1743SMark F. Adams /*@
9938263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
994676e1743SMark F. Adams 
995676e1743SMark F. Adams    Collective on PC
996676e1743SMark F. Adams 
997676e1743SMark F. Adams    Input Parameters:
998676e1743SMark F. Adams .  pc - the preconditioner context
999676e1743SMark F. Adams 
1000676e1743SMark F. Adams 
1001676e1743SMark F. Adams    Options Database Key:
10028263b398SMark F. Adams .  -pc_gamg_repartition
1003676e1743SMark F. Adams 
1004676e1743SMark F. Adams    Level: intermediate
1005676e1743SMark F. Adams 
1006676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1007676e1743SMark F. Adams 
1008676e1743SMark F. Adams .seealso: ()
1009676e1743SMark F. Adams @*/
10108263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1011676e1743SMark F. Adams {
1012676e1743SMark F. Adams   PetscErrorCode ierr;
1013676e1743SMark F. Adams 
1014676e1743SMark F. Adams   PetscFunctionBegin;
1015676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10168263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1017676e1743SMark F. Adams   PetscFunctionReturn(0);
1018676e1743SMark F. Adams }
1019676e1743SMark F. Adams 
1020676e1743SMark F. Adams EXTERN_C_BEGIN
1021676e1743SMark F. Adams #undef __FUNCT__
10228263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
10238263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1024676e1743SMark F. Adams {
1025c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1026c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1027676e1743SMark F. Adams 
1028676e1743SMark F. Adams   PetscFunctionBegin;
10299d5b6da9SMark F. Adams   pc_gamg->repart = n;
1030676e1743SMark F. Adams   PetscFunctionReturn(0);
1031676e1743SMark F. Adams }
1032676e1743SMark F. Adams EXTERN_C_END
1033676e1743SMark F. Adams 
1034676e1743SMark F. Adams #undef __FUNCT__
1035ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1036ffc955d6SMark F. Adams /*@
1037ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1038ffc955d6SMark F. Adams 
1039ffc955d6SMark F. Adams    Collective on PC
1040ffc955d6SMark F. Adams 
1041ffc955d6SMark F. Adams    Input Parameters:
1042ffc955d6SMark F. Adams .  pc - the preconditioner context
1043ffc955d6SMark F. Adams 
1044ffc955d6SMark F. Adams 
1045ffc955d6SMark F. Adams    Options Database Key:
1046ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1047ffc955d6SMark F. Adams 
1048ffc955d6SMark F. Adams    Level: intermediate
1049ffc955d6SMark F. Adams 
1050ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1051ffc955d6SMark F. Adams 
1052ffc955d6SMark F. Adams .seealso: ()
1053ffc955d6SMark F. Adams @*/
1054ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1055ffc955d6SMark F. Adams {
1056ffc955d6SMark F. Adams   PetscErrorCode ierr;
1057ffc955d6SMark F. Adams 
1058ffc955d6SMark F. Adams   PetscFunctionBegin;
1059ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1060ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1061ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1062ffc955d6SMark F. Adams }
1063ffc955d6SMark F. Adams 
1064ffc955d6SMark F. Adams EXTERN_C_BEGIN
1065ffc955d6SMark F. Adams #undef __FUNCT__
1066ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
1067ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1068ffc955d6SMark F. Adams {
1069ffc955d6SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1070ffc955d6SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1071ffc955d6SMark F. Adams 
1072ffc955d6SMark F. Adams   PetscFunctionBegin;
1073ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1074ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1075ffc955d6SMark F. Adams }
1076ffc955d6SMark F. Adams EXTERN_C_END
1077ffc955d6SMark F. Adams 
1078ffc955d6SMark F. Adams #undef __FUNCT__
10794ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
10804ef23d27SMark F. Adams /*@
10814ef23d27SMark F. Adams    PCGAMGSetNlevels -
10824ef23d27SMark F. Adams 
10834ef23d27SMark F. Adams    Not collective on PC
10844ef23d27SMark F. Adams 
10854ef23d27SMark F. Adams    Input Parameters:
10864ef23d27SMark F. Adams .  pc - the preconditioner context
10874ef23d27SMark F. Adams 
10884ef23d27SMark F. Adams 
10894ef23d27SMark F. Adams    Options Database Key:
10904ef23d27SMark F. Adams .  -pc_mg_levels
10914ef23d27SMark F. Adams 
10924ef23d27SMark F. Adams    Level: intermediate
10934ef23d27SMark F. Adams 
10944ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
10954ef23d27SMark F. Adams 
10964ef23d27SMark F. Adams .seealso: ()
10974ef23d27SMark F. Adams @*/
10984ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
10994ef23d27SMark F. Adams {
11004ef23d27SMark F. Adams   PetscErrorCode ierr;
11014ef23d27SMark F. Adams 
11024ef23d27SMark F. Adams   PetscFunctionBegin;
11034ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11044ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
11054ef23d27SMark F. Adams   PetscFunctionReturn(0);
11064ef23d27SMark F. Adams }
11074ef23d27SMark F. Adams 
11084ef23d27SMark F. Adams EXTERN_C_BEGIN
11094ef23d27SMark F. Adams #undef __FUNCT__
11104ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
11114ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
11124ef23d27SMark F. Adams {
11134ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
11144ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
11154ef23d27SMark F. Adams 
11164ef23d27SMark F. Adams   PetscFunctionBegin;
11179d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
11184ef23d27SMark F. Adams   PetscFunctionReturn(0);
11194ef23d27SMark F. Adams }
11204ef23d27SMark F. Adams EXTERN_C_END
11214ef23d27SMark F. Adams 
11224ef23d27SMark F. Adams #undef __FUNCT__
11233542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
11243542efc5SMark F. Adams /*@
11253542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
11263542efc5SMark F. Adams 
11273542efc5SMark F. Adams    Not collective on PC
11283542efc5SMark F. Adams 
11293542efc5SMark F. Adams    Input Parameters:
11303542efc5SMark F. Adams .  pc - the preconditioner context
11313542efc5SMark F. Adams 
11323542efc5SMark F. Adams 
11333542efc5SMark F. Adams    Options Database Key:
11343542efc5SMark F. Adams .  -pc_gamg_threshold
11353542efc5SMark F. Adams 
11363542efc5SMark F. Adams    Level: intermediate
11373542efc5SMark F. Adams 
11383542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
11393542efc5SMark F. Adams 
11403542efc5SMark F. Adams .seealso: ()
11413542efc5SMark F. Adams @*/
11423542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
11433542efc5SMark F. Adams {
11443542efc5SMark F. Adams   PetscErrorCode ierr;
11453542efc5SMark F. Adams 
11463542efc5SMark F. Adams   PetscFunctionBegin;
11473542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11483542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
11493542efc5SMark F. Adams   PetscFunctionReturn(0);
11503542efc5SMark F. Adams }
11513542efc5SMark F. Adams 
11523542efc5SMark F. Adams EXTERN_C_BEGIN
11533542efc5SMark F. Adams #undef __FUNCT__
11543542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
11553542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
11563542efc5SMark F. Adams {
1157c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1158c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
11593542efc5SMark F. Adams 
11603542efc5SMark F. Adams   PetscFunctionBegin;
11619d5b6da9SMark F. Adams   pc_gamg->threshold = n;
11623542efc5SMark F. Adams   PetscFunctionReturn(0);
11633542efc5SMark F. Adams }
11643542efc5SMark F. Adams EXTERN_C_END
11653542efc5SMark F. Adams 
11663542efc5SMark F. Adams #undef __FUNCT__
11679d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1168676e1743SMark F. Adams /*@
11699d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1170676e1743SMark F. Adams 
1171676e1743SMark F. Adams    Collective on PC
1172676e1743SMark F. Adams 
1173676e1743SMark F. Adams    Input Parameters:
1174676e1743SMark F. Adams .  pc - the preconditioner context
1175676e1743SMark F. Adams 
1176676e1743SMark F. Adams 
1177676e1743SMark F. Adams    Options Database Key:
11783542efc5SMark F. Adams .  -pc_gamg_type
1179676e1743SMark F. Adams 
1180676e1743SMark F. Adams    Level: intermediate
1181676e1743SMark F. Adams 
1182676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1183676e1743SMark F. Adams 
1184676e1743SMark F. Adams .seealso: ()
1185676e1743SMark F. Adams @*/
11860202ef74SSatish Balay PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type )
1187676e1743SMark F. Adams {
1188676e1743SMark F. Adams   PetscErrorCode ierr;
1189676e1743SMark F. Adams 
1190676e1743SMark F. Adams   PetscFunctionBegin;
1191676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11920202ef74SSatish Balay   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type));
1193676e1743SMark F. Adams   CHKERRQ(ierr);
1194676e1743SMark F. Adams   PetscFunctionReturn(0);
1195676e1743SMark F. Adams }
1196676e1743SMark F. Adams 
1197676e1743SMark F. Adams EXTERN_C_BEGIN
1198676e1743SMark F. Adams #undef __FUNCT__
11999d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
12000202ef74SSatish Balay PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type )
1201676e1743SMark F. Adams {
12029d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
1203676e1743SMark F. Adams 
1204676e1743SMark F. Adams   PetscFunctionBegin;
12059d5b6da9SMark F. Adams   ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r);
12069d5b6da9SMark F. Adams   CHKERRQ(ierr);
12079d5b6da9SMark F. Adams 
12089d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
12099d5b6da9SMark F. Adams 
12109d5b6da9SMark F. Adams   /* call sub create method */
12119d5b6da9SMark F. Adams   ierr = (*r)(pc); CHKERRQ(ierr);
12129d5b6da9SMark F. Adams 
1213676e1743SMark F. Adams   PetscFunctionReturn(0);
1214676e1743SMark F. Adams }
1215676e1743SMark F. Adams EXTERN_C_END
1216676e1743SMark F. Adams 
12175b89ad90SMark F. Adams #undef __FUNCT__
12185b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
12195b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG( PC pc )
12205b89ad90SMark F. Adams {
1221676e1743SMark F. Adams   PetscErrorCode  ierr;
1222676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1223676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1224676e1743SMark F. Adams   PetscBool        flag;
1225b43b03e9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
12265b89ad90SMark F. Adams 
12275b89ad90SMark F. Adams   PetscFunctionBegin;
1228676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1229676e1743SMark F. Adams   {
123075b74bdaSMark F. Adams     /* -pc_gamg_verbose */
12319d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
12329d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
12339d5b6da9SMark F. Adams                            &pc_gamg->verbose, PETSC_NULL );
12349d5b6da9SMark F. Adams     CHKERRQ(ierr);
123575b74bdaSMark F. Adams 
12368263b398SMark F. Adams     /* -pc_gamg_repartition */
12378263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
12388263b398SMark F. Adams                             "Repartion coarse grids (false)",
12398263b398SMark F. Adams                             "PCGAMGRepartitioning",
12409d5b6da9SMark F. Adams                             pc_gamg->repart,
12419d5b6da9SMark F. Adams                             &pc_gamg->repart,
1242676e1743SMark F. Adams                             &flag);
1243676e1743SMark F. Adams     CHKERRQ(ierr);
1244676e1743SMark F. Adams 
1245ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
1246ffc955d6SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1247ffc955d6SMark F. Adams                             "Use aggregation agragates for GASM smoother (false)",
1248ffc955d6SMark F. Adams                             "PCGAMGUseASMAggs",
1249ffc955d6SMark F. Adams                             pc_gamg->use_aggs_in_gasm,
1250ffc955d6SMark F. Adams                             &pc_gamg->use_aggs_in_gasm,
1251ffc955d6SMark F. Adams                             &flag);
1252ffc955d6SMark F. Adams     CHKERRQ(ierr);
1253ffc955d6SMark F. Adams 
1254c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1255676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1256676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1257676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
12589d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
12599d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1260676e1743SMark F. Adams                            &flag );
1261676e1743SMark F. Adams     CHKERRQ(ierr);
12623542efc5SMark F. Adams 
1263389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1264389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1265389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1266389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
12679d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
12689d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1269389730f3SMark F. Adams                            &flag );
1270389730f3SMark F. Adams     CHKERRQ(ierr);
1271389730f3SMark F. Adams 
1272c20e4228SMark F. Adams     /* -pc_gamg_threshold */
12733542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
12743542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
12753542efc5SMark F. Adams                             "PCGAMGSetThreshold",
12769d5b6da9SMark F. Adams                             pc_gamg->threshold,
12779d5b6da9SMark F. Adams                             &pc_gamg->threshold,
12783542efc5SMark F. Adams                             &flag );
12793542efc5SMark F. Adams     CHKERRQ(ierr);
1280b43b03e9SMark F. Adams     if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);
12814ef23d27SMark F. Adams 
12824ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
12834ef23d27SMark F. Adams                            "Set number of MG levels",
12844ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
12859d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
12869d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
12874ef23d27SMark F. Adams                            &flag );
1288676e1743SMark F. Adams   }
1289676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1290676e1743SMark F. Adams 
12915b89ad90SMark F. Adams   PetscFunctionReturn(0);
12925b89ad90SMark F. Adams }
12935b89ad90SMark F. Adams 
12945b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
12955b89ad90SMark F. Adams /*MC
1296856830bfSJed Brown      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
12979d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
12985b89ad90SMark F. Adams 
1299280d9858SJed Brown    Options Database Keys:
13005b89ad90SMark F. Adams    Multigrid options(inherited)
1301280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1302280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1303280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1304280d9858SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
13055b89ad90SMark F. Adams 
13065b89ad90SMark F. Adams   Level: intermediate
1307280d9858SJed Brown 
13085b89ad90SMark F. Adams   Concepts: multigrid
13095b89ad90SMark F. Adams 
13105b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1311280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
13125b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
13135b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
13145b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
13155b89ad90SMark F. Adams M*/
13165b89ad90SMark F. Adams EXTERN_C_BEGIN
13175b89ad90SMark F. Adams #undef __FUNCT__
13185b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
13195b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG( PC pc )
13205b89ad90SMark F. Adams {
13215b89ad90SMark F. Adams   PetscErrorCode  ierr;
13225b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
13235b89ad90SMark F. Adams   PC_MG           *mg;
13240cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13259d5b6da9SMark F. Adams   static long count = 0;
13265ee9c036SSatish Balay #endif
13275b89ad90SMark F. Adams 
13285b89ad90SMark F. Adams   PetscFunctionBegin;
132960a6d8f8SMark F. Adams 
13305b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
13315b89ad90SMark F. Adams   ierr = PCSetType( pc, PCMG );  CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
13325b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr);
13335b89ad90SMark F. Adams 
13345b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
13355b89ad90SMark F. Adams   ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr);
13365b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
1337ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
13385b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
13395b89ad90SMark F. Adams 
13409d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
13419d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
13429d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
13439d5b6da9SMark F. Adams   pc_gamg->data = 0;
13445b89ad90SMark F. Adams 
13459d5b6da9SMark F. Adams   /* register AMG type */
13469d5b6da9SMark F. Adams   if( !GAMGList ){
13479d5b6da9SMark F. Adams     ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
13489d5b6da9SMark F. Adams     ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
13499d5b6da9SMark F. Adams   }
13509d5b6da9SMark F. Adams 
13519d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
13525b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
13535b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
13545b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
13555b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
13565b89ad90SMark F. Adams 
13575b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1358676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1359676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1360676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1361676e1743SMark F. Adams   CHKERRQ(ierr);
1362676e1743SMark F. Adams 
1363676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1364389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_C",
1365389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_GAMG",
1366389730f3SMark F. Adams 					    PCGAMGSetCoarseEqLim_GAMG);
1367389730f3SMark F. Adams   CHKERRQ(ierr);
1368389730f3SMark F. Adams 
1369389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
13708263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_C",
13718263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_GAMG",
13728263b398SMark F. Adams 					    PCGAMGSetRepartitioning_GAMG);
1373676e1743SMark F. Adams   CHKERRQ(ierr);
1374676e1743SMark F. Adams 
1375676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1376ffc955d6SMark F. Adams 					    "PCGAMGSetUseASMAggs_C",
1377ffc955d6SMark F. Adams 					    "PCGAMGSetUseASMAggs_GAMG",
1378ffc955d6SMark F. Adams 					    PCGAMGSetUseASMAggs_GAMG);
1379ffc955d6SMark F. Adams   CHKERRQ(ierr);
1380ffc955d6SMark F. Adams 
1381ffc955d6SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
13823542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
13833542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
13843542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
13853542efc5SMark F. Adams   CHKERRQ(ierr);
13863542efc5SMark F. Adams 
13873542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
13889d5b6da9SMark F. Adams 					    "PCGAMGSetType_C",
13899d5b6da9SMark F. Adams 					    "PCGAMGSetType_GAMG",
13909d5b6da9SMark F. Adams 					    PCGAMGSetType_GAMG);
1391676e1743SMark F. Adams   CHKERRQ(ierr);
1392c97e1df0SMark F. Adams 
13939d5b6da9SMark F. Adams   pc_gamg->repart = PETSC_FALSE;
1394ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
13959d5b6da9SMark F. Adams   pc_gamg->min_eq_proc = 100;
1396c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit = 800;
13979d5b6da9SMark F. Adams   pc_gamg->threshold = 0.05;
13989d5b6da9SMark F. Adams   pc_gamg->Nlevels = GAMG_MAXLEVELS;
13999d5b6da9SMark F. Adams   pc_gamg->verbose = 0;
14009d5b6da9SMark F. Adams   pc_gamg->emax_id = -1;
14019d5b6da9SMark F. Adams   pc_gamg->col_bs_id = -1;
14029d5b6da9SMark F. Adams 
14030cbbd2e1SMark F. Adams   /* private events */
14040cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1405785cba28SMark F. Adams   if( count++ == 0 ) {
14060cbbd2e1SMark F. Adams     PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
14070cbbd2e1SMark F. Adams     PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
14080cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
14090cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
14100cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
14110cbbd2e1SMark F. Adams     PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
14120cbbd2e1SMark F. Adams     PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
14130cbbd2e1SMark F. Adams     PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
14140cbbd2e1SMark F. Adams     PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
141560a6d8f8SMark F. Adams     PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
14160cbbd2e1SMark F. Adams     PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
14170cbbd2e1SMark F. Adams     PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
14180cbbd2e1SMark F. Adams     PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
14190cbbd2e1SMark F. Adams     PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
14200cbbd2e1SMark F. Adams     PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
14210cbbd2e1SMark F. Adams     PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
14220cbbd2e1SMark F. Adams     PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);
1423f852f58cSMark F. Adams 
14240cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
14250cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
14260cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1427b4fbaa2aSMark F. Adams     /* create timer stages */
1428b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1429b4fbaa2aSMark F. Adams     {
1430b4fbaa2aSMark F. Adams       char str[32];
1431b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1432b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1433b4fbaa2aSMark F. Adams       PetscInt lidx;
1434b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1435b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1436b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1437b4fbaa2aSMark F. Adams       }
1438b4fbaa2aSMark F. Adams     }
1439b4fbaa2aSMark F. Adams #endif
1440b4fbaa2aSMark F. Adams   }
1441b4fbaa2aSMark F. Adams #endif
14420cbbd2e1SMark F. Adams   /* general events */
14430cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
144460a6d8f8SMark F. Adams   PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);
14450cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);
14460cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
14470cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
14480cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
14490cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
145060a6d8f8SMark F. Adams   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);
14510cbbd2e1SMark F. Adams #endif
14520cbbd2e1SMark F. Adams 
145360a6d8f8SMark F. Adams   /* instantiate derived type */
145460a6d8f8SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
145560a6d8f8SMark F. Adams   {
145660a6d8f8SMark F. Adams     char tname[256] = GAMGAGG;
145760a6d8f8SMark F. Adams     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",
145860a6d8f8SMark F. Adams                             GAMGList, tname, tname, sizeof(tname), PETSC_NULL );
145960a6d8f8SMark F. Adams     CHKERRQ(ierr);
146060a6d8f8SMark F. Adams     ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr);
146160a6d8f8SMark F. Adams   }
146260a6d8f8SMark F. Adams   ierr = PetscOptionsTail();   CHKERRQ(ierr);
146360a6d8f8SMark F. Adams 
14645b89ad90SMark F. Adams   PetscFunctionReturn(0);
14655b89ad90SMark F. Adams }
14665b89ad90SMark F. Adams EXTERN_C_END
1467