xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision a147abb01e48e1f017940485f0dcd6ad7c6cfce1)
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 /*
125a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
126a147abb0SMark F. Adams      of active processors.
1275b89ad90SMark F. Adams 
1285b89ad90SMark F. Adams    Input Parameter:
1298263b398SMark F. Adams    . a_pc - parameters
1303530afc2SMark F. Adams    . a_Amat_fine - matrix on this fine (k) level
131d3d6bff4SMark F. Adams    . a_ndata_rows - size of data to move (coarse grid)
132d3d6bff4SMark F. Adams    . a_ndata_cols - size of data to move (coarse grid)
1338263b398SMark F. Adams    . a_cbs - coarse block size
1348263b398SMark F. Adams    . a_isLast -
1353530afc2SMark F. Adams    In/Output Parameter:
1363530afc2SMark F. Adams    . a_P_inout - prolongation operator to the next level (k-1)
137eb07cef2SMark F. Adams    . a_coarse_data - data that need to be moved
138afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
13911e60469SMark F. Adams    Output Parameter:
1403530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1415b89ad90SMark F. Adams */
1425cb416c2SMark F. Adams 
1435b89ad90SMark F. Adams #undef __FUNCT__
1448263b398SMark F. Adams #define __FUNCT__ "createLevel"
1458263b398SMark F. Adams PetscErrorCode createLevel( const PC a_pc,
1468263b398SMark F. Adams                             const Mat a_Amat_fine,
1478263b398SMark F. Adams                             const PetscInt a_ndata_rows,
1488263b398SMark F. Adams                             const PetscInt a_ndata_cols,
1498263b398SMark F. Adams                             const PetscInt a_cbs,
1508263b398SMark F. Adams                             const PetscBool a_isLast,
1513530afc2SMark F. Adams                             Mat *a_P_inout,
152eb07cef2SMark F. Adams                             PetscReal **a_coarse_data,
153afc97cdcSMark F. Adams                             PetscMPIInt *a_nactive_proc,
154389730f3SMark F. Adams                             Mat *a_Amat_crs
15511e60469SMark F. Adams                             )
1565b89ad90SMark F. Adams {
1578263b398SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
158486a8d0bSJed Brown   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1598263b398SMark F. Adams   const PetscBool  repart = pc_gamg->m_repart;
160389730f3SMark F. Adams   const PetscInt   min_eq_proc = pc_gamg->m_min_eq_proc, coarse_max = pc_gamg->m_coarse_eq_limit;
1615b89ad90SMark F. Adams   PetscErrorCode   ierr;
162038e3b61SMark F. Adams   Mat              Cmat,Pnew,Pold=*a_P_inout;
16311e60469SMark F. Adams   IS               new_indices,isnum;
1643530afc2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_Amat_fine)->comm;
165fa31c87dSMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive;
166fa31c87dSMark F. Adams   PetscInt         neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0;
1675b89ad90SMark F. Adams 
1685b89ad90SMark F. Adams   PetscFunctionBegin;
16911e60469SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
17011e60469SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
171f6536408SMark F. Adams 
17211e60469SMark F. Adams   /* RAP */
17396434bc1SMark F. Adams #ifdef USE_R
17496434bc1SMark F. Adams   /* make R wih brute force for now */
17596434bc1SMark F. Adams   ierr = MatTranspose( Pold, Pnew );
17696434bc1SMark F. Adams   ierr = MatDestroy( &Pold );  CHKERRQ(ierr);
17796434bc1SMark F. Adams   ierr = MatRARt( a_Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
17896434bc1SMark F. Adams   Pold = Pnew;
17996434bc1SMark F. Adams #else
180038e3b61SMark F. Adams   ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
18196434bc1SMark F. Adams #endif
182038e3b61SMark F. Adams   ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
183038e3b61SMark F. Adams   ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr);
184038e3b61SMark F. Adams   ncrs0 = (Iend0-Istart0)/a_cbs; assert((Iend0-Istart0)%a_cbs == 0);
185038e3b61SMark F. Adams 
186f852f58cSMark F. Adams   /* get number of PEs to make active, reduce */
187f852f58cSMark F. Adams   ierr = MatGetSize( Cmat, &neq, &NN );  CHKERRQ(ierr);
188f852f58cSMark F. Adams   new_npe = neq/min_eq_proc; /* hardwire min. number of eq/proc */
189f852f58cSMark F. Adams   if( new_npe == 0 || neq < coarse_max ) new_npe = 1;
190f852f58cSMark F. Adams   else if (new_npe >= *a_nactive_proc ) new_npe = *a_nactive_proc; /* no change, rare */
1918263b398SMark F. Adams   if( a_isLast ) new_npe = 1;
192f852f58cSMark F. Adams 
1938263b398SMark F. Adams   if( !repart && !(new_npe == 1 && *a_nactive_proc != 1) ) {
194a147abb0SMark F. Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction */
19522063be5SMark F. Adams   }
19622063be5SMark F. Adams   else {
19722063be5SMark F. Adams     /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */
1985ef31b24SMark F. Adams     Mat              adj;
19922063be5SMark F. Adams     const PetscInt *idx,data_sz=a_ndata_rows*a_ndata_cols;
200b4fbaa2aSMark F. Adams     const PetscInt  stride0=ncrs0*a_ndata_rows;
2018263b398SMark F. Adams     PetscInt        is_sz,*newproc_idx,ii,jj,kk,strideNew,*tidx;
202c9a0b8beSMark F. Adams     /* create sub communicator  */
20371d60426SBarry Smith     MPI_Comm        cm;
204afc97cdcSMark F. Adams     MPI_Group       wg, g2;
205a147abb0SMark F. Adams     PetscInt       *counts,inpe,targetPE;
206fd3c6afaSMark F. Adams     PetscMPIInt    *ranks;
20722063be5SMark F. Adams     IS              isscat;
20822063be5SMark F. Adams     PetscScalar    *array;
20922063be5SMark F. Adams     Vec             src_crd, dest_crd;
21022063be5SMark F. Adams     PetscReal      *data = *a_coarse_data;
21122063be5SMark F. Adams     VecScatter      vecscat;
212b4fbaa2aSMark F. Adams     IS              isnewproc;
213e33ef3b1SMark F. Adams 
21492a756f0SMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscMPIInt), &ranks ); CHKERRQ(ierr);
215fd3c6afaSMark F. Adams     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
21692a756f0SMark F. Adams 
217fd3c6afaSMark F. Adams     ierr = MPI_Allgather( &ncrs0, 1, MPIU_INT, counts, 1, MPIU_INT, wcomm ); CHKERRQ(ierr);
218afc97cdcSMark F. Adams     assert(counts[mype]==ncrs0);
219afc97cdcSMark F. Adams     /* count real active pes */
220a147abb0SMark F. Adams     for( nactive = jj = 0, targetPE = -1 ; jj < npe ; jj++) {
221afc97cdcSMark F. Adams       if( counts[jj] != 0 ) {
22222063be5SMark F. Adams 	ranks[nactive++] = jj;
223a147abb0SMark F. Adams         if( targetPE == -1 ) targetPE = jj; /* find lowest active PE */
224afc97cdcSMark F. Adams       }
225afc97cdcSMark F. Adams     }
2263303bcf9Sadams 
227a147abb0SMark F. Adams     if (nactive < new_npe) {
228a147abb0SMark F. Adams       if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s reducing number of target PEs: %d --> %d (?)\n",mype,__FUNCT__,new_npe,nactive);
229a147abb0SMark F. Adams       new_npe = nactive; /* this can happen with empty input procs */
230a147abb0SMark F. Adams     }
23122063be5SMark F. Adams 
232bed7c2b9SMark F. Adams     if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s npe (active): %d --> %d. new npe = %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,nactive,new_npe,neq);
23358471d46SMark F. Adams 
23422063be5SMark F. Adams     *a_nactive_proc = new_npe; /* output */
2352f03bc48SMark F. Adams 
236afc97cdcSMark F. Adams     ierr = MPI_Comm_group( wcomm, &wg ); CHKERRQ(ierr);
23722063be5SMark F. Adams     ierr = MPI_Group_incl( wg, nactive, ranks, &g2 ); CHKERRQ(ierr);
238afc97cdcSMark F. Adams     ierr = MPI_Comm_create( wcomm, g2, &cm ); CHKERRQ(ierr);
239b4fbaa2aSMark F. Adams 
240afc97cdcSMark F. Adams     ierr = MPI_Group_free( &wg );                            CHKERRQ(ierr);
241afc97cdcSMark F. Adams     ierr = MPI_Group_free( &g2 );                            CHKERRQ(ierr);
242c9a0b8beSMark F. Adams 
2435ef31b24SMark F. Adams     /* MatPartitioningApply call MatConvert, which is collective */
244b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
245b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
246b4fbaa2aSMark F. Adams #endif
247038e3b61SMark F. Adams     if( a_cbs == 1) {
248038e3b61SMark F. Adams       ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
249eb07cef2SMark F. Adams     }
250eb07cef2SMark F. Adams     else{
251038e3b61SMark F. Adams       /* make a scalar matrix to partition */
252eb07cef2SMark F. Adams       Mat tMat;
25358471d46SMark F. Adams       PetscInt ncols,jj,Ii;
254b4fbaa2aSMark F. Adams       const PetscScalar *vals;
255b4fbaa2aSMark F. Adams       const PetscInt *idx;
25658471d46SMark F. Adams       PetscInt *d_nnz;
2575ee9c036SSatish Balay       static int llev = 0;
258b4fbaa2aSMark F. Adams 
25958471d46SMark F. Adams       ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
26058471d46SMark F. Adams       for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += a_cbs, jj++ ) {
26158471d46SMark F. Adams         ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
26258471d46SMark F. Adams         d_nnz[jj] = ncols/a_cbs;
26358471d46SMark F. Adams         if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0;
26458471d46SMark F. Adams         ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
26558471d46SMark F. Adams       }
2666876a03eSMark F. Adams 
267eb07cef2SMark F. Adams       ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0,
268eb07cef2SMark F. Adams                               PETSC_DETERMINE, PETSC_DETERMINE,
26958471d46SMark F. Adams                               0, d_nnz, 0, d_nnz,
270eb07cef2SMark F. Adams                               &tMat );
2716876a03eSMark F. Adams       CHKERRQ(ierr);
27258471d46SMark F. Adams       ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
273eb07cef2SMark F. Adams 
27422063be5SMark F. Adams       for ( ii = Istart0; ii < Iend0; ii++ ) {
27522063be5SMark F. Adams         PetscInt dest_row = ii/a_cbs;
27622063be5SMark F. Adams         ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
277eb07cef2SMark F. Adams         for( jj = 0 ; jj < ncols ; jj++ ){
278038e3b61SMark F. Adams           PetscInt dest_col = idx[jj]/a_cbs;
279eb07cef2SMark F. Adams           PetscScalar v = 1.0;
280eb07cef2SMark F. Adams           ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
281eb07cef2SMark F. Adams         }
28222063be5SMark F. Adams         ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
283eb07cef2SMark F. Adams       }
284eb07cef2SMark F. Adams       ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
285eb07cef2SMark F. Adams       ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
286eb07cef2SMark F. Adams 
287b4fbaa2aSMark F. Adams       if( llev++ == -1 ) {
288b4fbaa2aSMark F. Adams 	PetscViewer viewer; char fname[32];
289b4fbaa2aSMark F. Adams 	sprintf(fname,"part_mat_%d.mat",llev);
290b4fbaa2aSMark F. Adams 	PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
291b4fbaa2aSMark F. Adams 	ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
292b4fbaa2aSMark F. Adams 	ierr = PetscViewerDestroy( &viewer );
293b4fbaa2aSMark F. Adams       }
294b4fbaa2aSMark F. Adams 
295eb07cef2SMark F. Adams       ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
296eb07cef2SMark F. Adams 
297eb07cef2SMark F. Adams       ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
298eb07cef2SMark F. Adams     }
299f150b916SMark F. Adams 
300afc97cdcSMark F. Adams     if( ncrs0 != 0 ){
301b4fbaa2aSMark F. Adams       const PetscInt *is_idx;
302b4fbaa2aSMark F. Adams       MatPartitioning  mpart;
303a4b7d37bSMark F. Adams       IS proc_is;
3045ef31b24SMark F. Adams       /* hack to fix global data that pmetis.c uses in 'adj' */
30522063be5SMark F. Adams       for( nactive = jj = 0 ; jj < npe ; jj++ ) {
306afc97cdcSMark F. Adams 	if( counts[jj] != 0 ) {
30722063be5SMark F. Adams 	  adj->rmap->range[nactive++] = adj->rmap->range[jj];
308afc97cdcSMark F. Adams 	}
3095ef31b24SMark F. Adams       }
31022063be5SMark F. Adams       adj->rmap->range[nactive] = adj->rmap->range[npe];
3112f03bc48SMark F. Adams 
3128263b398SMark F. Adams       if( new_npe == 1 ) {
313a147abb0SMark F. Adams         /* simply dump on one proc w/o partitioner for final reduction */
3148263b398SMark F. Adams         ierr = MatGetLocalSize( adj, &is_sz, &ii );  CHKERRQ(ierr);
315a147abb0SMark F. Adams         ierr = ISCreateStride( cm, is_sz, targetPE, 0, &proc_is );  CHKERRQ(ierr);
31659a0be82SJed Brown       } else {
31759a0be82SJed Brown         char prefix[256];
31859a0be82SJed Brown         const char *pcpre;
31971d60426SBarry Smith         ierr = MatPartitioningCreate( cm, &mpart ); CHKERRQ(ierr);
3205ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
32159a0be82SJed Brown         ierr = PCGetOptionsPrefix(a_pc,&pcpre);CHKERRQ(ierr);
32259a0be82SJed Brown         ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
32359a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
32411e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
3255ef31b24SMark F. Adams         ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
326a4b7d37bSMark F. Adams         ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr);
32711e60469SMark F. Adams         ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
3288263b398SMark F. Adams       }
3295ef31b24SMark F. Adams       /* collect IS info */
330a4b7d37bSMark F. Adams       ierr = ISGetLocalSize( proc_is, &is_sz );       CHKERRQ(ierr);
3318263b398SMark F. Adams       ierr = PetscMalloc( a_cbs*is_sz*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr);
332a4b7d37bSMark F. Adams       ierr = ISGetIndices( proc_is, &is_idx );        CHKERRQ(ierr);
333a147abb0SMark F. Adams 
334a147abb0SMark F. Adams       NN = 1; /* bring to "front" of machine */
335a147abb0SMark F. Adams       /*NN = npe/new_npe;*/ /* spread partitioning across machine */
336b4fbaa2aSMark F. Adams       for( kk = jj = 0 ; kk < is_sz ; kk++ ){
337afc97cdcSMark F. Adams         for( ii = 0 ; ii < a_cbs ; ii++, jj++ ) {
3388263b398SMark F. Adams           newproc_idx[jj] = is_idx[kk] * NN; /* distribution */
339eb07cef2SMark F. Adams         }
3405ef31b24SMark F. Adams       }
341a4b7d37bSMark F. Adams       ierr = ISRestoreIndices( proc_is, &is_idx );     CHKERRQ(ierr);
342a4b7d37bSMark F. Adams       ierr = ISDestroy( &proc_is );                    CHKERRQ(ierr);
34371d60426SBarry Smith       ierr = MPI_Comm_free( &cm );              CHKERRQ(ierr);
344b4fbaa2aSMark F. Adams 
345b4fbaa2aSMark F. Adams       is_sz *= a_cbs;
3465ef31b24SMark F. Adams     }
3475ef31b24SMark F. Adams     else{
3488263b398SMark F. Adams       newproc_idx = 0;
3495ef31b24SMark F. Adams       is_sz = 0;
3505ef31b24SMark F. Adams     }
351f150b916SMark F. Adams 
3525ef31b24SMark F. Adams     ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
3538263b398SMark F. Adams     ierr = ISCreateGeneral( wcomm, is_sz, newproc_idx, PETSC_COPY_VALUES, &isnewproc );
3548263b398SMark F. Adams     if( newproc_idx != 0 ) {
3558263b398SMark F. Adams       ierr = PetscFree( newproc_idx );  CHKERRQ(ierr);
3565ef31b24SMark F. Adams     }
357e33ef3b1SMark F. Adams 
35811e60469SMark F. Adams     /*
35911e60469SMark F. Adams      Create an index set from the isnewproc index set to indicate the mapping TO
36011e60469SMark F. Adams      */
36111e60469SMark F. Adams     ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr);
36211e60469SMark F. Adams     /*
36311e60469SMark F. Adams      Determine how many elements are assigned to each processor
36411e60469SMark F. Adams      */
365fd3c6afaSMark F. Adams     inpe = npe;
366fd3c6afaSMark F. Adams     ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr);
36711e60469SMark F. Adams     ierr = ISDestroy( &isnewproc );                       CHKERRQ(ierr);
368038e3b61SMark F. Adams     ncrs_new = counts[mype]/a_cbs;
36922063be5SMark F. Adams     strideNew = ncrs_new*a_ndata_rows;
370b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
371b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
372b4fbaa2aSMark F. Adams #endif
37322063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
37411e60469SMark F. Adams     ierr = VecCreate( wcomm, &dest_crd );
375d3d6bff4SMark F. Adams     ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr);
376470e26f8SMark F. Adams     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
37711e60469SMark F. Adams     /*
378d3d6bff4SMark F. Adams      There are 'a_ndata_rows*a_ndata_cols' data items per node, (one can think of the vectors of having
379d3d6bff4SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'a_cbs'.
38011e60469SMark F. Adams      */
38192a756f0SMark F. Adams     ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
38211e60469SMark F. Adams     ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr);
383038e3b61SMark F. Adams     for(ii=0,jj=0; ii<ncrs0 ; ii++) {
384d3d6bff4SMark F. Adams       PetscInt id = idx[ii*a_cbs]/a_cbs; /* get node back */
385d3d6bff4SMark F. Adams       for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk;
38611e60469SMark F. Adams     }
387038e3b61SMark F. Adams     ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr);
388d3d6bff4SMark F. Adams     ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat );
38911e60469SMark F. Adams     CHKERRQ(ierr);
39092a756f0SMark F. Adams     ierr = PetscFree( tidx );  CHKERRQ(ierr);
39111e60469SMark F. Adams     /*
39211e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
39311e60469SMark F. Adams      */
394d3d6bff4SMark F. Adams     ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr);
395d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
396d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs0 ; ii++) {
397d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
398d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*a_ndata_cols + jj;
399676e1743SMark F. Adams           PetscScalar tt = (PetscScalar)data[ix];
400676e1743SMark F. Adams 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
401d3d6bff4SMark F. Adams 	}
402038e3b61SMark F. Adams       }
403eb07cef2SMark F. Adams     }
404eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
405eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
40611e60469SMark F. Adams     /*
40711e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
40811e60469SMark F. Adams       to the correct processor
40911e60469SMark F. Adams     */
41011e60469SMark F. Adams     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
41111e60469SMark F. Adams     CHKERRQ(ierr);
41211e60469SMark F. Adams     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
41311e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
41411e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
41511e60469SMark F. Adams     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
41611e60469SMark F. Adams     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
41711e60469SMark F. Adams     /*
41811e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
41911e60469SMark F. Adams     */
420eb07cef2SMark F. Adams     ierr = PetscFree( *a_coarse_data );    CHKERRQ(ierr);
421d3d6bff4SMark F. Adams     ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data );    CHKERRQ(ierr);
422eb07cef2SMark F. Adams 
42311e60469SMark F. Adams     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
424eb07cef2SMark F. Adams     data = *a_coarse_data;
425d3d6bff4SMark F. Adams     for( jj=0; jj<a_ndata_cols ; jj++ ) {
426d3d6bff4SMark F. Adams       for( ii=0 ; ii<ncrs_new ; ii++) {
427d3d6bff4SMark F. Adams 	for( kk=0; kk<a_ndata_rows ; kk++ ) {
428d3d6bff4SMark F. Adams 	  PetscInt ix = ii*a_ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*a_ndata_cols + jj;
429d3d6bff4SMark F. Adams 	  data[ix] = PetscRealPart(array[jx]);
430d3d6bff4SMark F. Adams 	  array[jx] = 1.e300;
431d3d6bff4SMark F. Adams 	}
432038e3b61SMark F. Adams       }
433038e3b61SMark F. Adams     }
43411e60469SMark F. Adams     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
43511e60469SMark F. Adams     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
43611e60469SMark F. Adams     /*
43711e60469SMark F. Adams       Invert for MatGetSubMatrix
43811e60469SMark F. Adams     */
439038e3b61SMark F. Adams     ierr = ISInvertPermutation( isnum, ncrs_new*a_cbs, &new_indices ); CHKERRQ(ierr);
44011e60469SMark F. Adams     ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */
44111e60469SMark F. Adams     ierr = ISDestroy( &isnum ); CHKERRQ(ierr);
44211e60469SMark F. Adams     /* A_crs output */
443038e3b61SMark F. Adams     ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs );
44411e60469SMark F. Adams     CHKERRQ(ierr);
445eb07cef2SMark F. Adams 
446038e3b61SMark F. Adams     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
447e33ef3b1SMark F. Adams     Cmat = *a_Amat_crs; /* output */
448038e3b61SMark F. Adams     ierr = MatSetBlockSize( Cmat, a_cbs );      CHKERRQ(ierr);
449eb07cef2SMark F. Adams 
45011e60469SMark F. Adams     /* prolongator */
45111e60469SMark F. Adams     ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
45211e60469SMark F. Adams     {
45311e60469SMark F. Adams       IS findices;
45411e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
45596434bc1SMark F. Adams #ifdef USE_R
45696434bc1SMark F. Adams       ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew );
45796434bc1SMark F. Adams #else
45811e60469SMark F. Adams       ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew );
45996434bc1SMark F. Adams #endif
46011e60469SMark F. Adams       CHKERRQ(ierr);
461d61a3a7dSMark F. Adams 
46211e60469SMark F. Adams       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
46311e60469SMark F. Adams     }
464d61a3a7dSMark F. Adams 
4653530afc2SMark F. Adams     ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
466a147abb0SMark F. Adams     *a_P_inout = Pnew; /* output - repartitioned */
46792a756f0SMark F. Adams 
46811e60469SMark F. Adams     ierr = ISDestroy( &new_indices ); CHKERRQ(ierr);
46992a756f0SMark F. Adams     ierr = PetscFree( counts );  CHKERRQ(ierr);
47092a756f0SMark F. Adams     ierr = PetscFree( ranks );  CHKERRQ(ierr);
471e33ef3b1SMark F. Adams   }
4725b89ad90SMark F. Adams 
4735b89ad90SMark F. Adams   PetscFunctionReturn(0);
4745b89ad90SMark F. Adams }
4755b89ad90SMark F. Adams 
4765b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4775b89ad90SMark F. Adams /*
4785b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4795b89ad90SMark F. Adams                     by setting data structures and options.
4805b89ad90SMark F. Adams 
4815b89ad90SMark F. Adams    Input Parameter:
4825b89ad90SMark F. Adams .  pc - the preconditioner context
4835b89ad90SMark F. Adams 
4845b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4855b89ad90SMark F. Adams 
4865b89ad90SMark F. Adams    Notes:
4875b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4885b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4895b89ad90SMark F. Adams */
4905b89ad90SMark F. Adams #undef __FUNCT__
4915b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
492eb07cef2SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC a_pc )
4935b89ad90SMark F. Adams {
4945b89ad90SMark F. Adams   PetscErrorCode  ierr;
495eb07cef2SMark F. Adams   PC_MG           *mg = (PC_MG*)a_pc->data;
4965b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
49758471d46SMark F. Adams   PC_MG_Levels   **mglevels = mg->levels;
498eb07cef2SMark F. Adams   Mat              Amat = a_pc->mat, Pmat = a_pc->pmat;
499d3d6bff4SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, nloc, lidx, Istart, Iend;
500eb07cef2SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)a_pc)->comm;
5013530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
502038e3b61SMark F. Adams   PetscBool        isOK;
503587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS];
504587fa25dSMark F. Adams   PetscReal       *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS];
505737a81a9SMark F. Adams   MatInfo          info;
5065ef31b24SMark F. Adams 
5075b89ad90SMark F. Adams   PetscFunctionBegin;
50803a628feSMark F. Adams   pc_gamg->m_count++;
509fca9ed99SMark F. Adams 
51058471d46SMark F. Adams   if( a_pc->setupcalled > 0 ) {
51103a628feSMark F. Adams     /* just do Galerkin grids */
51258471d46SMark F. Adams     Mat B,dA,dB;
51358471d46SMark F. Adams 
51458471d46SMark F. Adams     /* PCSetUp_MG seems to insists on setting this to GMRES */
51558471d46SMark F. Adams     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
51658471d46SMark F. Adams 
51758471d46SMark F. Adams     /* currently only handle case where mat and pmat are the same on coarser levels */
51858471d46SMark F. Adams     ierr = KSPGetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
51958471d46SMark F. Adams     /* (re)set to get dirty flag */
52058471d46SMark F. Adams     ierr = KSPSetOperators(mglevels[pc_gamg->m_Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
52158471d46SMark F. Adams     ierr = KSPSetUp( mglevels[pc_gamg->m_Nlevels-1]->smoothd ); CHKERRQ(ierr);
52258471d46SMark F. Adams 
52358471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-2; level>-1; level--) {
52458471d46SMark F. Adams       ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
52503a628feSMark F. Adams       /* the first time through the matrix structure has changed from repartitioning */
52603a628feSMark F. Adams       if( pc_gamg->m_count == 2 ) {
52703a628feSMark F. Adams         ierr = MatDestroy( &B );  CHKERRQ(ierr);
52803a628feSMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
52903a628feSMark F. Adams         mglevels[level]->A = B;
53003a628feSMark F. Adams       }
53103a628feSMark F. Adams       else {
53258471d46SMark F. Adams         ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
53303a628feSMark F. Adams       }
53458471d46SMark F. Adams       ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
53558471d46SMark F. Adams       dB = B;
53658471d46SMark F. Adams       /* setup KSP/PC */
53758471d46SMark F. Adams       ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr);
53858471d46SMark F. Adams     }
53958471d46SMark F. Adams 
540f6536408SMark F. Adams #define PRINT_MATS PETSC_FALSE
54158471d46SMark F. Adams     /* plot levels - A */
54258471d46SMark F. Adams     if( PRINT_MATS ) {
54358471d46SMark F. Adams       for (lidx=0, level=pc_gamg->m_Nlevels-1; level>0 ; level--,lidx++){
54458471d46SMark F. Adams         PetscViewer viewer;
54558471d46SMark F. Adams         char fname[32]; KSP smoother; Mat Tmat, TTm;
54658471d46SMark F. Adams         ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
54758471d46SMark F. Adams         ierr = KSPGetOperators( smoother, &Tmat, &TTm, 0 ); CHKERRQ(ierr);
54803a628feSMark F. Adams         sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
54958471d46SMark F. Adams         ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
55058471d46SMark F. Adams         ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
55158471d46SMark F. Adams         ierr = MatView( Tmat, viewer ); CHKERRQ(ierr);
55258471d46SMark F. Adams         ierr = PetscViewerDestroy( &viewer );
55358471d46SMark F. Adams       }
55458471d46SMark F. Adams     }
55558471d46SMark F. Adams 
55603a628feSMark F. Adams     a_pc->setupcalled = 2;
55703a628feSMark F. Adams 
55858471d46SMark F. Adams     PetscFunctionReturn(0);
559eb07cef2SMark F. Adams   }
560f6536408SMark F. Adams 
561baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
562baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
5635b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
5645b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
5653530afc2SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
566eb07cef2SMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
567eb07cef2SMark F. Adams   nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0);
568eb07cef2SMark F. Adams 
569038e3b61SMark F. Adams   /* get data of not around */
570038e3b61SMark F. Adams   if( pc_gamg->m_data == 0 && nloc > 0 ) {
571038e3b61SMark F. Adams     ierr  = PCSetCoordinates_GAMG( a_pc, -1, 0 );    CHKERRQ( ierr );
572eb07cef2SMark F. Adams   }
573eb07cef2SMark F. Adams   data = pc_gamg->m_data;
574038e3b61SMark F. Adams 
575eb07cef2SMark F. Adams   /* Get A_i and R_i */
576737a81a9SMark F. Adams   ierr = MatGetInfo(Amat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
577bed7c2b9SMark 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",
5782c047e2dSMark F. Adams 	      mype,__FUNCT__,0,N,pc_gamg->m_data_rows,pc_gamg->m_data_cols,
5792c047e2dSMark F. Adams 	      (int)(info.nz_used/(PetscReal)N),npe);
5808f4b7eb5SMark F. Adams   for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */
5814ef23d27SMark F. Adams         level < (pc_gamg->m_Nlevels-1) && (level==0 || M>pc_gamg->m_coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
5820205a208SMark F. Adams         level++ ){
5835b89ad90SMark F. Adams     level1 = level + 1;
584b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
585b4fbaa2aSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
586b4fbaa2aSMark F. Adams #endif
587b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
588b4fbaa2aSMark F. Adams     ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
589b4fbaa2aSMark F. Adams #endif
590389730f3SMark F. Adams     ierr = createProlongation(Aarr[level], data, level, &bs, &Parr[level1], &coarse_data, &isOK, &emaxs[level], pc_gamg );
5915b89ad90SMark F. Adams     CHKERRQ(ierr);
592d3d6bff4SMark F. Adams     ierr = PetscFree( data ); CHKERRQ( ierr );
593b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
594b4fbaa2aSMark F. Adams     ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
5955c7a89a6SBarry Smith 
596b4fbaa2aSMark F. Adams #endif
597baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
598baa4e9faSMark F. Adams     if( isOK ) {
599b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
600b4fbaa2aSMark F. Adams       ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
601b4fbaa2aSMark F. Adams #endif
6028263b398SMark F. Adams       ierr = createLevel( a_pc, Aarr[level], (pc_gamg->m_method != 0) ? bs : 1, pc_gamg->m_data_cols, bs,
6038263b398SMark F. Adams                           (PetscBool)(level == pc_gamg->m_Nlevels-2),
604389730f3SMark F. Adams                           &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] );
6053530afc2SMark F. Adams       CHKERRQ(ierr);
606b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
607b4fbaa2aSMark F. Adams       ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
608b4fbaa2aSMark F. Adams #endif
6093530afc2SMark F. Adams       ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr);
610737a81a9SMark F. Adams       ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
611bed7c2b9SMark 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",
6122c047e2dSMark F. Adams 		  mype,__FUNCT__,(int)level1,N,pc_gamg->m_data_cols,
6132c047e2dSMark F. Adams 		  (int)(info.nz_used/(PetscReal)N),nactivepe);
614e33ef3b1SMark F. Adams       /* coarse grids with SA can have zero row/cols from singleton aggregates */
61558471d46SMark F. Adams       /* aggregation method should gaurrentee this does not happen! */
616737a81a9SMark F. Adams 
617bed7c2b9SMark F. Adams       if (pc_gamg->m_verbose) {
618785cba28SMark F. Adams         Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id;
619785cba28SMark F. Adams         v = 1.e-10; /* LU factor has hard wired numbers for small diags so this needs to match (yuk) */
620e33ef3b1SMark F. Adams         ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr);
621785cba28SMark F. Adams         nloceq = Iend-Istart;
622e33ef3b1SMark F. Adams         ierr = MatGetVecs( Aarr[level1], &diag, 0 );    CHKERRQ(ierr);
623e33ef3b1SMark F. Adams         ierr = MatGetDiagonal( Aarr[level1], diag );    CHKERRQ(ierr);
624e33ef3b1SMark F. Adams         ierr = VecGetArray( diag, &data_arr );   CHKERRQ(ierr);
625785cba28SMark F. Adams         for(kk=0;kk<nloceq;kk++){
626e33ef3b1SMark F. Adams           if(data_arr[kk]==0.0) {
627785cba28SMark F. Adams             id = kk + Istart;
628e33ef3b1SMark F. Adams             ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES);
629e33ef3b1SMark F. Adams             CHKERRQ(ierr);
630ecd57171SMark F. Adams             PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1);
631e33ef3b1SMark F. Adams           }
632e33ef3b1SMark F. Adams         }
633e33ef3b1SMark F. Adams         ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr);
634e33ef3b1SMark F. Adams         ierr = VecDestroy( &diag );                CHKERRQ(ierr);
635e33ef3b1SMark F. Adams         ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
636e33ef3b1SMark F. Adams         ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
637e33ef3b1SMark F. Adams       }
638486a8d0bSJed Brown     } else {
639be544d3cSMark F. Adams       coarse_data = 0;
640baa4e9faSMark F. Adams       break;
641baa4e9faSMark F. Adams     }
642eb07cef2SMark F. Adams     data = coarse_data;
643b4fbaa2aSMark F. Adams 
644b4fbaa2aSMark F. Adams #if (defined PETSC_USE_LOG && defined GAMG_STAGES)
645b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop(); CHKERRQ( ierr );
646b4fbaa2aSMark F. Adams #endif
6475b89ad90SMark F. Adams   }
648be544d3cSMark F. Adams   if( coarse_data ) {
649eb07cef2SMark F. Adams     ierr = PetscFree( coarse_data ); CHKERRQ( ierr );
650be544d3cSMark F. Adams   }
651bed7c2b9SMark F. Adams   if (pc_gamg->m_verbose) PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
6525b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
6535b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
6545b89ad90SMark F. Adams   fine_level = level;
655eb07cef2SMark F. Adams   ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
6565b89ad90SMark F. Adams 
657fc4362bfSMark F. Adams   /* set default smoothers */
658587fa25dSMark F. Adams   for ( lidx=1, level = pc_gamg->m_Nlevels-2;
659587fa25dSMark F. Adams         lidx <= fine_level;
660587fa25dSMark F. Adams         lidx++, level--) {
6615745e0f5SMark F. Adams     PetscReal emax, emin;
6625b89ad90SMark F. Adams     KSP smoother; PC subpc;
663f6536408SMark F. Adams     PetscBool isCheb;
664f6536408SMark F. Adams     /* set defaults */
665587fa25dSMark F. Adams     ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr);
6665b89ad90SMark F. Adams     ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr);
667f6536408SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
668f6536408SMark F. Adams     /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */
669f6536408SMark F. Adams     ierr = PCSetType( subpc, PETSC_GAMG_SMOOTHER ); CHKERRQ(ierr);
670f6536408SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
671f6536408SMark F. Adams     /* overide defaults with input parameters */
672f6536408SMark F. Adams     ierr = KSPSetFromOptions( smoother ); CHKERRQ(ierr);
673f6536408SMark F. Adams 
674f6536408SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
675f6536408SMark F. Adams     /* do my own cheby */
676f6536408SMark F. Adams     ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr);
67741139db5SMark F. Adams 
678f6536408SMark F. Adams     if( isCheb ) {
679f6536408SMark F. Adams       ierr = PetscTypeCompare( (PetscObject)subpc, PETSC_GAMG_SMOOTHER, &isCheb ); CHKERRQ(ierr);
680f6536408SMark F. Adams       if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
681587fa25dSMark F. Adams       else{ /* eigen estimate 'emax' */
682587fa25dSMark F. Adams         KSP eksp; Mat Lmat = Aarr[level];
683fc4362bfSMark F. Adams         Vec bb, xx; PC pc;
684f6536408SMark F. Adams         const PCType type;
685038e3b61SMark F. Adams 
686f6536408SMark F. Adams         ierr = PCGetType( subpc, &type );   CHKERRQ(ierr);
6875745e0f5SMark F. Adams         ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
6885745e0f5SMark F. Adams         ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
689fc4362bfSMark F. Adams         {
690fc4362bfSMark F. Adams           PetscRandom    rctx;
691fc4362bfSMark F. Adams           ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
692fc4362bfSMark F. Adams           ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
693fc4362bfSMark F. Adams           ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
694fc4362bfSMark F. Adams           ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
6955b89ad90SMark F. Adams         }
696fc4362bfSMark F. Adams         ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
697e33ef3b1SMark F. Adams         ierr = KSPSetType( eksp, KSPCG );                      CHKERRQ(ierr);
698038e3b61SMark F. Adams         ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
699fc4362bfSMark F. Adams         CHKERRQ(ierr);
700fc4362bfSMark F. Adams         ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
701fc4362bfSMark F. Adams 
702f6536408SMark F. Adams         ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
703f6536408SMark F. Adams         ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
704f6536408SMark F. Adams 
705f6536408SMark F. Adams         ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
706f6536408SMark F. Adams         ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
707f6536408SMark F. Adams         ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr );
708f6536408SMark F. Adams         ierr = PCSetType( pc, type ); CHKERRQ(ierr); /* should be same as eigen estimates op. */
709f6536408SMark F. Adams 
710fc4362bfSMark F. Adams         ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
711fc4362bfSMark F. Adams         ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
712fc4362bfSMark F. Adams         ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
713fc4362bfSMark F. Adams         ierr = VecDestroy( &xx );       CHKERRQ(ierr);
714fc4362bfSMark F. Adams         ierr = VecDestroy( &bb );       CHKERRQ(ierr);
715fc4362bfSMark F. Adams         ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
716f6536408SMark F. Adams 
717f6536408SMark F. Adams         if (pc_gamg->m_verbose) {
71841139db5SMark F. Adams           PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s PC setup max eigen=%e min=%e\n",
71941139db5SMark F. Adams                       __FUNCT__,emax,emin);
720f6536408SMark F. Adams         }
721fc4362bfSMark F. Adams       }
722038e3b61SMark F. Adams       {
723038e3b61SMark F. Adams         PetscInt N1, N0, tt;
724038e3b61SMark F. Adams         ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
725038e3b61SMark F. Adams         ierr = MatGetSize( Aarr[level+1], &N0, &tt );       CHKERRQ(ierr);
726f6536408SMark F. Adams         /* heuristic - is this crap? */
727f6536408SMark F. Adams         emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
728038e3b61SMark F. Adams         emax *= 1.05;
729038e3b61SMark F. Adams       }
730fc4362bfSMark F. Adams       ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
731f6536408SMark F. Adams     }
7325745e0f5SMark F. Adams   }
7335745e0f5SMark F. Adams   {
734ecd57171SMark F. Adams     /* coarse grid */
735ecd57171SMark F. Adams     KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
7365745e0f5SMark F. Adams     Mat Lmat = Aarr[pc_gamg->m_Nlevels-1];
737eb07cef2SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
73858471d46SMark F. Adams     ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
7395745e0f5SMark F. Adams     ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
740ecd57171SMark F. Adams     ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
741ecd57171SMark F. Adams     ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
742ecd57171SMark F. Adams     ierr = PCSetUp( subpc ); CHKERRQ(ierr);
743ecd57171SMark F. Adams     ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
744ecd57171SMark F. Adams     assert(ii==1);
745ecd57171SMark F. Adams     ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
746ecd57171SMark F. Adams     ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
747fc4362bfSMark F. Adams   }
748737a81a9SMark F. Adams 
749fc4362bfSMark F. Adams   /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
750eb07cef2SMark F. Adams   ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr);
7515b89ad90SMark F. Adams   {
7525b89ad90SMark F. Adams     PetscBool galerkin;
753eb07cef2SMark F. Adams     ierr = PCMGGetGalerkin( a_pc,  &galerkin); CHKERRQ(ierr);
7545b89ad90SMark F. Adams     if(galerkin){
7555b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
7565b89ad90SMark F. Adams     }
7575b89ad90SMark F. Adams   }
7585745e0f5SMark F. Adams 
75958471d46SMark F. Adams   /* plot levels - R/P */
76058471d46SMark F. Adams   if( PRINT_MATS ) {
76158471d46SMark F. Adams     for (level=pc_gamg->m_Nlevels-1;level>0;level--){
76258471d46SMark F. Adams       PetscViewer viewer;
76358471d46SMark F. Adams       char fname[32];
76403a628feSMark F. Adams       sprintf(fname,"Pmat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
76511e60469SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
7665b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
767038e3b61SMark F. Adams       ierr = MatView( Parr[level], viewer ); CHKERRQ(ierr);
7685b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
76903a628feSMark F. Adams       sprintf(fname,"Amat_%d_%d.m",(int)pc_gamg->m_count,(int)level);
770e33ef3b1SMark F. Adams       ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer );  CHKERRQ(ierr);
771e33ef3b1SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
772e33ef3b1SMark F. Adams       ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr);
773e33ef3b1SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
7745b89ad90SMark F. Adams     }
77558471d46SMark F. Adams   }
77658471d46SMark F. Adams 
77758471d46SMark F. Adams   /* set interpolation between the levels, clean up */
77858471d46SMark F. Adams   for (lidx=0,level=pc_gamg->m_Nlevels-1;
77958471d46SMark F. Adams        lidx<fine_level;
78058471d46SMark F. Adams        lidx++, level--){
78158471d46SMark F. Adams     ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr);
782587fa25dSMark F. Adams     ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
783587fa25dSMark F. Adams     ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
7845b89ad90SMark F. Adams   }
7855745e0f5SMark F. Adams 
7865b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
787eb07cef2SMark F. Adams   a_pc->setupcalled = 0;
788eb07cef2SMark F. Adams   ierr = PCSetUp_MG( a_pc );CHKERRQ( ierr );
78903a628feSMark F. Adams   a_pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */
79058471d46SMark F. Adams 
79158471d46SMark F. Adams   {
79258471d46SMark F. Adams     KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
79358471d46SMark F. Adams     ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr);
79458471d46SMark F. Adams     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
79558471d46SMark F. Adams     ierr = KSPSetUp( smoother ); CHKERRQ(ierr);
79658471d46SMark F. Adams   }
797e33ef3b1SMark F. Adams 
7985b89ad90SMark F. Adams   PetscFunctionReturn(0);
7995b89ad90SMark F. Adams }
8005b89ad90SMark F. Adams 
801eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8025b89ad90SMark F. Adams /*
8035b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8045b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8055b89ad90SMark F. Adams 
8065b89ad90SMark F. Adams    Input Parameter:
8075b89ad90SMark F. Adams .  pc - the preconditioner context
8085b89ad90SMark F. Adams 
8095b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8105b89ad90SMark F. Adams */
8115b89ad90SMark F. Adams #undef __FUNCT__
8125b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
8135b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8145b89ad90SMark F. Adams {
8155b89ad90SMark F. Adams   PetscErrorCode  ierr;
8165b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
8175b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
8185b89ad90SMark F. Adams 
8195b89ad90SMark F. Adams   PetscFunctionBegin;
8205b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8215b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8225b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8235b89ad90SMark F. Adams   PetscFunctionReturn(0);
8245b89ad90SMark F. Adams }
8255b89ad90SMark F. Adams 
826676e1743SMark F. Adams 
827676e1743SMark F. Adams #undef __FUNCT__
828676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
829676e1743SMark F. Adams /*@
830676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
831676e1743SMark F. Adams    processor reduction.
832676e1743SMark F. Adams 
833676e1743SMark F. Adams    Not Collective on PC
834676e1743SMark F. Adams 
835676e1743SMark F. Adams    Input Parameters:
836676e1743SMark F. Adams .  pc - the preconditioner context
837676e1743SMark F. Adams 
838676e1743SMark F. Adams 
839676e1743SMark F. Adams    Options Database Key:
840676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
841676e1743SMark F. Adams 
842676e1743SMark F. Adams    Level: intermediate
843676e1743SMark F. Adams 
844676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
845676e1743SMark F. Adams 
846676e1743SMark F. Adams .seealso: ()
847676e1743SMark F. Adams @*/
848676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
849676e1743SMark F. Adams {
850676e1743SMark F. Adams   PetscErrorCode ierr;
851676e1743SMark F. Adams 
852676e1743SMark F. Adams   PetscFunctionBegin;
853676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
854676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
855676e1743SMark F. Adams   PetscFunctionReturn(0);
856676e1743SMark F. Adams }
857676e1743SMark F. Adams 
858676e1743SMark F. Adams EXTERN_C_BEGIN
859676e1743SMark F. Adams #undef __FUNCT__
860676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
861676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
862676e1743SMark F. Adams {
863c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
864c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
865676e1743SMark F. Adams 
866676e1743SMark F. Adams   PetscFunctionBegin;
867c20e4228SMark F. Adams   if(n>0) pc_gamg->m_min_eq_proc = n;
868676e1743SMark F. Adams   PetscFunctionReturn(0);
869676e1743SMark F. Adams }
870676e1743SMark F. Adams EXTERN_C_END
871676e1743SMark F. Adams 
872676e1743SMark F. Adams #undef __FUNCT__
873389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
874389730f3SMark F. Adams /*@
875389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
876389730f3SMark F. Adams 
877389730f3SMark F. Adams  Collective on PC
878389730f3SMark F. Adams 
879389730f3SMark F. Adams    Input Parameters:
880389730f3SMark F. Adams .  pc - the preconditioner context
881389730f3SMark F. Adams 
882389730f3SMark F. Adams 
883389730f3SMark F. Adams    Options Database Key:
884389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
885389730f3SMark F. Adams 
886389730f3SMark F. Adams    Level: intermediate
887389730f3SMark F. Adams 
888389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
889389730f3SMark F. Adams 
890389730f3SMark F. Adams .seealso: ()
891389730f3SMark F. Adams  @*/
892389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
893389730f3SMark F. Adams {
894389730f3SMark F. Adams   PetscErrorCode ierr;
895389730f3SMark F. Adams 
896389730f3SMark F. Adams   PetscFunctionBegin;
897389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
898389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
899389730f3SMark F. Adams   PetscFunctionReturn(0);
900389730f3SMark F. Adams }
901389730f3SMark F. Adams 
902389730f3SMark F. Adams EXTERN_C_BEGIN
903389730f3SMark F. Adams #undef __FUNCT__
904389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
905389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
906389730f3SMark F. Adams {
907389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
908389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
909389730f3SMark F. Adams 
910389730f3SMark F. Adams   PetscFunctionBegin;
911389730f3SMark F. Adams   if(n>0) pc_gamg->m_coarse_eq_limit = n;
912389730f3SMark F. Adams   PetscFunctionReturn(0);
913389730f3SMark F. Adams }
914389730f3SMark F. Adams EXTERN_C_END
915389730f3SMark F. Adams 
916389730f3SMark F. Adams #undef __FUNCT__
9178263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
918676e1743SMark F. Adams /*@
9198263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
920676e1743SMark F. Adams 
921676e1743SMark F. Adams    Collective on PC
922676e1743SMark F. Adams 
923676e1743SMark F. Adams    Input Parameters:
924676e1743SMark F. Adams .  pc - the preconditioner context
925676e1743SMark F. Adams 
926676e1743SMark F. Adams 
927676e1743SMark F. Adams    Options Database Key:
9288263b398SMark F. Adams .  -pc_gamg_repartition
929676e1743SMark F. Adams 
930676e1743SMark F. Adams    Level: intermediate
931676e1743SMark F. Adams 
932676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
933676e1743SMark F. Adams 
934676e1743SMark F. Adams .seealso: ()
935676e1743SMark F. Adams @*/
9368263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
937676e1743SMark F. Adams {
938676e1743SMark F. Adams   PetscErrorCode ierr;
939676e1743SMark F. Adams 
940676e1743SMark F. Adams   PetscFunctionBegin;
941676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9428263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
943676e1743SMark F. Adams   PetscFunctionReturn(0);
944676e1743SMark F. Adams }
945676e1743SMark F. Adams 
946676e1743SMark F. Adams EXTERN_C_BEGIN
947676e1743SMark F. Adams #undef __FUNCT__
9488263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
9498263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
950676e1743SMark F. Adams {
951c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
952c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
953676e1743SMark F. Adams 
954676e1743SMark F. Adams   PetscFunctionBegin;
9558263b398SMark F. Adams   pc_gamg->m_repart = n;
956676e1743SMark F. Adams   PetscFunctionReturn(0);
957676e1743SMark F. Adams }
958676e1743SMark F. Adams EXTERN_C_END
959676e1743SMark F. Adams 
960676e1743SMark F. Adams #undef __FUNCT__
9614ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
9624ef23d27SMark F. Adams /*@
9634ef23d27SMark F. Adams    PCGAMGSetNlevels -
9644ef23d27SMark F. Adams 
9654ef23d27SMark F. Adams    Not collective on PC
9664ef23d27SMark F. Adams 
9674ef23d27SMark F. Adams    Input Parameters:
9684ef23d27SMark F. Adams .  pc - the preconditioner context
9694ef23d27SMark F. Adams 
9704ef23d27SMark F. Adams 
9714ef23d27SMark F. Adams    Options Database Key:
9724ef23d27SMark F. Adams .  -pc_mg_levels
9734ef23d27SMark F. Adams 
9744ef23d27SMark F. Adams    Level: intermediate
9754ef23d27SMark F. Adams 
9764ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
9774ef23d27SMark F. Adams 
9784ef23d27SMark F. Adams .seealso: ()
9794ef23d27SMark F. Adams @*/
9804ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
9814ef23d27SMark F. Adams {
9824ef23d27SMark F. Adams   PetscErrorCode ierr;
9834ef23d27SMark F. Adams 
9844ef23d27SMark F. Adams   PetscFunctionBegin;
9854ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9864ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
9874ef23d27SMark F. Adams   PetscFunctionReturn(0);
9884ef23d27SMark F. Adams }
9894ef23d27SMark F. Adams 
9904ef23d27SMark F. Adams EXTERN_C_BEGIN
9914ef23d27SMark F. Adams #undef __FUNCT__
9924ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
9934ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
9944ef23d27SMark F. Adams {
9954ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
9964ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
9974ef23d27SMark F. Adams 
9984ef23d27SMark F. Adams   PetscFunctionBegin;
9994ef23d27SMark F. Adams   pc_gamg->m_Nlevels = n;
10004ef23d27SMark F. Adams   PetscFunctionReturn(0);
10014ef23d27SMark F. Adams }
10024ef23d27SMark F. Adams EXTERN_C_END
10034ef23d27SMark F. Adams 
10044ef23d27SMark F. Adams #undef __FUNCT__
10053542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
10063542efc5SMark F. Adams /*@
10073542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
10083542efc5SMark F. Adams 
10093542efc5SMark F. Adams    Not collective on PC
10103542efc5SMark F. Adams 
10113542efc5SMark F. Adams    Input Parameters:
10123542efc5SMark F. Adams .  pc - the preconditioner context
10133542efc5SMark F. Adams 
10143542efc5SMark F. Adams 
10153542efc5SMark F. Adams    Options Database Key:
10163542efc5SMark F. Adams .  -pc_gamg_threshold
10173542efc5SMark F. Adams 
10183542efc5SMark F. Adams    Level: intermediate
10193542efc5SMark F. Adams 
10203542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
10213542efc5SMark F. Adams 
10223542efc5SMark F. Adams .seealso: ()
10233542efc5SMark F. Adams @*/
10243542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
10253542efc5SMark F. Adams {
10263542efc5SMark F. Adams   PetscErrorCode ierr;
10273542efc5SMark F. Adams 
10283542efc5SMark F. Adams   PetscFunctionBegin;
10293542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10303542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
10313542efc5SMark F. Adams   PetscFunctionReturn(0);
10323542efc5SMark F. Adams }
10333542efc5SMark F. Adams 
10343542efc5SMark F. Adams EXTERN_C_BEGIN
10353542efc5SMark F. Adams #undef __FUNCT__
10363542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
10373542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
10383542efc5SMark F. Adams {
1039c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1040c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
10413542efc5SMark F. Adams 
10423542efc5SMark F. Adams   PetscFunctionBegin;
1043c20e4228SMark F. Adams   pc_gamg->m_threshold = n;
10443542efc5SMark F. Adams   PetscFunctionReturn(0);
10453542efc5SMark F. Adams }
10463542efc5SMark F. Adams EXTERN_C_END
10473542efc5SMark F. Adams 
10483542efc5SMark F. Adams #undef __FUNCT__
1049676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType"
1050676e1743SMark F. Adams /*@
1051676e1743SMark F. Adams    PCGAMGSetSolverType - Set solution method.
1052676e1743SMark F. Adams 
1053676e1743SMark F. Adams    Collective on PC
1054676e1743SMark F. Adams 
1055676e1743SMark F. Adams    Input Parameters:
1056676e1743SMark F. Adams .  pc - the preconditioner context
1057676e1743SMark F. Adams 
1058676e1743SMark F. Adams 
1059676e1743SMark F. Adams    Options Database Key:
10603542efc5SMark F. Adams .  -pc_gamg_type
1061676e1743SMark F. Adams 
1062676e1743SMark F. Adams    Level: intermediate
1063676e1743SMark F. Adams 
1064676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1065676e1743SMark F. Adams 
1066676e1743SMark F. Adams .seealso: ()
1067676e1743SMark F. Adams @*/
1068676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType(PC pc, char str[], PetscInt sz )
1069676e1743SMark F. Adams {
1070676e1743SMark F. Adams   PetscErrorCode ierr;
1071676e1743SMark F. Adams 
1072676e1743SMark F. Adams   PetscFunctionBegin;
1073676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1074676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSolverType_C",(PC,char[],PetscInt),(pc,str,sz));
1075676e1743SMark F. Adams   CHKERRQ(ierr);
1076676e1743SMark F. Adams   PetscFunctionReturn(0);
1077676e1743SMark F. Adams }
1078676e1743SMark F. Adams 
1079676e1743SMark F. Adams EXTERN_C_BEGIN
1080676e1743SMark F. Adams #undef __FUNCT__
1081676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetSolverType_GAMG"
1082676e1743SMark F. Adams PetscErrorCode PCGAMGSetSolverType_GAMG(PC pc, char str[], PetscInt sz )
1083676e1743SMark F. Adams {
1084676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1085676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1086676e1743SMark F. Adams 
1087676e1743SMark F. Adams   PetscFunctionBegin;
1088676e1743SMark F. Adams   if(sz < 64) strcpy(pc_gamg->m_type,str);
1089676e1743SMark F. Adams   PetscFunctionReturn(0);
1090676e1743SMark F. Adams }
1091676e1743SMark F. Adams EXTERN_C_END
1092676e1743SMark F. Adams 
10935b89ad90SMark F. Adams #undef __FUNCT__
10945b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
10955b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
10965b89ad90SMark F. Adams {
1097676e1743SMark F. Adams   PetscErrorCode  ierr;
1098676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1099676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1100676e1743SMark F. Adams   PetscBool flag;
11015b89ad90SMark F. Adams 
11025b89ad90SMark F. Adams   PetscFunctionBegin;
1103676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1104676e1743SMark F. Adams   {
1105676e1743SMark F. Adams     ierr = PetscOptionsString("-pc_gamg_type",
1106470e26f8SMark F. Adams                               "Solver type: plane aggregation ('pa'), smoothed aggregation ('sa') or geometric multigrid (default)",
1107676e1743SMark F. Adams                               "PCGAMGSetSolverType",
1108676e1743SMark F. Adams                               pc_gamg->m_type,
1109676e1743SMark F. Adams                               pc_gamg->m_type,
1110676e1743SMark F. Adams                               64,
1111676e1743SMark F. Adams                               &flag );
1112676e1743SMark F. Adams     CHKERRQ(ierr);
1113*d8c859f0SMark F. Adams     if( flag && pc_gamg->m_data != 0 ) {
1114*d8c859f0SMark F. Adams       if( (strcmp(pc_gamg->m_type,"sa")==0 && pc_gamg->m_method != 2) ||
1115*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"pa")==0 && pc_gamg->m_method != 1) ||
1116*d8c859f0SMark F. Adams           (strcmp(pc_gamg->m_type,"geo")==0 && pc_gamg->m_method != 0) ) {
1117*d8c859f0SMark F. Adams         SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "PCSetFromOptions called of PCSetCoordinates (with new method, after data was created)");
1118*d8c859f0SMark F. Adams       }
1119*d8c859f0SMark F. Adams     }
1120bed7c2b9SMark F. Adams 
1121bed7c2b9SMark F. Adams     /* -pc_gamg_verbose */
1122bed7c2b9SMark 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);
11232a44bfbaSMark F. Adams 
1124f6536408SMark F. Adams     pc_gamg->m_method = 1; /* default to plane aggregation */
1125f6536408SMark F. Adams     if (flag ) {
1126f6536408SMark F. Adams       if( strcmp(pc_gamg->m_type,"sa") == 0) pc_gamg->m_method = 2;
1127f6536408SMark F. Adams       else if( strcmp(pc_gamg->m_type,"pa") == 0) pc_gamg->m_method = 1;
1128f6536408SMark F. Adams       else if( strcmp(pc_gamg->m_type,"geo") == 0) pc_gamg->m_method = 0;
1129f6536408SMark F. Adams       else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Invalid gamg type: %s",pc_gamg->m_type);
1130f6536408SMark F. Adams     }
11318263b398SMark F. Adams     /* -pc_gamg_repartition */
11328263b398SMark F. Adams     pc_gamg->m_repart = PETSC_FALSE;
11338263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
11348263b398SMark F. Adams                             "Repartion coarse grids (false)",
11358263b398SMark F. Adams                             "PCGAMGRepartitioning",
11368263b398SMark F. Adams                             pc_gamg->m_repart,
11378263b398SMark F. Adams                             &pc_gamg->m_repart,
1138676e1743SMark F. Adams                             &flag);
1139676e1743SMark F. Adams     CHKERRQ(ierr);
1140676e1743SMark F. Adams 
1141c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1142c20e4228SMark F. Adams     pc_gamg->m_min_eq_proc = 600;
1143676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1144676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1145676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
1146c20e4228SMark F. Adams                            pc_gamg->m_min_eq_proc,
1147c20e4228SMark F. Adams                            &pc_gamg->m_min_eq_proc,
1148676e1743SMark F. Adams                            &flag );
1149676e1743SMark F. Adams     CHKERRQ(ierr);
11503542efc5SMark F. Adams 
1151389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1152389730f3SMark F. Adams     pc_gamg->m_coarse_eq_limit = 1500;
1153389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1154389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1155389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
1156389730f3SMark F. Adams                            pc_gamg->m_coarse_eq_limit,
1157389730f3SMark F. Adams                            &pc_gamg->m_coarse_eq_limit,
1158389730f3SMark F. Adams                            &flag );
1159389730f3SMark F. Adams     CHKERRQ(ierr);
1160389730f3SMark F. Adams 
1161c20e4228SMark F. Adams     /* -pc_gamg_threshold */
1162c20e4228SMark F. Adams     pc_gamg->m_threshold = 0.05;
11633542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
11643542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
11653542efc5SMark F. Adams                             "PCGAMGSetThreshold",
1166c20e4228SMark F. Adams                             pc_gamg->m_threshold,
1167c20e4228SMark F. Adams                             &pc_gamg->m_threshold,
11683542efc5SMark F. Adams                             &flag );
11693542efc5SMark F. Adams     CHKERRQ(ierr);
11704ef23d27SMark F. Adams 
11714ef23d27SMark F. Adams     pc_gamg->m_Nlevels = GAMG_MAXLEVELS;
11724ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
11734ef23d27SMark F. Adams                            "Set number of MG levels",
11744ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
11754ef23d27SMark F. Adams                            pc_gamg->m_Nlevels,
11764ef23d27SMark F. Adams                            &pc_gamg->m_Nlevels,
11774ef23d27SMark F. Adams                            &flag );
1178676e1743SMark F. Adams   }
1179676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1180676e1743SMark F. Adams 
11815b89ad90SMark F. Adams   PetscFunctionReturn(0);
11825b89ad90SMark F. Adams }
11835b89ad90SMark F. Adams 
11845b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
11855b89ad90SMark F. Adams /*MC
1186dc76db92SMark F. Adams      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning. This preconditioner currently has two
1187dc76db92SMark F. Adams            AMG methods: 1) an unstructured geometric method, which requires that you provide coordinates for each
1188dc76db92SMark F. Adams            vertex, and 2) smoothed aggregation.  Smoothed aggregation (SA) can work without coordinates but it
1189dc76db92SMark F. Adams            will generate some common non-trivial null spaces if coordinates are provided.  The input fine grid matrix
1190dc76db92SMark F. Adams            must have the block size set for 'system' problems (with multiple dof per vertex/cell) to work properly.
1191dc76db92SMark F. Adams            SA will generate rotational rigid body mode null space vectors, in addition to the trivial translational
1192dc76db92SMark F. Adams            modes, when coordinates are provide in 2D and 3D.
11935b89ad90SMark F. Adams 
1194280d9858SJed Brown    Options Database Keys:
11955b89ad90SMark F. Adams    Multigrid options(inherited)
1196280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1197280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1198280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1199280d9858SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
12005b89ad90SMark F. Adams 
12015b89ad90SMark F. Adams   Level: intermediate
1202280d9858SJed Brown 
12035b89ad90SMark F. Adams   Concepts: multigrid
12045b89ad90SMark F. Adams 
12055b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1206280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
12075b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
12085b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12095b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12105b89ad90SMark F. Adams M*/
12115b89ad90SMark F. Adams 
12125b89ad90SMark F. Adams EXTERN_C_BEGIN
12135b89ad90SMark F. Adams #undef __FUNCT__
12145b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
12155b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
12165b89ad90SMark F. Adams {
12175b89ad90SMark F. Adams   PetscErrorCode  ierr;
12185b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
12195b89ad90SMark F. Adams   PC_MG           *mg;
12205ef31b24SMark F. Adams   PetscClassId     cookie;
12215ee9c036SSatish Balay #if defined PETSC_USE_LOG
12225ee9c036SSatish Balay   static int count = 0;
12235ee9c036SSatish Balay #endif
12245b89ad90SMark F. Adams 
12255b89ad90SMark F. Adams   PetscFunctionBegin;
12265b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
12275b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
12285b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
12295b89ad90SMark F. Adams 
12305b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
12315b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
123203a628feSMark F. Adams   pc_gamg->m_data_sz = 0; pc_gamg->m_data = 0; pc_gamg->m_count = 0;
12335b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
12345b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
12355b89ad90SMark F. Adams 
12365b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
12375b89ad90SMark F. Adams 
12385b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
12395b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
12405b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
12415b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
12425b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
12435b89ad90SMark F. Adams 
12445b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12455b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
12465b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
1247c97e1df0SMark F. Adams 					    PCSetCoordinates_GAMG);
1248c97e1df0SMark F. Adams   CHKERRQ(ierr);
1249c97e1df0SMark F. Adams 
1250676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1251676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1252676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1253676e1743SMark F. Adams 					    PCGAMGSetProcEqLim_GAMG);
1254676e1743SMark F. Adams   CHKERRQ(ierr);
1255676e1743SMark F. Adams 
1256676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1257389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_C",
1258389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_GAMG",
1259389730f3SMark F. Adams 					    PCGAMGSetCoarseEqLim_GAMG);
1260389730f3SMark F. Adams   CHKERRQ(ierr);
1261389730f3SMark F. Adams 
1262389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12638263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_C",
12648263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_GAMG",
12658263b398SMark F. Adams 					    PCGAMGSetRepartitioning_GAMG);
1266676e1743SMark F. Adams   CHKERRQ(ierr);
1267676e1743SMark F. Adams 
1268676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
12693542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
12703542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
12713542efc5SMark F. Adams 					    PCGAMGSetThreshold_GAMG);
12723542efc5SMark F. Adams   CHKERRQ(ierr);
12733542efc5SMark F. Adams 
12743542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1275676e1743SMark F. Adams 					    "PCGAMGSetSolverType_C",
1276676e1743SMark F. Adams 					    "PCGAMGSetSolverType_GAMG",
1277676e1743SMark F. Adams 					    PCGAMGSetSolverType_GAMG);
1278676e1743SMark F. Adams   CHKERRQ(ierr);
1279c97e1df0SMark F. Adams 
1280b4fbaa2aSMark F. Adams #if defined PETSC_USE_LOG
1281785cba28SMark F. Adams   if( count++ == 0 ) {
12825ef31b24SMark F. Adams     PetscClassIdRegister("GAMG Setup",&cookie);
1283b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]);
12842a44bfbaSMark F. Adams     PetscLogEventRegister("  Graph", cookie, &gamg_setup_events[GRAPH]);
12852a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]);
12862a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]);
12872a44bfbaSMark F. Adams     PetscLogEventRegister("    G.Square", cookie, &gamg_setup_events[GRAPH_SQR]);
1288b4fbaa2aSMark F. Adams     PetscLogEventRegister("  MIS/Agg", cookie, &gamg_setup_events[SET4]);
1289b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: growSupp", cookie, &gamg_setup_events[SET5]);
1290b4fbaa2aSMark F. Adams     PetscLogEventRegister("  geo: triangle", cookie, &gamg_setup_events[SET6]);
1291b4fbaa2aSMark F. Adams     PetscLogEventRegister("    search&set", cookie, &gamg_setup_events[FIND_V]);
1292b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: init", cookie, &gamg_setup_events[SET7]);
1293b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("  SA: frmProl0", cookie, &gamg_setup_events[SET8]); */
1294b4fbaa2aSMark F. Adams     PetscLogEventRegister("  SA: smooth", cookie, &gamg_setup_events[SET9]);
1295b4fbaa2aSMark F. Adams     PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]);
1296b4fbaa2aSMark F. Adams     PetscLogEventRegister("  PL repartition", cookie, &gamg_setup_events[SET12]);
1297f852f58cSMark F. Adams 
1298fca9ed99SMark F. Adams     /* PetscLogEventRegister("  PL 1", cookie, &gamg_setup_events[SET13]); */
1299fca9ed99SMark F. Adams     /* PetscLogEventRegister("  PL 2", cookie, &gamg_setup_events[SET14]); */
1300fca9ed99SMark F. Adams     /* PetscLogEventRegister("  PL 3", cookie, &gamg_setup_events[SET15]); */
1301f852f58cSMark F. Adams 
1302b4fbaa2aSMark F. Adams     /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */
1303b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */
1304b4fbaa2aSMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */
13055ef31b24SMark F. Adams 
1306b4fbaa2aSMark F. Adams     /* create timer stages */
1307b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1308b4fbaa2aSMark F. Adams     {
1309b4fbaa2aSMark F. Adams       char str[32];
1310b4fbaa2aSMark F. Adams       sprintf(str,"MG Level %d (finest)",0);
1311b4fbaa2aSMark F. Adams       PetscLogStageRegister(str, &gamg_stages[0]);
1312b4fbaa2aSMark F. Adams       PetscInt lidx;
1313b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1314b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1315b4fbaa2aSMark F. Adams 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1316b4fbaa2aSMark F. Adams       }
1317b4fbaa2aSMark F. Adams     }
1318b4fbaa2aSMark F. Adams #endif
1319b4fbaa2aSMark F. Adams   }
1320b4fbaa2aSMark F. Adams #endif
13215b89ad90SMark F. Adams   PetscFunctionReturn(0);
13225b89ad90SMark F. Adams }
13235b89ad90SMark F. Adams EXTERN_C_END
1324