xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 4ef23d270073dc44b8ba5ee8ea19d38f441779a5)
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);
350b4fbaa2aSMark F. Adams #endif
35122063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
35211e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
353d3d6bff4SMark F. Adams     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
354470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
35511e60469SMark F. Adams     /*
356d3d6bff4SMark F. Adams      There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
357d3d6bff4SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
35811e60469SMark F. Adams      */
35992a756f0SMark F. Adams     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
36011e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
361038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
362d3d6bff4SMark F. Adams       PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
363d3d6bff4SMark F. Adams       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
36411e60469SMark F. Adams     }
365038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
366d3d6bff4SMark F. Adams     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
36711e60469SMark F. Adams     CHKERRQ(ierr);
36892a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
36911e60469SMark F. Adams     /*
37011e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
37111e60469SMark F. Adams      */
372d3d6bff4SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
373d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
374d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs0 ; ii++) {
375d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
376d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
377676e1743SMark F. Adams           PetscScalar tt = (PetscScalar)data[ix];
378676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
379d3d6bff4SMark F. Adams 	}
380038e3b61SMark F. Adams       }
381eb07cef2SMark F. Adams     }
382eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
383eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
38411e60469SMark F. Adams     /*
38511e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
38611e60469SMark F. Adams       to the correct processor
38711e60469SMark F. Adams     */
38811e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
38911e60469SMark F. Adams     CHKERRQ(ierr);
39011e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
39111e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39211e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39311e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
39411e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
39511e60469SMark F. Adams     /*
39611e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
39711e60469SMark F. Adams     */
398eb07cef2SMark F. Adams     ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
399d3d6bff4SMark F. Adams     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
400eb07cef2SMark F. Adams 
40111e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
402eb07cef2SMark F. Adams     data = *a_coarse_data;
403d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
404d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs_new ; ii++) {
405d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
406d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
407d3d6bff4SMark F. Adams 	  data[ix] = PetscRealPart(array[jx]);
408d3d6bff4SMark F. Adams 	  array[jx] = 1.e300;
409d3d6bff4SMark F. Adams 	}
410038e3b61SMark F. Adams       }
411038e3b61SMark F. Adams     }
41211e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
41311e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
41411e60469SMark F. Adams     /*
41511e60469SMark F. Adams       Invert for MatGetSubMatrix
41611e60469SMark F. Adams     */
417038e3b61SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
41811e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
41911e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
42011e60469SMark F. Adams     /* A_crs output */
421038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
42211e60469SMark F. Adams     CHKERRQ(ierr);
423eb07cef2SMark F. Adams 
424038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
425e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
426038e3b61SMark F. Adams     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
427eb07cef2SMark F. Adams 
42811e60469SMark F. Adams     /* prolongator */
42911e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
43011e60469SMark F. Adams     {
43111e60469SMark F. Adams       IS findices;
43211e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
43396434bc1SMark F. Adams #ifdef USE_R
43496434bc1SMark F. Adams       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
43596434bc1SMark F. Adams #else
43611e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
43796434bc1SMark F. Adams #endif
43811e60469SMark F. Adams       CHKERRQ(ierr);
439d61a3a7dSMark F. Adams 
44011e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
44111e60469SMark F. Adams     }
442d61a3a7dSMark F. Adams 
4433530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
4443530afc2SMark F. Adams     *a_P_inout = Pnew; /* output */
44592a756f0SMark F. Adams 
44611e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
44792a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
44892a756f0SMark F. Adams     ierr = PetscFree( ranks );  CHKERRQ(ierr);
449e33ef3b1SMark F. Adams   }
4505b89ad90SMark F. Adams 
4515b89ad90SMark F. Adams   PetscFunctionReturn(0);
4525b89ad90SMark F. Adams }
4535b89ad90SMark F. Adams 
4545b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4555b89ad90SMark F. Adams /*
4565b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4575b89ad90SMark F. Adams                     by setting data structures and options.
4585b89ad90SMark F. Adams 
4595b89ad90SMark F. Adams    Input Parameter:
4605b89ad90SMark F. Adams .  pc - the preconditioner context
4615b89ad90SMark F. Adams 
4625b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4635b89ad90SMark F. Adams 
4645b89ad90SMark F. Adams    Notes:
4655b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4665b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4675b89ad90SMark F. Adams */
4685b89ad90SMark F. Adams #undef __FUNCT__
4695b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
470eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
4715b89ad90SMark F. Adams {
4725b89ad90SMark F. Adams   PetscErrorCode  ierr;
473eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
4745b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
47558471d46SMark F. Adams   PC_MG_Levels   **mglevels = mg->levels;
476eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
477d3d6bff4SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
478eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
4793530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
480038e3b61SMark F. Adams   PetscBool        isOK;
481587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
482587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
483737a81a9SMark F. Adams   MatInfo          info;
4845ef31b24SMark F. Adams 
4855b89ad90SMark F. Adams   PetscFunctionBegin;
48603a628feSMark F. Adams   pc_gamg->m_count++;
487fca9ed99SMark F. Adams 
48858471d46SMark F. Adams   if( a_pc->setupcalled > 0 ) {
48903a628feSMark F. Adams     /* just do Galerkin grids */
49058471d46SMark F. Adams     Mat B,dA,dB;
49158471d46SMark F. Adams 
49258471d46SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
49358471d46SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
49458471d46SMark F. Adams 
49558471d46SMark F. Adams     /* currently only handle case where mat and pmat are the same on coarser levels */
49658471d46SMark F. Adams     ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
49758471d46SMark F. Adams     /* (re)set to get dirty flag */
49858471d46SMark F. Adams     ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
49958471d46SMark F. Adams     ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr);
50058471d46SMark F. Adams 
50158471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-2; level>-1; level--) {
50258471d46SMark F. Adams       ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
50303a628feSMark F. Adams       /* the first time through the matrix structure has changed from repartitioning */
50403a628feSMark F. Adams       if( pc_gamg->m_count == 2 ) {
50503a628feSMark F. Adams         ierr = MatDestroy( &B );  CHKERRQ(ierr);
50603a628feSMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
50703a628feSMark F. Adams         mglevels[level]->A = B;
50803a628feSMark F. Adams       }
50903a628feSMark F. Adams       else {
51058471d46SMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
51103a628feSMark F. Adams       }
51258471d46SMark F. Adams       ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
51358471d46SMark F. Adams       dB = B;
51458471d46SMark F. Adams       /* setup KSP/PC */
51558471d46SMark F. Adams       ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
51658471d46SMark F. Adams     }
51758471d46SMark F. Adams 
518f6536408SMark F. Adams #define PRINT_MATS PETSC_FALSE
51958471d46SMark F. Adams     /* plot levels - A */
52058471d46SMark F. Adams     if( PRINT_MATS ) {
52158471d46SMark F. Adams       for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){
52258471d46SMark F. Adams         PetscViewer viewer;
52358471d46SMark F. Adams         char fname[32]; KSP smoother; Mat Tmat, TTm;
52458471d46SMark F. Adams         ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
52558471d46SMark F. Adams         ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr);
52603a628feSMark F. Adams         sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
52758471d46SMark F. Adams         ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
52858471d46SMark F. Adams         ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
52958471d46SMark F. Adams         ierr = MatView( Tmat, viewer ); CHKERRQ(ierr);
53058471d46SMark F. Adams         ierr = PetscViewerDestroy( &viewer );
53158471d46SMark F. Adams       }
53258471d46SMark F. Adams     }
53358471d46SMark F. Adams 
53403a628feSMark F. Adams     a_pc->setupcalled = 2;
53503a628feSMark F. Adams 
53658471d46SMark F. Adams     PetscFunctionReturn(0);
537eb07cef2SMark F. Adams   }
538f6536408SMark F. Adams 
539baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
540baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
5415b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
5425b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
5433530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
544eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
545eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
546eb07cef2SMark F. Adams 
547038e3b61SMark F. Adams   /* get data of not around */
548038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
549038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
550eb07cef2SMark F. Adams   }
551eb07cef2SMark F. Adams   data = pc_gamg->m_data;
552038e3b61SMark F. Adams 
553eb07cef2SMark F. Adams   /* Get A_i and R_i */
554737a81a9SMark F. Adams   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
555bed7c2b9SMark 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",
5562c047e2dSMark F. Adams 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,
5572c047e2dSMark F. Adams 	      (int)(info.nz_used/(PetscReal)N),npe);
5588f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
5594ef23d27SMark F. Adams         level < (pc_gamg->m_Nlevels-1) && (level==0 || M>pc_gamg->m_coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
5600205a208SMark F. Adams         level++ ){
5615b89ad90SMark F. Adams     level1 = level + 1;
562b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
563b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
564b4fbaa2aSMark F. Adams #endif
565b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
566b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
567b4fbaa2aSMark F. Adams #endif
568389730f3SMark F. Adams     ierr = createProlongation(Aarr[level], data, level, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level], pc_gamg );
5695b89ad90SMark F. Adams     CHKERRQ(ierr);
570d3d6bff4SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
571b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
572b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
573b4fbaa2aSMark F. Adams #endif
574baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
575baa4e9faSMark F. Adams     if( isOK ) {
576b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
577b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
578b4fbaa2aSMark F. Adams #endif
579486a8d0bSJed Brown       ierr = PCGAMGPartitionLevel(a_pc, Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs,
580389730f3SMark F. Adams                              &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
5813530afc2SMark F. Adams       CHKERRQ(ierr);
582b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
583b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
584b4fbaa2aSMark F. Adams #endif
5853530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
586737a81a9SMark F. Adams       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
587bed7c2b9SMark 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",
5882c047e2dSMark F. Adams 		  mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols,
5892c047e2dSMark F. Adams 		  (int)(info.nz_used/(PetscReal)N),nactivepe);
590e33ef3b1SMark F. Adams       /* coarse grids with SA can have zero row/cols from singleton aggregates */
59158471d46SMark F. Adams       /* aggregation method should gaurrentee this does not happen! */
592737a81a9SMark F. Adams 
593bed7c2b9SMark F. Adams       if (pc_gamg->m_verbose) {
594785cba28SMark F. Adams         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
595785cba28SMark F. Adams         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
596e33ef3b1SMark F. Adams         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
597785cba28SMark F. Adams         nloceq = Iend-Istart;
598e33ef3b1SMark F. Adams         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
599e33ef3b1SMark F. Adams         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
600e33ef3b1SMark F. Adams         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
601785cba28SMark F. Adams         for(kk=0;kk<nloceq;kk++){
602e33ef3b1SMark F. Adams           if(data_arr[kk]==0.0) {
603785cba28SMark F. Adams             id = kk + Istart;
604e33ef3b1SMark F. Adams             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
605e33ef3b1SMark F. Adams             CHKERRQ(ierr);
606ecd57171SMark F. Adams             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1);
607e33ef3b1SMark F. Adams           }
608e33ef3b1SMark F. Adams         }
609e33ef3b1SMark F. Adams         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
610e33ef3b1SMark F. Adams         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
611e33ef3b1SMark F. Adams         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
612e33ef3b1SMark F. Adams         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
613e33ef3b1SMark F. Adams       }
614486a8d0bSJed Brown     } else {
615be544d3cSMark F. Adams       coarse_data = 0;
616baa4e9faSMark F. Adams       break;
617baa4e9faSMark F. Adams     }
618eb07cef2SMark F. Adams     data = coarse_data;
619b4fbaa2aSMark F. Adams 
620b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
621b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
622b4fbaa2aSMark F. Adams #endif
6235b89ad90SMark F. Adams   }
624be544d3cSMark F. Adams   if( coarse_data ) {
625eb07cef2SMark F. Adams     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
626be544d3cSMark F. Adams   }
627bed7c2b9SMark F. Adams   if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
6285b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
6295b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
6305b89ad90SMark F. Adams   fine_level = level;
631eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
6325b89ad90SMark F. Adams 
633fc4362bfSMark F. Adams   /* set default smoothers */
634587fa25dSMark F. Adams   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
635587fa25dSMark F. Adams         lidx <= fine_level;
636587fa25dSMark F. Adams         lidx++, level--) {
6375745e0f5SMark F. Adams     PetscReal emax, emin;
6385b89ad90SMark F. Adams     KSP smoother; PC subpc;
639f6536408SMark F. Adams     PetscBool isCheb;
640f6536408SMark F. Adams     /* set defaults */
641587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
6425b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
643f6536408SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
644f6536408SMark F. Adams     /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */
645f6536408SMark F. Adams     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
646f6536408SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
647f6536408SMark F. Adams     /* overide defaults with input parameters */
648f6536408SMark F. Adams     ierr = KSPSetFromOptions( smoother ); CHKERRQ(ierr);
649f6536408SMark F. Adams 
650f6536408SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
651f6536408SMark F. Adams     /* do my own cheby */
652f6536408SMark F. Adams     ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr);
653f6536408SMark F. Adams     if( isCheb ) {
654f6536408SMark F. Adams       ierr = PetscTypeCompare( (PetscObject)subpc, PETSC_GAMG_SMOOTHER, &isCheb ); CHKERRQ(ierr);
655f6536408SMark F. Adams       if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
656587fa25dSMark F. Adams       else{ /* eigen estimate 'emax' */
657587fa25dSMark F. Adams         KSP eksp; Mat Lmat = Aarr[level];
658fc4362bfSMark F. Adams         Vec bb, xx; PC pc;
659f6536408SMark F. Adams         const PCType type;
660038e3b61SMark F. Adams 
661f6536408SMark F. Adams         ierr = PCGetType( subpc, &type );   CHKERRQ(ierr);
6625745e0f5SMark F. Adams         ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
6635745e0f5SMark F. Adams         ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
664fc4362bfSMark F. Adams         {
665fc4362bfSMark F. Adams           PetscRandom    rctx;
666fc4362bfSMark F. Adams           ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
667fc4362bfSMark F. Adams           ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
668fc4362bfSMark F. Adams           ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
669fc4362bfSMark F. Adams           ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
6705b89ad90SMark F. Adams         }
671fc4362bfSMark F. Adams         ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
672e33ef3b1SMark F. Adams         ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr);
673038e3b61SMark F. Adams         ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
674fc4362bfSMark F. Adams         CHKERRQ(ierr);
675fc4362bfSMark F. Adams         ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
676fc4362bfSMark F. Adams 
677f6536408SMark F. Adams         ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
678f6536408SMark F. Adams         ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
679f6536408SMark F. Adams 
680f6536408SMark F. Adams         ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
681f6536408SMark F. Adams         ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
682f6536408SMark F. Adams         ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
683f6536408SMark F. Adams         ierr = PCSetType( pc, type ); CHKERRQ(ierr); /* should be same as eigen estimates op. */
684f6536408SMark F. Adams 
685fc4362bfSMark F. Adams         ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
686fc4362bfSMark F. Adams         ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
687fc4362bfSMark F. Adams         ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
688fc4362bfSMark F. Adams         ierr = VecDestroy( &xx );       CHKERRQ(ierr);
689fc4362bfSMark F. Adams         ierr = VecDestroy( &bb );       CHKERRQ(ierr);
690fc4362bfSMark F. Adams         ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
691f6536408SMark F. Adams 
692f6536408SMark F. Adams         if (pc_gamg->m_verbose) {
693f6536408SMark F. Adams           PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s PC setup max eigen=%e min=%e PC=%s\n",
694f6536408SMark F. Adams                       __FUNCT__,emax,emin,PETSC_GAMG_SMOOTHER);
695f6536408SMark F. Adams         }
696fc4362bfSMark F. Adams       }
697038e3b61SMark F. Adams       {
698038e3b61SMark F. Adams         PetscInt N1, N0, tt;
699038e3b61SMark F. Adams         ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
700038e3b61SMark F. Adams         ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
701f6536408SMark F. Adams         /* heuristic - is this crap? */
702f6536408SMark F. Adams         emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
703038e3b61SMark F. Adams         emax *= 1.05;
704038e3b61SMark F. Adams       }
705fc4362bfSMark F. Adams       ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
706f6536408SMark F. Adams     }
7075745e0f5SMark F. Adams   }
7085745e0f5SMark F. Adams   {
709ecd57171SMark F. Adams     /* coarse grid */
710ecd57171SMark F. Adams     KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
7115745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
712eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
71358471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
7145745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
715ecd57171SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
716ecd57171SMark F. Adams     ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
717ecd57171SMark F. Adams     ierr = PCSetUp( subpc ); CHKERRQ(ierr);
718ecd57171SMark F. Adams     ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
719ecd57171SMark F. Adams     assert(ii==1);
720ecd57171SMark F. Adams     ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
721ecd57171SMark F. Adams     ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
722fc4362bfSMark F. Adams   }
723737a81a9SMark F. Adams 
724fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
725eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
7265b89ad90SMark F. Adams   {
7275b89ad90SMark F. Adams     PetscBool galerkin;
728eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
7295b89ad90SMark F. Adams     if(galerkin){
7305b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
7315b89ad90SMark F. Adams     }
7325b89ad90SMark F. Adams   }
7335745e0f5SMark F. Adams 
73458471d46SMark F. Adams   /* plot levels - R/P */
73558471d46SMark F. Adams   if( PRINT_MATS ) {
73658471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-1;level>0;level--){
73758471d46SMark F. Adams       PetscViewer viewer;
73858471d46SMark F. Adams       char fname[32];
73903a628feSMark F. Adams       sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
74011e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
7415b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
742038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
7435b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
74403a628feSMark F. Adams       sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
745e33ef3b1SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
746e33ef3b1SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
747e33ef3b1SMark F. Adams       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
748e33ef3b1SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
7495b89ad90SMark F. Adams     }
75058471d46SMark F. Adams   }
75158471d46SMark F. Adams 
75258471d46SMark F. Adams   /* set interpolation between the levels, clean up */
75358471d46SMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
75458471d46SMark F. Adams        lidx<fine_level;
75558471d46SMark F. Adams        lidx++, level--){
75658471d46SMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
757587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
758587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
7595b89ad90SMark F. Adams   }
7605745e0f5SMark F. Adams 
7615b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
762eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
763eb07cef2SMark F. Adams   ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
76403a628feSMark F. Adams   a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
76558471d46SMark F. Adams 
76658471d46SMark F. Adams   {
76758471d46SMark F. Adams     KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
76858471d46SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
76958471d46SMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
77058471d46SMark F. Adams     ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
77158471d46SMark F. Adams   }
772e33ef3b1SMark F. Adams 
7735b89ad90SMark F. Adams   PetscFunctionReturn(0);
7745b89ad90SMark F. Adams }
7755b89ad90SMark F. Adams 
776eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7775b89ad90SMark F. Adams /*
7785b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7795b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7805b89ad90SMark F. Adams 
7815b89ad90SMark F. Adams    Input Parameter:
7825b89ad90SMark F. Adams .  pc - the preconditioner context
7835b89ad90SMark F. Adams 
7845b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7855b89ad90SMark F. Adams */
7865b89ad90SMark F. Adams #undef __FUNCT__
7875b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
7885b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
7895b89ad90SMark F. Adams {
7905b89ad90SMark F. Adams   PetscErrorCode  ierr;
7915b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
7925b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
7935b89ad90SMark F. Adams 
7945b89ad90SMark F. Adams   PetscFunctionBegin;
7955b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7965b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7975b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7985b89ad90SMark F. Adams   PetscFunctionReturn(0);
7995b89ad90SMark F. Adams }
8005b89ad90SMark F. Adams 
801676e1743SMark F. Adams 
802676e1743SMark F. Adams #undef __FUNCT__
803676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
804676e1743SMark F. Adams /*@
805676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
806676e1743SMark F. Adams    processor reduction.
807676e1743SMark F. Adams 
808676e1743SMark F. Adams    Not Collective on PC
809676e1743SMark F. Adams 
810676e1743SMark F. Adams    Input Parameters:
811676e1743SMark F. Adams .  pc - the preconditioner context
812676e1743SMark F. Adams 
813676e1743SMark F. Adams 
814676e1743SMark F. Adams    Options Database Key:
815676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
816676e1743SMark F. Adams 
817676e1743SMark F. Adams    Level: intermediate
818676e1743SMark F. Adams 
819676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
820676e1743SMark F. Adams 
821676e1743SMark F. Adams .seealso: ()
822676e1743SMark F. Adams @*/
823676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
824676e1743SMark F. Adams {
825676e1743SMark F. Adams   PetscErrorCode ierr;
826676e1743SMark F. Adams 
827676e1743SMark F. Adams   PetscFunctionBegin;
828676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
829676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
830676e1743SMark F. Adams   PetscFunctionReturn(0);
831676e1743SMark F. Adams }
832676e1743SMark F. Adams 
833676e1743SMark F. Adams EXTERN_C_BEGIN
834676e1743SMark F. Adams #undef __FUNCT__
835676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
836676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
837676e1743SMark F. Adams {
838c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
839c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
840676e1743SMark F. Adams 
841676e1743SMark F. Adams   PetscFunctionBegin;
842c20e4228SMark F. Adams   if(n>0) pc_gamg->m_min_eq_proc = n;
843676e1743SMark F. Adams   PetscFunctionReturn(0);
844676e1743SMark F. Adams }
845676e1743SMark F. Adams EXTERN_C_END
846676e1743SMark F. Adams 
847676e1743SMark F. Adams #undef __FUNCT__
848389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
849389730f3SMark F. Adams /*@
850389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
851389730f3SMark F. Adams 
852389730f3SMark F. Adams  Collective on PC
853389730f3SMark F. Adams 
854389730f3SMark F. Adams    Input Parameters:
855389730f3SMark F. Adams .  pc - the preconditioner context
856389730f3SMark F. Adams 
857389730f3SMark F. Adams 
858389730f3SMark F. Adams    Options Database Key:
859389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
860389730f3SMark F. Adams 
861389730f3SMark F. Adams    Level: intermediate
862389730f3SMark F. Adams 
863389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
864389730f3SMark F. Adams 
865389730f3SMark F. Adams .seealso: ()
866389730f3SMark F. Adams  @*/
867389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
868389730f3SMark F. Adams {
869389730f3SMark F. Adams   PetscErrorCode ierr;
870389730f3SMark F. Adams 
871389730f3SMark F. Adams   PetscFunctionBegin;
872389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
873389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
874389730f3SMark F. Adams   PetscFunctionReturn(0);
875389730f3SMark F. Adams }
876389730f3SMark F. Adams 
877389730f3SMark F. Adams EXTERN_C_BEGIN
878389730f3SMark F. Adams #undef __FUNCT__
879389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
880389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
881389730f3SMark F. Adams {
882389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
883389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
884389730f3SMark F. Adams 
885389730f3SMark F. Adams   PetscFunctionBegin;
886389730f3SMark F. Adams   if(n>0) pc_gamg->m_coarse_eq_limit = n;
887389730f3SMark F. Adams   PetscFunctionReturn(0);
888389730f3SMark F. Adams }
889389730f3SMark F. Adams EXTERN_C_END
890389730f3SMark F. Adams 
891389730f3SMark F. Adams #undef __FUNCT__
892676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning"
893676e1743SMark F. Adams /*@
894676e1743SMark F. Adams    PCGAMGAvoidRepartitioning - Do not repartition the coarse grids
895676e1743SMark F. Adams 
896676e1743SMark F. Adams    Collective on PC
897676e1743SMark F. Adams 
898676e1743SMark F. Adams    Input Parameters:
899676e1743SMark F. Adams .  pc - the preconditioner context
900676e1743SMark F. Adams 
901676e1743SMark F. Adams 
902676e1743SMark F. Adams    Options Database Key:
903676e1743SMark F. Adams .  -pc_gamg_avoid_repartitioning
904676e1743SMark F. Adams 
905676e1743SMark F. Adams    Level: intermediate
906676e1743SMark F. Adams 
907676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
908676e1743SMark F. Adams 
909676e1743SMark F. Adams .seealso: ()
910676e1743SMark F. Adams @*/
911676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning(PC pc, PetscBool n)
912676e1743SMark F. Adams {
913676e1743SMark F. Adams   PetscErrorCode ierr;
914676e1743SMark F. Adams 
915676e1743SMark F. Adams   PetscFunctionBegin;
916676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
917676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGAvoidRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
918676e1743SMark F. Adams   PetscFunctionReturn(0);
919676e1743SMark F. Adams }
920676e1743SMark F. Adams 
921676e1743SMark F. Adams EXTERN_C_BEGIN
922676e1743SMark F. Adams #undef __FUNCT__
923676e1743SMark F. Adams #define __FUNCT__ "PCGAMGAvoidRepartitioning_GAMG"
924676e1743SMark F. Adams PetscErrorCode PCGAMGAvoidRepartitioning_GAMG(PC pc, PetscBool n)
925676e1743SMark F. Adams {
926c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
927c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
928676e1743SMark F. Adams 
929676e1743SMark F. Adams   PetscFunctionBegin;
930c20e4228SMark F. Adams   pc_gamg->m_avoid_repart = n;
931676e1743SMark F. Adams   PetscFunctionReturn(0);
932676e1743SMark F. Adams }
933676e1743SMark F. Adams EXTERN_C_END
934676e1743SMark F. Adams 
935676e1743SMark F. Adams #undef __FUNCT__
9364ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
9374ef23d27SMark F. Adams /*@
9384ef23d27SMark F. Adams    PCGAMGSetNlevels -
9394ef23d27SMark F. Adams 
9404ef23d27SMark F. Adams    Not collective on PC
9414ef23d27SMark F. Adams 
9424ef23d27SMark F. Adams    Input Parameters:
9434ef23d27SMark F. Adams .  pc - the preconditioner context
9444ef23d27SMark F. Adams 
9454ef23d27SMark F. Adams 
9464ef23d27SMark F. Adams    Options Database Key:
9474ef23d27SMark F. Adams .  -pc_mg_levels
9484ef23d27SMark F. Adams 
9494ef23d27SMark F. Adams    Level: intermediate
9504ef23d27SMark F. Adams 
9514ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
9524ef23d27SMark F. Adams 
9534ef23d27SMark F. Adams .seealso: ()
9544ef23d27SMark F. Adams @*/
9554ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
9564ef23d27SMark F. Adams {
9574ef23d27SMark F. Adams   PetscErrorCode ierr;
9584ef23d27SMark F. Adams 
9594ef23d27SMark F. Adams   PetscFunctionBegin;
9604ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9614ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
9624ef23d27SMark F. Adams   PetscFunctionReturn(0);
9634ef23d27SMark F. Adams }
9644ef23d27SMark F. Adams 
9654ef23d27SMark F. Adams EXTERN_C_BEGIN
9664ef23d27SMark F. Adams #undef __FUNCT__
9674ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
9684ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
9694ef23d27SMark F. Adams {
9704ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
9714ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
9724ef23d27SMark F. Adams 
9734ef23d27SMark F. Adams   PetscFunctionBegin;
9744ef23d27SMark F. Adams   pc_gamg->m_Nlevels = n;
9754ef23d27SMark F. Adams   PetscFunctionReturn(0);
9764ef23d27SMark F. Adams }
9774ef23d27SMark F. Adams EXTERN_C_END
9784ef23d27SMark F. Adams 
9794ef23d27SMark F. Adams #undef __FUNCT__
9803542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
9813542efc5SMark F. Adams /*@
9823542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
9833542efc5SMark F. Adams 
9843542efc5SMark F. Adams    Not collective on PC
9853542efc5SMark F. Adams 
9863542efc5SMark F. Adams    Input Parameters:
9873542efc5SMark F. Adams .  pc - the preconditioner context
9883542efc5SMark F. Adams 
9893542efc5SMark F. Adams 
9903542efc5SMark F. Adams    Options Database Key:
9913542efc5SMark F. Adams .  -pc_gamg_threshold
9923542efc5SMark F. Adams 
9933542efc5SMark F. Adams    Level: intermediate
9943542efc5SMark F. Adams 
9953542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
9963542efc5SMark F. Adams 
9973542efc5SMark F. Adams .seealso: ()
9983542efc5SMark F. Adams @*/
9993542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
10003542efc5SMark F. Adams {
10013542efc5SMark F. Adams   PetscErrorCode ierr;
10023542efc5SMark F. Adams 
10033542efc5SMark F. Adams   PetscFunctionBegin;
10043542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10053542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
10063542efc5SMark F. Adams   PetscFunctionReturn(0);
10073542efc5SMark F. Adams }
10083542efc5SMark F. Adams 
10093542efc5SMark F. Adams EXTERN_C_BEGIN
10103542efc5SMark F. Adams #undef __FUNCT__
10113542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
10123542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
10133542efc5SMark F. Adams {
1014c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1015c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
10163542efc5SMark F. Adams 
10173542efc5SMark F. Adams   PetscFunctionBegin;
1018c20e4228SMark F. Adams   pc_gamg->m_threshold = n;
10193542efc5SMark F. Adams   PetscFunctionReturn(0);
10203542efc5SMark F. Adams }
10213542efc5SMark F. Adams EXTERN_C_END
10223542efc5SMark F. Adams 
10233542efc5SMark F. Adams #undef __FUNCT__
1024676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType"
1025676e1743SMark F. Adams /*@
1026676e1743SMark F. Adams    PCGAMGSetSolverType - Set solution method.
1027676e1743SMark F. Adams 
1028676e1743SMark F. Adams    Collective on PC
1029676e1743SMark F. Adams 
1030676e1743SMark F. Adams    Input Parameters:
1031676e1743SMark F. Adams .  pc - the preconditioner context
1032676e1743SMark F. Adams 
1033676e1743SMark F. Adams 
1034676e1743SMark F. Adams    Options Database Key:
10353542efc5SMark F. Adams .  -pc_gamg_type
1036676e1743SMark F. Adams 
1037676e1743SMark F. Adams    Level: intermediate
1038676e1743SMark F. Adams 
1039676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1040676e1743SMark F. Adams 
1041676e1743SMark F. Adams .seealso: ()
1042676e1743SMark F. Adams @*/
1043676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz )
1044676e1743SMark F. Adams {
1045676e1743SMark F. Adams   PetscErrorCode ierr;
1046676e1743SMark F. Adams 
1047676e1743SMark F. Adams   PetscFunctionBegin;
1048676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1049676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz));
1050676e1743SMark F. Adams   CHKERRQ(ierr);
1051676e1743SMark F. Adams   PetscFunctionReturn(0);
1052676e1743SMark F. Adams }
1053676e1743SMark F. Adams 
1054676e1743SMark F. Adams EXTERN_C_BEGIN
1055676e1743SMark F. Adams #undef __FUNCT__
1056676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG"
1057676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz )
1058676e1743SMark F. Adams {
1059676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1060676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1061676e1743SMark F. Adams 
1062676e1743SMark F. Adams   PetscFunctionBegin;
1063676e1743SMark F. Adams   if(sz < 64) strcpy(pc_gamg->m_type,str);
1064676e1743SMark F. Adams   PetscFunctionReturn(0);
1065676e1743SMark F. Adams }
1066676e1743SMark F. Adams EXTERN_C_END
1067676e1743SMark F. Adams 
10685b89ad90SMark F. Adams #undef __FUNCT__
10695b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
10705b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
10715b89ad90SMark F. Adams {
1072676e1743SMark F. Adams   PetscErrorCode  ierr;
1073676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1074676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1075676e1743SMark F. Adams   PetscBool flag;
10765b89ad90SMark F. Adams 
10775b89ad90SMark F. Adams   PetscFunctionBegin;
1078676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1079676e1743SMark F. Adams   {
1080676e1743SMark F. Adams     ierr = PetscOptionsString("-pc_gamg_type",
1081470e26f8SMark F. Adams                               "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid (default)",
1082676e1743SMark F. Adams                               "PCGAMGSetSolverType",
1083676e1743SMark F. Adams                               pc_gamg->m_type,
1084676e1743SMark F. Adams                               pc_gamg->m_type,
1085676e1743SMark F. Adams                               64,
1086676e1743SMark F. Adams                               &flag );
1087676e1743SMark F. Adams     CHKERRQ(ierr);
1088*d8c859f0SMark F. Adams     if( flag && pc_gamg->m_data != 0 ) {
1089*d8c859f0SMark F. Adams       if( (strcmp(pc_gamg->m_type,"sa")==0 && pc_gamg->m_method != 2) ||
1090*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"pa")==0 && pc_gamg->m_method != 1) ||
1091*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"geo")==0 && pc_gamg->m_method != 0) ) {
1092*d8c859f0SMark F. Adams         SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "PCSetFromOptions called of PCSetCoordinates (with new method, after data was created)");
1093*d8c859f0SMark F. Adams       }
1094*d8c859f0SMark F. Adams     }
1095bed7c2b9SMark F. Adams 
1096bed7c2b9SMark F. Adams     /* -pc_gamg_verbose */
1097bed7c2b9SMark 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);
10982a44bfbaSMark F. Adams 
1099f6536408SMark F. Adams     pc_gamg->m_method = 1; /* default to plane aggregation */
1100f6536408SMark F. Adams     if (flag ) {
1101f6536408SMark F. Adams       if( strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2;
1102f6536408SMark F. Adams       else if( strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1;
1103f6536408SMark F. Adams       else if( strcmp(pc_gamg->m_type,"geo") == 0) pc_gamg->m_method = 0;
1104f6536408SMark F. Adams       else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Invalid gamg type: %s",pc_gamg->m_type);
1105f6536408SMark F. Adams     }
1106c20e4228SMark F. Adams     /* -pc_gamg_avoid_repartitioning */
1107c20e4228SMark F. Adams     pc_gamg->m_avoid_repart = PETSC_FALSE;
1108676e1743SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_avoid_repartitioning",
1109676e1743SMark F. Adams                             "Do not repartion coarse grids (false)",
1110676e1743SMark F. Adams                             "PCGAMGAvoidRepartitioning",
1111c20e4228SMark F. Adams                             pc_gamg->m_avoid_repart,
1112c20e4228SMark F. Adams                             &pc_gamg->m_avoid_repart,
1113676e1743SMark F. Adams                             &flag);
1114676e1743SMark F. Adams     CHKERRQ(ierr);
1115676e1743SMark F. Adams 
1116c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1117c20e4228SMark F. Adams     pc_gamg->m_min_eq_proc = 600;
1118676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1119676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1120676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
1121c20e4228SMark F. Adams                            pc_gamg->m_min_eq_proc,
1122c20e4228SMark F. Adams                            &pc_gamg->m_min_eq_proc,
1123676e1743SMark F. Adams                            &flag );
1124676e1743SMark F. Adams     CHKERRQ(ierr);
11253542efc5SMark F. Adams 
1126389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1127389730f3SMark F. Adams     pc_gamg->m_coarse_eq_limit = 1500;
1128389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1129389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1130389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
1131389730f3SMark F. Adams                            pc_gamg->m_coarse_eq_limit,
1132389730f3SMark F. Adams                            &pc_gamg->m_coarse_eq_limit,
1133389730f3SMark F. Adams                            &flag );
1134389730f3SMark F. Adams     CHKERRQ(ierr);
1135389730f3SMark F. Adams 
1136c20e4228SMark F. Adams     /* -pc_gamg_threshold */
1137c20e4228SMark F. Adams     pc_gamg->m_threshold = 0.05;
11383542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
11393542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
11403542efc5SMark F. Adams                             "PCGAMGSetThreshold",
1141c20e4228SMark F. Adams                             pc_gamg->m_threshold,
1142c20e4228SMark F. Adams                             &pc_gamg->m_threshold,
11433542efc5SMark F. Adams                             &flag );
11443542efc5SMark F. Adams     CHKERRQ(ierr);
11454ef23d27SMark F. Adams 
11464ef23d27SMark F. Adams     pc_gamg->m_Nlevels = GAMG_MAXLEVELS;
11474ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
11484ef23d27SMark F. Adams                            "Set number of MG levels",
11494ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
11504ef23d27SMark F. Adams                            pc_gamg->m_Nlevels,
11514ef23d27SMark F. Adams                            &pc_gamg->m_Nlevels,
11524ef23d27SMark F. Adams                            &flag );
1153676e1743SMark F. Adams   }
1154676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1155676e1743SMark F. Adams 
11565b89ad90SMark F. Adams   PetscFunctionReturn(0);
11575b89ad90SMark F. Adams }
11585b89ad90SMark F. Adams 
11595b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
11605b89ad90SMark F. Adams /*MC
1161dc76db92SMark F. Adams      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning. This preconditioner currently has two
1162dc76db92SMark F. Adams            AMG methods: 1) an unstructured geometric method, which requires that you provide coordinates for each
1163dc76db92SMark F. Adams            vertex, and 2) smoothed aggregation.  Smoothed aggregation (SA) can work without coordinates but it
1164dc76db92SMark F. Adams            will generate some common non-trivial null spaces if coordinates are provided.  The input fine grid matrix
1165dc76db92SMark F. Adams            must have the block size set for 'system' problems (with multiple dof per vertex/cell) to work properly.
1166dc76db92SMark F. Adams            SA will generate rotational rigid body mode null space vectors, in addition to the trivial translational
1167dc76db92SMark F. Adams            modes, when coordinates are provide in 2D and 3D.
11685b89ad90SMark F. Adams 
1169280d9858SJed Brown    Options Database Keys:
11705b89ad90SMark F. Adams    Multigrid options(inherited)
1171280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1172280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1173280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1174280d9858SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
11755b89ad90SMark F. Adams 
11765b89ad90SMark F. Adams   Level: intermediate
1177280d9858SJed Brown 
11785b89ad90SMark F. Adams   Concepts: multigrid
11795b89ad90SMark F. Adams 
11805b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1181280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
11825b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
11835b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
11845b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11855b89ad90SMark F. Adams M*/
11865b89ad90SMark F. Adams 
11875b89ad90SMark F. Adams EXTERN_C_BEGIN
11885b89ad90SMark F. Adams #undef __FUNCT__
11895b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
11905b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
11915b89ad90SMark F. Adams {
11925b89ad90SMark F. Adams   PetscErrorCode  ierr;
11935b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
11945b89ad90SMark F. Adams   PC_MG           *mg;
11955ef31b24SMark F. Adams   PetscClassId     cookie;
11965ee9c036SSatish Balay #if defined PETSC_USE_LOG
11975ee9c036SSatish Balay   static int count = 0;
11985ee9c036SSatish Balay #endif
11995b89ad90SMark F. Adams 
12005b89ad90SMark F. Adams   PetscFunctionBegin;
12015b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
12025b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
12035b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
12045b89ad90SMark F. Adams 
12055b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
12065b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
120703a628feSMark F. Adams   pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0;
12085b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
12095b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
12105b89ad90SMark F. Adams 
12115b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
12125b89ad90SMark F. Adams 
12135b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
12145b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
12155b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
12165b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
12175b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
12185b89ad90SMark F. Adams 
12195b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12205b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
12215b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
1222c97e1df0SMark F. Adams 					    PCSetCoordinates_GAMG);
1223c97e1df0SMark F. Adams   CHKERRQ(ierr);
1224c97e1df0SMark F. Adams 
1225676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1226676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1227676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1228676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1229676e1743SMark F. Adams   CHKERRQ(ierr);
1230676e1743SMark F. Adams 
1231676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1232389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_C",
1233389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_GAMG",
1234389730f3SMark F. Adams 					    PCGAMGSetCoarseEqLim_GAMG);
1235389730f3SMark F. Adams   CHKERRQ(ierr);
1236389730f3SMark F. Adams 
1237389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1238676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_C",
1239676e1743SMark F. Adams 					    "PCGAMGAvoidRepartitioning_GAMG",
1240676e1743SMark F. Adams 					    PCGAMGAvoidRepartitioning_GAMG);
1241676e1743SMark F. Adams   CHKERRQ(ierr);
1242676e1743SMark F. Adams 
1243676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12443542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
12453542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
12463542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
12473542efc5SMark F. Adams   CHKERRQ(ierr);
12483542efc5SMark F. Adams 
12493542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1250676e1743SMark F. Adams 					    "PCGAMGSetSolverType_C",
1251676e1743SMark F. Adams 					    "PCGAMGSetSolverType_GAMG",
1252676e1743SMark F. Adams 					    PCGAMGSetSolverType_GAMG);
1253676e1743SMark F. Adams   CHKERRQ(ierr);
1254c97e1df0SMark F. Adams 
1255b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
1256785cba28SMark F. Adams   if( count++ == 0 ) {
12575ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
1258b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
12592a44bfbaSMark F. Adams     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
12602a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]);
12612a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]);
12622a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]);
1263b4fbaa2aSMark F. Adams     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1264b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1265b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1266b4fbaa2aSMark F. Adams     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1267b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_events[SET7]);
1268b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]); */
1269b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1270b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1271b4fbaa2aSMark F. Adams     PetscLogEventRegister("  PL repartition", cookie, &gamg_setup_events[SET12]);
1272f852f58cSMark F. Adams 
1273fca9ed99SMark F. Adams     /* PetscLogEventRegister("  PL 1", cookie, &gamg_setup_events[SET13]); */
1274fca9ed99SMark F. Adams     /* PetscLogEventRegister("  PL 2", cookie, &gamg_setup_events[SET14]); */
1275fca9ed99SMark F. Adams     /* PetscLogEventRegister("  PL 3", cookie, &gamg_setup_events[SET15]); */
1276f852f58cSMark F. Adams 
1277b4fbaa2aSMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1278b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1279b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
12805ef31b24SMark F. Adams 
1281b4fbaa2aSMark F. Adams     /* create timer stages */
1282b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1283b4fbaa2aSMark F. Adams     {
1284b4fbaa2aSMark F. Adams       char str[32];
1285b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1286b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1287b4fbaa2aSMark F. Adams       PetscInt lidx;
1288b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1289b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1290b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1291b4fbaa2aSMark F. Adams       }
1292b4fbaa2aSMark F. Adams     }
1293b4fbaa2aSMark F. Adams #endif
1294b4fbaa2aSMark F. Adams   }
1295b4fbaa2aSMark F. Adams #endif
12965b89ad90SMark F. Adams   PetscFunctionReturn(0);
12975b89ad90SMark F. Adams }
12985b89ad90SMark F. Adams EXTERN_C_END
1299