xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 5ca49e89713766d4e1b32ebf01b36ec902b56c10)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4aaa7dc30SBarry Smith #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 /*
57c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: 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
653530afc2SMark F. Adams    In/Output Parameter:
66a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
67afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
6811e60469SMark F. Adams    Output Parameter:
693530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
705b89ad90SMark F. Adams */
715cb416c2SMark F. Adams 
725b89ad90SMark F. Adams #undef __FUNCT__
73c238b0ebSToby Isaac #define __FUNCT__ "PCGAMGCreateLevel_GAMG"
74b34066adSToby Isaac static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,
753cb8563fSToby Isaac                                   Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,
763cb8563fSToby Isaac                                   IS * Pcolumnperm)
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;
81a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
823b4367a7SBarry Smith   MPI_Comm        comm;
83c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
843ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
855b89ad90SMark F. Adams 
865b89ad90SMark F. Adams   PetscFunctionBegin;
873b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
883b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
893b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
90c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
9111e60469SMark F. Adams   /* RAP */
929d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
93038e3b61SMark F. Adams 
943ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
950298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
963ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
973ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
983ae0bb68SMark Adams   }
993ae0bb68SMark Adams   else {
1003ae0bb68SMark Adams     PetscInt  bs;
1013ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
1023ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
1033ae0bb68SMark Adams   }
104a2f3521dSMark F. Adams 
105c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
106a2f3521dSMark F. Adams   {
107472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
1080298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
109a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
110a90e85d9SMark Adams     if (new_size == 0) new_size = 1; /* not likely, posible? */
111c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
112a2f3521dSMark F. Adams   }
113f852f58cSMark F. Adams 
1143cb8563fSToby Isaac   if (Pcolumnperm) *Pcolumnperm = NULL;
1153cb8563fSToby Isaac 
116a90e85d9SMark Adams   if (!pc_gamg->repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1172fa5cd67SKarl Rupp   else {
1183ae0bb68SMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
119885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
120e33ef3b1SMark F. Adams 
12171959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
12271959b99SBarry 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);
1230cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1240cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
125b4fbaa2aSMark F. Adams #endif
126a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
127785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
128a90e85d9SMark Adams     if (pc_gamg->repart) {
129a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1305a9b9e01SMark F. Adams       Mat adj;
1315a9b9e01SMark F. Adams 
132a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
1333b4367a7SBarry 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);
134a2f3521dSMark F. Adams         else {
135a2f3521dSMark F. Adams           PetscInt n;
1363b4367a7SBarry Smith           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
1373b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n);
138a2f3521dSMark F. Adams         }
139a2f3521dSMark F. Adams       }
1405a9b9e01SMark F. Adams 
141a2f3521dSMark F. Adams       /* get 'adj' */
142c5bfad50SMark F. Adams       if (cr_bs == 1) {
143038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
144806fa848SBarry Smith       } else {
145a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
146eb07cef2SMark F. Adams         Mat               tMat;
147a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
148b4fbaa2aSMark F. Adams         const PetscScalar *vals;
149b4fbaa2aSMark F. Adams         const PetscInt    *idx;
150a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1519057884aSMark F. Adams         static PetscInt   llev = 0;
152b4fbaa2aSMark F. Adams 
153578f55a3SPeter Brune         ierr = PetscMalloc1(ncrs, &d_nnz);CHKERRQ(ierr);
154578f55a3SPeter Brune         ierr = PetscMalloc1(ncrs, &o_nnz);CHKERRQ(ierr);
155a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
156a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
157c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
15858471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
159c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
160c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
16158471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1623ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1633ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
16458471d46SMark F. Adams         }
1656876a03eSMark F. Adams 
1663b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1673ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
168a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
169a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
170a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
17158471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
1725f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
173eb07cef2SMark F. Adams 
174a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
175c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
17622063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
177eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
178c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
179eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
180eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
181eb07cef2SMark F. Adams           }
18222063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
183eb07cef2SMark F. Adams         }
184eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
185eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
186eb07cef2SMark F. Adams 
187b4fbaa2aSMark F. Adams         if (llev++ == -1) {
188b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
1898caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
1903b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
191b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
1923bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
193b4fbaa2aSMark F. Adams         }
194b4fbaa2aSMark F. Adams 
195eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
196eb07cef2SMark F. Adams 
197eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
198a2f3521dSMark F. Adams       } /* create 'adj' */
199f150b916SMark F. Adams 
200a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2015a9b9e01SMark F. Adams         char            prefix[256];
2025a9b9e01SMark F. Adams         const char      *pcpre;
203b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
204b4fbaa2aSMark F. Adams         MatPartitioning mpart;
205a4b7d37bSMark F. Adams         IS              proc_is;
206a2f3521dSMark F. Adams         PetscInt        targetPE;
2072f03bc48SMark F. Adams 
2083b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2095ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2109d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2118caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
21259a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
21311e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
214c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
215a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
21611e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2175a9b9e01SMark F. Adams 
2185ef31b24SMark F. Adams         /* collect IS info */
219785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
220a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
221a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
222c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
223a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
224c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
225a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
226eb07cef2SMark F. Adams           }
2275ef31b24SMark F. Adams         }
228a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
229a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2305ef31b24SMark F. Adams       }
2315ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2325a9b9e01SMark F. Adams 
2333b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2348263b398SMark F. Adams       if (newproc_idx != 0) {
2358263b398SMark F. Adams         ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2365ef31b24SMark F. Adams       }
237806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
238a2f3521dSMark F. Adams 
239a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
2405a9b9e01SMark F. Adams       /* find factor */
241c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2425a9b9e01SMark F. Adams       else {
2435a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2445a9b9e01SMark F. Adams         jj = -1;
245c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
246c5df96a5SBarry Smith           if (size%kk==0) { /* a candidate */
247c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2485a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2495a9b9e01SMark F. Adams             if (fact > best_fact) {
2505a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
2515a9b9e01SMark F. Adams             }
2525a9b9e01SMark F. Adams           }
2535a9b9e01SMark F. Adams         }
2545a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
255a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
2565a9b9e01SMark F. Adams       }
257c5df96a5SBarry Smith       new_size = size/rfactor;
2585a9b9e01SMark F. Adams 
259c5df96a5SBarry Smith       if (new_size==nactive) {
260a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2615a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
262a2f3521dSMark F. Adams         if (pc_gamg->verbose>0) {
2633b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq);
264a2f3521dSMark F. Adams         }
2655a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2665a9b9e01SMark F. Adams       }
2675a9b9e01SMark F. Adams 
2683b4367a7SBarry Smith       if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq);
269c5df96a5SBarry Smith       targetPE = rank/rfactor;
2703b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
271a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
272e33ef3b1SMark F. Adams 
27311e60469SMark F. Adams     /*
274a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
27511e60469SMark F. Adams      */
276a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2777700e67bSMark Adams     is_eq_num_prim = is_eq_num;
27811e60469SMark F. Adams     /*
279a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
28011e60469SMark F. Adams      */
281c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
282c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
283a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
2843ae0bb68SMark Adams     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */
285a2f3521dSMark F. Adams 
286a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2870cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2880cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
289b4fbaa2aSMark F. Adams #endif
290885364a3SMark Adams     /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxilary data into some complex abstracted thing */
291885364a3SMark Adams     {
292885364a3SMark Adams     Vec            src_crd, dest_crd;
293885364a3SMark 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;
294885364a3SMark Adams     VecScatter     vecscat;
295885364a3SMark Adams     PetscScalar    *array;
296885364a3SMark Adams     IS isscat;
297a2f3521dSMark F. Adams 
298a2f3521dSMark F. Adams     /* move data (for primal equations only) */
29922063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
3003b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
3013ae0bb68SMark Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
302c0dedaeaSBarry Smith     ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
30311e60469SMark F. Adams     /*
3049d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
305c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
30611e60469SMark F. Adams      */
307854ce69bSBarry Smith     ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
308a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3093ae0bb68SMark Adams     for (ii=0,jj=0; ii<ncrs; ii++) {
310c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
311a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
31211e60469SMark F. Adams     }
313a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3143ae0bb68SMark Adams     ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
31592a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
31611e60469SMark F. Adams     /*
31711e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
31811e60469SMark F. Adams      */
3193ae0bb68SMark Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
3209d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3213ae0bb68SMark Adams       const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3223ae0bb68SMark Adams       for (ii=0; ii<ncrs; ii++) {
3239d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
324a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
325c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
326676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
327d3d6bff4SMark F. Adams         }
328038e3b61SMark F. Adams       }
329eb07cef2SMark F. Adams     }
330eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
331eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
33211e60469SMark F. Adams     /*
33311e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
33411e60469SMark F. Adams       to the correct processor
33511e60469SMark F. Adams     */
3360298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
33711e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
33811e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
33911e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
34011e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
34111e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
34211e60469SMark F. Adams     /*
34311e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
34411e60469SMark F. Adams     */
345c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
346578f55a3SPeter Brune     ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3472fa5cd67SKarl Rupp 
3483ae0bb68SMark Adams     pc_gamg->data_sz = node_data_sz*ncrs_new;
3493ae0bb68SMark Adams     strideNew        = ncrs_new*ndata_rows;
3502fa5cd67SKarl Rupp 
35111e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3529d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3533ae0bb68SMark Adams       for (ii=0; ii<ncrs_new; ii++) {
3549d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
355a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
356c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
357d3d6bff4SMark F. Adams         }
358038e3b61SMark F. Adams       }
359038e3b61SMark F. Adams     }
36011e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
36111e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
362885364a3SMark Adams     }
363a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3640cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3650cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
366ed3f9983SMark F. Adams #endif
367a2f3521dSMark F. Adams 
36811e60469SMark F. Adams     /*
36911e60469SMark F. Adams       Invert for MatGetSubMatrix
37011e60469SMark F. Adams     */
371a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
372a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
373c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
374a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
375a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
376a2f3521dSMark F. Adams     }
3773cb8563fSToby Isaac     if (Pcolumnperm) {
3783cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3793cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3803cb8563fSToby Isaac     }
381a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3820cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3830cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3840cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
385ed3f9983SMark F. Adams #endif
386a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
387a2f3521dSMark F. Adams     {
388a2f3521dSMark F. Adams       Mat mat;
389806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
390a2f3521dSMark F. Adams       *a_Amat_crs = mat;
391c5bfad50SMark F. Adams 
392c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
393c5bfad50SMark F. Adams         PetscInt cbs, rbs;
394c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
395c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
396c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
397c5df96a5SBarry 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);
398c5bfad50SMark F. Adams       }
399a2f3521dSMark F. Adams     }
400038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
401a2f3521dSMark F. Adams 
4020cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4030cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
404ed3f9983SMark F. Adams #endif
40511e60469SMark F. Adams     /* prolongator */
40611e60469SMark F. Adams     {
40711e60469SMark F. Adams       IS       findices;
408a2f3521dSMark F. Adams       PetscInt Istart,Iend;
409a2f3521dSMark F. Adams       Mat      Pnew;
410a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4110cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4120cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
413ed3f9983SMark F. Adams #endif
4143b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
415c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
416806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
41711e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
418c5bfad50SMark F. Adams 
419c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
420c5bfad50SMark F. Adams         PetscInt cbs, rbs;
421c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
422c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
423c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
424c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
425c5bfad50SMark F. Adams       }
4260cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4270cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
428ed3f9983SMark F. Adams #endif
4293530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
430a2f3521dSMark F. Adams 
431a2f3521dSMark F. Adams       /* output - repartitioned */
432a2f3521dSMark F. Adams       *a_P_inout = Pnew;
433e33ef3b1SMark F. Adams     }
434a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4355b89ad90SMark F. Adams 
436c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
437a2f3521dSMark F. Adams   }
4385a9b9e01SMark F. Adams 
439a2f3521dSMark F. Adams   /* outout matrix data */
440c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
441c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
442c8b0795cSMark F. Adams     if (llev==0) {
443c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
4443b4367a7SBarry Smith       PetscViewerASCIIOpen(comm,fname,&viewer);
445c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
446c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
447c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
448c8b0795cSMark F. Adams     }
449c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
4503b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
451c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
452c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
453c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
454c8b0795cSMark F. Adams   }
4555b89ad90SMark F. Adams   PetscFunctionReturn(0);
4565b89ad90SMark F. Adams }
4575b89ad90SMark F. Adams 
4585b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4595b89ad90SMark F. Adams /*
4605b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4615b89ad90SMark F. Adams                     by setting data structures and options.
4625b89ad90SMark F. Adams 
4635b89ad90SMark F. Adams    Input Parameter:
4645b89ad90SMark F. Adams .  pc - the preconditioner context
4655b89ad90SMark F. Adams 
4665b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4675b89ad90SMark F. Adams 
4685b89ad90SMark F. Adams    Notes:
4695b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4705b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4715b89ad90SMark F. Adams */
4725b89ad90SMark F. Adams #undef __FUNCT__
4735b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4749d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4755b89ad90SMark F. Adams {
4765b89ad90SMark F. Adams   PetscErrorCode ierr;
4779d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4785b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4792adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
480a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
4813b4367a7SBarry Smith   MPI_Comm       comm;
482c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
483587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
484c8b0795cSMark F. Adams   PetscReal      emaxs[GAMG_MAXLEVELS];
485e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
486a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
487737a81a9SMark F. Adams   MatInfo        info;
4887700e67bSMark Adams   PetscBool      redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);
4895ef31b24SMark F. Adams 
4905b89ad90SMark F. Adams   PetscFunctionBegin;
4913b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4923b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4933b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
494dfd5c07aSMark F. Adams 
4953b4367a7SBarry 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);
496dfd5c07aSMark F. Adams 
49784d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
498878e152fSMark F. Adams     if (redo_mesh_setup) {
499878e152fSMark F. Adams       /* reset everything */
500878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
501878e152fSMark F. Adams       pc->setupcalled = 0;
502806fa848SBarry Smith     } else {
50384d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
50403a628feSMark F. Adams       /* just do Galerkin grids */
50558471d46SMark F. Adams       Mat          B,dA,dB;
50658471d46SMark F. Adams 
50771959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
5089d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
50958471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
51023ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
51158471d46SMark F. Adams         /* (re)set to get dirty flag */
51223ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
51358471d46SMark F. Adams 
5142fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
51503a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
5160a97e771SToby Isaac           if (pc_gamg->setup_count==2) {
51703a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
518084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
5192fa5cd67SKarl Rupp 
52003a628feSMark F. Adams             mglevels[level]->A = B;
521806fa848SBarry Smith           } else {
52223ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
52358471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
52403a628feSMark F. Adams           }
52523ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
52658471d46SMark F. Adams           dB   = B;
52758471d46SMark F. Adams         }
5285f8cf99dSMark F. Adams       }
529d5280255SMark F. Adams 
530d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
531d5280255SMark F. Adams 
53258471d46SMark F. Adams       PetscFunctionReturn(0);
533eb07cef2SMark F. Adams     }
534878e152fSMark F. Adams   }
535f6536408SMark F. Adams 
536878e152fSMark F. Adams   if (!pc_gamg->data) {
537878e152fSMark F. Adams     if (pc_gamg->orig_data) {
538878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5390298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5402fa5cd67SKarl Rupp 
541878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
542878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
543878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5442fa5cd67SKarl Rupp 
545785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
546878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
547806fa848SBarry Smith     } else {
5481ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5497700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5509d5b6da9SMark F. Adams     }
551878e152fSMark F. Adams   }
552878e152fSMark F. Adams 
553878e152fSMark F. Adams   /* cache original data for reuse */
554878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
555785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
556878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
557878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
558878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
559878e152fSMark F. Adams   }
560038e3b61SMark F. Adams 
561302f38e8SMark F. Adams   /* get basic dims */
562302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
563a2f3521dSMark F. Adams 
564a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
565c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
56684f9421dSMark F. Adams     PetscInt NN = M;
56784f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
56884f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
5693bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
570fe06f982SMark Adams       if (!NN) NN=1;
571806fa848SBarry Smith     } else {
572806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
57384f9421dSMark F. Adams     }
574b2a4f308SMark F. Adams     nnz0   = info.nz_used;
575b2a4f308SMark F. Adams     nnztot = info.nz_used;
5763b4367a7SBarry 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",
577c5df96a5SBarry Smith                          rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
578c5df96a5SBarry Smith                          (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr);
579c8b0795cSMark F. Adams   }
58084d3f75bSMark F. Adams 
581a2f3521dSMark F. Adams   /* Get A_i and R_i */
582c5df96a5SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
583a90e85d9SMark Adams        level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit);
5840205a208SMark F. Adams        level++) {
58557d29afaSToby Isaac     pc_gamg->firstCoarsen = (level ? PETSC_FALSE : PETSC_TRUE);
5865b89ad90SMark F. Adams     level1 = level + 1;
5870cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5880cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
589a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
590a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
591b4fbaa2aSMark F. Adams #endif
592a2f3521dSMark F. Adams #endif
593c8b0795cSMark F. Adams     { /* construct prolongator */
594725b86d8SJed Brown       Mat              Gmat;
5950cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5967700e67bSMark Adams       Mat              Prol11;
597c8b0795cSMark F. Adams 
5987700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5991ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
6007700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
601c8b0795cSMark F. Adams 
602a2f3521dSMark F. Adams       /* could have failed to create new level */
603a2f3521dSMark F. Adams       if (Prol11) {
6049d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6050298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
606a2f3521dSMark F. Adams 
6071ab5ffc9SJed Brown         if (pc_gamg->ops->optprol) {
608c8b0795cSMark F. Adams           /* smooth */
6097700e67bSMark Adams           ierr = pc_gamg->ops->optprol(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
610c8b0795cSMark F. Adams         }
611c8b0795cSMark F. Adams 
6127700e67bSMark Adams         Parr[level1] = Prol11;
6130298fd71SBarry Smith       } else Parr[level1] = NULL;
614ffc955d6SMark F. Adams 
615ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
6161b18a24aSMark Adams         PetscInt bs;
6171b18a24aSMark Adams         ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
618806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
619ffc955d6SMark F. Adams       }
620ffc955d6SMark F. Adams 
621a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
62241b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
623a2f3521dSMark F. Adams     } /* construct prolongator scope */
6240cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6250cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
626c8b0795cSMark F. Adams #endif
6279d5b6da9SMark F. Adams     /* cache eigen estimate */
6289d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
6299d5b6da9SMark F. Adams       PetscBool flag;
6307700e67bSMark Adams       ierr = PetscObjectComposedDataGetReal((PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
6319d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
632806fa848SBarry Smith     } else emaxs[level] = -1.;
6332adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
634c8b0795cSMark F. Adams     if (!Parr[level1]) {
635806fa848SBarry Smith       if (pc_gamg->verbose) {
6363b4367a7SBarry Smith         ierr =  PetscPrintf(comm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr);
637806fa848SBarry Smith       }
638a90e85d9SMark Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
639a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
640a90e85d9SMark Adams #endif
641c8b0795cSMark F. Adams       break;
642c8b0795cSMark F. Adams     }
6430cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6440cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
645b4fbaa2aSMark F. Adams #endif
646a2f3521dSMark F. Adams 
647c238b0ebSToby Isaac     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs,
6483cb8563fSToby Isaac                        &Parr[level1], &Aarr[level1], &nactivepe, NULL);CHKERRQ(ierr);
649a2f3521dSMark F. Adams 
6500cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6510cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
652b4fbaa2aSMark F. Adams #endif
653a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
654a2f3521dSMark F. Adams 
655a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0) {
6560cbbd2e1SMark F. Adams       PetscInt NN = M;
6570cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
658a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
6593bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
660fe06f982SMark Adams         if (!NN) NN=1;
661806fa848SBarry Smith       } else {
662806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
6630cbbd2e1SMark F. Adams       }
664a2f3521dSMark F. Adams 
6650cbbd2e1SMark F. Adams       nnztot += info.nz_used;
6663b4367a7SBarry Smith       ierr    = PetscPrintf(comm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
667c5df96a5SBarry Smith                             rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
668806fa848SBarry Smith                             (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
669c8b0795cSMark F. Adams     }
6700cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
671b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
672b4fbaa2aSMark F. Adams #endif
673a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
674a90e85d9SMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M < 2) ) {
675a90e85d9SMark Adams       if (pc_gamg->verbose) {
676a90e85d9SMark Adams         ierr =  PetscPrintf(comm,"\t[%d]%s HARD stop of coarsening ?????????, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr);
677a90e85d9SMark Adams       }
678a90e85d9SMark Adams       level++;
679a90e85d9SMark Adams       break;
680a90e85d9SMark Adams     }
681c8b0795cSMark F. Adams   } /* levels */
68257d29afaSToby Isaac   pc_gamg->firstCoarsen = PETSC_FALSE;
683c8b0795cSMark F. Adams 
684c8b0795cSMark F. Adams   if (pc_gamg->data) {
685c8b0795cSMark F. Adams     ierr          = PetscFree(pc_gamg->data);CHKERRQ(ierr);
6860298fd71SBarry Smith     pc_gamg->data = NULL;
6875b89ad90SMark F. Adams   }
688c8b0795cSMark F. Adams 
6893b4367a7SBarry Smith   if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
6909d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6915b89ad90SMark F. Adams   fine_level       = level;
6920298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
6935b89ad90SMark F. Adams 
69484d3f75bSMark F. Adams   /* simple setup */
69584d3f75bSMark F. Adams   if (!PETSC_TRUE) {
69684d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
69784d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
69884d3f75bSMark F. Adams          lidx<fine_level;
69984d3f75bSMark F. Adams          lidx++, level--) {
70084d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
70123ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level]);CHKERRQ(ierr);
70284d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
70384d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
70484d3f75bSMark F. Adams     }
70523ee1639SBarry Smith     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0]);CHKERRQ(ierr);
70684d3f75bSMark F. Adams 
70784d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
708806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
709d5280255SMark F. Adams     /* set default smoothers & set operators */
7109d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
711587fa25dSMark F. Adams          lidx <= fine_level;
712587fa25dSMark F. Adams          lidx++, level--) {
713ffc955d6SMark F. Adams       KSP smoother;
714ffc955d6SMark F. Adams       PC  subpc;
715a2f3521dSMark F. Adams 
7169d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
717f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
718ffc955d6SMark F. Adams 
719a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
720a2f3521dSMark F. Adams       /* set ops */
72123ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
722a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
723a2f3521dSMark F. Adams 
724a2f3521dSMark F. Adams       /* set defaults */
7256c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
726a2f3521dSMark F. Adams 
7271b18a24aSMark Adams       /* set blocks for GASM smoother that uses the 'aggregates' */
728ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
7292d3561bbSSatish Balay         PetscInt sz;
7302d3561bbSSatish Balay         IS       *is;
731a2f3521dSMark F. Adams 
7322d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
7332d3561bbSSatish Balay         is   = ASMLocalIDsArr[level];
734ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
7351b18a24aSMark Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
736ffc955d6SMark F. Adams         if (sz==0) {
737ffc955d6SMark F. Adams           IS       is;
738ffc955d6SMark F. Adams           PetscInt my0,kk;
739ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
740ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
7410298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr);
742a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
743806fa848SBarry Smith         } else {
744a94c3b12SMark F. Adams           PetscInt kk;
7450298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr);
746a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
747a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
748a94c3b12SMark F. Adams           }
749ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
750ffc955d6SMark F. Adams         }
7510298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
752ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
753ffc955d6SMark F. Adams         ierr                  = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
754806fa848SBarry Smith       } else {
755890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
756ffc955d6SMark F. Adams       }
757d5280255SMark F. Adams     }
758d5280255SMark F. Adams     {
759d5280255SMark F. Adams       /* coarse grid */
760d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
761d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
762d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
76323ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
764d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
765d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
766d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
767d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
76871959b99SBarry Smith       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
76971959b99SBarry Smith       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
770d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
771d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
7729dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
7732fb0b348SMark F. Adams       ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
7745b42dca8SJed Brown       /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
7755b42dca8SJed Brown        * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
7765b42dca8SJed Brown        * view every subdomain as though they were different. */
7775b42dca8SJed Brown       ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
778d5280255SMark F. Adams     }
779d5280255SMark F. Adams 
780d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
781d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
782e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
783d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
7843b4367a7SBarry Smith     if (mg->galerkin != 2) SETERRQ(comm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
785d5280255SMark F. Adams 
786d5280255SMark F. Adams     /* create cheby smoothers */
787d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
788d5280255SMark F. Adams          lidx <= fine_level;
789d5280255SMark F. Adams          lidx++, level--) {
790d5280255SMark F. Adams       KSP       smoother;
791890ffe84SMark Adams       PetscBool flag,flag2;
792d5280255SMark F. Adams       PC        subpc;
793d5280255SMark F. Adams 
794ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
795a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
796a2f3521dSMark F. Adams 
797ffc955d6SMark F. Adams       /* do my own cheby */
7986c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
799ffc955d6SMark F. Adams       if (flag) {
800ffc955d6SMark F. Adams         PetscReal emax, emin;
801251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
802890ffe84SMark Adams         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCSOR, &flag2);CHKERRQ(ierr);
803*5ca49e89SMark Adams         /* eigen estimate only for diagnal PC but lets acccept SOR because it is close and safe (always lower) */
804*5ca49e89SMark Adams         if ((flag||flag2) && (emax=emaxs[level]) > 0.0) {
805c5bfad50SMark F. Adams           PetscInt N1, N0;
806*5ca49e89SMark Adams           emax=emaxs[level];
8070298fd71SBarry Smith           ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr);
8080298fd71SBarry Smith           ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr);
8095e7c91beSJed Brown           emin  = emax * pc_gamg->eigtarget[0];
8105e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
8116c9de887SHong Zhang           ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
812*5ca49e89SMark Adams         }
813ffc955d6SMark F. Adams       } /* setup checby flag */
814ffc955d6SMark F. Adams     } /* non-coarse levels */
815737a81a9SMark F. Adams 
816d5280255SMark F. Adams     /* clean up */
817d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
818587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
819587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
8205b89ad90SMark F. Adams     }
8210cbbd2e1SMark F. Adams 
8220cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
823806fa848SBarry Smith   } else {
8245f8cf99dSMark F. Adams     KSP smoother;
8253b4367a7SBarry 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__);
8269d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
82723ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
8285f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
8299d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8305f8cf99dSMark F. Adams   }
8315b89ad90SMark F. Adams   PetscFunctionReturn(0);
8325b89ad90SMark F. Adams }
8335b89ad90SMark F. Adams 
834eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8355b89ad90SMark F. Adams /*
8365b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8375b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8385b89ad90SMark F. Adams 
8395b89ad90SMark F. Adams    Input Parameter:
8405b89ad90SMark F. Adams .  pc - the preconditioner context
8415b89ad90SMark F. Adams 
8425b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8435b89ad90SMark F. Adams */
8445b89ad90SMark F. Adams #undef __FUNCT__
8455b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
8465b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8475b89ad90SMark F. Adams {
8485b89ad90SMark F. Adams   PetscErrorCode ierr;
8495b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
8505b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
8515b89ad90SMark F. Adams 
8525b89ad90SMark F. Adams   PetscFunctionBegin;
8535b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8549b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
8559b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
8569b8ffb57SJed Brown   }
8571ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
8581ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
8595b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8605b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8615b89ad90SMark F. Adams   PetscFunctionReturn(0);
8625b89ad90SMark F. Adams }
8635b89ad90SMark F. Adams 
864676e1743SMark F. Adams 
865676e1743SMark F. Adams #undef __FUNCT__
866676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
867676e1743SMark F. Adams /*@
8681cc46a46SBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via processor reduction.
869676e1743SMark F. Adams 
8701cc46a46SBarry Smith    Logically Collective on PC
871676e1743SMark F. Adams 
872676e1743SMark F. Adams    Input Parameters:
8731cc46a46SBarry Smith +  pc - the preconditioner context
8741cc46a46SBarry Smith -  n - the number of equations
875676e1743SMark F. Adams 
876676e1743SMark F. Adams 
877676e1743SMark F. Adams    Options Database Key:
8781cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
879676e1743SMark F. Adams 
880676e1743SMark F. Adams    Level: intermediate
881676e1743SMark F. Adams 
882676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
883676e1743SMark F. Adams 
884676e1743SMark F. Adams .seealso: ()
885676e1743SMark F. Adams @*/
886676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
887676e1743SMark F. Adams {
888676e1743SMark F. Adams   PetscErrorCode ierr;
889676e1743SMark F. Adams 
890676e1743SMark F. Adams   PetscFunctionBegin;
891676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
892676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
893676e1743SMark F. Adams   PetscFunctionReturn(0);
894676e1743SMark F. Adams }
895676e1743SMark F. Adams 
896676e1743SMark F. Adams #undef __FUNCT__
897676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
8981e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
899676e1743SMark F. Adams {
900c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
901c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
902676e1743SMark F. Adams 
903676e1743SMark F. Adams   PetscFunctionBegin;
9049d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
905676e1743SMark F. Adams   PetscFunctionReturn(0);
906676e1743SMark F. Adams }
907676e1743SMark F. Adams 
908676e1743SMark F. Adams #undef __FUNCT__
909389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
910389730f3SMark F. Adams /*@
911389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
912389730f3SMark F. Adams 
913389730f3SMark F. Adams  Collective on PC
914389730f3SMark F. Adams 
915389730f3SMark F. Adams    Input Parameters:
9161cc46a46SBarry Smith +  pc - the preconditioner context
9171cc46a46SBarry Smith -  n - maximum number of equations to aim for
918389730f3SMark F. Adams 
919389730f3SMark F. Adams    Options Database Key:
9201cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
921389730f3SMark F. Adams 
922389730f3SMark F. Adams    Level: intermediate
923389730f3SMark F. Adams 
924389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
925389730f3SMark F. Adams 
926389730f3SMark F. Adams @*/
927389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
928389730f3SMark F. Adams {
929389730f3SMark F. Adams   PetscErrorCode ierr;
930389730f3SMark F. Adams 
931389730f3SMark F. Adams   PetscFunctionBegin;
932389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
933389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
934389730f3SMark F. Adams   PetscFunctionReturn(0);
935389730f3SMark F. Adams }
936389730f3SMark F. Adams 
937389730f3SMark F. Adams #undef __FUNCT__
938389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
9391e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
940389730f3SMark F. Adams {
941389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
942389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
943389730f3SMark F. Adams 
944389730f3SMark F. Adams   PetscFunctionBegin;
9459d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
946389730f3SMark F. Adams   PetscFunctionReturn(0);
947389730f3SMark F. Adams }
948389730f3SMark F. Adams 
949389730f3SMark F. Adams #undef __FUNCT__
9508263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
951676e1743SMark F. Adams /*@
9528263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
953676e1743SMark F. Adams 
954676e1743SMark F. Adams    Collective on PC
955676e1743SMark F. Adams 
956676e1743SMark F. Adams    Input Parameters:
9571cc46a46SBarry Smith +  pc - the preconditioner context
9581cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
959676e1743SMark F. Adams 
960676e1743SMark F. Adams    Options Database Key:
9611cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
962676e1743SMark F. Adams 
963676e1743SMark F. Adams    Level: intermediate
964676e1743SMark F. Adams 
965676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
966676e1743SMark F. Adams 
967676e1743SMark F. Adams .seealso: ()
968676e1743SMark F. Adams @*/
9698263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
970676e1743SMark F. Adams {
971676e1743SMark F. Adams   PetscErrorCode ierr;
972676e1743SMark F. Adams 
973676e1743SMark F. Adams   PetscFunctionBegin;
974676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9758263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
976676e1743SMark F. Adams   PetscFunctionReturn(0);
977676e1743SMark F. Adams }
978676e1743SMark F. Adams 
979676e1743SMark F. Adams #undef __FUNCT__
9808263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
9811e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
982676e1743SMark F. Adams {
983c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
984c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
985676e1743SMark F. Adams 
986676e1743SMark F. Adams   PetscFunctionBegin;
9879d5b6da9SMark F. Adams   pc_gamg->repart = n;
988676e1743SMark F. Adams   PetscFunctionReturn(0);
989676e1743SMark F. Adams }
990676e1743SMark F. Adams 
991676e1743SMark F. Adams #undef __FUNCT__
9921cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation"
993dfd5c07aSMark F. Adams /*@
9941cc46a46SBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner
995dfd5c07aSMark F. Adams 
996dfd5c07aSMark F. Adams    Collective on PC
997dfd5c07aSMark F. Adams 
998dfd5c07aSMark F. Adams    Input Parameters:
9991cc46a46SBarry Smith +  pc - the preconditioner context
10001cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1001dfd5c07aSMark F. Adams 
1002dfd5c07aSMark F. Adams    Options Database Key:
10031cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
1004dfd5c07aSMark F. Adams 
1005dfd5c07aSMark F. Adams    Level: intermediate
1006dfd5c07aSMark F. Adams 
1007dfd5c07aSMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1008dfd5c07aSMark F. Adams 
1009dfd5c07aSMark F. Adams .seealso: ()
1010dfd5c07aSMark F. Adams @*/
10111cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1012dfd5c07aSMark F. Adams {
1013dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1014dfd5c07aSMark F. Adams 
1015dfd5c07aSMark F. Adams   PetscFunctionBegin;
1016dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10171cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1018dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1019dfd5c07aSMark F. Adams }
1020dfd5c07aSMark F. Adams 
1021dfd5c07aSMark F. Adams #undef __FUNCT__
10221cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG"
10231cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1024dfd5c07aSMark F. Adams {
1025dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1026dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1027dfd5c07aSMark F. Adams 
1028dfd5c07aSMark F. Adams   PetscFunctionBegin;
1029dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1030dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1031dfd5c07aSMark F. Adams }
1032dfd5c07aSMark F. Adams 
1033dfd5c07aSMark F. Adams #undef __FUNCT__
1034ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1035ffc955d6SMark F. Adams /*@
1036ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1037ffc955d6SMark F. Adams 
1038ffc955d6SMark F. Adams    Collective on PC
1039ffc955d6SMark F. Adams 
1040ffc955d6SMark F. Adams    Input Parameters:
1041ffc955d6SMark F. Adams .  pc - the preconditioner context
1042ffc955d6SMark F. Adams 
1043ffc955d6SMark F. Adams 
1044ffc955d6SMark F. Adams    Options Database Key:
1045ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1046ffc955d6SMark F. Adams 
1047ffc955d6SMark F. Adams    Level: intermediate
1048ffc955d6SMark F. Adams 
1049ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1050ffc955d6SMark F. Adams 
1051ffc955d6SMark F. Adams .seealso: ()
1052ffc955d6SMark F. Adams @*/
1053ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1054ffc955d6SMark F. Adams {
1055ffc955d6SMark F. Adams   PetscErrorCode ierr;
1056ffc955d6SMark F. Adams 
1057ffc955d6SMark F. Adams   PetscFunctionBegin;
1058ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1059ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1060ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1061ffc955d6SMark F. Adams }
1062ffc955d6SMark F. Adams 
1063ffc955d6SMark F. Adams #undef __FUNCT__
1064ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
10651e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1066ffc955d6SMark F. Adams {
1067ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1068ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1069ffc955d6SMark F. Adams 
1070ffc955d6SMark F. Adams   PetscFunctionBegin;
1071ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1072ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1073ffc955d6SMark F. Adams }
1074ffc955d6SMark F. Adams 
1075ffc955d6SMark F. Adams #undef __FUNCT__
10764ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
10774ef23d27SMark F. Adams /*@
10781cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
10794ef23d27SMark F. Adams 
10804ef23d27SMark F. Adams    Not collective on PC
10814ef23d27SMark F. Adams 
10824ef23d27SMark F. Adams    Input Parameters:
10831cc46a46SBarry Smith +  pc - the preconditioner
10841cc46a46SBarry Smith -  n - the maximum number of levels to use
10854ef23d27SMark F. Adams 
10864ef23d27SMark F. Adams    Options Database Key:
10874ef23d27SMark F. Adams .  -pc_mg_levels
10884ef23d27SMark F. Adams 
10894ef23d27SMark F. Adams    Level: intermediate
10904ef23d27SMark F. Adams 
10914ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
10924ef23d27SMark F. Adams 
10934ef23d27SMark F. Adams .seealso: ()
10944ef23d27SMark F. Adams @*/
10954ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
10964ef23d27SMark F. Adams {
10974ef23d27SMark F. Adams   PetscErrorCode ierr;
10984ef23d27SMark F. Adams 
10994ef23d27SMark F. Adams   PetscFunctionBegin;
11004ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11014ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
11024ef23d27SMark F. Adams   PetscFunctionReturn(0);
11034ef23d27SMark F. Adams }
11044ef23d27SMark F. Adams 
11054ef23d27SMark F. Adams #undef __FUNCT__
11064ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
11071e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
11084ef23d27SMark F. Adams {
11094ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
11104ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
11114ef23d27SMark F. Adams 
11124ef23d27SMark F. Adams   PetscFunctionBegin;
11139d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
11144ef23d27SMark F. Adams   PetscFunctionReturn(0);
11154ef23d27SMark F. Adams }
11164ef23d27SMark F. Adams 
11174ef23d27SMark F. Adams #undef __FUNCT__
11183542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
11193542efc5SMark F. Adams /*@
11203542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
11213542efc5SMark F. Adams 
11223542efc5SMark F. Adams    Not collective on PC
11233542efc5SMark F. Adams 
11243542efc5SMark F. Adams    Input Parameters:
11251cc46a46SBarry Smith +  pc - the preconditioner context
11261cc46a46SBarry Smith -  threshold - the threshold value, 0.0 is the lowest possible
11273542efc5SMark F. Adams 
11283542efc5SMark F. Adams    Options Database Key:
11291cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
11303542efc5SMark F. Adams 
11313542efc5SMark F. Adams    Level: intermediate
11323542efc5SMark F. Adams 
11333542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
11343542efc5SMark F. Adams 
11353542efc5SMark F. Adams .seealso: ()
11363542efc5SMark F. Adams @*/
11373542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
11383542efc5SMark F. Adams {
11393542efc5SMark F. Adams   PetscErrorCode ierr;
11403542efc5SMark F. Adams 
11413542efc5SMark F. Adams   PetscFunctionBegin;
11423542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11433542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
11443542efc5SMark F. Adams   PetscFunctionReturn(0);
11453542efc5SMark F. Adams }
11463542efc5SMark F. Adams 
11473542efc5SMark F. Adams #undef __FUNCT__
11483542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
11491e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
11503542efc5SMark F. Adams {
1151c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1152c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
11533542efc5SMark F. Adams 
11543542efc5SMark F. Adams   PetscFunctionBegin;
11559d5b6da9SMark F. Adams   pc_gamg->threshold = n;
11563542efc5SMark F. Adams   PetscFunctionReturn(0);
11573542efc5SMark F. Adams }
11583542efc5SMark F. Adams 
11593542efc5SMark F. Adams #undef __FUNCT__
11609d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1161676e1743SMark F. Adams /*@
1162c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1163676e1743SMark F. Adams 
1164676e1743SMark F. Adams    Collective on PC
1165676e1743SMark F. Adams 
1166676e1743SMark F. Adams    Input Parameters:
1167c60c7ad4SBarry Smith +  pc - the preconditioner context
1168c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1169676e1743SMark F. Adams 
1170676e1743SMark F. Adams    Options Database Key:
1171c60c7ad4SBarry Smith .  -pc_gamg_type <agg,geo,classical>
1172676e1743SMark F. Adams 
1173676e1743SMark F. Adams    Level: intermediate
1174676e1743SMark F. Adams 
1175676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1176676e1743SMark F. Adams 
1177c60c7ad4SBarry Smith .seealso: PCGAMGGetType(), PCGAMG
1178676e1743SMark F. Adams @*/
117919fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1180676e1743SMark F. Adams {
1181676e1743SMark F. Adams   PetscErrorCode ierr;
1182676e1743SMark F. Adams 
1183676e1743SMark F. Adams   PetscFunctionBegin;
1184676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1185806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1186676e1743SMark F. Adams   PetscFunctionReturn(0);
1187676e1743SMark F. Adams }
1188676e1743SMark F. Adams 
1189676e1743SMark F. Adams #undef __FUNCT__
1190c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType"
1191c60c7ad4SBarry Smith /*@
1192c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1193c60c7ad4SBarry Smith 
1194c60c7ad4SBarry Smith    Collective on PC
1195c60c7ad4SBarry Smith 
1196c60c7ad4SBarry Smith    Input Parameter:
1197c60c7ad4SBarry Smith .  pc - the preconditioner context
1198c60c7ad4SBarry Smith 
1199c60c7ad4SBarry Smith    Output Parameter:
1200c60c7ad4SBarry Smith .  type - the type of algorithm used
1201c60c7ad4SBarry Smith 
1202c60c7ad4SBarry Smith    Level: intermediate
1203c60c7ad4SBarry Smith 
1204c60c7ad4SBarry Smith    Concepts: Unstructured multrigrid preconditioner
1205c60c7ad4SBarry Smith 
1206c60c7ad4SBarry Smith .seealso: PCGAMGSetType()
1207c60c7ad4SBarry Smith @*/
1208c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1209c60c7ad4SBarry Smith {
1210c60c7ad4SBarry Smith   PetscErrorCode ierr;
1211c60c7ad4SBarry Smith 
1212c60c7ad4SBarry Smith   PetscFunctionBegin;
1213c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1214c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1215c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1216c60c7ad4SBarry Smith }
1217c60c7ad4SBarry Smith 
1218c60c7ad4SBarry Smith #undef __FUNCT__
1219c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG"
1220c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1221c60c7ad4SBarry Smith {
1222c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1223c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1224c60c7ad4SBarry Smith 
1225c60c7ad4SBarry Smith   PetscFunctionBegin;
1226c60c7ad4SBarry Smith   *type = pc_gamg->type;
1227c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1228c60c7ad4SBarry Smith }
1229c60c7ad4SBarry Smith 
1230c60c7ad4SBarry Smith #undef __FUNCT__
12319d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
12321e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1233676e1743SMark F. Adams {
12349d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
12351ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
12361ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1237676e1743SMark F. Adams 
1238676e1743SMark F. Adams   PetscFunctionBegin;
1239c60c7ad4SBarry Smith   pc_gamg->type = type;
12401c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
12419d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
12421ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
12433ae0bb68SMark Adams     /* there was something here - kill it */
12441ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
12451ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1246e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
12473ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
12483ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
12493ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
12503ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
12513ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
12523ae0bb68SMark Adams     if (pc_gamg->data_sz) {
12533ae0bb68SMark Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
12543ae0bb68SMark Adams       pc_gamg->data_sz = 0;
1255c60c7ad4SBarry Smith     } else if (pc_gamg->data) {
12563ae0bb68SMark Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); /* can this happen ? */
12573ae0bb68SMark Adams     }
12581ab5ffc9SJed Brown   }
12591ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
12601ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
12619d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1262676e1743SMark F. Adams   PetscFunctionReturn(0);
1263676e1743SMark F. Adams }
1264676e1743SMark F. Adams 
12655b89ad90SMark F. Adams #undef __FUNCT__
12665b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
12678c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptions *PetscOptionsObject,PC pc)
12685b89ad90SMark F. Adams {
1269676e1743SMark F. Adams   PetscErrorCode ierr;
1270676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1271676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1272676e1743SMark F. Adams   PetscBool      flag;
12735e7c91beSJed Brown   PetscInt       two   = 2;
12743b4367a7SBarry Smith   MPI_Comm       comm;
12755b89ad90SMark F. Adams 
12765b89ad90SMark F. Adams   PetscFunctionBegin;
12773b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1278e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1279676e1743SMark F. Adams   {
1280b7cbab4eSMark Adams     /* -pc_gamg_type */
1281b7cbab4eSMark Adams     {
1282bd94a7aaSJed Brown       char tname[256];
12831a1c1e04SBarry Smith       ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1284bd94a7aaSJed Brown       if (flag) {
1285bd94a7aaSJed Brown         ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
12861ab5ffc9SJed Brown       }
1287b7cbab4eSMark Adams     }
128875b74bdaSMark F. Adams     /* -pc_gamg_verbose */
128994ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG","none", pc_gamg->verbose,&pc_gamg->verbose, NULL);CHKERRQ(ierr);
12908263b398SMark F. Adams     /* -pc_gamg_repartition */
129194ae4db5SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
1292dfd5c07aSMark F. Adams     /* -pc_gamg_reuse_interpolation */
12931cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1294ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
129594ae4db5SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm","Use aggregation agragates for GASM smoother","PCGAMGUseASMAggs",pc_gamg->use_aggs_in_gasm,&pc_gamg->use_aggs_in_gasm,NULL);CHKERRQ(ierr);
1296c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
129794ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit","Limit (goal) on number of equations per process on coarse grids","PCGAMGSetProcEqLim",pc_gamg->min_eq_proc,&pc_gamg->min_eq_proc,NULL);CHKERRQ(ierr);
1298389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
129994ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit","Limit on number of equations for the coarse grid","PCGAMGSetCoarseEqLim",pc_gamg->coarse_eq_limit,&pc_gamg->coarse_eq_limit,NULL);CHKERRQ(ierr);
1300c20e4228SMark F. Adams     /* -pc_gamg_threshold */
130194ae4db5SBarry Smith     ierr = PetscOptionsReal("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&pc_gamg->threshold,&flag);CHKERRQ(ierr);
1302806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
13033b4367a7SBarry Smith       ierr = PetscPrintf(comm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1304806fa848SBarry Smith     }
1305b7cbab4eSMark Adams     /* -pc_gamg_eigtarget */
13060298fd71SBarry Smith     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr);
130794ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1308b7cbab4eSMark Adams 
1309b7cbab4eSMark Adams     /* set options for subtype */
1310e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1311676e1743SMark F. Adams   }
1312676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
13135b89ad90SMark F. Adams   PetscFunctionReturn(0);
13145b89ad90SMark F. Adams }
13155b89ad90SMark F. Adams 
13165b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
13175b89ad90SMark F. Adams /*MC
13181cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
13195b89ad90SMark F. Adams 
1320280d9858SJed Brown    Options Database Keys:
13215b89ad90SMark F. Adams    Multigrid options(inherited)
13221cc46a46SBarry Smith +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1323280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1324280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
13258c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
13265b89ad90SMark F. Adams 
13271cc46a46SBarry Smith 
13281cc46a46SBarry Smith   Notes: In order to obtain good performance for PCGAMG for vector valued problems you must
13291cc46a46SBarry Smith $       Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
13301cc46a46SBarry Smith $       Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
13311cc46a46SBarry Smith $       See the Users Manual Chapter 4 for more details
13321cc46a46SBarry Smith 
13335b89ad90SMark F. Adams   Level: intermediate
1334280d9858SJed Brown 
13351cc46a46SBarry Smith   Concepts: algebraic multigrid
13365b89ad90SMark F. Adams 
13371cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
13381cc46a46SBarry Smith            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartitioning(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseASMAggs(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType()
13395b89ad90SMark F. Adams M*/
1340b2573a8aSBarry Smith 
13415b89ad90SMark F. Adams #undef __FUNCT__
13425b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
13438cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
13445b89ad90SMark F. Adams {
13455b89ad90SMark F. Adams   PetscErrorCode ierr;
13465b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
13475b89ad90SMark F. Adams   PC_MG          *mg;
13480cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13499d5b6da9SMark F. Adams   static long count = 0;
13505ee9c036SSatish Balay #endif
13515b89ad90SMark F. Adams 
13525b89ad90SMark F. Adams   PetscFunctionBegin;
13535b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
13545b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
13555b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
13565b89ad90SMark F. Adams 
13575b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1358b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
13595b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
1360ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
13615b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
13625b89ad90SMark F. Adams 
1363b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
13641ab5ffc9SJed Brown 
13659d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
13669d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
13679d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
13689d5b6da9SMark F. Adams   pc_gamg->data    = 0;
13695b89ad90SMark F. Adams 
13709d5b6da9SMark F. Adams   /* register AMG type */
13713e3471ccSMark Adams   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
13729d5b6da9SMark F. Adams 
13739d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
13745b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
13755b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
13765b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
13775b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
13785b89ad90SMark F. Adams 
1379bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1380bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1381bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
13821cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1383bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1384bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1385bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1386c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1387bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
13889d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1389d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
1390ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1391038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
139225a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1393d3042614SMark Adams   pc_gamg->threshold        = 0.;
13949d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
13959d5b6da9SMark F. Adams   pc_gamg->verbose          = 0;
13969d5b6da9SMark F. Adams   pc_gamg->emax_id          = -1;
139757d29afaSToby Isaac   pc_gamg->firstCoarsen     = PETSC_FALSE;
13985e7c91beSJed Brown   pc_gamg->eigtarget[0]     = 0.05;
13995e7c91beSJed Brown   pc_gamg->eigtarget[1]     = 1.05;
1400c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
14019d5b6da9SMark F. Adams 
14020cbbd2e1SMark F. Adams   /* private events */
14030cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1404785cba28SMark F. Adams   if (count++ == 0) {
1405806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1406806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
14070cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
14080cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
14090cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1410806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1411806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1412806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1413806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1414806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1415806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1416806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1417806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1418806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1419806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1420806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1421806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1422f852f58cSMark F. Adams 
14230cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
14240cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
14250cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1426b4fbaa2aSMark F. Adams     /* create timer stages */
1427b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1428b4fbaa2aSMark F. Adams     {
1429b4fbaa2aSMark F. Adams       char     str[32];
1430b4fbaa2aSMark F. Adams       PetscInt lidx;
1431806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1432806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1433b4fbaa2aSMark F. Adams       for (lidx=1; lidx<9; lidx++) {
1434b4fbaa2aSMark F. Adams         sprintf(str,"MG Level %d",lidx);
1435806fa848SBarry Smith         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1436b4fbaa2aSMark F. Adams       }
1437b4fbaa2aSMark F. Adams     }
1438b4fbaa2aSMark F. Adams #endif
1439b4fbaa2aSMark F. Adams   }
1440b4fbaa2aSMark F. Adams #endif
1441bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1442bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
14435b89ad90SMark F. Adams   PetscFunctionReturn(0);
14445b89ad90SMark F. Adams }
14453e3471ccSMark Adams 
14463e3471ccSMark Adams #undef __FUNCT__
14473e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage"
14483e3471ccSMark Adams /*@C
14493e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
14503e3471ccSMark Adams  from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
14513e3471ccSMark Adams  when using static libraries.
14523e3471ccSMark Adams 
14533e3471ccSMark Adams  Level: developer
14543e3471ccSMark Adams 
14553e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
14563e3471ccSMark Adams  .seealso: PetscInitialize()
14573e3471ccSMark Adams @*/
14583e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
14593e3471ccSMark Adams {
14603e3471ccSMark Adams   PetscErrorCode ierr;
14613e3471ccSMark Adams 
14623e3471ccSMark Adams   PetscFunctionBegin;
14633e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
14643e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
14653e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
14663e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
14678e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
14683e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1469c1c463dbSMark Adams 
1470c1c463dbSMark Adams   /* general events */
1471c1c463dbSMark Adams #if defined PETSC_USE_LOG
1472c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1473c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1474c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1475c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1476c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1477c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1478c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1479c1c463dbSMark Adams #endif
1480c1c463dbSMark Adams 
14813e3471ccSMark Adams   PetscFunctionReturn(0);
14823e3471ccSMark Adams }
14833e3471ccSMark Adams 
14843e3471ccSMark Adams #undef __FUNCT__
14853e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage"
14863e3471ccSMark Adams /*@C
14873e3471ccSMark Adams  PCGAMGFinalizePackage - This function destroys everything in the PCGAMG package. It is
14883e3471ccSMark Adams  called from PetscFinalize().
14893e3471ccSMark Adams 
14903e3471ccSMark Adams  Level: developer
14913e3471ccSMark Adams 
14923e3471ccSMark Adams  .keywords: Petsc, destroy, package
14933e3471ccSMark Adams  .seealso: PetscFinalize()
14943e3471ccSMark Adams @*/
14953e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
14963e3471ccSMark Adams {
14973e3471ccSMark Adams   PetscErrorCode ierr;
14983e3471ccSMark Adams 
14993e3471ccSMark Adams   PetscFunctionBegin;
15003e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
15013e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
15023e3471ccSMark Adams   PetscFunctionReturn(0);
15033e3471ccSMark Adams }
1504a36cf38bSToby Isaac 
1505a36cf38bSToby Isaac #undef __FUNCT__
1506a36cf38bSToby Isaac #define __FUNCT__ "PCGAMGRegister"
1507a36cf38bSToby Isaac /*@C
1508a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1509a36cf38bSToby Isaac 
1510a36cf38bSToby Isaac  Input Parameters:
1511a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1512a36cf38bSToby Isaac  - create - function for creating the gamg context.
1513a36cf38bSToby Isaac 
1514a36cf38bSToby Isaac   Level: advanced
1515a36cf38bSToby Isaac 
1516a36cf38bSToby Isaac  .seealso: PCGAMGGetContext(), PCGAMGSetContext()
1517a36cf38bSToby Isaac @*/
1518a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1519a36cf38bSToby Isaac {
1520a36cf38bSToby Isaac   PetscErrorCode ierr;
1521a36cf38bSToby Isaac 
1522a36cf38bSToby Isaac   PetscFunctionBegin;
1523a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1524a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1525a36cf38bSToby Isaac   PetscFunctionReturn(0);
1526a36cf38bSToby Isaac }
1527a36cf38bSToby Isaac 
1528