xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 0cbbd2e198edea5eef85858c2ecfc70cb9545f65)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4313a3e24SSatish Balay #include "private/matimpl.h"
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
65a9b9e01SMark F. Adams #include <private/kspimpl.h>
7f96513f1SMatthew G Knepley 
8*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
9*0cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
10b4fbaa2aSMark F. Adams #endif
11*0cbbd2e1SMark F. Adams 
12*0cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
13*0cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_AGG;
14*0cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_GEO;
15*0cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
16*0cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
17*0cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
18*0cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
19*0cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGOptprol_AGG;
20*0cbbd2e1SMark F. Adams #endif
21*0cbbd2e1SMark F. Adams 
22b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
23b4fbaa2aSMark F. Adams 
24b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
25*0cbbd2e1SMark 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"
67*0cbbd2e1SMark 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);
127*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
128*0cbbd2e1SMark 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];
183b4fbaa2aSMark F. Adams 	  sprintf(fname,"part_mat_%d.mat",llev);
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;
284*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
285*0cbbd2e1SMark 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);
349*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
350*0cbbd2e1SMark 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);
358*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
359*0cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
360*0cbbd2e1SMark 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);
369*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
370*0cbbd2e1SMark 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;
376*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
377*0cbbd2e1SMark 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);
387*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
388*0cbbd2e1SMark 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;
440d3d6bff4SMark F. Adams   PetscInt         fine_level,level,level1,M,N,bs,nloc,lidx,Istart,Iend;
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];
4450534fa05SHong Zhang   PetscLogDouble   nnz0=0,nnztot=0;
446737a81a9SMark F. Adams   MatInfo          info;
4475ef31b24SMark F. Adams 
4485b89ad90SMark F. Adams   PetscFunctionBegin;
44984d3f75bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
45084d3f75bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
45184d3f75bSMark F. Adams   if( pc_gamg->setup_count++ > 0 ) {
45284d3f75bSMark F. Adams     PC_MG_Levels   **mglevels = mg->levels;
45303a628feSMark F. Adams     /* just do Galerkin grids */
45458471d46SMark F. Adams     Mat B,dA,dB;
455d5280255SMark F. Adams     assert(pc->setupcalled);
45658471d46SMark F. Adams 
4579d5b6da9SMark F. Adams     if( pc_gamg->Nlevels > 1 ) {
45858471d46SMark F. Adams       /* currently only handle case where mat and pmat are the same on coarser levels */
4599d5b6da9SMark F. Adams       ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
46058471d46SMark F. Adams       /* (re)set to get dirty flag */
4619d5b6da9SMark F. Adams       ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
46258471d46SMark F. Adams 
4639d5b6da9SMark F. Adams       for (level=pc_gamg->Nlevels-2; level>-1; level--) {
46403a628feSMark F. Adams         /* the first time through the matrix structure has changed from repartitioning */
46584d3f75bSMark F. Adams         if( pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
46603a628feSMark F. Adams           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
467084a8fe3SJed Brown           ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
46803a628feSMark F. Adams           mglevels[level]->A = B;
46903a628feSMark F. Adams         }
47003a628feSMark F. Adams         else {
471084a8fe3SJed Brown           ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
47258471d46SMark F. Adams           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
47303a628feSMark F. Adams         }
47458471d46SMark F. Adams         ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
47558471d46SMark F. Adams         dB = B;
47658471d46SMark F. Adams         }
4775f8cf99dSMark F. Adams     }
478d5280255SMark F. Adams 
479d5280255SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
480d5280255SMark F. Adams 
481d5280255SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
482d5280255SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
483d5280255SMark F. Adams 
48458471d46SMark F. Adams     PetscFunctionReturn(0);
485eb07cef2SMark F. Adams   }
48684d3f75bSMark F. Adams   assert(pc->setupcalled == 0);
487f6536408SMark F. Adams 
4885b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
4892adcac29SMark F. Adams   ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr);
4902adcac29SMark F. Adams   ierr = MatGetSize( Pmat, &M, &N );CHKERRQ(ierr);
4912adcac29SMark F. Adams   ierr = MatGetOwnershipRange( Pmat, &Istart, &Iend ); CHKERRQ(ierr);
492eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
493eb07cef2SMark F. Adams 
4949d5b6da9SMark F. Adams   if( pc_gamg->data == 0 && nloc > 0 ) {
4959d5b6da9SMark F. Adams     if(!pc_gamg->createdefaultdata){
496c8b0795cSMark F. Adams       SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"'createdefaultdata' not set?!?! need to support NULL data!!!");
497eb07cef2SMark F. Adams     }
4989d5b6da9SMark F. Adams     ierr = pc_gamg->createdefaultdata( pc ); CHKERRQ(ierr);
4999d5b6da9SMark F. Adams   }
500038e3b61SMark F. Adams 
501eb07cef2SMark F. Adams   /* Get A_i and R_i */
502c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
5032adcac29SMark F. Adams     if(pc_gamg->verbose==1) ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);
5042adcac29SMark F. Adams     else ierr =  MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);
505b2a4f308SMark F. Adams     CHKERRQ(ierr);
506b2a4f308SMark F. Adams     nnz0 = info.nz_used;
507b2a4f308SMark F. Adams     nnztot = info.nz_used;
508b43b03e9SMark 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",
509c8b0795cSMark F. Adams                 mype,__FUNCT__,0,N,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
510b2a4f308SMark F. Adams                 (int)(nnz0/(PetscReal)N),npe);
511c8b0795cSMark F. Adams   }
51284d3f75bSMark F. Adams 
5138f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
5149d5b6da9SMark F. Adams         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
5150205a208SMark F. Adams         level++ ){
5165b89ad90SMark F. Adams     level1 = level + 1;
517*0cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
518b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
519b4fbaa2aSMark F. Adams #endif
520*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
521*0cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr);
522b4fbaa2aSMark F. Adams #endif
523c8b0795cSMark F. Adams     { /* construct prolongator */
524725b86d8SJed Brown       Mat Gmat;
525*0cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
526c8b0795cSMark F. Adams 
527b43b03e9SMark F. Adams       ierr = pc_gamg->graph( pc, Aarr[level], &Gmat ); CHKERRQ(ierr);
528*0cbbd2e1SMark F. Adams       ierr = pc_gamg->coarsen( pc, &Gmat, &agg_lists ); CHKERRQ(ierr);
529*0cbbd2e1SMark F. Adams       ierr = pc_gamg->prolongator( pc, Aarr[level], Gmat, agg_lists, &Parr[level1] );
5305b89ad90SMark F. Adams       CHKERRQ(ierr);
531c8b0795cSMark F. Adams 
532c8b0795cSMark F. Adams       if( Parr[level1] ){
5339d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
534c8b0795cSMark F. Adams         if( pc_gamg->col_bs_id != -1 ){
5359d5b6da9SMark F. Adams           PetscBool flag;
5369d5b6da9SMark F. Adams           ierr = PetscObjectComposedDataGetInt( (PetscObject)Parr[level1], pc_gamg->col_bs_id, bs, flag );
537c8b0795cSMark F. Adams           CHKERRQ( ierr ); assert(flag);
5389d5b6da9SMark F. Adams         }
539c8b0795cSMark F. Adams       }
540c8b0795cSMark F. Adams 
541c8b0795cSMark F. Adams       if( pc_gamg->optprol && Parr[level1] ){
542c8b0795cSMark F. Adams         /* smooth */
543c8b0795cSMark F. Adams         ierr = pc_gamg->optprol( pc, Aarr[level], &Parr[level1] ); CHKERRQ(ierr);
544c8b0795cSMark F. Adams       }
545c8b0795cSMark F. Adams 
546c8b0795cSMark F. Adams       ierr = MatDestroy( &Gmat );      CHKERRQ(ierr);
547*0cbbd2e1SMark F. Adams       ierr = AILDestroy( agg_lists );  CHKERRQ(ierr);
548c8b0795cSMark F. Adams     }
549*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
550*0cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
551c8b0795cSMark F. Adams #endif
5529d5b6da9SMark F. Adams     /* cache eigen estimate */
5539d5b6da9SMark F. Adams     if( pc_gamg->emax_id != -1 ){
5549d5b6da9SMark F. Adams       PetscBool flag;
555f73f8d2cSSatish Balay       ierr = PetscObjectComposedDataGetReal( (PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag );
5569d5b6da9SMark F. Adams       CHKERRQ( ierr );
5579d5b6da9SMark F. Adams       if( !flag ) emaxs[level] = -1.;
5589d5b6da9SMark F. Adams     }
5599d5b6da9SMark F. Adams     else emaxs[level] = -1.;
5602adcac29SMark F. Adams     if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
561c8b0795cSMark F. Adams     if( !Parr[level1] ) {
562b43b03e9SMark F. Adams       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);
563c8b0795cSMark F. Adams       break;
564c8b0795cSMark F. Adams     }
565*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
566*0cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
567b4fbaa2aSMark F. Adams #endif
568c8b0795cSMark F. Adams     ierr = createLevel( pc, Aarr[level], bs,
5699d5b6da9SMark F. Adams                         (PetscBool)(level==pc_gamg->Nlevels-2),
570c8b0795cSMark F. Adams                         &Parr[level1], &nactivepe, &Aarr[level1] );
5713530afc2SMark F. Adams     CHKERRQ(ierr);
572*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
573*0cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
574b4fbaa2aSMark F. Adams #endif
5753530afc2SMark F. Adams     ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
576c8b0795cSMark F. Adams     if (pc_gamg->verbose){
577*0cbbd2e1SMark F. Adams       PetscInt NN = M;
578*0cbbd2e1SMark F. Adams       if(pc_gamg->verbose==1) {
579*0cbbd2e1SMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);
580*0cbbd2e1SMark F. Adams         ierr = MatGetLocalSize( Aarr[level1], &NN, &N );CHKERRQ(ierr);
581*0cbbd2e1SMark F. Adams       }
582*0cbbd2e1SMark F. Adams       else ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info);
583*0cbbd2e1SMark F. Adams       CHKERRQ(ierr);
584*0cbbd2e1SMark F. Adams       nnztot += info.nz_used;
585b43b03e9SMark F. Adams       PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
586c8b0795cSMark F. Adams                   mype,__FUNCT__,(int)level1,N,pc_gamg->data_cell_cols,
587*0cbbd2e1SMark F. Adams                   (int)(info.nz_used/(PetscReal)NN), nactivepe );
588*0cbbd2e1SMark F. Adams       CHKERRQ(ierr);
589c8b0795cSMark F. Adams     }
59081708dccSMark F. Adams     /* stop if one node */
591c8b0795cSMark F. Adams     if( M/pc_gamg->data_cell_cols < 2 ) {
59281708dccSMark F. Adams       level++;
59381708dccSMark F. Adams       break;
59481708dccSMark F. Adams     }
595c8b0795cSMark F. Adams     /* Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; */
596c8b0795cSMark F. Adams     /* v = 1.e-10; /\* LU factor has hard wired numbers for small diags so this needs to match (yuk) *\/ */
597c8b0795cSMark F. Adams     /* ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); */
598c8b0795cSMark F. Adams     /* nloceq = Iend-Istart; */
599c8b0795cSMark F. Adams     /* ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr); */
600c8b0795cSMark F. Adams     /* ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr); */
601c8b0795cSMark F. Adams     /* ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr); */
602c8b0795cSMark F. Adams     /* for(kk=0;kk<nloceq;kk++){ */
603c8b0795cSMark F. Adams     /*   if(data_arr[kk]==0.0) { */
604c8b0795cSMark F. Adams     /*     id = kk + Istart; */
605c8b0795cSMark F. Adams     /*     ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); */
606c8b0795cSMark F. Adams     /*     CHKERRQ(ierr); */
607c8b0795cSMark F. Adams     /*     PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1); */
608c8b0795cSMark F. Adams     /*   } */
609c8b0795cSMark F. Adams     /* } */
610c8b0795cSMark F. Adams     /* ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); */
611c8b0795cSMark F. Adams     /* ierr = VecDestroy( &diag );                CHKERRQ(ierr); */
61284d3f75bSMark F. Adams     /* ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
61384d3f75bSMark F. Adams     /* ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
614b4fbaa2aSMark F. Adams 
615*0cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
616b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
617b4fbaa2aSMark F. Adams #endif
618c8b0795cSMark F. Adams   } /* levels */
619c8b0795cSMark F. Adams 
620c8b0795cSMark F. Adams   if( pc_gamg->data ) {
621c8b0795cSMark F. Adams     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
622c8b0795cSMark F. Adams     pc_gamg->data = 0;
6235b89ad90SMark F. Adams   }
624c8b0795cSMark F. Adams 
625*0cbbd2e1SMark F. Adams   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
6269d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6275b89ad90SMark F. Adams   fine_level = level;
6289d5b6da9SMark F. Adams   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
6295b89ad90SMark F. Adams 
63084d3f75bSMark F. Adams   /* simple setup */
63184d3f75bSMark F. Adams   if( !PETSC_TRUE ){
63284d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
63384d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
63484d3f75bSMark F. Adams          lidx<fine_level;
63584d3f75bSMark F. Adams          lidx++, level--){
63684d3f75bSMark F. Adams       ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr);
63784d3f75bSMark F. Adams       ierr = KSPSetOperators( mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
63884d3f75bSMark F. Adams       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
63984d3f75bSMark F. Adams       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
64084d3f75bSMark F. Adams     }
64184d3f75bSMark F. Adams     ierr = KSPSetOperators( mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
64284d3f75bSMark F. Adams 
64384d3f75bSMark F. Adams     ierr = PCSetUp_MG( pc );  CHKERRQ( ierr );
64484d3f75bSMark F. Adams   }
64584d3f75bSMark F. Adams   else if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if one level */
646d5280255SMark F. Adams     /* set default smoothers & set operators */
6479d5b6da9SMark F. Adams     for ( lidx = 1, level = pc_gamg->Nlevels-2;
648587fa25dSMark F. Adams           lidx <= fine_level;
649587fa25dSMark F. Adams           lidx++, level--) {
6505b89ad90SMark F. Adams       KSP smoother; PC subpc;
651f6536408SMark F. Adams       /* set defaults */
6529d5b6da9SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
6535b89ad90SMark F. Adams       ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
654f6536408SMark F. Adams       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
6559d5b6da9SMark F. Adams       ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr);
656f6536408SMark F. Adams       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
657d5280255SMark F. Adams       /* set ops */
658f6536408SMark F. Adams       ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
659d5280255SMark F. Adams       ierr = PCMGSetInterpolation( pc, lidx, Parr[level+1] );CHKERRQ(ierr);
660d5280255SMark F. Adams     }
661d5280255SMark F. Adams     {
662d5280255SMark F. Adams       /* coarse grid */
663d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
664d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
665d5280255SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
666d5280255SMark F. Adams       ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
667d5280255SMark F. Adams       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
668d5280255SMark F. Adams       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
669d5280255SMark F. Adams       ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
670d5280255SMark F. Adams       ierr = PCSetUp( subpc ); CHKERRQ(ierr);
671d5280255SMark F. Adams       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
672d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
673d5280255SMark F. Adams       ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
674d5280255SMark F. Adams     }
675d5280255SMark F. Adams 
676d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
677d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
678d5280255SMark F. Adams     ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr);
679d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
680d5280255SMark 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.");
681d5280255SMark F. Adams 
682d5280255SMark F. Adams     /* create cheby smoothers */
683d5280255SMark F. Adams     for ( lidx = 1, level = pc_gamg->Nlevels-2;
684d5280255SMark F. Adams           lidx <= fine_level;
685d5280255SMark F. Adams           lidx++, level--) {
686d5280255SMark F. Adams       KSP smoother;
687d5280255SMark F. Adams       PetscBool isCheb;
688d5280255SMark F. Adams 
689d5280255SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
690f6536408SMark F. Adams       /* do my own cheby */
691f6536408SMark F. Adams       ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr);
692f6536408SMark F. Adams       if( isCheb ) {
693d5280255SMark F. Adams         PetscReal emax, emin;
694d5280255SMark F. Adams         PC subpc;
695d5280255SMark F. Adams 
696d5280255SMark F. Adams         ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
6979d5b6da9SMark F. Adams         ierr = PetscTypeCompare( (PetscObject)subpc, PCJACOBI, &isCheb ); CHKERRQ(ierr);
698f6536408SMark F. Adams         if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
699587fa25dSMark F. Adams         else{ /* eigen estimate 'emax' */
700587fa25dSMark F. Adams           KSP eksp; Mat Lmat = Aarr[level];
701b2a4f308SMark F. Adams           Vec bb, xx;
702038e3b61SMark F. Adams 
7035745e0f5SMark F. Adams           ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
7045745e0f5SMark F. Adams           ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
705fc4362bfSMark F. Adams           {
706fc4362bfSMark F. Adams             PetscRandom    rctx;
707fc4362bfSMark F. Adams             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
708fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
709fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
710fc4362bfSMark F. Adams             ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
7115b89ad90SMark F. Adams           }
712fc4362bfSMark F. Adams           ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
713038e3b61SMark F. Adams           ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
714fc4362bfSMark F. Adams           CHKERRQ(ierr);
715fc4362bfSMark F. Adams           ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
716fc4362bfSMark F. Adams 
717f6536408SMark F. Adams           ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
718f6536408SMark F. Adams           ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
719f6536408SMark F. Adams 
720f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
721f6536408SMark F. Adams           ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
722fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
7235a9b9e01SMark F. Adams 
724b2a4f308SMark F. Adams           { /* set PC type to be same as smoother - does not get all parameters!!! */
725b2a4f308SMark F. Adams             const PCType type;    PC pc;
726b2a4f308SMark F. Adams             ierr = PCGetType( subpc, &type );   CHKERRQ(ierr);
727b2a4f308SMark F. Adams             ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
728b2a4f308SMark F. Adams             ierr = PCSetType( pc, type ); CHKERRQ(ierr);
729b2a4f308SMark F. Adams           }
730b2a4f308SMark F. Adams 
7315a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
7325a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
7335a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
734fc4362bfSMark F. Adams           ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
7355a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
7365a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
7375a9b9e01SMark F. Adams 
738fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
7395a9b9e01SMark F. Adams 
740fc4362bfSMark F. Adams           ierr = VecDestroy( &xx );       CHKERRQ(ierr);
741fc4362bfSMark F. Adams           ierr = VecDestroy( &bb );       CHKERRQ(ierr);
742fc4362bfSMark F. Adams           ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
743f6536408SMark F. Adams 
7449d5b6da9SMark F. Adams           if (pc_gamg->verbose) {
745b2a4f308SMark F. Adams             PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e\n",__FUNCT__,emax,emin);
746f6536408SMark F. Adams           }
747fc4362bfSMark F. Adams         }
748038e3b61SMark F. Adams         {
749038e3b61SMark F. Adams           PetscInt N1, N0, tt;
750038e3b61SMark F. Adams           ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
751038e3b61SMark F. Adams           ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
752f6536408SMark F. Adams           /* heuristic - is this crap? */
753f6536408SMark F. Adams           emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
754038e3b61SMark F. Adams           emax *= 1.05;
755038e3b61SMark F. Adams         }
756fc4362bfSMark F. Adams         ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
757f6536408SMark F. Adams       }
7585745e0f5SMark F. Adams     }
759737a81a9SMark F. Adams 
760d5280255SMark F. Adams     /* clean up */
761d5280255SMark F. Adams     for(level=1;level<pc_gamg->Nlevels;level++){
762587fa25dSMark F. Adams       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
763587fa25dSMark F. Adams       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
7645b89ad90SMark F. Adams     }
765*0cbbd2e1SMark F. Adams 
766*0cbbd2e1SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
767*0cbbd2e1SMark F. Adams 
76858471d46SMark F. Adams     {
76958471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
7709d5b6da9SMark F. Adams       ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
77158471d46SMark F. Adams       ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
77258471d46SMark F. Adams     }
7735f8cf99dSMark F. Adams   }
7745f8cf99dSMark F. Adams   else {
7755f8cf99dSMark F. Adams     KSP smoother;
776b43b03e9SMark 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__);
7779d5b6da9SMark F. Adams     ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
7785f8cf99dSMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
7795f8cf99dSMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
7809d5b6da9SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
7815f8cf99dSMark F. Adams   }
78284d3f75bSMark F. Adams 
7835b89ad90SMark F. Adams   PetscFunctionReturn(0);
7845b89ad90SMark F. Adams }
7855b89ad90SMark F. Adams 
786eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7875b89ad90SMark F. Adams /*
7885b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7895b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7905b89ad90SMark F. Adams 
7915b89ad90SMark F. Adams    Input Parameter:
7925b89ad90SMark F. Adams .  pc - the preconditioner context
7935b89ad90SMark F. Adams 
7945b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7955b89ad90SMark F. Adams */
7965b89ad90SMark F. Adams #undef __FUNCT__
7975b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
7985b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG( PC pc )
7995b89ad90SMark F. Adams {
8005b89ad90SMark F. Adams   PetscErrorCode  ierr;
8015b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
8025b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
8035b89ad90SMark F. Adams 
8045b89ad90SMark F. Adams   PetscFunctionBegin;
8055b89ad90SMark F. Adams   ierr = PCReset_GAMG( pc );CHKERRQ(ierr);
8065b89ad90SMark F. Adams   ierr = PetscFree( pc_gamg );CHKERRQ(ierr);
8075b89ad90SMark F. Adams   ierr = PCDestroy_MG( pc );CHKERRQ(ierr);
8085b89ad90SMark F. Adams   PetscFunctionReturn(0);
8095b89ad90SMark F. Adams }
8105b89ad90SMark F. Adams 
811676e1743SMark F. Adams 
812676e1743SMark F. Adams #undef __FUNCT__
813676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
814676e1743SMark F. Adams /*@
815676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
816676e1743SMark F. Adams    processor reduction.
817676e1743SMark F. Adams 
818676e1743SMark F. Adams    Not Collective on PC
819676e1743SMark F. Adams 
820676e1743SMark F. Adams    Input Parameters:
821676e1743SMark F. Adams .  pc - the preconditioner context
822676e1743SMark F. Adams 
823676e1743SMark F. Adams 
824676e1743SMark F. Adams    Options Database Key:
825676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
826676e1743SMark F. Adams 
827676e1743SMark F. Adams    Level: intermediate
828676e1743SMark F. Adams 
829676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
830676e1743SMark F. Adams 
831676e1743SMark F. Adams .seealso: ()
832676e1743SMark F. Adams @*/
833676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
834676e1743SMark F. Adams {
835676e1743SMark F. Adams   PetscErrorCode ierr;
836676e1743SMark F. Adams 
837676e1743SMark F. Adams   PetscFunctionBegin;
838676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
839676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
840676e1743SMark F. Adams   PetscFunctionReturn(0);
841676e1743SMark F. Adams }
842676e1743SMark F. Adams 
843676e1743SMark F. Adams EXTERN_C_BEGIN
844676e1743SMark F. Adams #undef __FUNCT__
845676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
846676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
847676e1743SMark F. Adams {
848c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
849c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
850676e1743SMark F. Adams 
851676e1743SMark F. Adams   PetscFunctionBegin;
8529d5b6da9SMark F. Adams   if(n>0) pc_gamg->min_eq_proc = n;
853676e1743SMark F. Adams   PetscFunctionReturn(0);
854676e1743SMark F. Adams }
855676e1743SMark F. Adams EXTERN_C_END
856676e1743SMark F. Adams 
857676e1743SMark F. Adams #undef __FUNCT__
858389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
859389730f3SMark F. Adams /*@
860389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
861389730f3SMark F. Adams 
862389730f3SMark F. Adams  Collective on PC
863389730f3SMark F. Adams 
864389730f3SMark F. Adams    Input Parameters:
865389730f3SMark F. Adams .  pc - the preconditioner context
866389730f3SMark F. Adams 
867389730f3SMark F. Adams 
868389730f3SMark F. Adams    Options Database Key:
869389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
870389730f3SMark F. Adams 
871389730f3SMark F. Adams    Level: intermediate
872389730f3SMark F. Adams 
873389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
874389730f3SMark F. Adams 
875389730f3SMark F. Adams .seealso: ()
876389730f3SMark F. Adams  @*/
877389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
878389730f3SMark F. Adams {
879389730f3SMark F. Adams   PetscErrorCode ierr;
880389730f3SMark F. Adams 
881389730f3SMark F. Adams   PetscFunctionBegin;
882389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
883389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
884389730f3SMark F. Adams   PetscFunctionReturn(0);
885389730f3SMark F. Adams }
886389730f3SMark F. Adams 
887389730f3SMark F. Adams EXTERN_C_BEGIN
888389730f3SMark F. Adams #undef __FUNCT__
889389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
890389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
891389730f3SMark F. Adams {
892389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
893389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
894389730f3SMark F. Adams 
895389730f3SMark F. Adams   PetscFunctionBegin;
8969d5b6da9SMark F. Adams   if(n>0) pc_gamg->coarse_eq_limit = n;
897389730f3SMark F. Adams   PetscFunctionReturn(0);
898389730f3SMark F. Adams }
899389730f3SMark F. Adams EXTERN_C_END
900389730f3SMark F. Adams 
901389730f3SMark F. Adams #undef __FUNCT__
9028263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
903676e1743SMark F. Adams /*@
9048263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
905676e1743SMark F. Adams 
906676e1743SMark F. Adams    Collective on PC
907676e1743SMark F. Adams 
908676e1743SMark F. Adams    Input Parameters:
909676e1743SMark F. Adams .  pc - the preconditioner context
910676e1743SMark F. Adams 
911676e1743SMark F. Adams 
912676e1743SMark F. Adams    Options Database Key:
9138263b398SMark F. Adams .  -pc_gamg_repartition
914676e1743SMark F. Adams 
915676e1743SMark F. Adams    Level: intermediate
916676e1743SMark F. Adams 
917676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
918676e1743SMark F. Adams 
919676e1743SMark F. Adams .seealso: ()
920676e1743SMark F. Adams @*/
9218263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
922676e1743SMark F. Adams {
923676e1743SMark F. Adams   PetscErrorCode ierr;
924676e1743SMark F. Adams 
925676e1743SMark F. Adams   PetscFunctionBegin;
926676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9278263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
928676e1743SMark F. Adams   PetscFunctionReturn(0);
929676e1743SMark F. Adams }
930676e1743SMark F. Adams 
931676e1743SMark F. Adams EXTERN_C_BEGIN
932676e1743SMark F. Adams #undef __FUNCT__
9338263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
9348263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
935676e1743SMark F. Adams {
936c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
937c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
938676e1743SMark F. Adams 
939676e1743SMark F. Adams   PetscFunctionBegin;
9409d5b6da9SMark F. Adams   pc_gamg->repart = n;
941676e1743SMark F. Adams   PetscFunctionReturn(0);
942676e1743SMark F. Adams }
943676e1743SMark F. Adams EXTERN_C_END
944676e1743SMark F. Adams 
945676e1743SMark F. Adams #undef __FUNCT__
9464ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
9474ef23d27SMark F. Adams /*@
9484ef23d27SMark F. Adams    PCGAMGSetNlevels -
9494ef23d27SMark F. Adams 
9504ef23d27SMark F. Adams    Not collective on PC
9514ef23d27SMark F. Adams 
9524ef23d27SMark F. Adams    Input Parameters:
9534ef23d27SMark F. Adams .  pc - the preconditioner context
9544ef23d27SMark F. Adams 
9554ef23d27SMark F. Adams 
9564ef23d27SMark F. Adams    Options Database Key:
9574ef23d27SMark F. Adams .  -pc_mg_levels
9584ef23d27SMark F. Adams 
9594ef23d27SMark F. Adams    Level: intermediate
9604ef23d27SMark F. Adams 
9614ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
9624ef23d27SMark F. Adams 
9634ef23d27SMark F. Adams .seealso: ()
9644ef23d27SMark F. Adams @*/
9654ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
9664ef23d27SMark F. Adams {
9674ef23d27SMark F. Adams   PetscErrorCode ierr;
9684ef23d27SMark F. Adams 
9694ef23d27SMark F. Adams   PetscFunctionBegin;
9704ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9714ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
9724ef23d27SMark F. Adams   PetscFunctionReturn(0);
9734ef23d27SMark F. Adams }
9744ef23d27SMark F. Adams 
9754ef23d27SMark F. Adams EXTERN_C_BEGIN
9764ef23d27SMark F. Adams #undef __FUNCT__
9774ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
9784ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
9794ef23d27SMark F. Adams {
9804ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
9814ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
9824ef23d27SMark F. Adams 
9834ef23d27SMark F. Adams   PetscFunctionBegin;
9849d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
9854ef23d27SMark F. Adams   PetscFunctionReturn(0);
9864ef23d27SMark F. Adams }
9874ef23d27SMark F. Adams EXTERN_C_END
9884ef23d27SMark F. Adams 
9894ef23d27SMark F. Adams #undef __FUNCT__
9903542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
9913542efc5SMark F. Adams /*@
9923542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
9933542efc5SMark F. Adams 
9943542efc5SMark F. Adams    Not collective on PC
9953542efc5SMark F. Adams 
9963542efc5SMark F. Adams    Input Parameters:
9973542efc5SMark F. Adams .  pc - the preconditioner context
9983542efc5SMark F. Adams 
9993542efc5SMark F. Adams 
10003542efc5SMark F. Adams    Options Database Key:
10013542efc5SMark F. Adams .  -pc_gamg_threshold
10023542efc5SMark F. Adams 
10033542efc5SMark F. Adams    Level: intermediate
10043542efc5SMark F. Adams 
10053542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
10063542efc5SMark F. Adams 
10073542efc5SMark F. Adams .seealso: ()
10083542efc5SMark F. Adams @*/
10093542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
10103542efc5SMark F. Adams {
10113542efc5SMark F. Adams   PetscErrorCode ierr;
10123542efc5SMark F. Adams 
10133542efc5SMark F. Adams   PetscFunctionBegin;
10143542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10153542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
10163542efc5SMark F. Adams   PetscFunctionReturn(0);
10173542efc5SMark F. Adams }
10183542efc5SMark F. Adams 
10193542efc5SMark F. Adams EXTERN_C_BEGIN
10203542efc5SMark F. Adams #undef __FUNCT__
10213542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
10223542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
10233542efc5SMark F. Adams {
1024c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1025c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
10263542efc5SMark F. Adams 
10273542efc5SMark F. Adams   PetscFunctionBegin;
10289d5b6da9SMark F. Adams   pc_gamg->threshold = n;
10293542efc5SMark F. Adams   PetscFunctionReturn(0);
10303542efc5SMark F. Adams }
10313542efc5SMark F. Adams EXTERN_C_END
10323542efc5SMark F. Adams 
10333542efc5SMark F. Adams #undef __FUNCT__
10349d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1035676e1743SMark F. Adams /*@
10369d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1037676e1743SMark F. Adams 
1038676e1743SMark F. Adams    Collective on PC
1039676e1743SMark F. Adams 
1040676e1743SMark F. Adams    Input Parameters:
1041676e1743SMark F. Adams .  pc - the preconditioner context
1042676e1743SMark F. Adams 
1043676e1743SMark F. Adams 
1044676e1743SMark F. Adams    Options Database Key:
10453542efc5SMark F. Adams .  -pc_gamg_type
1046676e1743SMark F. Adams 
1047676e1743SMark F. Adams    Level: intermediate
1048676e1743SMark F. Adams 
1049676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1050676e1743SMark F. Adams 
1051676e1743SMark F. Adams .seealso: ()
1052676e1743SMark F. Adams @*/
10530202ef74SSatish Balay PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type )
1054676e1743SMark F. Adams {
1055676e1743SMark F. Adams   PetscErrorCode ierr;
1056676e1743SMark F. Adams 
1057676e1743SMark F. Adams   PetscFunctionBegin;
1058676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10590202ef74SSatish Balay   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type));
1060676e1743SMark F. Adams   CHKERRQ(ierr);
1061676e1743SMark F. Adams   PetscFunctionReturn(0);
1062676e1743SMark F. Adams }
1063676e1743SMark F. Adams 
1064676e1743SMark F. Adams EXTERN_C_BEGIN
1065676e1743SMark F. Adams #undef __FUNCT__
10669d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
10670202ef74SSatish Balay PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type )
1068676e1743SMark F. Adams {
10699d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
1070676e1743SMark F. Adams 
1071676e1743SMark F. Adams   PetscFunctionBegin;
10729d5b6da9SMark F. Adams   ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r);
10739d5b6da9SMark F. Adams   CHKERRQ(ierr);
10749d5b6da9SMark F. Adams 
10759d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
10769d5b6da9SMark F. Adams 
10779d5b6da9SMark F. Adams   /* call sub create method */
10789d5b6da9SMark F. Adams   ierr = (*r)(pc); CHKERRQ(ierr);
10799d5b6da9SMark F. Adams 
1080676e1743SMark F. Adams   PetscFunctionReturn(0);
1081676e1743SMark F. Adams }
1082676e1743SMark F. Adams EXTERN_C_END
1083676e1743SMark F. Adams 
10845b89ad90SMark F. Adams #undef __FUNCT__
10855b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
10865b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG( PC pc )
10875b89ad90SMark F. Adams {
1088676e1743SMark F. Adams   PetscErrorCode  ierr;
1089676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1090676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1091676e1743SMark F. Adams   PetscBool        flag;
1092b43b03e9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
10935b89ad90SMark F. Adams 
10945b89ad90SMark F. Adams   PetscFunctionBegin;
1095676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1096676e1743SMark F. Adams   {
109775b74bdaSMark F. Adams     /* -pc_gamg_verbose */
10989d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
10999d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
11009d5b6da9SMark F. Adams                            &pc_gamg->verbose, PETSC_NULL );
11019d5b6da9SMark F. Adams     CHKERRQ(ierr);
110275b74bdaSMark F. Adams 
11038263b398SMark F. Adams     /* -pc_gamg_repartition */
11048263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
11058263b398SMark F. Adams                             "Repartion coarse grids (false)",
11068263b398SMark F. Adams                             "PCGAMGRepartitioning",
11079d5b6da9SMark F. Adams                             pc_gamg->repart,
11089d5b6da9SMark F. Adams                             &pc_gamg->repart,
1109676e1743SMark F. Adams                             &flag);
1110676e1743SMark F. Adams     CHKERRQ(ierr);
1111676e1743SMark F. Adams 
1112c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1113676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1114676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1115676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
11169d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
11179d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1118676e1743SMark F. Adams                            &flag );
1119676e1743SMark F. Adams     CHKERRQ(ierr);
11203542efc5SMark F. Adams 
1121389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1122389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1123389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1124389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
11259d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
11269d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1127389730f3SMark F. Adams                            &flag );
1128389730f3SMark F. Adams     CHKERRQ(ierr);
1129389730f3SMark F. Adams 
1130c20e4228SMark F. Adams     /* -pc_gamg_threshold */
11313542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
11323542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
11333542efc5SMark F. Adams                             "PCGAMGSetThreshold",
11349d5b6da9SMark F. Adams                             pc_gamg->threshold,
11359d5b6da9SMark F. Adams                             &pc_gamg->threshold,
11363542efc5SMark F. Adams                             &flag );
11373542efc5SMark F. Adams     CHKERRQ(ierr);
1138b43b03e9SMark F. Adams     if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);
11394ef23d27SMark F. Adams 
11404ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
11414ef23d27SMark F. Adams                            "Set number of MG levels",
11424ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
11439d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
11449d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
11454ef23d27SMark F. Adams                            &flag );
1146676e1743SMark F. Adams   }
1147676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1148676e1743SMark F. Adams 
11495b89ad90SMark F. Adams   PetscFunctionReturn(0);
11505b89ad90SMark F. Adams }
11515b89ad90SMark F. Adams 
11525b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
11535b89ad90SMark F. Adams /*MC
11549d5b6da9SMark F. Adams      PCCreate_GAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
11559d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
11565b89ad90SMark F. Adams 
1157280d9858SJed Brown    Options Database Keys:
11585b89ad90SMark F. Adams    Multigrid options(inherited)
1159280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1160280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1161280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1162280d9858SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
11635b89ad90SMark F. Adams 
11645b89ad90SMark F. Adams   Level: intermediate
1165280d9858SJed Brown 
11665b89ad90SMark F. Adams   Concepts: multigrid
11675b89ad90SMark F. Adams 
11685b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1169280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
11705b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
11715b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
11725b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11735b89ad90SMark F. Adams M*/
11745b89ad90SMark F. Adams EXTERN_C_BEGIN
11755b89ad90SMark F. Adams #undef __FUNCT__
11765b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
11775b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG( PC pc )
11785b89ad90SMark F. Adams {
11795b89ad90SMark F. Adams   PetscErrorCode  ierr;
11805b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
11815b89ad90SMark F. Adams   PC_MG           *mg;
11829d5b6da9SMark F. Adams 
1183*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11849d5b6da9SMark F. Adams   static long count = 0;
11855ee9c036SSatish Balay #endif
11865b89ad90SMark F. Adams 
11875b89ad90SMark F. Adams   PetscFunctionBegin;
11885b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
11895b89ad90SMark F. Adams   ierr = PCSetType( pc, PCMG );  CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
11905b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr);
11915b89ad90SMark F. Adams 
11925b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
11935b89ad90SMark F. Adams   ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr);
11945b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
1195ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
11965b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
11975b89ad90SMark F. Adams 
11989d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
11999d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
12009d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
12019d5b6da9SMark F. Adams   pc_gamg->data = 0;
12025b89ad90SMark F. Adams 
12039d5b6da9SMark F. Adams   /* register AMG type */
12049d5b6da9SMark F. Adams   if( !GAMGList ){
12059d5b6da9SMark F. Adams     ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
12069d5b6da9SMark F. Adams     ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
12079d5b6da9SMark F. Adams   }
12089d5b6da9SMark F. Adams 
12099d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
12105b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
12115b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
12125b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
12135b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
12145b89ad90SMark F. Adams 
12155b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1216676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1217676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1218676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1219676e1743SMark F. Adams   CHKERRQ(ierr);
1220676e1743SMark F. Adams 
1221676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1222389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_C",
1223389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_GAMG",
1224389730f3SMark F. Adams 					    PCGAMGSetCoarseEqLim_GAMG);
1225389730f3SMark F. Adams   CHKERRQ(ierr);
1226389730f3SMark F. Adams 
1227389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12288263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_C",
12298263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_GAMG",
12308263b398SMark F. Adams 					    PCGAMGSetRepartitioning_GAMG);
1231676e1743SMark F. Adams   CHKERRQ(ierr);
1232676e1743SMark F. Adams 
1233676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12343542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
12353542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
12363542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
12373542efc5SMark F. Adams   CHKERRQ(ierr);
12383542efc5SMark F. Adams 
12393542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12409d5b6da9SMark F. Adams 					    "PCGAMGSetType_C",
12419d5b6da9SMark F. Adams 					    "PCGAMGSetType_GAMG",
12429d5b6da9SMark F. Adams 					    PCGAMGSetType_GAMG);
1243676e1743SMark F. Adams   CHKERRQ(ierr);
1244c97e1df0SMark F. Adams 
12459d5b6da9SMark F. Adams   pc_gamg->repart = PETSC_FALSE;
12469d5b6da9SMark F. Adams   pc_gamg->min_eq_proc = 100;
1247c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit = 800;
12489d5b6da9SMark F. Adams   pc_gamg->threshold = 0.05;
12499d5b6da9SMark F. Adams   pc_gamg->Nlevels = GAMG_MAXLEVELS;
12509d5b6da9SMark F. Adams   pc_gamg->verbose = 0;
12519d5b6da9SMark F. Adams   pc_gamg->emax_id = -1;
12529d5b6da9SMark F. Adams   pc_gamg->col_bs_id = -1;
12539d5b6da9SMark F. Adams 
12549d5b6da9SMark F. Adams   /* instantiate derived type */
12559d5b6da9SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
12569d5b6da9SMark F. Adams   {
12579d5b6da9SMark F. Adams     char tname[256] = GAMGAGG;
12589d5b6da9SMark F. Adams     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",
12599d5b6da9SMark F. Adams                             GAMGList, tname, tname, sizeof(tname), PETSC_NULL );
12609d5b6da9SMark F. Adams     CHKERRQ(ierr);
12619d5b6da9SMark F. Adams     ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr);
12629d5b6da9SMark F. Adams   }
12639d5b6da9SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
12649d5b6da9SMark F. Adams 
1265*0cbbd2e1SMark F. Adams   /* private events */
1266*0cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1267785cba28SMark F. Adams   if( count++ == 0 ) {
1268*0cbbd2e1SMark F. Adams     PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1269*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1270*0cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1271*0cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1272*0cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1273*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1274*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1275*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1276*0cbbd2e1SMark F. Adams     PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1277*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  SA: colect data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1278*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1279*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1280*0cbbd2e1SMark F. Adams     PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1281*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1282*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1283*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1284*0cbbd2e1SMark F. Adams     PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);
1285f852f58cSMark F. Adams 
1286*0cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1287*0cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1288*0cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1289b4fbaa2aSMark F. Adams     /* create timer stages */
1290b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1291b4fbaa2aSMark F. Adams     {
1292b4fbaa2aSMark F. Adams       char str[32];
1293b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1294b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1295b4fbaa2aSMark F. Adams       PetscInt lidx;
1296b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1297b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1298b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1299b4fbaa2aSMark F. Adams       }
1300b4fbaa2aSMark F. Adams     }
1301b4fbaa2aSMark F. Adams #endif
1302b4fbaa2aSMark F. Adams   }
1303b4fbaa2aSMark F. Adams #endif
1304*0cbbd2e1SMark F. Adams   /* general events */
1305*0cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
1306*0cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGgraph_AGG", PC_CLASSID, &PC_GAMGGgraph_AGG);
1307*0cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);
1308*0cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1309*0cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1310*0cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1311*0cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1312*0cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGProlOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);
1313*0cbbd2e1SMark F. Adams #endif
1314*0cbbd2e1SMark F. Adams 
13155b89ad90SMark F. Adams   PetscFunctionReturn(0);
13165b89ad90SMark F. Adams }
13175b89ad90SMark F. Adams EXTERN_C_END
1318