xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 1b18a24a426d773f19f5ebdc15fb0694bc41a23f)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4b45d2f2cSJed Brown #include "petsc-private/matimpl.h"
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
6b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
75b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */
8f96513f1SMatthew G Knepley 
90cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
100cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
11b4fbaa2aSMark F. Adams #endif
120cbbd2e1SMark F. Adams 
130cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
140cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_AGG;
150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_GEO;
160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
200cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGOptprol_AGG;
210cbbd2e1SMark F. Adams #endif
220cbbd2e1SMark F. Adams 
23b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
24b4fbaa2aSMark F. Adams 
25b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
260cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
27b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
28b4fbaa2aSMark F. Adams #endif
29f96513f1SMatthew G Knepley 
30140e18c1SBarry Smith static PetscFunctionList GAMGList = 0;
313e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
329d5b6da9SMark F. Adams 
33d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
34d3d6bff4SMark F. Adams #undef __FUNCT__
35d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
36d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
37d3d6bff4SMark F. Adams {
38d3d6bff4SMark F. Adams   PetscErrorCode ierr;
39d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
40d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
41d3d6bff4SMark F. Adams 
42d3d6bff4SMark F. Adams   PetscFunctionBegin;
43a2f3521dSMark F. Adams   if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */
44ce94432eSBarry Smith     PetscPrintf(PetscObjectComm((PetscObject)pc),"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__);
459d5b6da9SMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
4658471d46SMark F. Adams   }
470298fd71SBarry Smith   pc_gamg->data = NULL; pc_gamg->data_sz = 0;
48878e152fSMark F. Adams 
49878e152fSMark F. Adams   if (pc_gamg->orig_data) {
50878e152fSMark F. Adams     ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
51878e152fSMark F. Adams   }
52a2f3521dSMark F. Adams   PetscFunctionReturn(0);
53a2f3521dSMark F. Adams }
54a2f3521dSMark F. Adams 
555b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
565b89ad90SMark F. Adams /*
57a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
58a147abb0SMark F. Adams      of active processors.
595b89ad90SMark F. Adams 
605b89ad90SMark F. Adams    Input Parameter:
61a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
62a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
639d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
64c5bfad50SMark F. Adams    . cr_bs - coarse block size
659d5b6da9SMark F. Adams    . isLast -
663530afc2SMark F. Adams    In/Output Parameter:
67a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
68afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
6911e60469SMark F. Adams    Output Parameter:
703530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
715b89ad90SMark F. Adams */
725cb416c2SMark F. Adams 
735b89ad90SMark F. Adams #undef __FUNCT__
748263b398SMark F. Adams #define __FUNCT__ "createLevel"
757700e67bSMark Adams static PetscErrorCode createLevel(const PC pc,const Mat Amat_fine,const PetscInt cr_bs,const PetscBool isLast,
767700e67bSMark Adams                                   Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc)
775b89ad90SMark F. Adams {
78a2f3521dSMark F. Adams   PetscErrorCode  ierr;
799d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
80486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
819d5b6da9SMark F. Adams   const PetscBool repart      = pc_gamg->repart;
829d5b6da9SMark F. Adams   const PetscInt  min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
83a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
843b4367a7SBarry Smith   MPI_Comm        comm;
85c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
86c5bfad50SMark F. Adams   PetscInt        ncrs_eq,ncrs_prim,f_bs;
875b89ad90SMark F. Adams 
885b89ad90SMark F. Adams   PetscFunctionBegin;
893b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
903b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
913b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
92c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
9311e60469SMark F. Adams   /* RAP */
949d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
95038e3b61SMark F. Adams 
96a2f3521dSMark F. Adams   /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/
97a2f3521dSMark F. Adams   ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
9871959b99SBarry Smith   if (pc_gamg->data_sz % (pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data_sz %D not divisible by (pc_gamg->data_cell_cols %D *pc_gamg->data_cell_rows %D)",pc_gamg->data_sz,pc_gamg->data_cell_cols,pc_gamg->data_cell_rows);
990298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
100a2f3521dSMark F. Adams 
101c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
102a2f3521dSMark F. Adams   {
103472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
1040298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
105c5df96a5SBarry Smith     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
106c5df96a5SBarry Smith     if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1;
107c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
108c5df96a5SBarry Smith     if (isLast) new_size = 1;
109a2f3521dSMark F. Adams   }
110f852f58cSMark F. Adams 
1112fa5cd67SKarl Rupp   if (!repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1122fa5cd67SKarl Rupp   else {
113a2f3521dSMark F. Adams     const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
114a2f3521dSMark F. Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
1157700e67bSMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
1165a9b9e01SMark F. Adams     VecScatter     vecscat;
11722063be5SMark F. Adams     PetscScalar    *array;
11822063be5SMark F. Adams     Vec            src_crd, dest_crd;
119e33ef3b1SMark F. Adams 
12071959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
12171959b99SBarry Smith     if (ncrs_eq % cr_bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs);
1220cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1230cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
124b4fbaa2aSMark F. Adams #endif
125a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
126c5df96a5SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt), &counts);CHKERRQ(ierr);
1277700e67bSMark Adams     if (repart) {
128a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1295a9b9e01SMark F. Adams       Mat adj;
1305a9b9e01SMark F. Adams 
131a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
1323b4367a7SBarry Smith         if (pc_gamg->verbose==1) PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,ncrs_eq);
133a2f3521dSMark F. Adams         else {
134a2f3521dSMark F. Adams           PetscInt n;
1353b4367a7SBarry Smith           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
1363b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n);
137a2f3521dSMark F. Adams         }
138a2f3521dSMark F. Adams       }
1395a9b9e01SMark F. Adams 
140a2f3521dSMark F. Adams       /* get 'adj' */
141c5bfad50SMark F. Adams       if (cr_bs == 1) {
142038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
143806fa848SBarry Smith       } else {
144a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
145eb07cef2SMark F. Adams         Mat               tMat;
146a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
147b4fbaa2aSMark F. Adams         const PetscScalar *vals;
148b4fbaa2aSMark F. Adams         const PetscInt    *idx;
149a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1509057884aSMark F. Adams         static PetscInt   llev = 0;
151b4fbaa2aSMark F. Adams 
152a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &d_nnz);CHKERRQ(ierr);
153a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &o_nnz);CHKERRQ(ierr);
154a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
155a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
156c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
15758471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
158c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
159c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
16058471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
161a2f3521dSMark F. Adams           if (d_nnz[jj] > ncrs_prim) d_nnz[jj] = ncrs_prim;
162c5bfad50SMark F. Adams           if (o_nnz[jj] > (M/cr_bs-ncrs_prim)) o_nnz[jj] = M/cr_bs-ncrs_prim;
16358471d46SMark F. Adams         }
1646876a03eSMark F. Adams 
1653b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
166806fa848SBarry Smith         ierr = MatSetSizes(tMat, ncrs_prim, ncrs_prim,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
167a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
168a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
169a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
17058471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
1715f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
172eb07cef2SMark F. Adams 
173a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
174c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
17522063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
176eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
177c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
178eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
179eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
180eb07cef2SMark F. Adams           }
18122063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
182eb07cef2SMark F. Adams         }
183eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
184eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
185eb07cef2SMark F. Adams 
186b4fbaa2aSMark F. Adams         if (llev++ == -1) {
187b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
1888caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
1893b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
190b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
1913bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
192b4fbaa2aSMark F. Adams         }
193b4fbaa2aSMark F. Adams 
194eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
195eb07cef2SMark F. Adams 
196eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
197a2f3521dSMark F. Adams       } /* create 'adj' */
198f150b916SMark F. Adams 
199a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2005a9b9e01SMark F. Adams         char            prefix[256];
2015a9b9e01SMark F. Adams         const char      *pcpre;
202b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
203b4fbaa2aSMark F. Adams         MatPartitioning mpart;
204a4b7d37bSMark F. Adams         IS              proc_is;
205a2f3521dSMark F. Adams         PetscInt        targetPE;
2062f03bc48SMark F. Adams 
2073b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2085ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2099d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2108caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
21159a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
21211e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
213c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
214a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
21511e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2165a9b9e01SMark F. Adams 
2175ef31b24SMark F. Adams         /* collect IS info */
218a2f3521dSMark F. Adams         ierr     = PetscMalloc(ncrs_eq*sizeof(PetscInt), &newproc_idx);CHKERRQ(ierr);
219a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
220a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
221c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
222a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
223c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
224a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
225eb07cef2SMark F. Adams           }
2265ef31b24SMark F. Adams         }
227a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
228a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2295ef31b24SMark F. Adams       }
2305ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2315a9b9e01SMark F. Adams 
2323b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2338263b398SMark F. Adams       if (newproc_idx != 0) {
2348263b398SMark F. Adams         ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2355ef31b24SMark F. Adams       }
236806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
237a2f3521dSMark F. Adams 
238a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
2395a9b9e01SMark F. Adams       /* find factor */
240c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2415a9b9e01SMark F. Adams       else {
2425a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2435a9b9e01SMark F. Adams         jj = -1;
244c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
245c5df96a5SBarry Smith           if (size%kk==0) { /* a candidate */
246c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2475a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2485a9b9e01SMark F. Adams             if (fact > best_fact) {
2495a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
2505a9b9e01SMark F. Adams             }
2515a9b9e01SMark F. Adams           }
2525a9b9e01SMark F. Adams         }
2535a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
254a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
2555a9b9e01SMark F. Adams       }
256c5df96a5SBarry Smith       new_size = size/rfactor;
2575a9b9e01SMark F. Adams 
258c5df96a5SBarry Smith       if (new_size==nactive) {
259a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2605a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
261a2f3521dSMark F. Adams         if (pc_gamg->verbose>0) {
2623b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq);
263a2f3521dSMark F. Adams         }
2645a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2655a9b9e01SMark F. Adams       }
2665a9b9e01SMark F. Adams 
2673b4367a7SBarry Smith       if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq);
268c5df96a5SBarry Smith       targetPE = rank/rfactor;
2693b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
270a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
271e33ef3b1SMark F. Adams 
27211e60469SMark F. Adams     /*
273a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
27411e60469SMark F. Adams      */
275a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2767700e67bSMark Adams     is_eq_num_prim = is_eq_num;
27711e60469SMark F. Adams     /*
278a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
27911e60469SMark F. Adams      */
280c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
281c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
282a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
2837700e67bSMark Adams     ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */
284a2f3521dSMark F. Adams 
285a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2860cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2870cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
288b4fbaa2aSMark F. Adams #endif
289a2f3521dSMark F. Adams 
290a2f3521dSMark F. Adams     /* move data (for primal equations only) */
29122063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
2923b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
293a2f3521dSMark F. Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE);CHKERRQ(ierr);
294c0dedaeaSBarry Smith     ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
29511e60469SMark F. Adams     /*
2969d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
297c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
29811e60469SMark F. Adams      */
299a2f3521dSMark F. Adams     ierr = PetscMalloc((ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx);CHKERRQ(ierr);
300a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
301a2f3521dSMark F. Adams     for (ii=0,jj=0; ii<ncrs_prim; ii++) {
302c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
303a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
30411e60469SMark F. Adams     }
305a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3063b4367a7SBarry Smith     ierr = ISCreateGeneral(comm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
30792a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
30811e60469SMark F. Adams     /*
30911e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
31011e60469SMark F. Adams      */
311a2f3521dSMark F. Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd);CHKERRQ(ierr);
3129d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
313a2f3521dSMark F. Adams       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
314a2f3521dSMark F. Adams       for (ii=0; ii<ncrs_prim; ii++) {
3159d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
316a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
317c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
318676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
319d3d6bff4SMark F. Adams         }
320038e3b61SMark F. Adams       }
321eb07cef2SMark F. Adams     }
322eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
323eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
32411e60469SMark F. Adams     /*
32511e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
32611e60469SMark F. Adams       to the correct processor
32711e60469SMark F. Adams     */
3280298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
32911e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
33011e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
33111e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
33211e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
33311e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
33411e60469SMark F. Adams     /*
33511e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
33611e60469SMark F. Adams     */
337c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
338a2f3521dSMark F. Adams     ierr = PetscMalloc(node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
3392fa5cd67SKarl Rupp 
340a2f3521dSMark F. Adams     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
341a2f3521dSMark F. Adams     strideNew        = ncrs_prim_new*ndata_rows;
3422fa5cd67SKarl Rupp 
34311e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3449d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
345a2f3521dSMark F. Adams       for (ii=0; ii<ncrs_prim_new; ii++) {
3469d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
347a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
348c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
349d3d6bff4SMark F. Adams         }
350038e3b61SMark F. Adams       }
351038e3b61SMark F. Adams     }
35211e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
35311e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
354a2f3521dSMark F. Adams 
355a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3560cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3570cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
358ed3f9983SMark F. Adams #endif
359a2f3521dSMark F. Adams 
36011e60469SMark F. Adams     /*
36111e60469SMark F. Adams       Invert for MatGetSubMatrix
36211e60469SMark F. Adams     */
363a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
364a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
365c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
366a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
367a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
368a2f3521dSMark F. Adams     }
369a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3700cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3710cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3720cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
373ed3f9983SMark F. Adams #endif
374a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
375a2f3521dSMark F. Adams     {
376a2f3521dSMark F. Adams       Mat mat;
377806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
378a2f3521dSMark F. Adams       *a_Amat_crs = mat;
379c5bfad50SMark F. Adams 
380c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
381c5bfad50SMark F. Adams         PetscInt cbs, rbs;
382c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
383c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
384c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
385c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",rank,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr);
386c5bfad50SMark F. Adams       }
387a2f3521dSMark F. Adams     }
388038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
389a2f3521dSMark F. Adams 
3900cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3910cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
392ed3f9983SMark F. Adams #endif
39311e60469SMark F. Adams     /* prolongator */
39411e60469SMark F. Adams     {
39511e60469SMark F. Adams       IS       findices;
396a2f3521dSMark F. Adams       PetscInt Istart,Iend;
397a2f3521dSMark F. Adams       Mat      Pnew;
398a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
3990cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4000cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
401ed3f9983SMark F. Adams #endif
4023b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
403c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
404806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
40511e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
406c5bfad50SMark F. Adams 
407c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
408c5bfad50SMark F. Adams         PetscInt cbs, rbs;
409c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
410c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
411c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
412c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
413c5bfad50SMark F. Adams       }
4140cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4150cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
416ed3f9983SMark F. Adams #endif
4173530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
418a2f3521dSMark F. Adams 
419a2f3521dSMark F. Adams       /* output - repartitioned */
420a2f3521dSMark F. Adams       *a_P_inout = Pnew;
421e33ef3b1SMark F. Adams     }
422a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4235b89ad90SMark F. Adams 
424c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
425a2f3521dSMark F. Adams   }
4265a9b9e01SMark F. Adams 
427a2f3521dSMark F. Adams   /* outout matrix data */
428c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
429c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
430c8b0795cSMark F. Adams     if (llev==0) {
431c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
4323b4367a7SBarry Smith       PetscViewerASCIIOpen(comm,fname,&viewer);
433c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
434c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
435c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
436c8b0795cSMark F. Adams     }
437c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
4383b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
439c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
440c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
441c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
442c8b0795cSMark F. Adams   }
4435b89ad90SMark F. Adams   PetscFunctionReturn(0);
4445b89ad90SMark F. Adams }
4455b89ad90SMark F. Adams 
4465b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4475b89ad90SMark F. Adams /*
4485b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4495b89ad90SMark F. Adams                     by setting data structures and options.
4505b89ad90SMark F. Adams 
4515b89ad90SMark F. Adams    Input Parameter:
4525b89ad90SMark F. Adams .  pc - the preconditioner context
4535b89ad90SMark F. Adams 
4545b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4555b89ad90SMark F. Adams 
4565b89ad90SMark F. Adams    Notes:
4575b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4585b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4595b89ad90SMark F. Adams */
4605b89ad90SMark F. Adams #undef __FUNCT__
4615b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4629d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4635b89ad90SMark F. Adams {
4645b89ad90SMark F. Adams   PetscErrorCode ierr;
4659d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4665b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4672adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
468a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
4693b4367a7SBarry Smith   MPI_Comm       comm;
470c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
471587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
472c8b0795cSMark F. Adams   PetscReal      emaxs[GAMG_MAXLEVELS];
473e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
474a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
475737a81a9SMark F. Adams   MatInfo        info;
4767700e67bSMark Adams   PetscBool      redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);
4775ef31b24SMark F. Adams 
4785b89ad90SMark F. Adams   PetscFunctionBegin;
4793b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4803b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4813b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
482dfd5c07aSMark F. Adams 
4833b4367a7SBarry Smith   if (pc_gamg->verbose>2) PetscPrintf(comm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",rank,__FUNCT__,pc_gamg->setup_count,pc->setupcalled);
484dfd5c07aSMark F. Adams 
48584d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
486878e152fSMark F. Adams     if (redo_mesh_setup) {
487878e152fSMark F. Adams       /* reset everything */
488878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
489878e152fSMark F. Adams       pc->setupcalled = 0;
490806fa848SBarry Smith     } else {
49184d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
49203a628feSMark F. Adams       /* just do Galerkin grids */
49358471d46SMark F. Adams       Mat          B,dA,dB;
49458471d46SMark F. Adams 
49571959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
4969d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
49758471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
4980298fd71SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,NULL);CHKERRQ(ierr);
49958471d46SMark F. Adams         /* (re)set to get dirty flag */
5009d5b6da9SMark F. Adams         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
50158471d46SMark F. Adams 
5022fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
50303a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
5040a97e771SToby Isaac           if (pc_gamg->setup_count==2) {
50503a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
506084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
5072fa5cd67SKarl Rupp 
50803a628feSMark F. Adams             mglevels[level]->A = B;
509806fa848SBarry Smith           } else {
5100298fd71SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B,NULL);CHKERRQ(ierr);
51158471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
51203a628feSMark F. Adams           }
51358471d46SMark F. Adams           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
51458471d46SMark F. Adams           dB   = B;
51558471d46SMark F. Adams         }
5165f8cf99dSMark F. Adams       }
517d5280255SMark F. Adams 
518d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
519d5280255SMark F. Adams 
520d5280255SMark F. Adams       /* PCSetUp_MG seems to insists on setting this to GMRES */
5211ac9965eSMark F. Adams       ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr);
52258471d46SMark F. Adams       PetscFunctionReturn(0);
523eb07cef2SMark F. Adams     }
524878e152fSMark F. Adams   }
525f6536408SMark F. Adams 
526878e152fSMark F. Adams   if (!pc_gamg->data) {
527878e152fSMark F. Adams     if (pc_gamg->orig_data) {
528878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5290298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5302fa5cd67SKarl Rupp 
531878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
532878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
533878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5342fa5cd67SKarl Rupp 
535878e152fSMark F. Adams       ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
536878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
537806fa848SBarry Smith     } else {
5381ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5397700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5409d5b6da9SMark F. Adams     }
541878e152fSMark F. Adams   }
542878e152fSMark F. Adams 
543878e152fSMark F. Adams   /* cache original data for reuse */
544878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
545878e152fSMark F. Adams     ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr);
546878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
547878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
548878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
549878e152fSMark F. Adams   }
550038e3b61SMark F. Adams 
551302f38e8SMark F. Adams   /* get basic dims */
552302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
553a2f3521dSMark F. Adams 
554a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
555c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
55684f9421dSMark F. Adams     PetscInt NN = M;
55784f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
55884f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
5593bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
560806fa848SBarry Smith     } else {
561806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
56284f9421dSMark F. Adams     }
563b2a4f308SMark F. Adams     nnz0   = info.nz_used;
564b2a4f308SMark F. Adams     nnztot = info.nz_used;
5653b4367a7SBarry Smith     ierr   = PetscPrintf(comm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
566c5df96a5SBarry Smith                          rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
567c5df96a5SBarry Smith                          (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr);
568c8b0795cSMark F. Adams   }
56984d3f75bSMark F. Adams 
570a2f3521dSMark F. Adams   /* Get A_i and R_i */
571c5df96a5SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
572c5df96a5SBarry Smith        level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit);  /* && (size==1 || nactivepe>1); */
5730205a208SMark F. Adams        level++) {
5745b89ad90SMark F. Adams     level1 = level + 1;
5750cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5760cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
577a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
578a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
579b4fbaa2aSMark F. Adams #endif
580a2f3521dSMark F. Adams #endif
581c8b0795cSMark F. Adams     { /* construct prolongator */
582725b86d8SJed Brown       Mat              Gmat;
5830cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5847700e67bSMark Adams       Mat              Prol11;
585c8b0795cSMark F. Adams 
5867700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5871ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5887700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
589c8b0795cSMark F. Adams 
590a2f3521dSMark F. Adams       /* could have failed to create new level */
591a2f3521dSMark F. Adams       if (Prol11) {
5929d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5930298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
594a2f3521dSMark F. Adams 
5951ab5ffc9SJed Brown         if (pc_gamg->ops->optprol) {
596c8b0795cSMark F. Adams           /* smooth */
5977700e67bSMark Adams           ierr = pc_gamg->ops->optprol(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
598c8b0795cSMark F. Adams         }
599c8b0795cSMark F. Adams 
6007700e67bSMark Adams         Parr[level1] = Prol11;
6010298fd71SBarry Smith       } else Parr[level1] = NULL;
602ffc955d6SMark F. Adams 
603ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
604*1b18a24aSMark Adams         PetscInt bs;
605*1b18a24aSMark Adams         ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
606806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
607ffc955d6SMark F. Adams       }
608ffc955d6SMark F. Adams 
609a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
61041b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
611a2f3521dSMark F. Adams     } /* construct prolongator scope */
6120cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6130cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
614c8b0795cSMark F. Adams #endif
6159d5b6da9SMark F. Adams     /* cache eigen estimate */
6169d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
6179d5b6da9SMark F. Adams       PetscBool flag;
6187700e67bSMark Adams       ierr = PetscObjectComposedDataGetReal((PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
6199d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
620806fa848SBarry Smith     } else emaxs[level] = -1.;
6212adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
622c8b0795cSMark F. Adams     if (!Parr[level1]) {
623806fa848SBarry Smith       if (pc_gamg->verbose) {
6243b4367a7SBarry Smith         ierr =  PetscPrintf(comm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr);
625806fa848SBarry Smith       }
626c8b0795cSMark F. Adams       break;
627c8b0795cSMark F. Adams     }
6280cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6290cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
630b4fbaa2aSMark F. Adams #endif
631a2f3521dSMark F. Adams 
632a2f3521dSMark F. Adams     ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
6337700e67bSMark Adams                        &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr);
634a2f3521dSMark F. Adams 
6350cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6360cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
637b4fbaa2aSMark F. Adams #endif
638a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
639a2f3521dSMark F. Adams 
640a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0) {
6410cbbd2e1SMark F. Adams       PetscInt NN = M;
6420cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
643a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
6443bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
645806fa848SBarry Smith       } else {
646806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
6470cbbd2e1SMark F. Adams       }
648a2f3521dSMark F. Adams 
6490cbbd2e1SMark F. Adams       nnztot += info.nz_used;
6503b4367a7SBarry Smith       ierr    = PetscPrintf(comm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
651c5df96a5SBarry Smith                             rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
652806fa848SBarry Smith                             (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
653c8b0795cSMark F. Adams     }
654a2f3521dSMark F. Adams 
655a2f3521dSMark F. Adams     /* stop if one node -- could pull back for singular problems */
656c8b0795cSMark F. Adams     if (M/pc_gamg->data_cell_cols < 2) {
65781708dccSMark F. Adams       level++;
65881708dccSMark F. Adams       break;
65981708dccSMark F. Adams     }
6600cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
661b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
662b4fbaa2aSMark F. Adams #endif
663c8b0795cSMark F. Adams   } /* levels */
664c8b0795cSMark F. Adams 
665c8b0795cSMark F. Adams   if (pc_gamg->data) {
666c8b0795cSMark F. Adams     ierr          = PetscFree(pc_gamg->data);CHKERRQ(ierr);
6670298fd71SBarry Smith     pc_gamg->data = NULL;
6685b89ad90SMark F. Adams   }
669c8b0795cSMark F. Adams 
6703b4367a7SBarry Smith   if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
6719d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6725b89ad90SMark F. Adams   fine_level       = level;
6730298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
6745b89ad90SMark F. Adams 
67584d3f75bSMark F. Adams   /* simple setup */
67684d3f75bSMark F. Adams   if (!PETSC_TRUE) {
67784d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
67884d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
67984d3f75bSMark F. Adams          lidx<fine_level;
68084d3f75bSMark F. Adams          lidx++, level--) {
68184d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
68284d3f75bSMark F. Adams       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
68384d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
68484d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
68584d3f75bSMark F. Adams     }
68684d3f75bSMark F. Adams     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
68784d3f75bSMark F. Adams 
68884d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
689806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
690d5280255SMark F. Adams     /* set default smoothers & set operators */
6919d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
692587fa25dSMark F. Adams          lidx <= fine_level;
693587fa25dSMark F. Adams          lidx++, level--) {
694ffc955d6SMark F. Adams       KSP smoother;
695ffc955d6SMark F. Adams       PC  subpc;
696a2f3521dSMark F. Adams 
6979d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
698f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
699ffc955d6SMark F. Adams 
700a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
701a2f3521dSMark F. Adams       /* set ops */
702a2f3521dSMark F. Adams       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
703a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
704a2f3521dSMark F. Adams 
705a2f3521dSMark F. Adams       /* set defaults */
7066c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
707a2f3521dSMark F. Adams 
708*1b18a24aSMark Adams       /* set blocks for GASM smoother that uses the 'aggregates' */
709ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
7102d3561bbSSatish Balay         PetscInt sz;
7112d3561bbSSatish Balay         IS       *is;
712a2f3521dSMark F. Adams 
7132d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
7142d3561bbSSatish Balay         is   = ASMLocalIDsArr[level];
715ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
716*1b18a24aSMark Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
717ffc955d6SMark F. Adams         if (sz==0) {
718ffc955d6SMark F. Adams           IS       is;
719ffc955d6SMark F. Adams           PetscInt my0,kk;
720ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
721ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
7220298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr);
723a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
724806fa848SBarry Smith         } else {
725a94c3b12SMark F. Adams           PetscInt kk;
7260298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr);
727a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
728a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
729a94c3b12SMark F. Adams           }
730ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
731ffc955d6SMark F. Adams         }
7320298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
733ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
734ffc955d6SMark F. Adams         ierr                  = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
735806fa848SBarry Smith       } else {
736890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
737ffc955d6SMark F. Adams       }
738d5280255SMark F. Adams     }
739d5280255SMark F. Adams     {
740d5280255SMark F. Adams       /* coarse grid */
741d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
742d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
743d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
744d5280255SMark F. Adams       ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
745d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
746d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
747d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
748d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
74971959b99SBarry Smith       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
75071959b99SBarry Smith       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
751d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
752d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
7539dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
7542fb0b348SMark F. Adams       ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
7555b42dca8SJed Brown       /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
7565b42dca8SJed Brown        * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
7575b42dca8SJed Brown        * view every subdomain as though they were different. */
7585b42dca8SJed Brown       ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
759d5280255SMark F. Adams     }
760d5280255SMark F. Adams 
761d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
762d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
763d5280255SMark F. Adams     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr);
764d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
7653b4367a7SBarry Smith     if (mg->galerkin != 2) SETERRQ(comm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
766d5280255SMark F. Adams 
767d5280255SMark F. Adams     /* create cheby smoothers */
768d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
769d5280255SMark F. Adams          lidx <= fine_level;
770d5280255SMark F. Adams          lidx++, level--) {
771d5280255SMark F. Adams       KSP       smoother;
772890ffe84SMark Adams       PetscBool flag,flag2;
773d5280255SMark F. Adams       PC        subpc;
774d5280255SMark F. Adams 
775ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
776a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
777a2f3521dSMark F. Adams 
778ffc955d6SMark F. Adams       /* do my own cheby */
7796c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
780ffc955d6SMark F. Adams       if (flag) {
781ffc955d6SMark F. Adams         PetscReal emax, emin;
782251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
783890ffe84SMark Adams         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCSOR, &flag2);CHKERRQ(ierr);
784890ffe84SMark Adams         if ((flag||flag2) && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC but lets acccept SOR because it is close and safe (always lower) */
785890ffe84SMark Adams         else { /* eigen estimate 'emax' -- this is done in cheby */
786e696ed0bSMark F. Adams           KSP eksp;
787e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
788b2a4f308SMark F. Adams           Vec bb, xx;
789038e3b61SMark F. Adams 
7905745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr);
7915745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr);
792fc4362bfSMark F. Adams           {
793fc4362bfSMark F. Adams             PetscRandom rctx;
7943b4367a7SBarry Smith             ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
795fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
796fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
797fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
7985b89ad90SMark F. Adams           }
799a2f3521dSMark F. Adams 
800e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
801e696ed0bSMark F. Adams           {
802e696ed0bSMark F. Adams             PetscInt    Istart,Iend,ncols,jj,Ii;
803e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
804e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
805e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
806e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
807e696ed0bSMark F. Adams               if (ncols <= 1) {
808e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
809a94c3b12SMark F. Adams               }
810e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
811a94c3b12SMark F. Adams             }
812a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
813a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
814a94c3b12SMark F. Adams           }
815e696ed0bSMark F. Adams 
8163b4367a7SBarry Smith           ierr = KSPCreate(comm, &eksp);CHKERRQ(ierr);
817806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
818fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
8191a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
8201a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
821f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
822f6536408SMark F. Adams 
823f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
824f6536408SMark F. Adams           ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
825fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
8265a9b9e01SMark F. Adams 
827d3d0db20SJed Brown           /* set PC type to be same as smoother */
828ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
829b2a4f308SMark F. Adams 
8305a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
8315a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
8325a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
833fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
8345a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
8355a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
8365a9b9e01SMark F. Adams 
837fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
8385a9b9e01SMark F. Adams 
839fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
840fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
841fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
842f6536408SMark F. Adams 
843ffc955d6SMark F. Adams           if (pc_gamg->verbose > 0) {
844a94c3b12SMark F. Adams             PetscInt N1, tt;
845a94c3b12SMark F. Adams             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
8463b4367a7SBarry Smith             PetscPrintf(comm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1);
847f6536408SMark F. Adams           }
848fc4362bfSMark F. Adams         }
849038e3b61SMark F. Adams         {
850c5bfad50SMark F. Adams           PetscInt N1, N0;
8510298fd71SBarry Smith           ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr);
8520298fd71SBarry Smith           ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr);
853f6536408SMark F. Adams           /* heuristic - is this crap? */
854b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
8555e7c91beSJed Brown           emin  = emax * pc_gamg->eigtarget[0];
8565e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
857038e3b61SMark F. Adams         }
8586c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
859ffc955d6SMark F. Adams       } /* setup checby flag */
860ffc955d6SMark F. Adams     } /* non-coarse levels */
861737a81a9SMark F. Adams 
862d5280255SMark F. Adams     /* clean up */
863d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
864587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
865587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
8665b89ad90SMark F. Adams     }
8670cbbd2e1SMark F. Adams 
8680cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8690cbbd2e1SMark F. Adams 
8701ac9965eSMark F. Adams     if (PETSC_TRUE) {
87158471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
8729d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
87358471d46SMark F. Adams       ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
87458471d46SMark F. Adams     }
875806fa848SBarry Smith   } else {
8765f8cf99dSMark F. Adams     KSP smoother;
8773b4367a7SBarry Smith     if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",rank,__FUNCT__);
8789d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
8795f8cf99dSMark F. Adams     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8805f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
8819d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8825f8cf99dSMark F. Adams   }
8835b89ad90SMark F. Adams   PetscFunctionReturn(0);
8845b89ad90SMark F. Adams }
8855b89ad90SMark F. Adams 
886eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8875b89ad90SMark F. Adams /*
8885b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8895b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8905b89ad90SMark F. Adams 
8915b89ad90SMark F. Adams    Input Parameter:
8925b89ad90SMark F. Adams .  pc - the preconditioner context
8935b89ad90SMark F. Adams 
8945b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8955b89ad90SMark F. Adams */
8965b89ad90SMark F. Adams #undef __FUNCT__
8975b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
8985b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8995b89ad90SMark F. Adams {
9005b89ad90SMark F. Adams   PetscErrorCode ierr;
9015b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
9025b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
9035b89ad90SMark F. Adams 
9045b89ad90SMark F. Adams   PetscFunctionBegin;
9055b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
9069b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
9079b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
9089b8ffb57SJed Brown   }
9091ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
9101ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
9115b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
9125b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
9135b89ad90SMark F. Adams   PetscFunctionReturn(0);
9145b89ad90SMark F. Adams }
9155b89ad90SMark F. Adams 
916676e1743SMark F. Adams 
917676e1743SMark F. Adams #undef __FUNCT__
918676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
919676e1743SMark F. Adams /*@
920676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
921676e1743SMark F. Adams    processor reduction.
922676e1743SMark F. Adams 
923676e1743SMark F. Adams    Not Collective on PC
924676e1743SMark F. Adams 
925676e1743SMark F. Adams    Input Parameters:
926676e1743SMark F. Adams .  pc - the preconditioner context
927676e1743SMark F. Adams 
928676e1743SMark F. Adams 
929676e1743SMark F. Adams    Options Database Key:
930676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
931676e1743SMark F. Adams 
932676e1743SMark F. Adams    Level: intermediate
933676e1743SMark F. Adams 
934676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
935676e1743SMark F. Adams 
936676e1743SMark F. Adams .seealso: ()
937676e1743SMark F. Adams @*/
938676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
939676e1743SMark F. Adams {
940676e1743SMark F. Adams   PetscErrorCode ierr;
941676e1743SMark F. Adams 
942676e1743SMark F. Adams   PetscFunctionBegin;
943676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
944676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
945676e1743SMark F. Adams   PetscFunctionReturn(0);
946676e1743SMark F. Adams }
947676e1743SMark F. Adams 
948676e1743SMark F. Adams #undef __FUNCT__
949676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
9501e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
951676e1743SMark F. Adams {
952c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
953c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
954676e1743SMark F. Adams 
955676e1743SMark F. Adams   PetscFunctionBegin;
9569d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
957676e1743SMark F. Adams   PetscFunctionReturn(0);
958676e1743SMark F. Adams }
959676e1743SMark F. Adams 
960676e1743SMark F. Adams #undef __FUNCT__
961389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
962389730f3SMark F. Adams /*@
963389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
964389730f3SMark F. Adams 
965389730f3SMark F. Adams  Collective on PC
966389730f3SMark F. Adams 
967389730f3SMark F. Adams    Input Parameters:
968389730f3SMark F. Adams .  pc - the preconditioner context
969389730f3SMark F. Adams 
970389730f3SMark F. Adams 
971389730f3SMark F. Adams    Options Database Key:
972389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
973389730f3SMark F. Adams 
974389730f3SMark F. Adams    Level: intermediate
975389730f3SMark F. Adams 
976389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
977389730f3SMark F. Adams 
978389730f3SMark F. Adams .seealso: ()
979389730f3SMark F. Adams  @*/
980389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
981389730f3SMark F. Adams {
982389730f3SMark F. Adams   PetscErrorCode ierr;
983389730f3SMark F. Adams 
984389730f3SMark F. Adams   PetscFunctionBegin;
985389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
986389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
987389730f3SMark F. Adams   PetscFunctionReturn(0);
988389730f3SMark F. Adams }
989389730f3SMark F. Adams 
990389730f3SMark F. Adams #undef __FUNCT__
991389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
9921e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
993389730f3SMark F. Adams {
994389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
995389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
996389730f3SMark F. Adams 
997389730f3SMark F. Adams   PetscFunctionBegin;
9989d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
999389730f3SMark F. Adams   PetscFunctionReturn(0);
1000389730f3SMark F. Adams }
1001389730f3SMark F. Adams 
1002389730f3SMark F. Adams #undef __FUNCT__
10038263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1004676e1743SMark F. Adams /*@
10058263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1006676e1743SMark F. Adams 
1007676e1743SMark F. Adams    Collective on PC
1008676e1743SMark F. Adams 
1009676e1743SMark F. Adams    Input Parameters:
1010676e1743SMark F. Adams .  pc - the preconditioner context
1011676e1743SMark F. Adams 
1012676e1743SMark F. Adams 
1013676e1743SMark F. Adams    Options Database Key:
10148263b398SMark F. Adams .  -pc_gamg_repartition
1015676e1743SMark F. Adams 
1016676e1743SMark F. Adams    Level: intermediate
1017676e1743SMark F. Adams 
1018676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1019676e1743SMark F. Adams 
1020676e1743SMark F. Adams .seealso: ()
1021676e1743SMark F. Adams @*/
10228263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1023676e1743SMark F. Adams {
1024676e1743SMark F. Adams   PetscErrorCode ierr;
1025676e1743SMark F. Adams 
1026676e1743SMark F. Adams   PetscFunctionBegin;
1027676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10288263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1029676e1743SMark F. Adams   PetscFunctionReturn(0);
1030676e1743SMark F. Adams }
1031676e1743SMark F. Adams 
1032676e1743SMark F. Adams #undef __FUNCT__
10338263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
10341e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1035676e1743SMark F. Adams {
1036c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1037c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1038676e1743SMark F. Adams 
1039676e1743SMark F. Adams   PetscFunctionBegin;
10409d5b6da9SMark F. Adams   pc_gamg->repart = n;
1041676e1743SMark F. Adams   PetscFunctionReturn(0);
1042676e1743SMark F. Adams }
1043676e1743SMark F. Adams 
1044676e1743SMark F. Adams #undef __FUNCT__
1045dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl"
1046dfd5c07aSMark F. Adams /*@
1047dfd5c07aSMark F. Adams    PCGAMGSetReuseProl - Reuse prlongation
1048dfd5c07aSMark F. Adams 
1049dfd5c07aSMark F. Adams    Collective on PC
1050dfd5c07aSMark F. Adams 
1051dfd5c07aSMark F. Adams    Input Parameters:
1052dfd5c07aSMark F. Adams .  pc - the preconditioner context
1053dfd5c07aSMark F. Adams 
1054dfd5c07aSMark F. Adams 
1055dfd5c07aSMark F. Adams    Options Database Key:
1056dfd5c07aSMark F. Adams .  -pc_gamg_reuse_interpolation
1057dfd5c07aSMark F. Adams 
1058dfd5c07aSMark F. Adams    Level: intermediate
1059dfd5c07aSMark F. Adams 
1060dfd5c07aSMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1061dfd5c07aSMark F. Adams 
1062dfd5c07aSMark F. Adams .seealso: ()
1063dfd5c07aSMark F. Adams @*/
1064dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n)
1065dfd5c07aSMark F. Adams {
1066dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1067dfd5c07aSMark F. Adams 
1068dfd5c07aSMark F. Adams   PetscFunctionBegin;
1069dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1070dfd5c07aSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1071dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1072dfd5c07aSMark F. Adams }
1073dfd5c07aSMark F. Adams 
1074dfd5c07aSMark F. Adams #undef __FUNCT__
1075dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl_GAMG"
10761e6b0712SBarry Smith static PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n)
1077dfd5c07aSMark F. Adams {
1078dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1079dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1080dfd5c07aSMark F. Adams 
1081dfd5c07aSMark F. Adams   PetscFunctionBegin;
1082dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1083dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1084dfd5c07aSMark F. Adams }
1085dfd5c07aSMark F. Adams 
1086dfd5c07aSMark F. Adams #undef __FUNCT__
1087ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1088ffc955d6SMark F. Adams /*@
1089ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1090ffc955d6SMark F. Adams 
1091ffc955d6SMark F. Adams    Collective on PC
1092ffc955d6SMark F. Adams 
1093ffc955d6SMark F. Adams    Input Parameters:
1094ffc955d6SMark F. Adams .  pc - the preconditioner context
1095ffc955d6SMark F. Adams 
1096ffc955d6SMark F. Adams 
1097ffc955d6SMark F. Adams    Options Database Key:
1098ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1099ffc955d6SMark F. Adams 
1100ffc955d6SMark F. Adams    Level: intermediate
1101ffc955d6SMark F. Adams 
1102ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1103ffc955d6SMark F. Adams 
1104ffc955d6SMark F. Adams .seealso: ()
1105ffc955d6SMark F. Adams @*/
1106ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1107ffc955d6SMark F. Adams {
1108ffc955d6SMark F. Adams   PetscErrorCode ierr;
1109ffc955d6SMark F. Adams 
1110ffc955d6SMark F. Adams   PetscFunctionBegin;
1111ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1112ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1113ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1114ffc955d6SMark F. Adams }
1115ffc955d6SMark F. Adams 
1116ffc955d6SMark F. Adams #undef __FUNCT__
1117ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
11181e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1119ffc955d6SMark F. Adams {
1120ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1121ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1122ffc955d6SMark F. Adams 
1123ffc955d6SMark F. Adams   PetscFunctionBegin;
1124ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1125ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1126ffc955d6SMark F. Adams }
1127ffc955d6SMark F. Adams 
1128ffc955d6SMark F. Adams #undef __FUNCT__
11294ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
11304ef23d27SMark F. Adams /*@
11314ef23d27SMark F. Adams    PCGAMGSetNlevels -
11324ef23d27SMark F. Adams 
11334ef23d27SMark F. Adams    Not collective on PC
11344ef23d27SMark F. Adams 
11354ef23d27SMark F. Adams    Input Parameters:
11364ef23d27SMark F. Adams .  pc - the preconditioner context
11374ef23d27SMark F. Adams 
11384ef23d27SMark F. Adams 
11394ef23d27SMark F. Adams    Options Database Key:
11404ef23d27SMark F. Adams .  -pc_mg_levels
11414ef23d27SMark F. Adams 
11424ef23d27SMark F. Adams    Level: intermediate
11434ef23d27SMark F. Adams 
11444ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
11454ef23d27SMark F. Adams 
11464ef23d27SMark F. Adams .seealso: ()
11474ef23d27SMark F. Adams @*/
11484ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
11494ef23d27SMark F. Adams {
11504ef23d27SMark F. Adams   PetscErrorCode ierr;
11514ef23d27SMark F. Adams 
11524ef23d27SMark F. Adams   PetscFunctionBegin;
11534ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11544ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
11554ef23d27SMark F. Adams   PetscFunctionReturn(0);
11564ef23d27SMark F. Adams }
11574ef23d27SMark F. Adams 
11584ef23d27SMark F. Adams #undef __FUNCT__
11594ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
11601e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
11614ef23d27SMark F. Adams {
11624ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
11634ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
11644ef23d27SMark F. Adams 
11654ef23d27SMark F. Adams   PetscFunctionBegin;
11669d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
11674ef23d27SMark F. Adams   PetscFunctionReturn(0);
11684ef23d27SMark F. Adams }
11694ef23d27SMark F. Adams 
11704ef23d27SMark F. Adams #undef __FUNCT__
11713542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
11723542efc5SMark F. Adams /*@
11733542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
11743542efc5SMark F. Adams 
11753542efc5SMark F. Adams    Not collective on PC
11763542efc5SMark F. Adams 
11773542efc5SMark F. Adams    Input Parameters:
11783542efc5SMark F. Adams .  pc - the preconditioner context
11793542efc5SMark F. Adams 
11803542efc5SMark F. Adams 
11813542efc5SMark F. Adams    Options Database Key:
11823542efc5SMark F. Adams .  -pc_gamg_threshold
11833542efc5SMark F. Adams 
11843542efc5SMark F. Adams    Level: intermediate
11853542efc5SMark F. Adams 
11863542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
11873542efc5SMark F. Adams 
11883542efc5SMark F. Adams .seealso: ()
11893542efc5SMark F. Adams @*/
11903542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
11913542efc5SMark F. Adams {
11923542efc5SMark F. Adams   PetscErrorCode ierr;
11933542efc5SMark F. Adams 
11943542efc5SMark F. Adams   PetscFunctionBegin;
11953542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11963542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
11973542efc5SMark F. Adams   PetscFunctionReturn(0);
11983542efc5SMark F. Adams }
11993542efc5SMark F. Adams 
12003542efc5SMark F. Adams #undef __FUNCT__
12013542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
12021e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
12033542efc5SMark F. Adams {
1204c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1205c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12063542efc5SMark F. Adams 
12073542efc5SMark F. Adams   PetscFunctionBegin;
12089d5b6da9SMark F. Adams   pc_gamg->threshold = n;
12093542efc5SMark F. Adams   PetscFunctionReturn(0);
12103542efc5SMark F. Adams }
12113542efc5SMark F. Adams 
12123542efc5SMark F. Adams #undef __FUNCT__
12139d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1214676e1743SMark F. Adams /*@
12159d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1216676e1743SMark F. Adams 
1217676e1743SMark F. Adams    Collective on PC
1218676e1743SMark F. Adams 
1219676e1743SMark F. Adams    Input Parameters:
1220676e1743SMark F. Adams .  pc - the preconditioner context
1221676e1743SMark F. Adams 
1222676e1743SMark F. Adams 
1223676e1743SMark F. Adams    Options Database Key:
12243542efc5SMark F. Adams .  -pc_gamg_type
1225676e1743SMark F. Adams 
1226676e1743SMark F. Adams    Level: intermediate
1227676e1743SMark F. Adams 
1228676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1229676e1743SMark F. Adams 
1230676e1743SMark F. Adams .seealso: ()
1231676e1743SMark F. Adams @*/
123219fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1233676e1743SMark F. Adams {
1234676e1743SMark F. Adams   PetscErrorCode ierr;
1235676e1743SMark F. Adams 
1236676e1743SMark F. Adams   PetscFunctionBegin;
1237676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1238806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1239676e1743SMark F. Adams   PetscFunctionReturn(0);
1240676e1743SMark F. Adams }
1241676e1743SMark F. Adams 
1242676e1743SMark F. Adams #undef __FUNCT__
12439d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
12441e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1245676e1743SMark F. Adams {
12469d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
12471ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
12481ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1249676e1743SMark F. Adams 
1250676e1743SMark F. Adams   PetscFunctionBegin;
12511c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
12529d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
12531ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
12541ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
12551ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
12561ab5ffc9SJed Brown   }
12571ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
12581ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
12599d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1260676e1743SMark F. Adams   PetscFunctionReturn(0);
1261676e1743SMark F. Adams }
1262676e1743SMark F. Adams 
12635b89ad90SMark F. Adams #undef __FUNCT__
12645b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
12655b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
12665b89ad90SMark F. Adams {
1267676e1743SMark F. Adams   PetscErrorCode ierr;
1268676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1269676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1270676e1743SMark F. Adams   PetscBool      flag;
12715e7c91beSJed Brown   PetscInt       two   = 2;
12723b4367a7SBarry Smith   MPI_Comm       comm;
12735b89ad90SMark F. Adams 
12745b89ad90SMark F. Adams   PetscFunctionBegin;
12753b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1276676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
1277676e1743SMark F. Adams   {
1278b7cbab4eSMark Adams     /* -pc_gamg_type */
1279b7cbab4eSMark Adams     {
1280bd94a7aaSJed Brown       char tname[256];
1281bd94a7aaSJed Brown       ierr = PetscOptionsList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1282b7cbab4eSMark Adams       /* call PCCreateGAMG_XYZ */
1283bd94a7aaSJed Brown       if (flag) {
1284bd94a7aaSJed Brown         ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
12851ab5ffc9SJed Brown       }
1286b7cbab4eSMark Adams     }
128775b74bdaSMark F. Adams     /* -pc_gamg_verbose */
12889d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
12899d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
12900298fd71SBarry Smith                            &pc_gamg->verbose, NULL);CHKERRQ(ierr);
12918263b398SMark F. Adams     /* -pc_gamg_repartition */
12928263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
12938263b398SMark F. Adams                             "Repartion coarse grids (false)",
12948263b398SMark F. Adams                             "PCGAMGRepartitioning",
12959d5b6da9SMark F. Adams                             pc_gamg->repart,
12969d5b6da9SMark F. Adams                             &pc_gamg->repart,
1297806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1298dfd5c07aSMark F. Adams     /* -pc_gamg_reuse_interpolation */
1299dfd5c07aSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation",
1300dfd5c07aSMark F. Adams                             "Reuse prolongation operator (true)",
1301dfd5c07aSMark F. Adams                             "PCGAMGReuseProl",
1302dfd5c07aSMark F. Adams                             pc_gamg->reuse_prol,
1303dfd5c07aSMark F. Adams                             &pc_gamg->reuse_prol,
1304dfd5c07aSMark F. Adams                             &flag);CHKERRQ(ierr);
1305ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
1306ffc955d6SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1307ffc955d6SMark F. Adams                             "Use aggregation agragates for GASM smoother (false)",
1308ffc955d6SMark F. Adams                             "PCGAMGUseASMAggs",
1309ffc955d6SMark F. Adams                             pc_gamg->use_aggs_in_gasm,
1310ffc955d6SMark F. Adams                             &pc_gamg->use_aggs_in_gasm,
1311806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1312c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1313676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1314676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1315676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
13169d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
13179d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1318806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1319389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1320389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1321389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1322389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
13239d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
13249d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1325806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1326c20e4228SMark F. Adams     /* -pc_gamg_threshold */
13273542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
13283542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
13293542efc5SMark F. Adams                             "PCGAMGSetThreshold",
13309d5b6da9SMark F. Adams                             pc_gamg->threshold,
13319d5b6da9SMark F. Adams                             &pc_gamg->threshold,
1332806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1333806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
13343b4367a7SBarry Smith       ierr = PetscPrintf(comm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1335806fa848SBarry Smith     }
1336b7cbab4eSMark Adams     /* -pc_gamg_eigtarget */
13370298fd71SBarry Smith     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr);
13384ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
13394ef23d27SMark F. Adams                            "Set number of MG levels",
13404ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
13419d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
13429d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
1343806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1344b7cbab4eSMark Adams 
1345b7cbab4eSMark Adams     /* set options for subtype */
13461ab5ffc9SJed Brown     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(pc);CHKERRQ(ierr);}
1347676e1743SMark F. Adams   }
1348676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
13495b89ad90SMark F. Adams   PetscFunctionReturn(0);
13505b89ad90SMark F. Adams }
13515b89ad90SMark F. Adams 
13525b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
13535b89ad90SMark F. Adams /*MC
1354856830bfSJed Brown      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
13559d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
13565b89ad90SMark F. Adams 
1357280d9858SJed Brown    Options Database Keys:
13585b89ad90SMark F. Adams    Multigrid options(inherited)
1359280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1360280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1361280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
13628c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
13635b89ad90SMark F. Adams 
13645b89ad90SMark F. Adams   Level: intermediate
1365280d9858SJed Brown 
13665b89ad90SMark F. Adams   Concepts: multigrid
13675b89ad90SMark F. Adams 
13685b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1369280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
13705b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
13715b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
13725b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
13735b89ad90SMark F. Adams M*/
1374b2573a8aSBarry Smith 
13755b89ad90SMark F. Adams #undef __FUNCT__
13765b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
13778cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
13785b89ad90SMark F. Adams {
13795b89ad90SMark F. Adams   PetscErrorCode ierr;
13805b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
13815b89ad90SMark F. Adams   PC_MG          *mg;
13820cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13839d5b6da9SMark F. Adams   static long count = 0;
13845ee9c036SSatish Balay #endif
13855b89ad90SMark F. Adams 
13865b89ad90SMark F. Adams   PetscFunctionBegin;
13875b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
13885b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
13895b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
13905b89ad90SMark F. Adams 
13915b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
13925b89ad90SMark F. Adams   ierr         = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr);
13935b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
1394ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
13955b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
13965b89ad90SMark F. Adams 
13971ab5ffc9SJed Brown   ierr = PetscNewLog(pc,struct _PCGAMGOps,&pc_gamg->ops);CHKERRQ(ierr);
13981ab5ffc9SJed Brown 
13999d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14009d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14019d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14029d5b6da9SMark F. Adams   pc_gamg->data    = 0;
14035b89ad90SMark F. Adams 
14049d5b6da9SMark F. Adams   /* register AMG type */
14053e3471ccSMark Adams   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
14069d5b6da9SMark F. Adams 
14079d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
14085b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
14095b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
14105b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
14115b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
14125b89ad90SMark F. Adams 
1413bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1414bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1415bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
1416bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseProl_C",PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr);
1417bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1418bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1419bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1420bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
14219d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1422d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
1423ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1424038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
1425c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit  = 800;
1426d3042614SMark Adams   pc_gamg->threshold        = 0.;
14279d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
14289d5b6da9SMark F. Adams   pc_gamg->verbose          = 0;
14299d5b6da9SMark F. Adams   pc_gamg->emax_id          = -1;
14305e7c91beSJed Brown   pc_gamg->eigtarget[0]     = 0.05;
14315e7c91beSJed Brown   pc_gamg->eigtarget[1]     = 1.05;
14329d5b6da9SMark F. Adams 
14330cbbd2e1SMark F. Adams   /* private events */
14340cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1435785cba28SMark F. Adams   if (count++ == 0) {
1436806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1437806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
14380cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
14390cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
14400cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1441806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1442806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1443806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1444806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1445806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1446806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1447806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1448806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1449806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1450806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1451806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1452806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1453f852f58cSMark F. Adams 
14540cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
14550cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
14560cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1457b4fbaa2aSMark F. Adams     /* create timer stages */
1458b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1459b4fbaa2aSMark F. Adams     {
1460b4fbaa2aSMark F. Adams       char     str[32];
1461b4fbaa2aSMark F. Adams       PetscInt lidx;
1462806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1463806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1464b4fbaa2aSMark F. Adams       for (lidx=1; lidx<9; lidx++) {
1465b4fbaa2aSMark F. Adams         sprintf(str,"MG Level %d",lidx);
1466806fa848SBarry Smith         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1467b4fbaa2aSMark F. Adams       }
1468b4fbaa2aSMark F. Adams     }
1469b4fbaa2aSMark F. Adams #endif
1470b4fbaa2aSMark F. Adams   }
1471b4fbaa2aSMark F. Adams #endif
1472bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1473bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
14745b89ad90SMark F. Adams   PetscFunctionReturn(0);
14755b89ad90SMark F. Adams }
14763e3471ccSMark Adams 
14773e3471ccSMark Adams #undef __FUNCT__
14783e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage"
14793e3471ccSMark Adams /*@C
14803e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
14813e3471ccSMark Adams  from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
14823e3471ccSMark Adams  when using static libraries.
14833e3471ccSMark Adams 
14843e3471ccSMark Adams  Level: developer
14853e3471ccSMark Adams 
14863e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
14873e3471ccSMark Adams  .seealso: PetscInitialize()
14883e3471ccSMark Adams @*/
14893e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
14903e3471ccSMark Adams {
14913e3471ccSMark Adams   PetscErrorCode ierr;
14923e3471ccSMark Adams 
14933e3471ccSMark Adams   PetscFunctionBegin;
14943e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
14953e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
14963e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
14973e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
14988e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
14993e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1500c1c463dbSMark Adams 
1501c1c463dbSMark Adams   /* general events */
1502c1c463dbSMark Adams #if defined PETSC_USE_LOG
1503c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1504c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1505c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1506c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1507c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1508c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1509c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1510c1c463dbSMark Adams #endif
1511c1c463dbSMark Adams 
15123e3471ccSMark Adams   PetscFunctionReturn(0);
15133e3471ccSMark Adams }
15143e3471ccSMark Adams 
15153e3471ccSMark Adams #undef __FUNCT__
15163e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage"
15173e3471ccSMark Adams /*@C
15183e3471ccSMark Adams  PCGAMGFinalizePackage - This function destroys everything in the PCGAMG package. It is
15193e3471ccSMark Adams  called from PetscFinalize().
15203e3471ccSMark Adams 
15213e3471ccSMark Adams  Level: developer
15223e3471ccSMark Adams 
15233e3471ccSMark Adams  .keywords: Petsc, destroy, package
15243e3471ccSMark Adams  .seealso: PetscFinalize()
15253e3471ccSMark Adams @*/
15263e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
15273e3471ccSMark Adams {
15283e3471ccSMark Adams   PetscErrorCode ierr;
15293e3471ccSMark Adams 
15303e3471ccSMark Adams   PetscFunctionBegin;
15313e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
15323e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
15333e3471ccSMark Adams   PetscFunctionReturn(0);
15343e3471ccSMark Adams }
1535