xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision a2f3521de2a9d7869700f4c17e26c23fcfeaa6f6)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4b45d2f2cSJed Brown #include "petsc-private/matimpl.h"
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
6b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
7f96513f1SMatthew G Knepley 
80cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
90cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
10b4fbaa2aSMark F. Adams #endif
110cbbd2e1SMark F. Adams 
120cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
130cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_AGG;
140cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_GEO;
150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGOptprol_AGG;
20*a2f3521dSMark F. Adams PetscLogEvent PC_GAMGKKTProl_AGG;
210cbbd2e1SMark F. Adams #endif
220cbbd2e1SMark F. Adams 
23b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
24b4fbaa2aSMark F. Adams 
25b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
260cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
27b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
28b4fbaa2aSMark F. Adams #endif
29f96513f1SMatthew G Knepley 
309d5b6da9SMark F. Adams static PetscFList GAMGList = 0;
319d5b6da9SMark F. Adams 
32d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
33d3d6bff4SMark F. Adams #undef __FUNCT__
34d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
35d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
36d3d6bff4SMark F. Adams {
37d3d6bff4SMark F. Adams   PetscErrorCode  ierr;
38d3d6bff4SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
39d3d6bff4SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
40d3d6bff4SMark F. Adams 
41d3d6bff4SMark F. Adams   PetscFunctionBegin;
42*a2f3521dSMark F. Adams   if( pc_gamg->data ) { /* this should not happen, cleaned up in SetUp */
439d5b6da9SMark F. Adams     ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr);
4458471d46SMark F. Adams   }
45*a2f3521dSMark F. Adams   pc_gamg->data = PETSC_NULL; pc_gamg->data_sz = 0;
46*a2f3521dSMark F. Adams   PetscFunctionReturn(0);
47*a2f3521dSMark F. Adams }
48*a2f3521dSMark F. Adams 
49*a2f3521dSMark F. Adams /* private 2x2 Mat Nest for Stokes */
50*a2f3521dSMark F. Adams typedef struct{
51*a2f3521dSMark F. Adams   Mat A11,A21,A12,Amat;
52*a2f3521dSMark F. Adams   IS  prim_is,constr_is;
53*a2f3521dSMark F. Adams }GAMGKKTMat;
54*a2f3521dSMark F. Adams 
55*a2f3521dSMark F. Adams #undef __FUNCT__
56*a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatCreate"
57*a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatCreate( Mat A, PetscBool iskkt, GAMGKKTMat *out )
58*a2f3521dSMark F. Adams {
59*a2f3521dSMark F. Adams   PetscFunctionBegin;
60*a2f3521dSMark F. Adams   out->Amat = A;
61*a2f3521dSMark F. Adams   if( iskkt ){
62*a2f3521dSMark F. Adams     PetscErrorCode   ierr;
63*a2f3521dSMark F. Adams     IS       is_constraint, is_prime;
64*a2f3521dSMark F. Adams     PetscInt nmin,nmax;
65*a2f3521dSMark F. Adams 
66*a2f3521dSMark F. Adams     ierr = MatGetOwnershipRange( A, &nmin, &nmax );   CHKERRQ(ierr);
67*a2f3521dSMark F. Adams     ierr = MatFindZeroDiagonals( A, &is_constraint ); CHKERRQ(ierr);
68*a2f3521dSMark F. Adams     ierr = ISComplement( is_constraint, nmin, nmax, &is_prime ); CHKERRQ(ierr);
69*a2f3521dSMark F. Adams     out->prim_is = is_prime;
70*a2f3521dSMark F. Adams     out->constr_is = is_constraint;
71*a2f3521dSMark F. Adams 
72*a2f3521dSMark F. Adams     ierr = MatGetSubMatrix( A, is_prime, is_prime,      MAT_INITIAL_MATRIX, &out->A11); CHKERRQ(ierr);
73*a2f3521dSMark F. Adams     ierr = MatGetSubMatrix( A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12); CHKERRQ(ierr);
74*a2f3521dSMark F. Adams     ierr = MatGetSubMatrix( A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21); CHKERRQ(ierr);
75*a2f3521dSMark F. Adams   }
76*a2f3521dSMark F. Adams   else {
77*a2f3521dSMark F. Adams     out->A11 = A;
78*a2f3521dSMark F. Adams     out->A21 = PETSC_NULL;
79*a2f3521dSMark F. Adams     out->A12 = PETSC_NULL;
80*a2f3521dSMark F. Adams     out->prim_is = PETSC_NULL;
81*a2f3521dSMark F. Adams     out->constr_is = PETSC_NULL;
82*a2f3521dSMark F. Adams   }
83*a2f3521dSMark F. Adams   PetscFunctionReturn(0);
84*a2f3521dSMark F. Adams }
85*a2f3521dSMark F. Adams 
86*a2f3521dSMark F. Adams #undef __FUNCT__
87*a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatDestroy"
88*a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatDestroy( GAMGKKTMat *mat )
89*a2f3521dSMark F. Adams {
90*a2f3521dSMark F. Adams   PetscErrorCode   ierr;
91*a2f3521dSMark F. Adams 
92*a2f3521dSMark F. Adams   PetscFunctionBegin;
93*a2f3521dSMark F. Adams   if( mat->A11 && mat->A11 != mat->Amat ) {
94*a2f3521dSMark F. Adams     ierr = MatDestroy( &mat->A11 );  CHKERRQ(ierr);
95*a2f3521dSMark F. Adams   }
96*a2f3521dSMark F. Adams   ierr = MatDestroy( &mat->A21 );  CHKERRQ(ierr);
97*a2f3521dSMark F. Adams   ierr = MatDestroy( &mat->A12 );  CHKERRQ(ierr);
98*a2f3521dSMark F. Adams 
99*a2f3521dSMark F. Adams   ierr = ISDestroy( &mat->prim_is );    CHKERRQ(ierr);
100*a2f3521dSMark F. Adams   ierr = ISDestroy( &mat->constr_is );    CHKERRQ(ierr);
101*a2f3521dSMark F. Adams 
102d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
103d3d6bff4SMark F. Adams }
104d3d6bff4SMark F. Adams 
1055b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1065b89ad90SMark F. Adams /*
107a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
108a147abb0SMark F. Adams      of active processors.
1095b89ad90SMark F. Adams 
1105b89ad90SMark F. Adams    Input Parameter:
111*a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
112*a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
1139d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
1149d5b6da9SMark F. Adams    . cbs - coarse block size
1159d5b6da9SMark F. Adams    . isLast -
116*a2f3521dSMark F. Adams    . stokes -
1173530afc2SMark F. Adams    In/Output Parameter:
118*a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
119afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
12011e60469SMark F. Adams    Output Parameter:
1213530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1225b89ad90SMark F. Adams */
1235cb416c2SMark F. Adams 
1245b89ad90SMark F. Adams #undef __FUNCT__
1258263b398SMark F. Adams #define __FUNCT__ "createLevel"
1260cbbd2e1SMark F. Adams static PetscErrorCode createLevel( const PC pc,
1279d5b6da9SMark F. Adams                                    const Mat Amat_fine,
1289d5b6da9SMark F. Adams                                    const PetscInt cbs,
1299d5b6da9SMark F. Adams                                    const PetscBool isLast,
130*a2f3521dSMark F. Adams                                    const PetscBool stokes,
1313530afc2SMark F. Adams                                    Mat *a_P_inout,
132*a2f3521dSMark F. Adams                                    Mat *a_Amat_crs,
133*a2f3521dSMark F. Adams                                    PetscMPIInt *a_nactive_proc
13411e60469SMark F. Adams                                    )
1355b89ad90SMark F. Adams {
136*a2f3521dSMark F. Adams   PetscErrorCode   ierr;
1379d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
138486a8d0bSJed Brown   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1399d5b6da9SMark F. Adams   const PetscBool  repart = pc_gamg->repart;
1409d5b6da9SMark F. Adams   const PetscInt   min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
141*a2f3521dSMark F. Adams   Mat              Cmat,Pold=*a_P_inout;
1429d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)Amat_fine)->comm;
1435a9b9e01SMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive=*a_nactive_proc;
144*a2f3521dSMark F. Adams   PetscInt         ncrs_eq,ncrs_prim;
1455b89ad90SMark F. Adams 
1465b89ad90SMark F. Adams   PetscFunctionBegin;
14711e60469SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
14811e60469SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
14911e60469SMark F. Adams   /* RAP */
1509d5b6da9SMark F. Adams   ierr = MatPtAP( Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
151038e3b61SMark F. Adams 
152*a2f3521dSMark F. Adams   /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/
153*a2f3521dSMark F. Adams   ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
154*a2f3521dSMark F. Adams   assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)==0);
155*a2f3521dSMark F. Adams   ierr = MatGetLocalSize( Cmat, &ncrs_eq, PETSC_NULL );  CHKERRQ(ierr);
156*a2f3521dSMark F. Adams 
157*a2f3521dSMark F. Adams   /* get number of PEs to make active 'new_npe', reduce, can be any integer 1-P */
158*a2f3521dSMark F. Adams   {
159*a2f3521dSMark F. Adams     PetscInt ncrs_eq_glob,ncrs_eq_ave;
160*a2f3521dSMark F. Adams     ierr = MatGetSize( Cmat, &ncrs_eq_glob, PETSC_NULL );  CHKERRQ(ierr);
161*a2f3521dSMark F. Adams     ncrs_eq_ave = ncrs_eq_glob/npe;
162*a2f3521dSMark F. Adams     new_npe = (PetscMPIInt)((float)ncrs_eq_ave/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
163*a2f3521dSMark F. Adams     if( new_npe == 0 || ncrs_eq_ave < coarse_max ) new_npe = 1;
1645a9b9e01SMark F. Adams     else if ( new_npe >= nactive ) new_npe = nactive; /* no change, rare */
1659d5b6da9SMark F. Adams     if( isLast ) new_npe = 1;
166*a2f3521dSMark F. Adams   }
167f852f58cSMark F. Adams 
1685a9b9e01SMark F. Adams   if( !repart && new_npe==nactive ) {
169*a2f3521dSMark F. Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
17022063be5SMark F. Adams     }
17122063be5SMark F. Adams   else {
172*a2f3521dSMark F. Adams     const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
173*a2f3521dSMark F. Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
174*a2f3521dSMark F. Adams     IS              is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
1755a9b9e01SMark F. Adams     VecScatter      vecscat;
17622063be5SMark F. Adams     PetscScalar    *array;
17722063be5SMark F. Adams     Vec             src_crd, dest_crd;
178e33ef3b1SMark F. Adams 
179*a2f3521dSMark F. Adams     nloc_old = ncrs_eq/cbs; assert(ncrs_eq%cbs==0);
1800cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1810cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
182b4fbaa2aSMark F. Adams #endif
183*a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
184*a2f3521dSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
185*a2f3521dSMark F. Adams     if( repart && !stokes ) {
186*a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1875a9b9e01SMark F. Adams       Mat             adj;
1885a9b9e01SMark F. Adams 
189*a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
190*a2f3521dSMark F. Adams         if (pc_gamg->verbose==1) PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,ncrs_eq);
191*a2f3521dSMark F. Adams         else {
192*a2f3521dSMark F. Adams           PetscInt n;
193*a2f3521dSMark F. Adams           ierr = MPI_Allreduce( &ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm );CHKERRQ(ierr);
194*a2f3521dSMark F. Adams           PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,n);
195*a2f3521dSMark F. Adams         }
196*a2f3521dSMark F. Adams       }
1975a9b9e01SMark F. Adams 
198*a2f3521dSMark F. Adams       /* get 'adj' */
1999d5b6da9SMark F. Adams       if( cbs == 1 ) {
200038e3b61SMark F. Adams 	ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
201eb07cef2SMark F. Adams       }
202eb07cef2SMark F. Adams       else{
203*a2f3521dSMark F. Adams 	/* make a scalar matrix to partition (no Stokes here) */
204eb07cef2SMark F. Adams 	Mat tMat;
205*a2f3521dSMark F. Adams 	PetscInt Istart_crs,Iend_crs,ncols,jj,Ii;
206b4fbaa2aSMark F. Adams 	const PetscScalar *vals;
207b4fbaa2aSMark F. Adams 	const PetscInt *idx;
208*a2f3521dSMark F. Adams 	PetscInt *d_nnz, *o_nnz, M, N;
2099057884aSMark F. Adams 	static PetscInt llev = 0;
210b4fbaa2aSMark F. Adams 
211*a2f3521dSMark F. Adams 	ierr = PetscMalloc( ncrs_prim*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
212*a2f3521dSMark F. Adams 	ierr = PetscMalloc( ncrs_prim*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr);
213*a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange( Cmat, &Istart_crs, &Iend_crs );    CHKERRQ(ierr);
214*a2f3521dSMark F. Adams         ierr = MatGetSize( Cmat, &M, &N );CHKERRQ(ierr);
215*a2f3521dSMark F. Adams 	for ( Ii = Istart_crs, jj = 0 ; Ii < Iend_crs ; Ii += cbs, jj++ ) {
21658471d46SMark F. Adams 	  ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
2179d5b6da9SMark F. Adams 	  d_nnz[jj] = ncols/cbs;
2189d5b6da9SMark F. Adams 	  o_nnz[jj] = ncols/cbs;
21958471d46SMark F. Adams 	  ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
220*a2f3521dSMark F. Adams 	  if( d_nnz[jj] > ncrs_prim ) d_nnz[jj] = ncrs_prim;
221*a2f3521dSMark F. Adams 	  if( o_nnz[jj] > (M/cbs-ncrs_prim) ) o_nnz[jj] = M/cbs-ncrs_prim;
22258471d46SMark F. Adams 	}
2236876a03eSMark F. Adams 
224*a2f3521dSMark F. Adams 	ierr = MatCreate( wcomm, &tMat ); CHKERRQ(ierr);
225*a2f3521dSMark F. Adams 	ierr = MatSetSizes( tMat, ncrs_prim, ncrs_prim,
226*a2f3521dSMark F. Adams                             PETSC_DETERMINE, PETSC_DETERMINE );
2276876a03eSMark F. Adams         CHKERRQ(ierr);
228*a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);   CHKERRQ(ierr);
229*a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
230*a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
23158471d46SMark F. Adams 	ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
2325f8cf99dSMark F. Adams 	ierr = PetscFree( o_nnz ); CHKERRQ(ierr);
233eb07cef2SMark F. Adams 
234*a2f3521dSMark F. Adams 	for ( ii = Istart_crs; ii < Iend_crs; ii++ ) {
2359d5b6da9SMark F. Adams 	  PetscInt dest_row = ii/cbs;
23622063be5SMark F. Adams 	  ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
237eb07cef2SMark F. Adams 	  for( jj = 0 ; jj < ncols ; jj++ ){
2389d5b6da9SMark F. Adams 	    PetscInt dest_col = idx[jj]/cbs;
239eb07cef2SMark F. Adams 	    PetscScalar v = 1.0;
240eb07cef2SMark F. Adams 	    ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
241eb07cef2SMark F. Adams 	  }
24222063be5SMark F. Adams 	  ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
243eb07cef2SMark F. Adams 	}
244eb07cef2SMark F. Adams 	ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
245eb07cef2SMark F. Adams 	ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
246eb07cef2SMark F. Adams 
247b4fbaa2aSMark F. Adams 	if( llev++ == -1 ) {
248b4fbaa2aSMark F. Adams 	  PetscViewer viewer; char fname[32];
2494704b1adSJed Brown 	  ierr = PetscSNPrintf(fname,sizeof fname,"part_mat_%D.mat",llev);CHKERRQ(ierr);
250b4fbaa2aSMark F. Adams 	  PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
251b4fbaa2aSMark F. Adams 	  ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
252b4fbaa2aSMark F. Adams 	  ierr = PetscViewerDestroy( &viewer );
253b4fbaa2aSMark F. Adams 	}
254b4fbaa2aSMark F. Adams 
255eb07cef2SMark F. Adams 	ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
256eb07cef2SMark F. Adams 
257eb07cef2SMark F. Adams 	ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
258*a2f3521dSMark F. Adams       } /* create 'adj' */
259f150b916SMark F. Adams 
260*a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2615a9b9e01SMark F. Adams 	char prefix[256];
2625a9b9e01SMark F. Adams 	const char *pcpre;
263b4fbaa2aSMark F. Adams 	const PetscInt *is_idx;
264b4fbaa2aSMark F. Adams 	MatPartitioning  mpart;
265a4b7d37bSMark F. Adams 	IS proc_is;
266*a2f3521dSMark F. Adams 	PetscInt targetPE;
2672f03bc48SMark F. Adams 
2685a9b9e01SMark F. Adams 	ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr);
2695ef31b24SMark F. Adams 	ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
2709d5b6da9SMark F. Adams 	ierr = PCGetOptionsPrefix( pc, &pcpre );CHKERRQ(ierr);
27159a0be82SJed Brown 	ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
27259a0be82SJed Brown 	ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
27311e60469SMark F. Adams 	ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
2745ef31b24SMark F. Adams 	ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
275a4b7d37bSMark F. Adams 	ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr);
27611e60469SMark F. Adams 	ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
2775a9b9e01SMark F. Adams 
2785ef31b24SMark F. Adams 	/* collect IS info */
279*a2f3521dSMark F. Adams 	ierr = PetscMalloc( ncrs_eq*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr);
280a4b7d37bSMark F. Adams 	ierr = ISGetIndices( proc_is, &is_idx );        CHKERRQ(ierr);
281*a2f3521dSMark F. Adams 	targetPE = 1; /* bring to "front" of machine */
282*a2f3521dSMark F. Adams 	/*targetPE = npe/new_npe;*/ /* spread partitioning across machine */
283*a2f3521dSMark F. Adams 	for( kk = jj = 0 ; kk < nloc_old ; kk++ ){
2849d5b6da9SMark F. Adams 	  for( ii = 0 ; ii < cbs ; ii++, jj++ ){
285*a2f3521dSMark F. Adams 	    newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
286eb07cef2SMark F. Adams 	  }
2875ef31b24SMark F. Adams 	}
288a4b7d37bSMark F. Adams 	ierr = ISRestoreIndices( proc_is, &is_idx );     CHKERRQ(ierr);
289a4b7d37bSMark F. Adams 	ierr = ISDestroy( &proc_is );                    CHKERRQ(ierr);
2905ef31b24SMark F. Adams       }
2915ef31b24SMark F. Adams       ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
2925a9b9e01SMark F. Adams 
293*a2f3521dSMark F. Adams       ierr = ISCreateGeneral( wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc );
2945a9b9e01SMark F. Adams       CHKERRQ(ierr);
2958263b398SMark F. Adams       if( newproc_idx != 0 ) {
2968263b398SMark F. Adams 	ierr = PetscFree( newproc_idx );  CHKERRQ(ierr);
2975ef31b24SMark F. Adams       }
298*a2f3521dSMark F. Adams     } /* repartitioning */
299*a2f3521dSMark F. Adams     else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
300*a2f3521dSMark F. Adams 
301*a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
3025a9b9e01SMark F. Adams       /* find factor */
3035a9b9e01SMark F. Adams       if( new_npe == 1 ) rfactor = npe; /* easy */
3045a9b9e01SMark F. Adams       else {
3055a9b9e01SMark F. Adams 	PetscReal best_fact = 0.;
3065a9b9e01SMark F. Adams 	jj = -1;
3075a9b9e01SMark F. Adams 	for( kk = 1 ; kk <= npe ; kk++ ){
3085a9b9e01SMark F. Adams 	  if( npe%kk==0 ) { /* a candidate */
3095a9b9e01SMark F. Adams 	    PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe;
3105a9b9e01SMark F. Adams 	    if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */
3115a9b9e01SMark F. Adams 	    if( fact > best_fact ) {
3125a9b9e01SMark F. Adams 	      best_fact = fact; jj = kk;
3135a9b9e01SMark F. Adams 	    }
3145a9b9e01SMark F. Adams 	  }
3155a9b9e01SMark F. Adams 	}
3165a9b9e01SMark F. Adams 	if( jj != -1 ) rfactor = jj;
317*a2f3521dSMark F. Adams 	else rfactor = 1; /* does this happen .. a prime */
3185a9b9e01SMark F. Adams       }
3195a9b9e01SMark F. Adams       new_npe = npe/rfactor;
3205a9b9e01SMark F. Adams 
3215a9b9e01SMark F. Adams       if( new_npe==nactive ) {
322*a2f3521dSMark F. Adams 	*a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
3235a9b9e01SMark F. Adams 	ierr = PetscFree( counts );  CHKERRQ(ierr);
324*a2f3521dSMark F. Adams 	if (pc_gamg->verbose>0){
325*a2f3521dSMark F. Adams           PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq(loc)=%d\n",mype,__FUNCT__,new_npe,ncrs_eq);
326*a2f3521dSMark F. Adams         }
3275a9b9e01SMark F. Adams 	PetscFunctionReturn(0);
3285a9b9e01SMark F. Adams       }
3295a9b9e01SMark F. Adams 
330*a2f3521dSMark F. Adams       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",mype,__FUNCT__,ncrs_eq);
3315a9b9e01SMark F. Adams       targetPE = mype/rfactor;
332*a2f3521dSMark F. Adams       ierr = ISCreateStride( wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc ); CHKERRQ(ierr);
3335a9b9e01SMark F. Adams 
334*a2f3521dSMark F. Adams       if( stokes ) {
335*a2f3521dSMark F. Adams         ierr = ISCreateStride( wcomm, ncrs_prim*cbs, targetPE, 0, &is_eq_newproc_prim ); CHKERRQ(ierr);
3365a9b9e01SMark F. Adams       }
337*a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
338e33ef3b1SMark F. Adams 
33911e60469SMark F. Adams     /*
340*a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
34111e60469SMark F. Adams      */
342*a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering( is_eq_newproc, &is_eq_num ); CHKERRQ(ierr);
343*a2f3521dSMark F. Adams     if( stokes ) {
344*a2f3521dSMark F. Adams       ierr = ISPartitioningToNumbering( is_eq_newproc_prim, &is_eq_num_prim ); CHKERRQ(ierr);
345*a2f3521dSMark F. Adams     }
346*a2f3521dSMark F. Adams     else is_eq_num_prim = is_eq_num;
34711e60469SMark F. Adams     /*
348*a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
34911e60469SMark F. Adams      */
350*a2f3521dSMark F. Adams     ierr = ISPartitioningCount( is_eq_newproc, npe, counts ); CHKERRQ(ierr);
351*a2f3521dSMark F. Adams     ncrs_eq_new = counts[mype];
352*a2f3521dSMark F. Adams     ierr = ISDestroy( &is_eq_newproc );                       CHKERRQ(ierr);
353*a2f3521dSMark F. Adams     if( stokes ) {
354*a2f3521dSMark F. Adams       ierr = ISPartitioningCount( is_eq_newproc_prim, npe, counts ); CHKERRQ(ierr);
355*a2f3521dSMark F. Adams       ierr = ISDestroy( &is_eq_newproc_prim );                       CHKERRQ(ierr);
356*a2f3521dSMark F. Adams       ncrs_prim_new = counts[mype]/cbs; /* nodes */
357*a2f3521dSMark F. Adams     }
358*a2f3521dSMark F. Adams     else ncrs_prim_new = ncrs_eq_new/cbs; /* eqs */
359*a2f3521dSMark F. Adams 
360*a2f3521dSMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
3610cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3620cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
363b4fbaa2aSMark F. Adams #endif
364*a2f3521dSMark F. Adams 
365*a2f3521dSMark F. Adams     /* move data (for primal equations only) */
36622063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
36711e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
368*a2f3521dSMark F. Adams     ierr = VecSetSizes( dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE ); CHKERRQ(ierr);
369470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
37011e60469SMark F. Adams     /*
3719d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
3729d5b6da9SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cbs'.
37311e60469SMark F. Adams      */
374*a2f3521dSMark F. Adams     ierr = PetscMalloc( (ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
375*a2f3521dSMark F. Adams     ierr = ISGetIndices( is_eq_num_prim, &idx ); CHKERRQ(ierr);
376*a2f3521dSMark F. Adams     for(ii=0,jj=0; ii<ncrs_prim ; ii++) {
3779d5b6da9SMark F. Adams       PetscInt id = idx[ii*cbs]/cbs; /* get node back */
378*a2f3521dSMark F. Adams       for( kk=0; kk<node_data_sz ; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
37911e60469SMark F. Adams     }
380*a2f3521dSMark F. Adams     ierr = ISRestoreIndices( is_eq_num_prim, &idx ); CHKERRQ(ierr);
381*a2f3521dSMark F. Adams     ierr = ISCreateGeneral( wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat );
38211e60469SMark F. Adams     CHKERRQ(ierr);
38392a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
38411e60469SMark F. Adams     /*
38511e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
38611e60469SMark F. Adams      */
387*a2f3521dSMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd ); CHKERRQ(ierr);
3889d5b6da9SMark F. Adams     for( jj=0; jj<ndata_cols ; jj++ ) {
389*a2f3521dSMark F. Adams       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
390*a2f3521dSMark F. Adams       for( ii=0 ; ii<ncrs_prim ; ii++) {
3919d5b6da9SMark F. Adams 	for( kk=0; kk<ndata_rows ; kk++ ) {
392*a2f3521dSMark F. Adams 	  PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
393c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
394676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
395d3d6bff4SMark F. Adams 	}
396038e3b61SMark F. Adams       }
397eb07cef2SMark F. Adams     }
398eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
399eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
40011e60469SMark F. Adams     /*
40111e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
40211e60469SMark F. Adams       to the correct processor
40311e60469SMark F. Adams     */
40411e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
40511e60469SMark F. Adams     CHKERRQ(ierr);
40611e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
40711e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40811e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40911e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
41011e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
41111e60469SMark F. Adams     /*
41211e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
41311e60469SMark F. Adams     */
414c8b0795cSMark F. Adams     ierr = PetscFree( pc_gamg->data );    CHKERRQ(ierr);
415*a2f3521dSMark F. Adams     ierr = PetscMalloc( node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data );    CHKERRQ(ierr);
416*a2f3521dSMark F. Adams     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
417*a2f3521dSMark F. Adams     strideNew = ncrs_prim_new*ndata_rows;
41811e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
4199d5b6da9SMark F. Adams     for( jj=0; jj<ndata_cols ; jj++ ) {
420*a2f3521dSMark F. Adams       for( ii=0 ; ii<ncrs_prim_new ; ii++) {
4219d5b6da9SMark F. Adams 	for( kk=0; kk<ndata_rows ; kk++ ) {
422*a2f3521dSMark F. Adams 	  PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
423c8b0795cSMark F. Adams 	  pc_gamg->data[ix] = PetscRealPart(array[jx]);
424d3d6bff4SMark F. Adams 	}
425038e3b61SMark F. Adams       }
426038e3b61SMark F. Adams     }
42711e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
42811e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
429*a2f3521dSMark F. Adams 
430*a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
4310cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4320cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
433ed3f9983SMark F. Adams #endif
434*a2f3521dSMark F. Adams 
43511e60469SMark F. Adams     /*
43611e60469SMark F. Adams       Invert for MatGetSubMatrix
43711e60469SMark F. Adams     */
438*a2f3521dSMark F. Adams     ierr = ISInvertPermutation( is_eq_num, ncrs_eq_new, &new_eq_indices ); CHKERRQ(ierr);
439*a2f3521dSMark F. Adams     ierr = ISSort( new_eq_indices ); CHKERRQ(ierr); /* is this needed? */
440*a2f3521dSMark F. Adams     if(is_eq_num != is_eq_num_prim) {
441*a2f3521dSMark F. Adams       ierr = ISDestroy( &is_eq_num_prim ); CHKERRQ(ierr); /* could be same as 'is_eq_num' */
442*a2f3521dSMark F. Adams     }
443*a2f3521dSMark F. Adams     ierr = ISDestroy( &is_eq_num ); CHKERRQ(ierr);
4440cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4450cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
4460cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
447ed3f9983SMark F. Adams #endif
448*a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
449*a2f3521dSMark F. Adams     {
450*a2f3521dSMark F. Adams       Mat mat;
451*a2f3521dSMark F. Adams       ierr = MatGetSubMatrix( Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat );
45211e60469SMark F. Adams       CHKERRQ(ierr);
453*a2f3521dSMark F. Adams       *a_Amat_crs = mat;
454*a2f3521dSMark F. Adams     }
455038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
456*a2f3521dSMark F. Adams 
4570cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4580cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
459ed3f9983SMark F. Adams #endif
46011e60469SMark F. Adams     /* prolongator */
46111e60469SMark F. Adams     {
46211e60469SMark F. Adams       IS findices;
463*a2f3521dSMark F. Adams       PetscInt Istart,Iend;
464*a2f3521dSMark F. Adams       Mat Pnew;
465*a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
4660cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4670cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
468ed3f9983SMark F. Adams #endif
46911e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
470*a2f3521dSMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew );
47111e60469SMark F. Adams       CHKERRQ(ierr);
47211e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
4730cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4740cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
475ed3f9983SMark F. Adams #endif
4763530afc2SMark F. Adams       ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
477*a2f3521dSMark F. Adams 
478*a2f3521dSMark F. Adams       /* output - repartitioned */
479*a2f3521dSMark F. Adams       *a_P_inout = Pnew;
480e33ef3b1SMark F. Adams     }
481*a2f3521dSMark F. Adams     ierr = ISDestroy( &new_eq_indices ); CHKERRQ(ierr);
4825b89ad90SMark F. Adams 
4835a9b9e01SMark F. Adams     *a_nactive_proc = new_npe; /* output */
484*a2f3521dSMark F. Adams   }
4855a9b9e01SMark F. Adams 
486*a2f3521dSMark F. Adams   /* outout matrix data */
487c8b0795cSMark F. Adams   if( !PETSC_TRUE ) {
488c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
489c8b0795cSMark F. Adams     if(llev==0) {
490c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
491c8b0795cSMark F. Adams       PetscViewerASCIIOpen(wcomm,fname,&viewer);
492c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
493c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer ); CHKERRQ(ierr);
494c8b0795cSMark F. Adams       ierr = PetscViewerDestroy( &viewer );
495c8b0795cSMark F. Adams     }
496c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
497c8b0795cSMark F. Adams     PetscViewerASCIIOpen(wcomm,fname,&viewer);
498c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
499c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer ); CHKERRQ(ierr);
500c8b0795cSMark F. Adams     ierr = PetscViewerDestroy( &viewer );
501c8b0795cSMark F. Adams   }
502c8b0795cSMark F. Adams 
5035b89ad90SMark F. Adams   PetscFunctionReturn(0);
5045b89ad90SMark F. Adams }
5055b89ad90SMark F. Adams 
5065b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5075b89ad90SMark F. Adams /*
5085b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5095b89ad90SMark F. Adams                     by setting data structures and options.
5105b89ad90SMark F. Adams 
5115b89ad90SMark F. Adams    Input Parameter:
5125b89ad90SMark F. Adams .  pc - the preconditioner context
5135b89ad90SMark F. Adams 
5145b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
5155b89ad90SMark F. Adams 
5165b89ad90SMark F. Adams    Notes:
5175b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
5185b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
5195b89ad90SMark F. Adams */
5205b89ad90SMark F. Adams #undef __FUNCT__
5215b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
5229d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC pc )
5235b89ad90SMark F. Adams {
5245b89ad90SMark F. Adams   PetscErrorCode  ierr;
5259d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
5265b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
5272adcac29SMark F. Adams   Mat              Pmat = pc->pmat;
528*a2f3521dSMark F. Adams   PetscInt         fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
5299d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
5303530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
531587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
532c8b0795cSMark F. Adams   PetscReal        emaxs[GAMG_MAXLEVELS];
533a94c3b12SMark F. Adams   IS              *ASMLocalIDsArr[GAMG_MAXLEVELS],removedEqs[GAMG_MAXLEVELS];
534a94c3b12SMark F. Adams   PetscInt         level_bs[GAMG_MAXLEVELS];
535*a2f3521dSMark F. Adams   GAMGKKTMat       kktMatsArr[GAMG_MAXLEVELS];
536*a2f3521dSMark F. Adams   PetscLogDouble   nnz0=0.,nnztot=0.;
537737a81a9SMark F. Adams   MatInfo          info;
538302f38e8SMark F. Adams   PetscBool        stokes = PETSC_FALSE;
5395ef31b24SMark F. Adams 
5405b89ad90SMark F. Adams   PetscFunctionBegin;
54184d3f75bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
54284d3f75bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
543*a2f3521dSMark F. Adams   if (pc_gamg->verbose>2) PetscPrintf(wcomm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",mype,__FUNCT__,pc_gamg->setup_count,pc->setupcalled);
54484d3f75bSMark F. Adams   if( pc_gamg->setup_count++ > 0 ) {
54584d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
54603a628feSMark F. Adams     /* just do Galerkin grids */
54758471d46SMark F. Adams     Mat B,dA,dB;
548d5280255SMark F. Adams     assert(pc->setupcalled);
54958471d46SMark F. Adams 
5509d5b6da9SMark F. Adams     if( pc_gamg->Nlevels > 1 ) {
55158471d46SMark F. Adams       /* currently only handle case where mat and pmat are the same on coarser levels */
5529d5b6da9SMark F. Adams       ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
55358471d46SMark F. Adams       /* (re)set to get dirty flag */
5549d5b6da9SMark F. Adams       ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
55558471d46SMark F. Adams 
5569d5b6da9SMark F. Adams       for (level=pc_gamg->Nlevels-2; level>-1; level--) {
55703a628feSMark F. Adams         /* the first time through the matrix structure has changed from repartitioning */
55884d3f75bSMark F. Adams         if( pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
55903a628feSMark F. Adams           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
560084a8fe3SJed Brown           ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
56103a628feSMark F. Adams           mglevels[level]->A = B;
56203a628feSMark F. Adams         }
56303a628feSMark F. Adams         else {
564084a8fe3SJed Brown           ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
56558471d46SMark F. Adams           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
56603a628feSMark F. Adams         }
56758471d46SMark F. Adams         ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
56858471d46SMark F. Adams         dB = B;
56958471d46SMark F. Adams       }
5705f8cf99dSMark F. Adams     }
571d5280255SMark F. Adams 
572d5280255SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
573d5280255SMark F. Adams 
574d5280255SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
575d5280255SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
576d5280255SMark F. Adams 
57758471d46SMark F. Adams     PetscFunctionReturn(0);
578eb07cef2SMark F. Adams   }
57984d3f75bSMark F. Adams   assert(pc->setupcalled == 0);
580f6536408SMark F. Adams 
581302f38e8SMark F. Adams   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
582eb07cef2SMark F. Adams 
583302f38e8SMark F. Adams   if( pc_gamg->data == 0 ) {
5849d5b6da9SMark F. Adams     if( !pc_gamg->createdefaultdata ){
585*a2f3521dSMark F. Adams       SETERRQ(wcomm,PETSC_ERR_LIB,"'createdefaultdata' not set(?) need to support NULL data");
586302f38e8SMark F. Adams     }
587302f38e8SMark F. Adams     if( stokes ) {
588302f38e8SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_LIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
589eb07cef2SMark F. Adams     }
590b8cd405aSMark F. Adams     ierr = pc_gamg->createdefaultdata( pc, Pmat ); CHKERRQ(ierr);
5919d5b6da9SMark F. Adams   }
592038e3b61SMark F. Adams 
593302f38e8SMark F. Adams   /* get basic dims */
594302f38e8SMark F. Adams   if( stokes ) {
595*a2f3521dSMark F. Adams     bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
596*a2f3521dSMark F. Adams     /* nloc = pc_gamg->data_sz/pc_gamg->data_cell_cols/bs; assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*bs)==0); */
597302f38e8SMark F. Adams   }
598302f38e8SMark F. Adams   else {
599*a2f3521dSMark F. Adams     /* PetscInt Istart,Iend,M,N; */
600302f38e8SMark F. Adams     ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr);
601*a2f3521dSMark F. Adams     /* ierr = MatGetSize( Pmat, &M, &N );CHKERRQ(ierr); */
602*a2f3521dSMark F. Adams     /* ierr = MatGetOwnershipRange( Pmat, &Istart, &Iend ); CHKERRQ(ierr); */
603*a2f3521dSMark F. Adams     /* nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); */
604302f38e8SMark F. Adams   }
605*a2f3521dSMark F. Adams 
606*a2f3521dSMark F. Adams   ierr = MatGetSize( Pmat, &M, &qq );CHKERRQ(ierr);
607c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
6082adcac29SMark F. Adams     if(pc_gamg->verbose==1) ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);
6092adcac29SMark F. Adams     else ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);
610b2a4f308SMark F. Adams     CHKERRQ(ierr);
611b2a4f308SMark F. Adams     nnz0 = info.nz_used;
612b2a4f308SMark F. Adams     nnztot = info.nz_used;
613b43b03e9SMark 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",
614*a2f3521dSMark F. Adams                 mype,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
615*a2f3521dSMark F. Adams                 (int)(nnz0/(PetscReal)M),npe);
616c8b0795cSMark F. Adams   }
61784d3f75bSMark F. Adams 
618*a2f3521dSMark F. Adams   /* Get A_i and R_i */
6198f4b7eb5SMark F. Adams   for ( level=0, Aarr[0]=Pmat, nactivepe = npe; /* hard wired stopping logic */
6209d5b6da9SMark F. Adams         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
6210205a208SMark F. Adams         level++ ){
6225b89ad90SMark F. Adams     level1 = level + 1;
6230cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6240cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr);
625*a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
626*a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
627b4fbaa2aSMark F. Adams #endif
628*a2f3521dSMark F. Adams #endif
629*a2f3521dSMark F. Adams     /* deal with Stokes, get sub matrices */
630*a2f3521dSMark F. Adams     ierr = GAMGKKTMatCreate( Aarr[level], stokes, &kktMatsArr[level] ); CHKERRQ(ierr);
631*a2f3521dSMark F. Adams 
632c8b0795cSMark F. Adams     { /* construct prolongator */
633725b86d8SJed Brown       Mat Gmat;
6340cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
635*a2f3521dSMark F. Adams       Mat Prol11,Prol22;
636c8b0795cSMark F. Adams 
637a94c3b12SMark F. Adams       level_bs[level] = bs;
638*a2f3521dSMark F. Adams       ierr = pc_gamg->graph( pc,kktMatsArr[level].A11, &Gmat ); CHKERRQ(ierr);
6390cbbd2e1SMark F. Adams       ierr = pc_gamg->coarsen( pc, &Gmat, &agg_lists ); CHKERRQ(ierr);
640*a2f3521dSMark F. Adams       ierr = pc_gamg->prolongator( pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11 ); CHKERRQ(ierr);
641c8b0795cSMark F. Adams 
642*a2f3521dSMark F. Adams       /* could have failed to create new level */
643*a2f3521dSMark F. Adams       if( Prol11 ){
6449d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
645c8b0795cSMark F. Adams         if( pc_gamg->col_bs_id != -1 ){
6469d5b6da9SMark F. Adams           PetscBool flag;
647*a2f3521dSMark F. Adams           ierr = PetscObjectComposedDataGetInt( (PetscObject)Prol11, pc_gamg->col_bs_id, bs, flag );
648c8b0795cSMark F. Adams           CHKERRQ( ierr ); assert(flag);
6499d5b6da9SMark F. Adams         }
650*a2f3521dSMark F. Adams 
651*a2f3521dSMark F. Adams         if( stokes ) {
652*a2f3521dSMark F. Adams           if(!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
653*a2f3521dSMark F. Adams           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
654*a2f3521dSMark F. Adams           ierr = pc_gamg->formkktprol( pc, Prol11, kktMatsArr[level].A21, &Prol22 ); CHKERRQ(ierr);
655c8b0795cSMark F. Adams         }
656c8b0795cSMark F. Adams 
657*a2f3521dSMark F. Adams         if( pc_gamg->optprol ){
658c8b0795cSMark F. Adams           /* smooth */
659*a2f3521dSMark F. Adams           ierr = pc_gamg->optprol( pc, kktMatsArr[level].A11, &Prol11 ); CHKERRQ(ierr);
660c8b0795cSMark F. Adams         }
661c8b0795cSMark F. Adams 
662*a2f3521dSMark F. Adams         if( stokes ) {
663*a2f3521dSMark F. Adams           IS is_row[2] = {kktMatsArr[level].prim_is,kktMatsArr[level].constr_is};
664*a2f3521dSMark F. Adams           Mat a[4] = {Prol11, PETSC_NULL, PETSC_NULL, Prol22 };
665*a2f3521dSMark F. Adams           ierr = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1] ); CHKERRQ(ierr);
666*a2f3521dSMark F. Adams         }
667*a2f3521dSMark F. Adams         else {
668*a2f3521dSMark F. Adams           Parr[level1] = Prol11;
669*a2f3521dSMark F. Adams         }
670*a2f3521dSMark F. Adams       }
671*a2f3521dSMark F. Adams       else Parr[level1] = PETSC_NULL;
672ffc955d6SMark F. Adams 
673ffc955d6SMark F. Adams       if ( pc_gamg->use_aggs_in_gasm ) {
674a94c3b12SMark F. Adams         ierr = PetscCDGetASMBlocks(agg_lists, level_bs[level], &nASMBlocksArr[level], &ASMLocalIDsArr[level] );  CHKERRQ(ierr);
675ffc955d6SMark F. Adams       }
676ffc955d6SMark F. Adams 
677a94c3b12SMark F. Adams       ierr = PetscCDGetRemovedIS( agg_lists, &removedEqs[level] );  CHKERRQ(ierr);
678a94c3b12SMark F. Adams 
679*a2f3521dSMark F. Adams       ierr = MatDestroy( &Gmat );      CHKERRQ(ierr);
68041b27cdeSMark F. Adams       ierr = PetscCDDestroy( agg_lists );  CHKERRQ(ierr);
681*a2f3521dSMark F. Adams     } /* construct prolongator scope */
6820cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6830cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
684c8b0795cSMark F. Adams #endif
6859d5b6da9SMark F. Adams     /* cache eigen estimate */
6869d5b6da9SMark F. Adams     if( pc_gamg->emax_id != -1 ){
6879d5b6da9SMark F. Adams       PetscBool flag;
688*a2f3521dSMark F. Adams       ierr = PetscObjectComposedDataGetReal( (PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag );
6899d5b6da9SMark F. Adams       CHKERRQ( ierr );
6909d5b6da9SMark F. Adams       if( !flag ) emaxs[level] = -1.;
6919d5b6da9SMark F. Adams     }
6929d5b6da9SMark F. Adams     else emaxs[level] = -1.;
6932adcac29SMark F. Adams     if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
694c8b0795cSMark F. Adams     if( !Parr[level1] ) {
695b43b03e9SMark F. Adams       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);
696c8b0795cSMark F. Adams       break;
697c8b0795cSMark F. Adams     }
6980cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6990cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
700b4fbaa2aSMark F. Adams #endif
701*a2f3521dSMark F. Adams 
702*a2f3521dSMark F. Adams     ierr = createLevel( pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
703*a2f3521dSMark F. Adams                         stokes, &Parr[level1], &Aarr[level1], &nactivepe );
7043530afc2SMark F. Adams     CHKERRQ(ierr);
705*a2f3521dSMark F. Adams 
7060cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7070cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
708b4fbaa2aSMark F. Adams #endif
709*a2f3521dSMark F. Adams     ierr = MatGetSize( Aarr[level1], &M, &qq );CHKERRQ(ierr);
710*a2f3521dSMark F. Adams 
711*a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0){
7120cbbd2e1SMark F. Adams       PetscInt NN = M;
7130cbbd2e1SMark F. Adams       if(pc_gamg->verbose==1) {
714*a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info); CHKERRQ(ierr);
715*a2f3521dSMark F. Adams         ierr = MatGetLocalSize( Aarr[level1], &NN, &qq );
7160cbbd2e1SMark F. Adams       }
7170cbbd2e1SMark F. Adams       else ierr = MatGetInfo( Aarr[level1], MAT_GLOBAL_SUM, &info );
718*a2f3521dSMark F. Adams 
7190cbbd2e1SMark F. Adams       CHKERRQ(ierr);
7200cbbd2e1SMark F. Adams       nnztot += info.nz_used;
721b43b03e9SMark F. Adams       PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
722*a2f3521dSMark F. Adams                   mype,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
7230cbbd2e1SMark F. Adams                   (int)(info.nz_used/(PetscReal)NN), nactivepe );
7240cbbd2e1SMark F. Adams       CHKERRQ(ierr);
725c8b0795cSMark F. Adams     }
726*a2f3521dSMark F. Adams 
727*a2f3521dSMark F. Adams     /* stop if one node -- could pull back for singular problems */
728c8b0795cSMark F. Adams     if( M/pc_gamg->data_cell_cols < 2 ) {
72981708dccSMark F. Adams       level++;
73081708dccSMark F. Adams       break;
73181708dccSMark F. Adams     }
7320cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
733b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
734b4fbaa2aSMark F. Adams #endif
735c8b0795cSMark F. Adams   } /* levels */
736c8b0795cSMark F. Adams 
737c8b0795cSMark F. Adams   if( pc_gamg->data ) {
738c8b0795cSMark F. Adams     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
739*a2f3521dSMark F. Adams     pc_gamg->data = PETSC_NULL;
7405b89ad90SMark F. Adams   }
741c8b0795cSMark F. Adams 
7420cbbd2e1SMark F. Adams   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
7439d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7445b89ad90SMark F. Adams   fine_level = level;
7459d5b6da9SMark F. Adams   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
7465b89ad90SMark F. Adams 
74784d3f75bSMark F. Adams   /* simple setup */
74884d3f75bSMark F. Adams   if( !PETSC_TRUE ){
74984d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
75084d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
75184d3f75bSMark F. Adams          lidx<fine_level;
75284d3f75bSMark F. Adams          lidx++, level--){
75384d3f75bSMark F. Adams       ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr);
75484d3f75bSMark F. Adams       ierr = KSPSetOperators( mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );CHKERRQ(ierr);
75584d3f75bSMark F. Adams       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
75684d3f75bSMark F. Adams       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
75784d3f75bSMark F. Adams     }
75884d3f75bSMark F. Adams     ierr = KSPSetOperators( mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
75984d3f75bSMark F. Adams 
76084d3f75bSMark F. Adams     ierr = PCSetUp_MG( pc );  CHKERRQ( ierr );
76184d3f75bSMark F. Adams   }
76284d3f75bSMark F. Adams   else if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if one level */
763d5280255SMark F. Adams     /* set default smoothers & set operators */
7649d5b6da9SMark F. Adams     for ( lidx = 1, level = pc_gamg->Nlevels-2;
765587fa25dSMark F. Adams           lidx <= fine_level;
766587fa25dSMark F. Adams           lidx++, level--) {
767ffc955d6SMark F. Adams       KSP smoother;
768ffc955d6SMark F. Adams       PC subpc;
769*a2f3521dSMark F. Adams 
7709d5b6da9SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
771f6536408SMark F. Adams       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
772ffc955d6SMark F. Adams 
773*a2f3521dSMark F. Adams       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
774*a2f3521dSMark F. Adams       /* set ops */
775*a2f3521dSMark F. Adams       ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
776*a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation( pc, lidx, Parr[level+1] );CHKERRQ(ierr);
777*a2f3521dSMark F. Adams 
778*a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
779*a2f3521dSMark F. Adams       if( stokes ) {
780*a2f3521dSMark F. Adams         KSP *ksps;
781*a2f3521dSMark F. Adams         PetscInt nn;
782*a2f3521dSMark F. Adams         ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);   CHKERRQ(ierr);
783*a2f3521dSMark F. Adams         ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is); CHKERRQ(ierr);
784*a2f3521dSMark F. Adams         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps); CHKERRQ(ierr);
785*a2f3521dSMark F. Adams         smoother = ksps[0];
786*a2f3521dSMark F. Adams         ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
787*a2f3521dSMark F. Adams         ierr = PetscFree( ksps ); CHKERRQ(ierr);
788*a2f3521dSMark F. Adams       }
789*a2f3521dSMark F. Adams       ierr = GAMGKKTMatDestroy( &kktMatsArr[level] ); CHKERRQ(ierr);
790*a2f3521dSMark F. Adams 
791*a2f3521dSMark F. Adams       /* set defaults */
792*a2f3521dSMark F. Adams       ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
793*a2f3521dSMark F. Adams 
794ffc955d6SMark F. Adams       /* override defaults and command line args (!) */
795ffc955d6SMark F. Adams       if ( pc_gamg->use_aggs_in_gasm ) {
7962d3561bbSSatish Balay         PetscInt sz;
7972d3561bbSSatish Balay         IS *is;
798*a2f3521dSMark F. Adams 
7992d3561bbSSatish Balay         sz = nASMBlocksArr[level];
8002d3561bbSSatish Balay         is = ASMLocalIDsArr[level];
801ffc955d6SMark F. Adams         ierr = PCSetType( subpc, PCGASM ); CHKERRQ(ierr);
802ffc955d6SMark F. Adams         if(sz==0){
803ffc955d6SMark F. Adams           IS is;
804ffc955d6SMark F. Adams           PetscInt my0,kk;
805ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange( Aarr[level], &my0, &kk ); CHKERRQ(ierr);
806ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is ); CHKERRQ(ierr);
807a94c3b12SMark F. Adams           ierr = PCGASMSetLocalSubdomains( subpc, 1, &is, PETSC_NULL ); CHKERRQ(ierr);
808a94c3b12SMark F. Adams           ierr = ISDestroy( &is ); CHKERRQ(ierr);
809ffc955d6SMark F. Adams         }
810ffc955d6SMark F. Adams         else {
811a94c3b12SMark F. Adams           PetscInt kk;
812a94c3b12SMark F. Adams           ierr = PCGASMSetLocalSubdomains( subpc, sz, is, PETSC_NULL ); CHKERRQ(ierr);
813a94c3b12SMark F. Adams           for(kk=0;kk<sz;kk++){
814a94c3b12SMark F. Adams             ierr = ISDestroy( &is[kk] ); CHKERRQ(ierr);
815a94c3b12SMark F. Adams           }
816ffc955d6SMark F. Adams           ierr = PetscFree( is ); CHKERRQ(ierr);
817ffc955d6SMark F. Adams         }
818ffc955d6SMark F. Adams         ierr = PCGASMSetOverlap( subpc, 0 ); CHKERRQ(ierr);
819ffc955d6SMark F. Adams 
820ffc955d6SMark F. Adams         ASMLocalIDsArr[level] = PETSC_NULL;
821ffc955d6SMark F. Adams         nASMBlocksArr[level] = 0;
822ffc955d6SMark F. Adams         ierr = PCGASMSetType( subpc, PC_GASM_BASIC ); CHKERRQ(ierr);
823ffc955d6SMark F. Adams       }
824ffc955d6SMark F. Adams       else {
8259d5b6da9SMark F. Adams         ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr);
826ffc955d6SMark F. Adams       }
827d5280255SMark F. Adams     }
828d5280255SMark F. Adams     {
829d5280255SMark F. Adams       /* coarse grid */
830d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
831d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
832d5280255SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
833d5280255SMark F. Adams       ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
834d5280255SMark F. Adams       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
835d5280255SMark F. Adams       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
836d5280255SMark F. Adams       ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
837d5280255SMark F. Adams       ierr = PCSetUp( subpc ); CHKERRQ(ierr);
838d5280255SMark F. Adams       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
839d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
840d5280255SMark F. Adams       ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
841d5280255SMark F. Adams     }
842d5280255SMark F. Adams 
843d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
844d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin( (PetscObject)pc );CHKERRQ(ierr);
845d5280255SMark F. Adams     ierr = PCSetFromOptions_MG( pc ); CHKERRQ(ierr);
846d5280255SMark F. Adams     ierr = PetscOptionsEnd();  CHKERRQ(ierr);
847d5280255SMark 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.");
848d5280255SMark F. Adams 
849d5280255SMark F. Adams     /* create cheby smoothers */
850d5280255SMark F. Adams     for ( lidx = 1, level = pc_gamg->Nlevels-2;
851d5280255SMark F. Adams           lidx <= fine_level;
852d5280255SMark F. Adams           lidx++, level--) {
853d5280255SMark F. Adams       KSP smoother;
854ffc955d6SMark F. Adams       PetscBool flag;
855d5280255SMark F. Adams       PC subpc;
856d5280255SMark F. Adams 
857ffc955d6SMark F. Adams       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
858*a2f3521dSMark F. Adams       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
859*a2f3521dSMark F. Adams 
860*a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
861*a2f3521dSMark F. Adams       if( stokes ) {
862*a2f3521dSMark F. Adams         KSP *ksps;
863*a2f3521dSMark F. Adams         PetscInt nn;
864*a2f3521dSMark F. Adams         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps); CHKERRQ(ierr);
865*a2f3521dSMark F. Adams         smoother = ksps[0];
866*a2f3521dSMark F. Adams         ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
867*a2f3521dSMark F. Adams         ierr = PetscFree( ksps );  CHKERRQ(ierr);
868*a2f3521dSMark F. Adams       }
869ffc955d6SMark F. Adams 
870ffc955d6SMark F. Adams       /* do my own cheby */
871251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &flag ); CHKERRQ(ierr);
872ffc955d6SMark F. Adams       if( flag ) {
873ffc955d6SMark F. Adams         PetscReal emax, emin;
874251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare( (PetscObject)subpc, PCJACOBI, &flag ); CHKERRQ(ierr);
875ffc955d6SMark F. Adams         if( flag && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
876587fa25dSMark F. Adams         else{ /* eigen estimate 'emax' */
877587fa25dSMark F. Adams           KSP eksp; Mat Lmat = Aarr[level];
878b2a4f308SMark F. Adams           Vec bb, xx;
879038e3b61SMark F. Adams 
8805745e0f5SMark F. Adams           ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
8815745e0f5SMark F. Adams           ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
882fc4362bfSMark F. Adams           {
883fc4362bfSMark F. Adams             PetscRandom    rctx;
884fc4362bfSMark F. Adams             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
885fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
886fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
887fc4362bfSMark F. Adams             ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
8885b89ad90SMark F. Adams           }
889*a2f3521dSMark F. Adams 
890a94c3b12SMark F. Adams           if( removedEqs[level] ) {
891302f38e8SMark F. Adams             /* being very careful - zeroing out BC rows (this is not done in agg.c estimates) */
892a94c3b12SMark F. Adams             PetscScalar *zeros;
893bea1fdf0SMark F. Adams             PetscInt ii,jj, *idx_bs, sz, bs=level_bs[level];
894bea1fdf0SMark F. Adams             const PetscInt *idx;
895a94c3b12SMark F. Adams             ierr = ISGetLocalSize( removedEqs[level], &sz ); CHKERRQ(ierr);
896a94c3b12SMark F. Adams             ierr = PetscMalloc( bs*sz*sizeof(PetscScalar), &zeros ); CHKERRQ(ierr);
897a94c3b12SMark F. Adams             for(ii=0;ii<bs*sz;ii++) zeros[ii] = 0.;
898a94c3b12SMark F. Adams             ierr = PetscMalloc( bs*sz*sizeof(PetscInt), &idx_bs ); CHKERRQ(ierr);
899a94c3b12SMark F. Adams             ierr = ISGetIndices( removedEqs[level], &idx); CHKERRQ(ierr);
900a94c3b12SMark F. Adams             for(ii=0;ii<sz;ii++) {
901a94c3b12SMark F. Adams               for(jj=0;jj<bs;jj++) {
902a94c3b12SMark F. Adams                 idx_bs[ii] = bs*idx[ii]+jj;
903a94c3b12SMark F. Adams               }
904a94c3b12SMark F. Adams             }
905a94c3b12SMark F. Adams             ierr = ISRestoreIndices( removedEqs[level], &idx ); CHKERRQ(ierr);
906a94c3b12SMark F. Adams             if( sz > 0 ) {
907a94c3b12SMark F. Adams               ierr = VecSetValues( bb, sz, idx_bs, zeros, INSERT_VALUES );  CHKERRQ(ierr);
908a94c3b12SMark F. Adams             }
909a94c3b12SMark F. Adams             ierr = PetscFree( idx_bs );  CHKERRQ(ierr);
910a94c3b12SMark F. Adams             ierr = PetscFree( zeros );  CHKERRQ(ierr);
911a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb); CHKERRQ(ierr);
912a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb); CHKERRQ(ierr);
913a94c3b12SMark F. Adams           }
914fc4362bfSMark F. Adams           ierr = KSPCreate( wcomm, &eksp );CHKERRQ(ierr);
915038e3b61SMark F. Adams           ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
916fc4362bfSMark F. Adams           CHKERRQ(ierr);
917fc4362bfSMark F. Adams           ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
9181a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
9191a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix( eksp, "gamg_est_");         CHKERRQ(ierr);
920f6536408SMark F. Adams           ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
921f6536408SMark F. Adams 
922f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
923f6536408SMark F. Adams           ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
924fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
9255a9b9e01SMark F. Adams 
926d3d0db20SJed Brown           /* set PC type to be same as smoother */
927ffc955d6SMark F. Adams           ierr = KSPSetPC( eksp, subpc ); CHKERRQ( ierr );
928b2a4f308SMark F. Adams 
9295a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
9305a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
9315a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
932fc4362bfSMark F. Adams           ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
9335a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
9345a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
9355a9b9e01SMark F. Adams 
936fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
9375a9b9e01SMark F. Adams 
938fc4362bfSMark F. Adams           ierr = VecDestroy( &xx );       CHKERRQ(ierr);
939fc4362bfSMark F. Adams           ierr = VecDestroy( &bb );       CHKERRQ(ierr);
940fc4362bfSMark F. Adams           ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
941f6536408SMark F. Adams 
942ffc955d6SMark F. Adams           if( pc_gamg->verbose > 0 ) {
943a94c3b12SMark F. Adams             PetscInt N1, tt;
944a94c3b12SMark F. Adams             ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
945a94c3b12SMark F. Adams             PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1);
946f6536408SMark F. Adams           }
947fc4362bfSMark F. Adams         }
948038e3b61SMark F. Adams         {
949038e3b61SMark F. Adams           PetscInt N1, N0, tt;
950038e3b61SMark F. Adams           ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
951038e3b61SMark F. Adams           ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
952f6536408SMark F. Adams           /* heuristic - is this crap? */
953f6536408SMark F. Adams           emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
954038e3b61SMark F. Adams           emax *= 1.05;
955038e3b61SMark F. Adams         }
956fc4362bfSMark F. Adams         ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
957ffc955d6SMark F. Adams       } /* setup checby flag */
958*a2f3521dSMark F. Adams 
959*a2f3521dSMark F. Adams       if( removedEqs[level] ) {
960*a2f3521dSMark F. Adams         ierr = ISDestroy( &removedEqs[level] );                    CHKERRQ(ierr);
961*a2f3521dSMark F. Adams       }
962ffc955d6SMark F. Adams     } /* non-coarse levels */
963737a81a9SMark F. Adams 
964d5280255SMark F. Adams     /* clean up */
965d5280255SMark F. Adams     for(level=1;level<pc_gamg->Nlevels;level++){
966587fa25dSMark F. Adams       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
967587fa25dSMark F. Adams       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
9685b89ad90SMark F. Adams     }
9690cbbd2e1SMark F. Adams 
9700cbbd2e1SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
9710cbbd2e1SMark F. Adams 
972a94c3b12SMark F. Adams     if( PETSC_FALSE ){
97358471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
9749d5b6da9SMark F. Adams       ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
97558471d46SMark F. Adams       ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
97658471d46SMark F. Adams     }
9775f8cf99dSMark F. Adams   }
9785f8cf99dSMark F. Adams   else {
9795f8cf99dSMark F. Adams     KSP smoother;
980b43b03e9SMark 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__);
9819d5b6da9SMark F. Adams     ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
9825f8cf99dSMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
9835f8cf99dSMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
9849d5b6da9SMark F. Adams     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
9855f8cf99dSMark F. Adams   }
98684d3f75bSMark F. Adams 
9875b89ad90SMark F. Adams   PetscFunctionReturn(0);
9885b89ad90SMark F. Adams }
9895b89ad90SMark F. Adams 
990eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
9915b89ad90SMark F. Adams /*
9925b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
9935b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
9945b89ad90SMark F. Adams 
9955b89ad90SMark F. Adams    Input Parameter:
9965b89ad90SMark F. Adams .  pc - the preconditioner context
9975b89ad90SMark F. Adams 
9985b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
9995b89ad90SMark F. Adams */
10005b89ad90SMark F. Adams #undef __FUNCT__
10015b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
10025b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG( PC pc )
10035b89ad90SMark F. Adams {
10045b89ad90SMark F. Adams   PetscErrorCode  ierr;
10055b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
10065b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
10075b89ad90SMark F. Adams 
10085b89ad90SMark F. Adams   PetscFunctionBegin;
10095b89ad90SMark F. Adams   ierr = PCReset_GAMG( pc );CHKERRQ(ierr);
10105b89ad90SMark F. Adams   ierr = PetscFree( pc_gamg );CHKERRQ(ierr);
10115b89ad90SMark F. Adams   ierr = PCDestroy_MG( pc );CHKERRQ(ierr);
10125b89ad90SMark F. Adams   PetscFunctionReturn(0);
10135b89ad90SMark F. Adams }
10145b89ad90SMark F. Adams 
1015676e1743SMark F. Adams 
1016676e1743SMark F. Adams #undef __FUNCT__
1017676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
1018676e1743SMark F. Adams /*@
1019676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
1020676e1743SMark F. Adams    processor reduction.
1021676e1743SMark F. Adams 
1022676e1743SMark F. Adams    Not Collective on PC
1023676e1743SMark F. Adams 
1024676e1743SMark F. Adams    Input Parameters:
1025676e1743SMark F. Adams .  pc - the preconditioner context
1026676e1743SMark F. Adams 
1027676e1743SMark F. Adams 
1028676e1743SMark F. Adams    Options Database Key:
1029676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
1030676e1743SMark F. Adams 
1031676e1743SMark F. Adams    Level: intermediate
1032676e1743SMark F. Adams 
1033676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1034676e1743SMark F. Adams 
1035676e1743SMark F. Adams .seealso: ()
1036676e1743SMark F. Adams @*/
1037676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
1038676e1743SMark F. Adams {
1039676e1743SMark F. Adams   PetscErrorCode ierr;
1040676e1743SMark F. Adams 
1041676e1743SMark F. Adams   PetscFunctionBegin;
1042676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1043676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1044676e1743SMark F. Adams   PetscFunctionReturn(0);
1045676e1743SMark F. Adams }
1046676e1743SMark F. Adams 
1047676e1743SMark F. Adams EXTERN_C_BEGIN
1048676e1743SMark F. Adams #undef __FUNCT__
1049676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
1050676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1051676e1743SMark F. Adams {
1052c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1053c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1054676e1743SMark F. Adams 
1055676e1743SMark F. Adams   PetscFunctionBegin;
10569d5b6da9SMark F. Adams   if(n>0) pc_gamg->min_eq_proc = n;
1057676e1743SMark F. Adams   PetscFunctionReturn(0);
1058676e1743SMark F. Adams }
1059676e1743SMark F. Adams EXTERN_C_END
1060676e1743SMark F. Adams 
1061676e1743SMark F. Adams #undef __FUNCT__
1062389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
1063389730f3SMark F. Adams /*@
1064389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
1065389730f3SMark F. Adams 
1066389730f3SMark F. Adams  Collective on PC
1067389730f3SMark F. Adams 
1068389730f3SMark F. Adams    Input Parameters:
1069389730f3SMark F. Adams .  pc - the preconditioner context
1070389730f3SMark F. Adams 
1071389730f3SMark F. Adams 
1072389730f3SMark F. Adams    Options Database Key:
1073389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
1074389730f3SMark F. Adams 
1075389730f3SMark F. Adams    Level: intermediate
1076389730f3SMark F. Adams 
1077389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1078389730f3SMark F. Adams 
1079389730f3SMark F. Adams .seealso: ()
1080389730f3SMark F. Adams  @*/
1081389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
1082389730f3SMark F. Adams {
1083389730f3SMark F. Adams   PetscErrorCode ierr;
1084389730f3SMark F. Adams 
1085389730f3SMark F. Adams   PetscFunctionBegin;
1086389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1087389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1088389730f3SMark F. Adams   PetscFunctionReturn(0);
1089389730f3SMark F. Adams }
1090389730f3SMark F. Adams 
1091389730f3SMark F. Adams EXTERN_C_BEGIN
1092389730f3SMark F. Adams #undef __FUNCT__
1093389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
1094389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1095389730f3SMark F. Adams {
1096389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1097389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1098389730f3SMark F. Adams 
1099389730f3SMark F. Adams   PetscFunctionBegin;
11009d5b6da9SMark F. Adams   if(n>0) pc_gamg->coarse_eq_limit = n;
1101389730f3SMark F. Adams   PetscFunctionReturn(0);
1102389730f3SMark F. Adams }
1103389730f3SMark F. Adams EXTERN_C_END
1104389730f3SMark F. Adams 
1105389730f3SMark F. Adams #undef __FUNCT__
11068263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1107676e1743SMark F. Adams /*@
11088263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1109676e1743SMark F. Adams 
1110676e1743SMark F. Adams    Collective on PC
1111676e1743SMark F. Adams 
1112676e1743SMark F. Adams    Input Parameters:
1113676e1743SMark F. Adams .  pc - the preconditioner context
1114676e1743SMark F. Adams 
1115676e1743SMark F. Adams 
1116676e1743SMark F. Adams    Options Database Key:
11178263b398SMark F. Adams .  -pc_gamg_repartition
1118676e1743SMark F. Adams 
1119676e1743SMark F. Adams    Level: intermediate
1120676e1743SMark F. Adams 
1121676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1122676e1743SMark F. Adams 
1123676e1743SMark F. Adams .seealso: ()
1124676e1743SMark F. Adams @*/
11258263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1126676e1743SMark F. Adams {
1127676e1743SMark F. Adams   PetscErrorCode ierr;
1128676e1743SMark F. Adams 
1129676e1743SMark F. Adams   PetscFunctionBegin;
1130676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11318263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1132676e1743SMark F. Adams   PetscFunctionReturn(0);
1133676e1743SMark F. Adams }
1134676e1743SMark F. Adams 
1135676e1743SMark F. Adams EXTERN_C_BEGIN
1136676e1743SMark F. Adams #undef __FUNCT__
11378263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
11388263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1139676e1743SMark F. Adams {
1140c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1141c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1142676e1743SMark F. Adams 
1143676e1743SMark F. Adams   PetscFunctionBegin;
11449d5b6da9SMark F. Adams   pc_gamg->repart = n;
1145676e1743SMark F. Adams   PetscFunctionReturn(0);
1146676e1743SMark F. Adams }
1147676e1743SMark F. Adams EXTERN_C_END
1148676e1743SMark F. Adams 
1149676e1743SMark F. Adams #undef __FUNCT__
1150ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1151ffc955d6SMark F. Adams /*@
1152ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1153ffc955d6SMark F. Adams 
1154ffc955d6SMark F. Adams    Collective on PC
1155ffc955d6SMark F. Adams 
1156ffc955d6SMark F. Adams    Input Parameters:
1157ffc955d6SMark F. Adams .  pc - the preconditioner context
1158ffc955d6SMark F. Adams 
1159ffc955d6SMark F. Adams 
1160ffc955d6SMark F. Adams    Options Database Key:
1161ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1162ffc955d6SMark F. Adams 
1163ffc955d6SMark F. Adams    Level: intermediate
1164ffc955d6SMark F. Adams 
1165ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1166ffc955d6SMark F. Adams 
1167ffc955d6SMark F. Adams .seealso: ()
1168ffc955d6SMark F. Adams @*/
1169ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1170ffc955d6SMark F. Adams {
1171ffc955d6SMark F. Adams   PetscErrorCode ierr;
1172ffc955d6SMark F. Adams 
1173ffc955d6SMark F. Adams   PetscFunctionBegin;
1174ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1175ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1176ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1177ffc955d6SMark F. Adams }
1178ffc955d6SMark F. Adams 
1179ffc955d6SMark F. Adams EXTERN_C_BEGIN
1180ffc955d6SMark F. Adams #undef __FUNCT__
1181ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
1182ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1183ffc955d6SMark F. Adams {
1184ffc955d6SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1185ffc955d6SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1186ffc955d6SMark F. Adams 
1187ffc955d6SMark F. Adams   PetscFunctionBegin;
1188ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1189ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1190ffc955d6SMark F. Adams }
1191ffc955d6SMark F. Adams EXTERN_C_END
1192ffc955d6SMark F. Adams 
1193ffc955d6SMark F. Adams #undef __FUNCT__
11944ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
11954ef23d27SMark F. Adams /*@
11964ef23d27SMark F. Adams    PCGAMGSetNlevels -
11974ef23d27SMark F. Adams 
11984ef23d27SMark F. Adams    Not collective on PC
11994ef23d27SMark F. Adams 
12004ef23d27SMark F. Adams    Input Parameters:
12014ef23d27SMark F. Adams .  pc - the preconditioner context
12024ef23d27SMark F. Adams 
12034ef23d27SMark F. Adams 
12044ef23d27SMark F. Adams    Options Database Key:
12054ef23d27SMark F. Adams .  -pc_mg_levels
12064ef23d27SMark F. Adams 
12074ef23d27SMark F. Adams    Level: intermediate
12084ef23d27SMark F. Adams 
12094ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12104ef23d27SMark F. Adams 
12114ef23d27SMark F. Adams .seealso: ()
12124ef23d27SMark F. Adams @*/
12134ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12144ef23d27SMark F. Adams {
12154ef23d27SMark F. Adams   PetscErrorCode ierr;
12164ef23d27SMark F. Adams 
12174ef23d27SMark F. Adams   PetscFunctionBegin;
12184ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12194ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
12204ef23d27SMark F. Adams   PetscFunctionReturn(0);
12214ef23d27SMark F. Adams }
12224ef23d27SMark F. Adams 
12234ef23d27SMark F. Adams EXTERN_C_BEGIN
12244ef23d27SMark F. Adams #undef __FUNCT__
12254ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
12264ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12274ef23d27SMark F. Adams {
12284ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
12294ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
12304ef23d27SMark F. Adams 
12314ef23d27SMark F. Adams   PetscFunctionBegin;
12329d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12334ef23d27SMark F. Adams   PetscFunctionReturn(0);
12344ef23d27SMark F. Adams }
12354ef23d27SMark F. Adams EXTERN_C_END
12364ef23d27SMark F. Adams 
12374ef23d27SMark F. Adams #undef __FUNCT__
12383542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
12393542efc5SMark F. Adams /*@
12403542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12413542efc5SMark F. Adams 
12423542efc5SMark F. Adams    Not collective on PC
12433542efc5SMark F. Adams 
12443542efc5SMark F. Adams    Input Parameters:
12453542efc5SMark F. Adams .  pc - the preconditioner context
12463542efc5SMark F. Adams 
12473542efc5SMark F. Adams 
12483542efc5SMark F. Adams    Options Database Key:
12493542efc5SMark F. Adams .  -pc_gamg_threshold
12503542efc5SMark F. Adams 
12513542efc5SMark F. Adams    Level: intermediate
12523542efc5SMark F. Adams 
12533542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12543542efc5SMark F. Adams 
12553542efc5SMark F. Adams .seealso: ()
12563542efc5SMark F. Adams @*/
12573542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
12583542efc5SMark F. Adams {
12593542efc5SMark F. Adams   PetscErrorCode ierr;
12603542efc5SMark F. Adams 
12613542efc5SMark F. Adams   PetscFunctionBegin;
12623542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12633542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
12643542efc5SMark F. Adams   PetscFunctionReturn(0);
12653542efc5SMark F. Adams }
12663542efc5SMark F. Adams 
12673542efc5SMark F. Adams EXTERN_C_BEGIN
12683542efc5SMark F. Adams #undef __FUNCT__
12693542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
12703542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
12713542efc5SMark F. Adams {
1272c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1273c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
12743542efc5SMark F. Adams 
12753542efc5SMark F. Adams   PetscFunctionBegin;
12769d5b6da9SMark F. Adams   pc_gamg->threshold = n;
12773542efc5SMark F. Adams   PetscFunctionReturn(0);
12783542efc5SMark F. Adams }
12793542efc5SMark F. Adams EXTERN_C_END
12803542efc5SMark F. Adams 
12813542efc5SMark F. Adams #undef __FUNCT__
12829d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1283676e1743SMark F. Adams /*@
12849d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1285676e1743SMark F. Adams 
1286676e1743SMark F. Adams    Collective on PC
1287676e1743SMark F. Adams 
1288676e1743SMark F. Adams    Input Parameters:
1289676e1743SMark F. Adams .  pc - the preconditioner context
1290676e1743SMark F. Adams 
1291676e1743SMark F. Adams 
1292676e1743SMark F. Adams    Options Database Key:
12933542efc5SMark F. Adams .  -pc_gamg_type
1294676e1743SMark F. Adams 
1295676e1743SMark F. Adams    Level: intermediate
1296676e1743SMark F. Adams 
1297676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1298676e1743SMark F. Adams 
1299676e1743SMark F. Adams .seealso: ()
1300676e1743SMark F. Adams @*/
13010202ef74SSatish Balay PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type )
1302676e1743SMark F. Adams {
1303676e1743SMark F. Adams   PetscErrorCode ierr;
1304676e1743SMark F. Adams 
1305676e1743SMark F. Adams   PetscFunctionBegin;
1306676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13070202ef74SSatish Balay   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type));
1308676e1743SMark F. Adams   CHKERRQ(ierr);
1309676e1743SMark F. Adams   PetscFunctionReturn(0);
1310676e1743SMark F. Adams }
1311676e1743SMark F. Adams 
1312676e1743SMark F. Adams EXTERN_C_BEGIN
1313676e1743SMark F. Adams #undef __FUNCT__
13149d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
13150202ef74SSatish Balay PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type )
1316676e1743SMark F. Adams {
13179d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
1318676e1743SMark F. Adams 
1319676e1743SMark F. Adams   PetscFunctionBegin;
13209d5b6da9SMark F. Adams   ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r);
13219d5b6da9SMark F. Adams   CHKERRQ(ierr);
13229d5b6da9SMark F. Adams 
13239d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
13249d5b6da9SMark F. Adams 
13259d5b6da9SMark F. Adams   /* call sub create method */
13269d5b6da9SMark F. Adams   ierr = (*r)(pc); CHKERRQ(ierr);
13279d5b6da9SMark F. Adams 
1328676e1743SMark F. Adams   PetscFunctionReturn(0);
1329676e1743SMark F. Adams }
1330676e1743SMark F. Adams EXTERN_C_END
1331676e1743SMark F. Adams 
13325b89ad90SMark F. Adams #undef __FUNCT__
13335b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
13345b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG( PC pc )
13355b89ad90SMark F. Adams {
1336676e1743SMark F. Adams   PetscErrorCode  ierr;
1337676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1338676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1339676e1743SMark F. Adams   PetscBool        flag;
1340b43b03e9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
13415b89ad90SMark F. Adams 
13425b89ad90SMark F. Adams   PetscFunctionBegin;
1343676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1344676e1743SMark F. Adams   {
134575b74bdaSMark F. Adams     /* -pc_gamg_verbose */
13469d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
13479d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
13489d5b6da9SMark F. Adams                            &pc_gamg->verbose, PETSC_NULL );
13499d5b6da9SMark F. Adams     CHKERRQ(ierr);
135075b74bdaSMark F. Adams 
13518263b398SMark F. Adams     /* -pc_gamg_repartition */
13528263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
13538263b398SMark F. Adams                             "Repartion coarse grids (false)",
13548263b398SMark F. Adams                             "PCGAMGRepartitioning",
13559d5b6da9SMark F. Adams                             pc_gamg->repart,
13569d5b6da9SMark F. Adams                             &pc_gamg->repart,
1357676e1743SMark F. Adams                             &flag);
1358676e1743SMark F. Adams     CHKERRQ(ierr);
1359676e1743SMark F. Adams 
1360ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
1361ffc955d6SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1362ffc955d6SMark F. Adams                             "Use aggregation agragates for GASM smoother (false)",
1363ffc955d6SMark F. Adams                             "PCGAMGUseASMAggs",
1364ffc955d6SMark F. Adams                             pc_gamg->use_aggs_in_gasm,
1365ffc955d6SMark F. Adams                             &pc_gamg->use_aggs_in_gasm,
1366ffc955d6SMark F. Adams                             &flag);
1367ffc955d6SMark F. Adams     CHKERRQ(ierr);
1368ffc955d6SMark F. Adams 
1369c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1370676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1371676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1372676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
13739d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
13749d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1375676e1743SMark F. Adams                            &flag );
1376676e1743SMark F. Adams     CHKERRQ(ierr);
13773542efc5SMark F. Adams 
1378389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1379389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1380389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1381389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
13829d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
13839d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1384389730f3SMark F. Adams                            &flag );
1385389730f3SMark F. Adams     CHKERRQ(ierr);
1386389730f3SMark F. Adams 
1387c20e4228SMark F. Adams     /* -pc_gamg_threshold */
13883542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
13893542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
13903542efc5SMark F. Adams                             "PCGAMGSetThreshold",
13919d5b6da9SMark F. Adams                             pc_gamg->threshold,
13929d5b6da9SMark F. Adams                             &pc_gamg->threshold,
13933542efc5SMark F. Adams                             &flag );
13943542efc5SMark F. Adams     CHKERRQ(ierr);
1395b43b03e9SMark F. Adams     if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);
13964ef23d27SMark F. Adams 
13974ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
13984ef23d27SMark F. Adams                            "Set number of MG levels",
13994ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
14009d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
14019d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
14024ef23d27SMark F. Adams                            &flag );
1403676e1743SMark F. Adams   }
1404676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1405676e1743SMark F. Adams 
14065b89ad90SMark F. Adams   PetscFunctionReturn(0);
14075b89ad90SMark F. Adams }
14085b89ad90SMark F. Adams 
14095b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
14105b89ad90SMark F. Adams /*MC
1411856830bfSJed Brown      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
14129d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
14135b89ad90SMark F. Adams 
1414280d9858SJed Brown    Options Database Keys:
14155b89ad90SMark F. Adams    Multigrid options(inherited)
1416280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1417280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1418280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1419280d9858SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
14205b89ad90SMark F. Adams 
14215b89ad90SMark F. Adams   Level: intermediate
1422280d9858SJed Brown 
14235b89ad90SMark F. Adams   Concepts: multigrid
14245b89ad90SMark F. Adams 
14255b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1426280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
14275b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
14285b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
14295b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
14305b89ad90SMark F. Adams M*/
14315b89ad90SMark F. Adams EXTERN_C_BEGIN
14325b89ad90SMark F. Adams #undef __FUNCT__
14335b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
14345b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG( PC pc )
14355b89ad90SMark F. Adams {
14365b89ad90SMark F. Adams   PetscErrorCode  ierr;
14375b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
14385b89ad90SMark F. Adams   PC_MG           *mg;
14390cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14409d5b6da9SMark F. Adams   static long count = 0;
14415ee9c036SSatish Balay #endif
14425b89ad90SMark F. Adams 
14435b89ad90SMark F. Adams   PetscFunctionBegin;
144460a6d8f8SMark F. Adams 
14455b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14465b89ad90SMark F. Adams   ierr = PCSetType( pc, PCMG );  CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
14475b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr);
14485b89ad90SMark F. Adams 
14495b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
14505b89ad90SMark F. Adams   ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr);
14515b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
1452ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
14535b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
14545b89ad90SMark F. Adams 
14559d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14569d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14579d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14589d5b6da9SMark F. Adams   pc_gamg->data = 0;
14595b89ad90SMark F. Adams 
14609d5b6da9SMark F. Adams   /* register AMG type */
14619d5b6da9SMark F. Adams   if( !GAMGList ){
14629d5b6da9SMark F. Adams     ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
14639d5b6da9SMark F. Adams     ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
14649d5b6da9SMark F. Adams   }
14659d5b6da9SMark F. Adams 
14669d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
14675b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
14685b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
14695b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
14705b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
14715b89ad90SMark F. Adams 
14725b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1473676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1474676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1475676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1476676e1743SMark F. Adams   CHKERRQ(ierr);
1477676e1743SMark F. Adams 
1478676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1479389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_C",
1480389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_GAMG",
1481389730f3SMark F. Adams 					    PCGAMGSetCoarseEqLim_GAMG);
1482389730f3SMark F. Adams   CHKERRQ(ierr);
1483389730f3SMark F. Adams 
1484389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
14858263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_C",
14868263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_GAMG",
14878263b398SMark F. Adams 					    PCGAMGSetRepartitioning_GAMG);
1488676e1743SMark F. Adams   CHKERRQ(ierr);
1489676e1743SMark F. Adams 
1490676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1491ffc955d6SMark F. Adams 					    "PCGAMGSetUseASMAggs_C",
1492ffc955d6SMark F. Adams 					    "PCGAMGSetUseASMAggs_GAMG",
1493ffc955d6SMark F. Adams 					    PCGAMGSetUseASMAggs_GAMG);
1494ffc955d6SMark F. Adams   CHKERRQ(ierr);
1495ffc955d6SMark F. Adams 
1496ffc955d6SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
14973542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
14983542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
14993542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
15003542efc5SMark F. Adams   CHKERRQ(ierr);
15013542efc5SMark F. Adams 
15023542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
15039d5b6da9SMark F. Adams 					    "PCGAMGSetType_C",
15049d5b6da9SMark F. Adams 					    "PCGAMGSetType_GAMG",
15059d5b6da9SMark F. Adams 					    PCGAMGSetType_GAMG);
1506676e1743SMark F. Adams   CHKERRQ(ierr);
1507c97e1df0SMark F. Adams 
15089d5b6da9SMark F. Adams   pc_gamg->repart = PETSC_FALSE;
1509ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
15109d5b6da9SMark F. Adams   pc_gamg->min_eq_proc = 100;
1511c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit = 800;
1512*a2f3521dSMark F. Adams   pc_gamg->threshold = 0.001;
15139d5b6da9SMark F. Adams   pc_gamg->Nlevels = GAMG_MAXLEVELS;
15149d5b6da9SMark F. Adams   pc_gamg->verbose = 0;
15159d5b6da9SMark F. Adams   pc_gamg->emax_id = -1;
15169d5b6da9SMark F. Adams   pc_gamg->col_bs_id = -1;
15179d5b6da9SMark F. Adams 
15180cbbd2e1SMark F. Adams   /* private events */
15190cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1520785cba28SMark F. Adams   if( count++ == 0 ) {
15210cbbd2e1SMark F. Adams     PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
15220cbbd2e1SMark F. Adams     PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
15230cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
15240cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
15250cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
15260cbbd2e1SMark F. Adams     PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
15270cbbd2e1SMark F. Adams     PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
15280cbbd2e1SMark F. Adams     PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
15290cbbd2e1SMark F. Adams     PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
153060a6d8f8SMark F. Adams     PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
15310cbbd2e1SMark F. Adams     PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
15320cbbd2e1SMark F. Adams     PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
15330cbbd2e1SMark F. Adams     PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
15340cbbd2e1SMark F. Adams     PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
15350cbbd2e1SMark F. Adams     PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
15360cbbd2e1SMark F. Adams     PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
15370cbbd2e1SMark F. Adams     PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);
1538f852f58cSMark F. Adams 
15390cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
15400cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
15410cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1542b4fbaa2aSMark F. Adams     /* create timer stages */
1543b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1544b4fbaa2aSMark F. Adams     {
1545b4fbaa2aSMark F. Adams       char str[32];
1546b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1547b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1548b4fbaa2aSMark F. Adams       PetscInt lidx;
1549b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1550b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1551b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1552b4fbaa2aSMark F. Adams       }
1553b4fbaa2aSMark F. Adams     }
1554b4fbaa2aSMark F. Adams #endif
1555b4fbaa2aSMark F. Adams   }
1556b4fbaa2aSMark F. Adams #endif
15570cbbd2e1SMark F. Adams   /* general events */
15580cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
155960a6d8f8SMark F. Adams   PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);
15600cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);
15610cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
15620cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
15630cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
15640cbbd2e1SMark F. Adams   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
156560a6d8f8SMark F. Adams   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);
1566*a2f3521dSMark F. Adams   PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);
15670cbbd2e1SMark F. Adams #endif
15680cbbd2e1SMark F. Adams 
156960a6d8f8SMark F. Adams   /* instantiate derived type */
157060a6d8f8SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
157160a6d8f8SMark F. Adams   {
157260a6d8f8SMark F. Adams     char tname[256] = GAMGAGG;
157360a6d8f8SMark F. Adams     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",
157460a6d8f8SMark F. Adams                             GAMGList, tname, tname, sizeof(tname), PETSC_NULL );
157560a6d8f8SMark F. Adams     CHKERRQ(ierr);
157660a6d8f8SMark F. Adams     ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr);
157760a6d8f8SMark F. Adams   }
157860a6d8f8SMark F. Adams   ierr = PetscOptionsTail();   CHKERRQ(ierr);
157960a6d8f8SMark F. Adams 
15805b89ad90SMark F. Adams   PetscFunctionReturn(0);
15815b89ad90SMark F. Adams }
15825b89ad90SMark F. Adams EXTERN_C_END
1583