xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 084a8fe3df532bd8fd2af12244978a79ac3b3f5a)
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 
8b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
9b4fbaa2aSMark F. Adams PetscLogEvent gamg_setup_events[NUM_SET];
10b4fbaa2aSMark F. Adams #endif
11b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
12b4fbaa2aSMark F. Adams 
13b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
14b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
15b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
16b4fbaa2aSMark F. Adams #endif
17f96513f1SMatthew G Knepley 
189d5b6da9SMark F. Adams static PetscFList GAMGList = 0;
199d5b6da9SMark F. Adams 
20d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
21d3d6bff4SMark F. Adams #undef __FUNCT__
22d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
23d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
24d3d6bff4SMark F. Adams {
25d3d6bff4SMark F. Adams   PetscErrorCode  ierr;
26d3d6bff4SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
27d3d6bff4SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
28d3d6bff4SMark F. Adams 
29d3d6bff4SMark F. Adams   PetscFunctionBegin;
309d5b6da9SMark F. Adams   if( pc_gamg->data != 0 ) { /* this should not happen, cleaned up in SetUp */
319d5b6da9SMark F. Adams     ierr = PetscFree(pc_gamg->data); CHKERRQ(ierr);
3258471d46SMark F. Adams   }
339d5b6da9SMark F. Adams   pc_gamg->data = 0; pc_gamg->data_sz = 0;
34d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
35d3d6bff4SMark F. Adams }
36d3d6bff4SMark F. Adams 
375b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
385b89ad90SMark F. Adams /*
39a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
40a147abb0SMark F. Adams      of active processors.
415b89ad90SMark F. Adams 
425b89ad90SMark F. Adams    Input Parameter:
439d5b6da9SMark F. Adams    . pc - parameters
449d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
459d5b6da9SMark F. Adams    . cbs - coarse block size
469d5b6da9SMark F. Adams    . isLast -
473530afc2SMark F. Adams    In/Output Parameter:
483530afc2SMark F. Adams    . a_P_inout - prolongation operator to the next level (k-1)
49afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
5011e60469SMark F. Adams    Output Parameter:
513530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
525b89ad90SMark F. Adams */
535cb416c2SMark F. Adams 
545b89ad90SMark F. Adams #undef __FUNCT__
558263b398SMark F. Adams #define __FUNCT__ "createLevel"
569d5b6da9SMark F. Adams PetscErrorCode createLevel( const PC pc,
579d5b6da9SMark F. Adams                             const Mat Amat_fine,
589d5b6da9SMark F. Adams                             const PetscInt cbs,
599d5b6da9SMark F. Adams                             const PetscBool isLast,
603530afc2SMark F. Adams                             Mat *a_P_inout,
61afc97cdcSMark F. Adams                             PetscMPIInt *a_nactive_proc,
62389730f3SMark F. Adams                             Mat *a_Amat_crs
6311e60469SMark F. Adams                             )
645b89ad90SMark F. Adams {
659d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
66486a8d0bSJed Brown   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
679d5b6da9SMark F. Adams   const PetscBool  repart = pc_gamg->repart;
689d5b6da9SMark F. Adams   const PetscInt   min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
695b89ad90SMark F. Adams   PetscErrorCode   ierr;
70038e3b61SMark F. Adams   Mat              Cmat,Pnew,Pold=*a_P_inout;
7111e60469SMark F. Adams   IS               new_indices,isnum;
729d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)Amat_fine)->comm;
735a9b9e01SMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive = *a_nactive_proc;
745a9b9e01SMark F. Adams   PetscInt         neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0,rfactor;
755b89ad90SMark F. Adams 
765b89ad90SMark F. Adams   PetscFunctionBegin;
7711e60469SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
7811e60469SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
79f6536408SMark F. Adams 
8011e60469SMark F. Adams   /* RAP */
8196434bc1SMark F. Adams #ifdef USE_R
8296434bc1SMark F. Adams   /* make R wih brute force for now */
8396434bc1SMark F. Adams   ierr = MatTranspose( Pold, Pnew );
8496434bc1SMark F. Adams   ierr = MatDestroy( &Pold );  CHKERRQ(ierr);
859d5b6da9SMark F. Adams   ierr = MatRARt( Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
8696434bc1SMark F. Adams   Pold = Pnew;
8796434bc1SMark F. Adams #else
889d5b6da9SMark F. Adams   ierr = MatPtAP( Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
8996434bc1SMark F. Adams #endif
909d5b6da9SMark F. Adams   ierr = MatSetBlockSize( Cmat, cbs );      CHKERRQ(ierr);
91038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
929d5b6da9SMark F. Adams   ncrs0 = (Iend0-Istart0)/cbs; assert((Iend0-Istart0)%cbs == 0);
93038e3b61SMark F. Adams 
94f852f58cSMark F. Adams   /* get number of PEs to make active, reduce */
95f852f58cSMark F. Adams   ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
9681708dccSMark F. Adams   new_npe = (PetscMPIInt)((float)neq/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
97f852f58cSMark F. Adams   if( new_npe == 0 || neq < coarse_max ) new_npe = 1;
985a9b9e01SMark F. Adams   else if (new_npe >= nactive ) new_npe = nactive; /* no change, rare */
999d5b6da9SMark F. Adams   if( isLast ) new_npe = 1;
100f852f58cSMark F. Adams 
1015a9b9e01SMark F. Adams   if( !repart && new_npe==nactive ) {
102a147abb0SMark F. Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction */
10322063be5SMark F. Adams   }
10422063be5SMark F. Adams   else {
10522063be5SMark F. Adams     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
106c8b0795cSMark 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;
107c8b0795cSMark F. Adams     const PetscInt  stride0=ncrs0*pc_gamg->data_cell_rows;
1085a9b9e01SMark F. Adams     PetscInt        *counts,is_sz,*newproc_idx,ii,jj,kk,strideNew,*tidx,inpe,targetPE;
1095a9b9e01SMark F. Adams     IS              isnewproc;
1105a9b9e01SMark F. Adams     VecScatter      vecscat;
11122063be5SMark F. Adams     PetscScalar    *array;
11222063be5SMark F. Adams     Vec             src_crd, dest_crd;
1135a9b9e01SMark F. Adams     IS              isscat;
114e33ef3b1SMark F. Adams 
115fd3c6afaSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
116b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
117b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
118b4fbaa2aSMark F. Adams #endif
1195a9b9e01SMark F. Adams     if( repart ) {
1205a9b9e01SMark F. Adams       /* create sub communicator  */
1215a9b9e01SMark F. Adams       Mat             adj;
1225a9b9e01SMark F. Adams 
123b43b03e9SMark 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);
1245a9b9e01SMark F. Adams 
1255a9b9e01SMark F. Adams       /* MatPartitioningApply call MatConvert, which is collective */
1269d5b6da9SMark F. Adams       if( cbs == 1 ) {
127038e3b61SMark F. Adams 	ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
128eb07cef2SMark F. Adams       }
129eb07cef2SMark F. Adams       else{
130038e3b61SMark F. Adams 	/* make a scalar matrix to partition */
131eb07cef2SMark F. Adams 	Mat tMat;
13258471d46SMark F. Adams 	PetscInt ncols,jj,Ii;
133b4fbaa2aSMark F. Adams 	const PetscScalar *vals;
134b4fbaa2aSMark F. Adams 	const PetscInt *idx;
1355f8cf99dSMark F. Adams 	PetscInt *d_nnz, *o_nnz;
1369057884aSMark F. Adams 	static PetscInt llev = 0;
137b4fbaa2aSMark F. Adams 
13858471d46SMark F. Adams 	ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
1395f8cf99dSMark F. Adams 	ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr);
1409d5b6da9SMark F. Adams 	for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += cbs, jj++ ) {
14158471d46SMark F. Adams 	  ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
1429d5b6da9SMark F. Adams 	  d_nnz[jj] = ncols/cbs;
1439d5b6da9SMark F. Adams 	  o_nnz[jj] = ncols/cbs;
14458471d46SMark F. Adams 	  ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
1455f8cf99dSMark F. Adams 	  if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
1469d5b6da9SMark F. Adams 	  if( o_nnz[jj] > (neq/cbs-ncrs0) ) o_nnz[jj] = neq/cbs-ncrs0;
14758471d46SMark F. Adams 	}
1486876a03eSMark F. Adams 
149eb07cef2SMark F. Adams 	ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
150eb07cef2SMark F. Adams 				PETSC_DETERMINE, PETSC_DETERMINE,
1515f8cf99dSMark F. Adams 				0, d_nnz, 0, o_nnz,
152eb07cef2SMark F. Adams 				&tMat );
1536876a03eSMark F. Adams 	CHKERRQ(ierr);
15458471d46SMark F. Adams 	ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
1555f8cf99dSMark F. Adams 	ierr = PetscFree( o_nnz ); CHKERRQ(ierr);
156eb07cef2SMark F. Adams 
15722063be5SMark F. Adams 	for ( ii = Istart0; ii < Iend0; ii++ ) {
1589d5b6da9SMark F. Adams 	  PetscInt dest_row = ii/cbs;
15922063be5SMark F. Adams 	  ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
160eb07cef2SMark F. Adams 	  for( jj = 0 ; jj < ncols ; jj++ ){
1619d5b6da9SMark F. Adams 	    PetscInt dest_col = idx[jj]/cbs;
162eb07cef2SMark F. Adams 	    PetscScalar v = 1.0;
163eb07cef2SMark F. Adams 	    ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
164eb07cef2SMark F. Adams 	  }
16522063be5SMark F. Adams 	  ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
166eb07cef2SMark F. Adams 	}
167eb07cef2SMark F. Adams 	ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168eb07cef2SMark F. Adams 	ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169eb07cef2SMark F. Adams 
170b4fbaa2aSMark F. Adams 	if( llev++ == -1 ) {
171b4fbaa2aSMark F. Adams 	  PetscViewer viewer; char fname[32];
172b4fbaa2aSMark F. Adams 	  sprintf(fname,"part_mat_%d.mat",llev);
173b4fbaa2aSMark F. Adams 	  PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
174b4fbaa2aSMark F. Adams 	  ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
175b4fbaa2aSMark F. Adams 	  ierr = PetscViewerDestroy( &viewer );
176b4fbaa2aSMark F. Adams 	}
177b4fbaa2aSMark F. Adams 
178eb07cef2SMark F. Adams 	ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
179eb07cef2SMark F. Adams 
180eb07cef2SMark F. Adams 	ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
181eb07cef2SMark F. Adams       }
182f150b916SMark F. Adams 
1835a9b9e01SMark F. Adams       { /* partition: get newproc_idx, set is_sz */
1845a9b9e01SMark F. Adams 	char prefix[256];
1855a9b9e01SMark F. Adams 	const char *pcpre;
186b4fbaa2aSMark F. Adams 	const PetscInt *is_idx;
187b4fbaa2aSMark F. Adams 	MatPartitioning  mpart;
188a4b7d37bSMark F. Adams 	IS proc_is;
1892f03bc48SMark F. Adams 
1905a9b9e01SMark F. Adams 	ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr);
1915ef31b24SMark F. Adams 	ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
1929d5b6da9SMark F. Adams 	ierr = PCGetOptionsPrefix( pc, &pcpre );CHKERRQ(ierr);
19359a0be82SJed Brown 	ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
19459a0be82SJed Brown 	ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
19511e60469SMark F. Adams 	ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
1965ef31b24SMark F. Adams 	ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
197a4b7d37bSMark F. Adams 	ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr);
19811e60469SMark F. Adams 	ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
1995a9b9e01SMark F. Adams 
2005ef31b24SMark F. Adams 	/* collect IS info */
201a4b7d37bSMark F. Adams 	ierr = ISGetLocalSize( proc_is, &is_sz );       CHKERRQ(ierr);
2029d5b6da9SMark F. Adams 	ierr = PetscMalloc( cbs*is_sz*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr);
203a4b7d37bSMark F. Adams 	ierr = ISGetIndices( proc_is, &is_idx );        CHKERRQ(ierr);
204a147abb0SMark F. Adams 	NN = 1; /* bring to "front" of machine */
205a147abb0SMark F. Adams 	/*NN = npe/new_npe;*/ /* spread partitioning across machine */
206b4fbaa2aSMark F. Adams 	for( kk = jj = 0 ; kk < is_sz ; kk++ ){
2079d5b6da9SMark F. Adams 	  for( ii = 0 ; ii < cbs ; ii++, jj++ ){
2088263b398SMark F. Adams 	    newproc_idx[jj] = is_idx[kk] * NN; /* distribution */
209eb07cef2SMark F. Adams 	  }
2105ef31b24SMark F. Adams 	}
211a4b7d37bSMark F. Adams 	ierr = ISRestoreIndices( proc_is, &is_idx );     CHKERRQ(ierr);
212a4b7d37bSMark F. Adams 	ierr = ISDestroy( &proc_is );                    CHKERRQ(ierr);
213b4fbaa2aSMark F. Adams 
2149d5b6da9SMark F. Adams 	is_sz *= cbs;
2155ef31b24SMark F. Adams       }
2165ef31b24SMark F. Adams       ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
2175a9b9e01SMark F. Adams 
2188263b398SMark F. Adams       ierr = ISCreateGeneral( wcomm, is_sz, newproc_idx, PETSC_COPY_VALUES, &isnewproc );
2195a9b9e01SMark F. Adams       CHKERRQ(ierr);
2208263b398SMark F. Adams       if( newproc_idx != 0 ) {
2218263b398SMark F. Adams 	ierr = PetscFree( newproc_idx );  CHKERRQ(ierr);
2225ef31b24SMark F. Adams       }
2235a9b9e01SMark F. Adams     }
2245a9b9e01SMark F. Adams     else { /* simple aggreagtion of parts */
2255a9b9e01SMark F. Adams       /* find factor */
2265a9b9e01SMark F. Adams       if( new_npe == 1 ) rfactor = npe; /* easy */
2275a9b9e01SMark F. Adams       else {
2285a9b9e01SMark F. Adams 	PetscReal best_fact = 0.;
2295a9b9e01SMark F. Adams 	jj = -1;
2305a9b9e01SMark F. Adams 	for( kk = 1 ; kk <= npe ; kk++ ){
2315a9b9e01SMark F. Adams 	  if( npe%kk==0 ) { /* a candidate */
2325a9b9e01SMark F. Adams 	    PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe;
2335a9b9e01SMark F. Adams 	    if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2345a9b9e01SMark F. Adams 	    if( fact > best_fact ) {
2355a9b9e01SMark F. Adams 	      best_fact = fact; jj = kk;
2365a9b9e01SMark F. Adams 	    }
2375a9b9e01SMark F. Adams 	  }
2385a9b9e01SMark F. Adams 	}
2395a9b9e01SMark F. Adams 	if(jj!=-1) rfactor = jj;
2405a9b9e01SMark F. Adams 	else rfactor = 1; /* prime? */
2415a9b9e01SMark F. Adams       }
2425a9b9e01SMark F. Adams       new_npe = npe/rfactor;
2435a9b9e01SMark F. Adams 
2445a9b9e01SMark F. Adams       if( new_npe==nactive ) {
2455a9b9e01SMark F. Adams 	*a_Amat_crs = Cmat; /* output - no repartitioning or reduction */
2465a9b9e01SMark F. Adams 	ierr = PetscFree( counts );  CHKERRQ(ierr);
2475a9b9e01SMark F. Adams 	*a_nactive_proc = new_npe; /* output */
248b43b03e9SMark 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);
2495a9b9e01SMark F. Adams 	PetscFunctionReturn(0);
2505a9b9e01SMark F. Adams       }
2515a9b9e01SMark F. Adams 
2525a9b9e01SMark F. Adams       /* reduce - isnewproc */
2535a9b9e01SMark F. Adams       targetPE = mype/rfactor;
2545a9b9e01SMark F. Adams 
255b43b03e9SMark 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);
2569d5b6da9SMark F. Adams       is_sz = ncrs0*cbs;
2575a9b9e01SMark F. Adams       ierr = ISCreateStride( wcomm, is_sz, targetPE, 0, &isnewproc );
2585a9b9e01SMark F. Adams       CHKERRQ(ierr);
2595a9b9e01SMark F. Adams     }
260e33ef3b1SMark F. Adams 
26111e60469SMark F. Adams     /*
26211e60469SMark F. Adams       Create an index set from the isnewproc index set to indicate the mapping TO
26311e60469SMark F. Adams     */
26411e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
26511e60469SMark F. Adams     /*
26611e60469SMark F. Adams       Determine how many elements are assigned to each processor
26711e60469SMark F. Adams     */
268fd3c6afaSMark F. Adams     inpe = npe;
269fd3c6afaSMark F. Adams     ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr);
27011e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
2719d5b6da9SMark F. Adams     ncrs_new = counts[mype]/cbs;
2729d5b6da9SMark F. Adams     strideNew = ncrs_new*ndata_rows;
273b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
274b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
275b4fbaa2aSMark F. Adams #endif
27622063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
27711e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
278d3d6bff4SMark F. Adams     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
279470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
28011e60469SMark F. Adams     /*
2819d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
2829d5b6da9SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cbs'.
28311e60469SMark F. Adams      */
28492a756f0SMark F. Adams     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
28511e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
286038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
2879d5b6da9SMark F. Adams       PetscInt id = idx[ii*cbs]/cbs; /* get node back */
288d3d6bff4SMark F. Adams       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
28911e60469SMark F. Adams     }
290038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
291d3d6bff4SMark F. Adams     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
29211e60469SMark F. Adams     CHKERRQ(ierr);
29392a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
29411e60469SMark F. Adams     /*
29511e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
29611e60469SMark F. Adams      */
297d3d6bff4SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
2989d5b6da9SMark F. Adams     for( jj=0; jj<ndata_cols ; jj++ ) {
299d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs0 ; ii++) {
3009d5b6da9SMark F. Adams 	for( kk=0; kk<ndata_rows ; kk++ ) {
3019d5b6da9SMark F. Adams 	  PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*ndata_cols + jj;
302c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
303676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
304d3d6bff4SMark F. Adams 	}
305038e3b61SMark F. Adams       }
306eb07cef2SMark F. Adams     }
307eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
308eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
30911e60469SMark F. Adams     /*
31011e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
31111e60469SMark F. Adams       to the correct processor
31211e60469SMark F. Adams     */
31311e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
31411e60469SMark F. Adams     CHKERRQ(ierr);
31511e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
31611e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31711e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31811e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
31911e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
32011e60469SMark F. Adams     /*
32111e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
32211e60469SMark F. Adams     */
323c8b0795cSMark F. Adams     ierr = PetscFree( pc_gamg->data );    CHKERRQ(ierr);
324c8b0795cSMark F. Adams     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), &pc_gamg->data );    CHKERRQ(ierr);
325c8b0795cSMark F. Adams     pc_gamg->data_sz = data_sz*ncrs_new;
32611e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
3279d5b6da9SMark F. Adams     for( jj=0; jj<ndata_cols ; jj++ ) {
328d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs_new ; ii++) {
3299d5b6da9SMark F. Adams 	for( kk=0; kk<ndata_rows ; kk++ ) {
3309d5b6da9SMark F. Adams 	  PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*ndata_cols + jj;
331c8b0795cSMark F. Adams 	  pc_gamg->data[ix] = PetscRealPart(array[jx]);
332d3d6bff4SMark F. Adams 	  array[jx] = 1.e300;
333d3d6bff4SMark F. Adams 	}
334038e3b61SMark F. Adams       }
335038e3b61SMark F. Adams     }
33611e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
33711e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
338ed3f9983SMark F. Adams #if defined PETSC_USE_LOG
339ed3f9983SMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
340ed3f9983SMark F. Adams #endif
34111e60469SMark F. Adams     /*
34211e60469SMark F. Adams       Invert for MatGetSubMatrix
34311e60469SMark F. Adams     */
3449d5b6da9SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*cbs, &new_indices ); CHKERRQ(ierr);
34511e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
34611e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
347ed3f9983SMark F. Adams #if defined PETSC_USE_LOG
348ed3f9983SMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
349ed3f9983SMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
350ed3f9983SMark F. Adams #endif
35111e60469SMark F. Adams     /* A_crs output */
352038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
35311e60469SMark F. Adams     CHKERRQ(ierr);
354eb07cef2SMark F. Adams 
355038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
356e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
3579d5b6da9SMark F. Adams     ierr = MatSetBlockSize( Cmat, cbs );      CHKERRQ(ierr);
358ed3f9983SMark F. Adams #if defined PETSC_USE_LOG
359ed3f9983SMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
360ed3f9983SMark F. Adams #endif
36111e60469SMark F. Adams     /* prolongator */
36211e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
36311e60469SMark F. Adams     {
36411e60469SMark F. Adams       IS findices;
365ed3f9983SMark F. Adams #if defined PETSC_USE_LOG
366ed3f9983SMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
367ed3f9983SMark F. Adams #endif
36811e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
36996434bc1SMark F. Adams #ifdef USE_R
37096434bc1SMark F. Adams       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
37196434bc1SMark F. Adams #else
37211e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
37396434bc1SMark F. Adams #endif
37411e60469SMark F. Adams       CHKERRQ(ierr);
37511e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
376ed3f9983SMark F. Adams #if defined PETSC_USE_LOG
377ed3f9983SMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
378ed3f9983SMark F. Adams #endif
37911e60469SMark F. Adams     }
3803530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
381a147abb0SMark F. Adams     *a_P_inout = Pnew; /* output - repartitioned */
38211e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
38392a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
384e33ef3b1SMark F. Adams   }
3855b89ad90SMark F. Adams 
3865a9b9e01SMark F. Adams   *a_nactive_proc = new_npe; /* output */
3875a9b9e01SMark F. Adams 
388c8b0795cSMark F. Adams   if( !PETSC_TRUE ) {
389c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
390c8b0795cSMark F. Adams     if(llev==0) {
391c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
392c8b0795cSMark F. Adams       PetscViewerASCIIOpen(wcomm,fname,&viewer);
393c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
394c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer ); CHKERRQ(ierr);
395c8b0795cSMark F. Adams       ierr = PetscViewerDestroy( &viewer );
396c8b0795cSMark F. Adams     }
397c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
398c8b0795cSMark F. Adams     PetscViewerASCIIOpen(wcomm,fname,&viewer);
399c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
400c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer ); CHKERRQ(ierr);
401c8b0795cSMark F. Adams     ierr = PetscViewerDestroy( &viewer );
402c8b0795cSMark F. Adams   }
403c8b0795cSMark F. Adams 
4045b89ad90SMark F. Adams   PetscFunctionReturn(0);
4055b89ad90SMark F. Adams }
4065b89ad90SMark F. Adams 
4075b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4085b89ad90SMark F. Adams /*
4095b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4105b89ad90SMark F. Adams                     by setting data structures and options.
4115b89ad90SMark F. Adams 
4125b89ad90SMark F. Adams    Input Parameter:
4135b89ad90SMark F. Adams .  pc - the preconditioner context
4145b89ad90SMark F. Adams 
4155b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4165b89ad90SMark F. Adams 
4175b89ad90SMark F. Adams    Notes:
4185b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4195b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4205b89ad90SMark F. Adams */
4215b89ad90SMark F. Adams #undef __FUNCT__
4225b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4239d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC pc )
4245b89ad90SMark F. Adams {
4255b89ad90SMark F. Adams   PetscErrorCode  ierr;
4269d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
4275b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
4282adcac29SMark F. Adams   Mat              Pmat = pc->pmat;
429d3d6bff4SMark F. Adams   PetscInt         fine_level,level,level1,M,N,bs,nloc,lidx,Istart,Iend;
4309d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
4313530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
432587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
433c8b0795cSMark F. Adams   PetscReal        emaxs[GAMG_MAXLEVELS];
434b2a4f308SMark F. Adams   PetscLogDouble   nnz0,nnztot;
435737a81a9SMark F. Adams   MatInfo          info;
4365ef31b24SMark F. Adams 
4375b89ad90SMark F. Adams   PetscFunctionBegin;
43884d3f75bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
43984d3f75bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
44084d3f75bSMark F. Adams   if( pc_gamg->setup_count++ > 0 ) {
44184d3f75bSMark F. Adams     PC_MG_Levels   **mglevels = mg->levels;
44203a628feSMark F. Adams     /* just do Galerkin grids */
44358471d46SMark F. Adams     Mat B,dA,dB;
444d5280255SMark F. Adams     assert(pc->setupcalled);
44558471d46SMark F. Adams 
4469d5b6da9SMark F. Adams     if( pc_gamg->Nlevels > 1 ) {
44758471d46SMark F. Adams       /* currently only handle case where mat and pmat are the same on coarser levels */
4489d5b6da9SMark F. Adams       ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
44958471d46SMark F. Adams       /* (re)set to get dirty flag */
4509d5b6da9SMark F. Adams       ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
4519d5b6da9SMark F. Adams       ierr = KSPSetUp( mglevels[pc_gamg->Nlevels-1]->smoothd ); CHKERRQ(ierr);
45258471d46SMark F. Adams 
4539d5b6da9SMark F. Adams       for (level=pc_gamg->Nlevels-2; level>-1; level--) {
45403a628feSMark F. Adams         /* the first time through the matrix structure has changed from repartitioning */
45584d3f75bSMark F. Adams         if( pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
45603a628feSMark F. Adams           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
457*084a8fe3SJed Brown           ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
45803a628feSMark F. Adams           mglevels[level]->A = B;
45903a628feSMark F. Adams         }
46003a628feSMark F. Adams         else {
461*084a8fe3SJed Brown           ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
46258471d46SMark F. Adams           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
46303a628feSMark F. Adams         }
46458471d46SMark F. Adams         ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
46558471d46SMark F. Adams         dB = B;
46658471d46SMark F. Adams         /* setup KSP/PC */
46758471d46SMark F. Adams         ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
46858471d46SMark F. Adams       }
4695f8cf99dSMark F. Adams     }
470d5280255SMark F. Adams 
471d5280255SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
472d5280255SMark F. Adams 
473d5280255SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
474d5280255SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
475d5280255SMark F. Adams 
47658471d46SMark F. Adams     PetscFunctionReturn(0);
477eb07cef2SMark F. Adams   }
47884d3f75bSMark F. Adams   assert(pc->setupcalled == 0);
479f6536408SMark F. Adams 
4805b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
4812adcac29SMark F. Adams   ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr);
4822adcac29SMark F. Adams   ierr = MatGetSize( Pmat, &M, &N );CHKERRQ(ierr);
4832adcac29SMark F. Adams   ierr = MatGetOwnershipRange( Pmat, &Istart, &Iend ); CHKERRQ(ierr);
484eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
485eb07cef2SMark F. Adams 
4869d5b6da9SMark F. Adams   if( pc_gamg->data == 0 && nloc > 0 ) {
4879d5b6da9SMark F. Adams     if(!pc_gamg->createdefaultdata){
488c8b0795cSMark F. Adams       SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"'createdefaultdata' not set?!?! need to support NULL data!!!");
489eb07cef2SMark F. Adams     }
4909d5b6da9SMark F. Adams     ierr = pc_gamg->createdefaultdata( pc ); CHKERRQ(ierr);
4919d5b6da9SMark F. Adams   }
492038e3b61SMark F. Adams 
493eb07cef2SMark F. Adams   /* Get A_i and R_i */
494c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
4952adcac29SMark F. Adams     if(pc_gamg->verbose==1) ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);
4962adcac29SMark F. Adams     else ierr =  MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);
497b2a4f308SMark F. Adams     CHKERRQ(ierr);
498b2a4f308SMark F. Adams     nnz0 = info.nz_used;
499b2a4f308SMark F. Adams     nnztot = info.nz_used;
500b43b03e9SMark 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",
501c8b0795cSMark F. Adams                 mype,__FUNCT__,0,N,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
502b2a4f308SMark F. Adams                 (int)(nnz0/(PetscReal)N),npe);
503c8b0795cSMark F. Adams   }
50484d3f75bSMark F. Adams 
5058f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
5069d5b6da9SMark F. Adams         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
5070205a208SMark F. Adams         level++ ){
5085b89ad90SMark F. Adams     level1 = level + 1;
509b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
510b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
511b4fbaa2aSMark F. Adams #endif
512b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
513b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr);
514b4fbaa2aSMark F. Adams #endif
515c8b0795cSMark F. Adams     { /* construct prolongator */
516725b86d8SJed Brown       Mat Gmat;
517725b86d8SJed Brown       IS selected, llist;
518725b86d8SJed Brown       assert(pc_gamg->graph);
519725b86d8SJed Brown       assert(pc_gamg->coarsen);
520c8b0795cSMark F. Adams 
521b43b03e9SMark F. Adams       ierr = pc_gamg->graph( pc, Aarr[level], &Gmat ); CHKERRQ(ierr);
522b43b03e9SMark F. Adams       ierr = pc_gamg->coarsen( pc, Gmat, &selected, &llist ); CHKERRQ(ierr);
523c8b0795cSMark F. Adams       ierr = pc_gamg->prolongator( pc, Aarr[level], Gmat, selected, llist, &Parr[level1] );
5245b89ad90SMark F. Adams       CHKERRQ(ierr);
525c8b0795cSMark F. Adams 
526c8b0795cSMark F. Adams       if( Parr[level1] ){
5279d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
528c8b0795cSMark F. Adams         if( pc_gamg->col_bs_id != -1 ){
5299d5b6da9SMark F. Adams           PetscBool flag;
5309d5b6da9SMark F. Adams           ierr = PetscObjectComposedDataGetInt( (PetscObject)Parr[level1], pc_gamg->col_bs_id, bs, flag );
531c8b0795cSMark F. Adams           CHKERRQ( ierr ); assert(flag);
5329d5b6da9SMark F. Adams         }
533c8b0795cSMark F. Adams       }
534c8b0795cSMark F. Adams 
535c8b0795cSMark F. Adams       if( pc_gamg->optprol && Parr[level1] ){
536c8b0795cSMark F. Adams         /* smooth */
537c8b0795cSMark F. Adams         ierr = pc_gamg->optprol( pc, Aarr[level], &Parr[level1] ); CHKERRQ(ierr);
538c8b0795cSMark F. Adams       }
539c8b0795cSMark F. Adams 
540c8b0795cSMark F. Adams       ierr = MatDestroy( &Gmat );  CHKERRQ(ierr);
541c8b0795cSMark F. Adams       ierr = ISDestroy( &selected );  CHKERRQ(ierr);
542c8b0795cSMark F. Adams       ierr = ISDestroy( &llist );  CHKERRQ(ierr);
543c8b0795cSMark F. Adams     }
544c8b0795cSMark F. Adams #if defined PETSC_USE_LOG
545c8b0795cSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
546c8b0795cSMark F. Adams #endif
5479d5b6da9SMark F. Adams     /* cache eigen estimate */
5489d5b6da9SMark F. Adams     if( pc_gamg->emax_id != -1 ){
5499d5b6da9SMark F. Adams       PetscBool flag;
550f73f8d2cSSatish Balay       ierr = PetscObjectComposedDataGetReal( (PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag );
5519d5b6da9SMark F. Adams       CHKERRQ( ierr );
5529d5b6da9SMark F. Adams       if( !flag ) emaxs[level] = -1.;
5539d5b6da9SMark F. Adams     }
5549d5b6da9SMark F. Adams     else emaxs[level] = -1.;
5552adcac29SMark F. Adams     if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
556c8b0795cSMark F. Adams     if( !Parr[level1] ) {
557b43b03e9SMark F. Adams       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);
558c8b0795cSMark F. Adams       break;
559c8b0795cSMark F. Adams     }
560b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
561b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
562b4fbaa2aSMark F. Adams #endif
563c8b0795cSMark F. Adams     ierr = createLevel( pc, Aarr[level], bs,
5649d5b6da9SMark F. Adams                         (PetscBool)(level==pc_gamg->Nlevels-2),
565c8b0795cSMark F. Adams                         &Parr[level1], &nactivepe, &Aarr[level1] );
5663530afc2SMark F. Adams     CHKERRQ(ierr);
567b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
568b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
569b4fbaa2aSMark F. Adams #endif
5703530afc2SMark F. Adams     ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
571c8b0795cSMark F. Adams     if (pc_gamg->verbose){
57284d3f75bSMark F. Adams       ierr = MatGetInfo( Aarr[level1], MAT_GLOBAL_SUM, &info ); CHKERRQ(ierr);
573b43b03e9SMark F. Adams       PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
574c8b0795cSMark F. Adams                   mype,__FUNCT__,(int)level1,N,pc_gamg->data_cell_cols,
5752c047e2dSMark F. Adams                   (int)(info.nz_used/(PetscReal)N),nactivepe);
576c8b0795cSMark F. Adams     }
57781708dccSMark F. Adams     /* stop if one node */
578c8b0795cSMark F. Adams     if( M/pc_gamg->data_cell_cols < 2 ) {
57981708dccSMark F. Adams       level++;
58081708dccSMark F. Adams       break;
58181708dccSMark F. Adams     }
582c8b0795cSMark F. Adams     /* Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; */
583c8b0795cSMark F. Adams     /* v = 1.e-10; /\* LU factor has hard wired numbers for small diags so this needs to match (yuk) *\/ */
584c8b0795cSMark F. Adams     /* ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); */
585c8b0795cSMark F. Adams     /* nloceq = Iend-Istart; */
586c8b0795cSMark F. Adams     /* ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr); */
587c8b0795cSMark F. Adams     /* ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr); */
588c8b0795cSMark F. Adams     /* ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr); */
589c8b0795cSMark F. Adams     /* for(kk=0;kk<nloceq;kk++){ */
590c8b0795cSMark F. Adams     /*   if(data_arr[kk]==0.0) { */
591c8b0795cSMark F. Adams     /*     id = kk + Istart; */
592c8b0795cSMark F. Adams     /*     ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); */
593c8b0795cSMark F. Adams     /*     CHKERRQ(ierr); */
594c8b0795cSMark F. Adams     /*     PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1); */
595c8b0795cSMark F. Adams     /*   } */
596c8b0795cSMark F. Adams     /* } */
597c8b0795cSMark F. Adams     /* ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); */
598c8b0795cSMark F. Adams     /* ierr = VecDestroy( &diag );                CHKERRQ(ierr); */
59984d3f75bSMark F. Adams     /* ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
60084d3f75bSMark F. Adams     /* ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
601b4fbaa2aSMark F. Adams 
602b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
603b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
604b4fbaa2aSMark F. Adams #endif
605b2a4f308SMark F. Adams     if(pc_gamg->verbose){
606b2a4f308SMark F. Adams       if(pc_gamg->verbose==1) ierr =  MatGetInfo(Aarr[level1],MAT_LOCAL,&info);
607b2a4f308SMark F. Adams       else ierr =  MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info);
608b2a4f308SMark F. Adams       CHKERRQ(ierr);
609b2a4f308SMark F. Adams       nnztot += info.nz_used;
610b2a4f308SMark F. Adams     }
611c8b0795cSMark F. Adams   } /* levels */
612c8b0795cSMark F. Adams 
613c8b0795cSMark F. Adams   if( pc_gamg->data ) {
614c8b0795cSMark F. Adams     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
615c8b0795cSMark F. Adams     pc_gamg->data = 0;
6165b89ad90SMark F. Adams   }
617c8b0795cSMark F. Adams 
618b2a4f308SMark F. Adams   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid compexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
6199d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6205b89ad90SMark F. Adams   fine_level = level;
6219d5b6da9SMark F. Adams   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
6225b89ad90SMark F. Adams 
62384d3f75bSMark F. Adams   /* simple setup */
62484d3f75bSMark F. Adams   if( !PETSC_TRUE ){
62584d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
62684d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
62784d3f75bSMark F. Adams          lidx<fine_level;
62884d3f75bSMark F. Adams          lidx++, level--){
62984d3f75bSMark F. Adams       ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr);
63084d3f75bSMark F. Adams       ierr = KSPSetOperators( mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
63184d3f75bSMark F. Adams       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
63284d3f75bSMark F. Adams       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
63384d3f75bSMark F. Adams     }
63484d3f75bSMark F. Adams     ierr = KSPSetOperators( mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
63584d3f75bSMark F. Adams 
63684d3f75bSMark F. Adams     ierr = PCSetUp_MG( pc );  CHKERRQ( ierr );
63784d3f75bSMark F. Adams   }
63884d3f75bSMark F. Adams   else if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if one level */
639d5280255SMark F. Adams     /* set default smoothers & set operators */
6409d5b6da9SMark F. Adams     for ( lidx = 1, level = pc_gamg->Nlevels-2;
641587fa25dSMark F. Adams           lidx <= fine_level;
642587fa25dSMark F. Adams           lidx++, level--) {
6435b89ad90SMark F. Adams       KSP smoother; PC subpc;
644f6536408SMark F. Adams       /* set defaults */
6459d5b6da9SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
6465b89ad90SMark F. Adams       ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
647f6536408SMark F. Adams       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
6489d5b6da9SMark F. Adams       ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr);
649f6536408SMark F. Adams       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
650d5280255SMark F. Adams       /* set ops */
651f6536408SMark F. Adams       ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
652d5280255SMark F. Adams       ierr = PCMGSetInterpolation( pc, lidx, Parr[level+1] );CHKERRQ(ierr);
653d5280255SMark F. Adams     }
654d5280255SMark F. Adams     {
655d5280255SMark F. Adams       /* coarse grid */
656d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
657d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
658d5280255SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
659d5280255SMark F. Adams       ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
660d5280255SMark F. Adams       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
661d5280255SMark F. Adams       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
662d5280255SMark F. Adams       ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
663d5280255SMark F. Adams       ierr = PCSetUp( subpc ); CHKERRQ(ierr);
664d5280255SMark F. Adams       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
665d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
666d5280255SMark F. Adams       ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
667d5280255SMark F. Adams       /* set ops */
668d5280255SMark F. Adams       ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
669d5280255SMark F. Adams     }
670d5280255SMark F. Adams 
671d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
672d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
673d5280255SMark F. Adams     ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr);
674d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
675d5280255SMark 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.");
676d5280255SMark F. Adams 
677d5280255SMark F. Adams     /* create cheby smoothers */
678d5280255SMark F. Adams     for ( lidx = 1, level = pc_gamg->Nlevels-2;
679d5280255SMark F. Adams           lidx <= fine_level;
680d5280255SMark F. Adams           lidx++, level--) {
681d5280255SMark F. Adams       KSP smoother;
682d5280255SMark F. Adams       PetscBool isCheb;
683d5280255SMark F. Adams 
684d5280255SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
685f6536408SMark F. Adams       /* do my own cheby */
686f6536408SMark F. Adams       ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr);
687f6536408SMark F. Adams       if( isCheb ) {
688d5280255SMark F. Adams         PetscReal emax, emin;
689d5280255SMark F. Adams         PC subpc;
690d5280255SMark F. Adams 
691d5280255SMark F. Adams         ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
6929d5b6da9SMark F. Adams         ierr = PetscTypeCompare( (PetscObject)subpc, PCJACOBI, &isCheb ); CHKERRQ(ierr);
693f6536408SMark F. Adams         if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
694587fa25dSMark F. Adams         else{ /* eigen estimate 'emax' */
695587fa25dSMark F. Adams           KSP eksp; Mat Lmat = Aarr[level];
696b2a4f308SMark F. Adams           Vec bb, xx;
697038e3b61SMark F. Adams 
6985745e0f5SMark F. Adams           ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
6995745e0f5SMark F. Adams           ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
700fc4362bfSMark F. Adams           {
701fc4362bfSMark F. Adams             PetscRandom    rctx;
702fc4362bfSMark F. Adams             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
703fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
704fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
705fc4362bfSMark F. Adams             ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
7065b89ad90SMark F. Adams           }
707fc4362bfSMark F. Adams           ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
708038e3b61SMark F. Adams           ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
709fc4362bfSMark F. Adams           CHKERRQ(ierr);
710fc4362bfSMark F. Adams           ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
711fc4362bfSMark F. Adams 
712f6536408SMark F. Adams           ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
713f6536408SMark F. Adams           ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
714f6536408SMark F. Adams 
715f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
716f6536408SMark F. Adams           ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
717fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
7185a9b9e01SMark F. Adams 
719b2a4f308SMark F. Adams           { /* set PC type to be same as smoother - does not get all parameters!!! */
720b2a4f308SMark F. Adams             const PCType type;    PC pc;
721b2a4f308SMark F. Adams             ierr = PCGetType( subpc, &type );   CHKERRQ(ierr);
722b2a4f308SMark F. Adams             ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
723b2a4f308SMark F. Adams             ierr = PCSetType( pc, type ); CHKERRQ(ierr);
724b2a4f308SMark F. Adams           }
725b2a4f308SMark F. Adams 
7265a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
7275a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
7285a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
729fc4362bfSMark F. Adams           ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
7305a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
7315a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
7325a9b9e01SMark F. Adams 
733fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
7345a9b9e01SMark F. Adams 
735fc4362bfSMark F. Adams           ierr = VecDestroy( &xx );       CHKERRQ(ierr);
736fc4362bfSMark F. Adams           ierr = VecDestroy( &bb );       CHKERRQ(ierr);
737fc4362bfSMark F. Adams           ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
738f6536408SMark F. Adams 
7399d5b6da9SMark F. Adams           if (pc_gamg->verbose) {
740b2a4f308SMark F. Adams             PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e\n",__FUNCT__,emax,emin);
741f6536408SMark F. Adams           }
742fc4362bfSMark F. Adams         }
743038e3b61SMark F. Adams         {
744038e3b61SMark F. Adams           PetscInt N1, N0, tt;
745038e3b61SMark F. Adams           ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
746038e3b61SMark F. Adams           ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
747f6536408SMark F. Adams           /* heuristic - is this crap? */
748f6536408SMark F. Adams           emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
749038e3b61SMark F. Adams           emax *= 1.05;
750038e3b61SMark F. Adams         }
751fc4362bfSMark F. Adams         ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
752f6536408SMark F. Adams       }
7535745e0f5SMark F. Adams     }
754737a81a9SMark F. Adams 
755d5280255SMark F. Adams     /* clean up */
756d5280255SMark F. Adams     for(level=1;level<pc_gamg->Nlevels;level++){
757587fa25dSMark F. Adams       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
758587fa25dSMark F. Adams       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
7595b89ad90SMark F. Adams     }
7605745e0f5SMark F. Adams 
7619d5b6da9SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
76258471d46SMark F. Adams 
76358471d46SMark F. Adams     {
76458471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
7659d5b6da9SMark F. Adams       ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
76658471d46SMark F. Adams       ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
767d5280255SMark F. Adams       //ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
76858471d46SMark F. Adams     }
7695f8cf99dSMark F. Adams   }
7705f8cf99dSMark F. Adams   else {
7715f8cf99dSMark F. Adams     KSP smoother;
772b43b03e9SMark 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__);
7739d5b6da9SMark F. Adams     ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
7745f8cf99dSMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
7755f8cf99dSMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
7769d5b6da9SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
7775f8cf99dSMark F. Adams   }
77884d3f75bSMark F. Adams 
7795b89ad90SMark F. Adams   PetscFunctionReturn(0);
7805b89ad90SMark F. Adams }
7815b89ad90SMark F. Adams 
782eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7835b89ad90SMark F. Adams /*
7845b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7855b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7865b89ad90SMark F. Adams 
7875b89ad90SMark F. Adams    Input Parameter:
7885b89ad90SMark F. Adams .  pc - the preconditioner context
7895b89ad90SMark F. Adams 
7905b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7915b89ad90SMark F. Adams */
7925b89ad90SMark F. Adams #undef __FUNCT__
7935b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
7945b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG( PC pc )
7955b89ad90SMark F. Adams {
7965b89ad90SMark F. Adams   PetscErrorCode  ierr;
7975b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
7985b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
7995b89ad90SMark F. Adams 
8005b89ad90SMark F. Adams   PetscFunctionBegin;
8015b89ad90SMark F. Adams   ierr = PCReset_GAMG( pc );CHKERRQ(ierr);
8025b89ad90SMark F. Adams   ierr = PetscFree( pc_gamg );CHKERRQ(ierr);
8035b89ad90SMark F. Adams   ierr = PCDestroy_MG( pc );CHKERRQ(ierr);
8045b89ad90SMark F. Adams   PetscFunctionReturn(0);
8055b89ad90SMark F. Adams }
8065b89ad90SMark F. Adams 
807676e1743SMark F. Adams 
808676e1743SMark F. Adams #undef __FUNCT__
809676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
810676e1743SMark F. Adams /*@
811676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
812676e1743SMark F. Adams    processor reduction.
813676e1743SMark F. Adams 
814676e1743SMark F. Adams    Not Collective on PC
815676e1743SMark F. Adams 
816676e1743SMark F. Adams    Input Parameters:
817676e1743SMark F. Adams .  pc - the preconditioner context
818676e1743SMark F. Adams 
819676e1743SMark F. Adams 
820676e1743SMark F. Adams    Options Database Key:
821676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
822676e1743SMark F. Adams 
823676e1743SMark F. Adams    Level: intermediate
824676e1743SMark F. Adams 
825676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
826676e1743SMark F. Adams 
827676e1743SMark F. Adams .seealso: ()
828676e1743SMark F. Adams @*/
829676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
830676e1743SMark F. Adams {
831676e1743SMark F. Adams   PetscErrorCode ierr;
832676e1743SMark F. Adams 
833676e1743SMark F. Adams   PetscFunctionBegin;
834676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
835676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
836676e1743SMark F. Adams   PetscFunctionReturn(0);
837676e1743SMark F. Adams }
838676e1743SMark F. Adams 
839676e1743SMark F. Adams EXTERN_C_BEGIN
840676e1743SMark F. Adams #undef __FUNCT__
841676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
842676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
843676e1743SMark F. Adams {
844c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
845c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
846676e1743SMark F. Adams 
847676e1743SMark F. Adams   PetscFunctionBegin;
8489d5b6da9SMark F. Adams   if(n>0) pc_gamg->min_eq_proc = n;
849676e1743SMark F. Adams   PetscFunctionReturn(0);
850676e1743SMark F. Adams }
851676e1743SMark F. Adams EXTERN_C_END
852676e1743SMark F. Adams 
853676e1743SMark F. Adams #undef __FUNCT__
854389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
855389730f3SMark F. Adams /*@
856389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
857389730f3SMark F. Adams 
858389730f3SMark F. Adams  Collective on PC
859389730f3SMark F. Adams 
860389730f3SMark F. Adams    Input Parameters:
861389730f3SMark F. Adams .  pc - the preconditioner context
862389730f3SMark F. Adams 
863389730f3SMark F. Adams 
864389730f3SMark F. Adams    Options Database Key:
865389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
866389730f3SMark F. Adams 
867389730f3SMark F. Adams    Level: intermediate
868389730f3SMark F. Adams 
869389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
870389730f3SMark F. Adams 
871389730f3SMark F. Adams .seealso: ()
872389730f3SMark F. Adams  @*/
873389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
874389730f3SMark F. Adams {
875389730f3SMark F. Adams   PetscErrorCode ierr;
876389730f3SMark F. Adams 
877389730f3SMark F. Adams   PetscFunctionBegin;
878389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
879389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
880389730f3SMark F. Adams   PetscFunctionReturn(0);
881389730f3SMark F. Adams }
882389730f3SMark F. Adams 
883389730f3SMark F. Adams EXTERN_C_BEGIN
884389730f3SMark F. Adams #undef __FUNCT__
885389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
886389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
887389730f3SMark F. Adams {
888389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
889389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
890389730f3SMark F. Adams 
891389730f3SMark F. Adams   PetscFunctionBegin;
8929d5b6da9SMark F. Adams   if(n>0) pc_gamg->coarse_eq_limit = n;
893389730f3SMark F. Adams   PetscFunctionReturn(0);
894389730f3SMark F. Adams }
895389730f3SMark F. Adams EXTERN_C_END
896389730f3SMark F. Adams 
897389730f3SMark F. Adams #undef __FUNCT__
8988263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
899676e1743SMark F. Adams /*@
9008263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
901676e1743SMark F. Adams 
902676e1743SMark F. Adams    Collective on PC
903676e1743SMark F. Adams 
904676e1743SMark F. Adams    Input Parameters:
905676e1743SMark F. Adams .  pc - the preconditioner context
906676e1743SMark F. Adams 
907676e1743SMark F. Adams 
908676e1743SMark F. Adams    Options Database Key:
9098263b398SMark F. Adams .  -pc_gamg_repartition
910676e1743SMark F. Adams 
911676e1743SMark F. Adams    Level: intermediate
912676e1743SMark F. Adams 
913676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
914676e1743SMark F. Adams 
915676e1743SMark F. Adams .seealso: ()
916676e1743SMark F. Adams @*/
9178263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
918676e1743SMark F. Adams {
919676e1743SMark F. Adams   PetscErrorCode ierr;
920676e1743SMark F. Adams 
921676e1743SMark F. Adams   PetscFunctionBegin;
922676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9238263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
924676e1743SMark F. Adams   PetscFunctionReturn(0);
925676e1743SMark F. Adams }
926676e1743SMark F. Adams 
927676e1743SMark F. Adams EXTERN_C_BEGIN
928676e1743SMark F. Adams #undef __FUNCT__
9298263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
9308263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
931676e1743SMark F. Adams {
932c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
933c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
934676e1743SMark F. Adams 
935676e1743SMark F. Adams   PetscFunctionBegin;
9369d5b6da9SMark F. Adams   pc_gamg->repart = n;
937676e1743SMark F. Adams   PetscFunctionReturn(0);
938676e1743SMark F. Adams }
939676e1743SMark F. Adams EXTERN_C_END
940676e1743SMark F. Adams 
941676e1743SMark F. Adams #undef __FUNCT__
9424ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
9434ef23d27SMark F. Adams /*@
9444ef23d27SMark F. Adams    PCGAMGSetNlevels -
9454ef23d27SMark F. Adams 
9464ef23d27SMark F. Adams    Not collective on PC
9474ef23d27SMark F. Adams 
9484ef23d27SMark F. Adams    Input Parameters:
9494ef23d27SMark F. Adams .  pc - the preconditioner context
9504ef23d27SMark F. Adams 
9514ef23d27SMark F. Adams 
9524ef23d27SMark F. Adams    Options Database Key:
9534ef23d27SMark F. Adams .  -pc_mg_levels
9544ef23d27SMark F. Adams 
9554ef23d27SMark F. Adams    Level: intermediate
9564ef23d27SMark F. Adams 
9574ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
9584ef23d27SMark F. Adams 
9594ef23d27SMark F. Adams .seealso: ()
9604ef23d27SMark F. Adams @*/
9614ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
9624ef23d27SMark F. Adams {
9634ef23d27SMark F. Adams   PetscErrorCode ierr;
9644ef23d27SMark F. Adams 
9654ef23d27SMark F. Adams   PetscFunctionBegin;
9664ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9674ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
9684ef23d27SMark F. Adams   PetscFunctionReturn(0);
9694ef23d27SMark F. Adams }
9704ef23d27SMark F. Adams 
9714ef23d27SMark F. Adams EXTERN_C_BEGIN
9724ef23d27SMark F. Adams #undef __FUNCT__
9734ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
9744ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
9754ef23d27SMark F. Adams {
9764ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
9774ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
9784ef23d27SMark F. Adams 
9794ef23d27SMark F. Adams   PetscFunctionBegin;
9809d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
9814ef23d27SMark F. Adams   PetscFunctionReturn(0);
9824ef23d27SMark F. Adams }
9834ef23d27SMark F. Adams EXTERN_C_END
9844ef23d27SMark F. Adams 
9854ef23d27SMark F. Adams #undef __FUNCT__
9863542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
9873542efc5SMark F. Adams /*@
9883542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
9893542efc5SMark F. Adams 
9903542efc5SMark F. Adams    Not collective on PC
9913542efc5SMark F. Adams 
9923542efc5SMark F. Adams    Input Parameters:
9933542efc5SMark F. Adams .  pc - the preconditioner context
9943542efc5SMark F. Adams 
9953542efc5SMark F. Adams 
9963542efc5SMark F. Adams    Options Database Key:
9973542efc5SMark F. Adams .  -pc_gamg_threshold
9983542efc5SMark F. Adams 
9993542efc5SMark F. Adams    Level: intermediate
10003542efc5SMark F. Adams 
10013542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
10023542efc5SMark F. Adams 
10033542efc5SMark F. Adams .seealso: ()
10043542efc5SMark F. Adams @*/
10053542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
10063542efc5SMark F. Adams {
10073542efc5SMark F. Adams   PetscErrorCode ierr;
10083542efc5SMark F. Adams 
10093542efc5SMark F. Adams   PetscFunctionBegin;
10103542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10113542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
10123542efc5SMark F. Adams   PetscFunctionReturn(0);
10133542efc5SMark F. Adams }
10143542efc5SMark F. Adams 
10153542efc5SMark F. Adams EXTERN_C_BEGIN
10163542efc5SMark F. Adams #undef __FUNCT__
10173542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
10183542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
10193542efc5SMark F. Adams {
1020c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1021c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
10223542efc5SMark F. Adams 
10233542efc5SMark F. Adams   PetscFunctionBegin;
10249d5b6da9SMark F. Adams   pc_gamg->threshold = n;
10253542efc5SMark F. Adams   PetscFunctionReturn(0);
10263542efc5SMark F. Adams }
10273542efc5SMark F. Adams EXTERN_C_END
10283542efc5SMark F. Adams 
10293542efc5SMark F. Adams #undef __FUNCT__
10309d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1031676e1743SMark F. Adams /*@
10329d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1033676e1743SMark F. Adams 
1034676e1743SMark F. Adams    Collective on PC
1035676e1743SMark F. Adams 
1036676e1743SMark F. Adams    Input Parameters:
1037676e1743SMark F. Adams .  pc - the preconditioner context
1038676e1743SMark F. Adams 
1039676e1743SMark F. Adams 
1040676e1743SMark F. Adams    Options Database Key:
10413542efc5SMark F. Adams .  -pc_gamg_type
1042676e1743SMark F. Adams 
1043676e1743SMark F. Adams    Level: intermediate
1044676e1743SMark F. Adams 
1045676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1046676e1743SMark F. Adams 
1047676e1743SMark F. Adams .seealso: ()
1048676e1743SMark F. Adams @*/
10490202ef74SSatish Balay PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type )
1050676e1743SMark F. Adams {
1051676e1743SMark F. Adams   PetscErrorCode ierr;
1052676e1743SMark F. Adams 
1053676e1743SMark F. Adams   PetscFunctionBegin;
1054676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10550202ef74SSatish Balay   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type));
1056676e1743SMark F. Adams   CHKERRQ(ierr);
1057676e1743SMark F. Adams   PetscFunctionReturn(0);
1058676e1743SMark F. Adams }
1059676e1743SMark F. Adams 
1060676e1743SMark F. Adams EXTERN_C_BEGIN
1061676e1743SMark F. Adams #undef __FUNCT__
10629d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
10630202ef74SSatish Balay PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type )
1064676e1743SMark F. Adams {
10659d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
1066676e1743SMark F. Adams 
1067676e1743SMark F. Adams   PetscFunctionBegin;
10689d5b6da9SMark F. Adams   ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r);
10699d5b6da9SMark F. Adams   CHKERRQ(ierr);
10709d5b6da9SMark F. Adams 
10719d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
10729d5b6da9SMark F. Adams 
10739d5b6da9SMark F. Adams   /* call sub create method */
10749d5b6da9SMark F. Adams   ierr = (*r)(pc); CHKERRQ(ierr);
10759d5b6da9SMark F. Adams 
1076676e1743SMark F. Adams   PetscFunctionReturn(0);
1077676e1743SMark F. Adams }
1078676e1743SMark F. Adams EXTERN_C_END
1079676e1743SMark F. Adams 
10805b89ad90SMark F. Adams #undef __FUNCT__
10815b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
10825b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG( PC pc )
10835b89ad90SMark F. Adams {
1084676e1743SMark F. Adams   PetscErrorCode  ierr;
1085676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1086676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1087676e1743SMark F. Adams   PetscBool        flag;
1088b43b03e9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
10895b89ad90SMark F. Adams 
10905b89ad90SMark F. Adams   PetscFunctionBegin;
1091676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1092676e1743SMark F. Adams   {
109375b74bdaSMark F. Adams     /* -pc_gamg_verbose */
10949d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
10959d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
10969d5b6da9SMark F. Adams                            &pc_gamg->verbose, PETSC_NULL );
10979d5b6da9SMark F. Adams     CHKERRQ(ierr);
109875b74bdaSMark F. Adams 
10998263b398SMark F. Adams     /* -pc_gamg_repartition */
11008263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
11018263b398SMark F. Adams                             "Repartion coarse grids (false)",
11028263b398SMark F. Adams                             "PCGAMGRepartitioning",
11039d5b6da9SMark F. Adams                             pc_gamg->repart,
11049d5b6da9SMark F. Adams                             &pc_gamg->repart,
1105676e1743SMark F. Adams                             &flag);
1106676e1743SMark F. Adams     CHKERRQ(ierr);
1107676e1743SMark F. Adams 
1108c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1109676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1110676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1111676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
11129d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
11139d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1114676e1743SMark F. Adams                            &flag );
1115676e1743SMark F. Adams     CHKERRQ(ierr);
11163542efc5SMark F. Adams 
1117389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1118389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1119389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1120389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
11219d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
11229d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1123389730f3SMark F. Adams                            &flag );
1124389730f3SMark F. Adams     CHKERRQ(ierr);
1125389730f3SMark F. Adams 
1126c20e4228SMark F. Adams     /* -pc_gamg_threshold */
11273542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
11283542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
11293542efc5SMark F. Adams                             "PCGAMGSetThreshold",
11309d5b6da9SMark F. Adams                             pc_gamg->threshold,
11319d5b6da9SMark F. Adams                             &pc_gamg->threshold,
11323542efc5SMark F. Adams                             &flag );
11333542efc5SMark F. Adams     CHKERRQ(ierr);
1134b43b03e9SMark F. Adams     if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);
11354ef23d27SMark F. Adams 
11364ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
11374ef23d27SMark F. Adams                            "Set number of MG levels",
11384ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
11399d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
11409d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
11414ef23d27SMark F. Adams                            &flag );
1142676e1743SMark F. Adams   }
1143676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1144676e1743SMark F. Adams 
11455b89ad90SMark F. Adams   PetscFunctionReturn(0);
11465b89ad90SMark F. Adams }
11475b89ad90SMark F. Adams 
11485b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
11495b89ad90SMark F. Adams /*MC
11509d5b6da9SMark F. Adams      PCCreate_GAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
11519d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
11525b89ad90SMark F. Adams 
1153280d9858SJed Brown    Options Database Keys:
11545b89ad90SMark F. Adams    Multigrid options(inherited)
1155280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1156280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1157280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1158280d9858SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
11595b89ad90SMark F. Adams 
11605b89ad90SMark F. Adams   Level: intermediate
1161280d9858SJed Brown 
11625b89ad90SMark F. Adams   Concepts: multigrid
11635b89ad90SMark F. Adams 
11645b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1165280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
11665b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
11675b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
11685b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11695b89ad90SMark F. Adams M*/
11705b89ad90SMark F. Adams EXTERN_C_BEGIN
11715b89ad90SMark F. Adams #undef __FUNCT__
11725b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
11735b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG( PC pc )
11745b89ad90SMark F. Adams {
11755b89ad90SMark F. Adams   PetscErrorCode  ierr;
11765b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
11775b89ad90SMark F. Adams   PC_MG           *mg;
11785ef31b24SMark F. Adams   PetscClassId     cookie;
11799d5b6da9SMark F. Adams 
11805ee9c036SSatish Balay #if defined PETSC_USE_LOG
11819d5b6da9SMark F. Adams   static long count = 0;
11825ee9c036SSatish Balay #endif
11835b89ad90SMark F. Adams 
11845b89ad90SMark F. Adams   PetscFunctionBegin;
11855b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
11865b89ad90SMark F. Adams   ierr = PCSetType( pc, PCMG );  CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
11875b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr);
11885b89ad90SMark F. Adams 
11895b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
11905b89ad90SMark F. Adams   ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr);
11915b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
1192ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
11935b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
11945b89ad90SMark F. Adams 
11959d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
11969d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
11979d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
11989d5b6da9SMark F. Adams   pc_gamg->data = 0;
11995b89ad90SMark F. Adams 
12009d5b6da9SMark F. Adams   /* register AMG type */
12019d5b6da9SMark F. Adams   if( !GAMGList ){
12029d5b6da9SMark F. Adams     ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
12039d5b6da9SMark F. Adams     ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
12049d5b6da9SMark F. Adams   }
12059d5b6da9SMark F. Adams 
12069d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
12075b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
12085b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
12095b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
12105b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
12115b89ad90SMark F. Adams 
12125b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1213676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1214676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1215676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1216676e1743SMark F. Adams   CHKERRQ(ierr);
1217676e1743SMark F. Adams 
1218676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1219389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_C",
1220389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_GAMG",
1221389730f3SMark F. Adams 					    PCGAMGSetCoarseEqLim_GAMG);
1222389730f3SMark F. Adams   CHKERRQ(ierr);
1223389730f3SMark F. Adams 
1224389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12258263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_C",
12268263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_GAMG",
12278263b398SMark F. Adams 					    PCGAMGSetRepartitioning_GAMG);
1228676e1743SMark F. Adams   CHKERRQ(ierr);
1229676e1743SMark F. Adams 
1230676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12313542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
12323542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
12333542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
12343542efc5SMark F. Adams   CHKERRQ(ierr);
12353542efc5SMark F. Adams 
12363542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12379d5b6da9SMark F. Adams 					    "PCGAMGSetType_C",
12389d5b6da9SMark F. Adams 					    "PCGAMGSetType_GAMG",
12399d5b6da9SMark F. Adams 					    PCGAMGSetType_GAMG);
1240676e1743SMark F. Adams   CHKERRQ(ierr);
1241c97e1df0SMark F. Adams 
12429d5b6da9SMark F. Adams   pc_gamg->repart = PETSC_FALSE;
12439d5b6da9SMark F. Adams   pc_gamg->min_eq_proc = 100;
1244c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit = 800;
12459d5b6da9SMark F. Adams   pc_gamg->threshold = 0.05;
12469d5b6da9SMark F. Adams   pc_gamg->Nlevels = GAMG_MAXLEVELS;
12479d5b6da9SMark F. Adams   pc_gamg->verbose = 0;
12489d5b6da9SMark F. Adams   pc_gamg->emax_id = -1;
12499d5b6da9SMark F. Adams   pc_gamg->col_bs_id = -1;
12509d5b6da9SMark F. Adams 
12519d5b6da9SMark F. Adams   /* instantiate derived type */
12529d5b6da9SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
12539d5b6da9SMark F. Adams   {
12549d5b6da9SMark F. Adams     char tname[256] = GAMGAGG;
12559d5b6da9SMark F. Adams     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",
12569d5b6da9SMark F. Adams                             GAMGList, tname, tname, sizeof(tname), PETSC_NULL );
12579d5b6da9SMark F. Adams     CHKERRQ(ierr);
12589d5b6da9SMark F. Adams     ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr);
12599d5b6da9SMark F. Adams   }
12609d5b6da9SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
12619d5b6da9SMark F. Adams 
1262b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
1263785cba28SMark F. Adams   if( count++ == 0 ) {
12645ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
1265b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
12662a44bfbaSMark F. Adams     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
1267c8b0795cSMark F. Adams     /* PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]); */
1268c8b0795cSMark F. Adams     /* PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]); */
1269c8b0795cSMark F. Adams     /* PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]); */
1270b4fbaa2aSMark F. Adams     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1271b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1272b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1273b4fbaa2aSMark F. Adams     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1274c8b0795cSMark F. Adams     PetscLogEventRegister("  SA: colect data", cookie, &gamg_setup_events[SET7]);
1275c8b0795cSMark F. Adams     PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]);
1276b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1277b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1278ed3f9983SMark F. Adams     PetscLogEventRegister("  repartition", cookie, &gamg_setup_events[SET12]);
1279ed3f9983SMark F. Adams     PetscLogEventRegister("  Invert-Sort", cookie, &gamg_setup_events[SET13]);
1280ed3f9983SMark F. Adams     PetscLogEventRegister("  Move A", cookie, &gamg_setup_events[SET14]);
1281ed3f9983SMark F. Adams     PetscLogEventRegister("  Move P", cookie, &gamg_setup_events[SET15]);
1282f852f58cSMark F. Adams 
1283b4fbaa2aSMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1284b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1285b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
1286b4fbaa2aSMark F. Adams     /* create timer stages */
1287b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1288b4fbaa2aSMark F. Adams     {
1289b4fbaa2aSMark F. Adams       char str[32];
1290b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1291b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1292b4fbaa2aSMark F. Adams       PetscInt lidx;
1293b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1294b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1295b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1296b4fbaa2aSMark F. Adams       }
1297b4fbaa2aSMark F. Adams     }
1298b4fbaa2aSMark F. Adams #endif
1299b4fbaa2aSMark F. Adams   }
1300b4fbaa2aSMark F. Adams #endif
13015b89ad90SMark F. Adams   PetscFunctionReturn(0);
13025b89ad90SMark F. Adams }
13035b89ad90SMark F. Adams EXTERN_C_END
1304