xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision f150b916196796664ac83cc9dc3e953cdfd6292b)
15b89ad90SMark F. Adams /*
25b89ad90SMark F. Adams  GAMG geometric-algebric multiogrid 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 /*
125486a8d0bSJed Brown    PCGAMGPartitionLevel
1265b89ad90SMark F. Adams 
1275b89ad90SMark F. Adams    Input Parameter:
1283530afc2SMark F. Adams    . a_Amat_fine - matrix on this fine (k) level
129d3d6bff4SMark F. Adams    . a_ndata_rows - size of data to move (coarse grid)
130d3d6bff4SMark F. Adams    . a_ndata_cols - size of data to move (coarse grid)
131389730f3SMark F. Adams    . a_pc_gamg - parameters
1323530afc2SMark F. Adams    In/Output Parameter:
1333530afc2SMark F. Adams    . a_P_inout - prolongation operator to the next level (k-1)
134eb07cef2SMark F. Adams    . a_coarse_data - data that need to be moved
135afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
13611e60469SMark F. Adams    Output Parameter:
1373530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1385b89ad90SMark F. Adams */
1395cb416c2SMark F. Adams 
1405b89ad90SMark F. Adams #undef __FUNCT__
141486a8d0bSJed Brown #define __FUNCT__ "PCGAMGPartitionLevel"
142486a8d0bSJed Brown PetscErrorCode PCGAMGPartitionLevel(PC pc, Mat a_Amat_fine,
143d3d6bff4SMark F. Adams                                     PetscInt a_ndata_rows,
144d3d6bff4SMark F. Adams                                     PetscInt a_ndata_cols,
145038e3b61SMark F. Adams                                     PetscInt a_cbs,
1463530afc2SMark F. Adams                                     Mat *a_P_inout,
147eb07cef2SMark F. Adams                                     PetscReal **a_coarse_data,
148afc97cdcSMark F. Adams                                     PetscMPIInt *a_nactive_proc,
149389730f3SMark F. Adams                                     Mat *a_Amat_crs
15011e60469SMark F. Adams                                     )
1515b89ad90SMark F. Adams {
152486a8d0bSJed Brown   PC_MG           *mg = (PC_MG*)pc->data;
153486a8d0bSJed Brown   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
154389730f3SMark F. Adams   const PetscBool  avoid_repart = pc_gamg->m_avoid_repart;
155389730f3SMark F. Adams   const PetscInt   min_eq_proc = pc_gamg->m_min_eq_proc, coarse_max = pc_gamg->m_coarse_eq_limit;
1565b89ad90SMark F. Adams   PetscErrorCode   ierr;
157038e3b61SMark F. Adams   Mat              Cmat,Pnew,Pold=*a_P_inout;
15811e60469SMark F. Adams   IS               new_indices,isnum;
1593530afc2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
160fa31c87dSMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive;
161fa31c87dSMark F. Adams   PetscInt         neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0;
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);
166f6536408SMark F. Adams 
16711e60469SMark F. Adams   /* RAP */
16896434bc1SMark F. Adams #ifdef USE_R
16996434bc1SMark F. Adams   /* make R wih brute force for now */
17096434bc1SMark F. Adams   ierr = MatTranspose( Pold, Pnew );
17196434bc1SMark F. Adams   ierr = MatDestroy( &Pold );  CHKERRQ(ierr);
17296434bc1SMark F. Adams   ierr = MatRARt( a_Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
17396434bc1SMark F. Adams   Pold = Pnew;
17496434bc1SMark F. Adams #else
175038e3b61SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
17696434bc1SMark F. Adams #endif
177038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
178038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
179038e3b61SMark F. Adams   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
180038e3b61SMark F. Adams 
181f852f58cSMark F. Adams   /* get number of PEs to make active, reduce */
182f852f58cSMark F. Adams   ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
183f852f58cSMark F. Adams   new_npe = neq/min_eq_proc; /* hardwire min. number of eq/proc */
184f852f58cSMark F. Adams   if( new_npe == 0 || neq < coarse_max ) new_npe = 1;
185f852f58cSMark F. Adams   else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */
186f852f58cSMark F. Adams 
187f852f58cSMark F. Adams   if( avoid_repart && !(new_npe == 1 && *a_nactive_proc != 1) ) {
18822063be5SMark F. Adams     *a_Amat_crs = Cmat; /* output */
18922063be5SMark F. Adams   }
19022063be5SMark F. Adams   else {
19122063be5SMark F. Adams     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
1925ef31b24SMark F. Adams     Mat              adj;
19322063be5SMark F. Adams     const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols;
194b4fbaa2aSMark F. Adams     const PetscInt  stride0=ncrs0*a_ndata_rows;
19592a756f0SMark F. Adams     PetscInt        is_sz,*isnewproc_idx,ii,jj,kk,strideNew,*tidx;
196c9a0b8beSMark F. Adams     /* create sub communicator  */
19771d60426SBarry Smith     MPI_Comm        cm;
198afc97cdcSMark F. Adams     MPI_Group       wg, g2;
199fd3c6afaSMark F. Adams     PetscInt       *counts,inpe;
200fd3c6afaSMark F. Adams     PetscMPIInt    *ranks;
20122063be5SMark F. Adams     IS              isscat;
20222063be5SMark F. Adams     PetscScalar    *array;
20322063be5SMark F. Adams     Vec             src_crd, dest_crd;
20422063be5SMark F. Adams     PetscReal      *data = *a_coarse_data;
20522063be5SMark F. Adams     VecScatter      vecscat;
206b4fbaa2aSMark F. Adams     IS  isnewproc;
207e33ef3b1SMark F. Adams 
20892a756f0SMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr);
209fd3c6afaSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
21092a756f0SMark F. Adams 
211fd3c6afaSMark F. Adams     ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr);
212afc97cdcSMark F. Adams     assert(counts[mype]==ncrs0);
213afc97cdcSMark F. Adams     /* count real active pes */
21422063be5SMark F. Adams     for( nactive = jj = 0 ; jj < npe ; jj++) {
215afc97cdcSMark F. Adams       if( counts[jj] != 0 ) {
21622063be5SMark F. Adams 	ranks[nactive++] = jj;
217afc97cdcSMark F. Adams         }
218afc97cdcSMark F. Adams     }
2193303bcf9Sadams 
2203303bcf9Sadams     if (nactive < new_npe) new_npe = nactive; /* this can happen with empty input procs */
22122063be5SMark F. Adams 
222bed7c2b9SMark 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);
22358471d46SMark F. Adams 
22422063be5SMark F. Adams     *a_nactive_proc = new_npe; /* output */
2252f03bc48SMark F. Adams 
226afc97cdcSMark F. Adams     ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr);
22722063be5SMark F. Adams     ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr);
228afc97cdcSMark F. Adams     ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr);
229b4fbaa2aSMark F. Adams 
230afc97cdcSMark F. Adams     ierr = MPI_Group_free( &wg );                            CHKERRQ(ierr);
231afc97cdcSMark F. Adams     ierr = MPI_Group_free( &g2 );                            CHKERRQ(ierr);
232c9a0b8beSMark F. Adams 
2335ef31b24SMark F. Adams     /* MatPartitioningApply call MatConvert, which is collective */
234b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
235b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
236b4fbaa2aSMark F. Adams #endif
237038e3b61SMark F. Adams     if( a_cbs == 1) {
238038e3b61SMark F. Adams       ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
239eb07cef2SMark F. Adams     }
240eb07cef2SMark F. Adams     else{
241038e3b61SMark F. Adams       /* make a scalar matrix to partition */
242eb07cef2SMark F. Adams       Mat tMat;
24358471d46SMark F. Adams       PetscInt ncols,jj,Ii;
244b4fbaa2aSMark F. Adams       const PetscScalar *vals;
245b4fbaa2aSMark F. Adams       const PetscInt *idx;
24658471d46SMark F. Adams       PetscInt *d_nnz;
2475ee9c036SSatish Balay       static int llev = 0;
248b4fbaa2aSMark F. Adams 
24958471d46SMark F. Adams       ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
25058471d46SMark F. Adams       for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) {
25158471d46SMark F. Adams         ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
25258471d46SMark F. Adams         d_nnz[jj] = ncols/a_cbs;
25358471d46SMark F. Adams         if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
25458471d46SMark F. Adams         ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
25558471d46SMark F. Adams       }
2566876a03eSMark F. Adams 
257eb07cef2SMark F. Adams       ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
258eb07cef2SMark F. Adams                               PETSC_DETERMINE, PETSC_DETERMINE,
25958471d46SMark F. Adams                               0, d_nnz, 0, d_nnz,
260eb07cef2SMark F. Adams                               &tMat );
2616876a03eSMark F. Adams       CHKERRQ(ierr);
26258471d46SMark F. Adams       ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
263eb07cef2SMark F. Adams 
26422063be5SMark F. Adams       for ( ii = Istart0; ii < Iend0; ii++ ) {
26522063be5SMark F. Adams         PetscInt dest_row = ii/a_cbs;
26622063be5SMark F. Adams         ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
267eb07cef2SMark F. Adams         for( jj = 0 ; jj < ncols ; jj++ ){
268038e3b61SMark F. Adams           PetscInt dest_col = idx[jj]/a_cbs;
269eb07cef2SMark F. Adams           PetscScalar v = 1.0;
270eb07cef2SMark F. Adams           ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
271eb07cef2SMark F. Adams         }
27222063be5SMark F. Adams         ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
273eb07cef2SMark F. Adams       }
274eb07cef2SMark F. Adams       ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
275eb07cef2SMark F. Adams       ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
276eb07cef2SMark F. Adams 
277b4fbaa2aSMark F. Adams       if( llev++ == -1 ) {
278b4fbaa2aSMark F. Adams 	PetscViewer viewer; char fname[32];
279b4fbaa2aSMark F. Adams 	sprintf(fname,"part_mat_%d.mat",llev);
280b4fbaa2aSMark F. Adams 	PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
281b4fbaa2aSMark F. Adams 	ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
282b4fbaa2aSMark F. Adams 	ierr = PetscViewerDestroy( &viewer );
283b4fbaa2aSMark F. Adams       }
284b4fbaa2aSMark F. Adams 
285eb07cef2SMark F. Adams       ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
286eb07cef2SMark F. Adams 
287eb07cef2SMark F. Adams       ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
288eb07cef2SMark F. Adams     }
289f150b916SMark F. Adams 
290afc97cdcSMark F. Adams     if( ncrs0 != 0 ){
291b4fbaa2aSMark F. Adams       const PetscInt *is_idx;
292b4fbaa2aSMark F. Adams       MatPartitioning  mpart;
2935ef31b24SMark F. Adams       /* hack to fix global data that pmetis.c uses in 'adj' */
29422063be5SMark F. Adams       for( nactive = jj = 0 ; jj < npe ; jj++ ) {
295afc97cdcSMark F. Adams 	if( counts[jj] != 0 ) {
29622063be5SMark F. Adams 	  adj->rmap->range[nactive++] = adj->rmap->range[jj];
297afc97cdcSMark F. Adams 	}
2985ef31b24SMark F. Adams       }
29922063be5SMark F. Adams       adj->rmap->range[nactive] = adj->rmap->range[npe];
3002f03bc48SMark F. Adams 
30171d60426SBarry Smith       ierr = MatPartitioningCreate( cm, &mpart ); CHKERRQ(ierr);
3025ef31b24SMark F. Adams       ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
30311e60469SMark F. Adams       ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
3045ef31b24SMark F. Adams       ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
30511e60469SMark F. Adams       ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr);
30611e60469SMark F. Adams       ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
30722063be5SMark F. Adams 
3085ef31b24SMark F. Adams       /* collect IS info */
3095ef31b24SMark F. Adams       ierr = ISGetLocalSize( isnewproc, &is_sz );       CHKERRQ(ierr);
310038e3b61SMark F. Adams       ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr);
3115ef31b24SMark F. Adams       ierr = ISGetIndices( isnewproc, &is_idx );        CHKERRQ(ierr);
312b4fbaa2aSMark F. Adams       /* spread partitioning across machine - best way ??? */
313ecd57171SMark F. Adams       NN = 1; /*npe/new_npe;*/
314b4fbaa2aSMark F. Adams       for( kk = jj = 0 ; kk < is_sz ; kk++ ){
315afc97cdcSMark F. Adams         for( ii = 0 ; ii < a_cbs ; ii++, jj++ ) {
31622063be5SMark F. Adams           isnewproc_idx[jj] = is_idx[kk] * NN; /* distribution */
317eb07cef2SMark F. Adams         }
3185ef31b24SMark F. Adams       }
3195ef31b24SMark F. Adams       ierr = ISRestoreIndices( isnewproc, &is_idx );     CHKERRQ(ierr);
3205ef31b24SMark F. Adams       ierr = ISDestroy( &isnewproc );                    CHKERRQ(ierr);
32171d60426SBarry Smith       ierr = MPI_Comm_free( &cm );              CHKERRQ(ierr);
322b4fbaa2aSMark F. Adams 
323b4fbaa2aSMark F. Adams       is_sz *= a_cbs;
3245ef31b24SMark F. Adams     }
3255ef31b24SMark F. Adams     else{
3265ef31b24SMark F. Adams       isnewproc_idx = 0;
3275ef31b24SMark F. Adams       is_sz = 0;
3285ef31b24SMark F. Adams     }
329f150b916SMark F. Adams 
3305ef31b24SMark F. Adams     ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
3315ef31b24SMark F. Adams     ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc );
332afc97cdcSMark F. Adams     if( isnewproc_idx != 0 ) {
3335ef31b24SMark F. Adams       ierr = PetscFree( isnewproc_idx );  CHKERRQ(ierr);
3345ef31b24SMark F. Adams     }
335e33ef3b1SMark F. Adams 
33611e60469SMark F. Adams     /*
33711e60469SMark F. Adams      Create an index set from the isnewproc index set to indicate the mapping TO
33811e60469SMark F. Adams      */
33911e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
34011e60469SMark F. Adams     /*
34111e60469SMark F. Adams      Determine how many elements are assigned to each processor
34211e60469SMark F. Adams      */
343fd3c6afaSMark F. Adams     inpe = npe;
344fd3c6afaSMark F. Adams     ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr);
34511e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
346038e3b61SMark F. Adams     ncrs_new = counts[mype]/a_cbs;
34722063be5SMark F. Adams     strideNew = ncrs_new*a_ndata_rows;
348b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
349b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
350f852f58cSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
351b4fbaa2aSMark F. Adams #endif
35222063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
35311e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
354d3d6bff4SMark F. Adams     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
355470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
35611e60469SMark F. Adams     /*
357d3d6bff4SMark F. Adams      There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
358d3d6bff4SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
35911e60469SMark F. Adams      */
36092a756f0SMark F. Adams     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
36111e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
362038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
363d3d6bff4SMark F. Adams       PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
364d3d6bff4SMark F. Adams       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
36511e60469SMark F. Adams     }
366038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
367d3d6bff4SMark F. Adams     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
36811e60469SMark F. Adams     CHKERRQ(ierr);
36992a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
37011e60469SMark F. Adams     /*
37111e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
37211e60469SMark F. Adams      */
373d3d6bff4SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
374d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
375d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs0 ; ii++) {
376d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
377d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
378676e1743SMark F. Adams           PetscScalar tt = (PetscScalar)data[ix];
379676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
380d3d6bff4SMark F. Adams 	}
381038e3b61SMark F. Adams       }
382eb07cef2SMark F. Adams     }
383eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
384eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
385f852f58cSMark F. Adams #if defined PETSC_USE_LOG
386f852f58cSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET13],0,0,0,0);   CHKERRQ(ierr);
387f852f58cSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
388f852f58cSMark F. Adams #endif
38911e60469SMark F. Adams     /*
39011e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
39111e60469SMark F. Adams       to the correct processor
39211e60469SMark F. Adams     */
39311e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
39411e60469SMark F. Adams     CHKERRQ(ierr);
39511e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
39611e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39711e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39811e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
39911e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
40011e60469SMark F. Adams     /*
40111e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
40211e60469SMark F. Adams     */
403eb07cef2SMark F. Adams     ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
404d3d6bff4SMark F. Adams     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
405eb07cef2SMark F. Adams 
40611e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
407eb07cef2SMark F. Adams     data = *a_coarse_data;
408d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
409d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs_new ; ii++) {
410d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
411d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
412d3d6bff4SMark F. Adams 	  data[ix] = PetscRealPart(array[jx]);
413d3d6bff4SMark F. Adams 	  array[jx] = 1.e300;
414d3d6bff4SMark F. Adams 	}
415038e3b61SMark F. Adams       }
416038e3b61SMark F. Adams     }
41711e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
41811e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
419f852f58cSMark F. Adams #if defined PETSC_USE_LOG
420f852f58cSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET14],0,0,0,0);   CHKERRQ(ierr);
421f852f58cSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
422f852f58cSMark F. Adams #endif
42311e60469SMark F. Adams     /*
42411e60469SMark F. Adams       Invert for MatGetSubMatrix
42511e60469SMark F. Adams     */
426038e3b61SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
42711e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
42811e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
42911e60469SMark F. Adams     /* A_crs output */
430038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
43111e60469SMark F. Adams     CHKERRQ(ierr);
432eb07cef2SMark F. Adams 
433038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
434e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
435038e3b61SMark F. Adams     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
436eb07cef2SMark F. Adams 
43711e60469SMark F. Adams     /* prolongator */
43811e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
43911e60469SMark F. Adams     {
44011e60469SMark F. Adams       IS findices;
44111e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
44296434bc1SMark F. Adams #ifdef USE_R
44396434bc1SMark F. Adams       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
44496434bc1SMark F. Adams #else
44511e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
44696434bc1SMark F. Adams #endif
44711e60469SMark F. Adams       CHKERRQ(ierr);
448d61a3a7dSMark F. Adams 
44911e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
45011e60469SMark F. Adams     }
451d61a3a7dSMark F. Adams 
4523530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
4533530afc2SMark F. Adams     *a_P_inout = Pnew; /* output */
45492a756f0SMark F. Adams 
45511e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
45692a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
45792a756f0SMark F. Adams     ierr = PetscFree( ranks );  CHKERRQ(ierr);
458f852f58cSMark F. Adams #if defined PETSC_USE_LOG
459f852f58cSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET15],0,0,0,0);   CHKERRQ(ierr);
460f852f58cSMark F. Adams #endif
461e33ef3b1SMark F. Adams   }
4625b89ad90SMark F. Adams 
4635b89ad90SMark F. Adams   PetscFunctionReturn(0);
4645b89ad90SMark F. Adams }
4655b89ad90SMark F. Adams 
4665b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4675b89ad90SMark F. Adams /*
4685b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4695b89ad90SMark F. Adams                     by setting data structures and options.
4705b89ad90SMark F. Adams 
4715b89ad90SMark F. Adams    Input Parameter:
4725b89ad90SMark F. Adams .  pc - the preconditioner context
4735b89ad90SMark F. Adams 
4745b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4755b89ad90SMark F. Adams 
4765b89ad90SMark F. Adams    Notes:
4775b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4785b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4795b89ad90SMark F. Adams */
4805b89ad90SMark F. Adams #undef __FUNCT__
4815b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
482eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
4835b89ad90SMark F. Adams {
4845b89ad90SMark F. Adams   PetscErrorCode  ierr;
485eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
4865b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
48758471d46SMark F. Adams   PC_MG_Levels   **mglevels = mg->levels;
488eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
489d3d6bff4SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
490eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
4913530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
492038e3b61SMark F. Adams   PetscBool        isOK;
493587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
494587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
495737a81a9SMark F. Adams   MatInfo          info;
4965ef31b24SMark F. Adams 
4975b89ad90SMark F. Adams   PetscFunctionBegin;
49803a628feSMark F. Adams   pc_gamg->m_count++;
49958471d46SMark F. Adams   if( a_pc->setupcalled > 0 ) {
50003a628feSMark F. Adams     /* just do Galerkin grids */
50158471d46SMark F. Adams     Mat B,dA,dB;
50258471d46SMark F. Adams 
50358471d46SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
50458471d46SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
50558471d46SMark F. Adams 
50658471d46SMark F. Adams     /* currently only handle case where mat and pmat are the same on coarser levels */
50758471d46SMark F. Adams     ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
50858471d46SMark F. Adams     /* (re)set to get dirty flag */
50958471d46SMark F. Adams     ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
51058471d46SMark F. Adams     ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr);
51158471d46SMark F. Adams 
51258471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-2; level>-1; level--) {
51358471d46SMark F. Adams       ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
51403a628feSMark F. Adams       /* the first time through the matrix structure has changed from repartitioning */
51503a628feSMark F. Adams       if( pc_gamg->m_count == 2 ) {
51603a628feSMark F. Adams         ierr = MatDestroy( &B );  CHKERRQ(ierr);
51703a628feSMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
51803a628feSMark F. Adams         mglevels[level]->A = B;
51903a628feSMark F. Adams       }
52003a628feSMark F. Adams       else {
52158471d46SMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
52203a628feSMark F. Adams       }
52358471d46SMark F. Adams       ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
52458471d46SMark F. Adams       dB = B;
52558471d46SMark F. Adams       /* setup KSP/PC */
52658471d46SMark F. Adams       ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
52758471d46SMark F. Adams     }
52858471d46SMark F. Adams 
529f6536408SMark F. Adams #define PRINT_MATS PETSC_FALSE
53058471d46SMark F. Adams     /* plot levels - A */
53158471d46SMark F. Adams     if( PRINT_MATS ) {
53258471d46SMark F. Adams       for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){
53358471d46SMark F. Adams         PetscViewer viewer;
53458471d46SMark F. Adams         char fname[32]; KSP smoother; Mat Tmat, TTm;
53558471d46SMark F. Adams         ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
53658471d46SMark F. Adams         ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr);
53703a628feSMark F. Adams         sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
53858471d46SMark F. Adams         ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
53958471d46SMark F. Adams         ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
54058471d46SMark F. Adams         ierr = MatView( Tmat, viewer ); CHKERRQ(ierr);
54158471d46SMark F. Adams         ierr = PetscViewerDestroy( &viewer );
54258471d46SMark F. Adams       }
54358471d46SMark F. Adams     }
54458471d46SMark F. Adams 
54503a628feSMark F. Adams     a_pc->setupcalled = 2;
54603a628feSMark F. Adams 
54758471d46SMark F. Adams     PetscFunctionReturn(0);
548eb07cef2SMark F. Adams   }
549f6536408SMark F. Adams 
550baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
551baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
5525b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
5535b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
5543530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
555eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
556eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
557eb07cef2SMark F. Adams 
558038e3b61SMark F. Adams   /* get data of not around */
559038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
560038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
561eb07cef2SMark F. Adams   }
562eb07cef2SMark F. Adams   data = pc_gamg->m_data;
563038e3b61SMark F. Adams 
564eb07cef2SMark F. Adams   /* Get A_i and R_i */
565737a81a9SMark F. Adams   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
566bed7c2b9SMark 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",
5672c047e2dSMark F. Adams 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,
5682c047e2dSMark F. Adams 	      (int)(info.nz_used/(PetscReal)N),npe);
5698f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
570389730f3SMark F. Adams         level < (GAMG_MAXLEVELS-1) && (level==0 || M>pc_gamg->m_coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
5710205a208SMark F. Adams         level++ ){
5725b89ad90SMark F. Adams     level1 = level + 1;
573b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
574b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
575b4fbaa2aSMark F. Adams #endif
576b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
577b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
578b4fbaa2aSMark F. Adams #endif
579389730f3SMark F. Adams     ierr = createProlongation(Aarr[level], data, level, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level], pc_gamg );
5805b89ad90SMark F. Adams     CHKERRQ(ierr);
581d3d6bff4SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
582b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
583b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
584b4fbaa2aSMark F. Adams #endif
585baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
586baa4e9faSMark F. Adams     if( isOK ) {
587b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
588b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
589b4fbaa2aSMark F. Adams #endif
590486a8d0bSJed Brown       ierr = PCGAMGPartitionLevel(a_pc, Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs,
591389730f3SMark F. Adams                              &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
5923530afc2SMark F. Adams       CHKERRQ(ierr);
593b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
594b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
595b4fbaa2aSMark F. Adams #endif
5963530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
597737a81a9SMark F. Adams       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
598bed7c2b9SMark 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",
5992c047e2dSMark F. Adams 		  mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols,
6002c047e2dSMark F. Adams 		  (int)(info.nz_used/(PetscReal)N),nactivepe);
601e33ef3b1SMark F. Adams       /* coarse grids with SA can have zero row/cols from singleton aggregates */
60258471d46SMark F. Adams       /* aggregation method should gaurrentee this does not happen! */
603737a81a9SMark F. Adams 
604bed7c2b9SMark F. Adams       if (pc_gamg->m_verbose) {
605785cba28SMark F. Adams         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
606785cba28SMark F. Adams         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
607e33ef3b1SMark F. Adams         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
608785cba28SMark F. Adams         nloceq = Iend-Istart;
609e33ef3b1SMark F. Adams         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
610e33ef3b1SMark F. Adams         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
611e33ef3b1SMark F. Adams         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
612785cba28SMark F. Adams         for(kk=0;kk<nloceq;kk++){
613e33ef3b1SMark F. Adams           if(data_arr[kk]==0.0) {
614785cba28SMark F. Adams             id = kk + Istart;
615e33ef3b1SMark F. Adams             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
616e33ef3b1SMark F. Adams             CHKERRQ(ierr);
617ecd57171SMark F. Adams             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1);
618e33ef3b1SMark F. Adams           }
619e33ef3b1SMark F. Adams         }
620e33ef3b1SMark F. Adams         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
621e33ef3b1SMark F. Adams         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
622e33ef3b1SMark F. Adams         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
623e33ef3b1SMark F. Adams         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
624e33ef3b1SMark F. Adams       }
625486a8d0bSJed Brown     } else {
626be544d3cSMark F. Adams       coarse_data = 0;
627baa4e9faSMark F. Adams       break;
628baa4e9faSMark F. Adams     }
629eb07cef2SMark F. Adams     data = coarse_data;
630b4fbaa2aSMark F. Adams 
631b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
632b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
633b4fbaa2aSMark F. Adams #endif
6345b89ad90SMark F. Adams   }
635be544d3cSMark F. Adams   if( coarse_data ) {
636eb07cef2SMark F. Adams     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
637be544d3cSMark F. Adams   }
638bed7c2b9SMark F. Adams   if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
6395b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
6405b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
6415b89ad90SMark F. Adams   fine_level = level;
642eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
6435b89ad90SMark F. Adams 
644fc4362bfSMark F. Adams   /* set default smoothers */
645587fa25dSMark F. Adams   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
646587fa25dSMark F. Adams         lidx <= fine_level;
647587fa25dSMark F. Adams         lidx++, level--) {
6485745e0f5SMark F. Adams     PetscReal emax, emin;
6495b89ad90SMark F. Adams     KSP smoother; PC subpc;
650f6536408SMark F. Adams     PetscBool isCheb;
651f6536408SMark F. Adams     /* set defaults */
652587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
6535b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
654f6536408SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
655f6536408SMark F. Adams     /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */
656f6536408SMark F. Adams     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
657f6536408SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
658f6536408SMark F. Adams     /* overide defaults with input parameters */
659f6536408SMark F. Adams     ierr = KSPSetFromOptions( smoother ); CHKERRQ(ierr);
660f6536408SMark F. Adams 
661f6536408SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
662f6536408SMark F. Adams     /* do my own cheby */
663f6536408SMark F. Adams     ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr);
664f6536408SMark F. Adams     if( isCheb ) {
665f6536408SMark F. Adams       ierr = PetscTypeCompare( (PetscObject)subpc, PETSC_GAMG_SMOOTHER, &isCheb ); CHKERRQ(ierr);
666f6536408SMark F. Adams       if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
667587fa25dSMark F. Adams       else{ /* eigen estimate 'emax' */
668587fa25dSMark F. Adams         KSP eksp; Mat Lmat = Aarr[level];
669fc4362bfSMark F. Adams         Vec bb, xx; PC pc;
670f6536408SMark F. Adams         const PCType type;
671038e3b61SMark F. Adams 
672f6536408SMark F. Adams         ierr = PCGetType( subpc, &type );   CHKERRQ(ierr);
6735745e0f5SMark F. Adams         ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
6745745e0f5SMark F. Adams         ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
675fc4362bfSMark F. Adams         {
676fc4362bfSMark F. Adams           PetscRandom    rctx;
677fc4362bfSMark F. Adams           ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
678fc4362bfSMark F. Adams           ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
679fc4362bfSMark F. Adams           ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
680fc4362bfSMark F. Adams           ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
6815b89ad90SMark F. Adams         }
682fc4362bfSMark F. Adams         ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
683e33ef3b1SMark F. Adams         ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr);
684038e3b61SMark F. Adams         ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
685fc4362bfSMark F. Adams         CHKERRQ(ierr);
686fc4362bfSMark F. Adams         ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
687fc4362bfSMark F. Adams 
688f6536408SMark F. Adams         ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
689f6536408SMark F. Adams         ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
690f6536408SMark F. Adams 
691f6536408SMark F. Adams         ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
692f6536408SMark F. Adams         ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
693f6536408SMark F. Adams         ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
694f6536408SMark F. Adams         ierr = PCSetType( pc, type ); CHKERRQ(ierr); /* should be same as eigen estimates op. */
695f6536408SMark F. Adams 
696fc4362bfSMark F. Adams         ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
697fc4362bfSMark F. Adams         ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
698fc4362bfSMark F. Adams         ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
699fc4362bfSMark F. Adams         ierr = VecDestroy( &xx );       CHKERRQ(ierr);
700fc4362bfSMark F. Adams         ierr = VecDestroy( &bb );       CHKERRQ(ierr);
701fc4362bfSMark F. Adams         ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
702f6536408SMark F. Adams 
703f6536408SMark F. Adams         if (pc_gamg->m_verbose) {
704f6536408SMark F. Adams           PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s PC setup max eigen=%e min=%e PC=%s\n",
705f6536408SMark F. Adams                       __FUNCT__,emax,emin,PETSC_GAMG_SMOOTHER);
706f6536408SMark F. Adams         }
707fc4362bfSMark F. Adams       }
708038e3b61SMark F. Adams       {
709038e3b61SMark F. Adams         PetscInt N1, N0, tt;
710038e3b61SMark F. Adams         ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
711038e3b61SMark F. Adams         ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
712f6536408SMark F. Adams         /* heuristic - is this crap? */
713f6536408SMark F. Adams         emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
714038e3b61SMark F. Adams         emax *= 1.05;
715038e3b61SMark F. Adams       }
716fc4362bfSMark F. Adams       ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
717f6536408SMark F. Adams     }
7185745e0f5SMark F. Adams   }
7195745e0f5SMark F. Adams   {
720ecd57171SMark F. Adams     /* coarse grid */
721ecd57171SMark F. Adams     KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
7225745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
723eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
72458471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
7255745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
726ecd57171SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
727ecd57171SMark F. Adams     ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
728ecd57171SMark F. Adams     ierr = PCSetUp( subpc ); CHKERRQ(ierr);
729ecd57171SMark F. Adams     ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
730ecd57171SMark F. Adams     assert(ii==1);
731ecd57171SMark F. Adams     ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
732ecd57171SMark F. Adams     ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
733fc4362bfSMark F. Adams   }
734737a81a9SMark F. Adams 
735fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
736eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
7375b89ad90SMark F. Adams   {
7385b89ad90SMark F. Adams     PetscBool galerkin;
739eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
7405b89ad90SMark F. Adams     if(galerkin){
7415b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
7425b89ad90SMark F. Adams     }
7435b89ad90SMark F. Adams   }
7445745e0f5SMark F. Adams 
74558471d46SMark F. Adams   /* plot levels - R/P */
74658471d46SMark F. Adams   if( PRINT_MATS ) {
74758471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-1;level>0;level--){
74858471d46SMark F. Adams       PetscViewer viewer;
74958471d46SMark F. Adams       char fname[32];
75003a628feSMark F. Adams       sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
75111e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
7525b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
753038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
7545b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
75503a628feSMark F. Adams       sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
756e33ef3b1SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
757e33ef3b1SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
758e33ef3b1SMark F. Adams       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
759e33ef3b1SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
7605b89ad90SMark F. Adams     }
76158471d46SMark F. Adams   }
76258471d46SMark F. Adams 
76358471d46SMark F. Adams   /* set interpolation between the levels, clean up */
76458471d46SMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
76558471d46SMark F. Adams        lidx<fine_level;
76658471d46SMark F. Adams        lidx++, level--){
76758471d46SMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
768587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
769587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
7705b89ad90SMark F. Adams   }
7715745e0f5SMark F. Adams 
7725b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
773eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
774eb07cef2SMark F. Adams   ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
77503a628feSMark F. Adams   a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
77658471d46SMark F. Adams 
77758471d46SMark F. Adams   {
77858471d46SMark F. Adams     KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
77958471d46SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
78058471d46SMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
78158471d46SMark F. Adams     ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
78258471d46SMark F. Adams   }
783e33ef3b1SMark F. Adams 
7845b89ad90SMark F. Adams   PetscFunctionReturn(0);
7855b89ad90SMark F. Adams }
7865b89ad90SMark F. Adams 
787eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7885b89ad90SMark F. Adams /*
7895b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7905b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7915b89ad90SMark F. Adams 
7925b89ad90SMark F. Adams    Input Parameter:
7935b89ad90SMark F. Adams .  pc - the preconditioner context
7945b89ad90SMark F. Adams 
7955b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7965b89ad90SMark F. Adams */
7975b89ad90SMark F. Adams #undef __FUNCT__
7985b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
7995b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8005b89ad90SMark F. Adams {
8015b89ad90SMark F. Adams   PetscErrorCode  ierr;
8025b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
8035b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
8045b89ad90SMark F. Adams 
8055b89ad90SMark F. Adams   PetscFunctionBegin;
8065b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8075b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8085b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8095b89ad90SMark F. Adams   PetscFunctionReturn(0);
8105b89ad90SMark F. Adams }
8115b89ad90SMark F. Adams 
812676e1743SMark F. Adams 
813676e1743SMark F. Adams #undef __FUNCT__
814676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
815676e1743SMark F. Adams /*@
816676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
817676e1743SMark F. Adams    processor reduction.
818676e1743SMark F. Adams 
819676e1743SMark F. Adams    Not Collective on PC
820676e1743SMark F. Adams 
821676e1743SMark F. Adams    Input Parameters:
822676e1743SMark F. Adams .  pc - the preconditioner context
823676e1743SMark F. Adams 
824676e1743SMark F. Adams 
825676e1743SMark F. Adams    Options Database Key:
826676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
827676e1743SMark F. Adams 
828676e1743SMark F. Adams    Level: intermediate
829676e1743SMark F. Adams 
830676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
831676e1743SMark F. Adams 
832676e1743SMark F. Adams .seealso: ()
833676e1743SMark F. Adams @*/
834676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
835676e1743SMark F. Adams {
836676e1743SMark F. Adams   PetscErrorCode ierr;
837676e1743SMark F. Adams 
838676e1743SMark F. Adams   PetscFunctionBegin;
839676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
840676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
841676e1743SMark F. Adams   PetscFunctionReturn(0);
842676e1743SMark F. Adams }
843676e1743SMark F. Adams 
844676e1743SMark F. Adams EXTERN_C_BEGIN
845676e1743SMark F. Adams #undef __FUNCT__
846676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
847676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
848676e1743SMark F. Adams {
849c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
850c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
851676e1743SMark F. Adams 
852676e1743SMark F. Adams   PetscFunctionBegin;
853c20e4228SMark F. Adams   if(n>0) pc_gamg->m_min_eq_proc = n;
854676e1743SMark F. Adams   PetscFunctionReturn(0);
855676e1743SMark F. Adams }
856676e1743SMark F. Adams EXTERN_C_END
857676e1743SMark F. Adams 
858676e1743SMark F. Adams #undef __FUNCT__
859389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
860389730f3SMark F. Adams /*@
861389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
862389730f3SMark F. Adams 
863389730f3SMark F. Adams  Collective on PC
864389730f3SMark F. Adams 
865389730f3SMark F. Adams    Input Parameters:
866389730f3SMark F. Adams .  pc - the preconditioner context
867389730f3SMark F. Adams 
868389730f3SMark F. Adams 
869389730f3SMark F. Adams    Options Database Key:
870389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
871389730f3SMark F. Adams 
872389730f3SMark F. Adams    Level: intermediate
873389730f3SMark F. Adams 
874389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
875389730f3SMark F. Adams 
876389730f3SMark F. Adams .seealso: ()
877389730f3SMark F. Adams  @*/
878389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
879389730f3SMark F. Adams {
880389730f3SMark F. Adams   PetscErrorCode ierr;
881389730f3SMark F. Adams 
882389730f3SMark F. Adams   PetscFunctionBegin;
883389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
884389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
885389730f3SMark F. Adams   PetscFunctionReturn(0);
886389730f3SMark F. Adams }
887389730f3SMark F. Adams 
888389730f3SMark F. Adams EXTERN_C_BEGIN
889389730f3SMark F. Adams #undef __FUNCT__
890389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
891389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
892389730f3SMark F. Adams {
893389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
894389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
895389730f3SMark F. Adams 
896389730f3SMark F. Adams   PetscFunctionBegin;
897389730f3SMark F. Adams   if(n>0) pc_gamg->m_coarse_eq_limit = n;
898389730f3SMark F. Adams   PetscFunctionReturn(0);
899389730f3SMark F. Adams }
900389730f3SMark F. Adams EXTERN_C_END
901389730f3SMark F. Adams 
902389730f3SMark F. Adams #undef __FUNCT__
903676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning"
904676e1743SMark F. Adams /*@
905676e1743SMark F. Adams    PCGAMGAvoidRepartitioning - Do not repartition the coarse grids
906676e1743SMark F. Adams 
907676e1743SMark F. Adams    Collective on PC
908676e1743SMark F. Adams 
909676e1743SMark F. Adams    Input Parameters:
910676e1743SMark F. Adams .  pc - the preconditioner context
911676e1743SMark F. Adams 
912676e1743SMark F. Adams 
913676e1743SMark F. Adams    Options Database Key:
914676e1743SMark F. Adams .  -pc_gamg_avoid_repartitioning
915676e1743SMark F. Adams 
916676e1743SMark F. Adams    Level: intermediate
917676e1743SMark F. Adams 
918676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
919676e1743SMark F. Adams 
920676e1743SMark F. Adams .seealso: ()
921676e1743SMark F. Adams @*/
922676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning(PC pc, PetscBool n)
923676e1743SMark F. Adams {
924676e1743SMark F. Adams   PetscErrorCode ierr;
925676e1743SMark F. Adams 
926676e1743SMark F. Adams   PetscFunctionBegin;
927676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
928676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGAvoidRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
929676e1743SMark F. Adams   PetscFunctionReturn(0);
930676e1743SMark F. Adams }
931676e1743SMark F. Adams 
932676e1743SMark F. Adams EXTERN_C_BEGIN
933676e1743SMark F. Adams #undef __FUNCT__
934676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning_GAMG"
935676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning_GAMG(PC pc, PetscBool n)
936676e1743SMark F. Adams {
937c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
938c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
939676e1743SMark F. Adams 
940676e1743SMark F. Adams   PetscFunctionBegin;
941c20e4228SMark F. Adams   pc_gamg->m_avoid_repart = n;
942676e1743SMark F. Adams   PetscFunctionReturn(0);
943676e1743SMark F. Adams }
944676e1743SMark F. Adams EXTERN_C_END
945676e1743SMark F. Adams 
946676e1743SMark F. Adams #undef __FUNCT__
9473542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
9483542efc5SMark F. Adams /*@
9493542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
9503542efc5SMark F. Adams 
9513542efc5SMark F. Adams    Not collective on PC
9523542efc5SMark F. Adams 
9533542efc5SMark F. Adams    Input Parameters:
9543542efc5SMark F. Adams .  pc - the preconditioner context
9553542efc5SMark F. Adams 
9563542efc5SMark F. Adams 
9573542efc5SMark F. Adams    Options Database Key:
9583542efc5SMark F. Adams .  -pc_gamg_threshold
9593542efc5SMark F. Adams 
9603542efc5SMark F. Adams    Level: intermediate
9613542efc5SMark F. Adams 
9623542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
9633542efc5SMark F. Adams 
9643542efc5SMark F. Adams .seealso: ()
9653542efc5SMark F. Adams @*/
9663542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
9673542efc5SMark F. Adams {
9683542efc5SMark F. Adams   PetscErrorCode ierr;
9693542efc5SMark F. Adams 
9703542efc5SMark F. Adams   PetscFunctionBegin;
9713542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9723542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
9733542efc5SMark F. Adams   PetscFunctionReturn(0);
9743542efc5SMark F. Adams }
9753542efc5SMark F. Adams 
9763542efc5SMark F. Adams EXTERN_C_BEGIN
9773542efc5SMark F. Adams #undef __FUNCT__
9783542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
9793542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
9803542efc5SMark F. Adams {
981c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
982c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
9833542efc5SMark F. Adams 
9843542efc5SMark F. Adams   PetscFunctionBegin;
985c20e4228SMark F. Adams   pc_gamg->m_threshold = n;
9863542efc5SMark F. Adams   PetscFunctionReturn(0);
9873542efc5SMark F. Adams }
9883542efc5SMark F. Adams EXTERN_C_END
9893542efc5SMark F. Adams 
9903542efc5SMark F. Adams #undef __FUNCT__
991676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType"
992676e1743SMark F. Adams /*@
993676e1743SMark F. Adams    PCGAMGSetSolverType - Set solution method.
994676e1743SMark F. Adams 
995676e1743SMark F. Adams    Collective on PC
996676e1743SMark F. Adams 
997676e1743SMark F. Adams    Input Parameters:
998676e1743SMark F. Adams .  pc - the preconditioner context
999676e1743SMark F. Adams 
1000676e1743SMark F. Adams 
1001676e1743SMark F. Adams    Options Database Key:
10023542efc5SMark F. Adams .  -pc_gamg_type
1003676e1743SMark F. Adams 
1004676e1743SMark F. Adams    Level: intermediate
1005676e1743SMark F. Adams 
1006676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1007676e1743SMark F. Adams 
1008676e1743SMark F. Adams .seealso: ()
1009676e1743SMark F. Adams @*/
1010676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz )
1011676e1743SMark F. Adams {
1012676e1743SMark F. Adams   PetscErrorCode ierr;
1013676e1743SMark F. Adams 
1014676e1743SMark F. Adams   PetscFunctionBegin;
1015676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1016676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz));
1017676e1743SMark F. Adams   CHKERRQ(ierr);
1018676e1743SMark F. Adams   PetscFunctionReturn(0);
1019676e1743SMark F. Adams }
1020676e1743SMark F. Adams 
1021676e1743SMark F. Adams EXTERN_C_BEGIN
1022676e1743SMark F. Adams #undef __FUNCT__
1023676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG"
1024676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz )
1025676e1743SMark F. Adams {
1026676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1027676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1028676e1743SMark F. Adams 
1029676e1743SMark F. Adams   PetscFunctionBegin;
1030676e1743SMark F. Adams   if(sz < 64) strcpy(pc_gamg->m_type,str);
1031676e1743SMark F. Adams   PetscFunctionReturn(0);
1032676e1743SMark F. Adams }
1033676e1743SMark F. Adams EXTERN_C_END
1034676e1743SMark F. Adams 
10355b89ad90SMark F. Adams #undef __FUNCT__
10365b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
10375b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
10385b89ad90SMark F. Adams {
1039676e1743SMark F. Adams   PetscErrorCode  ierr;
1040676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1041676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1042676e1743SMark F. Adams   PetscBool flag;
10435b89ad90SMark F. Adams 
10445b89ad90SMark F. Adams   PetscFunctionBegin;
1045676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1046676e1743SMark F. Adams   {
1047676e1743SMark F. Adams     ierr = PetscOptionsString("-pc_gamg_type",
1048470e26f8SMark F. Adams                               "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid (default)",
1049676e1743SMark F. Adams                               "PCGAMGSetSolverType",
1050676e1743SMark F. Adams                               pc_gamg->m_type,
1051676e1743SMark F. Adams                               pc_gamg->m_type,
1052676e1743SMark F. Adams                               64,
1053676e1743SMark F. Adams                               &flag );
1054676e1743SMark F. Adams     CHKERRQ(ierr);
1055*d8c859f0SMark F. Adams     if( flag && pc_gamg->m_data != 0 ) {
1056*d8c859f0SMark F. Adams       if( (strcmp(pc_gamg->m_type,"sa")==0 && pc_gamg->m_method != 2) ||
1057*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"pa")==0 && pc_gamg->m_method != 1) ||
1058*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"geo")==0 && pc_gamg->m_method != 0) ) {
1059*d8c859f0SMark F. Adams         SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "PCSetFromOptions called of PCSetCoordinates (with new method, after data was created)");
1060*d8c859f0SMark F. Adams       }
1061*d8c859f0SMark F. Adams     }
1062bed7c2b9SMark F. Adams 
1063bed7c2b9SMark F. Adams     /* -pc_gamg_verbose */
1064bed7c2b9SMark 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);
10652a44bfbaSMark F. Adams 
1066f6536408SMark F. Adams     pc_gamg->m_method = 1; /* default to plane aggregation */
1067f6536408SMark F. Adams     if (flag ) {
1068f6536408SMark F. Adams       if( strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2;
1069f6536408SMark F. Adams       else if( strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1;
1070f6536408SMark F. Adams       else if( strcmp(pc_gamg->m_type,"geo") == 0) pc_gamg->m_method = 0;
1071f6536408SMark F. Adams       else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Invalid gamg type: %s",pc_gamg->m_type);
1072f6536408SMark F. Adams     }
1073c20e4228SMark F. Adams     /* -pc_gamg_avoid_repartitioning */
1074c20e4228SMark F. Adams     pc_gamg->m_avoid_repart = PETSC_FALSE;
1075676e1743SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_avoid_repartitioning",
1076676e1743SMark F. Adams                             "Do not repartion coarse grids (false)",
1077676e1743SMark F. Adams                             "PCGAMGAvoidRepartitioning",
1078c20e4228SMark F. Adams                             pc_gamg->m_avoid_repart,
1079c20e4228SMark F. Adams                             &pc_gamg->m_avoid_repart,
1080676e1743SMark F. Adams                             &flag);
1081676e1743SMark F. Adams     CHKERRQ(ierr);
1082676e1743SMark F. Adams 
1083c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1084c20e4228SMark F. Adams     pc_gamg->m_min_eq_proc = 600;
1085676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1086676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1087676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
1088c20e4228SMark F. Adams                            pc_gamg->m_min_eq_proc,
1089c20e4228SMark F. Adams                            &pc_gamg->m_min_eq_proc,
1090676e1743SMark F. Adams                            &flag );
1091676e1743SMark F. Adams     CHKERRQ(ierr);
10923542efc5SMark F. Adams 
1093389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1094389730f3SMark F. Adams     pc_gamg->m_coarse_eq_limit = 1500;
1095389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1096389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1097389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
1098389730f3SMark F. Adams                            pc_gamg->m_coarse_eq_limit,
1099389730f3SMark F. Adams                            &pc_gamg->m_coarse_eq_limit,
1100389730f3SMark F. Adams                            &flag );
1101389730f3SMark F. Adams     CHKERRQ(ierr);
1102389730f3SMark F. Adams 
1103c20e4228SMark F. Adams     /* -pc_gamg_threshold */
1104c20e4228SMark F. Adams     pc_gamg->m_threshold = 0.05;
11053542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
11063542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
11073542efc5SMark F. Adams                             "PCGAMGSetThreshold",
1108c20e4228SMark F. Adams                             pc_gamg->m_threshold,
1109c20e4228SMark F. Adams                             &pc_gamg->m_threshold,
11103542efc5SMark F. Adams                             &flag );
11113542efc5SMark F. Adams     CHKERRQ(ierr);
1112676e1743SMark F. Adams   }
1113676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1114676e1743SMark F. Adams 
11155b89ad90SMark F. Adams   PetscFunctionReturn(0);
11165b89ad90SMark F. Adams }
11175b89ad90SMark F. Adams 
11185b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
11195b89ad90SMark F. Adams /*MC
1120dc76db92SMark F. Adams      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning. This preconditioner currently has two
1121dc76db92SMark F. Adams            AMG methods: 1) an unstructured geometric method, which requires that you provide coordinates for each
1122dc76db92SMark F. Adams            vertex, and 2) smoothed aggregation.  Smoothed aggregation (SA) can work without coordinates but it
1123dc76db92SMark F. Adams            will generate some common non-trivial null spaces if coordinates are provided.  The input fine grid matrix
1124dc76db92SMark F. Adams            must have the block size set for 'system' problems (with multiple dof per vertex/cell) to work properly.
1125dc76db92SMark F. Adams            SA will generate rotational rigid body mode null space vectors, in addition to the trivial translational
1126dc76db92SMark F. Adams            modes, when coordinates are provide in 2D and 3D.
11275b89ad90SMark F. Adams 
1128280d9858SJed Brown    Options Database Keys:
11295b89ad90SMark F. Adams    Multigrid options(inherited)
1130280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1131280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1132280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1133280d9858SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
11345b89ad90SMark F. Adams 
11355b89ad90SMark F. Adams   Level: intermediate
1136280d9858SJed Brown 
11375b89ad90SMark F. Adams   Concepts: multigrid
11385b89ad90SMark F. Adams 
11395b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1140280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
11415b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
11425b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
11435b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11445b89ad90SMark F. Adams M*/
11455b89ad90SMark F. Adams 
11465b89ad90SMark F. Adams EXTERN_C_BEGIN
11475b89ad90SMark F. Adams #undef __FUNCT__
11485b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
11495b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
11505b89ad90SMark F. Adams {
11515b89ad90SMark F. Adams   PetscErrorCode  ierr;
11525b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
11535b89ad90SMark F. Adams   PC_MG           *mg;
11545ef31b24SMark F. Adams   PetscClassId     cookie;
11555ee9c036SSatish Balay #if defined PETSC_USE_LOG
11565ee9c036SSatish Balay   static int count = 0;
11575ee9c036SSatish Balay #endif
11585b89ad90SMark F. Adams 
11595b89ad90SMark F. Adams   PetscFunctionBegin;
11605b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
11615b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
11625b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
11635b89ad90SMark F. Adams 
11645b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
11655b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
116603a628feSMark F. Adams   pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0;
11675b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
11685b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
11695b89ad90SMark F. Adams 
11705b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
11715b89ad90SMark F. Adams 
11725b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
11735b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
11745b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
11755b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
11765b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
11775b89ad90SMark F. Adams 
11785b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
11795b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
11805b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
1181c97e1df0SMark F. Adams 					    PCSetCoordinates_GAMG);
1182c97e1df0SMark F. Adams   CHKERRQ(ierr);
1183c97e1df0SMark F. Adams 
1184676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1185676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1186676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1187676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1188676e1743SMark F. Adams   CHKERRQ(ierr);
1189676e1743SMark F. Adams 
1190676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1191389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_C",
1192389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_GAMG",
1193389730f3SMark F. Adams 					    PCGAMGSetCoarseEqLim_GAMG);
1194389730f3SMark F. Adams   CHKERRQ(ierr);
1195389730f3SMark F. Adams 
1196389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1197676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_C",
1198676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_GAMG",
1199676e1743SMark F. Adams 					    PCGAMGAvoidRepartitioning_GAMG);
1200676e1743SMark F. Adams   CHKERRQ(ierr);
1201676e1743SMark F. Adams 
1202676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12033542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
12043542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
12053542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
12063542efc5SMark F. Adams   CHKERRQ(ierr);
12073542efc5SMark F. Adams 
12083542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1209676e1743SMark F. Adams 					    "PCGAMGSetSolverType_C",
1210676e1743SMark F. Adams 					    "PCGAMGSetSolverType_GAMG",
1211676e1743SMark F. Adams 					    PCGAMGSetSolverType_GAMG);
1212676e1743SMark F. Adams   CHKERRQ(ierr);
1213c97e1df0SMark F. Adams 
1214b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
1215785cba28SMark F. Adams   if( count++ == 0 ) {
12165ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
1217b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
12182a44bfbaSMark F. Adams     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
12192a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]);
12202a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]);
12212a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]);
1222b4fbaa2aSMark F. Adams     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1223b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1224b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1225b4fbaa2aSMark F. Adams     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1226b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_events[SET7]);
1227b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]); */
1228b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1229b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1230b4fbaa2aSMark F. Adams     PetscLogEventRegister("  PL repartition", cookie, &gamg_setup_events[SET12]);
1231f852f58cSMark F. Adams 
1232f852f58cSMark F. Adams     PetscLogEventRegister("  PL 1", cookie, &gamg_setup_events[SET13]);
1233f852f58cSMark F. Adams     PetscLogEventRegister("  PL 2", cookie, &gamg_setup_events[SET14]);
1234f852f58cSMark F. Adams     PetscLogEventRegister("  PL 3", cookie, &gamg_setup_events[SET15]);
1235f852f58cSMark F. Adams 
1236b4fbaa2aSMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1237b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1238b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
12395ef31b24SMark F. Adams 
1240b4fbaa2aSMark F. Adams     /* create timer stages */
1241b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1242b4fbaa2aSMark F. Adams     {
1243b4fbaa2aSMark F. Adams       char str[32];
1244b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1245b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1246b4fbaa2aSMark F. Adams       PetscInt lidx;
1247b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1248b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1249b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1250b4fbaa2aSMark F. Adams       }
1251b4fbaa2aSMark F. Adams     }
1252b4fbaa2aSMark F. Adams #endif
1253b4fbaa2aSMark F. Adams   }
1254b4fbaa2aSMark F. Adams #endif
12555b89ad90SMark F. Adams   PetscFunctionReturn(0);
12565b89ad90SMark F. Adams }
12575b89ad90SMark F. Adams EXTERN_C_END
1258