xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision cf8ae1d3cf0eec4bff4f7d4eff91937e3edf3150)
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 
22b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
23b4fbaa2aSMark F. Adams 
24b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
250cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
26b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
27b4fbaa2aSMark F. Adams #endif
28f96513f1SMatthew G Knepley 
29140e18c1SBarry Smith static PetscFunctionList GAMGList = 0;
303e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
319d5b6da9SMark F. Adams 
32d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
33d3d6bff4SMark F. Adams #undef __FUNCT__
34d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
35d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
36d3d6bff4SMark F. Adams {
37d3d6bff4SMark F. Adams   PetscErrorCode ierr;
38d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
39d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
40d3d6bff4SMark F. Adams 
41d3d6bff4SMark F. Adams   PetscFunctionBegin;
4262294041SBarry Smith   if (pc_gamg->data) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen, cleaned up in SetUp\n");
431c1aac46SBarry Smith   pc_gamg->data_sz = 0;
44878e152fSMark F. Adams   ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
45a2f3521dSMark F. Adams   PetscFunctionReturn(0);
46a2f3521dSMark F. Adams }
47a2f3521dSMark F. Adams 
485b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
495b89ad90SMark F. Adams /*
50c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
51a147abb0SMark F. Adams      of active processors.
525b89ad90SMark F. Adams 
535b89ad90SMark F. Adams    Input Parameter:
54a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
55a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
569d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
57c5bfad50SMark F. Adams    . cr_bs - coarse block size
583530afc2SMark F. Adams    In/Output Parameter:
59a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
60afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
6111e60469SMark F. Adams    Output Parameter:
623530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
635b89ad90SMark F. Adams */
645cb416c2SMark F. Adams 
655b89ad90SMark F. Adams #undef __FUNCT__
66c238b0ebSToby Isaac #define __FUNCT__ "PCGAMGCreateLevel_GAMG"
67171cca9aSMark 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)
685b89ad90SMark F. Adams {
69a2f3521dSMark F. Adams   PetscErrorCode  ierr;
709d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
71486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
72a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
733b4367a7SBarry Smith   MPI_Comm        comm;
74c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
753ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
765b89ad90SMark F. Adams 
775b89ad90SMark F. Adams   PetscFunctionBegin;
783b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
793b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
803b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
81c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
829d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
83038e3b61SMark F. Adams 
843ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
850298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
863ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
873ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
8873911c69SBarry Smith   } else {
893ae0bb68SMark Adams     PetscInt  bs;
903ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
913ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
923ae0bb68SMark Adams   }
93a2f3521dSMark F. Adams 
94c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
95171cca9aSMark Adams   if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) new_size = 1;
96171cca9aSMark Adams   else {
97472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
980298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
99a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
1007f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
101c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
102a2f3521dSMark F. Adams   }
103f852f58cSMark F. Adams 
1043cb8563fSToby Isaac   if (Pcolumnperm) *Pcolumnperm = NULL;
1053cb8563fSToby Isaac 
106a90e85d9SMark Adams   if (!pc_gamg->repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1072fa5cd67SKarl Rupp   else {
1083ae0bb68SMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
109885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
110e33ef3b1SMark F. Adams 
11171959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
11271959b99SBarry 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);
1130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1140cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
115b4fbaa2aSMark F. Adams #endif
116a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
117785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
118a90e85d9SMark Adams     if (pc_gamg->repart) {
119a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1205a9b9e01SMark F. Adams       Mat adj;
1215a9b9e01SMark F. Adams 
122*cf8ae1d3SMark Adams      ierr = PetscInfo3(pc,"Repartition: size (active): %D --> %D, %D local equations\n",*a_nactive_proc,new_size,ncrs_eq);CHKERRQ(ierr);
1235a9b9e01SMark F. Adams 
124a2f3521dSMark F. Adams       /* get 'adj' */
125c5bfad50SMark F. Adams       if (cr_bs == 1) {
126038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
127806fa848SBarry Smith       } else {
128a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
129eb07cef2SMark F. Adams         Mat               tMat;
130a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
131b4fbaa2aSMark F. Adams         const PetscScalar *vals;
132b4fbaa2aSMark F. Adams         const PetscInt    *idx;
133a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1349057884aSMark F. Adams         static PetscInt   llev = 0;
135d9558ea9SBarry Smith         MatType           mtype;
136b4fbaa2aSMark F. Adams 
137e632b94dSBarry Smith         ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr);
138a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
139a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
140c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
14158471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
142c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
143c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
14458471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1453ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1463ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
14758471d46SMark F. Adams         }
1486876a03eSMark F. Adams 
149d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
1503b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1513ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
152d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
153a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
154a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
155e632b94dSBarry Smith         ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
156eb07cef2SMark F. Adams 
157a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
158c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
15922063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
160eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
161c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
162eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
163eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
164eb07cef2SMark F. Adams           }
16522063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
166eb07cef2SMark F. Adams         }
167eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169eb07cef2SMark F. Adams 
170b4fbaa2aSMark F. Adams         if (llev++ == -1) {
171b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
1728caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
1733b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
174b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
1753bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
176b4fbaa2aSMark F. Adams         }
177eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
178eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
179a2f3521dSMark F. Adams       } /* create 'adj' */
180f150b916SMark F. Adams 
181a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
1825a9b9e01SMark F. Adams         char            prefix[256];
1835a9b9e01SMark F. Adams         const char      *pcpre;
184b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
185b4fbaa2aSMark F. Adams         MatPartitioning mpart;
186a4b7d37bSMark F. Adams         IS              proc_is;
187a2f3521dSMark F. Adams         PetscInt        targetPE;
1882f03bc48SMark F. Adams 
1893b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
1905ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
1919d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
1928caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
19359a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
19411e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
195c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
196a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
19711e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1985a9b9e01SMark F. Adams 
1995ef31b24SMark F. Adams         /* collect IS info */
200785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
201a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
202a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
203c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
204a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
205c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
206a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
207eb07cef2SMark F. Adams           }
2085ef31b24SMark F. Adams         }
209a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
210a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2115ef31b24SMark F. Adams       }
2125ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2135a9b9e01SMark F. Adams 
2143b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2158263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
216806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
217a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
21862294041SBarry Smith 
2195a9b9e01SMark F. Adams       /* find factor */
220c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2215a9b9e01SMark F. Adams       else {
2225a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2235a9b9e01SMark F. Adams         jj = -1;
224c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
2252a82295dSMark Adams           if (!(size%kk)) { /* a candidate */
226c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2275a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2285a9b9e01SMark F. Adams             if (fact > best_fact) {
2295a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
2305a9b9e01SMark F. Adams             }
2315a9b9e01SMark F. Adams           }
2325a9b9e01SMark F. Adams         }
2335a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
234a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
2355a9b9e01SMark F. Adams       }
236c5df96a5SBarry Smith       new_size = size/rfactor;
2375a9b9e01SMark F. Adams 
238c5df96a5SBarry Smith       if (new_size==nactive) {
239a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2405a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
241302440fdSBarry Smith         ierr = PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr);
242fe682e87SBarry Smith #if defined PETSC_GAMG_USE_LOG
243fe682e87SBarry Smith         ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
244fe682e87SBarry Smith #endif
2455a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2465a9b9e01SMark F. Adams       }
2475a9b9e01SMark F. Adams 
248302440fdSBarry Smith       ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr);
249c5df96a5SBarry Smith       targetPE = rank/rfactor;
2503b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
251a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
252e33ef3b1SMark F. Adams 
25311e60469SMark F. Adams     /*
254a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
25511e60469SMark F. Adams      */
256a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2577700e67bSMark Adams     is_eq_num_prim = is_eq_num;
25811e60469SMark F. Adams     /*
259a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
26011e60469SMark F. Adams      */
261c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
262c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
263a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
2643ae0bb68SMark Adams     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */
265a2f3521dSMark F. Adams 
266a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2670cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2680cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
269b4fbaa2aSMark F. Adams #endif
270885364a3SMark 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 */
271885364a3SMark Adams     {
272885364a3SMark Adams     Vec            src_crd, dest_crd;
273885364a3SMark 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;
274885364a3SMark Adams     VecScatter     vecscat;
275885364a3SMark Adams     PetscScalar    *array;
276885364a3SMark Adams     IS isscat;
277a2f3521dSMark F. Adams 
278a2f3521dSMark F. Adams     /* move data (for primal equations only) */
27922063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
2803b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
2813ae0bb68SMark Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
282c0dedaeaSBarry Smith     ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
28311e60469SMark F. Adams     /*
2849d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
285c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
28611e60469SMark F. Adams      */
287854ce69bSBarry Smith     ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
288a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2893ae0bb68SMark Adams     for (ii=0,jj=0; ii<ncrs; ii++) {
290c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
291a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
29211e60469SMark F. Adams     }
293a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2943ae0bb68SMark Adams     ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
29592a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
29611e60469SMark F. Adams     /*
29711e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
29811e60469SMark F. Adams      */
2993ae0bb68SMark Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
3009d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3013ae0bb68SMark Adams       const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3023ae0bb68SMark Adams       for (ii=0; ii<ncrs; ii++) {
3039d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
304a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
305c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
306676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
307d3d6bff4SMark F. Adams         }
308038e3b61SMark F. Adams       }
309eb07cef2SMark F. Adams     }
310eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
311eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
31211e60469SMark F. Adams     /*
31311e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
31411e60469SMark F. Adams       to the correct processor
31511e60469SMark F. Adams     */
3160298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
31711e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
31811e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31911e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32011e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
32111e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
32211e60469SMark F. Adams     /*
32311e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
32411e60469SMark F. Adams     */
325c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
326578f55a3SPeter Brune     ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3272fa5cd67SKarl Rupp 
3283ae0bb68SMark Adams     pc_gamg->data_sz = node_data_sz*ncrs_new;
3293ae0bb68SMark Adams     strideNew        = ncrs_new*ndata_rows;
3302fa5cd67SKarl Rupp 
33111e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3329d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3333ae0bb68SMark Adams       for (ii=0; ii<ncrs_new; ii++) {
3349d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
335a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
336c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
337d3d6bff4SMark F. Adams         }
338038e3b61SMark F. Adams       }
339038e3b61SMark F. Adams     }
34011e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
34111e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
342885364a3SMark Adams     }
343a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3440cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3450cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
346ed3f9983SMark F. Adams #endif
347a2f3521dSMark F. Adams 
34811e60469SMark F. Adams     /*
34911e60469SMark F. Adams       Invert for MatGetSubMatrix
35011e60469SMark F. Adams     */
351a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
352a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
353c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
354a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
355a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
356a2f3521dSMark F. Adams     }
3573cb8563fSToby Isaac     if (Pcolumnperm) {
3583cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3593cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3603cb8563fSToby Isaac     }
361a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3620cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3630cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3640cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
365ed3f9983SMark F. Adams #endif
366a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
367a2f3521dSMark F. Adams     {
368a2f3521dSMark F. Adams       Mat mat;
369806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
370a2f3521dSMark F. Adams       *a_Amat_crs = mat;
371a2f3521dSMark F. Adams     }
372038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
373a2f3521dSMark F. Adams 
3740cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3750cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
376ed3f9983SMark F. Adams #endif
37711e60469SMark F. Adams     /* prolongator */
37811e60469SMark F. Adams     {
37911e60469SMark F. Adams       IS       findices;
380a2f3521dSMark F. Adams       PetscInt Istart,Iend;
381a2f3521dSMark F. Adams       Mat      Pnew;
38262294041SBarry Smith 
383a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
3840cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3850cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
386ed3f9983SMark F. Adams #endif
3873b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
388c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
389806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
39011e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
391c5bfad50SMark F. Adams 
3920cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3930cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
394ed3f9983SMark F. Adams #endif
3953530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
396a2f3521dSMark F. Adams 
397a2f3521dSMark F. Adams       /* output - repartitioned */
398a2f3521dSMark F. Adams       *a_P_inout = Pnew;
399e33ef3b1SMark F. Adams     }
400a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4015b89ad90SMark F. Adams 
402c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
403a2f3521dSMark F. Adams   }
4045b89ad90SMark F. Adams   PetscFunctionReturn(0);
4055b89ad90SMark F. Adams }
4065b89ad90SMark F. Adams 
4075b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4085b89ad90SMark F. Adams /*
4095b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4105b89ad90SMark F. Adams                     by setting data structures and options.
4115b89ad90SMark F. Adams 
4125b89ad90SMark F. Adams    Input Parameter:
4135b89ad90SMark F. Adams .  pc - the preconditioner context
4145b89ad90SMark F. Adams 
4155b89ad90SMark F. Adams */
4165b89ad90SMark F. Adams #undef __FUNCT__
4175b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4189d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4195b89ad90SMark F. Adams {
4205b89ad90SMark F. Adams   PetscErrorCode ierr;
4219d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4225b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4232adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
424171cca9aSMark Adams   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
4253b4367a7SBarry Smith   MPI_Comm       comm;
426c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
427587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
428e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
429a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
430569f4572SMark Adams   MatInfo        info;
431171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
4325ef31b24SMark F. Adams 
4335b89ad90SMark F. Adams   PetscFunctionBegin;
4343b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4353b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4363b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
437dfd5c07aSMark F. Adams 
43884d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
4391c1aac46SBarry Smith     if ((PetscBool)(!pc_gamg->reuse_prol)) {
440878e152fSMark F. Adams       /* reset everything */
441878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
442878e152fSMark F. Adams       pc->setupcalled = 0;
443806fa848SBarry Smith     } else {
44484d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
44503a628feSMark F. Adams       /* just do Galerkin grids */
44658471d46SMark F. Adams       Mat          B,dA,dB;
44758471d46SMark F. Adams 
44871959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
4499d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
45058471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
45123ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
45258471d46SMark F. Adams         /* (re)set to get dirty flag */
45323ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
45458471d46SMark F. Adams 
4552fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
45603a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
4570a97e771SToby Isaac           if (pc_gamg->setup_count==2) {
45803a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
459084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
4602fa5cd67SKarl Rupp 
46103a628feSMark F. Adams             mglevels[level]->A = B;
462806fa848SBarry Smith           } else {
46323ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
46458471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
46503a628feSMark F. Adams           }
46623ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
46758471d46SMark F. Adams           dB   = B;
46858471d46SMark F. Adams         }
4695f8cf99dSMark F. Adams       }
470d5280255SMark F. Adams 
471d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
47258471d46SMark F. Adams       PetscFunctionReturn(0);
473eb07cef2SMark F. Adams     }
474878e152fSMark F. Adams   }
475f6536408SMark F. Adams 
476878e152fSMark F. Adams   if (!pc_gamg->data) {
477878e152fSMark F. Adams     if (pc_gamg->orig_data) {
478878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
4790298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
4802fa5cd67SKarl Rupp 
481878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
482878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
483878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
4842fa5cd67SKarl Rupp 
485785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
486878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
487806fa848SBarry Smith     } else {
4881ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
4897700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
4909d5b6da9SMark F. Adams     }
491878e152fSMark F. Adams   }
492878e152fSMark F. Adams 
493878e152fSMark F. Adams   /* cache original data for reuse */
4941c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
495785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
496878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
497878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
498878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
499878e152fSMark F. Adams   }
500038e3b61SMark F. Adams 
501302f38e8SMark F. Adams   /* get basic dims */
502302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
503171cca9aSMark Adams   ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr);
50484d3f75bSMark F. Adams 
505569f4572SMark Adams   ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
506569f4572SMark Adams   nnz0   = info.nz_used;
507569f4572SMark Adams   nnztot = info.nz_used;
50862294041SBarry 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);
509569f4572SMark Adams 
510a2f3521dSMark F. Adams   /* Get A_i and R_i */
51162294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
5129ab59c8bSMark Adams     pc_gamg->current_level = level;
5135b89ad90SMark F. Adams     level1 = level + 1;
5140cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5150cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
516a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
517a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
518b4fbaa2aSMark F. Adams #endif
519a2f3521dSMark F. Adams #endif
520c8b0795cSMark F. Adams     { /* construct prolongator */
521725b86d8SJed Brown       Mat              Gmat;
5220cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5237700e67bSMark Adams       Mat              Prol11;
524c8b0795cSMark F. Adams 
5257700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5261ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5277700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
528c8b0795cSMark F. Adams 
529a2f3521dSMark F. Adams       /* could have failed to create new level */
530a2f3521dSMark F. Adams       if (Prol11) {
5319d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5320298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
533a2f3521dSMark F. Adams 
534fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
535c8b0795cSMark F. Adams           /* smooth */
536fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
537c8b0795cSMark F. Adams         }
538c8b0795cSMark F. Adams 
5397700e67bSMark Adams         Parr[level1] = Prol11;
540171cca9aSMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
541ffc955d6SMark F. Adams 
5420c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
5431b18a24aSMark Adams         PetscInt bs;
5441b18a24aSMark Adams         ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
5450a3c815dSMark Adams         ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
546ffc955d6SMark F. Adams       }
547ffc955d6SMark F. Adams 
548a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
54941b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
550a2f3521dSMark F. Adams     } /* construct prolongator scope */
5510cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5520cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
553c8b0795cSMark F. Adams #endif
5547f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
555171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
556569f4572SMark Adams       ierr =  PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr);
55762294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
558a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
559a90e85d9SMark Adams #endif
560c8b0795cSMark F. Adams       break;
561c8b0795cSMark F. Adams     }
5620cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5630cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
564b4fbaa2aSMark F. Adams #endif
565171cca9aSMark Adams     ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */
566171cca9aSMark Adams     if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ????????");
567171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
568171cca9aSMark Adams     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr);
569a2f3521dSMark F. Adams 
5700cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5710cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
572b4fbaa2aSMark F. Adams #endif
573171cca9aSMark Adams     ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */
574569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
575569f4572SMark Adams     nnztot += info.nz_used;
5761d5b2942SMark 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);
577569f4572SMark Adams 
5780cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
579b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
580b4fbaa2aSMark F. Adams #endif
581a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
5829ab59c8bSMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) {
5839ab59c8bSMark Adams       ierr =  PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr);
584a90e85d9SMark Adams       level++;
585a90e85d9SMark Adams       break;
586a90e85d9SMark Adams     }
587c8b0795cSMark F. Adams   } /* levels */
588c8b0795cSMark F. Adams   ierr                  = PetscFree(pc_gamg->data);CHKERRQ(ierr);
589c8b0795cSMark F. Adams 
590569f4572SMark Adams   ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr);
5919d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
5925b89ad90SMark F. Adams   fine_level       = level;
5930298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
5945b89ad90SMark F. Adams 
59562294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
596d5280255SMark F. Adams     /* set default smoothers & set operators */
59762294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
598ffc955d6SMark F. Adams       KSP smoother;
599ffc955d6SMark F. Adams       PC  subpc;
600a2f3521dSMark F. Adams 
6019d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
602f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
603ffc955d6SMark F. Adams 
604a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
605a2f3521dSMark F. Adams       /* set ops */
60623ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
607a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
608a2f3521dSMark F. Adams 
609a2f3521dSMark F. Adams       /* set defaults */
6106c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
611a2f3521dSMark F. Adams 
6120c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
6130c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
6142d3561bbSSatish Balay         PetscInt sz;
6157a28f3e5SMark Adams         IS       *iss;
616a2f3521dSMark F. Adams 
6172d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
6187a28f3e5SMark Adams         iss   = ASMLocalIDsArr[level];
6190c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
6200a3c815dSMark Adams         ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr);
6210c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
6227f66b68fSMark Adams         if (!sz) {
623ffc955d6SMark F. Adams           IS       is;
6240a3c815dSMark Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
6257a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
626a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
627806fa848SBarry Smith         } else {
628a94c3b12SMark F. Adams           PetscInt kk;
6297a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr);
630a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
6317a28f3e5SMark Adams             ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr);
632a94c3b12SMark F. Adams           }
6337a28f3e5SMark Adams           ierr = PetscFree(iss);CHKERRQ(ierr);
634ffc955d6SMark F. Adams         }
6350298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
636ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
637806fa848SBarry Smith       } else {
638890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
639ffc955d6SMark F. Adams       }
640d5280255SMark F. Adams     }
641d5280255SMark F. Adams     {
642d5280255SMark F. Adams       /* coarse grid */
643d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
644d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
645d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
64623ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
647*cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
648d5280255SMark F. Adams         ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
649d5280255SMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
650d5280255SMark F. Adams         ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
651d5280255SMark F. Adams         ierr = PCSetUp(subpc);CHKERRQ(ierr);
65271959b99SBarry Smith         ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
65371959b99SBarry Smith         if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
654d5280255SMark F. Adams         ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
655d5280255SMark F. Adams         ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
6569dbfc187SHong Zhang         ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
6572fb0b348SMark F. Adams         ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
65808e36f19SMark Adams         ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
6595b42dca8SJed Brown         /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
6605b42dca8SJed Brown          * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
6615b42dca8SJed Brown          * view every subdomain as though they were different. */
6625b42dca8SJed Brown         ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
663d5280255SMark F. Adams       }
664*cf8ae1d3SMark Adams     }
665d5280255SMark F. Adams 
666d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
667d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
668e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
669d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
6701c1aac46SBarry Smith     if (!mg->galerkin) SETERRQ(comm,PETSC_ERR_USER,"PCGAMG must use Galerkin for coarse operators.");
6711c1aac46SBarry Smith     if (mg->galerkin == 1) mg->galerkin = 2;
672d5280255SMark F. Adams 
673d5280255SMark F. Adams     /* clean up */
674d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
675587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
676587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
6775b89ad90SMark F. Adams     }
6780cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
679806fa848SBarry Smith   } else {
6805f8cf99dSMark F. Adams     KSP smoother;
681302440fdSBarry Smith     ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
6829d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
68323ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
6845f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
6859d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
6865f8cf99dSMark F. Adams   }
6875b89ad90SMark F. Adams   PetscFunctionReturn(0);
6885b89ad90SMark F. Adams }
6895b89ad90SMark F. Adams 
690eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
6915b89ad90SMark F. Adams /*
6925b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
6935b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
6945b89ad90SMark F. Adams 
6955b89ad90SMark F. Adams    Input Parameter:
6965b89ad90SMark F. Adams .  pc - the preconditioner context
6975b89ad90SMark F. Adams 
6985b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
6995b89ad90SMark F. Adams */
7005b89ad90SMark F. Adams #undef __FUNCT__
7015b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
7025b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
7035b89ad90SMark F. Adams {
7045b89ad90SMark F. Adams   PetscErrorCode ierr;
7055b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
7065b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
7075b89ad90SMark F. Adams 
7085b89ad90SMark F. Adams   PetscFunctionBegin;
7095b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7109b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
7119b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
7129b8ffb57SJed Brown   }
71314a9496bSBarry Smith   ierr = PetscRandomDestroy(&pc_gamg->random);CHKERRQ(ierr);
7141ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
7151ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
7165b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7175b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7185b89ad90SMark F. Adams   PetscFunctionReturn(0);
7195b89ad90SMark F. Adams }
7205b89ad90SMark F. Adams 
721676e1743SMark F. Adams #undef __FUNCT__
722676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
723676e1743SMark F. Adams /*@
724cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
725676e1743SMark F. Adams 
7261cc46a46SBarry Smith    Logically Collective on PC
727676e1743SMark F. Adams 
728676e1743SMark F. Adams    Input Parameters:
7291cc46a46SBarry Smith +  pc - the preconditioner context
7301cc46a46SBarry Smith -  n - the number of equations
731676e1743SMark F. Adams 
732676e1743SMark F. Adams 
733676e1743SMark F. Adams    Options Database Key:
7341cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
735676e1743SMark F. Adams 
736cab9ed1eSBarry Smith    Notes: GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
737cab9ed1eSBarry Smith           that has degrees of freedom
738cab9ed1eSBarry Smith 
739676e1743SMark F. Adams    Level: intermediate
740676e1743SMark F. Adams 
7411c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
742676e1743SMark F. Adams 
7431c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim()
744676e1743SMark F. Adams @*/
745676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
746676e1743SMark F. Adams {
747676e1743SMark F. Adams   PetscErrorCode ierr;
748676e1743SMark F. Adams 
749676e1743SMark F. Adams   PetscFunctionBegin;
750676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
751676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
752676e1743SMark F. Adams   PetscFunctionReturn(0);
753676e1743SMark F. Adams }
754676e1743SMark F. Adams 
755676e1743SMark F. Adams #undef __FUNCT__
756676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
7571e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
758676e1743SMark F. Adams {
759c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
760c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
761676e1743SMark F. Adams 
762676e1743SMark F. Adams   PetscFunctionBegin;
7639d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
764676e1743SMark F. Adams   PetscFunctionReturn(0);
765676e1743SMark F. Adams }
766676e1743SMark F. Adams 
767676e1743SMark F. Adams #undef __FUNCT__
768389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
769389730f3SMark F. Adams /*@
770cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
771389730f3SMark F. Adams 
772389730f3SMark F. Adams  Collective on PC
773389730f3SMark F. Adams 
774389730f3SMark F. Adams    Input Parameters:
7751cc46a46SBarry Smith +  pc - the preconditioner context
7761cc46a46SBarry Smith -  n - maximum number of equations to aim for
777389730f3SMark F. Adams 
778389730f3SMark F. Adams    Options Database Key:
7791cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
780389730f3SMark F. Adams 
781389730f3SMark F. Adams    Level: intermediate
782389730f3SMark F. Adams 
7831c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
784389730f3SMark F. Adams 
7851c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim()
786389730f3SMark F. Adams @*/
787389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
788389730f3SMark F. Adams {
789389730f3SMark F. Adams   PetscErrorCode ierr;
790389730f3SMark F. Adams 
791389730f3SMark F. Adams   PetscFunctionBegin;
792389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
793389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
794389730f3SMark F. Adams   PetscFunctionReturn(0);
795389730f3SMark F. Adams }
796389730f3SMark F. Adams 
797389730f3SMark F. Adams #undef __FUNCT__
798389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
7991e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
800389730f3SMark F. Adams {
801389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
802389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
803389730f3SMark F. Adams 
804389730f3SMark F. Adams   PetscFunctionBegin;
8059d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
806389730f3SMark F. Adams   PetscFunctionReturn(0);
807389730f3SMark F. Adams }
808389730f3SMark F. Adams 
809389730f3SMark F. Adams #undef __FUNCT__
810cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGSetRepartition"
811676e1743SMark F. Adams /*@
812cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
813676e1743SMark F. Adams 
814676e1743SMark F. Adams    Collective on PC
815676e1743SMark F. Adams 
816676e1743SMark F. Adams    Input Parameters:
8171cc46a46SBarry Smith +  pc - the preconditioner context
8181cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
819676e1743SMark F. Adams 
820676e1743SMark F. Adams    Options Database Key:
8211cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
822676e1743SMark F. Adams 
823cab9ed1eSBarry Smith    Notes: this will generally improve the loading balancing of the work on each level
824cab9ed1eSBarry Smith 
825676e1743SMark F. Adams    Level: intermediate
826676e1743SMark F. Adams 
8271c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
828676e1743SMark F. Adams 
829676e1743SMark F. Adams .seealso: ()
830676e1743SMark F. Adams @*/
831cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
832676e1743SMark F. Adams {
833676e1743SMark F. Adams   PetscErrorCode ierr;
834676e1743SMark F. Adams 
835676e1743SMark F. Adams   PetscFunctionBegin;
836676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
837cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
838676e1743SMark F. Adams   PetscFunctionReturn(0);
839676e1743SMark F. Adams }
840676e1743SMark F. Adams 
841676e1743SMark F. Adams #undef __FUNCT__
842cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGSetRepartition_GAMG"
843cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
844676e1743SMark F. Adams {
845c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
846c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
847676e1743SMark F. Adams 
848676e1743SMark F. Adams   PetscFunctionBegin;
8499d5b6da9SMark F. Adams   pc_gamg->repart = n;
850676e1743SMark F. Adams   PetscFunctionReturn(0);
851676e1743SMark F. Adams }
852676e1743SMark F. Adams 
853676e1743SMark F. Adams #undef __FUNCT__
8541cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation"
855dfd5c07aSMark F. Adams /*@
856cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
857dfd5c07aSMark F. Adams 
858dfd5c07aSMark F. Adams    Collective on PC
859dfd5c07aSMark F. Adams 
860dfd5c07aSMark F. Adams    Input Parameters:
8611cc46a46SBarry Smith +  pc - the preconditioner context
8621cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
863dfd5c07aSMark F. Adams 
864dfd5c07aSMark F. Adams    Options Database Key:
8651cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
866dfd5c07aSMark F. Adams 
867dfd5c07aSMark F. Adams    Level: intermediate
868dfd5c07aSMark F. Adams 
869cab9ed1eSBarry Smith    Notes: this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
870cab9ed1eSBarry Smith           rebuilding the preconditioner quicker.
871cab9ed1eSBarry Smith 
8721c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
873dfd5c07aSMark F. Adams 
874dfd5c07aSMark F. Adams .seealso: ()
875dfd5c07aSMark F. Adams @*/
8761cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
877dfd5c07aSMark F. Adams {
878dfd5c07aSMark F. Adams   PetscErrorCode ierr;
879dfd5c07aSMark F. Adams 
880dfd5c07aSMark F. Adams   PetscFunctionBegin;
881dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8821cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
883dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
884dfd5c07aSMark F. Adams }
885dfd5c07aSMark F. Adams 
886dfd5c07aSMark F. Adams #undef __FUNCT__
8871cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG"
8881cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
889dfd5c07aSMark F. Adams {
890dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
891dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
892dfd5c07aSMark F. Adams 
893dfd5c07aSMark F. Adams   PetscFunctionBegin;
894dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
895dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
896dfd5c07aSMark F. Adams }
897dfd5c07aSMark F. Adams 
898dfd5c07aSMark F. Adams #undef __FUNCT__
899cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGASMSetUseAggs"
900ffc955d6SMark F. Adams /*@
901cab9ed1eSBarry 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.
902ffc955d6SMark F. Adams 
903ffc955d6SMark F. Adams    Collective on PC
904ffc955d6SMark F. Adams 
905ffc955d6SMark F. Adams    Input Parameters:
906cab9ed1eSBarry Smith +  pc - the preconditioner context
907cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
908ffc955d6SMark F. Adams 
909ffc955d6SMark F. Adams    Options Database Key:
910cab9ed1eSBarry Smith .  -pc_gamg_asm_use_agg
911ffc955d6SMark F. Adams 
912ffc955d6SMark F. Adams    Level: intermediate
913ffc955d6SMark F. Adams 
9141c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
915ffc955d6SMark F. Adams 
916ffc955d6SMark F. Adams .seealso: ()
917ffc955d6SMark F. Adams @*/
918cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
919ffc955d6SMark F. Adams {
920ffc955d6SMark F. Adams   PetscErrorCode ierr;
921ffc955d6SMark F. Adams 
922ffc955d6SMark F. Adams   PetscFunctionBegin;
923ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
924cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
925ffc955d6SMark F. Adams   PetscFunctionReturn(0);
926ffc955d6SMark F. Adams }
927ffc955d6SMark F. Adams 
928ffc955d6SMark F. Adams #undef __FUNCT__
929cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGASMSetUseAggs_GAMG"
930cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
931ffc955d6SMark F. Adams {
932ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
933ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
934ffc955d6SMark F. Adams 
935ffc955d6SMark F. Adams   PetscFunctionBegin;
936cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
937ffc955d6SMark F. Adams   PetscFunctionReturn(0);
938ffc955d6SMark F. Adams }
939ffc955d6SMark F. Adams 
940ffc955d6SMark F. Adams #undef __FUNCT__
941171cca9aSMark Adams #define __FUNCT__ "PCGAMGSetUseParallelCoarseGridSolve"
942171cca9aSMark Adams /*@
943*cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
944171cca9aSMark Adams 
945171cca9aSMark Adams    Collective on PC
946171cca9aSMark Adams 
947171cca9aSMark Adams    Input Parameters:
948171cca9aSMark Adams +  pc - the preconditioner context
949*cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
950171cca9aSMark Adams 
951171cca9aSMark Adams    Options Database Key:
952*cf8ae1d3SMark Adams .  -pc_gamg_use_parallel_coarse_grid_solver
953171cca9aSMark Adams 
954171cca9aSMark Adams    Level: intermediate
955171cca9aSMark Adams 
956171cca9aSMark Adams    Concepts: Unstructured multigrid preconditioner
957171cca9aSMark Adams 
958171cca9aSMark Adams .seealso: ()
959171cca9aSMark Adams @*/
960171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
961171cca9aSMark Adams {
962171cca9aSMark Adams   PetscErrorCode ierr;
963171cca9aSMark Adams 
964171cca9aSMark Adams   PetscFunctionBegin;
965171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
966171cca9aSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
967171cca9aSMark Adams   PetscFunctionReturn(0);
968171cca9aSMark Adams }
969171cca9aSMark Adams 
970171cca9aSMark Adams #undef __FUNCT__
971171cca9aSMark Adams #define __FUNCT__ "PCGAMGSetUseParallelCoarseGridSolve_GAMG"
972171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
973171cca9aSMark Adams {
974171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
975171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
976171cca9aSMark Adams 
977171cca9aSMark Adams   PetscFunctionBegin;
978171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
979171cca9aSMark Adams   PetscFunctionReturn(0);
980171cca9aSMark Adams }
981171cca9aSMark Adams 
982171cca9aSMark Adams #undef __FUNCT__
9834ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
9844ef23d27SMark F. Adams /*@
9851cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
9864ef23d27SMark F. Adams 
9874ef23d27SMark F. Adams    Not collective on PC
9884ef23d27SMark F. Adams 
9894ef23d27SMark F. Adams    Input Parameters:
9901cc46a46SBarry Smith +  pc - the preconditioner
9911cc46a46SBarry Smith -  n - the maximum number of levels to use
9924ef23d27SMark F. Adams 
9934ef23d27SMark F. Adams    Options Database Key:
9944ef23d27SMark F. Adams .  -pc_mg_levels
9954ef23d27SMark F. Adams 
9964ef23d27SMark F. Adams    Level: intermediate
9974ef23d27SMark F. Adams 
9981c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
9994ef23d27SMark F. Adams 
10004ef23d27SMark F. Adams .seealso: ()
10014ef23d27SMark F. Adams @*/
10024ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
10034ef23d27SMark F. Adams {
10044ef23d27SMark F. Adams   PetscErrorCode ierr;
10054ef23d27SMark F. Adams 
10064ef23d27SMark F. Adams   PetscFunctionBegin;
10074ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10084ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
10094ef23d27SMark F. Adams   PetscFunctionReturn(0);
10104ef23d27SMark F. Adams }
10114ef23d27SMark F. Adams 
10124ef23d27SMark F. Adams #undef __FUNCT__
10134ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
10141e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
10154ef23d27SMark F. Adams {
10164ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
10174ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
10184ef23d27SMark F. Adams 
10194ef23d27SMark F. Adams   PetscFunctionBegin;
10209d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
10214ef23d27SMark F. Adams   PetscFunctionReturn(0);
10224ef23d27SMark F. Adams }
10234ef23d27SMark F. Adams 
10244ef23d27SMark F. Adams #undef __FUNCT__
10253542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
10263542efc5SMark F. Adams /*@
10273542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
10283542efc5SMark F. Adams 
10293542efc5SMark F. Adams    Not collective on PC
10303542efc5SMark F. Adams 
10313542efc5SMark F. Adams    Input Parameters:
10321cc46a46SBarry Smith +  pc - the preconditioner context
1033b001cb0fSBarry Smith -  threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
10343542efc5SMark F. Adams 
10353542efc5SMark F. Adams    Options Database Key:
10361cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
10373542efc5SMark F. Adams 
1038cab9ed1eSBarry Smith    Notes: Before aggregating the graph GAMG will remove small values from the graph thus reducing the coupling in the graph and a different
1039cab9ed1eSBarry Smith     (perhaps better) coarser set of points.
1040cab9ed1eSBarry Smith 
10413542efc5SMark F. Adams    Level: intermediate
10423542efc5SMark F. Adams 
10431c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
10443542efc5SMark F. Adams 
10453542efc5SMark F. Adams .seealso: ()
10463542efc5SMark F. Adams @*/
10473542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
10483542efc5SMark F. Adams {
10493542efc5SMark F. Adams   PetscErrorCode ierr;
10503542efc5SMark F. Adams 
10513542efc5SMark F. Adams   PetscFunctionBegin;
10523542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10533542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
10543542efc5SMark F. Adams   PetscFunctionReturn(0);
10553542efc5SMark F. Adams }
10563542efc5SMark F. Adams 
10573542efc5SMark F. Adams #undef __FUNCT__
10583542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
10591e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
10603542efc5SMark F. Adams {
1061c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1062c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
10633542efc5SMark F. Adams 
10643542efc5SMark F. Adams   PetscFunctionBegin;
10659d5b6da9SMark F. Adams   pc_gamg->threshold = n;
10663542efc5SMark F. Adams   PetscFunctionReturn(0);
10673542efc5SMark F. Adams }
10683542efc5SMark F. Adams 
10693542efc5SMark F. Adams #undef __FUNCT__
10709d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1071676e1743SMark F. Adams /*@
1072c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1073676e1743SMark F. Adams 
1074676e1743SMark F. Adams    Collective on PC
1075676e1743SMark F. Adams 
1076676e1743SMark F. Adams    Input Parameters:
1077c60c7ad4SBarry Smith +  pc - the preconditioner context
1078c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1079676e1743SMark F. Adams 
1080676e1743SMark F. Adams    Options Database Key:
1081cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1082676e1743SMark F. Adams 
1083676e1743SMark F. Adams    Level: intermediate
1084676e1743SMark F. Adams 
10851c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1086676e1743SMark F. Adams 
1087cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1088676e1743SMark F. Adams @*/
108919fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1090676e1743SMark F. Adams {
1091676e1743SMark F. Adams   PetscErrorCode ierr;
1092676e1743SMark F. Adams 
1093676e1743SMark F. Adams   PetscFunctionBegin;
1094676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1095806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1096676e1743SMark F. Adams   PetscFunctionReturn(0);
1097676e1743SMark F. Adams }
1098676e1743SMark F. Adams 
1099676e1743SMark F. Adams #undef __FUNCT__
1100c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType"
1101c60c7ad4SBarry Smith /*@
1102c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1103c60c7ad4SBarry Smith 
1104c60c7ad4SBarry Smith    Collective on PC
1105c60c7ad4SBarry Smith 
1106c60c7ad4SBarry Smith    Input Parameter:
1107c60c7ad4SBarry Smith .  pc - the preconditioner context
1108c60c7ad4SBarry Smith 
1109c60c7ad4SBarry Smith    Output Parameter:
1110c60c7ad4SBarry Smith .  type - the type of algorithm used
1111c60c7ad4SBarry Smith 
1112c60c7ad4SBarry Smith    Level: intermediate
1113c60c7ad4SBarry Smith 
11141c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1115c60c7ad4SBarry Smith 
11161c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1117c60c7ad4SBarry Smith @*/
1118c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1119c60c7ad4SBarry Smith {
1120c60c7ad4SBarry Smith   PetscErrorCode ierr;
1121c60c7ad4SBarry Smith 
1122c60c7ad4SBarry Smith   PetscFunctionBegin;
1123c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1124c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1125c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1126c60c7ad4SBarry Smith }
1127c60c7ad4SBarry Smith 
1128c60c7ad4SBarry Smith #undef __FUNCT__
1129c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG"
1130c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1131c60c7ad4SBarry Smith {
1132c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1133c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1134c60c7ad4SBarry Smith 
1135c60c7ad4SBarry Smith   PetscFunctionBegin;
1136c60c7ad4SBarry Smith   *type = pc_gamg->type;
1137c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1138c60c7ad4SBarry Smith }
1139c60c7ad4SBarry Smith 
1140c60c7ad4SBarry Smith #undef __FUNCT__
11419d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
11421e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1143676e1743SMark F. Adams {
11449d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
11451ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
11461ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1147676e1743SMark F. Adams 
1148676e1743SMark F. Adams   PetscFunctionBegin;
1149c60c7ad4SBarry Smith   pc_gamg->type = type;
11501c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
11519d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
11521ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
11531ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
11541ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1155e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
11563ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
11573ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
11583ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
11593ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
11603ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
11613ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
11623ae0bb68SMark Adams     pc_gamg->data_sz = 0;
11631ab5ffc9SJed Brown   }
11641ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
11651ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
11669d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1167676e1743SMark F. Adams   PetscFunctionReturn(0);
1168676e1743SMark F. Adams }
1169676e1743SMark F. Adams 
11705b89ad90SMark F. Adams #undef __FUNCT__
11715adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG"
11725adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
11735adeb434SBarry Smith {
11745adeb434SBarry Smith   PetscErrorCode ierr;
11755adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
11765adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
11775adeb434SBarry Smith 
11785adeb434SBarry Smith   PetscFunctionBegin;
11795adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1180b001cb0fSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values from graph %g\n",(double)pc_gamg->threshold);CHKERRQ(ierr);
1181cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
1182cab9ed1eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr);
1183cab9ed1eSBarry Smith   }
1184171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
1185171cca9aSMark Adams     ierr = PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr);
1186171cca9aSMark Adams   }
11875adeb434SBarry Smith   if (pc_gamg->ops->view) {
11885adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
11895adeb434SBarry Smith   }
11905adeb434SBarry Smith   PetscFunctionReturn(0);
11915adeb434SBarry Smith }
11925adeb434SBarry Smith 
11935adeb434SBarry Smith #undef __FUNCT__
11945b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
11954416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
11965b89ad90SMark F. Adams {
1197676e1743SMark F. Adams   PetscErrorCode ierr;
1198676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1199676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1200676e1743SMark F. Adams   PetscBool      flag;
12013b4367a7SBarry Smith   MPI_Comm       comm;
120214a9496bSBarry Smith   char           prefix[256];
120314a9496bSBarry Smith   const char     *pcpre;
12045b89ad90SMark F. Adams 
12055b89ad90SMark F. Adams   PetscFunctionBegin;
12063b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1207e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1208676e1743SMark F. Adams   {
1209bd94a7aaSJed Brown     char tname[256];
12101a1c1e04SBarry Smith     ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1211bd94a7aaSJed Brown     if (flag) {
1212bd94a7aaSJed Brown       ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
12131ab5ffc9SJed Brown     }
1214cab9ed1eSBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
12151cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1216cab9ed1eSBarry Smith     ierr = PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation agragates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL);CHKERRQ(ierr);
1217*cf8ae1d3SMark 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);
121894ae4db5SBarry 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);
121994ae4db5SBarry 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);
122094ae4db5SBarry Smith     ierr = PetscOptionsReal("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&pc_gamg->threshold,&flag);CHKERRQ(ierr);
122194ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1222b7cbab4eSMark Adams 
1223b7cbab4eSMark Adams     /* set options for subtype */
1224e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1225676e1743SMark F. Adams   }
122614a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
122714a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
122814a9496bSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)pc_gamg->random,prefix);CHKERRQ(ierr);
122914a9496bSBarry Smith   ierr = PetscRandomSetFromOptions(pc_gamg->random);CHKERRQ(ierr);
1230676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
12315b89ad90SMark F. Adams   PetscFunctionReturn(0);
12325b89ad90SMark F. Adams }
12335b89ad90SMark F. Adams 
12345b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
12355b89ad90SMark F. Adams /*MC
12361cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
12375b89ad90SMark F. Adams 
1238280d9858SJed Brown    Options Database Keys:
1239cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1240cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1241cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1242cab9ed1eSBarry 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
1243cab9ed1eSBarry 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>
1244cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1245cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
1246cab9ed1eSBarry Smith -   -pc_gamg_threshold <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph thus reducing the coupling in the graph and a different
1247cab9ed1eSBarry Smith 
1248cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1249cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1250cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1251cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1252cab9ed1eSBarry Smith 
1253cab9ed1eSBarry Smith    Multigrid options(inherited):
12541cc46a46SBarry Smith +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1255280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1256280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1257cab9ed1eSBarry Smith .  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
1258cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
12595b89ad90SMark F. Adams 
12601cc46a46SBarry Smith 
12611cc46a46SBarry Smith   Notes: In order to obtain good performance for PCGAMG for vector valued problems you must
12621cc46a46SBarry Smith $       Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
12631cc46a46SBarry Smith $       Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
12641cc46a46SBarry Smith $       See the Users Manual Chapter 4 for more details
12651cc46a46SBarry Smith 
12665b89ad90SMark F. Adams   Level: intermediate
1267280d9858SJed Brown 
12681cc46a46SBarry Smith   Concepts: algebraic multigrid
12695b89ad90SMark F. Adams 
12701cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1271171cca9aSMark Adams            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation()
12725b89ad90SMark F. Adams M*/
1273b2573a8aSBarry Smith 
12745b89ad90SMark F. Adams #undef __FUNCT__
12755b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
12768cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
12775b89ad90SMark F. Adams {
12785b89ad90SMark F. Adams   PetscErrorCode ierr;
12795b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
12805b89ad90SMark F. Adams   PC_MG          *mg;
12815b89ad90SMark F. Adams 
12825b89ad90SMark F. Adams   PetscFunctionBegin;
12831c1aac46SBarry Smith    /* register AMG type */
12841c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
12851c1aac46SBarry Smith 
12865b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
12871c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
12885b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
12895b89ad90SMark F. Adams 
12905b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1291b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
12925b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
12931c1aac46SBarry Smith   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally from PCMG by GAMG code */
12945b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
12955b89ad90SMark F. Adams 
1296b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
12971ab5ffc9SJed Brown 
12989d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
12999d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
13009d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
13019d5b6da9SMark F. Adams   pc_gamg->data    = 0;
13025b89ad90SMark F. Adams 
13039d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
13045b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
13055b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
13065b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
13075b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
13085adeb434SBarry Smith   mg->view                = PCView_GAMG;
13095b89ad90SMark F. Adams 
1310bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1311bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1312cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr);
13131cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1314cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr);
1315171cca9aSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr);
1316bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1317bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1318c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1319bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
13209d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1321d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
13220c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1323171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1324038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
132525a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1326d3042614SMark Adams   pc_gamg->threshold        = 0.;
13279d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
13289ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
1329c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
13309d5b6da9SMark F. Adams 
133114a9496bSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)pc),&pc_gamg->random);CHKERRQ(ierr);
133214a9496bSBarry Smith 
1333bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1334bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
13355b89ad90SMark F. Adams   PetscFunctionReturn(0);
13365b89ad90SMark F. Adams }
13373e3471ccSMark Adams 
13383e3471ccSMark Adams #undef __FUNCT__
13393e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage"
13403e3471ccSMark Adams /*@C
13413e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
13423e3471ccSMark Adams     from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
13433e3471ccSMark Adams     when using static libraries.
13443e3471ccSMark Adams 
13453e3471ccSMark Adams  Level: developer
13463e3471ccSMark Adams 
13473e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
13483e3471ccSMark Adams  .seealso: PetscInitialize()
13493e3471ccSMark Adams @*/
13503e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
13513e3471ccSMark Adams {
13523e3471ccSMark Adams   PetscErrorCode ierr;
13533e3471ccSMark Adams 
13543e3471ccSMark Adams   PetscFunctionBegin;
13553e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
13563e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
13573e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
13583e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
13598e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
13603e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1361c1c463dbSMark Adams 
1362c1c463dbSMark Adams   /* general events */
1363fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1364fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1365fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1366fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1367c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1368c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1369fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1370c1c463dbSMark Adams 
13715b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13725b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
13735b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
13745b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
13755b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
13765b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
13775b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
13785b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
13795b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1380bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
13815b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
13825b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
13835b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
13845b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
13855b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
13865b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
13875b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
13885b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
13895b89ad90SMark F. Adams 
13905b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
13915b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
13925b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
13935b89ad90SMark F. Adams   /* create timer stages */
13945b89ad90SMark F. Adams #if defined GAMG_STAGES
13955b89ad90SMark F. Adams   {
13965b89ad90SMark F. Adams     char     str[32];
13975b89ad90SMark F. Adams     PetscInt lidx;
13985b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
13995b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
14005b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
14015b89ad90SMark F. Adams       sprintf(str,"MG Level %d",lidx);
14025b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
14035b89ad90SMark F. Adams     }
14045b89ad90SMark F. Adams   }
14055b89ad90SMark F. Adams #endif
14065b89ad90SMark F. Adams #endif
14073e3471ccSMark Adams   PetscFunctionReturn(0);
14083e3471ccSMark Adams }
14093e3471ccSMark Adams 
14103e3471ccSMark Adams #undef __FUNCT__
14113e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage"
14123e3471ccSMark Adams /*@C
14131c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
14141c1aac46SBarry Smith     called from PetscFinalize() automatically.
14153e3471ccSMark Adams 
14163e3471ccSMark Adams  Level: developer
14173e3471ccSMark Adams 
14183e3471ccSMark Adams  .keywords: Petsc, destroy, package
14193e3471ccSMark Adams  .seealso: PetscFinalize()
14203e3471ccSMark Adams @*/
14213e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
14223e3471ccSMark Adams {
14233e3471ccSMark Adams   PetscErrorCode ierr;
14243e3471ccSMark Adams 
14253e3471ccSMark Adams   PetscFunctionBegin;
14263e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
14273e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
14283e3471ccSMark Adams   PetscFunctionReturn(0);
14293e3471ccSMark Adams }
1430a36cf38bSToby Isaac 
1431a36cf38bSToby Isaac #undef __FUNCT__
1432a36cf38bSToby Isaac #define __FUNCT__ "PCGAMGRegister"
1433a36cf38bSToby Isaac /*@C
1434a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1435a36cf38bSToby Isaac 
1436a36cf38bSToby Isaac  Input Parameters:
1437a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1438a36cf38bSToby Isaac  - create - function for creating the gamg context.
1439a36cf38bSToby Isaac 
1440a36cf38bSToby Isaac   Level: advanced
1441a36cf38bSToby Isaac 
14421c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1443a36cf38bSToby Isaac @*/
1444a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1445a36cf38bSToby Isaac {
1446a36cf38bSToby Isaac   PetscErrorCode ierr;
1447a36cf38bSToby Isaac 
1448a36cf38bSToby Isaac   PetscFunctionBegin;
1449a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1450a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1451a36cf38bSToby Isaac   PetscFunctionReturn(0);
1452a36cf38bSToby Isaac }
1453a36cf38bSToby Isaac 
1454