xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 5f8cf99da6441de8d319d2c2cfa5792d12e1a83c)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4313a3e24SSatish Balay #include "private/matimpl.h"
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
6f96513f1SMatthew G Knepley 
7b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
8b4fbaa2aSMark F. Adams PetscLogEvent gamg_setup_events[NUM_SET];
9b4fbaa2aSMark F. Adams #endif
10b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
11b4fbaa2aSMark F. Adams 
12b8fd24d8SMark F. Adams /*#define GAMG_STAGES*/
13b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
14b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
15b4fbaa2aSMark F. Adams #endif
16f96513f1SMatthew G Knepley 
175b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
185b89ad90SMark F. Adams /*
195b89ad90SMark F. Adams    PCSetCoordinates_GAMG
205b89ad90SMark F. Adams 
215b89ad90SMark F. Adams    Input Parameter:
225b89ad90SMark F. Adams    .  pc - the preconditioner context
235b89ad90SMark F. Adams */
24a92563c5SMark F. Adams EXTERN_C_BEGIN
255b89ad90SMark F. Adams #undef __FUNCT__
265b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG"
27eb07cef2SMark F. Adams PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords )
285b89ad90SMark F. Adams {
29eb07cef2SMark F. Adams   PC_MG          *mg = (PC_MG*)a_pc->data;
305b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
316c237d78SBarry Smith   PetscErrorCode ierr;
32d3d6bff4SMark F. Adams   PetscInt       arrsz,bs,my0,kk,ii,jj,nloc,Iend;
33038e3b61SMark F. Adams   Mat            Amat = a_pc->pmat;
345b89ad90SMark F. Adams 
355b89ad90SMark F. Adams   PetscFunctionBegin;
3658471d46SMark F. Adams   PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 );
37038e3b61SMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
38d3d6bff4SMark F. Adams   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
39d3d6bff4SMark F. Adams   nloc = (Iend-my0)/bs;
40d3d6bff4SMark F. Adams   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
41038e3b61SMark F. Adams 
42d3d6bff4SMark F. Adams   pc_gamg->m_data_rows = 1;
43f6536408SMark F. Adams   if(a_coords==0 && pc_gamg->m_method==0) {
44f6536408SMark F. Adams     SETERRQ(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
45f6536408SMark F. Adams   }
46470e26f8SMark F. Adams   if( pc_gamg->m_method==0 ) pc_gamg->m_data_cols = a_ndm; /* coordinates */
47038e3b61SMark F. Adams   else{ /* SA: null space vectors */
48d3d6bff4SMark F. Adams     if(a_coords != 0 && bs==1 ) pc_gamg->m_data_cols = 1; /* scalar w/ coords and SA (not needed) */
49d3d6bff4SMark F. Adams     else if(a_coords != 0 ) pc_gamg->m_data_cols = (a_ndm==2 ? 3 : 6); /* elasticity */
50d3d6bff4SMark F. Adams     else pc_gamg->m_data_cols = bs; /* no data, force SA with constant null space vectors */
51d3d6bff4SMark F. Adams     pc_gamg->m_data_rows = bs;
52038e3b61SMark F. Adams   }
53d3d6bff4SMark F. Adams   arrsz = nloc*pc_gamg->m_data_rows*pc_gamg->m_data_cols;
545b89ad90SMark F. Adams 
55038e3b61SMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
566e3e101aSMark F. Adams   if (pc_gamg->m_data==0 || (pc_gamg->m_data_sz != arrsz)) {
575b89ad90SMark F. Adams     ierr = PetscFree( pc_gamg->m_data );  CHKERRQ(ierr);
5833a2b366SJed Brown     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->m_data ); CHKERRQ(ierr);
595b89ad90SMark F. Adams   }
60038e3b61SMark F. Adams   for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.;
61eb07cef2SMark F. Adams   pc_gamg->m_data[arrsz] = -99.;
62038e3b61SMark F. Adams   /* copy data in - column oriented */
63470e26f8SMark F. Adams   if( pc_gamg->m_method != 0 ) {
64d3d6bff4SMark F. Adams     const PetscInt M = Iend - my0;
65038e3b61SMark F. Adams     for(kk=0;kk<nloc;kk++){
66038e3b61SMark F. Adams       PetscReal *data = &pc_gamg->m_data[kk*bs];
67d3d6bff4SMark F. Adams       if( pc_gamg->m_data_cols==1 ) *data = 1.0;
68038e3b61SMark F. Adams       else {
69038e3b61SMark F. Adams         for(ii=0;ii<bs;ii++)
70038e3b61SMark F. Adams 	  for(jj=0;jj<bs;jj++)
71038e3b61SMark F. Adams 	    if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */
72038e3b61SMark F. Adams 	    else data[ii*M + jj] = 0.0;
73eb07cef2SMark F. Adams         if( a_coords != 0 ) {
74038e3b61SMark F. Adams           if( a_ndm == 2 ){ /* rotational modes */
75038e3b61SMark F. Adams             data += 2*M;
76038e3b61SMark F. Adams             data[0] = -a_coords[2*kk+1];
77038e3b61SMark F. Adams             data[1] =  a_coords[2*kk];
785b89ad90SMark F. Adams           }
79eb07cef2SMark F. Adams           else {
80038e3b61SMark F. Adams             data += 3*M;
81038e3b61SMark F. Adams             data[0] = 0.0;               data[M+0] =  a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1];
82038e3b61SMark F. Adams             data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0;               data[2*M+1] =  a_coords[3*kk];
83038e3b61SMark F. Adams             data[2] =  a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk];   data[2*M+2] = 0.0;
84038e3b61SMark F. Adams           }
85eb07cef2SMark F. Adams         }
86eb07cef2SMark F. Adams       }
87eb07cef2SMark F. Adams     }
88eb07cef2SMark F. Adams   }
89eb07cef2SMark F. Adams   else {
90038e3b61SMark F. Adams     for( kk = 0 ; kk < nloc ; kk++ ){
91038e3b61SMark F. Adams       for( ii = 0 ; ii < a_ndm ; ii++ ) {
92038e3b61SMark F. Adams         pc_gamg->m_data[ii*nloc + kk] =  a_coords[kk*a_ndm + ii];
93eb07cef2SMark F. Adams       }
94eb07cef2SMark F. Adams     }
95038e3b61SMark F. Adams   }
96eb07cef2SMark F. Adams   assert(pc_gamg->m_data[arrsz] == -99.);
97038e3b61SMark F. Adams 
985b89ad90SMark F. Adams   pc_gamg->m_data_sz = arrsz;
99eb07cef2SMark F. Adams   pc_gamg->m_dim = a_ndm;
1005b89ad90SMark F. Adams 
1015b89ad90SMark F. Adams   PetscFunctionReturn(0);
1025b89ad90SMark F. Adams }
103a92563c5SMark F. Adams EXTERN_C_END
1045b89ad90SMark F. Adams 
105d3d6bff4SMark F. Adams 
106d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
107d3d6bff4SMark F. Adams #undef __FUNCT__
108d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
109d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
110d3d6bff4SMark F. Adams {
111d3d6bff4SMark F. Adams   PetscErrorCode  ierr;
112d3d6bff4SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
113d3d6bff4SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
114d3d6bff4SMark F. Adams 
115d3d6bff4SMark F. Adams   PetscFunctionBegin;
11658471d46SMark F. Adams   if( pc_gamg->m_data != 0 ) { /* this should not happen, cleaned up in SetUp */
117d3d6bff4SMark F. Adams     ierr = PetscFree(pc_gamg->m_data); CHKERRQ(ierr);
11858471d46SMark F. Adams   }
119d3d6bff4SMark F. Adams   pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0;
120d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
121d3d6bff4SMark F. Adams }
122d3d6bff4SMark F. Adams 
1235b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1245b89ad90SMark F. Adams /*
125a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
126a147abb0SMark F. Adams      of active processors.
1275b89ad90SMark F. Adams 
1285b89ad90SMark F. Adams    Input Parameter:
1298263b398SMark F. Adams    . a_pc - parameters
1303530afc2SMark F. Adams    . a_Amat_fine - matrix on this fine (k) level
131d3d6bff4SMark F. Adams    . a_ndata_rows - size of data to move (coarse grid)
132d3d6bff4SMark F. Adams    . a_ndata_cols - size of data to move (coarse grid)
1338263b398SMark F. Adams    . a_cbs - coarse block size
1348263b398SMark F. Adams    . a_isLast -
1353530afc2SMark F. Adams    In/Output Parameter:
1363530afc2SMark F. Adams    . a_P_inout - prolongation operator to the next level (k-1)
137eb07cef2SMark F. Adams    . a_coarse_data - data that need to be moved
138afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
13911e60469SMark F. Adams    Output Parameter:
1403530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1415b89ad90SMark F. Adams */
1425cb416c2SMark F. Adams 
1435b89ad90SMark F. Adams #undef __FUNCT__
1448263b398SMark F. Adams #define __FUNCT__ "createLevel"
1458263b398SMark F. Adams PetscErrorCode createLevel( const PC a_pc,
1468263b398SMark F. Adams                             const Mat a_Amat_fine,
1478263b398SMark F. Adams                             const PetscInt a_ndata_rows,
1488263b398SMark F. Adams                             const PetscInt a_ndata_cols,
1498263b398SMark F. Adams                             const PetscInt a_cbs,
1508263b398SMark F. Adams                             const PetscBool a_isLast,
1513530afc2SMark F. Adams                             Mat *a_P_inout,
152eb07cef2SMark F. Adams                             PetscReal **a_coarse_data,
153afc97cdcSMark F. Adams                             PetscMPIInt *a_nactive_proc,
154389730f3SMark F. Adams                             Mat *a_Amat_crs
15511e60469SMark F. Adams                             )
1565b89ad90SMark F. Adams {
1578263b398SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
158486a8d0bSJed Brown   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1598263b398SMark F. Adams   const PetscBool  repart = pc_gamg->m_repart;
160389730f3SMark F. Adams   const PetscInt   min_eq_proc = pc_gamg->m_min_eq_proc, coarse_max = pc_gamg->m_coarse_eq_limit;
1615b89ad90SMark F. Adams   PetscErrorCode   ierr;
162038e3b61SMark F. Adams   Mat              Cmat,Pnew,Pold=*a_P_inout;
16311e60469SMark F. Adams   IS               new_indices,isnum;
1643530afc2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
165fa31c87dSMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive;
166fa31c87dSMark F. Adams   PetscInt         neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0;
1675b89ad90SMark F. Adams 
1685b89ad90SMark F. Adams   PetscFunctionBegin;
16911e60469SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
17011e60469SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
171f6536408SMark F. Adams 
17211e60469SMark F. Adams   /* RAP */
17396434bc1SMark F. Adams #ifdef USE_R
17496434bc1SMark F. Adams   /* make R wih brute force for now */
17596434bc1SMark F. Adams   ierr = MatTranspose( Pold, Pnew );
17696434bc1SMark F. Adams   ierr = MatDestroy( &Pold );  CHKERRQ(ierr);
17796434bc1SMark F. Adams   ierr = MatRARt( a_Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
17896434bc1SMark F. Adams   Pold = Pnew;
17996434bc1SMark F. Adams #else
180038e3b61SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
18196434bc1SMark F. Adams #endif
182038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
183038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
184038e3b61SMark F. Adams   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
185038e3b61SMark F. Adams 
186f852f58cSMark F. Adams   /* get number of PEs to make active, reduce */
187f852f58cSMark F. Adams   ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
18881708dccSMark F. Adams   new_npe = (PetscMPIInt)((float)neq/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
189f852f58cSMark F. Adams   if( new_npe == 0 || neq < coarse_max ) new_npe = 1;
190f852f58cSMark F. Adams   else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */
1918263b398SMark F. Adams   if( a_isLast ) new_npe = 1;
192f852f58cSMark F. Adams 
1938263b398SMark F. Adams   if( !repart && !(new_npe == 1 && *a_nactive_proc != 1) ) {
194a147abb0SMark F. Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction */
19522063be5SMark F. Adams   }
19622063be5SMark F. Adams   else {
19722063be5SMark F. Adams     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
1985ef31b24SMark F. Adams     Mat              adj;
19922063be5SMark F. Adams     const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols;
200b4fbaa2aSMark F. Adams     const PetscInt  stride0=ncrs0*a_ndata_rows;
2018263b398SMark F. Adams     PetscInt        is_sz,*newproc_idx,ii,jj,kk,strideNew,*tidx;
202c9a0b8beSMark F. Adams     /* create sub communicator  */
20371d60426SBarry Smith     MPI_Comm        cm;
204afc97cdcSMark F. Adams     MPI_Group       wg, g2;
205a147abb0SMark F. Adams     PetscInt       *counts,inpe,targetPE;
206fd3c6afaSMark F. Adams     PetscMPIInt    *ranks;
20722063be5SMark F. Adams     IS              isscat;
20822063be5SMark F. Adams     PetscScalar    *array;
20922063be5SMark F. Adams     Vec             src_crd, dest_crd;
21022063be5SMark F. Adams     PetscReal      *data = *a_coarse_data;
21122063be5SMark F. Adams     VecScatter      vecscat;
212b4fbaa2aSMark F. Adams     IS              isnewproc;
213e33ef3b1SMark F. Adams 
21492a756f0SMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr);
215fd3c6afaSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
21692a756f0SMark F. Adams 
217fd3c6afaSMark F. Adams     ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr);
218afc97cdcSMark F. Adams     assert(counts[mype]==ncrs0);
219afc97cdcSMark F. Adams     /* count real active pes */
220a147abb0SMark F. Adams     for( nactive = jj = 0, targetPE = -1 ; jj < npe ; jj++) {
221afc97cdcSMark F. Adams       if( counts[jj] != 0 ) {
22222063be5SMark F. Adams 	ranks[nactive++] = jj;
223a147abb0SMark F. Adams         if( targetPE == -1 ) targetPE = jj; /* find lowest active PE */
224afc97cdcSMark F. Adams       }
225afc97cdcSMark F. Adams     }
2263303bcf9Sadams 
227a147abb0SMark F. Adams     if (nactive < new_npe) {
228a147abb0SMark F. Adams       if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s reducing number of target PEs: %d --> %d (?)\n",mype,__FUNCT__,new_npe,nactive);
229a147abb0SMark F. Adams       new_npe = nactive; /* this can happen with empty input procs */
230a147abb0SMark F. Adams     }
23122063be5SMark F. Adams 
232bed7c2b9SMark F. Adams     if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s npe (active): %d --> %d. new npe = %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,nactive,new_npe,neq);
23358471d46SMark F. Adams 
23422063be5SMark F. Adams     *a_nactive_proc = new_npe; /* output */
2352f03bc48SMark F. Adams 
236afc97cdcSMark F. Adams     ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr);
23722063be5SMark F. Adams     ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr);
238afc97cdcSMark F. Adams     ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr);
239b4fbaa2aSMark F. Adams 
240afc97cdcSMark F. Adams     ierr = MPI_Group_free( &wg );                            CHKERRQ(ierr);
241afc97cdcSMark F. Adams     ierr = MPI_Group_free( &g2 );                            CHKERRQ(ierr);
242c9a0b8beSMark F. Adams 
2435ef31b24SMark F. Adams     /* MatPartitioningApply call MatConvert, which is collective */
244b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
245b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
246b4fbaa2aSMark F. Adams #endif
247038e3b61SMark F. Adams     if( a_cbs == 1) {
248038e3b61SMark F. Adams       ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
249eb07cef2SMark F. Adams     }
250eb07cef2SMark F. Adams     else{
251038e3b61SMark F. Adams       /* make a scalar matrix to partition */
252eb07cef2SMark F. Adams       Mat tMat;
25358471d46SMark F. Adams       PetscInt ncols,jj,Ii;
254b4fbaa2aSMark F. Adams       const PetscScalar *vals;
255b4fbaa2aSMark F. Adams       const PetscInt *idx;
256*5f8cf99dSMark F. Adams       PetscInt *d_nnz, *o_nnz;
2575ee9c036SSatish Balay       static int llev = 0;
258b4fbaa2aSMark F. Adams 
25958471d46SMark F. Adams       ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
260*5f8cf99dSMark F. Adams       ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr);
26158471d46SMark F. Adams       for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) {
26258471d46SMark F. Adams         ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
26358471d46SMark F. Adams         d_nnz[jj] = ncols/a_cbs;
264*5f8cf99dSMark F. Adams         o_nnz[jj] = ncols/a_cbs;
26558471d46SMark F. Adams         ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
266*5f8cf99dSMark F. Adams         if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
267*5f8cf99dSMark F. Adams         if( o_nnz[jj] > (neq/a_cbs-ncrs0) ) o_nnz[jj] = neq/a_cbs-ncrs0;
26858471d46SMark F. Adams       }
2696876a03eSMark F. Adams 
270eb07cef2SMark F. Adams       ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
271eb07cef2SMark F. Adams                               PETSC_DETERMINE, PETSC_DETERMINE,
272*5f8cf99dSMark F. Adams                               0, d_nnz, 0, o_nnz,
273eb07cef2SMark F. Adams                               &tMat );
2746876a03eSMark F. Adams       CHKERRQ(ierr);
27558471d46SMark F. Adams       ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
276*5f8cf99dSMark F. Adams       ierr = PetscFree( o_nnz ); CHKERRQ(ierr);
277eb07cef2SMark F. Adams 
27822063be5SMark F. Adams       for ( ii = Istart0; ii < Iend0; ii++ ) {
27922063be5SMark F. Adams         PetscInt dest_row = ii/a_cbs;
28022063be5SMark F. Adams         ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
281eb07cef2SMark F. Adams         for( jj = 0 ; jj < ncols ; jj++ ){
282038e3b61SMark F. Adams           PetscInt dest_col = idx[jj]/a_cbs;
283eb07cef2SMark F. Adams           PetscScalar v = 1.0;
284eb07cef2SMark F. Adams           ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
285eb07cef2SMark F. Adams         }
28622063be5SMark F. Adams         ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
287eb07cef2SMark F. Adams       }
288eb07cef2SMark F. Adams       ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
289eb07cef2SMark F. Adams       ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
290eb07cef2SMark F. Adams 
291b4fbaa2aSMark F. Adams       if( llev++ == -1 ) {
292b4fbaa2aSMark F. Adams 	PetscViewer viewer; char fname[32];
293b4fbaa2aSMark F. Adams 	sprintf(fname,"part_mat_%d.mat",llev);
294b4fbaa2aSMark F. Adams 	PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
295b4fbaa2aSMark F. Adams 	ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
296b4fbaa2aSMark F. Adams 	ierr = PetscViewerDestroy( &viewer );
297b4fbaa2aSMark F. Adams       }
298b4fbaa2aSMark F. Adams 
299eb07cef2SMark F. Adams       ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
300eb07cef2SMark F. Adams 
301eb07cef2SMark F. Adams       ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
302eb07cef2SMark F. Adams     }
303f150b916SMark F. Adams 
304afc97cdcSMark F. Adams     if( ncrs0 != 0 ){
305b4fbaa2aSMark F. Adams       const PetscInt *is_idx;
306b4fbaa2aSMark F. Adams       MatPartitioning  mpart;
307a4b7d37bSMark F. Adams       IS proc_is;
3085ef31b24SMark F. Adams       /* hack to fix global data that pmetis.c uses in 'adj' */
30922063be5SMark F. Adams       for( nactive = jj = 0 ; jj < npe ; jj++ ) {
310afc97cdcSMark F. Adams 	if( counts[jj] != 0 ) {
31122063be5SMark F. Adams 	  adj->rmap->range[nactive++] = adj->rmap->range[jj];
312afc97cdcSMark F. Adams 	}
3135ef31b24SMark F. Adams       }
31422063be5SMark F. Adams       adj->rmap->range[nactive] = adj->rmap->range[npe];
3152f03bc48SMark F. Adams 
3168263b398SMark F. Adams       if( new_npe == 1 ) {
317a147abb0SMark F. Adams         /* simply dump on one proc w/o partitioner for final reduction */
3188263b398SMark F. Adams         ierr = MatGetLocalSize( adj, &is_sz, &ii );  CHKERRQ(ierr);
319a147abb0SMark F. Adams         ierr = ISCreateStride( cm, is_sz, targetPE, 0, &proc_is );  CHKERRQ(ierr);
32059a0be82SJed Brown       } else {
32159a0be82SJed Brown         char prefix[256];
32259a0be82SJed Brown         const char *pcpre;
32371d60426SBarry Smith         ierr = MatPartitioningCreate( cm, &mpart ); CHKERRQ(ierr);
3245ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
32559a0be82SJed Brown         ierr = PCGetOptionsPrefix(a_pc,&pcpre);CHKERRQ(ierr);
32659a0be82SJed Brown         ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
32759a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
32811e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
3295ef31b24SMark F. Adams         ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
330a4b7d37bSMark F. Adams         ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr);
33111e60469SMark F. Adams         ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
3328263b398SMark F. Adams       }
3335ef31b24SMark F. Adams       /* collect IS info */
334a4b7d37bSMark F. Adams       ierr = ISGetLocalSize( proc_is, &is_sz );       CHKERRQ(ierr);
3358263b398SMark F. Adams       ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr);
336a4b7d37bSMark F. Adams       ierr = ISGetIndices( proc_is, &is_idx );        CHKERRQ(ierr);
337a147abb0SMark F. Adams 
338a147abb0SMark F. Adams       NN = 1; /* bring to "front" of machine */
339a147abb0SMark F. Adams       /*NN = npe/new_npe;*/ /* spread partitioning across machine */
340b4fbaa2aSMark F. Adams       for( kk = jj = 0 ; kk < is_sz ; kk++ ){
341afc97cdcSMark F. Adams         for( ii = 0 ; ii < a_cbs ; ii++, jj++ ) {
3428263b398SMark F. Adams           newproc_idx[jj] = is_idx[kk] * NN; /* distribution */
343eb07cef2SMark F. Adams         }
3445ef31b24SMark F. Adams       }
345a4b7d37bSMark F. Adams       ierr = ISRestoreIndices( proc_is, &is_idx );     CHKERRQ(ierr);
346a4b7d37bSMark F. Adams       ierr = ISDestroy( &proc_is );                    CHKERRQ(ierr);
34771d60426SBarry Smith       ierr = MPI_Comm_free( &cm );              CHKERRQ(ierr);
348b4fbaa2aSMark F. Adams 
349b4fbaa2aSMark F. Adams       is_sz *= a_cbs;
3505ef31b24SMark F. Adams     }
3515ef31b24SMark F. Adams     else{
3528263b398SMark F. Adams       newproc_idx = 0;
3535ef31b24SMark F. Adams       is_sz = 0;
3545ef31b24SMark F. Adams     }
355f150b916SMark F. Adams 
3565ef31b24SMark F. Adams     ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
3578263b398SMark F. Adams     ierr = ISCreateGeneral( wcomm, is_sz, newproc_idx, PETSC_COPY_VALUES, &isnewproc );
3588263b398SMark F. Adams     if( newproc_idx != 0 ) {
3598263b398SMark F. Adams       ierr = PetscFree( newproc_idx );  CHKERRQ(ierr);
3605ef31b24SMark F. Adams     }
361e33ef3b1SMark F. Adams 
36211e60469SMark F. Adams     /*
36311e60469SMark F. Adams      Create an index set from the isnewproc index set to indicate the mapping TO
36411e60469SMark F. Adams      */
36511e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
36611e60469SMark F. Adams     /*
36711e60469SMark F. Adams      Determine how many elements are assigned to each processor
36811e60469SMark F. Adams      */
369fd3c6afaSMark F. Adams     inpe = npe;
370fd3c6afaSMark F. Adams     ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr);
37111e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
372038e3b61SMark F. Adams     ncrs_new = counts[mype]/a_cbs;
37322063be5SMark F. Adams     strideNew = ncrs_new*a_ndata_rows;
374b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
375b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
376b4fbaa2aSMark F. Adams #endif
37722063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
37811e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
379d3d6bff4SMark F. Adams     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
380470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
38111e60469SMark F. Adams     /*
382d3d6bff4SMark F. Adams      There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
383d3d6bff4SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
38411e60469SMark F. Adams      */
38592a756f0SMark F. Adams     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
38611e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
387038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
388d3d6bff4SMark F. Adams       PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
389d3d6bff4SMark F. Adams       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
39011e60469SMark F. Adams     }
391038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
392d3d6bff4SMark F. Adams     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
39311e60469SMark F. Adams     CHKERRQ(ierr);
39492a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
39511e60469SMark F. Adams     /*
39611e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
39711e60469SMark F. Adams      */
398d3d6bff4SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
399d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
400d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs0 ; ii++) {
401d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
402d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
403676e1743SMark F. Adams           PetscScalar tt = (PetscScalar)data[ix];
404676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
405d3d6bff4SMark F. Adams 	}
406038e3b61SMark F. Adams       }
407eb07cef2SMark F. Adams     }
408eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
409eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
41011e60469SMark F. Adams     /*
41111e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
41211e60469SMark F. Adams       to the correct processor
41311e60469SMark F. Adams     */
41411e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
41511e60469SMark F. Adams     CHKERRQ(ierr);
41611e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
41711e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
41811e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
41911e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
42011e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
42111e60469SMark F. Adams     /*
42211e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
42311e60469SMark F. Adams     */
424eb07cef2SMark F. Adams     ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
425d3d6bff4SMark F. Adams     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
426eb07cef2SMark F. Adams 
42711e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
428eb07cef2SMark F. Adams     data = *a_coarse_data;
429d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
430d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs_new ; ii++) {
431d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
432d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
433d3d6bff4SMark F. Adams 	  data[ix] = PetscRealPart(array[jx]);
434d3d6bff4SMark F. Adams 	  array[jx] = 1.e300;
435d3d6bff4SMark F. Adams 	}
436038e3b61SMark F. Adams       }
437038e3b61SMark F. Adams     }
43811e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
43911e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
44011e60469SMark F. Adams     /*
44111e60469SMark F. Adams       Invert for MatGetSubMatrix
44211e60469SMark F. Adams     */
443038e3b61SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
44411e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
44511e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
44611e60469SMark F. Adams     /* A_crs output */
447038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
44811e60469SMark F. Adams     CHKERRQ(ierr);
449eb07cef2SMark F. Adams 
450038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
451e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
452038e3b61SMark F. Adams     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
453eb07cef2SMark F. Adams 
45411e60469SMark F. Adams     /* prolongator */
45511e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
45611e60469SMark F. Adams     {
45711e60469SMark F. Adams       IS findices;
45811e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
45996434bc1SMark F. Adams #ifdef USE_R
46096434bc1SMark F. Adams       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
46196434bc1SMark F. Adams #else
46211e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
46396434bc1SMark F. Adams #endif
46411e60469SMark F. Adams       CHKERRQ(ierr);
465d61a3a7dSMark F. Adams 
46611e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
46711e60469SMark F. Adams     }
468d61a3a7dSMark F. Adams 
4693530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
470a147abb0SMark F. Adams     *a_P_inout = Pnew; /* output - repartitioned */
47192a756f0SMark F. Adams 
47211e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
47392a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
47492a756f0SMark F. Adams     ierr = PetscFree( ranks );  CHKERRQ(ierr);
475e33ef3b1SMark F. Adams   }
4765b89ad90SMark F. Adams 
4775b89ad90SMark F. Adams   PetscFunctionReturn(0);
4785b89ad90SMark F. Adams }
4795b89ad90SMark F. Adams 
4805b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4815b89ad90SMark F. Adams /*
4825b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4835b89ad90SMark F. Adams                     by setting data structures and options.
4845b89ad90SMark F. Adams 
4855b89ad90SMark F. Adams    Input Parameter:
4865b89ad90SMark F. Adams .  pc - the preconditioner context
4875b89ad90SMark F. Adams 
4885b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4895b89ad90SMark F. Adams 
4905b89ad90SMark F. Adams    Notes:
4915b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4925b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4935b89ad90SMark F. Adams */
4945b89ad90SMark F. Adams #undef __FUNCT__
4955b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
496eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
4975b89ad90SMark F. Adams {
4985b89ad90SMark F. Adams   PetscErrorCode  ierr;
499eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
5005b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
50158471d46SMark F. Adams   PC_MG_Levels   **mglevels = mg->levels;
502eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
503d3d6bff4SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
504eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
5053530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
506038e3b61SMark F. Adams   PetscBool        isOK;
507587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
508587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
509737a81a9SMark F. Adams   MatInfo          info;
5105ef31b24SMark F. Adams 
5115b89ad90SMark F. Adams   PetscFunctionBegin;
51203a628feSMark F. Adams   pc_gamg->m_count++;
513fca9ed99SMark F. Adams 
51458471d46SMark F. Adams   if( a_pc->setupcalled > 0 ) {
51503a628feSMark F. Adams     /* just do Galerkin grids */
51658471d46SMark F. Adams     Mat B,dA,dB;
51758471d46SMark F. Adams 
51858471d46SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
51958471d46SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
52058471d46SMark F. Adams 
521*5f8cf99dSMark F. Adams     if( pc_gamg->m_Nlevels > 1 ) {
52258471d46SMark F. Adams     /* currently only handle case where mat and pmat are the same on coarser levels */
52358471d46SMark F. Adams     ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
52458471d46SMark F. Adams     /* (re)set to get dirty flag */
52558471d46SMark F. Adams     ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
52658471d46SMark F. Adams     ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr);
52758471d46SMark F. Adams 
52858471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-2; level>-1; level--) {
52958471d46SMark F. Adams       ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
53003a628feSMark F. Adams       /* the first time through the matrix structure has changed from repartitioning */
53103a628feSMark F. Adams       if( pc_gamg->m_count == 2 ) {
53203a628feSMark F. Adams         ierr = MatDestroy( &B );  CHKERRQ(ierr);
53303a628feSMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
53403a628feSMark F. Adams         mglevels[level]->A = B;
53503a628feSMark F. Adams       }
53603a628feSMark F. Adams       else {
53758471d46SMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
53803a628feSMark F. Adams       }
53958471d46SMark F. Adams       ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
54058471d46SMark F. Adams       dB = B;
54158471d46SMark F. Adams       /* setup KSP/PC */
54258471d46SMark F. Adams       ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
54358471d46SMark F. Adams     }
544*5f8cf99dSMark F. Adams     }
545f6536408SMark F. Adams #define PRINT_MATS PETSC_FALSE
54658471d46SMark F. Adams     /* plot levels - A */
54758471d46SMark F. Adams     if( PRINT_MATS ) {
54858471d46SMark F. Adams       for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){
54958471d46SMark F. Adams         PetscViewer viewer;
55058471d46SMark F. Adams         char fname[32]; KSP smoother; Mat Tmat, TTm;
55158471d46SMark F. Adams         ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
55258471d46SMark F. Adams         ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr);
55303a628feSMark F. Adams         sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
55458471d46SMark F. Adams         ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
55558471d46SMark F. Adams         ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
55658471d46SMark F. Adams         ierr = MatView( Tmat, viewer ); CHKERRQ(ierr);
55758471d46SMark F. Adams         ierr = PetscViewerDestroy( &viewer );
55858471d46SMark F. Adams       }
55958471d46SMark F. Adams     }
56058471d46SMark F. Adams 
56103a628feSMark F. Adams     a_pc->setupcalled = 2;
56203a628feSMark F. Adams 
56358471d46SMark F. Adams     PetscFunctionReturn(0);
564eb07cef2SMark F. Adams   }
565f6536408SMark F. Adams 
566baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
567baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
5685b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
5695b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
5703530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
571eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
572eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
573eb07cef2SMark F. Adams 
574038e3b61SMark F. Adams   /* get data of not around */
575038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
576038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
577eb07cef2SMark F. Adams   }
578eb07cef2SMark F. Adams   data = pc_gamg->m_data;
579038e3b61SMark F. Adams 
580eb07cef2SMark F. Adams   /* Get A_i and R_i */
581737a81a9SMark F. Adams   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
582bed7c2b9SMark F. Adams   if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
5832c047e2dSMark F. Adams 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,
5842c047e2dSMark F. Adams 	      (int)(info.nz_used/(PetscReal)N),npe);
5858f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
5864ef23d27SMark F. Adams         level < (pc_gamg->m_Nlevels-1) && (level==0 || M>pc_gamg->m_coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
5870205a208SMark F. Adams         level++ ){
5885b89ad90SMark F. Adams     level1 = level + 1;
589b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
590b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
591b4fbaa2aSMark F. Adams #endif
592b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
593b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
594b4fbaa2aSMark F. Adams #endif
595389730f3SMark F. Adams     ierr = createProlongation(Aarr[level], data, level, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level], pc_gamg );
5965b89ad90SMark F. Adams     CHKERRQ(ierr);
597d3d6bff4SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
598b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
599b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
6005c7a89a6SBarry Smith 
601b4fbaa2aSMark F. Adams #endif
602baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
603baa4e9faSMark F. Adams     if( isOK ) {
604b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
605b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
606b4fbaa2aSMark F. Adams #endif
6078263b398SMark F. Adams       ierr = createLevel( a_pc, Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs,
6088263b398SMark F. Adams                           (PetscBool)(level == pc_gamg->m_Nlevels-2),
609389730f3SMark F. Adams                           &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
6103530afc2SMark F. Adams       CHKERRQ(ierr);
611b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
612b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
613b4fbaa2aSMark F. Adams #endif
6143530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
615737a81a9SMark F. Adams       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
616bed7c2b9SMark F. Adams       if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
6172c047e2dSMark F. Adams 		  mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols,
6182c047e2dSMark F. Adams 		  (int)(info.nz_used/(PetscReal)N),nactivepe);
619e33ef3b1SMark F. Adams       /* coarse grids with SA can have zero row/cols from singleton aggregates */
620737a81a9SMark F. Adams 
62181708dccSMark F. Adams       /* stop if one node */
62281708dccSMark F. Adams       if( M/pc_gamg->m_data_cols < 2 ) {
62381708dccSMark F. Adams         level++;
62481708dccSMark F. Adams         break;
62581708dccSMark F. Adams       }
62681708dccSMark F. Adams 
62781708dccSMark F. Adams       if (PETSC_FALSE) {
628785cba28SMark F. Adams         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
629785cba28SMark F. Adams         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
630e33ef3b1SMark F. Adams         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
631785cba28SMark F. Adams         nloceq = Iend-Istart;
632e33ef3b1SMark F. Adams         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
633e33ef3b1SMark F. Adams         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
634e33ef3b1SMark F. Adams         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
635785cba28SMark F. Adams         for(kk=0;kk<nloceq;kk++){
636e33ef3b1SMark F. Adams           if(data_arr[kk]==0.0) {
637785cba28SMark F. Adams             id = kk + Istart;
638e33ef3b1SMark F. Adams             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
639e33ef3b1SMark F. Adams             CHKERRQ(ierr);
640ecd57171SMark F. Adams             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1);
641e33ef3b1SMark F. Adams           }
642e33ef3b1SMark F. Adams         }
643e33ef3b1SMark F. Adams         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
644e33ef3b1SMark F. Adams         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
645e33ef3b1SMark F. Adams         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
646e33ef3b1SMark F. Adams         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
647e33ef3b1SMark F. Adams       }
648486a8d0bSJed Brown     } else {
649be544d3cSMark F. Adams       coarse_data = 0;
650baa4e9faSMark F. Adams       break;
651baa4e9faSMark F. Adams     }
652eb07cef2SMark F. Adams     data = coarse_data;
653b4fbaa2aSMark F. Adams 
654b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
655b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
656b4fbaa2aSMark F. Adams #endif
6575b89ad90SMark F. Adams   }
658be544d3cSMark F. Adams   if( coarse_data ) {
659eb07cef2SMark F. Adams     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
660be544d3cSMark F. Adams   }
661bed7c2b9SMark F. Adams   if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
6625b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
6635b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
6645b89ad90SMark F. Adams   fine_level = level;
665eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
6665b89ad90SMark F. Adams 
667*5f8cf99dSMark F. Adams   if( pc_gamg->m_Nlevels > 1 ) { /* don't setup MG if  */
668fc4362bfSMark F. Adams   /* set default smoothers */
669587fa25dSMark F. Adams   for ( lidx = 1, level = pc_gamg->m_Nlevels-2;
670587fa25dSMark F. Adams         lidx <= fine_level;
671587fa25dSMark F. Adams         lidx++, level--) {
6725745e0f5SMark F. Adams     PetscReal emax, emin;
6735b89ad90SMark F. Adams     KSP smoother; PC subpc;
674f6536408SMark F. Adams     PetscBool isCheb;
675f6536408SMark F. Adams     /* set defaults */
676587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
6775b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
678f6536408SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
679f6536408SMark F. Adams     /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */
680f6536408SMark F. Adams     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
681f6536408SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
682f6536408SMark F. Adams     /* overide defaults with input parameters */
683f6536408SMark F. Adams     ierr = KSPSetFromOptions( smoother ); CHKERRQ(ierr);
684f6536408SMark F. Adams 
685f6536408SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
686f6536408SMark F. Adams     /* do my own cheby */
687f6536408SMark F. Adams     ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr);
68841139db5SMark F. Adams 
689f6536408SMark F. Adams     if( isCheb ) {
690f6536408SMark F. Adams       ierr = PetscTypeCompare( (PetscObject)subpc, PETSC_GAMG_SMOOTHER, &isCheb ); CHKERRQ(ierr);
691f6536408SMark F. Adams       if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
692587fa25dSMark F. Adams       else{ /* eigen estimate 'emax' */
693587fa25dSMark F. Adams         KSP eksp; Mat Lmat = Aarr[level];
694fc4362bfSMark F. Adams         Vec bb, xx; PC pc;
695f6536408SMark F. Adams         const PCType type;
696038e3b61SMark F. Adams 
697f6536408SMark F. Adams         ierr = PCGetType( subpc, &type );   CHKERRQ(ierr);
6985745e0f5SMark F. Adams         ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
6995745e0f5SMark F. Adams         ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
700fc4362bfSMark F. Adams         {
701fc4362bfSMark F. Adams           PetscRandom    rctx;
702fc4362bfSMark F. Adams           ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
703fc4362bfSMark F. Adams           ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
704fc4362bfSMark F. Adams           ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
705fc4362bfSMark F. Adams           ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
7065b89ad90SMark F. Adams         }
707fc4362bfSMark F. Adams         ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
70874050f66SMark F. Adams         /* ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr); */
709038e3b61SMark F. Adams         ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
710fc4362bfSMark F. Adams         CHKERRQ(ierr);
711fc4362bfSMark F. Adams         ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
712fc4362bfSMark F. Adams 
713f6536408SMark F. Adams         ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
714f6536408SMark F. Adams         ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
715f6536408SMark F. Adams 
716f6536408SMark F. Adams         ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
717f6536408SMark F. Adams         ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
718f6536408SMark F. Adams         ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
719f6536408SMark F. Adams         ierr = PCSetType( pc, type ); CHKERRQ(ierr); /* should be same as eigen estimates op. */
720f6536408SMark F. Adams 
721fc4362bfSMark F. Adams         ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
722fc4362bfSMark F. Adams         ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
723fc4362bfSMark F. Adams         ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
724fc4362bfSMark F. Adams         ierr = VecDestroy( &xx );       CHKERRQ(ierr);
725fc4362bfSMark F. Adams         ierr = VecDestroy( &bb );       CHKERRQ(ierr);
726fc4362bfSMark F. Adams         ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
727f6536408SMark F. Adams 
728f6536408SMark F. Adams         if (pc_gamg->m_verbose) {
72941139db5SMark F. Adams           PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s PC setup max eigen=%e min=%e\n",
73041139db5SMark F. Adams                       __FUNCT__,emax,emin);
731f6536408SMark F. Adams         }
732fc4362bfSMark F. Adams       }
733038e3b61SMark F. Adams       {
734038e3b61SMark F. Adams         PetscInt N1, N0, tt;
735038e3b61SMark F. Adams         ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
736038e3b61SMark F. Adams         ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
737f6536408SMark F. Adams         /* heuristic - is this crap? */
738f6536408SMark F. Adams         emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
739038e3b61SMark F. Adams         emax *= 1.05;
740038e3b61SMark F. Adams       }
741fc4362bfSMark F. Adams       ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
742f6536408SMark F. Adams     }
7435745e0f5SMark F. Adams   }
7445745e0f5SMark F. Adams   {
745ecd57171SMark F. Adams     /* coarse grid */
746ecd57171SMark F. Adams     KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
7475745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
748eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
74958471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
7505745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
751ecd57171SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
752ecd57171SMark F. Adams     ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
753ecd57171SMark F. Adams     ierr = PCSetUp( subpc ); CHKERRQ(ierr);
754ecd57171SMark F. Adams     ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
755ecd57171SMark F. Adams     assert(ii==1);
756ecd57171SMark F. Adams     ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
757ecd57171SMark F. Adams     ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
758fc4362bfSMark F. Adams   }
759737a81a9SMark F. Adams 
760fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
761eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
7625b89ad90SMark F. Adams   {
7635b89ad90SMark F. Adams     PetscBool galerkin;
764eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
7655b89ad90SMark F. Adams     if(galerkin){
7665b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
7675b89ad90SMark F. Adams     }
7685b89ad90SMark F. Adams   }
7695745e0f5SMark F. Adams 
77058471d46SMark F. Adams   /* plot levels - R/P */
77158471d46SMark F. Adams   if( PRINT_MATS ) {
77258471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-1;level>0;level--){
77358471d46SMark F. Adams       PetscViewer viewer;
77458471d46SMark F. Adams       char fname[32];
77503a628feSMark F. Adams       sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
77611e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
7775b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
778038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
7795b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
78003a628feSMark F. Adams       sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
781e33ef3b1SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
782e33ef3b1SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
783e33ef3b1SMark F. Adams       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
784e33ef3b1SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
7855b89ad90SMark F. Adams     }
78658471d46SMark F. Adams   }
78758471d46SMark F. Adams 
78858471d46SMark F. Adams   /* set interpolation between the levels, clean up */
78958471d46SMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
79058471d46SMark F. Adams        lidx<fine_level;
79158471d46SMark F. Adams        lidx++, level--){
79258471d46SMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
793587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
794587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
7955b89ad90SMark F. Adams   }
7965745e0f5SMark F. Adams 
7975b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
798eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
799eb07cef2SMark F. Adams   ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
80003a628feSMark F. Adams   a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
80158471d46SMark F. Adams 
80258471d46SMark F. Adams   {
80358471d46SMark F. Adams     KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
80458471d46SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
80558471d46SMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
80658471d46SMark F. Adams     ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
80758471d46SMark F. Adams   }
808*5f8cf99dSMark F. Adams   }
809*5f8cf99dSMark F. Adams   else {
810*5f8cf99dSMark F. Adams     KSP smoother;
811*5f8cf99dSMark F. Adams     if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__);
812*5f8cf99dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
813*5f8cf99dSMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
814*5f8cf99dSMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
815*5f8cf99dSMark F. Adams     /* setupcalled is set to 0 so that MG is setup from scratch */
816*5f8cf99dSMark F. Adams     a_pc->setupcalled = 0;
817*5f8cf99dSMark F. Adams     ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
818*5f8cf99dSMark F. Adams     a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
819*5f8cf99dSMark F. Adams   }
8205b89ad90SMark F. Adams   PetscFunctionReturn(0);
8215b89ad90SMark F. Adams }
8225b89ad90SMark F. Adams 
823eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8245b89ad90SMark F. Adams /*
8255b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8265b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8275b89ad90SMark F. Adams 
8285b89ad90SMark F. Adams    Input Parameter:
8295b89ad90SMark F. Adams .  pc - the preconditioner context
8305b89ad90SMark F. Adams 
8315b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8325b89ad90SMark F. Adams */
8335b89ad90SMark F. Adams #undef __FUNCT__
8345b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
8355b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8365b89ad90SMark F. Adams {
8375b89ad90SMark F. Adams   PetscErrorCode  ierr;
8385b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
8395b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
8405b89ad90SMark F. Adams 
8415b89ad90SMark F. Adams   PetscFunctionBegin;
8425b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8435b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8445b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8455b89ad90SMark F. Adams   PetscFunctionReturn(0);
8465b89ad90SMark F. Adams }
8475b89ad90SMark F. Adams 
848676e1743SMark F. Adams 
849676e1743SMark F. Adams #undef __FUNCT__
850676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
851676e1743SMark F. Adams /*@
852676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
853676e1743SMark F. Adams    processor reduction.
854676e1743SMark F. Adams 
855676e1743SMark F. Adams    Not Collective on PC
856676e1743SMark F. Adams 
857676e1743SMark F. Adams    Input Parameters:
858676e1743SMark F. Adams .  pc - the preconditioner context
859676e1743SMark F. Adams 
860676e1743SMark F. Adams 
861676e1743SMark F. Adams    Options Database Key:
862676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
863676e1743SMark F. Adams 
864676e1743SMark F. Adams    Level: intermediate
865676e1743SMark F. Adams 
866676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
867676e1743SMark F. Adams 
868676e1743SMark F. Adams .seealso: ()
869676e1743SMark F. Adams @*/
870676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
871676e1743SMark F. Adams {
872676e1743SMark F. Adams   PetscErrorCode ierr;
873676e1743SMark F. Adams 
874676e1743SMark F. Adams   PetscFunctionBegin;
875676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
876676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
877676e1743SMark F. Adams   PetscFunctionReturn(0);
878676e1743SMark F. Adams }
879676e1743SMark F. Adams 
880676e1743SMark F. Adams EXTERN_C_BEGIN
881676e1743SMark F. Adams #undef __FUNCT__
882676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
883676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
884676e1743SMark F. Adams {
885c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
886c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
887676e1743SMark F. Adams 
888676e1743SMark F. Adams   PetscFunctionBegin;
889c20e4228SMark F. Adams   if(n>0) pc_gamg->m_min_eq_proc = n;
890676e1743SMark F. Adams   PetscFunctionReturn(0);
891676e1743SMark F. Adams }
892676e1743SMark F. Adams EXTERN_C_END
893676e1743SMark F. Adams 
894676e1743SMark F. Adams #undef __FUNCT__
895389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
896389730f3SMark F. Adams /*@
897389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
898389730f3SMark F. Adams 
899389730f3SMark F. Adams  Collective on PC
900389730f3SMark F. Adams 
901389730f3SMark F. Adams    Input Parameters:
902389730f3SMark F. Adams .  pc - the preconditioner context
903389730f3SMark F. Adams 
904389730f3SMark F. Adams 
905389730f3SMark F. Adams    Options Database Key:
906389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
907389730f3SMark F. Adams 
908389730f3SMark F. Adams    Level: intermediate
909389730f3SMark F. Adams 
910389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
911389730f3SMark F. Adams 
912389730f3SMark F. Adams .seealso: ()
913389730f3SMark F. Adams  @*/
914389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
915389730f3SMark F. Adams {
916389730f3SMark F. Adams   PetscErrorCode ierr;
917389730f3SMark F. Adams 
918389730f3SMark F. Adams   PetscFunctionBegin;
919389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
920389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
921389730f3SMark F. Adams   PetscFunctionReturn(0);
922389730f3SMark F. Adams }
923389730f3SMark F. Adams 
924389730f3SMark F. Adams EXTERN_C_BEGIN
925389730f3SMark F. Adams #undef __FUNCT__
926389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
927389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
928389730f3SMark F. Adams {
929389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
930389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
931389730f3SMark F. Adams 
932389730f3SMark F. Adams   PetscFunctionBegin;
933389730f3SMark F. Adams   if(n>0) pc_gamg->m_coarse_eq_limit = n;
934389730f3SMark F. Adams   PetscFunctionReturn(0);
935389730f3SMark F. Adams }
936389730f3SMark F. Adams EXTERN_C_END
937389730f3SMark F. Adams 
938389730f3SMark F. Adams #undef __FUNCT__
9398263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
940676e1743SMark F. Adams /*@
9418263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
942676e1743SMark F. Adams 
943676e1743SMark F. Adams    Collective on PC
944676e1743SMark F. Adams 
945676e1743SMark F. Adams    Input Parameters:
946676e1743SMark F. Adams .  pc - the preconditioner context
947676e1743SMark F. Adams 
948676e1743SMark F. Adams 
949676e1743SMark F. Adams    Options Database Key:
9508263b398SMark F. Adams .  -pc_gamg_repartition
951676e1743SMark F. Adams 
952676e1743SMark F. Adams    Level: intermediate
953676e1743SMark F. Adams 
954676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
955676e1743SMark F. Adams 
956676e1743SMark F. Adams .seealso: ()
957676e1743SMark F. Adams @*/
9588263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
959676e1743SMark F. Adams {
960676e1743SMark F. Adams   PetscErrorCode ierr;
961676e1743SMark F. Adams 
962676e1743SMark F. Adams   PetscFunctionBegin;
963676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9648263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
965676e1743SMark F. Adams   PetscFunctionReturn(0);
966676e1743SMark F. Adams }
967676e1743SMark F. Adams 
968676e1743SMark F. Adams EXTERN_C_BEGIN
969676e1743SMark F. Adams #undef __FUNCT__
9708263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
9718263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
972676e1743SMark F. Adams {
973c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
974c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
975676e1743SMark F. Adams 
976676e1743SMark F. Adams   PetscFunctionBegin;
9778263b398SMark F. Adams   pc_gamg->m_repart = n;
978676e1743SMark F. Adams   PetscFunctionReturn(0);
979676e1743SMark F. Adams }
980676e1743SMark F. Adams EXTERN_C_END
981676e1743SMark F. Adams 
982676e1743SMark F. Adams #undef __FUNCT__
9834ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
9844ef23d27SMark F. Adams /*@
9854ef23d27SMark F. Adams    PCGAMGSetNlevels -
9864ef23d27SMark F. Adams 
9874ef23d27SMark F. Adams    Not collective on PC
9884ef23d27SMark F. Adams 
9894ef23d27SMark F. Adams    Input Parameters:
9904ef23d27SMark F. Adams .  pc - the preconditioner context
9914ef23d27SMark F. Adams 
9924ef23d27SMark F. Adams 
9934ef23d27SMark F. Adams    Options Database Key:
9944ef23d27SMark F. Adams .  -pc_mg_levels
9954ef23d27SMark F. Adams 
9964ef23d27SMark F. Adams    Level: intermediate
9974ef23d27SMark F. Adams 
9984ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
9994ef23d27SMark F. Adams 
10004ef23d27SMark F. Adams .seealso: ()
10014ef23d27SMark F. Adams @*/
10024ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
10034ef23d27SMark F. Adams {
10044ef23d27SMark F. Adams   PetscErrorCode ierr;
10054ef23d27SMark F. Adams 
10064ef23d27SMark F. Adams   PetscFunctionBegin;
10074ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10084ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
10094ef23d27SMark F. Adams   PetscFunctionReturn(0);
10104ef23d27SMark F. Adams }
10114ef23d27SMark F. Adams 
10124ef23d27SMark F. Adams EXTERN_C_BEGIN
10134ef23d27SMark F. Adams #undef __FUNCT__
10144ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
10154ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
10164ef23d27SMark F. Adams {
10174ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
10184ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
10194ef23d27SMark F. Adams 
10204ef23d27SMark F. Adams   PetscFunctionBegin;
10214ef23d27SMark F. Adams   pc_gamg->m_Nlevels = n;
10224ef23d27SMark F. Adams   PetscFunctionReturn(0);
10234ef23d27SMark F. Adams }
10244ef23d27SMark F. Adams EXTERN_C_END
10254ef23d27SMark F. Adams 
10264ef23d27SMark F. Adams #undef __FUNCT__
10273542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
10283542efc5SMark F. Adams /*@
10293542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
10303542efc5SMark F. Adams 
10313542efc5SMark F. Adams    Not collective on PC
10323542efc5SMark F. Adams 
10333542efc5SMark F. Adams    Input Parameters:
10343542efc5SMark F. Adams .  pc - the preconditioner context
10353542efc5SMark F. Adams 
10363542efc5SMark F. Adams 
10373542efc5SMark F. Adams    Options Database Key:
10383542efc5SMark F. Adams .  -pc_gamg_threshold
10393542efc5SMark F. Adams 
10403542efc5SMark F. Adams    Level: intermediate
10413542efc5SMark F. Adams 
10423542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
10433542efc5SMark F. Adams 
10443542efc5SMark F. Adams .seealso: ()
10453542efc5SMark F. Adams @*/
10463542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
10473542efc5SMark F. Adams {
10483542efc5SMark F. Adams   PetscErrorCode ierr;
10493542efc5SMark F. Adams 
10503542efc5SMark F. Adams   PetscFunctionBegin;
10513542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10523542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
10533542efc5SMark F. Adams   PetscFunctionReturn(0);
10543542efc5SMark F. Adams }
10553542efc5SMark F. Adams 
10563542efc5SMark F. Adams EXTERN_C_BEGIN
10573542efc5SMark F. Adams #undef __FUNCT__
10583542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
10593542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
10603542efc5SMark F. Adams {
1061c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1062c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
10633542efc5SMark F. Adams 
10643542efc5SMark F. Adams   PetscFunctionBegin;
1065c20e4228SMark F. Adams   pc_gamg->m_threshold = n;
10663542efc5SMark F. Adams   PetscFunctionReturn(0);
10673542efc5SMark F. Adams }
10683542efc5SMark F. Adams EXTERN_C_END
10693542efc5SMark F. Adams 
10703542efc5SMark F. Adams #undef __FUNCT__
1071676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType"
1072676e1743SMark F. Adams /*@
1073676e1743SMark F. Adams    PCGAMGSetSolverType - Set solution method.
1074676e1743SMark F. Adams 
1075676e1743SMark F. Adams    Collective on PC
1076676e1743SMark F. Adams 
1077676e1743SMark F. Adams    Input Parameters:
1078676e1743SMark F. Adams .  pc - the preconditioner context
1079676e1743SMark F. Adams 
1080676e1743SMark F. Adams 
1081676e1743SMark F. Adams    Options Database Key:
10823542efc5SMark F. Adams .  -pc_gamg_type
1083676e1743SMark F. Adams 
1084676e1743SMark F. Adams    Level: intermediate
1085676e1743SMark F. Adams 
1086676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1087676e1743SMark F. Adams 
1088676e1743SMark F. Adams .seealso: ()
1089676e1743SMark F. Adams @*/
1090676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz )
1091676e1743SMark F. Adams {
1092676e1743SMark F. Adams   PetscErrorCode ierr;
1093676e1743SMark F. Adams 
1094676e1743SMark F. Adams   PetscFunctionBegin;
1095676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1096676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz));
1097676e1743SMark F. Adams   CHKERRQ(ierr);
1098676e1743SMark F. Adams   PetscFunctionReturn(0);
1099676e1743SMark F. Adams }
1100676e1743SMark F. Adams 
1101676e1743SMark F. Adams EXTERN_C_BEGIN
1102676e1743SMark F. Adams #undef __FUNCT__
1103676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG"
1104676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz )
1105676e1743SMark F. Adams {
1106676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1107676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1108676e1743SMark F. Adams 
1109676e1743SMark F. Adams   PetscFunctionBegin;
1110676e1743SMark F. Adams   if(sz < 64) strcpy(pc_gamg->m_type,str);
1111676e1743SMark F. Adams   PetscFunctionReturn(0);
1112676e1743SMark F. Adams }
1113676e1743SMark F. Adams EXTERN_C_END
1114676e1743SMark F. Adams 
11155b89ad90SMark F. Adams #undef __FUNCT__
11165b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
11175b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
11185b89ad90SMark F. Adams {
1119676e1743SMark F. Adams   PetscErrorCode  ierr;
1120676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1121676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1122676e1743SMark F. Adams   PetscBool flag;
11235b89ad90SMark F. Adams 
11245b89ad90SMark F. Adams   PetscFunctionBegin;
1125676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1126676e1743SMark F. Adams   {
112775b74bdaSMark F. Adams 
112875b74bdaSMark F. Adams     pc_gamg->m_method = 1; /* default to plane aggregation */
1129676e1743SMark F. Adams     ierr = PetscOptionsString("-pc_gamg_type",
113075b74bdaSMark F. Adams                               "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid ('geo')",
1131676e1743SMark F. Adams                               "PCGAMGSetSolverType",
113275b74bdaSMark F. Adams                               "pa",
1133676e1743SMark F. Adams                               pc_gamg->m_type,
1134676e1743SMark F. Adams                               64,
1135676e1743SMark F. Adams                               &flag );
1136676e1743SMark F. Adams     CHKERRQ(ierr);
1137d8c859f0SMark F. Adams     if( flag && pc_gamg->m_data != 0 ) {
1138d8c859f0SMark F. Adams       if( (strcmp(pc_gamg->m_type,"sa")==0 && pc_gamg->m_method != 2) ||
1139d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"pa")==0 && pc_gamg->m_method != 1) ||
1140d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"geo")==0 && pc_gamg->m_method != 0) ) {
1141d8c859f0SMark F. Adams         SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "PCSetFromOptions called of PCSetCoordinates (with new method, after data was created)");
1142d8c859f0SMark F. Adams       }
1143d8c859f0SMark F. Adams     }
1144bed7c2b9SMark F. Adams 
1145f6536408SMark F. Adams     if (flag ) {
1146f6536408SMark F. Adams       if( strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2;
1147f6536408SMark F. Adams       else if( strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1;
1148f6536408SMark F. Adams       else if( strcmp(pc_gamg->m_type,"geo") == 0) pc_gamg->m_method = 0;
1149f6536408SMark F. Adams       else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Invalid gamg type: %s",pc_gamg->m_type);
1150f6536408SMark F. Adams     }
115175b74bdaSMark F. Adams 
115275b74bdaSMark F. Adams     /* -pc_gamg_verbose */
115375b74bdaSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG","none",pc_gamg->m_verbose,&pc_gamg->m_verbose,PETSC_NULL);CHKERRQ(ierr);
115475b74bdaSMark F. Adams 
11558263b398SMark F. Adams     /* -pc_gamg_repartition */
11568263b398SMark F. Adams     pc_gamg->m_repart = PETSC_FALSE;
11578263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
11588263b398SMark F. Adams                             "Repartion coarse grids (false)",
11598263b398SMark F. Adams                             "PCGAMGRepartitioning",
11608263b398SMark F. Adams                             pc_gamg->m_repart,
11618263b398SMark F. Adams                             &pc_gamg->m_repart,
1162676e1743SMark F. Adams                             &flag);
1163676e1743SMark F. Adams     CHKERRQ(ierr);
1164676e1743SMark F. Adams 
1165c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1166c20e4228SMark F. Adams     pc_gamg->m_min_eq_proc = 600;
1167676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1168676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1169676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
1170c20e4228SMark F. Adams                            pc_gamg->m_min_eq_proc,
1171c20e4228SMark F. Adams                            &pc_gamg->m_min_eq_proc,
1172676e1743SMark F. Adams                            &flag );
1173676e1743SMark F. Adams     CHKERRQ(ierr);
11743542efc5SMark F. Adams 
1175389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1176389730f3SMark F. Adams     pc_gamg->m_coarse_eq_limit = 1500;
1177389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1178389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1179389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
1180389730f3SMark F. Adams                            pc_gamg->m_coarse_eq_limit,
1181389730f3SMark F. Adams                            &pc_gamg->m_coarse_eq_limit,
1182389730f3SMark F. Adams                            &flag );
1183389730f3SMark F. Adams     CHKERRQ(ierr);
1184389730f3SMark F. Adams 
1185c20e4228SMark F. Adams     /* -pc_gamg_threshold */
1186c20e4228SMark F. Adams     pc_gamg->m_threshold = 0.05;
11873542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
11883542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
11893542efc5SMark F. Adams                             "PCGAMGSetThreshold",
1190c20e4228SMark F. Adams                             pc_gamg->m_threshold,
1191c20e4228SMark F. Adams                             &pc_gamg->m_threshold,
11923542efc5SMark F. Adams                             &flag );
11933542efc5SMark F. Adams     CHKERRQ(ierr);
11944ef23d27SMark F. Adams 
11954ef23d27SMark F. Adams     pc_gamg->m_Nlevels = GAMG_MAXLEVELS;
11964ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
11974ef23d27SMark F. Adams                            "Set number of MG levels",
11984ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
11994ef23d27SMark F. Adams                            pc_gamg->m_Nlevels,
12004ef23d27SMark F. Adams                            &pc_gamg->m_Nlevels,
12014ef23d27SMark F. Adams                            &flag );
1202676e1743SMark F. Adams   }
1203676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1204676e1743SMark F. Adams 
12055b89ad90SMark F. Adams   PetscFunctionReturn(0);
12065b89ad90SMark F. Adams }
12075b89ad90SMark F. Adams 
12085b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
12095b89ad90SMark F. Adams /*MC
1210dc76db92SMark F. Adams      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning. This preconditioner currently has two
1211dc76db92SMark F. Adams            AMG methods: 1) an unstructured geometric method, which requires that you provide coordinates for each
1212dc76db92SMark F. Adams            vertex, and 2) smoothed aggregation.  Smoothed aggregation (SA) can work without coordinates but it
1213dc76db92SMark F. Adams            will generate some common non-trivial null spaces if coordinates are provided.  The input fine grid matrix
1214dc76db92SMark F. Adams            must have the block size set for 'system' problems (with multiple dof per vertex/cell) to work properly.
1215dc76db92SMark F. Adams            SA will generate rotational rigid body mode null space vectors, in addition to the trivial translational
1216dc76db92SMark F. Adams            modes, when coordinates are provide in 2D and 3D.
12175b89ad90SMark F. Adams 
1218280d9858SJed Brown    Options Database Keys:
12195b89ad90SMark F. Adams    Multigrid options(inherited)
1220280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1221280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1222280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1223280d9858SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
12245b89ad90SMark F. Adams 
12255b89ad90SMark F. Adams   Level: intermediate
1226280d9858SJed Brown 
12275b89ad90SMark F. Adams   Concepts: multigrid
12285b89ad90SMark F. Adams 
12295b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1230280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
12315b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
12325b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12335b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12345b89ad90SMark F. Adams M*/
12355b89ad90SMark F. Adams 
12365b89ad90SMark F. Adams EXTERN_C_BEGIN
12375b89ad90SMark F. Adams #undef __FUNCT__
12385b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
12395b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
12405b89ad90SMark F. Adams {
12415b89ad90SMark F. Adams   PetscErrorCode  ierr;
12425b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
12435b89ad90SMark F. Adams   PC_MG           *mg;
12445ef31b24SMark F. Adams   PetscClassId     cookie;
12455ee9c036SSatish Balay #if defined PETSC_USE_LOG
12465ee9c036SSatish Balay   static int count = 0;
12475ee9c036SSatish Balay #endif
12485b89ad90SMark F. Adams 
12495b89ad90SMark F. Adams   PetscFunctionBegin;
12505b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
12515b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
12525b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
12535b89ad90SMark F. Adams 
12545b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
12555b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
125603a628feSMark F. Adams   pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0;
12575b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
12585b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
12595b89ad90SMark F. Adams 
12605b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
12615b89ad90SMark F. Adams 
12625b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
12635b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
12645b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
12655b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
12665b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
12675b89ad90SMark F. Adams 
12685b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12695b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
12705b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
1271c97e1df0SMark F. Adams 					    PCSetCoordinates_GAMG);
1272c97e1df0SMark F. Adams   CHKERRQ(ierr);
1273c97e1df0SMark F. Adams 
1274676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1275676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1276676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1277676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1278676e1743SMark F. Adams   CHKERRQ(ierr);
1279676e1743SMark F. Adams 
1280676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1281389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_C",
1282389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_GAMG",
1283389730f3SMark F. Adams 					    PCGAMGSetCoarseEqLim_GAMG);
1284389730f3SMark F. Adams   CHKERRQ(ierr);
1285389730f3SMark F. Adams 
1286389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12878263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_C",
12888263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_GAMG",
12898263b398SMark F. Adams 					    PCGAMGSetRepartitioning_GAMG);
1290676e1743SMark F. Adams   CHKERRQ(ierr);
1291676e1743SMark F. Adams 
1292676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12933542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
12943542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
12953542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
12963542efc5SMark F. Adams   CHKERRQ(ierr);
12973542efc5SMark F. Adams 
12983542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1299676e1743SMark F. Adams 					    "PCGAMGSetSolverType_C",
1300676e1743SMark F. Adams 					    "PCGAMGSetSolverType_GAMG",
1301676e1743SMark F. Adams 					    PCGAMGSetSolverType_GAMG);
1302676e1743SMark F. Adams   CHKERRQ(ierr);
1303c97e1df0SMark F. Adams 
1304b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
1305785cba28SMark F. Adams   if( count++ == 0 ) {
13065ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
1307b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
13082a44bfbaSMark F. Adams     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
13092a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]);
13102a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]);
13112a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]);
1312b4fbaa2aSMark F. Adams     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1313b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1314b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1315b4fbaa2aSMark F. Adams     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1316b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_events[SET7]);
1317b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]); */
1318b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1319b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1320b4fbaa2aSMark F. Adams     PetscLogEventRegister("  PL repartition", cookie, &gamg_setup_events[SET12]);
1321f852f58cSMark F. Adams 
1322fca9ed99SMark F. Adams     /* PetscLogEventRegister("  PL 1", cookie, &gamg_setup_events[SET13]); */
1323fca9ed99SMark F. Adams     /* PetscLogEventRegister("  PL 2", cookie, &gamg_setup_events[SET14]); */
1324fca9ed99SMark F. Adams     /* PetscLogEventRegister("  PL 3", cookie, &gamg_setup_events[SET15]); */
1325f852f58cSMark F. Adams 
1326b4fbaa2aSMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1327b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1328b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
13295ef31b24SMark F. Adams 
1330b4fbaa2aSMark F. Adams     /* create timer stages */
1331b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1332b4fbaa2aSMark F. Adams     {
1333b4fbaa2aSMark F. Adams       char str[32];
1334b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1335b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1336b4fbaa2aSMark F. Adams       PetscInt lidx;
1337b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1338b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1339b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1340b4fbaa2aSMark F. Adams       }
1341b4fbaa2aSMark F. Adams     }
1342b4fbaa2aSMark F. Adams #endif
1343b4fbaa2aSMark F. Adams   }
1344b4fbaa2aSMark F. Adams #endif
13455b89ad90SMark F. Adams   PetscFunctionReturn(0);
13465b89ad90SMark F. Adams }
13475b89ad90SMark F. Adams EXTERN_C_END
1348