xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 0c3bc534e6b4d730785c434f036dde260507c8bc)
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"
6762294041SBarry Smith 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)
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 */
95a2f3521dSMark F. Adams   {
96472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
970298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
98a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
997f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
100c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
101a2f3521dSMark F. Adams   }
102f852f58cSMark F. Adams 
1033cb8563fSToby Isaac   if (Pcolumnperm) *Pcolumnperm = NULL;
1043cb8563fSToby Isaac 
105a90e85d9SMark Adams   if (!pc_gamg->repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1062fa5cd67SKarl Rupp   else {
1073ae0bb68SMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
108885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
109e33ef3b1SMark F. Adams 
11071959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
11171959b99SBarry 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);
1120cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1130cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
114b4fbaa2aSMark F. Adams #endif
115a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
116785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
117a90e85d9SMark Adams     if (pc_gamg->repart) {
118a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1195a9b9e01SMark F. Adams       Mat adj;
1205a9b9e01SMark F. Adams 
121302440fdSBarry Smith      ierr = PetscInfo3(pc,"Repartition: size (active): %D --> %D, neq = %D\n",*a_nactive_proc,new_size,ncrs_eq);CHKERRQ(ierr);
1225a9b9e01SMark F. Adams 
123a2f3521dSMark F. Adams       /* get 'adj' */
124c5bfad50SMark F. Adams       if (cr_bs == 1) {
125038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
126806fa848SBarry Smith       } else {
127a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
128eb07cef2SMark F. Adams         Mat               tMat;
129a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
130b4fbaa2aSMark F. Adams         const PetscScalar *vals;
131b4fbaa2aSMark F. Adams         const PetscInt    *idx;
132a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1339057884aSMark F. Adams         static PetscInt   llev = 0;
134d9558ea9SBarry Smith         MatType           mtype;
135b4fbaa2aSMark F. Adams 
136e632b94dSBarry Smith         ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr);
137a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
138a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
139c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
14058471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
141c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
142c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
14358471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1443ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1453ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
14658471d46SMark F. Adams         }
1476876a03eSMark F. Adams 
148d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
1493b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1503ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
151d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
152a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
153a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
154e632b94dSBarry Smith         ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
155eb07cef2SMark F. Adams 
156a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
157c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
15822063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
159eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
160c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
161eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
162eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
163eb07cef2SMark F. Adams           }
16422063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
165eb07cef2SMark F. Adams         }
166eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168eb07cef2SMark F. Adams 
169b4fbaa2aSMark F. Adams         if (llev++ == -1) {
170b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
1718caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
1723b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
173b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
1743bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
175b4fbaa2aSMark F. Adams         }
176eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
177eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
178a2f3521dSMark F. Adams       } /* create 'adj' */
179f150b916SMark F. Adams 
180a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
1815a9b9e01SMark F. Adams         char            prefix[256];
1825a9b9e01SMark F. Adams         const char      *pcpre;
183b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
184b4fbaa2aSMark F. Adams         MatPartitioning mpart;
185a4b7d37bSMark F. Adams         IS              proc_is;
186a2f3521dSMark F. Adams         PetscInt        targetPE;
1872f03bc48SMark F. Adams 
1883b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
1895ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
1909d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
1918caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
19259a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
19311e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
194c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
195a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
19611e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1975a9b9e01SMark F. Adams 
1985ef31b24SMark F. Adams         /* collect IS info */
199785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
200a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
201a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
202c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
203a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
204c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
205a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
206eb07cef2SMark F. Adams           }
2075ef31b24SMark F. Adams         }
208a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
209a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2105ef31b24SMark F. Adams       }
2115ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2125a9b9e01SMark F. Adams 
2133b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2148263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
215806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
216a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
21762294041SBarry Smith 
2185a9b9e01SMark F. Adams       /* find factor */
219c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2205a9b9e01SMark F. Adams       else {
2215a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2225a9b9e01SMark F. Adams         jj = -1;
223c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
2242a82295dSMark Adams           if (!(size%kk)) { /* a candidate */
225c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2265a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2275a9b9e01SMark F. Adams             if (fact > best_fact) {
2285a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
2295a9b9e01SMark F. Adams             }
2305a9b9e01SMark F. Adams           }
2315a9b9e01SMark F. Adams         }
2325a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
233a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
2345a9b9e01SMark F. Adams       }
235c5df96a5SBarry Smith       new_size = size/rfactor;
2365a9b9e01SMark F. Adams 
237c5df96a5SBarry Smith       if (new_size==nactive) {
238a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2395a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
240302440fdSBarry Smith         ierr = PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr);
241fe682e87SBarry Smith #if defined PETSC_GAMG_USE_LOG
242fe682e87SBarry Smith         ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
243fe682e87SBarry Smith #endif
2445a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2455a9b9e01SMark F. Adams       }
2465a9b9e01SMark F. Adams 
247302440fdSBarry Smith       ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr);
248c5df96a5SBarry Smith       targetPE = rank/rfactor;
2493b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
250a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
251e33ef3b1SMark F. Adams 
25211e60469SMark F. Adams     /*
253a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
25411e60469SMark F. Adams      */
255a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2567700e67bSMark Adams     is_eq_num_prim = is_eq_num;
25711e60469SMark F. Adams     /*
258a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
25911e60469SMark F. Adams      */
260c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
261c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
262a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
2633ae0bb68SMark Adams     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */
264a2f3521dSMark F. Adams 
265a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2660cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2670cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
268b4fbaa2aSMark F. Adams #endif
269885364a3SMark 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 */
270885364a3SMark Adams     {
271885364a3SMark Adams     Vec            src_crd, dest_crd;
272885364a3SMark 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;
273885364a3SMark Adams     VecScatter     vecscat;
274885364a3SMark Adams     PetscScalar    *array;
275885364a3SMark Adams     IS isscat;
276a2f3521dSMark F. Adams 
277a2f3521dSMark F. Adams     /* move data (for primal equations only) */
27822063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
2793b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
2803ae0bb68SMark Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
281c0dedaeaSBarry Smith     ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
28211e60469SMark F. Adams     /*
2839d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
284c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
28511e60469SMark F. Adams      */
286854ce69bSBarry Smith     ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
287a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2883ae0bb68SMark Adams     for (ii=0,jj=0; ii<ncrs; ii++) {
289c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
290a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
29111e60469SMark F. Adams     }
292a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2933ae0bb68SMark Adams     ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
29492a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
29511e60469SMark F. Adams     /*
29611e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
29711e60469SMark F. Adams      */
2983ae0bb68SMark Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
2999d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3003ae0bb68SMark Adams       const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3013ae0bb68SMark Adams       for (ii=0; ii<ncrs; ii++) {
3029d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
303a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
304c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
305676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
306d3d6bff4SMark F. Adams         }
307038e3b61SMark F. Adams       }
308eb07cef2SMark F. Adams     }
309eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
310eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
31111e60469SMark F. Adams     /*
31211e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
31311e60469SMark F. Adams       to the correct processor
31411e60469SMark F. Adams     */
3150298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
31611e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
31711e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31811e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31911e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
32011e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
32111e60469SMark F. Adams     /*
32211e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
32311e60469SMark F. Adams     */
324c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
325578f55a3SPeter Brune     ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3262fa5cd67SKarl Rupp 
3273ae0bb68SMark Adams     pc_gamg->data_sz = node_data_sz*ncrs_new;
3283ae0bb68SMark Adams     strideNew        = ncrs_new*ndata_rows;
3292fa5cd67SKarl Rupp 
33011e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3319d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3323ae0bb68SMark Adams       for (ii=0; ii<ncrs_new; ii++) {
3339d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
334a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
335c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
336d3d6bff4SMark F. Adams         }
337038e3b61SMark F. Adams       }
338038e3b61SMark F. Adams     }
33911e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
34011e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
341885364a3SMark Adams     }
342a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3430cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3440cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
345ed3f9983SMark F. Adams #endif
346a2f3521dSMark F. Adams 
34711e60469SMark F. Adams     /*
34811e60469SMark F. Adams       Invert for MatGetSubMatrix
34911e60469SMark F. Adams     */
350a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
351a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
352c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
353a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
354a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
355a2f3521dSMark F. Adams     }
3563cb8563fSToby Isaac     if (Pcolumnperm) {
3573cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3583cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3593cb8563fSToby Isaac     }
360a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3610cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3620cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3630cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
364ed3f9983SMark F. Adams #endif
365a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
366a2f3521dSMark F. Adams     {
367a2f3521dSMark F. Adams       Mat mat;
368806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
369a2f3521dSMark F. Adams       *a_Amat_crs = mat;
370a2f3521dSMark F. Adams     }
371038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
372a2f3521dSMark F. Adams 
3730cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3740cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
375ed3f9983SMark F. Adams #endif
37611e60469SMark F. Adams     /* prolongator */
37711e60469SMark F. Adams     {
37811e60469SMark F. Adams       IS       findices;
379a2f3521dSMark F. Adams       PetscInt Istart,Iend;
380a2f3521dSMark F. Adams       Mat      Pnew;
38162294041SBarry Smith 
382a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
3830cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3840cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
385ed3f9983SMark F. Adams #endif
3863b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
387c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
388806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
38911e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
390c5bfad50SMark F. Adams 
3910cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3920cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
393ed3f9983SMark F. Adams #endif
3943530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
395a2f3521dSMark F. Adams 
396a2f3521dSMark F. Adams       /* output - repartitioned */
397a2f3521dSMark F. Adams       *a_P_inout = Pnew;
398e33ef3b1SMark F. Adams     }
399a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4005b89ad90SMark F. Adams 
401c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
402a2f3521dSMark F. Adams   }
4035b89ad90SMark F. Adams   PetscFunctionReturn(0);
4045b89ad90SMark F. Adams }
4055b89ad90SMark F. Adams 
4065b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4075b89ad90SMark F. Adams /*
4085b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4095b89ad90SMark F. Adams                     by setting data structures and options.
4105b89ad90SMark F. Adams 
4115b89ad90SMark F. Adams    Input Parameter:
4125b89ad90SMark F. Adams .  pc - the preconditioner context
4135b89ad90SMark F. Adams 
4145b89ad90SMark F. Adams */
4155b89ad90SMark F. Adams #undef __FUNCT__
4165b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4179d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4185b89ad90SMark F. Adams {
4195b89ad90SMark F. Adams   PetscErrorCode ierr;
4209d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4215b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4222adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
423a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
4243b4367a7SBarry Smith   MPI_Comm       comm;
425c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
426587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
427e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
428a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
429569f4572SMark Adams   MatInfo        info;
4305ef31b24SMark F. Adams 
4315b89ad90SMark F. Adams   PetscFunctionBegin;
4323b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4333b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4343b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
435dfd5c07aSMark F. Adams 
43684d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
4371c1aac46SBarry Smith     if ((PetscBool)(!pc_gamg->reuse_prol)) {
438878e152fSMark F. Adams       /* reset everything */
439878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
440878e152fSMark F. Adams       pc->setupcalled = 0;
441806fa848SBarry Smith     } else {
44284d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
44303a628feSMark F. Adams       /* just do Galerkin grids */
44458471d46SMark F. Adams       Mat          B,dA,dB;
44558471d46SMark F. Adams 
44671959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
4479d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
44858471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
44923ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
45058471d46SMark F. Adams         /* (re)set to get dirty flag */
45123ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
45258471d46SMark F. Adams 
4532fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
45403a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
4550a97e771SToby Isaac           if (pc_gamg->setup_count==2) {
45603a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
457084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
4582fa5cd67SKarl Rupp 
45903a628feSMark F. Adams             mglevels[level]->A = B;
460806fa848SBarry Smith           } else {
46123ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
46258471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
46303a628feSMark F. Adams           }
46423ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
46558471d46SMark F. Adams           dB   = B;
46658471d46SMark F. Adams         }
4675f8cf99dSMark F. Adams       }
468d5280255SMark F. Adams 
469d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
47058471d46SMark F. Adams       PetscFunctionReturn(0);
471eb07cef2SMark F. Adams     }
472878e152fSMark F. Adams   }
473f6536408SMark F. Adams 
474878e152fSMark F. Adams   if (!pc_gamg->data) {
475878e152fSMark F. Adams     if (pc_gamg->orig_data) {
476878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
4770298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
4782fa5cd67SKarl Rupp 
479878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
480878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
481878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
4822fa5cd67SKarl Rupp 
483785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
484878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
485806fa848SBarry Smith     } else {
4861ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
4877700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
4889d5b6da9SMark F. Adams     }
489878e152fSMark F. Adams   }
490878e152fSMark F. Adams 
491878e152fSMark F. Adams   /* cache original data for reuse */
4921c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
493785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
494878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
495878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
496878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
497878e152fSMark F. Adams   }
498038e3b61SMark F. Adams 
499302f38e8SMark F. Adams   /* get basic dims */
500302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
501a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
50284d3f75bSMark F. Adams 
503569f4572SMark Adams   ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
504569f4572SMark Adams   nnz0   = info.nz_used;
505569f4572SMark Adams   nnztot = info.nz_used;
50662294041SBarry 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);
507569f4572SMark Adams 
508a2f3521dSMark F. Adams   /* Get A_i and R_i */
50962294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
5109ab59c8bSMark Adams     pc_gamg->current_level = level;
5115b89ad90SMark F. Adams     level1 = level + 1;
5120cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5130cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
514a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
515a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
516b4fbaa2aSMark F. Adams #endif
517a2f3521dSMark F. Adams #endif
518c8b0795cSMark F. Adams     { /* construct prolongator */
519725b86d8SJed Brown       Mat              Gmat;
5200cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5217700e67bSMark Adams       Mat              Prol11;
522c8b0795cSMark F. Adams 
5237700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5241ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5257700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
526c8b0795cSMark F. Adams 
527a2f3521dSMark F. Adams       /* could have failed to create new level */
528a2f3521dSMark F. Adams       if (Prol11) {
5299d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5300298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
531a2f3521dSMark F. Adams 
532fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
533c8b0795cSMark F. Adams           /* smooth */
534fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
535c8b0795cSMark F. Adams         }
536c8b0795cSMark F. Adams 
5377700e67bSMark Adams         Parr[level1] = Prol11;
5380298fd71SBarry Smith       } else Parr[level1] = NULL;
539ffc955d6SMark F. Adams 
540*0c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
5411b18a24aSMark Adams         PetscInt bs;
5421b18a24aSMark Adams         ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
543806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
544ffc955d6SMark F. Adams       }
545ffc955d6SMark F. Adams 
546a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
54741b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
548a2f3521dSMark F. Adams     } /* construct prolongator scope */
5490cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5500cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
551c8b0795cSMark F. Adams #endif
5527f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
553c8b0795cSMark F. Adams     if (!Parr[level1]) {
554569f4572SMark Adams       ierr =  PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr);
55562294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
556a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
557a90e85d9SMark Adams #endif
558c8b0795cSMark F. Adams       break;
559c8b0795cSMark F. Adams     }
5600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5610cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
562b4fbaa2aSMark F. Adams #endif
563a2f3521dSMark F. Adams 
5641c1aac46SBarry Smith     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs,&Parr[level1], &Aarr[level1], &nactivepe, NULL);CHKERRQ(ierr);
565a2f3521dSMark F. Adams 
5660cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5670cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
568b4fbaa2aSMark F. Adams #endif
569a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
570569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
571569f4572SMark Adams     nnztot += info.nz_used;
5721d5b2942SMark 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);
573569f4572SMark Adams 
5740cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
575b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
576b4fbaa2aSMark F. Adams #endif
577a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
5789ab59c8bSMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) {
5799ab59c8bSMark Adams       ierr =  PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr);
580a90e85d9SMark Adams       level++;
581a90e85d9SMark Adams       break;
582a90e85d9SMark Adams     }
583c8b0795cSMark F. Adams   } /* levels */
584c8b0795cSMark F. Adams   ierr                  = PetscFree(pc_gamg->data);CHKERRQ(ierr);
585c8b0795cSMark F. Adams 
586569f4572SMark Adams   ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr);
5879d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
5885b89ad90SMark F. Adams   fine_level       = level;
5890298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
5905b89ad90SMark F. Adams 
59162294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
592d5280255SMark F. Adams     /* set default smoothers & set operators */
59362294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
594ffc955d6SMark F. Adams       KSP smoother;
595ffc955d6SMark F. Adams       PC  subpc;
596a2f3521dSMark F. Adams 
5979d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
598f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
599ffc955d6SMark F. Adams 
600a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
601a2f3521dSMark F. Adams       /* set ops */
60223ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
603a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
604a2f3521dSMark F. Adams 
605a2f3521dSMark F. Adams       /* set defaults */
6066c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
607a2f3521dSMark F. Adams 
608*0c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
609*0c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
6102d3561bbSSatish Balay         PetscInt sz;
6112d3561bbSSatish Balay         IS       *is;
612a2f3521dSMark F. Adams 
6132d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
6142d3561bbSSatish Balay         is   = ASMLocalIDsArr[level];
615*0c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
616*0c3bc534SBarry Smith         ierr = PCASMSetOverlap(subpc, 1);CHKERRQ(ierr);
617*0c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
6187f66b68fSMark Adams         if (!sz) {
619ffc955d6SMark F. Adams           IS       is;
620ffc955d6SMark F. Adams           PetscInt my0,kk;
621ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
622ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
623*0c3bc534SBarry Smith           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
624a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
625806fa848SBarry Smith         } else {
626a94c3b12SMark F. Adams           PetscInt kk;
627*0c3bc534SBarry Smith           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, is);CHKERRQ(ierr);
628a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
629a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
630a94c3b12SMark F. Adams           }
631ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
632ffc955d6SMark F. Adams         }
6330298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
634ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
635806fa848SBarry Smith       } else {
636890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
637ffc955d6SMark F. Adams       }
638d5280255SMark F. Adams     }
639d5280255SMark F. Adams     {
640d5280255SMark F. Adams       /* coarse grid */
641d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
642d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
643d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
64423ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
645d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
646d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
647d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
648d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
64971959b99SBarry Smith       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
65071959b99SBarry Smith       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
651d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
652d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
6539dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
6542fb0b348SMark F. Adams       ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
65508e36f19SMark Adams       ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
6565b42dca8SJed Brown       /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
6575b42dca8SJed Brown        * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
6585b42dca8SJed Brown        * view every subdomain as though they were different. */
6595b42dca8SJed Brown       ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
660d5280255SMark F. Adams     }
661d5280255SMark F. Adams 
662d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
663d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
664e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
665d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
6661c1aac46SBarry Smith     if (!mg->galerkin) SETERRQ(comm,PETSC_ERR_USER,"PCGAMG must use Galerkin for coarse operators.");
6671c1aac46SBarry Smith     if (mg->galerkin == 1) mg->galerkin = 2;
668d5280255SMark F. Adams 
669d5280255SMark F. Adams     /* clean up */
670d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
671587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
672587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
6735b89ad90SMark F. Adams     }
6740cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
675806fa848SBarry Smith   } else {
6765f8cf99dSMark F. Adams     KSP smoother;
677302440fdSBarry Smith     ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
6789d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
67923ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
6805f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
6819d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
6825f8cf99dSMark F. Adams   }
6835b89ad90SMark F. Adams   PetscFunctionReturn(0);
6845b89ad90SMark F. Adams }
6855b89ad90SMark F. Adams 
686eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
6875b89ad90SMark F. Adams /*
6885b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
6895b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
6905b89ad90SMark F. Adams 
6915b89ad90SMark F. Adams    Input Parameter:
6925b89ad90SMark F. Adams .  pc - the preconditioner context
6935b89ad90SMark F. Adams 
6945b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
6955b89ad90SMark F. Adams */
6965b89ad90SMark F. Adams #undef __FUNCT__
6975b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
6985b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
6995b89ad90SMark F. Adams {
7005b89ad90SMark F. Adams   PetscErrorCode ierr;
7015b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
7025b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
7035b89ad90SMark F. Adams 
7045b89ad90SMark F. Adams   PetscFunctionBegin;
7055b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7069b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
7079b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
7089b8ffb57SJed Brown   }
70914a9496bSBarry Smith   ierr = PetscRandomDestroy(&pc_gamg->random);CHKERRQ(ierr);
7101ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
7111ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
7125b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7135b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7145b89ad90SMark F. Adams   PetscFunctionReturn(0);
7155b89ad90SMark F. Adams }
7165b89ad90SMark F. Adams 
717676e1743SMark F. Adams #undef __FUNCT__
718676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
719676e1743SMark F. Adams /*@
7201cc46a46SBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via processor reduction.
721676e1743SMark F. Adams 
7221cc46a46SBarry Smith    Logically Collective on PC
723676e1743SMark F. Adams 
724676e1743SMark F. Adams    Input Parameters:
7251cc46a46SBarry Smith +  pc - the preconditioner context
7261cc46a46SBarry Smith -  n - the number of equations
727676e1743SMark F. Adams 
728676e1743SMark F. Adams 
729676e1743SMark F. Adams    Options Database Key:
7301cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
731676e1743SMark F. Adams 
732676e1743SMark F. Adams    Level: intermediate
733676e1743SMark F. Adams 
7341c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
735676e1743SMark F. Adams 
7361c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim()
737676e1743SMark F. Adams @*/
738676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
739676e1743SMark F. Adams {
740676e1743SMark F. Adams   PetscErrorCode ierr;
741676e1743SMark F. Adams 
742676e1743SMark F. Adams   PetscFunctionBegin;
743676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
744676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
745676e1743SMark F. Adams   PetscFunctionReturn(0);
746676e1743SMark F. Adams }
747676e1743SMark F. Adams 
748676e1743SMark F. Adams #undef __FUNCT__
749676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
7501e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
751676e1743SMark F. Adams {
752c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
753c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
754676e1743SMark F. Adams 
755676e1743SMark F. Adams   PetscFunctionBegin;
7569d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
757676e1743SMark F. Adams   PetscFunctionReturn(0);
758676e1743SMark F. Adams }
759676e1743SMark F. Adams 
760676e1743SMark F. Adams #undef __FUNCT__
761389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
762389730f3SMark F. Adams /*@
763389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
764389730f3SMark F. Adams 
765389730f3SMark F. Adams  Collective on PC
766389730f3SMark F. Adams 
767389730f3SMark F. Adams    Input Parameters:
7681cc46a46SBarry Smith +  pc - the preconditioner context
7691cc46a46SBarry Smith -  n - maximum number of equations to aim for
770389730f3SMark F. Adams 
771389730f3SMark F. Adams    Options Database Key:
7721cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
773389730f3SMark F. Adams 
774389730f3SMark F. Adams    Level: intermediate
775389730f3SMark F. Adams 
7761c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
777389730f3SMark F. Adams 
7781c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim()
779389730f3SMark F. Adams @*/
780389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
781389730f3SMark F. Adams {
782389730f3SMark F. Adams   PetscErrorCode ierr;
783389730f3SMark F. Adams 
784389730f3SMark F. Adams   PetscFunctionBegin;
785389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
786389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
787389730f3SMark F. Adams   PetscFunctionReturn(0);
788389730f3SMark F. Adams }
789389730f3SMark F. Adams 
790389730f3SMark F. Adams #undef __FUNCT__
791389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
7921e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
793389730f3SMark F. Adams {
794389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
795389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
796389730f3SMark F. Adams 
797389730f3SMark F. Adams   PetscFunctionBegin;
7989d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
799389730f3SMark F. Adams   PetscFunctionReturn(0);
800389730f3SMark F. Adams }
801389730f3SMark F. Adams 
802389730f3SMark F. Adams #undef __FUNCT__
8038263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
804676e1743SMark F. Adams /*@
8058263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
806676e1743SMark F. Adams 
807676e1743SMark F. Adams    Collective on PC
808676e1743SMark F. Adams 
809676e1743SMark F. Adams    Input Parameters:
8101cc46a46SBarry Smith +  pc - the preconditioner context
8111cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
812676e1743SMark F. Adams 
813676e1743SMark F. Adams    Options Database Key:
8141cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
815676e1743SMark F. Adams 
816676e1743SMark F. Adams    Level: intermediate
817676e1743SMark F. Adams 
8181c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
819676e1743SMark F. Adams 
820676e1743SMark F. Adams .seealso: ()
821676e1743SMark F. Adams @*/
8228263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
823676e1743SMark F. Adams {
824676e1743SMark F. Adams   PetscErrorCode ierr;
825676e1743SMark F. Adams 
826676e1743SMark F. Adams   PetscFunctionBegin;
827676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8288263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
829676e1743SMark F. Adams   PetscFunctionReturn(0);
830676e1743SMark F. Adams }
831676e1743SMark F. Adams 
832676e1743SMark F. Adams #undef __FUNCT__
8338263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
8341e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
835676e1743SMark F. Adams {
836c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
837c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
838676e1743SMark F. Adams 
839676e1743SMark F. Adams   PetscFunctionBegin;
8409d5b6da9SMark F. Adams   pc_gamg->repart = n;
841676e1743SMark F. Adams   PetscFunctionReturn(0);
842676e1743SMark F. Adams }
843676e1743SMark F. Adams 
844676e1743SMark F. Adams #undef __FUNCT__
8451cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation"
846dfd5c07aSMark F. Adams /*@
8471cc46a46SBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner
848dfd5c07aSMark F. Adams 
849dfd5c07aSMark F. Adams    Collective on PC
850dfd5c07aSMark F. Adams 
851dfd5c07aSMark F. Adams    Input Parameters:
8521cc46a46SBarry Smith +  pc - the preconditioner context
8531cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
854dfd5c07aSMark F. Adams 
855dfd5c07aSMark F. Adams    Options Database Key:
8561cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
857dfd5c07aSMark F. Adams 
858dfd5c07aSMark F. Adams    Level: intermediate
859dfd5c07aSMark F. Adams 
8601c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
861dfd5c07aSMark F. Adams 
862dfd5c07aSMark F. Adams .seealso: ()
863dfd5c07aSMark F. Adams @*/
8641cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
865dfd5c07aSMark F. Adams {
866dfd5c07aSMark F. Adams   PetscErrorCode ierr;
867dfd5c07aSMark F. Adams 
868dfd5c07aSMark F. Adams   PetscFunctionBegin;
869dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8701cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
871dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
872dfd5c07aSMark F. Adams }
873dfd5c07aSMark F. Adams 
874dfd5c07aSMark F. Adams #undef __FUNCT__
8751cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG"
8761cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
877dfd5c07aSMark F. Adams {
878dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
879dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
880dfd5c07aSMark F. Adams 
881dfd5c07aSMark F. Adams   PetscFunctionBegin;
882dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
883dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
884dfd5c07aSMark F. Adams }
885dfd5c07aSMark F. Adams 
886dfd5c07aSMark F. Adams #undef __FUNCT__
887ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
888ffc955d6SMark F. Adams /*@
889ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
890ffc955d6SMark F. Adams 
891ffc955d6SMark F. Adams    Collective on PC
892ffc955d6SMark F. Adams 
893ffc955d6SMark F. Adams    Input Parameters:
894ffc955d6SMark F. Adams .  pc - the preconditioner context
895ffc955d6SMark F. Adams 
896ffc955d6SMark F. Adams 
897ffc955d6SMark F. Adams    Options Database Key:
898*0c3bc534SBarry Smith .  -pc_gamg_use_agg_asm
899ffc955d6SMark F. Adams 
900ffc955d6SMark F. Adams    Level: intermediate
901ffc955d6SMark F. Adams 
9021c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
903ffc955d6SMark F. Adams 
904ffc955d6SMark F. Adams .seealso: ()
905ffc955d6SMark F. Adams @*/
906ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
907ffc955d6SMark F. Adams {
908ffc955d6SMark F. Adams   PetscErrorCode ierr;
909ffc955d6SMark F. Adams 
910ffc955d6SMark F. Adams   PetscFunctionBegin;
911ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
912ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
913ffc955d6SMark F. Adams   PetscFunctionReturn(0);
914ffc955d6SMark F. Adams }
915ffc955d6SMark F. Adams 
916ffc955d6SMark F. Adams #undef __FUNCT__
917ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
9181e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
919ffc955d6SMark F. Adams {
920ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
921ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
922ffc955d6SMark F. Adams 
923ffc955d6SMark F. Adams   PetscFunctionBegin;
924*0c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm = n;
925ffc955d6SMark F. Adams   PetscFunctionReturn(0);
926ffc955d6SMark F. Adams }
927ffc955d6SMark F. Adams 
928ffc955d6SMark F. Adams #undef __FUNCT__
9294ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
9304ef23d27SMark F. Adams /*@
9311cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
9324ef23d27SMark F. Adams 
9334ef23d27SMark F. Adams    Not collective on PC
9344ef23d27SMark F. Adams 
9354ef23d27SMark F. Adams    Input Parameters:
9361cc46a46SBarry Smith +  pc - the preconditioner
9371cc46a46SBarry Smith -  n - the maximum number of levels to use
9384ef23d27SMark F. Adams 
9394ef23d27SMark F. Adams    Options Database Key:
9404ef23d27SMark F. Adams .  -pc_mg_levels
9414ef23d27SMark F. Adams 
9424ef23d27SMark F. Adams    Level: intermediate
9434ef23d27SMark F. Adams 
9441c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
9454ef23d27SMark F. Adams 
9464ef23d27SMark F. Adams .seealso: ()
9474ef23d27SMark F. Adams @*/
9484ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
9494ef23d27SMark F. Adams {
9504ef23d27SMark F. Adams   PetscErrorCode ierr;
9514ef23d27SMark F. Adams 
9524ef23d27SMark F. Adams   PetscFunctionBegin;
9534ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9544ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
9554ef23d27SMark F. Adams   PetscFunctionReturn(0);
9564ef23d27SMark F. Adams }
9574ef23d27SMark F. Adams 
9584ef23d27SMark F. Adams #undef __FUNCT__
9594ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
9601e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
9614ef23d27SMark F. Adams {
9624ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
9634ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
9644ef23d27SMark F. Adams 
9654ef23d27SMark F. Adams   PetscFunctionBegin;
9669d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
9674ef23d27SMark F. Adams   PetscFunctionReturn(0);
9684ef23d27SMark F. Adams }
9694ef23d27SMark F. Adams 
9704ef23d27SMark F. Adams #undef __FUNCT__
9713542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
9723542efc5SMark F. Adams /*@
9733542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
9743542efc5SMark F. Adams 
9753542efc5SMark F. Adams    Not collective on PC
9763542efc5SMark F. Adams 
9773542efc5SMark F. Adams    Input Parameters:
9781cc46a46SBarry Smith +  pc - the preconditioner context
979b001cb0fSBarry Smith -  threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
9803542efc5SMark F. Adams 
9813542efc5SMark F. Adams    Options Database Key:
9821cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
9833542efc5SMark F. Adams 
9843542efc5SMark F. Adams    Level: intermediate
9853542efc5SMark F. Adams 
9861c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
9873542efc5SMark F. Adams 
9883542efc5SMark F. Adams .seealso: ()
9893542efc5SMark F. Adams @*/
9903542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
9913542efc5SMark F. Adams {
9923542efc5SMark F. Adams   PetscErrorCode ierr;
9933542efc5SMark F. Adams 
9943542efc5SMark F. Adams   PetscFunctionBegin;
9953542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9963542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
9973542efc5SMark F. Adams   PetscFunctionReturn(0);
9983542efc5SMark F. Adams }
9993542efc5SMark F. Adams 
10003542efc5SMark F. Adams #undef __FUNCT__
10013542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
10021e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
10033542efc5SMark F. Adams {
1004c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1005c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
10063542efc5SMark F. Adams 
10073542efc5SMark F. Adams   PetscFunctionBegin;
10089d5b6da9SMark F. Adams   pc_gamg->threshold = n;
10093542efc5SMark F. Adams   PetscFunctionReturn(0);
10103542efc5SMark F. Adams }
10113542efc5SMark F. Adams 
10123542efc5SMark F. Adams #undef __FUNCT__
10139d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1014676e1743SMark F. Adams /*@
1015c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1016676e1743SMark F. Adams 
1017676e1743SMark F. Adams    Collective on PC
1018676e1743SMark F. Adams 
1019676e1743SMark F. Adams    Input Parameters:
1020c60c7ad4SBarry Smith +  pc - the preconditioner context
1021c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1022676e1743SMark F. Adams 
1023676e1743SMark F. Adams    Options Database Key:
1024c60c7ad4SBarry Smith .  -pc_gamg_type <agg,geo,classical>
1025676e1743SMark F. Adams 
1026676e1743SMark F. Adams    Level: intermediate
1027676e1743SMark F. Adams 
10281c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1029676e1743SMark F. Adams 
1030c60c7ad4SBarry Smith .seealso: PCGAMGGetType(), PCGAMG
1031676e1743SMark F. Adams @*/
103219fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1033676e1743SMark F. Adams {
1034676e1743SMark F. Adams   PetscErrorCode ierr;
1035676e1743SMark F. Adams 
1036676e1743SMark F. Adams   PetscFunctionBegin;
1037676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1038806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1039676e1743SMark F. Adams   PetscFunctionReturn(0);
1040676e1743SMark F. Adams }
1041676e1743SMark F. Adams 
1042676e1743SMark F. Adams #undef __FUNCT__
1043c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType"
1044c60c7ad4SBarry Smith /*@
1045c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1046c60c7ad4SBarry Smith 
1047c60c7ad4SBarry Smith    Collective on PC
1048c60c7ad4SBarry Smith 
1049c60c7ad4SBarry Smith    Input Parameter:
1050c60c7ad4SBarry Smith .  pc - the preconditioner context
1051c60c7ad4SBarry Smith 
1052c60c7ad4SBarry Smith    Output Parameter:
1053c60c7ad4SBarry Smith .  type - the type of algorithm used
1054c60c7ad4SBarry Smith 
1055c60c7ad4SBarry Smith    Level: intermediate
1056c60c7ad4SBarry Smith 
10571c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1058c60c7ad4SBarry Smith 
10591c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1060c60c7ad4SBarry Smith @*/
1061c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1062c60c7ad4SBarry Smith {
1063c60c7ad4SBarry Smith   PetscErrorCode ierr;
1064c60c7ad4SBarry Smith 
1065c60c7ad4SBarry Smith   PetscFunctionBegin;
1066c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1067c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1068c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1069c60c7ad4SBarry Smith }
1070c60c7ad4SBarry Smith 
1071c60c7ad4SBarry Smith #undef __FUNCT__
1072c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG"
1073c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1074c60c7ad4SBarry Smith {
1075c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1076c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1077c60c7ad4SBarry Smith 
1078c60c7ad4SBarry Smith   PetscFunctionBegin;
1079c60c7ad4SBarry Smith   *type = pc_gamg->type;
1080c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1081c60c7ad4SBarry Smith }
1082c60c7ad4SBarry Smith 
1083c60c7ad4SBarry Smith #undef __FUNCT__
10849d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
10851e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1086676e1743SMark F. Adams {
10879d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
10881ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
10891ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1090676e1743SMark F. Adams 
1091676e1743SMark F. Adams   PetscFunctionBegin;
1092c60c7ad4SBarry Smith   pc_gamg->type = type;
10931c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
10949d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
10951ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
10961ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
10971ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1098e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
10993ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
11003ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
11013ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
11023ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
11033ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
11043ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
11053ae0bb68SMark Adams     pc_gamg->data_sz = 0;
11061ab5ffc9SJed Brown   }
11071ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
11081ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
11099d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1110676e1743SMark F. Adams   PetscFunctionReturn(0);
1111676e1743SMark F. Adams }
1112676e1743SMark F. Adams 
11135b89ad90SMark F. Adams #undef __FUNCT__
11145adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG"
11155adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
11165adeb434SBarry Smith {
11175adeb434SBarry Smith   PetscErrorCode ierr;
11185adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
11195adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
11205adeb434SBarry Smith 
11215adeb434SBarry Smith   PetscFunctionBegin;
11225adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1123b001cb0fSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values from graph %g\n",(double)pc_gamg->threshold);CHKERRQ(ierr);
11245adeb434SBarry Smith   if (pc_gamg->ops->view) {
11255adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
11265adeb434SBarry Smith   }
11275adeb434SBarry Smith   PetscFunctionReturn(0);
11285adeb434SBarry Smith }
11295adeb434SBarry Smith 
11305adeb434SBarry Smith #undef __FUNCT__
11315b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
11324416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
11335b89ad90SMark F. Adams {
1134676e1743SMark F. Adams   PetscErrorCode ierr;
1135676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1136676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1137676e1743SMark F. Adams   PetscBool      flag;
11383b4367a7SBarry Smith   MPI_Comm       comm;
113914a9496bSBarry Smith   char           prefix[256];
114014a9496bSBarry Smith   const char     *pcpre;
11415b89ad90SMark F. Adams 
11425b89ad90SMark F. Adams   PetscFunctionBegin;
11433b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1144e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1145676e1743SMark F. Adams   {
1146bd94a7aaSJed Brown     char tname[256];
11471a1c1e04SBarry Smith     ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1148bd94a7aaSJed Brown     if (flag) {
1149bd94a7aaSJed Brown       ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
11501ab5ffc9SJed Brown     }
115194ae4db5SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
11521cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1153*0c3bc534SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_use_agg_asm","Use aggregation agragates for ASM smoother","PCGAMGUseASMAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL);CHKERRQ(ierr);
115494ae4db5SBarry 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);
115594ae4db5SBarry 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);
115694ae4db5SBarry 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);
115794ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1158b7cbab4eSMark Adams 
1159b7cbab4eSMark Adams     /* set options for subtype */
1160e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1161676e1743SMark F. Adams   }
116214a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
116314a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
116414a9496bSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)pc_gamg->random,prefix);CHKERRQ(ierr);
116514a9496bSBarry Smith   ierr = PetscRandomSetFromOptions(pc_gamg->random);CHKERRQ(ierr);
1166676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
11675b89ad90SMark F. Adams   PetscFunctionReturn(0);
11685b89ad90SMark F. Adams }
11695b89ad90SMark F. Adams 
11705b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
11715b89ad90SMark F. Adams /*MC
11721cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
11735b89ad90SMark F. Adams 
1174280d9858SJed Brown    Options Database Keys:
11755b89ad90SMark F. Adams    Multigrid options(inherited)
11761cc46a46SBarry Smith +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1177280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1178280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
11798c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
11805b89ad90SMark F. Adams 
11811cc46a46SBarry Smith 
11821cc46a46SBarry Smith   Notes: In order to obtain good performance for PCGAMG for vector valued problems you must
11831cc46a46SBarry Smith $       Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
11841cc46a46SBarry Smith $       Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
11851cc46a46SBarry Smith $       See the Users Manual Chapter 4 for more details
11861cc46a46SBarry Smith 
11875b89ad90SMark F. Adams   Level: intermediate
1188280d9858SJed Brown 
11891cc46a46SBarry Smith   Concepts: algebraic multigrid
11905b89ad90SMark F. Adams 
11911cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
11921cc46a46SBarry Smith            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartitioning(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseASMAggs(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType()
11935b89ad90SMark F. Adams M*/
1194b2573a8aSBarry Smith 
11955b89ad90SMark F. Adams #undef __FUNCT__
11965b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
11978cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
11985b89ad90SMark F. Adams {
11995b89ad90SMark F. Adams   PetscErrorCode ierr;
12005b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
12015b89ad90SMark F. Adams   PC_MG          *mg;
12025b89ad90SMark F. Adams 
12035b89ad90SMark F. Adams   PetscFunctionBegin;
12041c1aac46SBarry Smith   /* register AMG type */
12051c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
12061c1aac46SBarry Smith 
12075b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
12081c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
12095b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
12105b89ad90SMark F. Adams 
12115b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1212b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
12135b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
12141c1aac46SBarry Smith   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally from PCMG by GAMG code */
12155b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
12165b89ad90SMark F. Adams 
1217b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
12181ab5ffc9SJed Brown 
12199d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
12209d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
12219d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
12229d5b6da9SMark F. Adams   pc_gamg->data    = 0;
12235b89ad90SMark F. Adams 
12249d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
12255b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
12265b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
12275b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
12285b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
12295adeb434SBarry Smith   mg->view                = PCView_GAMG;
12305b89ad90SMark F. Adams 
1231bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1232bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1233bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
12341cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1235bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1236bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1237bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1238c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1239bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
12409d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1241d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
1242*0c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1243038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
124425a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1245d3042614SMark Adams   pc_gamg->threshold        = 0.;
12469d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
12479ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
1248c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
12499d5b6da9SMark F. Adams 
125014a9496bSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)pc),&pc_gamg->random);CHKERRQ(ierr);
125114a9496bSBarry Smith 
1252bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1253bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
12545b89ad90SMark F. Adams   PetscFunctionReturn(0);
12555b89ad90SMark F. Adams }
12563e3471ccSMark Adams 
12573e3471ccSMark Adams #undef __FUNCT__
12583e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage"
12593e3471ccSMark Adams /*@C
12603e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
12613e3471ccSMark Adams     from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
12623e3471ccSMark Adams     when using static libraries.
12633e3471ccSMark Adams 
12643e3471ccSMark Adams  Level: developer
12653e3471ccSMark Adams 
12663e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
12673e3471ccSMark Adams  .seealso: PetscInitialize()
12683e3471ccSMark Adams @*/
12693e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
12703e3471ccSMark Adams {
12713e3471ccSMark Adams   PetscErrorCode ierr;
12723e3471ccSMark Adams 
12733e3471ccSMark Adams   PetscFunctionBegin;
12743e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
12753e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
12763e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
12773e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
12788e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
12793e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1280c1c463dbSMark Adams 
1281c1c463dbSMark Adams   /* general events */
1282fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1283fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1284fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1285fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1286c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1287c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1288fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1289c1c463dbSMark Adams 
12905b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12915b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
12925b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
12935b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
12945b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
12955b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
12965b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
12975b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
12985b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1299bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
13005b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
13015b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
13025b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
13035b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
13045b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
13055b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
13065b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
13075b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
13085b89ad90SMark F. Adams 
13095b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
13105b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
13115b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
13125b89ad90SMark F. Adams   /* create timer stages */
13135b89ad90SMark F. Adams #if defined GAMG_STAGES
13145b89ad90SMark F. Adams   {
13155b89ad90SMark F. Adams     char     str[32];
13165b89ad90SMark F. Adams     PetscInt lidx;
13175b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
13185b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
13195b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
13205b89ad90SMark F. Adams       sprintf(str,"MG Level %d",lidx);
13215b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
13225b89ad90SMark F. Adams     }
13235b89ad90SMark F. Adams   }
13245b89ad90SMark F. Adams #endif
13255b89ad90SMark F. Adams #endif
13263e3471ccSMark Adams   PetscFunctionReturn(0);
13273e3471ccSMark Adams }
13283e3471ccSMark Adams 
13293e3471ccSMark Adams #undef __FUNCT__
13303e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage"
13313e3471ccSMark Adams /*@C
13321c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
13331c1aac46SBarry Smith     called from PetscFinalize() automatically.
13343e3471ccSMark Adams 
13353e3471ccSMark Adams  Level: developer
13363e3471ccSMark Adams 
13373e3471ccSMark Adams  .keywords: Petsc, destroy, package
13383e3471ccSMark Adams  .seealso: PetscFinalize()
13393e3471ccSMark Adams @*/
13403e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
13413e3471ccSMark Adams {
13423e3471ccSMark Adams   PetscErrorCode ierr;
13433e3471ccSMark Adams 
13443e3471ccSMark Adams   PetscFunctionBegin;
13453e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
13463e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
13473e3471ccSMark Adams   PetscFunctionReturn(0);
13483e3471ccSMark Adams }
1349a36cf38bSToby Isaac 
1350a36cf38bSToby Isaac #undef __FUNCT__
1351a36cf38bSToby Isaac #define __FUNCT__ "PCGAMGRegister"
1352a36cf38bSToby Isaac /*@C
1353a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1354a36cf38bSToby Isaac 
1355a36cf38bSToby Isaac  Input Parameters:
1356a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1357a36cf38bSToby Isaac  - create - function for creating the gamg context.
1358a36cf38bSToby Isaac 
1359a36cf38bSToby Isaac   Level: advanced
1360a36cf38bSToby Isaac 
13611c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1362a36cf38bSToby Isaac @*/
1363a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1364a36cf38bSToby Isaac {
1365a36cf38bSToby Isaac   PetscErrorCode ierr;
1366a36cf38bSToby Isaac 
1367a36cf38bSToby Isaac   PetscFunctionBegin;
1368a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1369a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1370a36cf38bSToby Isaac   PetscFunctionReturn(0);
1371a36cf38bSToby Isaac }
1372a36cf38bSToby Isaac 
1373