xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision ef3f0257682f750e4ac8ac8ff986e80e09216830)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4af0996ceSBarry Smith #include <petsc/private/matimpl.h>
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
65b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */
7f96513f1SMatthew G Knepley 
80cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
90cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
10b4fbaa2aSMark F. Adams #endif
110cbbd2e1SMark F. Adams 
120cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
13fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG;
14fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO;
150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
19fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG;
200cbbd2e1SMark F. Adams #endif
210cbbd2e1SMark F. Adams 
22b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
230cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
24c1eae691SMark Adams static PetscLogStage gamg_stages[PETSC_GAMG_MAXLEVELS];
25b4fbaa2aSMark F. Adams #endif
26f96513f1SMatthew G Knepley 
27140e18c1SBarry Smith static PetscFunctionList GAMGList = 0;
283e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
299d5b6da9SMark F. Adams 
30d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
31d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
32d3d6bff4SMark F. Adams {
33d3d6bff4SMark F. Adams   PetscErrorCode ierr;
34d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
35d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
36d3d6bff4SMark F. Adams 
37d3d6bff4SMark F. Adams   PetscFunctionBegin;
3862294041SBarry Smith   if (pc_gamg->data) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen, cleaned up in SetUp\n");
391c1aac46SBarry Smith   pc_gamg->data_sz = 0;
40878e152fSMark F. Adams   ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
41a2f3521dSMark F. Adams   PetscFunctionReturn(0);
42a2f3521dSMark F. Adams }
43a2f3521dSMark F. Adams 
445b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
455b89ad90SMark F. Adams /*
46c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
47a147abb0SMark F. Adams      of active processors.
485b89ad90SMark F. Adams 
495b89ad90SMark F. Adams    Input Parameter:
50a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
51a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
529d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
53c5bfad50SMark F. Adams    . cr_bs - coarse block size
543530afc2SMark F. Adams    In/Output Parameter:
55a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
56afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
5711e60469SMark F. Adams    Output Parameter:
583530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
595b89ad90SMark F. Adams */
605cb416c2SMark F. Adams 
61171cca9aSMark Adams static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,IS * Pcolumnperm, PetscBool is_last)
625b89ad90SMark F. Adams {
63a2f3521dSMark F. Adams   PetscErrorCode  ierr;
649d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
65486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
66a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
673b4367a7SBarry Smith   MPI_Comm        comm;
68c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
693ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
705b89ad90SMark F. Adams 
715b89ad90SMark F. Adams   PetscFunctionBegin;
723b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
733b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
743b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
75c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
769d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
77038e3b61SMark F. Adams 
783ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
790298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
803ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
813ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
8273911c69SBarry Smith   } else {
833ae0bb68SMark Adams     PetscInt  bs;
843ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
853ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
863ae0bb68SMark Adams   }
87a2f3521dSMark F. Adams 
88c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
89171cca9aSMark Adams   if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) new_size = 1;
90171cca9aSMark Adams   else {
91472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
920298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
93a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
947f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
95c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
96a2f3521dSMark F. Adams   }
97f852f58cSMark F. Adams 
983cb8563fSToby Isaac   if (Pcolumnperm) *Pcolumnperm = NULL;
993cb8563fSToby Isaac 
100*ef3f0257SMark Adams   if (!pc_gamg->repart && new_size==nactive) {
101*ef3f0257SMark Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
102*ef3f0257SMark Adams     /* we know that the grid structure can be reused in MatPtAP */
103*ef3f0257SMark Adams   } else {
104*ef3f0257SMark Adams     /* we know that the grid structure can NOT be reused in MatPtAP */
1053ae0bb68SMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
106885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
107e33ef3b1SMark F. Adams 
10871959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
10971959b99SBarry 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);
1100cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1110cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
112b4fbaa2aSMark F. Adams #endif
113a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
114785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
115a90e85d9SMark Adams     if (pc_gamg->repart) {
116a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1175a9b9e01SMark F. Adams       Mat adj;
1185a9b9e01SMark F. Adams 
119cf8ae1d3SMark Adams      ierr = PetscInfo3(pc,"Repartition: size (active): %D --> %D, %D local equations\n",*a_nactive_proc,new_size,ncrs_eq);CHKERRQ(ierr);
1205a9b9e01SMark F. Adams 
121a2f3521dSMark F. Adams       /* get 'adj' */
122c5bfad50SMark F. Adams       if (cr_bs == 1) {
123038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
124806fa848SBarry Smith       } else {
125a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
126eb07cef2SMark F. Adams         Mat               tMat;
127a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
128b4fbaa2aSMark F. Adams         const PetscScalar *vals;
129b4fbaa2aSMark F. Adams         const PetscInt    *idx;
130a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1319057884aSMark F. Adams         static PetscInt   llev = 0;
132d9558ea9SBarry Smith         MatType           mtype;
133b4fbaa2aSMark F. Adams 
134e632b94dSBarry Smith         ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr);
135a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
136a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
137c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
13858471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
139c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
140c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
14158471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1423ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1433ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
14458471d46SMark F. Adams         }
1456876a03eSMark F. Adams 
146d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
1473b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1483ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
149d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
150a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
151a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
152e632b94dSBarry Smith         ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
153eb07cef2SMark F. Adams 
154a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
155c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
15622063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
157eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
158c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
159eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
160eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
161eb07cef2SMark F. Adams           }
16222063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
163eb07cef2SMark F. Adams         }
164eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166eb07cef2SMark F. Adams 
167b4fbaa2aSMark F. Adams         if (llev++ == -1) {
168b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
1698caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
1703b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
171b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
1723bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
173b4fbaa2aSMark F. Adams         }
174eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
175eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
176a2f3521dSMark F. Adams       } /* create 'adj' */
177f150b916SMark F. Adams 
178a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
1795a9b9e01SMark F. Adams         char            prefix[256];
1805a9b9e01SMark F. Adams         const char      *pcpre;
181b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
182b4fbaa2aSMark F. Adams         MatPartitioning mpart;
183a4b7d37bSMark F. Adams         IS              proc_is;
184a2f3521dSMark F. Adams         PetscInt        targetPE;
1852f03bc48SMark F. Adams 
1863b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
1875ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
1889d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
1898caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
19059a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
19111e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
192c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
193a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
19411e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1955a9b9e01SMark F. Adams 
1965ef31b24SMark F. Adams         /* collect IS info */
197785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
198a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
199a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
200c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
201a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
202c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
203a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
204eb07cef2SMark F. Adams           }
2055ef31b24SMark F. Adams         }
206a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
207a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2085ef31b24SMark F. Adams       }
2095ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2105a9b9e01SMark F. Adams 
2113b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2128263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
213806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
214a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
21562294041SBarry Smith 
2165a9b9e01SMark F. Adams       /* find factor */
217c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2185a9b9e01SMark F. Adams       else {
2195a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2205a9b9e01SMark F. Adams         jj = -1;
221c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
2222a82295dSMark Adams           if (!(size%kk)) { /* a candidate */
223c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2245a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2255a9b9e01SMark F. Adams             if (fact > best_fact) {
2265a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
2275a9b9e01SMark F. Adams             }
2285a9b9e01SMark F. Adams           }
2295a9b9e01SMark F. Adams         }
2305a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
231a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
2325a9b9e01SMark F. Adams       }
233c5df96a5SBarry Smith       new_size = size/rfactor;
2345a9b9e01SMark F. Adams 
235c5df96a5SBarry Smith       if (new_size==nactive) {
236a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2375a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
238302440fdSBarry Smith         ierr = PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr);
239fe682e87SBarry Smith #if defined PETSC_GAMG_USE_LOG
240fe682e87SBarry Smith         ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
241fe682e87SBarry Smith #endif
2425a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2435a9b9e01SMark F. Adams       }
2445a9b9e01SMark F. Adams 
245302440fdSBarry Smith       ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr);
246c5df96a5SBarry Smith       targetPE = rank/rfactor;
2473b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
248a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
249e33ef3b1SMark F. Adams 
25011e60469SMark F. Adams     /*
251a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
25211e60469SMark F. Adams      */
253a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2547700e67bSMark Adams     is_eq_num_prim = is_eq_num;
25511e60469SMark F. Adams     /*
256a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
25711e60469SMark F. Adams      */
258c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
259c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
260a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
2613ae0bb68SMark Adams     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */
262a2f3521dSMark F. Adams 
263a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2640cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2650cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
266b4fbaa2aSMark F. Adams #endif
267885364a3SMark 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 */
268885364a3SMark Adams     {
269885364a3SMark Adams     Vec            src_crd, dest_crd;
270885364a3SMark 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;
271885364a3SMark Adams     VecScatter     vecscat;
272885364a3SMark Adams     PetscScalar    *array;
273885364a3SMark Adams     IS isscat;
274a2f3521dSMark F. Adams 
275a2f3521dSMark F. Adams     /* move data (for primal equations only) */
27622063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
2773b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
2783ae0bb68SMark Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
279c0dedaeaSBarry Smith     ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
28011e60469SMark F. Adams     /*
2819d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
282c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
28311e60469SMark F. Adams      */
284854ce69bSBarry Smith     ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
285a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2863ae0bb68SMark Adams     for (ii=0,jj=0; ii<ncrs; ii++) {
287c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
288a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
28911e60469SMark F. Adams     }
290a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2913ae0bb68SMark Adams     ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
29292a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
29311e60469SMark F. Adams     /*
29411e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
29511e60469SMark F. Adams      */
2963ae0bb68SMark Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
2979d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
2983ae0bb68SMark Adams       const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
2993ae0bb68SMark Adams       for (ii=0; ii<ncrs; ii++) {
3009d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
301a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
302c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
303676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
304d3d6bff4SMark F. Adams         }
305038e3b61SMark F. Adams       }
306eb07cef2SMark F. Adams     }
307eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
308eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
30911e60469SMark F. Adams     /*
31011e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
31111e60469SMark F. Adams       to the correct processor
31211e60469SMark F. Adams     */
3130298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
31411e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
31511e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31611e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31711e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
31811e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
31911e60469SMark F. Adams     /*
32011e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
32111e60469SMark F. Adams     */
322c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
323578f55a3SPeter Brune     ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3242fa5cd67SKarl Rupp 
3253ae0bb68SMark Adams     pc_gamg->data_sz = node_data_sz*ncrs_new;
3263ae0bb68SMark Adams     strideNew        = ncrs_new*ndata_rows;
3272fa5cd67SKarl Rupp 
32811e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3299d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3303ae0bb68SMark Adams       for (ii=0; ii<ncrs_new; ii++) {
3319d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
332a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
333c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
334d3d6bff4SMark F. Adams         }
335038e3b61SMark F. Adams       }
336038e3b61SMark F. Adams     }
33711e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
33811e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
339885364a3SMark Adams     }
340a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3410cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3420cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
343ed3f9983SMark F. Adams #endif
344a2f3521dSMark F. Adams 
34511e60469SMark F. Adams     /*
3467dae84e0SHong Zhang       Invert for MatCreateSubMatrix
34711e60469SMark F. Adams     */
348a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
349a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
350c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
351a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
352a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
353a2f3521dSMark F. Adams     }
3543cb8563fSToby Isaac     if (Pcolumnperm) {
3553cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3563cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3573cb8563fSToby Isaac     }
358a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3590cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3600cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3610cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
362ed3f9983SMark F. Adams #endif
363a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
364a2f3521dSMark F. Adams     {
365a2f3521dSMark F. Adams       Mat mat;
3667dae84e0SHong Zhang       ierr        = MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
367a2f3521dSMark F. Adams       *a_Amat_crs = mat;
368a2f3521dSMark F. Adams     }
369038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
370a2f3521dSMark F. Adams 
3710cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3720cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
373ed3f9983SMark F. Adams #endif
37411e60469SMark F. Adams     /* prolongator */
37511e60469SMark F. Adams     {
37611e60469SMark F. Adams       IS       findices;
377a2f3521dSMark F. Adams       PetscInt Istart,Iend;
378a2f3521dSMark F. Adams       Mat      Pnew;
37962294041SBarry Smith 
380a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
3810cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3820cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
383ed3f9983SMark F. Adams #endif
3843b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
385c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
3867dae84e0SHong Zhang       ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
38711e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
388c5bfad50SMark F. Adams 
3890cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3900cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
391ed3f9983SMark F. Adams #endif
3923530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
393a2f3521dSMark F. Adams 
394a2f3521dSMark F. Adams       /* output - repartitioned */
395a2f3521dSMark F. Adams       *a_P_inout = Pnew;
396e33ef3b1SMark F. Adams     }
397a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
3985b89ad90SMark F. Adams 
399c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
400a2f3521dSMark F. Adams   }
4015b89ad90SMark F. Adams   PetscFunctionReturn(0);
4025b89ad90SMark F. Adams }
4035b89ad90SMark F. Adams 
4045b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4055b89ad90SMark F. Adams /*
4065b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4075b89ad90SMark F. Adams                     by setting data structures and options.
4085b89ad90SMark F. Adams 
4095b89ad90SMark F. Adams    Input Parameter:
4105b89ad90SMark F. Adams .  pc - the preconditioner context
4115b89ad90SMark F. Adams 
4125b89ad90SMark F. Adams */
4139d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4145b89ad90SMark F. Adams {
4155b89ad90SMark F. Adams   PetscErrorCode ierr;
4169d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4175b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4182adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
419c1eae691SMark Adams   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_GAMG_MAXLEVELS];
4203b4367a7SBarry Smith   MPI_Comm       comm;
421c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
422c1eae691SMark Adams   Mat            Aarr[PETSC_GAMG_MAXLEVELS],Parr[PETSC_GAMG_MAXLEVELS];
423c1eae691SMark Adams   IS             *ASMLocalIDsArr[PETSC_GAMG_MAXLEVELS];
424a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
425569f4572SMark Adams   MatInfo        info;
426171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
4275ef31b24SMark F. Adams 
4285b89ad90SMark F. Adams   PetscFunctionBegin;
4293b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4303b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4313b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
432dfd5c07aSMark F. Adams 
43384d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
4341c1aac46SBarry Smith     if ((PetscBool)(!pc_gamg->reuse_prol)) {
435878e152fSMark F. Adams       /* reset everything */
436878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
437878e152fSMark F. Adams       pc->setupcalled = 0;
438806fa848SBarry Smith     } else {
43984d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
44003a628feSMark F. Adams       /* just do Galerkin grids */
44158471d46SMark F. Adams       Mat          B,dA,dB;
44258471d46SMark F. Adams 
44371959b99SBarry Smith       if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
4449d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
44558471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
44623ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
44758471d46SMark F. Adams         /* (re)set to get dirty flag */
44823ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
44958471d46SMark F. Adams 
4502fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
451*ef3f0257SMark Adams           /* 2nd solve, matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
452*ef3f0257SMark Adams           if (pc_gamg->setup_count==2 /* && pc_gamg->repart||reduction */) {
453*ef3f0257SMark Adams             ierr = PetscInfo2(pc,"new RAP after first solve level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr);
454*ef3f0257SMark Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,2.0,&B);CHKERRQ(ierr);
455084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
45603a628feSMark F. Adams             mglevels[level]->A = B;
457806fa848SBarry Smith           } else {
458*ef3f0257SMark Adams             ierr = PetscInfo2(pc,"RAP after first solve reusing matrix level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr);
45923ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
46058471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
46103a628feSMark F. Adams           }
46223ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
46358471d46SMark F. Adams           dB   = B;
46458471d46SMark F. Adams         }
4655f8cf99dSMark F. Adams       }
466d5280255SMark F. Adams 
467d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
46858471d46SMark F. Adams       PetscFunctionReturn(0);
469eb07cef2SMark F. Adams     }
470878e152fSMark F. Adams   }
471f6536408SMark F. Adams 
472878e152fSMark F. Adams   if (!pc_gamg->data) {
473878e152fSMark F. Adams     if (pc_gamg->orig_data) {
474878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
4750298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
4762fa5cd67SKarl Rupp 
477878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
478878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
479878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
4802fa5cd67SKarl Rupp 
481785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
482878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
483806fa848SBarry Smith     } else {
4841ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
4857700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
4869d5b6da9SMark F. Adams     }
487878e152fSMark F. Adams   }
488878e152fSMark F. Adams 
489878e152fSMark F. Adams   /* cache original data for reuse */
4901c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
491785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
492878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
493878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
494878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
495878e152fSMark F. Adams   }
496038e3b61SMark F. Adams 
497302f38e8SMark F. Adams   /* get basic dims */
498302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
499171cca9aSMark Adams   ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr);
50084d3f75bSMark F. Adams 
501569f4572SMark Adams   ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
502569f4572SMark Adams   nnz0   = info.nz_used;
503569f4572SMark Adams   nnztot = info.nz_used;
50462294041SBarry Smith   ierr = PetscInfo6(pc,"level %d) N=%D, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(int)(nnz0/(PetscReal)M+0.5),size);CHKERRQ(ierr);
505569f4572SMark Adams 
506a2f3521dSMark F. Adams   /* Get A_i and R_i */
50762294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
5089ab59c8bSMark Adams     pc_gamg->current_level = level;
5095b89ad90SMark F. Adams     level1 = level + 1;
5100cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5110cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
512a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
513a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
514b4fbaa2aSMark F. Adams #endif
515a2f3521dSMark F. Adams #endif
516c8b0795cSMark F. Adams     { /* construct prolongator */
517725b86d8SJed Brown       Mat              Gmat;
5180cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5197700e67bSMark Adams       Mat              Prol11;
520c8b0795cSMark F. Adams 
5217700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5221ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5237700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
524c8b0795cSMark F. Adams 
525a2f3521dSMark F. Adams       /* could have failed to create new level */
526a2f3521dSMark F. Adams       if (Prol11) {
5279d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5280298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
529a2f3521dSMark F. Adams 
530fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
531c8b0795cSMark F. Adams           /* smooth */
532fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
533c8b0795cSMark F. Adams         }
534c8b0795cSMark F. Adams 
5357700e67bSMark Adams         Parr[level1] = Prol11;
536171cca9aSMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
537ffc955d6SMark F. Adams 
5380c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
5391b18a24aSMark Adams         PetscInt bs;
5401b18a24aSMark Adams         ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
5410a3c815dSMark Adams         ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
542ffc955d6SMark F. Adams       }
543ffc955d6SMark F. Adams 
544a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
54541b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
546a2f3521dSMark F. Adams     } /* construct prolongator scope */
5470cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5480cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
549c8b0795cSMark F. Adams #endif
5507f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
551171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
552569f4572SMark Adams       ierr =  PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr);
55362294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
554a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
555a90e85d9SMark Adams #endif
556c8b0795cSMark F. Adams       break;
557c8b0795cSMark F. Adams     }
5580cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5590cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
560b4fbaa2aSMark F. Adams #endif
561171cca9aSMark Adams     ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */
562171cca9aSMark Adams     if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ????????");
563171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
5640e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
565171cca9aSMark Adams     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr);
566a2f3521dSMark F. Adams 
5670cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5680cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
569b4fbaa2aSMark F. Adams #endif
570171cca9aSMark Adams     ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */
571569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
572569f4572SMark Adams     nnztot += info.nz_used;
5731d5b2942SMark Adams     ierr = PetscInfo5(pc,"%d) N=%D, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",level1,M,pc_gamg->data_cell_cols,(int)(info.nz_used/(PetscReal)M),nactivepe);CHKERRQ(ierr);
574569f4572SMark Adams 
5750cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
576b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
577b4fbaa2aSMark F. Adams #endif
578a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
5799ab59c8bSMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) {
5809ab59c8bSMark Adams       ierr =  PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr);
581a90e85d9SMark Adams       level++;
582a90e85d9SMark Adams       break;
583a90e85d9SMark Adams     }
584c8b0795cSMark F. Adams   } /* levels */
585c8b0795cSMark F. Adams   ierr                  = PetscFree(pc_gamg->data);CHKERRQ(ierr);
586c8b0795cSMark F. Adams 
587569f4572SMark Adams   ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr);
5889d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
5895b89ad90SMark F. Adams   fine_level       = level;
5900298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
5915b89ad90SMark F. Adams 
59262294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
593d5280255SMark F. Adams     /* set default smoothers & set operators */
59462294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
595ffc955d6SMark F. Adams       KSP smoother;
596ffc955d6SMark F. Adams       PC  subpc;
597a2f3521dSMark F. Adams 
5989d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
599f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
600ffc955d6SMark F. Adams 
601a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
602a2f3521dSMark F. Adams       /* set ops */
60323ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
604a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
605a2f3521dSMark F. Adams 
606a2f3521dSMark F. Adams       /* set defaults */
6076c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
608a2f3521dSMark F. Adams 
6090c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
6100c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
6112d3561bbSSatish Balay         PetscInt sz;
6127a28f3e5SMark Adams         IS       *iss;
613a2f3521dSMark F. Adams 
6142d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
6157a28f3e5SMark Adams         iss   = ASMLocalIDsArr[level];
6160c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
6170a3c815dSMark Adams         ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr);
6180c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
6197f66b68fSMark Adams         if (!sz) {
620ffc955d6SMark F. Adams           IS       is;
6210a3c815dSMark Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
6227a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
623a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
624806fa848SBarry Smith         } else {
625a94c3b12SMark F. Adams           PetscInt kk;
6267a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr);
627a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
6287a28f3e5SMark Adams             ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr);
629a94c3b12SMark F. Adams           }
6307a28f3e5SMark Adams           ierr = PetscFree(iss);CHKERRQ(ierr);
631ffc955d6SMark F. Adams         }
6320298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
633ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
634806fa848SBarry Smith       } else {
635890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
636ffc955d6SMark F. Adams       }
637d5280255SMark F. Adams     }
638d5280255SMark F. Adams     {
639d5280255SMark F. Adams       /* coarse grid */
640d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
641d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
642d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
64323ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
644cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
645d5280255SMark F. Adams         ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
646d5280255SMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
647d5280255SMark F. Adams         ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
648d5280255SMark F. Adams         ierr = PCSetUp(subpc);CHKERRQ(ierr);
64971959b99SBarry Smith         ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
65071959b99SBarry Smith         if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
651d5280255SMark F. Adams         ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
652d5280255SMark F. Adams         ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
6539dbfc187SHong Zhang         ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
6542fb0b348SMark F. Adams         ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
65508e36f19SMark Adams         ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
6565b42dca8SJed Brown         /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
6575b42dca8SJed Brown          * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
6585b42dca8SJed Brown          * view every subdomain as though they were different. */
6595b42dca8SJed Brown         ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
660d5280255SMark F. Adams       }
661cf8ae1d3SMark Adams     }
662d5280255SMark F. Adams 
663d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
664d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
665e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
666d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
66769aca0b8SBarry Smith     ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
668d5280255SMark F. Adams 
669d5280255SMark F. Adams     /* clean up */
670d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
671587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
672587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
6735b89ad90SMark F. Adams     }
6740cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
675806fa848SBarry Smith   } else {
6765f8cf99dSMark F. Adams     KSP smoother;
677302440fdSBarry Smith     ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
6789d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
67923ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
6805f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
6819d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
6825f8cf99dSMark F. Adams   }
6835b89ad90SMark F. Adams   PetscFunctionReturn(0);
6845b89ad90SMark F. Adams }
6855b89ad90SMark F. Adams 
686eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
6875b89ad90SMark F. Adams /*
6885b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
6895b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
6905b89ad90SMark F. Adams 
6915b89ad90SMark F. Adams    Input Parameter:
6925b89ad90SMark F. Adams .  pc - the preconditioner context
6935b89ad90SMark F. Adams 
6945b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
6955b89ad90SMark F. Adams */
6965b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
6975b89ad90SMark F. Adams {
6985b89ad90SMark F. Adams   PetscErrorCode ierr;
6995b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
7005b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
7015b89ad90SMark F. Adams 
7025b89ad90SMark F. Adams   PetscFunctionBegin;
7035b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7049b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
7059b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
7069b8ffb57SJed Brown   }
7071ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
7081ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
7095b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7105b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7115b89ad90SMark F. Adams   PetscFunctionReturn(0);
7125b89ad90SMark F. Adams }
7135b89ad90SMark F. Adams 
714676e1743SMark F. Adams /*@
715cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
716676e1743SMark F. Adams 
7171cc46a46SBarry Smith    Logically Collective on PC
718676e1743SMark F. Adams 
719676e1743SMark F. Adams    Input Parameters:
7201cc46a46SBarry Smith +  pc - the preconditioner context
7211cc46a46SBarry Smith -  n - the number of equations
722676e1743SMark F. Adams 
723676e1743SMark F. Adams 
724676e1743SMark F. Adams    Options Database Key:
7251cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
726676e1743SMark F. Adams 
72795452b02SPatrick Sanan    Notes:
72895452b02SPatrick Sanan     GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
729cab9ed1eSBarry Smith           that has degrees of freedom
730cab9ed1eSBarry Smith 
731676e1743SMark F. Adams    Level: intermediate
732676e1743SMark F. Adams 
7331c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
734676e1743SMark F. Adams 
7351c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim()
736676e1743SMark F. Adams @*/
737676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
738676e1743SMark F. Adams {
739676e1743SMark F. Adams   PetscErrorCode ierr;
740676e1743SMark F. Adams 
741676e1743SMark F. Adams   PetscFunctionBegin;
742676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
743676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
744676e1743SMark F. Adams   PetscFunctionReturn(0);
745676e1743SMark F. Adams }
746676e1743SMark F. Adams 
7471e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
748676e1743SMark F. Adams {
749c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
750c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
751676e1743SMark F. Adams 
752676e1743SMark F. Adams   PetscFunctionBegin;
7539d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
754676e1743SMark F. Adams   PetscFunctionReturn(0);
755676e1743SMark F. Adams }
756676e1743SMark F. Adams 
757389730f3SMark F. Adams /*@
758cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
759389730f3SMark F. Adams 
760389730f3SMark F. Adams  Collective on PC
761389730f3SMark F. Adams 
762389730f3SMark F. Adams    Input Parameters:
7631cc46a46SBarry Smith +  pc - the preconditioner context
7641cc46a46SBarry Smith -  n - maximum number of equations to aim for
765389730f3SMark F. Adams 
766389730f3SMark F. Adams    Options Database Key:
7671cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
768389730f3SMark F. Adams 
769389730f3SMark F. Adams    Level: intermediate
770389730f3SMark F. Adams 
7711c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
772389730f3SMark F. Adams 
7731c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim()
774389730f3SMark F. Adams @*/
775389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
776389730f3SMark F. Adams {
777389730f3SMark F. Adams   PetscErrorCode ierr;
778389730f3SMark F. Adams 
779389730f3SMark F. Adams   PetscFunctionBegin;
780389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
781389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
782389730f3SMark F. Adams   PetscFunctionReturn(0);
783389730f3SMark F. Adams }
784389730f3SMark F. Adams 
7851e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
786389730f3SMark F. Adams {
787389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
788389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
789389730f3SMark F. Adams 
790389730f3SMark F. Adams   PetscFunctionBegin;
7919d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
792389730f3SMark F. Adams   PetscFunctionReturn(0);
793389730f3SMark F. Adams }
794389730f3SMark F. Adams 
795676e1743SMark F. Adams /*@
796cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
797676e1743SMark F. Adams 
798676e1743SMark F. Adams    Collective on PC
799676e1743SMark F. Adams 
800676e1743SMark F. Adams    Input Parameters:
8011cc46a46SBarry Smith +  pc - the preconditioner context
8021cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
803676e1743SMark F. Adams 
804676e1743SMark F. Adams    Options Database Key:
8051cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
806676e1743SMark F. Adams 
80795452b02SPatrick Sanan    Notes:
80895452b02SPatrick Sanan     this will generally improve the loading balancing of the work on each level
809cab9ed1eSBarry Smith 
810676e1743SMark F. Adams    Level: intermediate
811676e1743SMark F. Adams 
8121c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
813676e1743SMark F. Adams 
814676e1743SMark F. Adams .seealso: ()
815676e1743SMark F. Adams @*/
816cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
817676e1743SMark F. Adams {
818676e1743SMark F. Adams   PetscErrorCode ierr;
819676e1743SMark F. Adams 
820676e1743SMark F. Adams   PetscFunctionBegin;
821676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
822cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
823676e1743SMark F. Adams   PetscFunctionReturn(0);
824676e1743SMark F. Adams }
825676e1743SMark F. Adams 
826cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
827676e1743SMark F. Adams {
828c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
829c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
830676e1743SMark F. Adams 
831676e1743SMark F. Adams   PetscFunctionBegin;
8329d5b6da9SMark F. Adams   pc_gamg->repart = n;
833676e1743SMark F. Adams   PetscFunctionReturn(0);
834676e1743SMark F. Adams }
835676e1743SMark F. Adams 
836dfd5c07aSMark F. Adams /*@
837cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
838dfd5c07aSMark F. Adams 
839dfd5c07aSMark F. Adams    Collective on PC
840dfd5c07aSMark F. Adams 
841dfd5c07aSMark F. Adams    Input Parameters:
8421cc46a46SBarry Smith +  pc - the preconditioner context
8431cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
844dfd5c07aSMark F. Adams 
845dfd5c07aSMark F. Adams    Options Database Key:
8461cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
847dfd5c07aSMark F. Adams 
848dfd5c07aSMark F. Adams    Level: intermediate
849dfd5c07aSMark F. Adams 
85095452b02SPatrick Sanan    Notes:
85195452b02SPatrick Sanan     this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
852cab9ed1eSBarry Smith           rebuilding the preconditioner quicker.
853cab9ed1eSBarry Smith 
8541c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
855dfd5c07aSMark F. Adams 
856dfd5c07aSMark F. Adams .seealso: ()
857dfd5c07aSMark F. Adams @*/
8581cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
859dfd5c07aSMark F. Adams {
860dfd5c07aSMark F. Adams   PetscErrorCode ierr;
861dfd5c07aSMark F. Adams 
862dfd5c07aSMark F. Adams   PetscFunctionBegin;
863dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8641cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
865dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
866dfd5c07aSMark F. Adams }
867dfd5c07aSMark F. Adams 
8681cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
869dfd5c07aSMark F. Adams {
870dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
871dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
872dfd5c07aSMark F. Adams 
873dfd5c07aSMark F. Adams   PetscFunctionBegin;
874dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
875dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
876dfd5c07aSMark F. Adams }
877dfd5c07aSMark F. Adams 
878ffc955d6SMark F. Adams /*@
879cab9ed1eSBarry Smith    PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner.
880ffc955d6SMark F. Adams 
881ffc955d6SMark F. Adams    Collective on PC
882ffc955d6SMark F. Adams 
883ffc955d6SMark F. Adams    Input Parameters:
884cab9ed1eSBarry Smith +  pc - the preconditioner context
885cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
886ffc955d6SMark F. Adams 
887ffc955d6SMark F. Adams    Options Database Key:
888cab9ed1eSBarry Smith .  -pc_gamg_asm_use_agg
889ffc955d6SMark F. Adams 
890ffc955d6SMark F. Adams    Level: intermediate
891ffc955d6SMark F. Adams 
8921c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
893ffc955d6SMark F. Adams 
894ffc955d6SMark F. Adams .seealso: ()
895ffc955d6SMark F. Adams @*/
896cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
897ffc955d6SMark F. Adams {
898ffc955d6SMark F. Adams   PetscErrorCode ierr;
899ffc955d6SMark F. Adams 
900ffc955d6SMark F. Adams   PetscFunctionBegin;
901ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
902cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
903ffc955d6SMark F. Adams   PetscFunctionReturn(0);
904ffc955d6SMark F. Adams }
905ffc955d6SMark F. Adams 
906cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
907ffc955d6SMark F. Adams {
908ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
909ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
910ffc955d6SMark F. Adams 
911ffc955d6SMark F. Adams   PetscFunctionBegin;
912cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
913ffc955d6SMark F. Adams   PetscFunctionReturn(0);
914ffc955d6SMark F. Adams }
915ffc955d6SMark F. Adams 
916171cca9aSMark Adams /*@
917cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
918171cca9aSMark Adams 
919171cca9aSMark Adams    Collective on PC
920171cca9aSMark Adams 
921171cca9aSMark Adams    Input Parameters:
922171cca9aSMark Adams +  pc - the preconditioner context
923cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
924171cca9aSMark Adams 
925171cca9aSMark Adams    Options Database Key:
926cf8ae1d3SMark Adams .  -pc_gamg_use_parallel_coarse_grid_solver
927171cca9aSMark Adams 
928171cca9aSMark Adams    Level: intermediate
929171cca9aSMark Adams 
930171cca9aSMark Adams    Concepts: Unstructured multigrid preconditioner
931171cca9aSMark Adams 
932171cca9aSMark Adams .seealso: ()
933171cca9aSMark Adams @*/
934171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
935171cca9aSMark Adams {
936171cca9aSMark Adams   PetscErrorCode ierr;
937171cca9aSMark Adams 
938171cca9aSMark Adams   PetscFunctionBegin;
939171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
940171cca9aSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
941171cca9aSMark Adams   PetscFunctionReturn(0);
942171cca9aSMark Adams }
943171cca9aSMark Adams 
944171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
945171cca9aSMark Adams {
946171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
947171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
948171cca9aSMark Adams 
949171cca9aSMark Adams   PetscFunctionBegin;
950171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
951ffc955d6SMark F. Adams   PetscFunctionReturn(0);
952ffc955d6SMark F. Adams }
953ffc955d6SMark F. Adams 
9544ef23d27SMark F. Adams /*@
9551cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
9564ef23d27SMark F. Adams 
9574ef23d27SMark F. Adams    Not collective on PC
9584ef23d27SMark F. Adams 
9594ef23d27SMark F. Adams    Input Parameters:
9601cc46a46SBarry Smith +  pc - the preconditioner
9611cc46a46SBarry Smith -  n - the maximum number of levels to use
9624ef23d27SMark F. Adams 
9634ef23d27SMark F. Adams    Options Database Key:
9644ef23d27SMark F. Adams .  -pc_mg_levels
9654ef23d27SMark F. Adams 
9664ef23d27SMark F. Adams    Level: intermediate
9674ef23d27SMark F. Adams 
9681c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
9694ef23d27SMark F. Adams 
9704ef23d27SMark F. Adams .seealso: ()
9714ef23d27SMark F. Adams @*/
9724ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
9734ef23d27SMark F. Adams {
9744ef23d27SMark F. Adams   PetscErrorCode ierr;
9754ef23d27SMark F. Adams 
9764ef23d27SMark F. Adams   PetscFunctionBegin;
9774ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9784ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
9794ef23d27SMark F. Adams   PetscFunctionReturn(0);
9804ef23d27SMark F. Adams }
9814ef23d27SMark F. Adams 
9821e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
9834ef23d27SMark F. Adams {
9844ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
9854ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
9864ef23d27SMark F. Adams 
9874ef23d27SMark F. Adams   PetscFunctionBegin;
9889d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
9894ef23d27SMark F. Adams   PetscFunctionReturn(0);
9904ef23d27SMark F. Adams }
9914ef23d27SMark F. Adams 
9923542efc5SMark F. Adams /*@
9933542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
9943542efc5SMark F. Adams 
9953542efc5SMark F. Adams    Not collective on PC
9963542efc5SMark F. Adams 
9973542efc5SMark F. Adams    Input Parameters:
9981cc46a46SBarry Smith +  pc - the preconditioner context
999b001cb0fSBarry Smith -  threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
10003542efc5SMark F. Adams 
10013542efc5SMark F. Adams    Options Database Key:
10021cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
10033542efc5SMark F. Adams 
100495452b02SPatrick Sanan    Notes:
100595452b02SPatrick Sanan     Before aggregating the graph GAMG will remove small values from the graph thus reducing the coupling in the graph and a different
1006cab9ed1eSBarry Smith     (perhaps better) coarser set of points.
1007cab9ed1eSBarry Smith 
10083542efc5SMark F. Adams    Level: intermediate
10093542efc5SMark F. Adams 
10101c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
10113542efc5SMark F. Adams 
1012a37438d7SBarry Smith .seealso: PCGAMGFilterGraph()
10133542efc5SMark F. Adams @*/
1014c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
10153542efc5SMark F. Adams {
10163542efc5SMark F. Adams   PetscErrorCode ierr;
10173542efc5SMark F. Adams 
10183542efc5SMark F. Adams   PetscFunctionBegin;
10193542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1020c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr);
10213542efc5SMark F. Adams   PetscFunctionReturn(0);
10223542efc5SMark F. Adams }
10233542efc5SMark F. Adams 
1024c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
10253542efc5SMark F. Adams {
1026c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1027c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1028c1eae691SMark Adams   PetscInt i;
1029c1eae691SMark Adams   PetscFunctionBegin;
1030c1eae691SMark Adams   for (i=0;i<n;i++) pc_gamg->threshold[i] = v[i];
1031c1eae691SMark Adams   do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1032c1eae691SMark Adams   PetscFunctionReturn(0);
1033c1eae691SMark Adams }
1034c1eae691SMark Adams 
1035c1eae691SMark Adams /*@
1036c1eae691SMark Adams    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1037c1eae691SMark Adams 
1038c1eae691SMark Adams    Not collective on PC
1039c1eae691SMark Adams 
1040c1eae691SMark Adams    Input Parameters:
1041c1eae691SMark Adams +  pc - the preconditioner context
1042c1eae691SMark Adams -  scale - the threshold value reduction, ussually < 1.0
1043c1eae691SMark Adams 
1044c1eae691SMark Adams    Options Database Key:
1045c1eae691SMark Adams .  -pc_gamg_threshold_scale <v>
1046c1eae691SMark Adams 
1047c1eae691SMark Adams    Level: advanced
1048c1eae691SMark Adams 
1049c1eae691SMark Adams    Concepts: Unstructured multigrid preconditioner
1050c1eae691SMark Adams 
1051c1eae691SMark Adams .seealso: ()
1052c1eae691SMark Adams @*/
1053c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1054c1eae691SMark Adams {
1055c1eae691SMark Adams   PetscErrorCode ierr;
10563542efc5SMark F. Adams 
10573542efc5SMark F. Adams   PetscFunctionBegin;
1058c1eae691SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1059c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr);
1060c1eae691SMark Adams   PetscFunctionReturn(0);
1061c1eae691SMark Adams }
1062c1eae691SMark Adams 
1063c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1064c1eae691SMark Adams {
1065c1eae691SMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1066c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1067c1eae691SMark Adams   PetscFunctionBegin;
1068c1eae691SMark Adams   pc_gamg->threshold_scale = v;
10693542efc5SMark F. Adams   PetscFunctionReturn(0);
10703542efc5SMark F. Adams }
10713542efc5SMark F. Adams 
1072e20c40e8SBarry Smith /*@C
1073c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1074676e1743SMark F. Adams 
1075676e1743SMark F. Adams    Collective on PC
1076676e1743SMark F. Adams 
1077676e1743SMark F. Adams    Input Parameters:
1078c60c7ad4SBarry Smith +  pc - the preconditioner context
1079c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1080676e1743SMark F. Adams 
1081676e1743SMark F. Adams    Options Database Key:
1082cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1083676e1743SMark F. Adams 
1084676e1743SMark F. Adams    Level: intermediate
1085676e1743SMark F. Adams 
10861c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1087676e1743SMark F. Adams 
1088cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1089676e1743SMark F. Adams @*/
109019fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1091676e1743SMark F. Adams {
1092676e1743SMark F. Adams   PetscErrorCode ierr;
1093676e1743SMark F. Adams 
1094676e1743SMark F. Adams   PetscFunctionBegin;
1095676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1096806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1097676e1743SMark F. Adams   PetscFunctionReturn(0);
1098676e1743SMark F. Adams }
1099676e1743SMark F. Adams 
1100e20c40e8SBarry Smith /*@C
1101c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1102c60c7ad4SBarry Smith 
1103c60c7ad4SBarry Smith    Collective on PC
1104c60c7ad4SBarry Smith 
1105c60c7ad4SBarry Smith    Input Parameter:
1106c60c7ad4SBarry Smith .  pc - the preconditioner context
1107c60c7ad4SBarry Smith 
1108c60c7ad4SBarry Smith    Output Parameter:
1109c60c7ad4SBarry Smith .  type - the type of algorithm used
1110c60c7ad4SBarry Smith 
1111c60c7ad4SBarry Smith    Level: intermediate
1112c60c7ad4SBarry Smith 
11131c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1114c60c7ad4SBarry Smith 
11151c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1116c60c7ad4SBarry Smith @*/
1117c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1118c60c7ad4SBarry Smith {
1119c60c7ad4SBarry Smith   PetscErrorCode ierr;
1120c60c7ad4SBarry Smith 
1121c60c7ad4SBarry Smith   PetscFunctionBegin;
1122c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1123c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1124c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1125c60c7ad4SBarry Smith }
1126c60c7ad4SBarry Smith 
1127c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1128c60c7ad4SBarry Smith {
1129c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1130c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1131c60c7ad4SBarry Smith 
1132c60c7ad4SBarry Smith   PetscFunctionBegin;
1133c60c7ad4SBarry Smith   *type = pc_gamg->type;
1134c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1135c60c7ad4SBarry Smith }
1136c60c7ad4SBarry Smith 
11371e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1138676e1743SMark F. Adams {
11399d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
11401ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
11411ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1142676e1743SMark F. Adams 
1143676e1743SMark F. Adams   PetscFunctionBegin;
1144c60c7ad4SBarry Smith   pc_gamg->type = type;
11451c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
11469d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
11471ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
11481ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
11491ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1150e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
11513ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
11523ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
11533ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
11543ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
11553ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
11563ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
11573ae0bb68SMark Adams     pc_gamg->data_sz = 0;
11581ab5ffc9SJed Brown   }
11591ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
11601ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
11619d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1162676e1743SMark F. Adams   PetscFunctionReturn(0);
1163676e1743SMark F. Adams }
1164676e1743SMark F. Adams 
11655adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
11665adeb434SBarry Smith {
1167c1eae691SMark Adams   PetscErrorCode ierr,i;
11685adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
11695adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
11705adeb434SBarry Smith 
11715adeb434SBarry Smith   PetscFunctionBegin;
11725adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1173459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");CHKERRQ(ierr);
1174c1eae691SMark Adams   for (i=0;i<pc_gamg->current_level;i++) {
1175c1eae691SMark Adams     ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr);
1176c1eae691SMark Adams   }
1177459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1178459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);CHKERRQ(ierr);
1179cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
1180cab9ed1eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr);
1181cab9ed1eSBarry Smith   }
1182171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
1183171cca9aSMark Adams     ierr = PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr);
1184171cca9aSMark Adams   }
11855adeb434SBarry Smith   if (pc_gamg->ops->view) {
11865adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
11875adeb434SBarry Smith   }
11885adeb434SBarry Smith   PetscFunctionReturn(0);
11895adeb434SBarry Smith }
11905adeb434SBarry Smith 
11914416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
11925b89ad90SMark F. Adams {
1193676e1743SMark F. Adams   PetscErrorCode ierr;
1194676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1195676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1196676e1743SMark F. Adams   PetscBool      flag;
11973b4367a7SBarry Smith   MPI_Comm       comm;
119814a9496bSBarry Smith   char           prefix[256];
1199c1eae691SMark Adams   PetscInt       i,n;
120014a9496bSBarry Smith   const char     *pcpre;
12015b89ad90SMark F. Adams 
12025b89ad90SMark F. Adams   PetscFunctionBegin;
12033b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1204e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1205676e1743SMark F. Adams   {
1206bd94a7aaSJed Brown     char tname[256];
12071a1c1e04SBarry Smith     ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1208bd94a7aaSJed Brown     if (flag) {
1209bd94a7aaSJed Brown       ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
12101ab5ffc9SJed Brown     }
1211cab9ed1eSBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
12121cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1213a303c832SJed Brown     ierr = PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation aggregates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL);CHKERRQ(ierr);
1214cf8ae1d3SMark Adams     ierr = PetscOptionsBool("-pc_gamg_use_parallel_coarse_grid_solver","Use parallel coarse grid solver (otherwise put last grid on one process)","PCGAMGSetUseParallelCoarseGridSolve",pc_gamg->use_parallel_coarse_grid_solver,&pc_gamg->use_parallel_coarse_grid_solver,NULL);CHKERRQ(ierr);
121594ae4db5SBarry 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);
121694ae4db5SBarry 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);
1217a303c832SJed Brown     ierr = PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL);CHKERRQ(ierr);
1218c1eae691SMark Adams     n = PETSC_GAMG_MAXLEVELS;
1219c1eae691SMark Adams     ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr);
1220efd3c5ceSMark Adams     if (!flag || n < PETSC_GAMG_MAXLEVELS) {
1221efd3c5ceSMark Adams       if (!flag) n = 1;
1222c1eae691SMark Adams       i = n;
1223c1eae691SMark Adams       do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1224c1eae691SMark Adams     }
122594ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1226b7cbab4eSMark Adams 
1227b7cbab4eSMark Adams     /* set options for subtype */
1228e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1229676e1743SMark F. Adams   }
123014a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
123114a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
1232676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
12335b89ad90SMark F. Adams   PetscFunctionReturn(0);
12345b89ad90SMark F. Adams }
12355b89ad90SMark F. Adams 
12365b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
12375b89ad90SMark F. Adams /*MC
12381cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
12395b89ad90SMark F. Adams 
1240280d9858SJed Brown    Options Database Keys:
1241cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1242cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1243cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1244cab9ed1eSBarry Smith .   -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother
1245cab9ed1eSBarry Smith .   -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit>
1246cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1247cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
12486008e27bSRichard Tran Mills .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1249c1eae691SMark Adams -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1250cab9ed1eSBarry Smith 
1251cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1252cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1253cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1254cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1255cab9ed1eSBarry Smith 
1256db9745e2SBarry Smith    Multigrid options:
1257db9745e2SBarry Smith +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1258db9745e2SBarry Smith .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1259db9745e2SBarry Smith .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1260cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
12615b89ad90SMark F. Adams 
12621cc46a46SBarry Smith 
126395452b02SPatrick Sanan   Notes:
126495452b02SPatrick Sanan     In order to obtain good performance for PCGAMG for vector valued problems you must
1265db9745e2SBarry Smith        Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1266db9745e2SBarry Smith        Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1267db9745e2SBarry Smith        See the Users Manual Chapter 4 for more details
12681cc46a46SBarry Smith 
12695b89ad90SMark F. Adams   Level: intermediate
1270280d9858SJed Brown 
12711cc46a46SBarry Smith   Concepts: algebraic multigrid
12725b89ad90SMark F. Adams 
12731cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1274171cca9aSMark Adams            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation()
12755b89ad90SMark F. Adams M*/
1276b2573a8aSBarry Smith 
12778cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
12785b89ad90SMark F. Adams {
1279c1eae691SMark Adams   PetscErrorCode ierr,i;
12805b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
12815b89ad90SMark F. Adams   PC_MG          *mg;
12825b89ad90SMark F. Adams 
12835b89ad90SMark F. Adams   PetscFunctionBegin;
12841c1aac46SBarry Smith    /* register AMG type */
12851c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
12861c1aac46SBarry Smith 
12875b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
12881c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
12895b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
12905b89ad90SMark F. Adams 
12915b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1292b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
129369aca0b8SBarry Smith   ierr         = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
12945b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
12955b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
12965b89ad90SMark F. Adams 
1297b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
12981ab5ffc9SJed Brown 
12999d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
13009d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
13019d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
13029d5b6da9SMark F. Adams   pc_gamg->data    = 0;
13035b89ad90SMark F. Adams 
13049d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
13055b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
13065b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
13075b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
13085b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
13095adeb434SBarry Smith   mg->view                = PCView_GAMG;
13105b89ad90SMark F. Adams 
1311bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1312bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1313cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr);
13141cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1315cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr);
1316171cca9aSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr);
1317bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1318c1eae691SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr);
1319bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1320c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1321bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
13229d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1323d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
13240c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1325171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1326038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
132725a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1328c1eae691SMark Adams   for (i=0;i<PETSC_GAMG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1329c1eae691SMark Adams   pc_gamg->threshold_scale = 1.;
1330c1eae691SMark Adams   pc_gamg->Nlevels          = PETSC_GAMG_MAXLEVELS;
13319ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
1332c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
13339d5b6da9SMark F. Adams 
1334bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1335bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
13365b89ad90SMark F. Adams   PetscFunctionReturn(0);
13375b89ad90SMark F. Adams }
13383e3471ccSMark Adams 
13393e3471ccSMark Adams /*@C
13403e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
13413e3471ccSMark Adams     from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
13423e3471ccSMark Adams     when using static libraries.
13433e3471ccSMark Adams 
13443e3471ccSMark Adams  Level: developer
13453e3471ccSMark Adams 
13463e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
13473e3471ccSMark Adams  .seealso: PetscInitialize()
13483e3471ccSMark Adams @*/
13493e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
13503e3471ccSMark Adams {
13513e3471ccSMark Adams   PetscErrorCode ierr;
13523e3471ccSMark Adams 
13533e3471ccSMark Adams   PetscFunctionBegin;
13543e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
13553e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
13563e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
13573e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
13588e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
13593e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1360c1c463dbSMark Adams 
1361c1c463dbSMark Adams   /* general events */
1362fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1363fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1364fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1365fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1366c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1367c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1368fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1369c1c463dbSMark Adams 
13705b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13715b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
13725b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
13735b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
13745b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
13755b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
13765b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
13775b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
13785b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1379bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
13805b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
13815b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
13825b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
13835b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
13845b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
13855b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
13865b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
13875b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
13885b89ad90SMark F. Adams 
13895b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
13905b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
13915b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
13925b89ad90SMark F. Adams   /* create timer stages */
13935b89ad90SMark F. Adams #if defined GAMG_STAGES
13945b89ad90SMark F. Adams   {
13955b89ad90SMark F. Adams     char     str[32];
13965b89ad90SMark F. Adams     PetscInt lidx;
13975b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
13985b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
13995b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
14005b89ad90SMark F. Adams       sprintf(str,"MG Level %d",lidx);
14015b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
14025b89ad90SMark F. Adams     }
14035b89ad90SMark F. Adams   }
14045b89ad90SMark F. Adams #endif
14055b89ad90SMark F. Adams #endif
14063e3471ccSMark Adams   PetscFunctionReturn(0);
14073e3471ccSMark Adams }
14083e3471ccSMark Adams 
14093e3471ccSMark Adams /*@C
14101c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
14111c1aac46SBarry Smith     called from PetscFinalize() automatically.
14123e3471ccSMark Adams 
14133e3471ccSMark Adams  Level: developer
14143e3471ccSMark Adams 
14153e3471ccSMark Adams  .keywords: Petsc, destroy, package
14163e3471ccSMark Adams  .seealso: PetscFinalize()
14173e3471ccSMark Adams @*/
14183e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
14193e3471ccSMark Adams {
14203e3471ccSMark Adams   PetscErrorCode ierr;
14213e3471ccSMark Adams 
14223e3471ccSMark Adams   PetscFunctionBegin;
14233e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
14243e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
14253e3471ccSMark Adams   PetscFunctionReturn(0);
14263e3471ccSMark Adams }
1427a36cf38bSToby Isaac 
1428a36cf38bSToby Isaac /*@C
1429a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1430a36cf38bSToby Isaac 
1431a36cf38bSToby Isaac  Input Parameters:
1432a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1433a36cf38bSToby Isaac  - create - function for creating the gamg context.
1434a36cf38bSToby Isaac 
1435a36cf38bSToby Isaac   Level: advanced
1436a36cf38bSToby Isaac 
14371c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1438a36cf38bSToby Isaac @*/
1439a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1440a36cf38bSToby Isaac {
1441a36cf38bSToby Isaac   PetscErrorCode ierr;
1442a36cf38bSToby Isaac 
1443a36cf38bSToby Isaac   PetscFunctionBegin;
1444a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1445a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1446a36cf38bSToby Isaac   PetscFunctionReturn(0);
1447a36cf38bSToby Isaac }
1448a36cf38bSToby Isaac 
1449