xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision ed3f998397a0c8edb2d5630e649db309894f19e5)
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 
185b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
195b89ad90SMark F. Adams /*
205b89ad90SMark F. Adams    PCSetCoordinates_GAMG
215b89ad90SMark F. Adams 
225b89ad90SMark F. Adams    Input Parameter:
235b89ad90SMark F. Adams    .  pc - the preconditioner context
245b89ad90SMark F. Adams */
25a92563c5SMark F. Adams EXTERN_C_BEGIN
265b89ad90SMark F. Adams #undef __FUNCT__
275b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG"
28eb07cef2SMark F. Adams PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords )
295b89ad90SMark F. Adams {
30eb07cef2SMark F. Adams   PC_MG          *mg = (PC_MG*)a_pc->data;
315b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
326c237d78SBarry Smith   PetscErrorCode ierr;
33d3d6bff4SMark F. Adams   PetscInt       arrsz,bs,my0,kk,ii,jj,nloc,Iend;
34038e3b61SMark F. Adams   Mat            Amat = a_pc->pmat;
355b89ad90SMark F. Adams 
365b89ad90SMark F. Adams   PetscFunctionBegin;
3758471d46SMark F. Adams   PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 );
38038e3b61SMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
39d3d6bff4SMark F. Adams   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
40d3d6bff4SMark F. Adams   nloc = (Iend-my0)/bs;
41d3d6bff4SMark F. Adams   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
42038e3b61SMark F. Adams 
43d3d6bff4SMark F. Adams   pc_gamg->m_data_rows = 1;
44f6536408SMark F. Adams   if(a_coords==0 && pc_gamg->m_method==0) {
45f6536408SMark F. Adams     SETERRQ(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
46f6536408SMark F. Adams   }
47470e26f8SMark F. Adams   if( pc_gamg->m_method==0 ) pc_gamg->m_data_cols = a_ndm; /* coordinates */
48038e3b61SMark F. Adams   else{ /* SA: null space vectors */
49d3d6bff4SMark F. Adams     if(a_coords != 0 && bs==1 ) pc_gamg->m_data_cols = 1; /* scalar w/ coords and SA (not needed) */
50d3d6bff4SMark F. Adams     else if(a_coords != 0 ) pc_gamg->m_data_cols = (a_ndm==2 ? 3 : 6); /* elasticity */
51d3d6bff4SMark F. Adams     else pc_gamg->m_data_cols = bs; /* no data, force SA with constant null space vectors */
52d3d6bff4SMark F. Adams     pc_gamg->m_data_rows = bs;
53038e3b61SMark F. Adams   }
54d3d6bff4SMark F. Adams   arrsz = nloc*pc_gamg->m_data_rows*pc_gamg->m_data_cols;
555b89ad90SMark F. Adams 
56038e3b61SMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
576e3e101aSMark F. Adams   if (pc_gamg->m_data==0 || (pc_gamg->m_data_sz != arrsz)) {
585b89ad90SMark F. Adams     ierr = PetscFree( pc_gamg->m_data );  CHKERRQ(ierr);
5933a2b366SJed Brown     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->m_data ); CHKERRQ(ierr);
605b89ad90SMark F. Adams   }
61038e3b61SMark F. Adams   for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.;
62eb07cef2SMark F. Adams   pc_gamg->m_data[arrsz] = -99.;
63038e3b61SMark F. Adams   /* copy data in - column oriented */
64470e26f8SMark F. Adams   if( pc_gamg->m_method != 0 ) {
65d3d6bff4SMark F. Adams     const PetscInt M = Iend - my0;
66038e3b61SMark F. Adams     for(kk=0;kk<nloc;kk++){
67038e3b61SMark F. Adams       PetscReal *data = &pc_gamg->m_data[kk*bs];
68d3d6bff4SMark F. Adams       if( pc_gamg->m_data_cols==1 ) *data = 1.0;
69038e3b61SMark F. Adams       else {
70038e3b61SMark F. Adams         for(ii=0;ii<bs;ii++)
71038e3b61SMark F. Adams 	  for(jj=0;jj<bs;jj++)
72038e3b61SMark F. Adams 	    if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */
73038e3b61SMark F. Adams 	    else data[ii*M + jj] = 0.0;
74eb07cef2SMark F. Adams         if( a_coords != 0 ) {
75038e3b61SMark F. Adams           if( a_ndm == 2 ){ /* rotational modes */
76038e3b61SMark F. Adams             data += 2*M;
77038e3b61SMark F. Adams             data[0] = -a_coords[2*kk+1];
78038e3b61SMark F. Adams             data[1] =  a_coords[2*kk];
795b89ad90SMark F. Adams           }
80eb07cef2SMark F. Adams           else {
81038e3b61SMark F. Adams             data += 3*M;
82038e3b61SMark F. Adams             data[0] = 0.0;               data[M+0] =  a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1];
83038e3b61SMark F. Adams             data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0;               data[2*M+1] =  a_coords[3*kk];
84038e3b61SMark F. Adams             data[2] =  a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk];   data[2*M+2] = 0.0;
85038e3b61SMark F. Adams           }
86eb07cef2SMark F. Adams         }
87eb07cef2SMark F. Adams       }
88eb07cef2SMark F. Adams     }
89eb07cef2SMark F. Adams   }
90eb07cef2SMark F. Adams   else {
91038e3b61SMark F. Adams     for( kk = 0 ; kk < nloc ; kk++ ){
92038e3b61SMark F. Adams       for( ii = 0 ; ii < a_ndm ; ii++ ) {
93038e3b61SMark F. Adams         pc_gamg->m_data[ii*nloc + kk] =  a_coords[kk*a_ndm + ii];
94eb07cef2SMark F. Adams       }
95eb07cef2SMark F. Adams     }
96038e3b61SMark F. Adams   }
97eb07cef2SMark F. Adams   assert(pc_gamg->m_data[arrsz] == -99.);
98038e3b61SMark F. Adams 
995b89ad90SMark F. Adams   pc_gamg->m_data_sz = arrsz;
100eb07cef2SMark F. Adams   pc_gamg->m_dim = a_ndm;
1015b89ad90SMark F. Adams 
1025b89ad90SMark F. Adams   PetscFunctionReturn(0);
1035b89ad90SMark F. Adams }
104a92563c5SMark F. Adams EXTERN_C_END
1055b89ad90SMark F. Adams 
106d3d6bff4SMark F. Adams 
107d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
108d3d6bff4SMark F. Adams #undef __FUNCT__
109d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
110d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
111d3d6bff4SMark F. Adams {
112d3d6bff4SMark F. Adams   PetscErrorCode  ierr;
113d3d6bff4SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
114d3d6bff4SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
115d3d6bff4SMark F. Adams 
116d3d6bff4SMark F. Adams   PetscFunctionBegin;
11758471d46SMark F. Adams   if( pc_gamg->m_data != 0 ) { /* this should not happen, cleaned up in SetUp */
118d3d6bff4SMark F. Adams     ierr = PetscFree(pc_gamg->m_data); CHKERRQ(ierr);
11958471d46SMark F. Adams   }
120d3d6bff4SMark F. Adams   pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0;
121d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
122d3d6bff4SMark F. Adams }
123d3d6bff4SMark F. Adams 
1245b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1255b89ad90SMark F. Adams /*
126a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
127a147abb0SMark F. Adams      of active processors.
1285b89ad90SMark F. Adams 
1295b89ad90SMark F. Adams    Input Parameter:
1308263b398SMark F. Adams    . a_pc - parameters
1313530afc2SMark F. Adams    . a_Amat_fine - matrix on this fine (k) level
132d3d6bff4SMark F. Adams    . a_ndata_rows - size of data to move (coarse grid)
133d3d6bff4SMark F. Adams    . a_ndata_cols - size of data to move (coarse grid)
1348263b398SMark F. Adams    . a_cbs - coarse block size
1358263b398SMark F. Adams    . a_isLast -
1363530afc2SMark F. Adams    In/Output Parameter:
1373530afc2SMark F. Adams    . a_P_inout - prolongation operator to the next level (k-1)
138eb07cef2SMark F. Adams    . a_coarse_data - data that need to be moved
139afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
14011e60469SMark F. Adams    Output Parameter:
1413530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1425b89ad90SMark F. Adams */
1435cb416c2SMark F. Adams 
1445b89ad90SMark F. Adams #undef __FUNCT__
1458263b398SMark F. Adams #define __FUNCT__ "createLevel"
1468263b398SMark F. Adams PetscErrorCode createLevel( const PC a_pc,
1478263b398SMark F. Adams                             const Mat a_Amat_fine,
1488263b398SMark F. Adams                             const PetscInt a_ndata_rows,
1498263b398SMark F. Adams                             const PetscInt a_ndata_cols,
1508263b398SMark F. Adams                             const PetscInt a_cbs,
1518263b398SMark F. Adams                             const PetscBool a_isLast,
1523530afc2SMark F. Adams                             Mat *a_P_inout,
153eb07cef2SMark F. Adams                             PetscReal **a_coarse_data,
154afc97cdcSMark F. Adams                             PetscMPIInt *a_nactive_proc,
155389730f3SMark F. Adams                             Mat *a_Amat_crs
15611e60469SMark F. Adams                             )
1575b89ad90SMark F. Adams {
1588263b398SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
159486a8d0bSJed Brown   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1608263b398SMark F. Adams   const PetscBool  repart = pc_gamg->m_repart;
161389730f3SMark F. Adams   const PetscInt   min_eq_proc = pc_gamg->m_min_eq_proc, coarse_max = pc_gamg->m_coarse_eq_limit;
1625b89ad90SMark F. Adams   PetscErrorCode   ierr;
163038e3b61SMark F. Adams   Mat              Cmat,Pnew,Pold=*a_P_inout;
16411e60469SMark F. Adams   IS               new_indices,isnum;
1653530afc2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
1665a9b9e01SMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive = *a_nactive_proc;
1675a9b9e01SMark F. Adams   PetscInt         neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0,rfactor;
1685b89ad90SMark F. Adams 
1695b89ad90SMark F. Adams   PetscFunctionBegin;
17011e60469SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
17111e60469SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
172f6536408SMark F. Adams 
17311e60469SMark F. Adams   /* RAP */
17496434bc1SMark F. Adams #ifdef USE_R
17596434bc1SMark F. Adams   /* make R wih brute force for now */
17696434bc1SMark F. Adams   ierr = MatTranspose( Pold, Pnew );
17796434bc1SMark F. Adams   ierr = MatDestroy( &Pold );  CHKERRQ(ierr);
17896434bc1SMark F. Adams   ierr = MatRARt( a_Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
17996434bc1SMark F. Adams   Pold = Pnew;
18096434bc1SMark F. Adams #else
181038e3b61SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
18296434bc1SMark F. Adams #endif
183038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
184038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
185038e3b61SMark F. Adams   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
186038e3b61SMark F. Adams 
187f852f58cSMark F. Adams   /* get number of PEs to make active, reduce */
188f852f58cSMark F. Adams   ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
18981708dccSMark F. Adams   new_npe = (PetscMPIInt)((float)neq/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
190f852f58cSMark F. Adams   if( new_npe == 0 || neq < coarse_max ) new_npe = 1;
1915a9b9e01SMark F. Adams   else if (new_npe >= nactive ) new_npe = nactive; /* no change, rare */
1928263b398SMark F. Adams   if( a_isLast ) new_npe = 1;
193f852f58cSMark F. Adams 
1945a9b9e01SMark F. Adams   if( !repart && new_npe==nactive ) {
195a147abb0SMark F. Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction */
19622063be5SMark F. Adams   }
19722063be5SMark F. Adams   else {
19822063be5SMark F. Adams     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
19922063be5SMark F. Adams     const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols;
200b4fbaa2aSMark F. Adams     const PetscInt  stride0=ncrs0*a_ndata_rows;
2015a9b9e01SMark F. Adams     PetscInt        *counts,is_sz,*newproc_idx,ii,jj,kk,strideNew,*tidx,inpe,targetPE;
2025a9b9e01SMark F. Adams     IS              isnewproc;
2035a9b9e01SMark F. Adams     VecScatter      vecscat;
20422063be5SMark F. Adams     PetscScalar    *array;
20522063be5SMark F. Adams     Vec             src_crd, dest_crd;
20622063be5SMark F. Adams     PetscReal      *data = *a_coarse_data;
2075a9b9e01SMark F. Adams     IS              isscat;
208e33ef3b1SMark F. Adams 
209fd3c6afaSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
21092a756f0SMark F. Adams 
211b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
212b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
213b4fbaa2aSMark F. Adams #endif
2145a9b9e01SMark F. Adams     if( repart ) {
2155a9b9e01SMark F. Adams       /* create sub communicator  */
2165a9b9e01SMark F. Adams       Mat             adj;
2175a9b9e01SMark F. Adams 
2185a9b9e01SMark F. Adams       if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,neq);
2195a9b9e01SMark F. Adams 
2205a9b9e01SMark F. Adams       /* MatPartitioningApply call MatConvert, which is collective */
221038e3b61SMark F. Adams       if( a_cbs == 1 ) {
222038e3b61SMark F. Adams 	ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
223eb07cef2SMark F. Adams       }
224eb07cef2SMark F. Adams       else{
225038e3b61SMark F. Adams 	/* make a scalar matrix to partition */
226eb07cef2SMark F. Adams 	Mat tMat;
22758471d46SMark F. Adams 	PetscInt ncols,jj,Ii;
228b4fbaa2aSMark F. Adams 	const PetscScalar *vals;
229b4fbaa2aSMark F. Adams 	const PetscInt *idx;
2305f8cf99dSMark F. Adams 	PetscInt *d_nnz, *o_nnz;
2315ee9c036SSatish Balay 	static int llev = 0;
232b4fbaa2aSMark F. Adams 
23358471d46SMark F. Adams 	ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
2345f8cf99dSMark F. Adams 	ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr);
23558471d46SMark F. Adams 	for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) {
23658471d46SMark F. Adams 	  ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
23758471d46SMark F. Adams 	  d_nnz[jj] = ncols/a_cbs;
2385f8cf99dSMark F. Adams 	  o_nnz[jj] = ncols/a_cbs;
23958471d46SMark F. Adams 	  ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
2405f8cf99dSMark F. Adams 	  if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
2415f8cf99dSMark F. Adams 	  if( o_nnz[jj] > (neq/a_cbs-ncrs0) ) o_nnz[jj] = neq/a_cbs-ncrs0;
24258471d46SMark F. Adams 	}
2436876a03eSMark F. Adams 
244eb07cef2SMark F. Adams 	ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
245eb07cef2SMark F. Adams 				PETSC_DETERMINE, PETSC_DETERMINE,
2465f8cf99dSMark F. Adams 				0, d_nnz, 0, o_nnz,
247eb07cef2SMark F. Adams 				&tMat );
2486876a03eSMark F. Adams 	CHKERRQ(ierr);
24958471d46SMark F. Adams 	ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
2505f8cf99dSMark F. Adams 	ierr = PetscFree( o_nnz ); CHKERRQ(ierr);
251eb07cef2SMark F. Adams 
25222063be5SMark F. Adams 	for ( ii = Istart0; ii < Iend0; ii++ ) {
25322063be5SMark F. Adams 	  PetscInt dest_row = ii/a_cbs;
25422063be5SMark F. Adams 	  ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
255eb07cef2SMark F. Adams 	  for( jj = 0 ; jj < ncols ; jj++ ){
256038e3b61SMark F. Adams 	    PetscInt dest_col = idx[jj]/a_cbs;
257eb07cef2SMark F. Adams 	    PetscScalar v = 1.0;
258eb07cef2SMark F. Adams 	    ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
259eb07cef2SMark F. Adams 	  }
26022063be5SMark F. Adams 	  ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
261eb07cef2SMark F. Adams 	}
262eb07cef2SMark F. Adams 	ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
263eb07cef2SMark F. Adams 	ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
264eb07cef2SMark F. Adams 
265b4fbaa2aSMark F. Adams 	if( llev++ == -1 ) {
266b4fbaa2aSMark F. Adams 	  PetscViewer viewer; char fname[32];
267b4fbaa2aSMark F. Adams 	  sprintf(fname,"part_mat_%d.mat",llev);
268b4fbaa2aSMark F. Adams 	  PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
269b4fbaa2aSMark F. Adams 	  ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
270b4fbaa2aSMark F. Adams 	  ierr = PetscViewerDestroy( &viewer );
271b4fbaa2aSMark F. Adams 	}
272b4fbaa2aSMark F. Adams 
273eb07cef2SMark F. Adams 	ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
274eb07cef2SMark F. Adams 
275eb07cef2SMark F. Adams 	ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
276eb07cef2SMark F. Adams       }
277f150b916SMark F. Adams 
2785a9b9e01SMark F. Adams       { /* partition: get newproc_idx, set is_sz */
2795a9b9e01SMark F. Adams 	char prefix[256];
2805a9b9e01SMark F. Adams 	const char *pcpre;
281b4fbaa2aSMark F. Adams 	const PetscInt *is_idx;
282b4fbaa2aSMark F. Adams 	MatPartitioning  mpart;
283a4b7d37bSMark F. Adams 	IS proc_is;
2842f03bc48SMark F. Adams 
2855a9b9e01SMark F. Adams 	ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr);
2865ef31b24SMark F. Adams 	ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
28759a0be82SJed Brown 	ierr = PCGetOptionsPrefix(a_pc,&pcpre);CHKERRQ(ierr);
28859a0be82SJed Brown 	ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
28959a0be82SJed Brown 	ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
29011e60469SMark F. Adams 	ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
2915ef31b24SMark F. Adams 	ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
292a4b7d37bSMark F. Adams 	ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr);
29311e60469SMark F. Adams 	ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
2945a9b9e01SMark F. Adams 
2955ef31b24SMark F. Adams 	/* collect IS info */
296a4b7d37bSMark F. Adams 	ierr = ISGetLocalSize( proc_is, &is_sz );       CHKERRQ(ierr);
2978263b398SMark F. Adams 	ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr);
298a4b7d37bSMark F. Adams 	ierr = ISGetIndices( proc_is, &is_idx );        CHKERRQ(ierr);
299a147abb0SMark F. Adams 	NN = 1; /* bring to "front" of machine */
300a147abb0SMark F. Adams 	/*NN = npe/new_npe;*/ /* spread partitioning across machine */
301b4fbaa2aSMark F. Adams 	for( kk = jj = 0 ; kk < is_sz ; kk++ ){
302afc97cdcSMark F. Adams 	  for( ii = 0 ; ii < a_cbs ; ii++, jj++ ){
3038263b398SMark F. Adams 	    newproc_idx[jj] = is_idx[kk] * NN; /* distribution */
304eb07cef2SMark F. Adams 	  }
3055ef31b24SMark F. Adams 	}
306a4b7d37bSMark F. Adams 	ierr = ISRestoreIndices( proc_is, &is_idx );     CHKERRQ(ierr);
307a4b7d37bSMark F. Adams 	ierr = ISDestroy( &proc_is );                    CHKERRQ(ierr);
308b4fbaa2aSMark F. Adams 
309b4fbaa2aSMark F. Adams 	is_sz *= a_cbs;
3105ef31b24SMark F. Adams       }
3115ef31b24SMark F. Adams       ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
3125a9b9e01SMark F. Adams 
3138263b398SMark F. Adams       ierr = ISCreateGeneral( wcomm, is_sz, newproc_idx, PETSC_COPY_VALUES, &isnewproc );
3145a9b9e01SMark F. Adams       CHKERRQ(ierr);
3158263b398SMark F. Adams       if( newproc_idx != 0 ) {
3168263b398SMark F. Adams 	ierr = PetscFree( newproc_idx );  CHKERRQ(ierr);
3175ef31b24SMark F. Adams       }
3185a9b9e01SMark F. Adams     }
3195a9b9e01SMark F. Adams     else { /* simple aggreagtion of parts */
3205a9b9e01SMark F. Adams       /* find factor */
3215a9b9e01SMark F. Adams       if( new_npe == 1 ) rfactor = npe; /* easy */
3225a9b9e01SMark F. Adams       else {
3235a9b9e01SMark F. Adams 	PetscReal best_fact = 0.;
3245a9b9e01SMark F. Adams 	jj = -1;
3255a9b9e01SMark F. Adams 	for( kk = 1 ; kk <= npe ; kk++ ){
3265a9b9e01SMark F. Adams 	  if( npe%kk==0 ) { /* a candidate */
3275a9b9e01SMark F. Adams 	    PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe;
3285a9b9e01SMark F. Adams 	    if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */
3295a9b9e01SMark F. Adams 	    if( fact > best_fact ) {
3305a9b9e01SMark F. Adams 	      best_fact = fact; jj = kk;
3315a9b9e01SMark F. Adams 	    }
3325a9b9e01SMark F. Adams 	  }
3335a9b9e01SMark F. Adams 	}
3345a9b9e01SMark F. Adams 	if(jj!=-1) rfactor = jj;
3355a9b9e01SMark F. Adams 	else rfactor = 1; /* prime? */
3365a9b9e01SMark F. Adams       }
3375a9b9e01SMark F. Adams       new_npe = npe/rfactor;
3385a9b9e01SMark F. Adams 
3395a9b9e01SMark F. Adams       if( new_npe==nactive ) {
3405a9b9e01SMark F. Adams 	*a_Amat_crs = Cmat; /* output - no repartitioning or reduction */
3415a9b9e01SMark F. Adams 	ierr = PetscFree( counts );  CHKERRQ(ierr);
3425a9b9e01SMark F. Adams 	*a_nactive_proc = new_npe; /* output */
3435a9b9e01SMark F. Adams 	if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s aggregate processors noop: new_npe=%d, neq=%d\n",mype,__FUNCT__,new_npe,neq);
3445a9b9e01SMark F. Adams 	PetscFunctionReturn(0);
3455a9b9e01SMark F. Adams       }
3465a9b9e01SMark F. Adams 
3475a9b9e01SMark F. Adams       /* reduce - isnewproc */
3485a9b9e01SMark F. Adams       targetPE = mype/rfactor;
3495a9b9e01SMark F. Adams 
3505a9b9e01SMark F. Adams       if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s aggregate processors: npe: %d --> %d, neq=%d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,neq);
3515a9b9e01SMark F. Adams       is_sz = ncrs0*a_cbs;
3525a9b9e01SMark F. Adams       ierr = ISCreateStride( wcomm, is_sz, targetPE, 0, &isnewproc );
3535a9b9e01SMark F. Adams       CHKERRQ(ierr);
3545a9b9e01SMark F. Adams     }
355e33ef3b1SMark F. Adams 
35611e60469SMark F. Adams     /*
35711e60469SMark F. Adams       Create an index set from the isnewproc index set to indicate the mapping TO
35811e60469SMark F. Adams     */
35911e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
36011e60469SMark F. Adams     /*
36111e60469SMark F. Adams       Determine how many elements are assigned to each processor
36211e60469SMark F. Adams     */
363fd3c6afaSMark F. Adams     inpe = npe;
364fd3c6afaSMark F. Adams     ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr);
36511e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
366038e3b61SMark F. Adams     ncrs_new = counts[mype]/a_cbs;
36722063be5SMark F. Adams     strideNew = ncrs_new*a_ndata_rows;
368b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
369b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
370b4fbaa2aSMark F. Adams #endif
37122063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
37211e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
373d3d6bff4SMark F. Adams     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
374470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
37511e60469SMark F. Adams     /*
376d3d6bff4SMark F. Adams      There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
377d3d6bff4SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
37811e60469SMark F. Adams      */
37992a756f0SMark F. Adams     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
38011e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
381038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
382d3d6bff4SMark F. Adams       PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
383d3d6bff4SMark F. Adams       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
38411e60469SMark F. Adams     }
385038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
386d3d6bff4SMark F. Adams     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
38711e60469SMark F. Adams     CHKERRQ(ierr);
38892a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
38911e60469SMark F. Adams     /*
39011e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
39111e60469SMark F. Adams      */
392d3d6bff4SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
393d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
394d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs0 ; ii++) {
395d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
396d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
397676e1743SMark F. Adams           PetscScalar tt = (PetscScalar)data[ix];
398676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
399d3d6bff4SMark F. Adams 	}
400038e3b61SMark F. Adams       }
401eb07cef2SMark F. Adams     }
402eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
403eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
40411e60469SMark F. Adams     /*
40511e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
40611e60469SMark F. Adams       to the correct processor
40711e60469SMark F. Adams     */
40811e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
40911e60469SMark F. Adams     CHKERRQ(ierr);
41011e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
41111e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
41211e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
41311e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
41411e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
41511e60469SMark F. Adams     /*
41611e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
41711e60469SMark F. Adams     */
418eb07cef2SMark F. Adams     ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
419d3d6bff4SMark F. Adams     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
420eb07cef2SMark F. Adams 
42111e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
422eb07cef2SMark F. Adams     data = *a_coarse_data;
423d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
424d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs_new ; ii++) {
425d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
426d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
427d3d6bff4SMark F. Adams 	  data[ix] = PetscRealPart(array[jx]);
428d3d6bff4SMark F. Adams 	  array[jx] = 1.e300;
429d3d6bff4SMark F. Adams 	}
430038e3b61SMark F. Adams       }
431038e3b61SMark F. Adams     }
43211e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
43311e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
434*ed3f9983SMark F. Adams #if defined PETSC_USE_LOG
435*ed3f9983SMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
436*ed3f9983SMark F. Adams #endif
43711e60469SMark F. Adams     /*
43811e60469SMark F. Adams       Invert for MatGetSubMatrix
43911e60469SMark F. Adams     */
440038e3b61SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
44111e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
44211e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
443*ed3f9983SMark F. Adams #if defined PETSC_USE_LOG
444*ed3f9983SMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
445*ed3f9983SMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
446*ed3f9983SMark F. Adams #endif
44711e60469SMark F. Adams     /* A_crs output */
448038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
44911e60469SMark F. Adams     CHKERRQ(ierr);
450eb07cef2SMark F. Adams 
451038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
452e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
453038e3b61SMark F. Adams     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
454*ed3f9983SMark F. Adams #if defined PETSC_USE_LOG
455*ed3f9983SMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
456*ed3f9983SMark F. Adams #endif
45711e60469SMark F. Adams     /* prolongator */
45811e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
45911e60469SMark F. Adams     {
46011e60469SMark F. Adams       IS findices;
461*ed3f9983SMark F. Adams #if defined PETSC_USE_LOG
462*ed3f9983SMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
463*ed3f9983SMark F. Adams #endif
46411e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
46596434bc1SMark F. Adams #ifdef USE_R
46696434bc1SMark F. Adams       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
46796434bc1SMark F. Adams #else
46811e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
46996434bc1SMark F. Adams #endif
47011e60469SMark F. Adams       CHKERRQ(ierr);
47111e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
472*ed3f9983SMark F. Adams #if defined PETSC_USE_LOG
473*ed3f9983SMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
474*ed3f9983SMark F. Adams #endif
47511e60469SMark F. Adams     }
4763530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
477a147abb0SMark F. Adams     *a_P_inout = Pnew; /* output - repartitioned */
47811e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
47992a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
480e33ef3b1SMark F. Adams   }
4815b89ad90SMark F. Adams 
4825a9b9e01SMark F. Adams   *a_nactive_proc = new_npe; /* output */
4835a9b9e01SMark F. Adams 
4845b89ad90SMark F. Adams   PetscFunctionReturn(0);
4855b89ad90SMark F. Adams }
4865b89ad90SMark F. Adams 
4875b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4885b89ad90SMark F. Adams /*
4895b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4905b89ad90SMark F. Adams                     by setting data structures and options.
4915b89ad90SMark F. Adams 
4925b89ad90SMark F. Adams    Input Parameter:
4935b89ad90SMark F. Adams .  pc - the preconditioner context
4945b89ad90SMark F. Adams 
4955b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4965b89ad90SMark F. Adams 
4975b89ad90SMark F. Adams    Notes:
4985b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4995b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
5005b89ad90SMark F. Adams */
5015b89ad90SMark F. Adams #undef __FUNCT__
5025b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
503eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
5045b89ad90SMark F. Adams {
5055b89ad90SMark F. Adams   PetscErrorCode  ierr;
506eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
5075b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
50858471d46SMark F. Adams   PC_MG_Levels   **mglevels = mg->levels;
509eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
510d3d6bff4SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
511eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
5123530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
513038e3b61SMark F. Adams   PetscBool        isOK;
514587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
515587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
516737a81a9SMark F. Adams   MatInfo          info;
5175ef31b24SMark F. Adams 
5185b89ad90SMark F. Adams   PetscFunctionBegin;
51903a628feSMark F. Adams   pc_gamg->m_count++;
520fca9ed99SMark F. Adams 
52158471d46SMark F. Adams   if( a_pc->setupcalled > 0 ) {
52203a628feSMark F. Adams     /* just do Galerkin grids */
52358471d46SMark F. Adams     Mat B,dA,dB;
52458471d46SMark F. Adams 
52558471d46SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
52658471d46SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
52758471d46SMark F. Adams 
5285f8cf99dSMark F. Adams     if( pc_gamg->m_Nlevels > 1 ) {
52958471d46SMark F. Adams     /* currently only handle case where mat and pmat are the same on coarser levels */
53058471d46SMark F. Adams     ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
53158471d46SMark F. Adams     /* (re)set to get dirty flag */
53258471d46SMark F. Adams     ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
53358471d46SMark F. Adams     ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr);
53458471d46SMark F. Adams 
53558471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-2; level>-1; level--) {
53658471d46SMark F. Adams       ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
53703a628feSMark F. Adams       /* the first time through the matrix structure has changed from repartitioning */
53803a628feSMark F. Adams       if( pc_gamg->m_count == 2 ) {
53903a628feSMark F. Adams         ierr = MatDestroy( &B );  CHKERRQ(ierr);
54003a628feSMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
54103a628feSMark F. Adams         mglevels[level]->A = B;
54203a628feSMark F. Adams       }
54303a628feSMark F. Adams       else {
54458471d46SMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
54503a628feSMark F. Adams       }
54658471d46SMark F. Adams       ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
54758471d46SMark F. Adams       dB = B;
54858471d46SMark F. Adams       /* setup KSP/PC */
54958471d46SMark F. Adams       ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
55058471d46SMark F. Adams     }
5515f8cf99dSMark F. Adams     }
552f6536408SMark F. Adams #define PRINT_MATS PETSC_FALSE
55358471d46SMark F. Adams     /* plot levels - A */
55458471d46SMark F. Adams     if( PRINT_MATS ) {
55558471d46SMark F. Adams       for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){
55658471d46SMark F. Adams         PetscViewer viewer;
55758471d46SMark F. Adams         char fname[32]; KSP smoother; Mat Tmat, TTm;
55858471d46SMark F. Adams         ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
55958471d46SMark F. Adams         ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr);
56003a628feSMark F. Adams         sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
56158471d46SMark F. Adams         ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
56258471d46SMark F. Adams         ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
56358471d46SMark F. Adams         ierr = MatView( Tmat, viewer ); CHKERRQ(ierr);
56458471d46SMark F. Adams         ierr = PetscViewerDestroy( &viewer );
56558471d46SMark F. Adams       }
56658471d46SMark F. Adams     }
56758471d46SMark F. Adams 
56803a628feSMark F. Adams     a_pc->setupcalled = 2;
56903a628feSMark F. Adams 
57058471d46SMark F. Adams     PetscFunctionReturn(0);
571eb07cef2SMark F. Adams   }
572f6536408SMark F. Adams 
573baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
574baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
5755b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
5765b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
5773530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
578eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
579eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
580eb07cef2SMark F. Adams 
581038e3b61SMark F. Adams   /* get data of not around */
582038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
583038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
584eb07cef2SMark F. Adams   }
585eb07cef2SMark F. Adams   data = pc_gamg->m_data;
586038e3b61SMark F. Adams 
587eb07cef2SMark F. Adams   /* Get A_i and R_i */
588737a81a9SMark F. Adams   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
589bed7c2b9SMark F. Adams   if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
5902c047e2dSMark F. Adams 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,
5912c047e2dSMark F. Adams 	      (int)(info.nz_used/(PetscReal)N),npe);
5928f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
5934ef23d27SMark F. Adams         level < (pc_gamg->m_Nlevels-1) && (level==0 || M>pc_gamg->m_coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
5940205a208SMark F. Adams         level++ ){
5955b89ad90SMark F. Adams     level1 = level + 1;
596b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
597b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
598b4fbaa2aSMark F. Adams #endif
599b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
600b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
601b4fbaa2aSMark F. Adams #endif
602389730f3SMark F. Adams     ierr = createProlongation(Aarr[level], data, level, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level], pc_gamg );
6035b89ad90SMark F. Adams     CHKERRQ(ierr);
604d3d6bff4SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
605b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
606b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
6075c7a89a6SBarry Smith 
608b4fbaa2aSMark F. Adams #endif
609baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
610baa4e9faSMark F. Adams     if( isOK ) {
611b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
612b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
613b4fbaa2aSMark F. Adams #endif
6148263b398SMark F. Adams       ierr = createLevel( a_pc, Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs,
6158263b398SMark F. Adams                           (PetscBool)(level == pc_gamg->m_Nlevels-2),
616389730f3SMark F. Adams                           &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
6173530afc2SMark F. Adams       CHKERRQ(ierr);
618b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
619b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
620b4fbaa2aSMark F. Adams #endif
6213530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
622737a81a9SMark F. Adams       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
623bed7c2b9SMark F. Adams       if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
6242c047e2dSMark F. Adams 		  mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols,
6252c047e2dSMark F. Adams 		  (int)(info.nz_used/(PetscReal)N),nactivepe);
626e33ef3b1SMark F. Adams       /* coarse grids with SA can have zero row/cols from singleton aggregates */
627737a81a9SMark F. Adams 
62881708dccSMark F. Adams       /* stop if one node */
62981708dccSMark F. Adams       if( M/pc_gamg->m_data_cols < 2 ) {
63081708dccSMark F. Adams         level++;
63181708dccSMark F. Adams         break;
63281708dccSMark F. Adams       }
63381708dccSMark F. Adams 
63481708dccSMark F. Adams       if (PETSC_FALSE) {
635785cba28SMark F. Adams         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
636785cba28SMark F. Adams         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
637e33ef3b1SMark F. Adams         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
638785cba28SMark F. Adams         nloceq = Iend-Istart;
639e33ef3b1SMark F. Adams         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
640e33ef3b1SMark F. Adams         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
641e33ef3b1SMark F. Adams         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
642785cba28SMark F. Adams         for(kk=0;kk<nloceq;kk++){
643e33ef3b1SMark F. Adams           if(data_arr[kk]==0.0) {
644785cba28SMark F. Adams             id = kk + Istart;
645e33ef3b1SMark F. Adams             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
646e33ef3b1SMark F. Adams             CHKERRQ(ierr);
647ecd57171SMark F. Adams             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1);
648e33ef3b1SMark F. Adams           }
649e33ef3b1SMark F. Adams         }
650e33ef3b1SMark F. Adams         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
651e33ef3b1SMark F. Adams         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
652e33ef3b1SMark F. Adams         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
653e33ef3b1SMark F. Adams         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
654e33ef3b1SMark F. Adams       }
655486a8d0bSJed Brown     } else {
656be544d3cSMark F. Adams       coarse_data = 0;
657baa4e9faSMark F. Adams       break;
658baa4e9faSMark F. Adams     }
659eb07cef2SMark F. Adams     data = coarse_data;
660b4fbaa2aSMark F. Adams 
661b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
662b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
663b4fbaa2aSMark F. Adams #endif
6645b89ad90SMark F. Adams   }
665be544d3cSMark F. Adams   if( coarse_data ) {
666eb07cef2SMark F. Adams     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
667be544d3cSMark F. Adams   }
668bed7c2b9SMark F. Adams   if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
6695b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
6705b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
6715b89ad90SMark F. Adams   fine_level = level;
672eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
6735b89ad90SMark F. Adams 
6745f8cf99dSMark F. Adams   if( pc_gamg->m_Nlevels > 1 ) { /* don't setup MG if  */
675fc4362bfSMark F. Adams   /* set default smoothers */
676587fa25dSMark F. Adams   for ( lidx = 1, level = pc_gamg->m_Nlevels-2;
677587fa25dSMark F. Adams         lidx <= fine_level;
678587fa25dSMark F. Adams         lidx++, level--) {
6795745e0f5SMark F. Adams     PetscReal emax, emin;
6805b89ad90SMark F. Adams     KSP smoother; PC subpc;
681f6536408SMark F. Adams     PetscBool isCheb;
682f6536408SMark F. Adams     /* set defaults */
683587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
6845b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
685f6536408SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
686f6536408SMark F. Adams     /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */
687f6536408SMark F. Adams     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
688f6536408SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
689f6536408SMark F. Adams     /* overide defaults with input parameters */
690f6536408SMark F. Adams     ierr = KSPSetFromOptions( smoother ); CHKERRQ(ierr);
691f6536408SMark F. Adams 
692f6536408SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
693f6536408SMark F. Adams     /* do my own cheby */
694f6536408SMark F. Adams     ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr);
69541139db5SMark F. Adams 
696f6536408SMark F. Adams     if( isCheb ) {
697f6536408SMark F. Adams       ierr = PetscTypeCompare( (PetscObject)subpc, PETSC_GAMG_SMOOTHER, &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];
701fc4362bfSMark F. Adams         Vec bb, xx; PC pc;
702f6536408SMark F. Adams         const PCType type;
703038e3b61SMark F. Adams 
704f6536408SMark F. Adams         ierr = PCGetType( subpc, &type );   CHKERRQ(ierr);
7055745e0f5SMark F. Adams         ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
7065745e0f5SMark F. Adams         ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
707fc4362bfSMark F. Adams         {
708fc4362bfSMark F. Adams           PetscRandom    rctx;
709fc4362bfSMark F. Adams           ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
710fc4362bfSMark F. Adams           ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
711fc4362bfSMark F. Adams           ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
712fc4362bfSMark F. Adams           ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
7135b89ad90SMark F. Adams         }
714fc4362bfSMark F. Adams         ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
71574050f66SMark F. Adams         /* ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr); */
716038e3b61SMark F. Adams         ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
717fc4362bfSMark F. Adams         CHKERRQ(ierr);
718fc4362bfSMark F. Adams         ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
719fc4362bfSMark F. Adams 
720f6536408SMark F. Adams         ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
721f6536408SMark F. Adams         ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
722f6536408SMark F. Adams 
723f6536408SMark F. Adams         ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
724f6536408SMark F. Adams         ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
725f6536408SMark F. Adams         ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
726f6536408SMark F. Adams         ierr = PCSetType( pc, type ); CHKERRQ(ierr); /* should be same as eigen estimates op. */
727f6536408SMark F. Adams 
728fc4362bfSMark F. Adams         ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
7295a9b9e01SMark F. Adams 
7305a9b9e01SMark F. Adams 	/* solve - keep stuff out of logging */
7315a9b9e01SMark F. Adams 	ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
7325a9b9e01SMark F. Adams 	ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
733fc4362bfSMark F. Adams         ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
7345a9b9e01SMark F. Adams 	ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
7355a9b9e01SMark F. Adams 	ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
7365a9b9e01SMark F. Adams 
737fc4362bfSMark F. Adams         ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
7385a9b9e01SMark F. Adams 
739fc4362bfSMark F. Adams         ierr = VecDestroy( &xx );       CHKERRQ(ierr);
740fc4362bfSMark F. Adams         ierr = VecDestroy( &bb );       CHKERRQ(ierr);
741fc4362bfSMark F. Adams         ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
742f6536408SMark F. Adams 
743f6536408SMark F. Adams         if (pc_gamg->m_verbose) {
74441139db5SMark F. Adams           PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s PC setup max eigen=%e min=%e\n",
74541139db5SMark F. Adams                       __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   }
7595745e0f5SMark F. Adams   {
760ecd57171SMark F. Adams     /* coarse grid */
761ecd57171SMark F. Adams     KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
7625745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
763eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
76458471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
7655745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
766ecd57171SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
767ecd57171SMark F. Adams     ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
768ecd57171SMark F. Adams     ierr = PCSetUp( subpc ); CHKERRQ(ierr);
769ecd57171SMark F. Adams     ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
770ecd57171SMark F. Adams     assert(ii==1);
771ecd57171SMark F. Adams     ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
772ecd57171SMark F. Adams     ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
773fc4362bfSMark F. Adams   }
774737a81a9SMark F. Adams 
775fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
776eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
7775b89ad90SMark F. Adams   {
7785b89ad90SMark F. Adams     PetscBool galerkin;
779eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
7805b89ad90SMark F. Adams     if(galerkin){
7815b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
7825b89ad90SMark F. Adams     }
7835b89ad90SMark F. Adams   }
7845745e0f5SMark F. Adams 
78558471d46SMark F. Adams   /* plot levels - R/P */
78658471d46SMark F. Adams   if( PRINT_MATS ) {
78758471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-1;level>0;level--){
78858471d46SMark F. Adams       PetscViewer viewer;
78958471d46SMark F. Adams       char fname[32];
79003a628feSMark F. Adams       sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
79111e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
7925b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
793038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
7945b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
79503a628feSMark F. Adams       sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
796e33ef3b1SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
797e33ef3b1SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
798e33ef3b1SMark F. Adams       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
799e33ef3b1SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
8005b89ad90SMark F. Adams     }
80158471d46SMark F. Adams   }
80258471d46SMark F. Adams 
80358471d46SMark F. Adams   /* set interpolation between the levels, clean up */
80458471d46SMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
80558471d46SMark F. Adams        lidx<fine_level;
80658471d46SMark F. Adams        lidx++, level--){
80758471d46SMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
808587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
809587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
8105b89ad90SMark F. Adams   }
8115745e0f5SMark F. Adams 
8125b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
813eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
814eb07cef2SMark F. Adams   ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
81503a628feSMark F. Adams   a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
81658471d46SMark F. Adams 
81758471d46SMark F. Adams   {
81858471d46SMark F. Adams     KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
81958471d46SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
82058471d46SMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
82158471d46SMark F. Adams     ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
82258471d46SMark F. Adams   }
8235f8cf99dSMark F. Adams   }
8245f8cf99dSMark F. Adams   else {
8255f8cf99dSMark F. Adams     KSP smoother;
8265f8cf99dSMark F. Adams     if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__);
8275f8cf99dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
8285f8cf99dSMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
8295f8cf99dSMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
8305f8cf99dSMark F. Adams     /* setupcalled is set to 0 so that MG is setup from scratch */
8315f8cf99dSMark F. Adams     a_pc->setupcalled = 0;
8325f8cf99dSMark F. Adams     ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
8335f8cf99dSMark F. Adams     a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
8345f8cf99dSMark F. Adams   }
8355b89ad90SMark F. Adams   PetscFunctionReturn(0);
8365b89ad90SMark F. Adams }
8375b89ad90SMark F. Adams 
838eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8395b89ad90SMark F. Adams /*
8405b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8415b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8425b89ad90SMark F. Adams 
8435b89ad90SMark F. Adams    Input Parameter:
8445b89ad90SMark F. Adams .  pc - the preconditioner context
8455b89ad90SMark F. Adams 
8465b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8475b89ad90SMark F. Adams */
8485b89ad90SMark F. Adams #undef __FUNCT__
8495b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
8505b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8515b89ad90SMark F. Adams {
8525b89ad90SMark F. Adams   PetscErrorCode  ierr;
8535b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
8545b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
8555b89ad90SMark F. Adams 
8565b89ad90SMark F. Adams   PetscFunctionBegin;
8575b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8585b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8595b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8605b89ad90SMark F. Adams   PetscFunctionReturn(0);
8615b89ad90SMark F. Adams }
8625b89ad90SMark F. Adams 
863676e1743SMark F. Adams 
864676e1743SMark F. Adams #undef __FUNCT__
865676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
866676e1743SMark F. Adams /*@
867676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
868676e1743SMark F. Adams    processor reduction.
869676e1743SMark F. Adams 
870676e1743SMark F. Adams    Not Collective on PC
871676e1743SMark F. Adams 
872676e1743SMark F. Adams    Input Parameters:
873676e1743SMark F. Adams .  pc - the preconditioner context
874676e1743SMark F. Adams 
875676e1743SMark F. Adams 
876676e1743SMark F. Adams    Options Database Key:
877676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
878676e1743SMark F. Adams 
879676e1743SMark F. Adams    Level: intermediate
880676e1743SMark F. Adams 
881676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
882676e1743SMark F. Adams 
883676e1743SMark F. Adams .seealso: ()
884676e1743SMark F. Adams @*/
885676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
886676e1743SMark F. Adams {
887676e1743SMark F. Adams   PetscErrorCode ierr;
888676e1743SMark F. Adams 
889676e1743SMark F. Adams   PetscFunctionBegin;
890676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
891676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
892676e1743SMark F. Adams   PetscFunctionReturn(0);
893676e1743SMark F. Adams }
894676e1743SMark F. Adams 
895676e1743SMark F. Adams EXTERN_C_BEGIN
896676e1743SMark F. Adams #undef __FUNCT__
897676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
898676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
899676e1743SMark F. Adams {
900c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
901c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
902676e1743SMark F. Adams 
903676e1743SMark F. Adams   PetscFunctionBegin;
904c20e4228SMark F. Adams   if(n>0) pc_gamg->m_min_eq_proc = n;
905676e1743SMark F. Adams   PetscFunctionReturn(0);
906676e1743SMark F. Adams }
907676e1743SMark F. Adams EXTERN_C_END
908676e1743SMark F. Adams 
909676e1743SMark F. Adams #undef __FUNCT__
910389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
911389730f3SMark F. Adams /*@
912389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
913389730f3SMark F. Adams 
914389730f3SMark F. Adams  Collective on PC
915389730f3SMark F. Adams 
916389730f3SMark F. Adams    Input Parameters:
917389730f3SMark F. Adams .  pc - the preconditioner context
918389730f3SMark F. Adams 
919389730f3SMark F. Adams 
920389730f3SMark F. Adams    Options Database Key:
921389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
922389730f3SMark F. Adams 
923389730f3SMark F. Adams    Level: intermediate
924389730f3SMark F. Adams 
925389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
926389730f3SMark F. Adams 
927389730f3SMark F. Adams .seealso: ()
928389730f3SMark F. Adams  @*/
929389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
930389730f3SMark F. Adams {
931389730f3SMark F. Adams   PetscErrorCode ierr;
932389730f3SMark F. Adams 
933389730f3SMark F. Adams   PetscFunctionBegin;
934389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
935389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
936389730f3SMark F. Adams   PetscFunctionReturn(0);
937389730f3SMark F. Adams }
938389730f3SMark F. Adams 
939389730f3SMark F. Adams EXTERN_C_BEGIN
940389730f3SMark F. Adams #undef __FUNCT__
941389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
942389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
943389730f3SMark F. Adams {
944389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
945389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
946389730f3SMark F. Adams 
947389730f3SMark F. Adams   PetscFunctionBegin;
948389730f3SMark F. Adams   if(n>0) pc_gamg->m_coarse_eq_limit = n;
949389730f3SMark F. Adams   PetscFunctionReturn(0);
950389730f3SMark F. Adams }
951389730f3SMark F. Adams EXTERN_C_END
952389730f3SMark F. Adams 
953389730f3SMark F. Adams #undef __FUNCT__
9548263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
955676e1743SMark F. Adams /*@
9568263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
957676e1743SMark F. Adams 
958676e1743SMark F. Adams    Collective on PC
959676e1743SMark F. Adams 
960676e1743SMark F. Adams    Input Parameters:
961676e1743SMark F. Adams .  pc - the preconditioner context
962676e1743SMark F. Adams 
963676e1743SMark F. Adams 
964676e1743SMark F. Adams    Options Database Key:
9658263b398SMark F. Adams .  -pc_gamg_repartition
966676e1743SMark F. Adams 
967676e1743SMark F. Adams    Level: intermediate
968676e1743SMark F. Adams 
969676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
970676e1743SMark F. Adams 
971676e1743SMark F. Adams .seealso: ()
972676e1743SMark F. Adams @*/
9738263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
974676e1743SMark F. Adams {
975676e1743SMark F. Adams   PetscErrorCode ierr;
976676e1743SMark F. Adams 
977676e1743SMark F. Adams   PetscFunctionBegin;
978676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9798263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
980676e1743SMark F. Adams   PetscFunctionReturn(0);
981676e1743SMark F. Adams }
982676e1743SMark F. Adams 
983676e1743SMark F. Adams EXTERN_C_BEGIN
984676e1743SMark F. Adams #undef __FUNCT__
9858263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
9868263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
987676e1743SMark F. Adams {
988c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
989c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
990676e1743SMark F. Adams 
991676e1743SMark F. Adams   PetscFunctionBegin;
9928263b398SMark F. Adams   pc_gamg->m_repart = n;
993676e1743SMark F. Adams   PetscFunctionReturn(0);
994676e1743SMark F. Adams }
995676e1743SMark F. Adams EXTERN_C_END
996676e1743SMark F. Adams 
997676e1743SMark F. Adams #undef __FUNCT__
9984ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
9994ef23d27SMark F. Adams /*@
10004ef23d27SMark F. Adams    PCGAMGSetNlevels -
10014ef23d27SMark F. Adams 
10024ef23d27SMark F. Adams    Not collective on PC
10034ef23d27SMark F. Adams 
10044ef23d27SMark F. Adams    Input Parameters:
10054ef23d27SMark F. Adams .  pc - the preconditioner context
10064ef23d27SMark F. Adams 
10074ef23d27SMark F. Adams 
10084ef23d27SMark F. Adams    Options Database Key:
10094ef23d27SMark F. Adams .  -pc_mg_levels
10104ef23d27SMark F. Adams 
10114ef23d27SMark F. Adams    Level: intermediate
10124ef23d27SMark F. Adams 
10134ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
10144ef23d27SMark F. Adams 
10154ef23d27SMark F. Adams .seealso: ()
10164ef23d27SMark F. Adams @*/
10174ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
10184ef23d27SMark F. Adams {
10194ef23d27SMark F. Adams   PetscErrorCode ierr;
10204ef23d27SMark F. Adams 
10214ef23d27SMark F. Adams   PetscFunctionBegin;
10224ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10234ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
10244ef23d27SMark F. Adams   PetscFunctionReturn(0);
10254ef23d27SMark F. Adams }
10264ef23d27SMark F. Adams 
10274ef23d27SMark F. Adams EXTERN_C_BEGIN
10284ef23d27SMark F. Adams #undef __FUNCT__
10294ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
10304ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
10314ef23d27SMark F. Adams {
10324ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
10334ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
10344ef23d27SMark F. Adams 
10354ef23d27SMark F. Adams   PetscFunctionBegin;
10364ef23d27SMark F. Adams   pc_gamg->m_Nlevels = n;
10374ef23d27SMark F. Adams   PetscFunctionReturn(0);
10384ef23d27SMark F. Adams }
10394ef23d27SMark F. Adams EXTERN_C_END
10404ef23d27SMark F. Adams 
10414ef23d27SMark F. Adams #undef __FUNCT__
10423542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
10433542efc5SMark F. Adams /*@
10443542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
10453542efc5SMark F. Adams 
10463542efc5SMark F. Adams    Not collective on PC
10473542efc5SMark F. Adams 
10483542efc5SMark F. Adams    Input Parameters:
10493542efc5SMark F. Adams .  pc - the preconditioner context
10503542efc5SMark F. Adams 
10513542efc5SMark F. Adams 
10523542efc5SMark F. Adams    Options Database Key:
10533542efc5SMark F. Adams .  -pc_gamg_threshold
10543542efc5SMark F. Adams 
10553542efc5SMark F. Adams    Level: intermediate
10563542efc5SMark F. Adams 
10573542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
10583542efc5SMark F. Adams 
10593542efc5SMark F. Adams .seealso: ()
10603542efc5SMark F. Adams @*/
10613542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
10623542efc5SMark F. Adams {
10633542efc5SMark F. Adams   PetscErrorCode ierr;
10643542efc5SMark F. Adams 
10653542efc5SMark F. Adams   PetscFunctionBegin;
10663542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10673542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
10683542efc5SMark F. Adams   PetscFunctionReturn(0);
10693542efc5SMark F. Adams }
10703542efc5SMark F. Adams 
10713542efc5SMark F. Adams EXTERN_C_BEGIN
10723542efc5SMark F. Adams #undef __FUNCT__
10733542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
10743542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
10753542efc5SMark F. Adams {
1076c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1077c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
10783542efc5SMark F. Adams 
10793542efc5SMark F. Adams   PetscFunctionBegin;
1080c20e4228SMark F. Adams   pc_gamg->m_threshold = n;
10813542efc5SMark F. Adams   PetscFunctionReturn(0);
10823542efc5SMark F. Adams }
10833542efc5SMark F. Adams EXTERN_C_END
10843542efc5SMark F. Adams 
10853542efc5SMark F. Adams #undef __FUNCT__
1086676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType"
1087676e1743SMark F. Adams /*@
1088676e1743SMark F. Adams    PCGAMGSetSolverType - Set solution method.
1089676e1743SMark F. Adams 
1090676e1743SMark F. Adams    Collective on PC
1091676e1743SMark F. Adams 
1092676e1743SMark F. Adams    Input Parameters:
1093676e1743SMark F. Adams .  pc - the preconditioner context
1094676e1743SMark F. Adams 
1095676e1743SMark F. Adams 
1096676e1743SMark F. Adams    Options Database Key:
10973542efc5SMark F. Adams .  -pc_gamg_type
1098676e1743SMark F. Adams 
1099676e1743SMark F. Adams    Level: intermediate
1100676e1743SMark F. Adams 
1101676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1102676e1743SMark F. Adams 
1103676e1743SMark F. Adams .seealso: ()
1104676e1743SMark F. Adams @*/
1105676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz )
1106676e1743SMark F. Adams {
1107676e1743SMark F. Adams   PetscErrorCode ierr;
1108676e1743SMark F. Adams 
1109676e1743SMark F. Adams   PetscFunctionBegin;
1110676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1111676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz));
1112676e1743SMark F. Adams   CHKERRQ(ierr);
1113676e1743SMark F. Adams   PetscFunctionReturn(0);
1114676e1743SMark F. Adams }
1115676e1743SMark F. Adams 
1116676e1743SMark F. Adams EXTERN_C_BEGIN
1117676e1743SMark F. Adams #undef __FUNCT__
1118676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG"
1119676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz )
1120676e1743SMark F. Adams {
1121676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1122676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1123676e1743SMark F. Adams 
1124676e1743SMark F. Adams   PetscFunctionBegin;
1125676e1743SMark F. Adams   if(sz < 64) strcpy(pc_gamg->m_type,str);
1126676e1743SMark F. Adams   PetscFunctionReturn(0);
1127676e1743SMark F. Adams }
1128676e1743SMark F. Adams EXTERN_C_END
1129676e1743SMark F. Adams 
11305b89ad90SMark F. Adams #undef __FUNCT__
11315b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
11325b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
11335b89ad90SMark F. Adams {
1134676e1743SMark F. Adams   PetscErrorCode  ierr;
1135676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1136676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1137676e1743SMark F. Adams   PetscBool flag;
11385b89ad90SMark F. Adams 
11395b89ad90SMark F. Adams   PetscFunctionBegin;
1140676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1141676e1743SMark F. Adams   {
114275b74bdaSMark F. Adams 
114375b74bdaSMark F. Adams     pc_gamg->m_method = 1; /* default to plane aggregation */
1144676e1743SMark F. Adams     ierr = PetscOptionsString("-pc_gamg_type",
114575b74bdaSMark F. Adams                               "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid ('geo')",
1146676e1743SMark F. Adams                               "PCGAMGSetSolverType",
114775b74bdaSMark F. Adams                               "pa",
1148676e1743SMark F. Adams                               pc_gamg->m_type,
1149676e1743SMark F. Adams                               64,
1150676e1743SMark F. Adams                               &flag );
1151676e1743SMark F. Adams     CHKERRQ(ierr);
1152d8c859f0SMark F. Adams     if( flag && pc_gamg->m_data != 0 ) {
1153d8c859f0SMark F. Adams       if( (strcmp(pc_gamg->m_type,"sa")==0 && pc_gamg->m_method != 2) ||
1154d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"pa")==0 && pc_gamg->m_method != 1) ||
1155d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"geo")==0 && pc_gamg->m_method != 0) ) {
1156d8c859f0SMark F. Adams         SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "PCSetFromOptions called of PCSetCoordinates (with new method, after data was created)");
1157d8c859f0SMark F. Adams       }
1158d8c859f0SMark F. Adams     }
1159bed7c2b9SMark F. Adams 
1160f6536408SMark F. Adams     if (flag ) {
1161f6536408SMark F. Adams       if( strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2;
1162f6536408SMark F. Adams       else if( strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1;
1163f6536408SMark F. Adams       else if( strcmp(pc_gamg->m_type,"geo") == 0) pc_gamg->m_method = 0;
1164f6536408SMark F. Adams       else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Invalid gamg type: %s",pc_gamg->m_type);
1165f6536408SMark F. Adams     }
116675b74bdaSMark F. Adams 
116775b74bdaSMark F. Adams     /* -pc_gamg_verbose */
116875b74bdaSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG","none",pc_gamg->m_verbose,&pc_gamg->m_verbose,PETSC_NULL);CHKERRQ(ierr);
116975b74bdaSMark F. Adams 
11708263b398SMark F. Adams     /* -pc_gamg_repartition */
11718263b398SMark F. Adams     pc_gamg->m_repart = PETSC_FALSE;
11728263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
11738263b398SMark F. Adams                             "Repartion coarse grids (false)",
11748263b398SMark F. Adams                             "PCGAMGRepartitioning",
11758263b398SMark F. Adams                             pc_gamg->m_repart,
11768263b398SMark F. Adams                             &pc_gamg->m_repart,
1177676e1743SMark F. Adams                             &flag);
1178676e1743SMark F. Adams     CHKERRQ(ierr);
1179676e1743SMark F. Adams 
1180c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1181c20e4228SMark F. Adams     pc_gamg->m_min_eq_proc = 600;
1182676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1183676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1184676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
1185c20e4228SMark F. Adams                            pc_gamg->m_min_eq_proc,
1186c20e4228SMark F. Adams                            &pc_gamg->m_min_eq_proc,
1187676e1743SMark F. Adams                            &flag );
1188676e1743SMark F. Adams     CHKERRQ(ierr);
11893542efc5SMark F. Adams 
1190389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1191389730f3SMark F. Adams     pc_gamg->m_coarse_eq_limit = 1500;
1192389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1193389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1194389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
1195389730f3SMark F. Adams                            pc_gamg->m_coarse_eq_limit,
1196389730f3SMark F. Adams                            &pc_gamg->m_coarse_eq_limit,
1197389730f3SMark F. Adams                            &flag );
1198389730f3SMark F. Adams     CHKERRQ(ierr);
1199389730f3SMark F. Adams 
1200c20e4228SMark F. Adams     /* -pc_gamg_threshold */
1201c20e4228SMark F. Adams     pc_gamg->m_threshold = 0.05;
12023542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
12033542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
12043542efc5SMark F. Adams                             "PCGAMGSetThreshold",
1205c20e4228SMark F. Adams                             pc_gamg->m_threshold,
1206c20e4228SMark F. Adams                             &pc_gamg->m_threshold,
12073542efc5SMark F. Adams                             &flag );
12083542efc5SMark F. Adams     CHKERRQ(ierr);
12094ef23d27SMark F. Adams 
12104ef23d27SMark F. Adams     pc_gamg->m_Nlevels = GAMG_MAXLEVELS;
12114ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
12124ef23d27SMark F. Adams                            "Set number of MG levels",
12134ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
12144ef23d27SMark F. Adams                            pc_gamg->m_Nlevels,
12154ef23d27SMark F. Adams                            &pc_gamg->m_Nlevels,
12164ef23d27SMark F. Adams                            &flag );
1217676e1743SMark F. Adams   }
1218676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1219676e1743SMark F. Adams 
12205b89ad90SMark F. Adams   PetscFunctionReturn(0);
12215b89ad90SMark F. Adams }
12225b89ad90SMark F. Adams 
12235b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
12245b89ad90SMark F. Adams /*MC
1225dc76db92SMark F. Adams      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning. This preconditioner currently has two
1226dc76db92SMark F. Adams            AMG methods: 1) an unstructured geometric method, which requires that you provide coordinates for each
1227dc76db92SMark F. Adams            vertex, and 2) smoothed aggregation.  Smoothed aggregation (SA) can work without coordinates but it
1228dc76db92SMark F. Adams            will generate some common non-trivial null spaces if coordinates are provided.  The input fine grid matrix
1229dc76db92SMark F. Adams            must have the block size set for 'system' problems (with multiple dof per vertex/cell) to work properly.
1230dc76db92SMark F. Adams            SA will generate rotational rigid body mode null space vectors, in addition to the trivial translational
1231dc76db92SMark F. Adams            modes, when coordinates are provide in 2D and 3D.
12325b89ad90SMark F. Adams 
1233280d9858SJed Brown    Options Database Keys:
12345b89ad90SMark F. Adams    Multigrid options(inherited)
1235280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1236280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1237280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1238280d9858SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
12395b89ad90SMark F. Adams 
12405b89ad90SMark F. Adams   Level: intermediate
1241280d9858SJed Brown 
12425b89ad90SMark F. Adams   Concepts: multigrid
12435b89ad90SMark F. Adams 
12445b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1245280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
12465b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
12475b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12485b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12495b89ad90SMark F. Adams M*/
12505b89ad90SMark F. Adams 
12515b89ad90SMark F. Adams EXTERN_C_BEGIN
12525b89ad90SMark F. Adams #undef __FUNCT__
12535b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
12545b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
12555b89ad90SMark F. Adams {
12565b89ad90SMark F. Adams   PetscErrorCode  ierr;
12575b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
12585b89ad90SMark F. Adams   PC_MG           *mg;
12595ef31b24SMark F. Adams   PetscClassId     cookie;
12605ee9c036SSatish Balay #if defined PETSC_USE_LOG
12615ee9c036SSatish Balay   static int count = 0;
12625ee9c036SSatish Balay #endif
12635b89ad90SMark F. Adams 
12645b89ad90SMark F. Adams   PetscFunctionBegin;
12655b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
12665b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
12675b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
12685b89ad90SMark F. Adams 
12695b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
12705b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
127103a628feSMark F. Adams   pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0;
12725b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
12735b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
12745b89ad90SMark F. Adams 
12755b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
12765b89ad90SMark F. Adams 
12775b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
12785b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
12795b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
12805b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
12815b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
12825b89ad90SMark F. Adams 
12835b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12845b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
12855b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
1286c97e1df0SMark F. Adams 					    PCSetCoordinates_GAMG);
1287c97e1df0SMark F. Adams   CHKERRQ(ierr);
1288c97e1df0SMark F. Adams 
1289676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1290676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1291676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1292676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1293676e1743SMark F. Adams   CHKERRQ(ierr);
1294676e1743SMark F. Adams 
1295676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1296389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_C",
1297389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_GAMG",
1298389730f3SMark F. Adams 					    PCGAMGSetCoarseEqLim_GAMG);
1299389730f3SMark F. Adams   CHKERRQ(ierr);
1300389730f3SMark F. Adams 
1301389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
13028263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_C",
13038263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_GAMG",
13048263b398SMark F. Adams 					    PCGAMGSetRepartitioning_GAMG);
1305676e1743SMark F. Adams   CHKERRQ(ierr);
1306676e1743SMark F. Adams 
1307676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
13083542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
13093542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
13103542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
13113542efc5SMark F. Adams   CHKERRQ(ierr);
13123542efc5SMark F. Adams 
13133542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1314676e1743SMark F. Adams 					    "PCGAMGSetSolverType_C",
1315676e1743SMark F. Adams 					    "PCGAMGSetSolverType_GAMG",
1316676e1743SMark F. Adams 					    PCGAMGSetSolverType_GAMG);
1317676e1743SMark F. Adams   CHKERRQ(ierr);
1318c97e1df0SMark F. Adams 
1319b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
1320785cba28SMark F. Adams   if( count++ == 0 ) {
13215ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
1322b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
13232a44bfbaSMark F. Adams     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
13242a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]);
13252a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]);
13262a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]);
1327b4fbaa2aSMark F. Adams     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1328b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1329b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1330b4fbaa2aSMark F. Adams     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1331b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_events[SET7]);
1332b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]); */
1333b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1334b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1335*ed3f9983SMark F. Adams     PetscLogEventRegister("  repartition", cookie, &gamg_setup_events[SET12]);
1336*ed3f9983SMark F. Adams     PetscLogEventRegister("  Invert-Sort", cookie, &gamg_setup_events[SET13]);
1337*ed3f9983SMark F. Adams     PetscLogEventRegister("  Move A", cookie, &gamg_setup_events[SET14]);
1338*ed3f9983SMark F. Adams     PetscLogEventRegister("  Move P", cookie, &gamg_setup_events[SET15]);
1339f852f58cSMark F. Adams 
1340b4fbaa2aSMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1341b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1342b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
13435ef31b24SMark F. Adams 
1344b4fbaa2aSMark F. Adams     /* create timer stages */
1345b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1346b4fbaa2aSMark F. Adams     {
1347b4fbaa2aSMark F. Adams       char str[32];
1348b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1349b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1350b4fbaa2aSMark F. Adams       PetscInt lidx;
1351b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1352b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1353b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1354b4fbaa2aSMark F. Adams       }
1355b4fbaa2aSMark F. Adams     }
1356b4fbaa2aSMark F. Adams #endif
1357b4fbaa2aSMark F. Adams   }
1358b4fbaa2aSMark F. Adams #endif
13595b89ad90SMark F. Adams   PetscFunctionReturn(0);
13605b89ad90SMark F. Adams }
13615b89ad90SMark F. Adams EXTERN_C_END
1362