xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 0e1b4bd6f46fbc30cc04dd3bbd1f4e6bbef12688)
15b89ad90SMark F. Adams /*
25b89ad90SMark F. Adams  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
45ef31b24SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>
5f96513f1SMatthew G Knepley 
6f96513f1SMatthew G Knepley PetscLogEvent gamg_setup_stages[NUM_SET];
7f96513f1SMatthew G Knepley 
85b89ad90SMark F. Adams /* Private context for the GAMG preconditioner */
95b89ad90SMark F. Adams typedef struct gamg_TAG {
105b89ad90SMark F. Adams   PetscInt       m_dim;
115b89ad90SMark F. Adams   PetscInt       m_Nlevels;
125b89ad90SMark F. Adams   PetscInt       m_data_sz;
13d3d6bff4SMark F. Adams   PetscInt       m_data_rows;
14d3d6bff4SMark F. Adams   PetscInt       m_data_cols;
15d3d6bff4SMark F. Adams   PetscBool      m_useSA;
165b89ad90SMark F. Adams   PetscReal     *m_data; /* blocked vector of vertex data on fine grid (coordinates) */
175b89ad90SMark F. Adams } PC_GAMG;
185b89ad90SMark F. Adams 
195b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
205b89ad90SMark F. Adams /*
215b89ad90SMark F. Adams    PCSetCoordinates_GAMG
225b89ad90SMark F. Adams 
235b89ad90SMark F. Adams    Input Parameter:
245b89ad90SMark F. Adams    .  pc - the preconditioner context
255b89ad90SMark F. Adams */
26a92563c5SMark F. Adams EXTERN_C_BEGIN
275b89ad90SMark F. Adams #undef __FUNCT__
285b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG"
29eb07cef2SMark F. Adams PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords )
305b89ad90SMark F. Adams {
31eb07cef2SMark F. Adams   PC_MG          *mg = (PC_MG*)a_pc->data;
325b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
336c237d78SBarry Smith   PetscErrorCode ierr;
34d3d6bff4SMark F. Adams   PetscInt       arrsz,bs,my0,kk,ii,jj,nloc,Iend;
35038e3b61SMark F. Adams   Mat            Amat = a_pc->pmat;
36d3d6bff4SMark F. Adams   PetscBool      flag;
37785cba28SMark F. Adams   char           str[64];
385b89ad90SMark F. Adams 
395b89ad90SMark F. Adams   PetscFunctionBegin;
40038e3b61SMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
41d3d6bff4SMark F. Adams   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
42d3d6bff4SMark F. Adams   nloc = (Iend-my0)/bs;
43d3d6bff4SMark F. Adams   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
44038e3b61SMark F. Adams 
45785cba28SMark F. Adams   ierr  = PetscOptionsGetString(PETSC_NULL,"-pc_gamg_type",str,64,&flag);    CHKERRQ( ierr );
46d3d6bff4SMark F. Adams   pc_gamg->m_useSA = (PetscBool)(flag && strcmp(str,"sa") == 0);
47038e3b61SMark F. Adams 
48d3d6bff4SMark F. Adams   pc_gamg->m_data_rows = 1;
49d3d6bff4SMark F. Adams   if(a_coords == 0) pc_gamg->m_useSA = PETSC_TRUE; /* use SA if no data */
50d3d6bff4SMark F. Adams   if( !pc_gamg->m_useSA ) pc_gamg->m_data_cols = a_ndm; /* coordinates */
51038e3b61SMark F. Adams   else{ /* SA: null space vectors */
52d3d6bff4SMark F. Adams     if(a_coords != 0 && bs==1 ) pc_gamg->m_data_cols = 1; /* scalar w/ coords and SA (not needed) */
53d3d6bff4SMark F. Adams     else if(a_coords != 0) pc_gamg->m_data_cols = (a_ndm==2 ? 3 : 6); /* elasticity */
54d3d6bff4SMark F. Adams     else pc_gamg->m_data_cols = bs; /* no data, force SA with constant null space vectors */
55d3d6bff4SMark F. Adams     pc_gamg->m_data_rows = bs;
56038e3b61SMark F. Adams   }
57d3d6bff4SMark F. Adams   arrsz = nloc*pc_gamg->m_data_rows*pc_gamg->m_data_cols;
585b89ad90SMark F. Adams 
59038e3b61SMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
606c237d78SBarry Smith   if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) {
615b89ad90SMark F. Adams     ierr = PetscFree( pc_gamg->m_data );  CHKERRQ(ierr);
62eb07cef2SMark F. Adams     ierr = PetscMalloc((arrsz+1)*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr);
635b89ad90SMark F. Adams   }
64038e3b61SMark F. Adams   for(kk=0;kk<arrsz;kk++)pc_gamg->m_data[kk] = -999.;
65eb07cef2SMark F. Adams   pc_gamg->m_data[arrsz] = -99.;
66038e3b61SMark F. Adams   /* copy data in - column oriented */
67d3d6bff4SMark F. Adams   if( pc_gamg->m_useSA ) {
68d3d6bff4SMark F. Adams     const PetscInt M = Iend - my0;
69038e3b61SMark F. Adams     for(kk=0;kk<nloc;kk++){
70038e3b61SMark F. Adams       PetscReal *data = &pc_gamg->m_data[kk*bs];
71d3d6bff4SMark F. Adams       if( pc_gamg->m_data_cols==1 ) *data = 1.0;
72038e3b61SMark F. Adams       else {
73038e3b61SMark F. Adams         for(ii=0;ii<bs;ii++)
74038e3b61SMark F. Adams 	  for(jj=0;jj<bs;jj++)
75038e3b61SMark F. Adams 	    if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */
76038e3b61SMark F. Adams 	    else data[ii*M + jj] = 0.0;
77eb07cef2SMark F. Adams         if( a_coords != 0 ) {
78038e3b61SMark F. Adams           if( a_ndm == 2 ){ /* rotational modes */
79038e3b61SMark F. Adams             data += 2*M;
80038e3b61SMark F. Adams             data[0] = -a_coords[2*kk+1];
81038e3b61SMark F. Adams             data[1] =  a_coords[2*kk];
825b89ad90SMark F. Adams           }
83eb07cef2SMark F. Adams           else {
84038e3b61SMark F. Adams             data += 3*M;
85038e3b61SMark F. Adams             data[0] = 0.0;               data[M+0] =  a_coords[3*kk+2]; data[2*M+0] = -a_coords[3*kk+1];
86038e3b61SMark F. Adams             data[1] = -a_coords[3*kk+2]; data[M+1] = 0.0;               data[2*M+1] =  a_coords[3*kk];
87038e3b61SMark F. Adams             data[2] =  a_coords[3*kk+1]; data[M+2] = -a_coords[3*kk];   data[2*M+2] = 0.0;
88038e3b61SMark F. Adams           }
89eb07cef2SMark F. Adams         }
90eb07cef2SMark F. Adams       }
91eb07cef2SMark F. Adams     }
92eb07cef2SMark F. Adams   }
93eb07cef2SMark F. Adams   else {
94038e3b61SMark F. Adams     for( kk = 0 ; kk < nloc ; kk++ ){
95038e3b61SMark F. Adams       for( ii = 0 ; ii < a_ndm ; ii++ ) {
96038e3b61SMark F. Adams         pc_gamg->m_data[ii*nloc + kk] =  a_coords[kk*a_ndm + ii];
97eb07cef2SMark F. Adams       }
98eb07cef2SMark F. Adams     }
99038e3b61SMark F. Adams   }
100eb07cef2SMark F. Adams   assert(pc_gamg->m_data[arrsz] == -99.);
101038e3b61SMark F. Adams 
1025b89ad90SMark F. Adams   pc_gamg->m_data_sz = arrsz;
103eb07cef2SMark F. Adams   pc_gamg->m_dim = a_ndm;
1045b89ad90SMark F. Adams 
1055b89ad90SMark F. Adams   PetscFunctionReturn(0);
1065b89ad90SMark F. Adams }
107a92563c5SMark F. Adams EXTERN_C_END
1085b89ad90SMark F. Adams 
109d3d6bff4SMark F. Adams 
110d3d6bff4SMark F. Adams /* -----------------------------------------------------------------------------*/
111d3d6bff4SMark F. Adams #undef __FUNCT__
112d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
113d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
114d3d6bff4SMark F. Adams {
115d3d6bff4SMark F. Adams   PetscErrorCode  ierr;
116d3d6bff4SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
117d3d6bff4SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
118d3d6bff4SMark F. Adams 
119d3d6bff4SMark F. Adams   PetscFunctionBegin;
120d3d6bff4SMark F. Adams   ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr);
121d3d6bff4SMark F. Adams   pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0;
122d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
123d3d6bff4SMark F. Adams }
124d3d6bff4SMark F. Adams 
1255b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1265b89ad90SMark F. Adams /*
12711e60469SMark F. Adams    partitionLevel
1285b89ad90SMark F. Adams 
1295b89ad90SMark F. Adams    Input Parameter:
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)
1333530afc2SMark F. Adams    In/Output Parameter:
1343530afc2SMark F. Adams    . a_P_inout - prolongation operator to the next level (k-1)
135eb07cef2SMark F. Adams    . a_coarse_data - data that need to be moved
1363530afc2SMark F. Adams    . a_active_proc - number of active procs
13711e60469SMark F. Adams    Output Parameter:
1383530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1395b89ad90SMark F. Adams */
1405cb416c2SMark F. Adams 
1415b89ad90SMark F. Adams #undef __FUNCT__
14211e60469SMark F. Adams #define __FUNCT__ "partitionLevel"
1433530afc2SMark F. Adams PetscErrorCode partitionLevel( Mat a_Amat_fine,
144d3d6bff4SMark F. Adams                                PetscInt a_ndata_rows,
145d3d6bff4SMark F. Adams                                PetscInt a_ndata_cols,
146038e3b61SMark F. Adams 			       PetscInt a_cbs,
1473530afc2SMark F. Adams                                Mat *a_P_inout,
148eb07cef2SMark F. Adams                                PetscReal **a_coarse_data,
1493530afc2SMark F. Adams                                PetscMPIInt *a_active_proc,
1503530afc2SMark F. Adams                                Mat *a_Amat_crs
15111e60469SMark F. Adams                                )
1525b89ad90SMark F. Adams {
1535b89ad90SMark F. Adams   PetscErrorCode   ierr;
154038e3b61SMark F. Adams   Mat              Cmat,Pnew,Pold=*a_P_inout;
15511e60469SMark F. Adams   IS               new_indices,isnum;
1563530afc2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
1575ef31b24SMark F. Adams   PetscMPIInt      nactive_procs,mype,npe;
158038e3b61SMark F. Adams   PetscInt         Istart,Iend,Istart0,Iend0,ncrs0,ncrs_new,fbs;
159a6828334SMark F. Adams   PetscInt         neq,NN;
160e33ef3b1SMark F. Adams   PetscMPIInt      new_npe,targ_npe;
161737a81a9SMark F. Adams   PetscBool        flag = PETSC_FALSE;
1625b89ad90SMark F. Adams 
1635b89ad90SMark F. Adams   PetscFunctionBegin;
16411e60469SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
16511e60469SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
166038e3b61SMark F. Adams   ierr = MatGetBlockSize( a_Amat_fine, &fbs ); CHKERRQ(ierr);
16711e60469SMark F. Adams   /* RAP */
168038e3b61SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
169737a81a9SMark F. Adams 
170038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
171acadaa72SMark F. Adams 
172038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
173038e3b61SMark F. Adams   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
174038e3b61SMark F. Adams 
175038e3b61SMark F. Adams   /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
176038e3b61SMark F. Adams   ierr = MatGetSize( Cmat, &neq, &NN );CHKERRQ(ierr);
177*0e1b4bd6SMark F. Adams #define MIN_EQ_PROC 2000
1785ef31b24SMark F. Adams   nactive_procs = *a_active_proc;
1793530afc2SMark F. Adams   targ_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */
180*0e1b4bd6SMark F. Adams #define TOP_GRID_LIM 2000
1815ef31b24SMark F. Adams   if( targ_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1; /* chop coarsest grid */
1825ef31b24SMark F. Adams   else if (targ_npe >= nactive_procs ) new_npe = nactive_procs; /* no change */
1833530afc2SMark F. Adams   else {
1845ef31b24SMark F. Adams     PetscMPIInt factstart,fact;
1853530afc2SMark F. Adams     new_npe = -9999;
1865ef31b24SMark F. Adams     factstart = nactive_procs;
1873530afc2SMark F. Adams     for(fact=factstart;fact>0;fact--){ /* try to find a better number of procs */
1885ef31b24SMark F. Adams       if( nactive_procs%fact==0 && neq/(nactive_procs/fact) > MIN_EQ_PROC ) {
1895ef31b24SMark F. Adams         new_npe = nactive_procs/fact;
1903530afc2SMark F. Adams       }
1913530afc2SMark F. Adams     }
192038e3b61SMark F. Adams     if(new_npe == -9999) new_npe = nactive_procs;
1933530afc2SMark F. Adams   }
194e33ef3b1SMark F. Adams 
195737a81a9SMark F. Adams   ierr  = PetscOptionsHasName(PETSC_NULL,"-pc_gamg_avoid_repartitioning",&flag);
196737a81a9SMark F. Adams   CHKERRQ( ierr );
197737a81a9SMark F. Adams   if( !flag ) { /* re-partition */
19811e60469SMark F. Adams     MatPartitioning  mpart;
1995ef31b24SMark F. Adams     Mat              adj;
2005ef31b24SMark F. Adams     const PetscInt  *is_idx;
201d3d6bff4SMark F. Adams     PetscInt         is_sz,kk,jj,ii,old_fact=(npe/nactive_procs), *isnewproc_idx;
202c9a0b8beSMark F. Adams     /* create sub communicator  */
203c9a0b8beSMark F. Adams     MPI_Comm         cm,new_comm;
204d3d6bff4SMark F. Adams     PetscInt         membershipKey = mype%old_fact, new_fact=(npe/new_npe),counts[npe];
205e33ef3b1SMark F. Adams     IS               isnewproc;
206e33ef3b1SMark F. Adams 
207d3d6bff4SMark F. Adams     *a_active_proc = new_npe; /* output for next level */
208c9a0b8beSMark F. Adams     ierr = MPI_Comm_split(wcomm, membershipKey, mype, &cm); CHKERRQ(ierr);
209c9a0b8beSMark F. Adams     ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr);
210c9a0b8beSMark F. Adams     ierr = MPI_Comm_free( &cm );                            CHKERRQ(ierr);
211c9a0b8beSMark F. Adams 
2125ef31b24SMark F. Adams     /* MatPartitioningApply call MatConvert, which is collective */
213737a81a9SMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_stages[SET12],0,0,0,0);CHKERRQ(ierr);
214038e3b61SMark F. Adams     if( a_cbs==1) {
215038e3b61SMark F. Adams       ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
216eb07cef2SMark F. Adams     }
217eb07cef2SMark F. Adams     else {
218038e3b61SMark F. Adams       /* make a scalar matrix to partition */
219eb07cef2SMark F. Adams       Mat tMat;
220eb07cef2SMark F. Adams       PetscInt Ii,ncols; const PetscScalar *vals; const PetscInt *idx;
2216876a03eSMark F. Adams       MatInfo info;
2226876a03eSMark F. Adams       ierr = MatGetInfo(Cmat,MAT_LOCAL,&info); CHKERRQ(ierr);
22355ea7f60SMark F. Adams       ncols = (PetscInt)info.nz_used/((ncrs0+1)*a_cbs*a_cbs)+1;
2246876a03eSMark F. Adams 
225eb07cef2SMark F. Adams       ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
226eb07cef2SMark F. Adams                               PETSC_DETERMINE, PETSC_DETERMINE,
2276876a03eSMark F. Adams                               2*ncols, PETSC_NULL, ncols, PETSC_NULL,
228eb07cef2SMark F. Adams                               &tMat );
2296876a03eSMark F. Adams       CHKERRQ(ierr);
230eb07cef2SMark F. Adams 
231eb07cef2SMark F. Adams       for ( Ii = Istart0; Ii < Iend0; Ii++ ) {
232038e3b61SMark F. Adams         PetscInt dest_row = Ii/a_cbs;
233038e3b61SMark F. Adams         ierr = MatGetRow(Cmat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr);
234eb07cef2SMark F. Adams         for( jj = 0 ; jj < ncols ; jj++ ){
235038e3b61SMark F. Adams           PetscInt dest_col = idx[jj]/a_cbs;
236eb07cef2SMark F. Adams           PetscScalar v = 1.0;
237eb07cef2SMark F. Adams           ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
238eb07cef2SMark F. Adams         }
239038e3b61SMark F. Adams         ierr = MatRestoreRow(Cmat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr);
240eb07cef2SMark F. Adams       }
241eb07cef2SMark F. Adams       ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
242eb07cef2SMark F. Adams       ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
243eb07cef2SMark F. Adams 
244eb07cef2SMark F. Adams       ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
245eb07cef2SMark F. Adams 
246eb07cef2SMark F. Adams       ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
247eb07cef2SMark F. Adams     }
2485ef31b24SMark F. Adams     if( membershipKey == 0 ) {
249737a81a9SMark F. Adams       if(ncrs0==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"zero local nodes -- increase min");
2505ef31b24SMark F. Adams       /* hack to fix global data that pmetis.c uses in 'adj' */
2515ef31b24SMark F. Adams       for( kk=0 , jj=0 ; kk<=npe ; jj++, kk += old_fact ) {
2525ef31b24SMark F. Adams         adj->rmap->range[jj] = adj->rmap->range[kk];
2535ef31b24SMark F. Adams       }
2545ef31b24SMark F. Adams       ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr);
2555ef31b24SMark F. Adams       ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
25611e60469SMark F. Adams       ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr);
2575ef31b24SMark F. Adams       ierr = MatPartitioningSetNParts( mpart, new_npe ); CHKERRQ(ierr);
25811e60469SMark F. Adams       ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr);
25911e60469SMark F. Adams       ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr);
2605ef31b24SMark F. Adams       /* collect IS info */
2615ef31b24SMark F. Adams       ierr = ISGetLocalSize( isnewproc, &is_sz );        CHKERRQ(ierr);
262038e3b61SMark F. Adams       ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr);
2635ef31b24SMark F. Adams       ierr = ISGetIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
264eb07cef2SMark F. Adams       /* spread partitioning across machine - probably the right thing to do but machine spec. */
265eb07cef2SMark F. Adams       for(kk=0,jj=0;kk<is_sz;kk++){
266038e3b61SMark F. Adams         for(ii=0 ; ii<a_cbs ; ii++, jj++ ) { /* expand for equation level by 'bs' */
267eb07cef2SMark F. Adams           isnewproc_idx[jj] = is_idx[kk] * new_fact;
268eb07cef2SMark F. Adams         }
2695ef31b24SMark F. Adams       }
2705ef31b24SMark F. Adams       ierr = ISRestoreIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
2715ef31b24SMark F. Adams       ierr = ISDestroy( &isnewproc );                    CHKERRQ(ierr);
272038e3b61SMark F. Adams       is_sz *= a_cbs;
2735ef31b24SMark F. Adams     }
2745ef31b24SMark F. Adams     else {
2755ef31b24SMark F. Adams       isnewproc_idx = 0;
2765ef31b24SMark F. Adams       is_sz = 0;
2775ef31b24SMark F. Adams     }
2785ef31b24SMark F. Adams     ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
2795ef31b24SMark F. Adams     ierr = MPI_Comm_free( &new_comm );    CHKERRQ(ierr);
2805ef31b24SMark F. Adams     ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc );
2815ef31b24SMark F. Adams     if( membershipKey == 0 ) {
2825ef31b24SMark F. Adams       ierr = PetscFree( isnewproc_idx );  CHKERRQ(ierr);
2835ef31b24SMark F. Adams     }
284e33ef3b1SMark F. Adams 
28511e60469SMark F. Adams     /*
28611e60469SMark F. Adams      Create an index set from the isnewproc index set to indicate the mapping TO
28711e60469SMark F. Adams      */
28811e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
28911e60469SMark F. Adams     /*
29011e60469SMark F. Adams      Determine how many elements are assigned to each processor
29111e60469SMark F. Adams      */
29211e60469SMark F. Adams     ierr = ISPartitioningCount( isnewproc, npe, counts ); CHKERRQ(ierr);
29311e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
294038e3b61SMark F. Adams     ncrs_new = counts[mype]/a_cbs;
295737a81a9SMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_stages[SET12],0,0,0,0);   CHKERRQ(ierr);
2965ef31b24SMark F. Adams 
29711e60469SMark F. Adams     { /* Create a vector to contain the newly ordered element information */
298d3d6bff4SMark F. Adams       const PetscInt *idx, data_sz=a_ndata_rows*a_ndata_cols;
299d3d6bff4SMark F. Adams       const PetscInt  stride0=ncrs0*a_ndata_rows,strideNew=ncrs_new*a_ndata_rows;
300038e3b61SMark F. Adams       PetscInt        ii,jj,kk;
30111e60469SMark F. Adams       IS              isscat;
30211e60469SMark F. Adams       PetscScalar    *array;
30311e60469SMark F. Adams       Vec             src_crd, dest_crd;
304eb07cef2SMark F. Adams       PetscReal      *data = *a_coarse_data;
30511e60469SMark F. Adams       VecScatter      vecscat;
306d3d6bff4SMark F. Adams       PetscInt        tidx[ncrs0*data_sz];
307acadaa72SMark F. Adams 
30811e60469SMark F. Adams       ierr = VecCreate( wcomm, &dest_crd );
309d3d6bff4SMark F. Adams       ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
31011e60469SMark F. Adams       ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/
31111e60469SMark F. Adams       /*
312d3d6bff4SMark F. Adams        There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
313d3d6bff4SMark F. Adams        a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
31411e60469SMark F. Adams        */
31511e60469SMark F. Adams       ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
316038e3b61SMark F. Adams       for(ii=0,jj=0; ii<ncrs0 ; ii++) {
317d3d6bff4SMark F. Adams         PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
318d3d6bff4SMark F. Adams         for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
31911e60469SMark F. Adams       }
320038e3b61SMark F. Adams       ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
321d3d6bff4SMark F. Adams       ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
32211e60469SMark F. Adams       CHKERRQ(ierr);
32311e60469SMark F. Adams       /*
32411e60469SMark F. Adams        Create a vector to contain the original vertex information for each element
32511e60469SMark F. Adams        */
326d3d6bff4SMark F. Adams       ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
327d3d6bff4SMark F. Adams       for( jj=0; jj<a_ndata_cols ; jj++ ) {
328d3d6bff4SMark F. Adams         for( ii=0 ; ii<ncrs0 ; ii++) {
329d3d6bff4SMark F. Adams           for( kk=0; kk<a_ndata_rows ; kk++ ) {
330d3d6bff4SMark F. Adams             PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
331d3d6bff4SMark F. Adams             ierr = VecSetValues( src_crd, 1, &jx, &data[ix], INSERT_VALUES );  CHKERRQ(ierr);
332d3d6bff4SMark F. Adams           }
333038e3b61SMark F. Adams         }
334eb07cef2SMark F. Adams       }
335eb07cef2SMark F. Adams       ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
336eb07cef2SMark F. Adams       ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
33711e60469SMark F. Adams       /*
33811e60469SMark F. Adams        Scatter the element vertex information (still in the original vertex ordering)
33911e60469SMark F. Adams        to the correct processor
34011e60469SMark F. Adams        */
34111e60469SMark F. Adams       ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
34211e60469SMark F. Adams       CHKERRQ(ierr);
34311e60469SMark F. Adams       ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
34411e60469SMark F. Adams       ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
34511e60469SMark F. Adams       ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
34611e60469SMark F. Adams       ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
34711e60469SMark F. Adams       ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
34811e60469SMark F. Adams       /*
34911e60469SMark F. Adams        Put the element vertex data into a new allocation of the gdata->ele
35011e60469SMark F. Adams        */
351eb07cef2SMark F. Adams       ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
352d3d6bff4SMark F. Adams       ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
353eb07cef2SMark F. Adams 
35411e60469SMark F. Adams       ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
355eb07cef2SMark F. Adams       data = *a_coarse_data;
356d3d6bff4SMark F. Adams       for( jj=0; jj<a_ndata_cols ; jj++ ) {
357d3d6bff4SMark F. Adams         for( ii=0 ; ii<ncrs_new ; ii++) {
358d3d6bff4SMark F. Adams           for( kk=0; kk<a_ndata_rows ; kk++ ) {
359d3d6bff4SMark F. Adams             PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
360d3d6bff4SMark F. Adams             data[ix] = PetscRealPart(array[jx]);
361d3d6bff4SMark F. Adams             array[jx] = 1.e300;
362d3d6bff4SMark F. Adams           }
363038e3b61SMark F. Adams         }
364038e3b61SMark F. Adams       }
36511e60469SMark F. Adams       ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
36611e60469SMark F. Adams       ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
36711e60469SMark F. Adams     }
368737a81a9SMark F. Adams 
36911e60469SMark F. Adams     /*
37011e60469SMark F. Adams      Invert for MatGetSubMatrix
37111e60469SMark F. Adams      */
372038e3b61SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
37311e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
37411e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
37511e60469SMark F. Adams     /* A_crs output */
376038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
37711e60469SMark F. Adams     CHKERRQ(ierr);
378eb07cef2SMark F. Adams 
379038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
380e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
381038e3b61SMark F. Adams     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
382eb07cef2SMark F. Adams 
38311e60469SMark F. Adams     /* prolongator */
38411e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
38511e60469SMark F. Adams     {
38611e60469SMark F. Adams       IS findices;
38711e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
38811e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
38911e60469SMark F. Adams       CHKERRQ(ierr);
39011e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
39111e60469SMark F. Adams     }
3923530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
3933530afc2SMark F. Adams     *a_P_inout = Pnew; /* output */
39411e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
395e33ef3b1SMark F. Adams   }
396e33ef3b1SMark F. Adams   else {
397e33ef3b1SMark F. Adams     *a_Amat_crs = Cmat; /* output */
398e33ef3b1SMark F. Adams   }
3995b89ad90SMark F. Adams 
4005b89ad90SMark F. Adams   PetscFunctionReturn(0);
4015b89ad90SMark F. Adams }
4025b89ad90SMark F. Adams 
4031de87581SMark F. Adams #define GAMG_MAXLEVELS 30
4045ef31b24SMark F. Adams #if defined(PETSC_USE_LOG)
4055ef31b24SMark F. Adams PetscLogStage  gamg_stages[20];
4065ef31b24SMark F. Adams #endif
4075b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4085b89ad90SMark F. Adams /*
4095b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4105b89ad90SMark F. Adams                     by setting data structures and options.
4115b89ad90SMark F. Adams 
4125b89ad90SMark F. Adams    Input Parameter:
4135b89ad90SMark F. Adams .  pc - the preconditioner context
4145b89ad90SMark F. Adams 
4155b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4165b89ad90SMark F. Adams 
4175b89ad90SMark F. Adams    Notes:
4185b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4195b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4205b89ad90SMark F. Adams */
4215b89ad90SMark F. Adams #undef __FUNCT__
4225b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
423eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
4245b89ad90SMark F. Adams {
4255b89ad90SMark F. Adams   PetscErrorCode  ierr;
426eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
4275b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
428eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
429d3d6bff4SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
430eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
4313530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
432038e3b61SMark F. Adams   PetscBool        isOK;
433587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
434587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
435737a81a9SMark F. Adams   MatInfo          info;
4365ef31b24SMark F. Adams 
4375b89ad90SMark F. Adams   PetscFunctionBegin;
438eb07cef2SMark F. Adams   if( a_pc->setupcalled ) {
439eb07cef2SMark F. Adams     /* no state data in GAMG to destroy */
440eb07cef2SMark F. Adams     ierr = PCReset_MG( a_pc ); CHKERRQ(ierr);
441eb07cef2SMark F. Adams   }
442baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
443baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
4445b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
4455b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
4463530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
447eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
448eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
449eb07cef2SMark F. Adams 
450038e3b61SMark F. Adams   /* get data of not around */
451038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
452038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
453eb07cef2SMark F. Adams   }
454eb07cef2SMark F. Adams   data = pc_gamg->m_data;
455038e3b61SMark F. Adams 
456eb07cef2SMark F. Adams   /* Get A_i and R_i */
457737a81a9SMark F. Adams   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
458737a81a9SMark F. Adams   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",
45904afdd00SMark F. Adams 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,(PetscInt)(info.nz_used/(PetscReal)N),npe);
4608f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
461785cba28SMark F. Adams         level < GAMG_MAXLEVELS-1 && (level==0 || M>TOP_GRID_LIM) && (npe==1 || nactivepe>1);
4620205a208SMark F. Adams         level++ ){
4635b89ad90SMark F. Adams     level1 = level + 1;
4645ef31b24SMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr);
465d3d6bff4SMark F. Adams     ierr = createProlongation(Aarr[level], data, pc_gamg->m_dim, pc_gamg->m_data_cols, pc_gamg->m_useSA,
466785cba28SMark F. Adams                               level, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level] );
4675b89ad90SMark F. Adams     CHKERRQ(ierr);
468d3d6bff4SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
4695ef31b24SMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr);
4705ef31b24SMark F. Adams 
471baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
472d3d6bff4SMark F. Adams 
473baa4e9faSMark F. Adams     if( isOK ) {
4745ef31b24SMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr);
475d3d6bff4SMark F. Adams       ierr = partitionLevel( Aarr[level], pc_gamg->m_useSA ? bs : 1, pc_gamg->m_data_cols, bs,
476eb07cef2SMark F. Adams                              &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
4773530afc2SMark F. Adams       CHKERRQ(ierr);
4785ef31b24SMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr);
4793530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
480737a81a9SMark F. Adams       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
481737a81a9SMark F. Adams       PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, bs=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
48204afdd00SMark F. Adams 		  mype,__FUNCT__,level1,N,bs,pc_gamg->m_data_cols,(PetscInt)(info.nz_used/(PetscReal)N),nactivepe);
483e33ef3b1SMark F. Adams       /* coarse grids with SA can have zero row/cols from singleton aggregates */
484e33ef3b1SMark F. Adams       /* aggregation method can probably gaurrentee this does not happen! - be safe for now */
485737a81a9SMark F. Adams 
486e33ef3b1SMark F. Adams       if( PETSC_TRUE ){
487785cba28SMark F. Adams         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
488785cba28SMark F. Adams         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
489e33ef3b1SMark F. Adams         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
490785cba28SMark F. Adams         nloceq = Iend-Istart;
491e33ef3b1SMark F. Adams         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
492e33ef3b1SMark F. Adams         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
493e33ef3b1SMark F. Adams         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
494785cba28SMark F. Adams         for(kk=0;kk<nloceq;kk++){
495e33ef3b1SMark F. Adams           if(data_arr[kk]==0.0) {
496785cba28SMark F. Adams             id = kk + Istart;
497e33ef3b1SMark F. Adams             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
498e33ef3b1SMark F. Adams             CHKERRQ(ierr);
499785cba28SMark F. Adams             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added diag to zero (%d) on level %d \n",mype,__FUNCT__,id,level);
500e33ef3b1SMark F. Adams           }
501e33ef3b1SMark F. Adams         }
502e33ef3b1SMark F. Adams         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
503e33ef3b1SMark F. Adams         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
504e33ef3b1SMark F. Adams         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
505e33ef3b1SMark F. Adams         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
506e33ef3b1SMark F. Adams       }
507baa4e9faSMark F. Adams     }
508baa4e9faSMark F. Adams     else{
509be544d3cSMark F. Adams       coarse_data = 0;
510baa4e9faSMark F. Adams       break;
511baa4e9faSMark F. Adams     }
512eb07cef2SMark F. Adams     data = coarse_data;
5135b89ad90SMark F. Adams   }
514be544d3cSMark F. Adams   if( coarse_data ) {
515eb07cef2SMark F. Adams     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
516be544d3cSMark F. Adams   }
517baa4e9faSMark F. Adams   PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
5185b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
5195b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
5205b89ad90SMark F. Adams   fine_level = level;
521eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
5225b89ad90SMark F. Adams 
523fc4362bfSMark F. Adams   /* set default smoothers */
524587fa25dSMark F. Adams   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
525587fa25dSMark F. Adams         lidx <= fine_level;
526587fa25dSMark F. Adams         lidx++, level--) {
5275745e0f5SMark F. Adams     PetscReal emax, emin;
5285b89ad90SMark F. Adams     KSP smoother; PC subpc;
529587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
5305b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
531038e3b61SMark F. Adams     if( emaxs[level] > 0.0 ) emax = emaxs[level];
532587fa25dSMark F. Adams     else{ /* eigen estimate 'emax' */
533587fa25dSMark F. Adams       KSP eksp; Mat Lmat = Aarr[level];
534fc4362bfSMark F. Adams       Vec bb, xx; PC pc;
535038e3b61SMark F. Adams 
5365745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
5375745e0f5SMark F. Adams       ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
538fc4362bfSMark F. Adams       {
539fc4362bfSMark F. Adams 	PetscRandom    rctx;
540fc4362bfSMark F. Adams 	ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
541fc4362bfSMark F. Adams 	ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
542fc4362bfSMark F. Adams 	ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
543fc4362bfSMark F. Adams 	ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
5445b89ad90SMark F. Adams       }
545fc4362bfSMark F. Adams       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
546e33ef3b1SMark F. Adams       ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr);
547fc4362bfSMark F. Adams       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
5485745e0f5SMark F. Adams       ierr = KSPSetOperators( eksp, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN ); CHKERRQ( ierr );
549fc4362bfSMark F. Adams       ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
550e33ef3b1SMark F. Adams       ierr = PCSetType( pc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr); /* should be same as above */
551038e3b61SMark F. Adams       ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
552fc4362bfSMark F. Adams       CHKERRQ(ierr);
553e33ef3b1SMark F. Adams       //ierr = KSPSetConvergenceTest( eksp, KSPSkipConverged, 0, 0 ); CHKERRQ(ierr);
554fc4362bfSMark F. Adams       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
555fc4362bfSMark F. Adams 
556fc4362bfSMark F. Adams       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
557fc4362bfSMark F. Adams       ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
558fc4362bfSMark F. Adams       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
559fc4362bfSMark F. Adams       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
560fc4362bfSMark F. Adams       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
561fc4362bfSMark F. Adams       ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
562fc4362bfSMark F. Adams     }
563038e3b61SMark F. Adams     {
564038e3b61SMark F. Adams       PetscInt N1, N0, tt;
565038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
566038e3b61SMark F. Adams       ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
567785cba28SMark F. Adams       emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */
568038e3b61SMark F. Adams       emax *= 1.05;
569785cba28SMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"%s max eigen = %e min = %e\n",__FUNCT__,emax,emin);
570038e3b61SMark F. Adams     }
571038e3b61SMark F. Adams 
572587fa25dSMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], DIFFERENT_NONZERO_PATTERN );
573fc4362bfSMark F. Adams     ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
574*0e1b4bd6SMark F. Adams     /*ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr);*/
5755745e0f5SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
576e33ef3b1SMark F. Adams     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
5775745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
5785745e0f5SMark F. Adams   }
5795745e0f5SMark F. Adams   {
5805745e0f5SMark F. Adams     KSP smoother; /* coarse grid */
5815745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
582eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
5835745e0f5SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN );
5845745e0f5SMark F. Adams     CHKERRQ(ierr);
5855745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
586fc4362bfSMark F. Adams   }
587737a81a9SMark F. Adams 
588fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
589eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
5905b89ad90SMark F. Adams   {
5915b89ad90SMark F. Adams     PetscBool galerkin;
592eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
5935b89ad90SMark F. Adams     if(galerkin){
5945b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
5955b89ad90SMark F. Adams     }
5965b89ad90SMark F. Adams   }
5975745e0f5SMark F. Adams 
5985745e0f5SMark F. Adams   /* set interpolation between the levels, create timer stages, clean up */
5998f4b7eb5SMark F. Adams   if( PETSC_FALSE ) {
6005ef31b24SMark F. Adams     char str[32];
6015ef31b24SMark F. Adams     sprintf(str,"MG Level %d (%d)",0,pc_gamg->m_Nlevels-1);
6025ef31b24SMark F. Adams     PetscLogStageRegister(str, &gamg_stages[fine_level]);
6035ef31b24SMark F. Adams   }
604587fa25dSMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
605587fa25dSMark F. Adams        lidx<fine_level;
606587fa25dSMark F. Adams        lidx++, level--){
607587fa25dSMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
6086c237d78SBarry Smith     if( !PETSC_TRUE ) {
60911e60469SMark F. Adams       PetscViewer viewer; char fname[32];
610038e3b61SMark F. Adams       sprintf(fname,"Pmat_%d.m",level);
61111e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
6125b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
613038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
6145b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
615e33ef3b1SMark F. Adams       sprintf(fname,"Amat_%d.m",level);
616e33ef3b1SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
617e33ef3b1SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
618e33ef3b1SMark F. Adams       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
619e33ef3b1SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
6205b89ad90SMark F. Adams     }
621587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
622587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
6238f4b7eb5SMark F. Adams     if( PETSC_FALSE ) {
6245ef31b24SMark F. Adams       char str[32];
625587fa25dSMark F. Adams       sprintf(str,"MG Level %d (%d)",lidx+1,level-1);
626587fa25dSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[level-1]);
627a92563c5SMark F. Adams     }
6285b89ad90SMark F. Adams   }
6295745e0f5SMark F. Adams 
6305b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
631eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
632eb07cef2SMark F. Adams   ierr = PCSetUp_MG(a_pc);CHKERRQ(ierr);
633e33ef3b1SMark F. Adams 
6345b89ad90SMark F. Adams   PetscFunctionReturn(0);
6355b89ad90SMark F. Adams }
6365b89ad90SMark F. Adams 
637eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
6385b89ad90SMark F. Adams /*
6395b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
6405b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
6415b89ad90SMark F. Adams 
6425b89ad90SMark F. Adams    Input Parameter:
6435b89ad90SMark F. Adams .  pc - the preconditioner context
6445b89ad90SMark F. Adams 
6455b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
6465b89ad90SMark F. Adams */
6475b89ad90SMark F. Adams #undef __FUNCT__
6485b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
6495b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
6505b89ad90SMark F. Adams {
6515b89ad90SMark F. Adams   PetscErrorCode  ierr;
6525b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
6535b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
6545b89ad90SMark F. Adams 
6555b89ad90SMark F. Adams   PetscFunctionBegin;
6565b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
6575b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
6585b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
6595b89ad90SMark F. Adams   PetscFunctionReturn(0);
6605b89ad90SMark F. Adams }
6615b89ad90SMark F. Adams 
6625b89ad90SMark F. Adams #undef __FUNCT__
6635b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
6645b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
6655b89ad90SMark F. Adams {
6665b89ad90SMark F. Adams   /* PetscErrorCode  ierr; */
6675b89ad90SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
6685b89ad90SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
6695b89ad90SMark F. Adams   /* MPI_Comm        comm = ((PetscObject)pc)->comm; */
6705b89ad90SMark F. Adams 
6715b89ad90SMark F. Adams   PetscFunctionBegin;
6725b89ad90SMark F. Adams   PetscFunctionReturn(0);
6735b89ad90SMark F. Adams }
6745b89ad90SMark F. Adams 
6755b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
6765b89ad90SMark F. Adams /*
6775b89ad90SMark F. Adams  PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG
6785b89ad90SMark F. Adams 
6795b89ad90SMark F. Adams    Input Parameter:
6805b89ad90SMark F. Adams .  pc - the preconditioner context
6815b89ad90SMark F. Adams 
6825b89ad90SMark F. Adams    Application Interface Routine: PCCreate()
6835b89ad90SMark F. Adams 
6845b89ad90SMark F. Adams   */
6855b89ad90SMark F. Adams  /* MC
6865b89ad90SMark F. Adams      PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide
6875b89ad90SMark F. Adams        fine grid discretization matrix and coordinates on the fine grid.
6885b89ad90SMark F. Adams 
6895b89ad90SMark F. Adams    Options Database Key:
6905b89ad90SMark F. Adams    Multigrid options(inherited)
6915b89ad90SMark F. Adams +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
6925b89ad90SMark F. Adams .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
6935b89ad90SMark F. Adams .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
6945b89ad90SMark F. Adams    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
6955b89ad90SMark F. Adams    GAMG options:
6965b89ad90SMark F. Adams 
6975b89ad90SMark F. Adams    Level: intermediate
6985b89ad90SMark F. Adams   Concepts: multigrid
6995b89ad90SMark F. Adams 
7005b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
7015b89ad90SMark F. Adams            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
7025b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
7035b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
7045b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
7055b89ad90SMark F. Adams M */
7065b89ad90SMark F. Adams 
7075b89ad90SMark F. Adams EXTERN_C_BEGIN
7085b89ad90SMark F. Adams #undef __FUNCT__
7095b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
7105b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
7115b89ad90SMark F. Adams {
7125b89ad90SMark F. Adams   PetscErrorCode  ierr;
7135b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
7145b89ad90SMark F. Adams   PC_MG           *mg;
7155ef31b24SMark F. Adams   PetscClassId     cookie;
7165b89ad90SMark F. Adams 
7175b89ad90SMark F. Adams   PetscFunctionBegin;
7185b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
7195b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
7205b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
7215b89ad90SMark F. Adams 
7225b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
7235b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
7245b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
7255b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
7265b89ad90SMark F. Adams 
7275b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
7285b89ad90SMark F. Adams 
7295b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
7305b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
7315b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
7325b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
7335b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
7345b89ad90SMark F. Adams 
7355b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
7365b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
7375b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
7385b89ad90SMark F. Adams 					    PCSetCoordinates_GAMG);CHKERRQ(ierr);
739785cba28SMark F. Adams   static int count = 0;
740785cba28SMark F. Adams   if( count++ == 0 ) {
7415ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
742737a81a9SMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_stages[SET1]);
743737a81a9SMark F. Adams     PetscLogEventRegister(" make graph", cookie, &gamg_setup_stages[SET3]);
744737a81a9SMark F. Adams     PetscLogEventRegister(" MIS/Agg", cookie, &gamg_setup_stages[SET4]);
745737a81a9SMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_stages[SET5]);
746737a81a9SMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_stages[SET6]);
747737a81a9SMark F. Adams     PetscLogEventRegister("   search & set", cookie, &gamg_setup_stages[FIND_V]);
748737a81a9SMark F. Adams     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_stages[SET7]);
749737a81a9SMark F. Adams     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_stages[SET8]); */
750737a81a9SMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_stages[SET9]);
751737a81a9SMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_stages[SET2]);
752737a81a9SMark F. Adams     PetscLogEventRegister(" PL repartition", cookie, &gamg_setup_stages[SET12]);
753737a81a9SMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_stages[SET13]); */
754737a81a9SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_stages[SET10]); */
755737a81a9SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_stages[SET11]); */
756785cba28SMark F. Adams   }
7575ef31b24SMark F. Adams 
7585b89ad90SMark F. Adams   PetscFunctionReturn(0);
7595b89ad90SMark F. Adams }
7605b89ad90SMark F. Adams EXTERN_C_END
761