xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision ce7c7f2f1fe04cf563f92500d5a0ea512c738ad7)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4af0996ceSBarry Smith #include <petsc/private/matimpl.h>
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
65b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */
7f96513f1SMatthew G Knepley 
80cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
90cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
10b4fbaa2aSMark F. Adams #endif
110cbbd2e1SMark F. Adams 
120cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
13fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG;
14fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO;
150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
19fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG;
200cbbd2e1SMark F. Adams #endif
210cbbd2e1SMark F. Adams 
22b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
230cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
24c1eae691SMark Adams static PetscLogStage gamg_stages[PETSC_GAMG_MAXLEVELS];
25b4fbaa2aSMark F. Adams #endif
26f96513f1SMatthew G Knepley 
27140e18c1SBarry Smith static PetscFunctionList GAMGList = 0;
283e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
299d5b6da9SMark F. Adams 
30d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
31d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
32d3d6bff4SMark F. Adams {
33d3d6bff4SMark F. Adams   PetscErrorCode ierr;
34d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
35d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
36d3d6bff4SMark F. Adams 
37d3d6bff4SMark F. Adams   PetscFunctionBegin;
3822a233eaSStefano Zampini   ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
391c1aac46SBarry Smith   pc_gamg->data_sz = 0;
40878e152fSMark F. Adams   ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
41a2f3521dSMark F. Adams   PetscFunctionReturn(0);
42a2f3521dSMark F. Adams }
43a2f3521dSMark F. Adams 
445b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
455b89ad90SMark F. Adams /*
46c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
47a147abb0SMark F. Adams      of active processors.
485b89ad90SMark F. Adams 
495b89ad90SMark F. Adams    Input Parameter:
50a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
51a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
529d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
53c5bfad50SMark F. Adams    . cr_bs - coarse block size
543530afc2SMark F. Adams    In/Output Parameter:
55a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
56afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
5711e60469SMark F. Adams    Output Parameter:
583530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
595b89ad90SMark F. Adams */
605cb416c2SMark F. Adams 
61171cca9aSMark Adams static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,IS * Pcolumnperm, PetscBool is_last)
625b89ad90SMark F. Adams {
63a2f3521dSMark F. Adams   PetscErrorCode  ierr;
649d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
65486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
66a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
673b4367a7SBarry Smith   MPI_Comm        comm;
68c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
693ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
705b89ad90SMark F. Adams 
715b89ad90SMark F. Adams   PetscFunctionBegin;
723b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
733b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
743b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
75c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
769d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
77038e3b61SMark F. Adams 
78*ce7c7f2fSMark Adams   if (Pcolumnperm) *Pcolumnperm = NULL;
79*ce7c7f2fSMark Adams 
803ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
810298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
823ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
833ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
8473911c69SBarry Smith   } else {
853ae0bb68SMark Adams     PetscInt  bs;
863ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
873ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
883ae0bb68SMark Adams   }
89c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
90171cca9aSMark Adams   if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) new_size = 1;
91171cca9aSMark Adams   else {
92472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
930298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
94a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
957f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
96c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
97a2f3521dSMark F. Adams   }
98f852f58cSMark F. Adams 
992e3501ffSMark Adams   if (new_size==nactive) {
100ef3f0257SMark Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
101*ce7c7f2fSMark Adams     if (new_size < size) {
102*ce7c7f2fSMark Adams       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
103*ce7c7f2fSMark Adams       ierr = PetscInfo1(pc,"reduced grid using same number of processors (%D) as last grid (use larger coarse grid)\n",nactive);CHKERRQ(ierr);
104*ce7c7f2fSMark Adams       if (pc_gamg->cpu_pin_coarse_grids) {
105*ce7c7f2fSMark Adams         ierr = MatPinToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr);
106*ce7c7f2fSMark Adams         ierr = MatPinToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr);
107*ce7c7f2fSMark Adams       }
108*ce7c7f2fSMark Adams     }
109ef3f0257SMark Adams     /* we know that the grid structure can be reused in MatPtAP */
1102e3501ffSMark Adams   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
111*ce7c7f2fSMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor;
112885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
11371959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
11471959b99SBarry 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);
1150cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1160cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
117b4fbaa2aSMark F. Adams #endif
118*ce7c7f2fSMark Adams     /* get new_size and rfactor */
119*ce7c7f2fSMark Adams     if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
120*ce7c7f2fSMark Adams       /* find factor */
121*ce7c7f2fSMark Adams       if (new_size == 1) rfactor = size; /* don't modify */
122*ce7c7f2fSMark Adams       else {
123*ce7c7f2fSMark Adams         PetscReal best_fact = 0.;
124*ce7c7f2fSMark Adams         jj = -1;
125*ce7c7f2fSMark Adams         for (kk = 1 ; kk <= size ; kk++) {
126*ce7c7f2fSMark Adams           if (!(size%kk)) { /* a candidate */
127*ce7c7f2fSMark Adams             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
128*ce7c7f2fSMark Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
129*ce7c7f2fSMark Adams             if (fact > best_fact) {
130*ce7c7f2fSMark Adams               best_fact = fact; jj = kk;
131*ce7c7f2fSMark Adams             }
132*ce7c7f2fSMark Adams           }
133*ce7c7f2fSMark Adams         }
134*ce7c7f2fSMark Adams         if (jj != -1) rfactor = jj;
135*ce7c7f2fSMark Adams         else rfactor = 1; /* a prime */
136*ce7c7f2fSMark Adams         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
137*ce7c7f2fSMark Adams         else expand_factor = rfactor;
138*ce7c7f2fSMark Adams       }
139*ce7c7f2fSMark Adams       new_size = size/rfactor; /* make new size one that is factor */
140*ce7c7f2fSMark Adams       if (new_size==nactive) {
141*ce7c7f2fSMark Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
142*ce7c7f2fSMark Adams         ierr = PetscInfo2(pc,"Finding factorable processor set stopped reduction: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr);
143*ce7c7f2fSMark Adams #if defined PETSC_GAMG_USE_LOG
144*ce7c7f2fSMark Adams         ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
145*ce7c7f2fSMark Adams #endif
146*ce7c7f2fSMark Adams         PetscFunctionReturn(0);
147*ce7c7f2fSMark Adams       }
148*ce7c7f2fSMark Adams     }
149a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
150785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
1512e3501ffSMark Adams     if (pc_gamg->repart) {
152a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1535a9b9e01SMark F. Adams       Mat      adj;
154*ce7c7f2fSMark Adams       ierr = PetscInfo4(pc,"Repartition: size (active): %D --> %D, %D local equations, using %s process layout\n",*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread");CHKERRQ(ierr);
155a2f3521dSMark F. Adams       /* get 'adj' */
156c5bfad50SMark F. Adams       if (cr_bs == 1) {
157038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
158806fa848SBarry Smith       } else {
159a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
160eb07cef2SMark F. Adams         Mat               tMat;
161a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
162b4fbaa2aSMark F. Adams         const PetscScalar *vals;
163b4fbaa2aSMark F. Adams         const PetscInt    *idx;
164a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1659057884aSMark F. Adams         static PetscInt   llev = 0;
166d9558ea9SBarry Smith         MatType           mtype;
167b4fbaa2aSMark F. Adams 
168e632b94dSBarry Smith         ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr);
169a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
170a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
171c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
17258471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
173c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
174c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
17558471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1763ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1773ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
17858471d46SMark F. Adams         }
1796876a03eSMark F. Adams 
180d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
1813b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1823ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
183d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
184a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
185a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
186e632b94dSBarry Smith         ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
187eb07cef2SMark F. Adams 
188a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
189c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
19022063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
191eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
192c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
193eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
194eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
195eb07cef2SMark F. Adams           }
19622063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
197eb07cef2SMark F. Adams         }
198eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
199eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
200eb07cef2SMark F. Adams 
201b4fbaa2aSMark F. Adams         if (llev++ == -1) {
202b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2038caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
2043b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
205b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2063bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
207b4fbaa2aSMark F. Adams         }
208eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
209eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
210a2f3521dSMark F. Adams       } /* create 'adj' */
211f150b916SMark F. Adams 
212a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2135a9b9e01SMark F. Adams         char            prefix[256];
2145a9b9e01SMark F. Adams         const char      *pcpre;
215b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
216b4fbaa2aSMark F. Adams         MatPartitioning mpart;
217a4b7d37bSMark F. Adams         IS              proc_is;
2182f03bc48SMark F. Adams 
2193b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2205ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2219d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2228caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
22359a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
22411e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
225c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
226a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
22711e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2285a9b9e01SMark F. Adams 
2295ef31b24SMark F. Adams         /* collect IS info */
230785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
231a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
232a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
233c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
234*ce7c7f2fSMark Adams             newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
235eb07cef2SMark F. Adams           }
2365ef31b24SMark F. Adams         }
237a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
238a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2395ef31b24SMark F. Adams       }
2405ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2415a9b9e01SMark F. Adams 
2423b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2438263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
244806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
245*ce7c7f2fSMark Adams       PetscInt targetPE;
246c5df96a5SBarry Smith       if (new_size==nactive) {
247a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2485a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
249302440fdSBarry Smith         ierr = PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr);
250fe682e87SBarry Smith #if defined PETSC_GAMG_USE_LOG
251fe682e87SBarry Smith         ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
252fe682e87SBarry Smith #endif
2535a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2545a9b9e01SMark F. Adams       }
255302440fdSBarry Smith       ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr);
256*ce7c7f2fSMark Adams       targetPE = (rank/rfactor)*expand_factor;
2573b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
258a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
259e33ef3b1SMark F. Adams 
26011e60469SMark F. Adams     /*
261a2f3521dSMark F. Adams       Create an index set from the is_eq_newproc index set to indicate the mapping TO
26211e60469SMark F. Adams     */
263a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2647700e67bSMark Adams     is_eq_num_prim = is_eq_num;
26511e60469SMark F. Adams     /*
266a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
26711e60469SMark F. Adams     */
268c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
269c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
270a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
271*ce7c7f2fSMark Adams     ncrs_new = ncrs_eq_new/cr_bs;
272a2f3521dSMark F. Adams 
273a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2740cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2750cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
276b4fbaa2aSMark F. Adams #endif
277885364a3SMark 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 */
278885364a3SMark Adams     {
279885364a3SMark Adams       Vec            src_crd, dest_crd;
280885364a3SMark 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;
281885364a3SMark Adams       VecScatter     vecscat;
282885364a3SMark Adams       PetscScalar    *array;
283885364a3SMark Adams       IS isscat;
284a2f3521dSMark F. Adams       /* move data (for primal equations only) */
28522063be5SMark F. Adams       /* Create a vector to contain the newly ordered element information */
2863b4367a7SBarry Smith       ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
2873ae0bb68SMark Adams       ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
288c0dedaeaSBarry Smith       ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
28911e60469SMark F. Adams       /*
2909d5b6da9SMark F. Adams 	There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
291c5bfad50SMark F. Adams 	a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
29211e60469SMark F. Adams       */
293854ce69bSBarry Smith       ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
294a2f3521dSMark F. Adams       ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2953ae0bb68SMark Adams       for (ii=0,jj=0; ii<ncrs; ii++) {
296c5bfad50SMark F. Adams 	PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
297a2f3521dSMark F. Adams 	for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
29811e60469SMark F. Adams       }
299a2f3521dSMark F. Adams       ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3003ae0bb68SMark Adams       ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
30192a756f0SMark F. Adams       ierr = PetscFree(tidx);CHKERRQ(ierr);
30211e60469SMark F. Adams       /*
30311e60469SMark F. Adams 	Create a vector to contain the original vertex information for each element
30411e60469SMark F. Adams       */
3053ae0bb68SMark Adams       ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
3069d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3073ae0bb68SMark Adams 	const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3083ae0bb68SMark Adams 	for (ii=0; ii<ncrs; ii++) {
3099d5b6da9SMark F. Adams 	  for (kk=0; kk<ndata_rows; kk++) {
310a2f3521dSMark F. Adams 	    PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
311c8b0795cSMark F. Adams 	    PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
312676e1743SMark F. Adams 	    ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
313d3d6bff4SMark F. Adams 	  }
314038e3b61SMark F. Adams 	}
315eb07cef2SMark F. Adams       }
316eb07cef2SMark F. Adams       ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
317eb07cef2SMark F. Adams       ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
31811e60469SMark F. Adams       /*
31911e60469SMark F. Adams 	Scatter the element vertex information (still in the original vertex ordering)
32011e60469SMark F. Adams 	to the correct processor
32111e60469SMark F. Adams       */
3229448b7f1SJunchao Zhang       ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
32311e60469SMark F. Adams       ierr = ISDestroy(&isscat);CHKERRQ(ierr);
32411e60469SMark F. Adams       ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32511e60469SMark F. Adams       ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32611e60469SMark F. Adams       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
32711e60469SMark F. Adams       ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
32811e60469SMark F. Adams       /*
32911e60469SMark F. Adams 	Put the element vertex data into a new allocation of the gdata->ele
33011e60469SMark F. Adams       */
331c8b0795cSMark F. Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
332578f55a3SPeter Brune       ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3332fa5cd67SKarl Rupp 
3343ae0bb68SMark Adams       pc_gamg->data_sz = node_data_sz*ncrs_new;
3353ae0bb68SMark Adams       strideNew        = ncrs_new*ndata_rows;
3362fa5cd67SKarl Rupp 
33711e60469SMark F. Adams       ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3389d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3393ae0bb68SMark Adams 	for (ii=0; ii<ncrs_new; ii++) {
3409d5b6da9SMark F. Adams 	  for (kk=0; kk<ndata_rows; kk++) {
341a2f3521dSMark F. Adams 	    PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
342c8b0795cSMark F. Adams 	    pc_gamg->data[ix] = PetscRealPart(array[jx]);
343d3d6bff4SMark F. Adams 	  }
344038e3b61SMark F. Adams 	}
345038e3b61SMark F. Adams       }
34611e60469SMark F. Adams       ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
34711e60469SMark F. Adams       ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
348885364a3SMark Adams     }
349a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3500cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3510cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
352ed3f9983SMark F. Adams #endif
35311e60469SMark F. Adams     /*
3547dae84e0SHong Zhang       Invert for MatCreateSubMatrix
35511e60469SMark F. Adams     */
356a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
357a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
358c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
359a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
360a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
361a2f3521dSMark F. Adams     }
3623cb8563fSToby Isaac     if (Pcolumnperm) {
3633cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3643cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3653cb8563fSToby Isaac     }
366a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3670cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3680cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3690cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
370ed3f9983SMark F. Adams #endif
371a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
372a2f3521dSMark F. Adams     {
373a2f3521dSMark F. Adams       Mat mat;
3747dae84e0SHong Zhang       ierr        = MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
375a2f3521dSMark F. Adams       *a_Amat_crs = mat;
376a2f3521dSMark F. Adams     }
377038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
378a2f3521dSMark F. Adams 
3790cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3800cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
381ed3f9983SMark F. Adams #endif
38211e60469SMark F. Adams     /* prolongator */
38311e60469SMark F. Adams     {
38411e60469SMark F. Adams       IS       findices;
385a2f3521dSMark F. Adams       PetscInt Istart,Iend;
386a2f3521dSMark F. Adams       Mat      Pnew;
38762294041SBarry Smith 
388a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
3890cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3900cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
391ed3f9983SMark F. Adams #endif
3923b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
393c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
3947dae84e0SHong Zhang       ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
39511e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
396c5bfad50SMark F. Adams 
3970cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3980cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
399ed3f9983SMark F. Adams #endif
4003530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
401a2f3521dSMark F. Adams 
402a2f3521dSMark F. Adams       /* output - repartitioned */
403a2f3521dSMark F. Adams       *a_P_inout = Pnew;
404e33ef3b1SMark F. Adams     }
405a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4065b89ad90SMark F. Adams 
407c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
408*ce7c7f2fSMark Adams 
409*ce7c7f2fSMark Adams     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
410*ce7c7f2fSMark Adams     if (pc_gamg->cpu_pin_coarse_grids) {
411*ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
412*ce7c7f2fSMark Adams       ierr = PetscInfo(pc,"Pinning level to the CPU\n");CHKERRQ(ierr);
413*ce7c7f2fSMark Adams #endif
414*ce7c7f2fSMark Adams       ierr = MatPinToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr);
415*ce7c7f2fSMark Adams       ierr = MatPinToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr);
416*ce7c7f2fSMark Adams       if (1) { /* lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
417*ce7c7f2fSMark Adams         Mat         A = *a_Amat_crs, P = *a_P_inout;
418*ce7c7f2fSMark Adams         PetscMPIInt size;
419*ce7c7f2fSMark Adams         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
420*ce7c7f2fSMark Adams         if (size > 1) {
421*ce7c7f2fSMark Adams           Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data;
422*ce7c7f2fSMark Adams           ierr = VecPinToCPU(a->lvec,PETSC_TRUE);CHKERRQ(ierr);
423*ce7c7f2fSMark Adams           ierr = VecPinToCPU(p->lvec,PETSC_TRUE);CHKERRQ(ierr);
424*ce7c7f2fSMark Adams         }
425*ce7c7f2fSMark Adams       }
426*ce7c7f2fSMark Adams     }
427a2f3521dSMark F. Adams   }
4285b89ad90SMark F. Adams   PetscFunctionReturn(0);
4295b89ad90SMark F. Adams }
4305b89ad90SMark F. Adams 
4315b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4325b89ad90SMark F. Adams /*
4335b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4345b89ad90SMark F. Adams                     by setting data structures and options.
4355b89ad90SMark F. Adams 
4365b89ad90SMark F. Adams    Input Parameter:
4375b89ad90SMark F. Adams .  pc - the preconditioner context
4385b89ad90SMark F. Adams 
4395b89ad90SMark F. Adams */
4409d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4415b89ad90SMark F. Adams {
4425b89ad90SMark F. Adams   PetscErrorCode ierr;
4439d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4445b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4452adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
446c1eae691SMark Adams   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_GAMG_MAXLEVELS];
4473b4367a7SBarry Smith   MPI_Comm       comm;
448c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
449c1eae691SMark Adams   Mat            Aarr[PETSC_GAMG_MAXLEVELS],Parr[PETSC_GAMG_MAXLEVELS];
450c1eae691SMark Adams   IS             *ASMLocalIDsArr[PETSC_GAMG_MAXLEVELS];
451a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
452569f4572SMark Adams   MatInfo        info;
453171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
4545ef31b24SMark F. Adams 
4555b89ad90SMark F. Adams   PetscFunctionBegin;
4563b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4573b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4583b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
459dfd5c07aSMark F. Adams 
46084d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
4611c1aac46SBarry Smith     if ((PetscBool)(!pc_gamg->reuse_prol)) {
462878e152fSMark F. Adams       /* reset everything */
463878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
464878e152fSMark F. Adams       pc->setupcalled = 0;
465806fa848SBarry Smith     } else {
46684d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
46703a628feSMark F. Adams       /* just do Galerkin grids */
46858471d46SMark F. Adams       Mat          B,dA,dB;
46958471d46SMark F. Adams 
47071959b99SBarry Smith       if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
4719d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
47258471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
47323ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
47458471d46SMark F. Adams         /* (re)set to get dirty flag */
47523ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
47658471d46SMark F. Adams 
4772fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
478ef3f0257SMark Adams           /* 2nd solve, matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
479ef3f0257SMark Adams           if (pc_gamg->setup_count==2 /* && pc_gamg->repart||reduction */) {
480ef3f0257SMark Adams             ierr = PetscInfo2(pc,"new RAP after first solve level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr);
481ef3f0257SMark Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,2.0,&B);CHKERRQ(ierr);
482084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
48303a628feSMark F. Adams             mglevels[level]->A = B;
484806fa848SBarry Smith           } else {
485ef3f0257SMark Adams             ierr = PetscInfo2(pc,"RAP after first solve reusing matrix level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr);
48623ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
48758471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
48803a628feSMark F. Adams           }
48923ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
49058471d46SMark F. Adams           dB   = B;
49158471d46SMark F. Adams         }
4925f8cf99dSMark F. Adams       }
493d5280255SMark F. Adams 
494d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
49558471d46SMark F. Adams       PetscFunctionReturn(0);
496eb07cef2SMark F. Adams     }
497878e152fSMark F. Adams   }
498f6536408SMark F. Adams 
499878e152fSMark F. Adams   if (!pc_gamg->data) {
500878e152fSMark F. Adams     if (pc_gamg->orig_data) {
501878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5020298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5032fa5cd67SKarl Rupp 
504878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
505878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
506878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5072fa5cd67SKarl Rupp 
508785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
509878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
510806fa848SBarry Smith     } else {
5111ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5127700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5139d5b6da9SMark F. Adams     }
514878e152fSMark F. Adams   }
515878e152fSMark F. Adams 
516878e152fSMark F. Adams   /* cache original data for reuse */
5171c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
518785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
519878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
520878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
521878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
522878e152fSMark F. Adams   }
523038e3b61SMark F. Adams 
524302f38e8SMark F. Adams   /* get basic dims */
525302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
526171cca9aSMark Adams   ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr);
52784d3f75bSMark F. Adams 
528569f4572SMark Adams   ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
529569f4572SMark Adams   nnz0   = info.nz_used;
530569f4572SMark Adams   nnztot = info.nz_used;
53162294041SBarry 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);
532569f4572SMark Adams 
533a2f3521dSMark F. Adams   /* Get A_i and R_i */
53462294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
5359ab59c8bSMark Adams     pc_gamg->current_level = level;
5365b89ad90SMark F. Adams     level1 = level + 1;
5370cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5380cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
539a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
540a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
541b4fbaa2aSMark F. Adams #endif
542a2f3521dSMark F. Adams #endif
543c8b0795cSMark F. Adams     { /* construct prolongator */
544725b86d8SJed Brown       Mat              Gmat;
5450cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5467700e67bSMark Adams       Mat              Prol11;
547c8b0795cSMark F. Adams 
5487700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5491ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5507700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
551c8b0795cSMark F. Adams 
552a2f3521dSMark F. Adams       /* could have failed to create new level */
553a2f3521dSMark F. Adams       if (Prol11) {
5549d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5550298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
556a2f3521dSMark F. Adams 
557fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
558c8b0795cSMark F. Adams           /* smooth */
559fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
560c8b0795cSMark F. Adams         }
561c8b0795cSMark F. Adams 
5620c3bc534SBarry Smith         if (pc_gamg->use_aggs_in_asm) {
5631b18a24aSMark Adams           PetscInt bs;
5641b18a24aSMark Adams           ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
5650a3c815dSMark Adams           ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
566ffc955d6SMark F. Adams         }
567ffc955d6SMark F. Adams 
5684bde40a0SMark Adams         Parr[level1] = Prol11;
5694bde40a0SMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
5704bde40a0SMark Adams 
571a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
57241b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
573a2f3521dSMark F. Adams     } /* construct prolongator scope */
5740cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5750cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
576c8b0795cSMark F. Adams #endif
5777f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
578171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
579569f4572SMark Adams       ierr =  PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr);
58062294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
581a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
582a90e85d9SMark Adams #endif
583c8b0795cSMark F. Adams       break;
584c8b0795cSMark F. Adams     }
5850cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5860cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
587b4fbaa2aSMark F. Adams #endif
588171cca9aSMark Adams     ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */
589171cca9aSMark Adams     if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ????????");
590171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
5910e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
592171cca9aSMark Adams     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr);
593a2f3521dSMark F. Adams 
5940cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5950cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
596b4fbaa2aSMark F. Adams #endif
597171cca9aSMark Adams     ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */
598569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
599569f4572SMark Adams     nnztot += info.nz_used;
6001d5b2942SMark 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);
601569f4572SMark Adams 
6020cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
603b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
604b4fbaa2aSMark F. Adams #endif
605a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
6069ab59c8bSMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) {
6079ab59c8bSMark Adams       ierr =  PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr);
608a90e85d9SMark Adams       level++;
609a90e85d9SMark Adams       break;
610a90e85d9SMark Adams     }
611c8b0795cSMark F. Adams   } /* levels */
612c8b0795cSMark F. Adams   ierr                  = PetscFree(pc_gamg->data);CHKERRQ(ierr);
613c8b0795cSMark F. Adams 
614569f4572SMark Adams   ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr);
6159d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6165b89ad90SMark F. Adams   fine_level       = level;
6170298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
6185b89ad90SMark F. Adams 
61962294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
620d5280255SMark F. Adams     /* set default smoothers & set operators */
62162294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
622ffc955d6SMark F. Adams       KSP smoother;
623ffc955d6SMark F. Adams       PC  subpc;
624a2f3521dSMark F. Adams 
6259d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
626f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
627ffc955d6SMark F. Adams 
628a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
629a2f3521dSMark F. Adams       /* set ops */
63023ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
631a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
632a2f3521dSMark F. Adams 
633a2f3521dSMark F. Adams       /* set defaults */
6346c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
635a2f3521dSMark F. Adams 
6360c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
6370c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
6382d3561bbSSatish Balay         PetscInt sz;
6397a28f3e5SMark Adams         IS       *iss;
640a2f3521dSMark F. Adams 
6412d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
6427a28f3e5SMark Adams         iss  = ASMLocalIDsArr[level];
6430c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
6440a3c815dSMark Adams         ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr);
6450c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
6467f66b68fSMark Adams         if (!sz) {
647ffc955d6SMark F. Adams           IS       is;
6480a3c815dSMark Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
6497a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
650a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
651806fa848SBarry Smith         } else {
652a94c3b12SMark F. Adams           PetscInt kk;
6537a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr);
654a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
6557a28f3e5SMark Adams             ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr);
656a94c3b12SMark F. Adams           }
6577a28f3e5SMark Adams           ierr = PetscFree(iss);CHKERRQ(ierr);
658ffc955d6SMark F. Adams         }
6590298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
660ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
661806fa848SBarry Smith       } else {
662890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
663ffc955d6SMark F. Adams       }
664d5280255SMark F. Adams     }
665d5280255SMark F. Adams     {
666d5280255SMark F. Adams       /* coarse grid */
667d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
668d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
669d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
67023ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
671cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
672d5280255SMark F. Adams         ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
673d5280255SMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
674d5280255SMark F. Adams         ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
675d5280255SMark F. Adams         ierr = PCSetUp(subpc);CHKERRQ(ierr);
67671959b99SBarry Smith         ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
67771959b99SBarry Smith         if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
678d5280255SMark F. Adams         ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
679d5280255SMark F. Adams         ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
6809dbfc187SHong Zhang         ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
6812fb0b348SMark F. Adams         ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
68208e36f19SMark Adams         ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
6835b42dca8SJed Brown         /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
6845b42dca8SJed Brown          * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
6855b42dca8SJed Brown          * view every subdomain as though they were different. */
6865b42dca8SJed Brown         ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
687d5280255SMark F. Adams       }
688cf8ae1d3SMark Adams     }
689d5280255SMark F. Adams 
690d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
691d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
692e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
693d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
69469aca0b8SBarry Smith     ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
695d5280255SMark F. Adams 
696d5280255SMark F. Adams     /* clean up */
697d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
698587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
699587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
7005b89ad90SMark F. Adams     }
7010cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
702806fa848SBarry Smith   } else {
7035f8cf99dSMark F. Adams     KSP smoother;
704302440fdSBarry Smith     ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
7059d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
70623ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
7075f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
7089d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
7095f8cf99dSMark F. Adams   }
7105b89ad90SMark F. Adams   PetscFunctionReturn(0);
7115b89ad90SMark F. Adams }
7125b89ad90SMark F. Adams 
713eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7145b89ad90SMark F. Adams /*
7155b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7165b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7175b89ad90SMark F. Adams 
7185b89ad90SMark F. Adams    Input Parameter:
7195b89ad90SMark F. Adams .  pc - the preconditioner context
7205b89ad90SMark F. Adams 
7215b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7225b89ad90SMark F. Adams */
7235b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
7245b89ad90SMark F. Adams {
7255b89ad90SMark F. Adams   PetscErrorCode ierr;
7265b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
7275b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
7285b89ad90SMark F. Adams 
7295b89ad90SMark F. Adams   PetscFunctionBegin;
7305b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7319b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
7329b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
7339b8ffb57SJed Brown   }
7341ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
7351ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
7365b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7375b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7385b89ad90SMark F. Adams   PetscFunctionReturn(0);
7395b89ad90SMark F. Adams }
7405b89ad90SMark F. Adams 
741676e1743SMark F. Adams /*@
742cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
743676e1743SMark F. Adams 
7441cc46a46SBarry Smith    Logically Collective on PC
745676e1743SMark F. Adams 
746676e1743SMark F. Adams    Input Parameters:
7471cc46a46SBarry Smith +  pc - the preconditioner context
7481cc46a46SBarry Smith -  n - the number of equations
749676e1743SMark F. Adams 
750676e1743SMark F. Adams 
751676e1743SMark F. Adams    Options Database Key:
7521cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
753676e1743SMark F. Adams 
75495452b02SPatrick Sanan    Notes:
75595452b02SPatrick Sanan     GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
756cab9ed1eSBarry Smith           that has degrees of freedom
757cab9ed1eSBarry Smith 
758676e1743SMark F. Adams    Level: intermediate
759676e1743SMark F. Adams 
7601c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim()
761676e1743SMark F. Adams @*/
762676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
763676e1743SMark F. Adams {
764676e1743SMark F. Adams   PetscErrorCode ierr;
765676e1743SMark F. Adams 
766676e1743SMark F. Adams   PetscFunctionBegin;
767676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
768676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
769676e1743SMark F. Adams   PetscFunctionReturn(0);
770676e1743SMark F. Adams }
771676e1743SMark F. Adams 
7721e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
773676e1743SMark F. Adams {
774c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
775c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
776676e1743SMark F. Adams 
777676e1743SMark F. Adams   PetscFunctionBegin;
7789d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
779676e1743SMark F. Adams   PetscFunctionReturn(0);
780676e1743SMark F. Adams }
781676e1743SMark F. Adams 
782389730f3SMark F. Adams /*@
783cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
784389730f3SMark F. Adams 
785389730f3SMark F. Adams  Collective on PC
786389730f3SMark F. Adams 
787389730f3SMark F. Adams    Input Parameters:
7881cc46a46SBarry Smith +  pc - the preconditioner context
7891cc46a46SBarry Smith -  n - maximum number of equations to aim for
790389730f3SMark F. Adams 
791389730f3SMark F. Adams    Options Database Key:
7921cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
793389730f3SMark F. Adams 
79474329af1SBarry Smith    Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
79574329af1SBarry Smith      has less than 1000 unknowns.
79674329af1SBarry Smith 
797389730f3SMark F. Adams    Level: intermediate
798389730f3SMark F. Adams 
7991c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim()
800389730f3SMark F. Adams @*/
801389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
802389730f3SMark F. Adams {
803389730f3SMark F. Adams   PetscErrorCode ierr;
804389730f3SMark F. Adams 
805389730f3SMark F. Adams   PetscFunctionBegin;
806389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
807389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
808389730f3SMark F. Adams   PetscFunctionReturn(0);
809389730f3SMark F. Adams }
810389730f3SMark F. Adams 
8111e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
812389730f3SMark F. Adams {
813389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
814389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
815389730f3SMark F. Adams 
816389730f3SMark F. Adams   PetscFunctionBegin;
8179d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
818389730f3SMark F. Adams   PetscFunctionReturn(0);
819389730f3SMark F. Adams }
820389730f3SMark F. Adams 
821676e1743SMark F. Adams /*@
822cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
823676e1743SMark F. Adams 
824676e1743SMark F. Adams    Collective on PC
825676e1743SMark F. Adams 
826676e1743SMark F. Adams    Input Parameters:
8271cc46a46SBarry Smith +  pc - the preconditioner context
8281cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
829676e1743SMark F. Adams 
830676e1743SMark F. Adams    Options Database Key:
8311cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
832676e1743SMark F. Adams 
83395452b02SPatrick Sanan    Notes:
83495452b02SPatrick Sanan     this will generally improve the loading balancing of the work on each level
835cab9ed1eSBarry Smith 
836676e1743SMark F. Adams    Level: intermediate
837676e1743SMark F. Adams 
838676e1743SMark F. Adams .seealso: ()
839676e1743SMark F. Adams @*/
840cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
841676e1743SMark F. Adams {
842676e1743SMark F. Adams   PetscErrorCode ierr;
843676e1743SMark F. Adams 
844676e1743SMark F. Adams   PetscFunctionBegin;
845676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
846cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
847676e1743SMark F. Adams   PetscFunctionReturn(0);
848676e1743SMark F. Adams }
849676e1743SMark F. Adams 
850cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
851676e1743SMark F. Adams {
852c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
853c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
854676e1743SMark F. Adams 
855676e1743SMark F. Adams   PetscFunctionBegin;
8569d5b6da9SMark F. Adams   pc_gamg->repart = n;
857676e1743SMark F. Adams   PetscFunctionReturn(0);
858676e1743SMark F. Adams }
859676e1743SMark F. Adams 
860dfd5c07aSMark F. Adams /*@
861cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
862dfd5c07aSMark F. Adams 
863dfd5c07aSMark F. Adams    Collective on PC
864dfd5c07aSMark F. Adams 
865dfd5c07aSMark F. Adams    Input Parameters:
8661cc46a46SBarry Smith +  pc - the preconditioner context
8671cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
868dfd5c07aSMark F. Adams 
869dfd5c07aSMark F. Adams    Options Database Key:
8701cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
871dfd5c07aSMark F. Adams 
872dfd5c07aSMark F. Adams    Level: intermediate
873dfd5c07aSMark F. Adams 
87495452b02SPatrick Sanan    Notes:
87595452b02SPatrick Sanan     this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
876cab9ed1eSBarry Smith           rebuilding the preconditioner quicker.
877cab9ed1eSBarry Smith 
878dfd5c07aSMark F. Adams .seealso: ()
879dfd5c07aSMark F. Adams @*/
8801cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
881dfd5c07aSMark F. Adams {
882dfd5c07aSMark F. Adams   PetscErrorCode ierr;
883dfd5c07aSMark F. Adams 
884dfd5c07aSMark F. Adams   PetscFunctionBegin;
885dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8861cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
887dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
888dfd5c07aSMark F. Adams }
889dfd5c07aSMark F. Adams 
8901cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
891dfd5c07aSMark F. Adams {
892dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
893dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
894dfd5c07aSMark F. Adams 
895dfd5c07aSMark F. Adams   PetscFunctionBegin;
896dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
897dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
898dfd5c07aSMark F. Adams }
899dfd5c07aSMark F. Adams 
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 
914ffc955d6SMark F. Adams .seealso: ()
915ffc955d6SMark F. Adams @*/
916cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
917ffc955d6SMark F. Adams {
918ffc955d6SMark F. Adams   PetscErrorCode ierr;
919ffc955d6SMark F. Adams 
920ffc955d6SMark F. Adams   PetscFunctionBegin;
921ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
922cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
923ffc955d6SMark F. Adams   PetscFunctionReturn(0);
924ffc955d6SMark F. Adams }
925ffc955d6SMark F. Adams 
926cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
927ffc955d6SMark F. Adams {
928ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
929ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
930ffc955d6SMark F. Adams 
931ffc955d6SMark F. Adams   PetscFunctionBegin;
932cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
933ffc955d6SMark F. Adams   PetscFunctionReturn(0);
934ffc955d6SMark F. Adams }
935ffc955d6SMark F. Adams 
936171cca9aSMark Adams /*@
937cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
938171cca9aSMark Adams 
939171cca9aSMark Adams    Collective on PC
940171cca9aSMark Adams 
941171cca9aSMark Adams    Input Parameters:
942171cca9aSMark Adams +  pc - the preconditioner context
943cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
944171cca9aSMark Adams 
945171cca9aSMark Adams    Options Database Key:
946cf8ae1d3SMark Adams .  -pc_gamg_use_parallel_coarse_grid_solver
947171cca9aSMark Adams 
948171cca9aSMark Adams    Level: intermediate
949171cca9aSMark Adams 
950171cca9aSMark Adams .seealso: ()
951171cca9aSMark Adams @*/
952171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
953171cca9aSMark Adams {
954171cca9aSMark Adams   PetscErrorCode ierr;
955171cca9aSMark Adams 
956171cca9aSMark Adams   PetscFunctionBegin;
957171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
958171cca9aSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
959171cca9aSMark Adams   PetscFunctionReturn(0);
960171cca9aSMark Adams }
961171cca9aSMark Adams 
962171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
963171cca9aSMark Adams {
964171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
965171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
966171cca9aSMark Adams 
967171cca9aSMark Adams   PetscFunctionBegin;
968171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
969ffc955d6SMark F. Adams   PetscFunctionReturn(0);
970ffc955d6SMark F. Adams }
971ffc955d6SMark F. Adams 
9724ef23d27SMark F. Adams /*@
973*ce7c7f2fSMark Adams    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU
974*ce7c7f2fSMark Adams 
975*ce7c7f2fSMark Adams    Collective on PC
976*ce7c7f2fSMark Adams 
977*ce7c7f2fSMark Adams    Input Parameters:
978*ce7c7f2fSMark Adams +  pc - the preconditioner context
979*ce7c7f2fSMark Adams -  flg - PETSC_TRUE to pin coarse grids to CPU
980*ce7c7f2fSMark Adams 
981*ce7c7f2fSMark Adams    Options Database Key:
982*ce7c7f2fSMark Adams .  -pc_gamg_cpu_pin_coarse_grids
983*ce7c7f2fSMark Adams 
984*ce7c7f2fSMark Adams    Level: intermediate
985*ce7c7f2fSMark Adams 
986*ce7c7f2fSMark Adams .seealso: ()
987*ce7c7f2fSMark Adams @*/
988*ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
989*ce7c7f2fSMark Adams {
990*ce7c7f2fSMark Adams   PetscErrorCode ierr;
991*ce7c7f2fSMark Adams 
992*ce7c7f2fSMark Adams   PetscFunctionBegin;
993*ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
994*ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
995*ce7c7f2fSMark Adams   PetscFunctionReturn(0);
996*ce7c7f2fSMark Adams }
997*ce7c7f2fSMark Adams 
998*ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
999*ce7c7f2fSMark Adams {
1000*ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1001*ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1002*ce7c7f2fSMark Adams 
1003*ce7c7f2fSMark Adams   PetscFunctionBegin;
1004*ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = flg;
1005*ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1006*ce7c7f2fSMark Adams }
1007*ce7c7f2fSMark Adams 
1008*ce7c7f2fSMark Adams /*@
1009*ce7c7f2fSMark Adams    PCGAMGSetCoarseGridLayoutType - place reduce grids on processors with natural order (compact type)
1010*ce7c7f2fSMark Adams 
1011*ce7c7f2fSMark Adams    Collective on PC
1012*ce7c7f2fSMark Adams 
1013*ce7c7f2fSMark Adams    Input Parameters:
1014*ce7c7f2fSMark Adams +  pc - the preconditioner context
1015*ce7c7f2fSMark Adams -  flg - Layout type
1016*ce7c7f2fSMark Adams 
1017*ce7c7f2fSMark Adams    Options Database Key:
1018*ce7c7f2fSMark Adams .  -pc_gamg_coarse_grid_layout_type
1019*ce7c7f2fSMark Adams 
1020*ce7c7f2fSMark Adams    Level: intermediate
1021*ce7c7f2fSMark Adams 
1022*ce7c7f2fSMark Adams .seealso: ()
1023*ce7c7f2fSMark Adams @*/
1024*ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1025*ce7c7f2fSMark Adams {
1026*ce7c7f2fSMark Adams   PetscErrorCode ierr;
1027*ce7c7f2fSMark Adams 
1028*ce7c7f2fSMark Adams   PetscFunctionBegin;
1029*ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1030*ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));CHKERRQ(ierr);
1031*ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1032*ce7c7f2fSMark Adams }
1033*ce7c7f2fSMark Adams 
1034*ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1035*ce7c7f2fSMark Adams {
1036*ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1037*ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1038*ce7c7f2fSMark Adams 
1039*ce7c7f2fSMark Adams   PetscFunctionBegin;
1040*ce7c7f2fSMark Adams   pc_gamg->layout_type = flg;
1041*ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1042*ce7c7f2fSMark Adams }
1043*ce7c7f2fSMark Adams 
1044*ce7c7f2fSMark Adams /*@
10451cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
10464ef23d27SMark F. Adams 
10474ef23d27SMark F. Adams    Not collective on PC
10484ef23d27SMark F. Adams 
10494ef23d27SMark F. Adams    Input Parameters:
10501cc46a46SBarry Smith +  pc - the preconditioner
10511cc46a46SBarry Smith -  n - the maximum number of levels to use
10524ef23d27SMark F. Adams 
10534ef23d27SMark F. Adams    Options Database Key:
10544ef23d27SMark F. Adams .  -pc_mg_levels
10554ef23d27SMark F. Adams 
10564ef23d27SMark F. Adams    Level: intermediate
10574ef23d27SMark F. Adams 
10584ef23d27SMark F. Adams .seealso: ()
10594ef23d27SMark F. Adams @*/
10604ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
10614ef23d27SMark F. Adams {
10624ef23d27SMark F. Adams   PetscErrorCode ierr;
10634ef23d27SMark F. Adams 
10644ef23d27SMark F. Adams   PetscFunctionBegin;
10654ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10664ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
10674ef23d27SMark F. Adams   PetscFunctionReturn(0);
10684ef23d27SMark F. Adams }
10694ef23d27SMark F. Adams 
10701e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
10714ef23d27SMark F. Adams {
10724ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
10734ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
10744ef23d27SMark F. Adams 
10754ef23d27SMark F. Adams   PetscFunctionBegin;
10769d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
10774ef23d27SMark F. Adams   PetscFunctionReturn(0);
10784ef23d27SMark F. Adams }
10794ef23d27SMark F. Adams 
10803542efc5SMark F. Adams /*@
10813542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
10823542efc5SMark F. Adams 
10833542efc5SMark F. Adams    Not collective on PC
10843542efc5SMark F. Adams 
10853542efc5SMark F. Adams    Input Parameters:
10861cc46a46SBarry Smith +  pc - the preconditioner context
1087b001cb0fSBarry Smith -  threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
10883542efc5SMark F. Adams 
10893542efc5SMark F. Adams    Options Database Key:
10901cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
10913542efc5SMark F. Adams 
109295452b02SPatrick Sanan    Notes:
1093af3c827dSMark Adams     Increasing the threshold decreases the rate of coarsening. Conversely reducing the threshold increases the rate of coarsening (aggressive coarsening) and thereby reduces the complexity of the coarse grids, and generally results in slower solver converge rates. Reducing coarse grid complexity reduced the complexity of Galerkin coarse grid construction considerably.
1094af3c827dSMark Adams     Before coarsening or aggregating the graph, GAMG removes small values from the graph with this threshold, and thus reducing the coupling in the graph and a different (perhaps better) coarser set of points.
1095cab9ed1eSBarry Smith 
10963542efc5SMark F. Adams    Level: intermediate
10973542efc5SMark F. Adams 
1098af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
10993542efc5SMark F. Adams @*/
1100c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
11013542efc5SMark F. Adams {
11023542efc5SMark F. Adams   PetscErrorCode ierr;
11033542efc5SMark F. Adams 
11043542efc5SMark F. Adams   PetscFunctionBegin;
11053542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1106c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr);
11073542efc5SMark F. Adams   PetscFunctionReturn(0);
11083542efc5SMark F. Adams }
11093542efc5SMark F. Adams 
1110c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
11113542efc5SMark F. Adams {
1112c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1113c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1114c1eae691SMark Adams   PetscInt i;
1115c1eae691SMark Adams   PetscFunctionBegin;
1116c1eae691SMark Adams   for (i=0;i<n;i++) pc_gamg->threshold[i] = v[i];
1117c1eae691SMark Adams   do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1118c1eae691SMark Adams   PetscFunctionReturn(0);
1119c1eae691SMark Adams }
1120c1eae691SMark Adams 
1121c1eae691SMark Adams /*@
1122c1eae691SMark Adams    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1123c1eae691SMark Adams 
1124c1eae691SMark Adams    Not collective on PC
1125c1eae691SMark Adams 
1126c1eae691SMark Adams    Input Parameters:
1127c1eae691SMark Adams +  pc - the preconditioner context
1128c1eae691SMark Adams -  scale - the threshold value reduction, ussually < 1.0
1129c1eae691SMark Adams 
1130c1eae691SMark Adams    Options Database Key:
1131c1eae691SMark Adams .  -pc_gamg_threshold_scale <v>
1132c1eae691SMark Adams 
1133c1eae691SMark Adams    Level: advanced
1134c1eae691SMark Adams 
1135c1eae691SMark Adams .seealso: ()
1136c1eae691SMark Adams @*/
1137c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1138c1eae691SMark Adams {
1139c1eae691SMark Adams   PetscErrorCode ierr;
11403542efc5SMark F. Adams 
11413542efc5SMark F. Adams   PetscFunctionBegin;
1142c1eae691SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1143c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr);
1144c1eae691SMark Adams   PetscFunctionReturn(0);
1145c1eae691SMark Adams }
1146c1eae691SMark Adams 
1147c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1148c1eae691SMark Adams {
1149c1eae691SMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1150c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1151c1eae691SMark Adams   PetscFunctionBegin;
1152c1eae691SMark Adams   pc_gamg->threshold_scale = v;
11533542efc5SMark F. Adams   PetscFunctionReturn(0);
11543542efc5SMark F. Adams }
11553542efc5SMark F. Adams 
1156e20c40e8SBarry Smith /*@C
1157c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1158676e1743SMark F. Adams 
1159676e1743SMark F. Adams    Collective on PC
1160676e1743SMark F. Adams 
1161676e1743SMark F. Adams    Input Parameters:
1162c60c7ad4SBarry Smith +  pc - the preconditioner context
1163c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1164676e1743SMark F. Adams 
1165676e1743SMark F. Adams    Options Database Key:
1166cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1167676e1743SMark F. Adams 
1168676e1743SMark F. Adams    Level: intermediate
1169676e1743SMark F. Adams 
1170cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1171676e1743SMark F. Adams @*/
117219fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1173676e1743SMark F. Adams {
1174676e1743SMark F. Adams   PetscErrorCode ierr;
1175676e1743SMark F. Adams 
1176676e1743SMark F. Adams   PetscFunctionBegin;
1177676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1178806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1179676e1743SMark F. Adams   PetscFunctionReturn(0);
1180676e1743SMark F. Adams }
1181676e1743SMark F. Adams 
1182e20c40e8SBarry Smith /*@C
1183c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1184c60c7ad4SBarry Smith 
1185c60c7ad4SBarry Smith    Collective on PC
1186c60c7ad4SBarry Smith 
1187c60c7ad4SBarry Smith    Input Parameter:
1188c60c7ad4SBarry Smith .  pc - the preconditioner context
1189c60c7ad4SBarry Smith 
1190c60c7ad4SBarry Smith    Output Parameter:
1191c60c7ad4SBarry Smith .  type - the type of algorithm used
1192c60c7ad4SBarry Smith 
1193c60c7ad4SBarry Smith    Level: intermediate
1194c60c7ad4SBarry Smith 
11951c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1196c60c7ad4SBarry Smith @*/
1197c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1198c60c7ad4SBarry Smith {
1199c60c7ad4SBarry Smith   PetscErrorCode ierr;
1200c60c7ad4SBarry Smith 
1201c60c7ad4SBarry Smith   PetscFunctionBegin;
1202c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1203c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1204c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1205c60c7ad4SBarry Smith }
1206c60c7ad4SBarry Smith 
1207c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1208c60c7ad4SBarry Smith {
1209c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1210c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1211c60c7ad4SBarry Smith 
1212c60c7ad4SBarry Smith   PetscFunctionBegin;
1213c60c7ad4SBarry Smith   *type = pc_gamg->type;
1214c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1215c60c7ad4SBarry Smith }
1216c60c7ad4SBarry Smith 
12171e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1218676e1743SMark F. Adams {
12199d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
12201ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
12211ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1222676e1743SMark F. Adams 
1223676e1743SMark F. Adams   PetscFunctionBegin;
1224c60c7ad4SBarry Smith   pc_gamg->type = type;
12251c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
12269d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
12271ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
12281ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
12291ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1230e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
12313ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
12323ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
12333ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
12343ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
12353ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
12363ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
12373ae0bb68SMark Adams     pc_gamg->data_sz = 0;
12381ab5ffc9SJed Brown   }
12391ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
12401ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
12419d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1242676e1743SMark F. Adams   PetscFunctionReturn(0);
1243676e1743SMark F. Adams }
1244676e1743SMark F. Adams 
12459d766c59SMark Adams /* -------------------------------------------------------------------------- */
12469d766c59SMark Adams /*
12479d766c59SMark Adams    PCMGGetGridComplexity - compute coarse grid complexity of MG hierarchy
12489d766c59SMark Adams 
12499d766c59SMark Adams    Input Parameter:
12509d766c59SMark Adams .  pc - the preconditioner context
12519d766c59SMark Adams 
12529d766c59SMark Adams    Output Parameter:
12539d766c59SMark Adams .  gc - grid complexity = sum_i(nnz_i) / nnz_0
12549d766c59SMark Adams 
12559d766c59SMark Adams    Level: advanced
12569d766c59SMark Adams */
12579d766c59SMark Adams static PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc)
12589d766c59SMark Adams {
12599d766c59SMark Adams   PetscErrorCode ierr;
12609d766c59SMark Adams   PC_MG          *mg      = (PC_MG*)pc->data;
12619d766c59SMark Adams   PC_MG_Levels   **mglevels = mg->levels;
12629d766c59SMark Adams   PetscInt       lev, nnz0 = -1;
12639d766c59SMark Adams   MatInfo        info;
12649d766c59SMark Adams   PetscFunctionBegin;
12659d766c59SMark Adams   if (!mg->nlevels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels");
12669d766c59SMark Adams   for (lev=0, *gc=0; lev<mg->nlevels; lev++) {
126762a6e064SMark Adams     Mat dB;
126862a6e064SMark Adams     ierr = KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB);CHKERRQ(ierr);
126962a6e064SMark Adams     ierr = MatGetInfo(dB,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
12709d766c59SMark Adams     *gc += (PetscReal)info.nz_used;
12719d766c59SMark Adams     if (lev==mg->nlevels-1) nnz0 = info.nz_used;
12729d766c59SMark Adams   }
12739d766c59SMark Adams   if (nnz0) *gc /= (PetscReal)nnz0;
12749d766c59SMark Adams   else *gc = 0;
12759d766c59SMark Adams   PetscFunctionReturn(0);
12769d766c59SMark Adams }
12779d766c59SMark Adams 
12785adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
12795adeb434SBarry Smith {
1280c1eae691SMark Adams   PetscErrorCode ierr,i;
12815adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
12825adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
128323b2d91dSMark Adams   PetscReal       gc=0;
12845adeb434SBarry Smith   PetscFunctionBegin;
12855adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1286459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");CHKERRQ(ierr);
1287c1eae691SMark Adams   for (i=0;i<pc_gamg->current_level;i++) {
1288c1eae691SMark Adams     ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr);
1289c1eae691SMark Adams   }
1290459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1291459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);CHKERRQ(ierr);
1292cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
1293cab9ed1eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr);
1294cab9ed1eSBarry Smith   }
1295171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
1296171cca9aSMark Adams     ierr = PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr);
1297171cca9aSMark Adams   }
1298*ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1299*ce7c7f2fSMark Adams   if (pc_gamg->cpu_pin_coarse_grids) {
1300*ce7c7f2fSMark Adams     /* ierr = PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n");CHKERRQ(ierr); */
1301*ce7c7f2fSMark Adams   }
1302*ce7c7f2fSMark Adams #endif
1303*ce7c7f2fSMark Adams   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
1304*ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n");CHKERRQ(ierr); */
1305*ce7c7f2fSMark Adams   /* } else { */
1306*ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n");CHKERRQ(ierr); */
1307*ce7c7f2fSMark Adams   /* } */
13085adeb434SBarry Smith   if (pc_gamg->ops->view) {
13095adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
13105adeb434SBarry Smith   }
13119d766c59SMark Adams   ierr = PCMGGetGridComplexity(pc,&gc);CHKERRQ(ierr);
13129d766c59SMark Adams   ierr = PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g\n",gc);CHKERRQ(ierr);
13135adeb434SBarry Smith   PetscFunctionReturn(0);
13145adeb434SBarry Smith }
13155adeb434SBarry Smith 
13164416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
13175b89ad90SMark F. Adams {
1318676e1743SMark F. Adams   PetscErrorCode ierr;
1319676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1320676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1321676e1743SMark F. Adams   PetscBool      flag;
13223b4367a7SBarry Smith   MPI_Comm       comm;
132314a9496bSBarry Smith   char           prefix[256];
1324c1eae691SMark Adams   PetscInt       i,n;
132514a9496bSBarry Smith   const char     *pcpre;
13265b89ad90SMark F. Adams   PetscFunctionBegin;
13273b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1328e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1329676e1743SMark F. Adams   {
1330bd94a7aaSJed Brown     char tname[256];
13311a1c1e04SBarry Smith     ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1332bd94a7aaSJed Brown     if (flag) {
1333bd94a7aaSJed Brown       ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
13341ab5ffc9SJed Brown     }
1335cab9ed1eSBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
13361cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1337a303c832SJed Brown     ierr = PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation aggregates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL);CHKERRQ(ierr);
1338cf8ae1d3SMark 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);
1339*ce7c7f2fSMark Adams     ierr = PetscOptionsBool("-pc_gamg_cpu_pin_coarse_grids","Pin coarse grids to the CPU","PCGAMGSetCpuPinCoarseGrids",pc_gamg->cpu_pin_coarse_grids,&pc_gamg->cpu_pin_coarse_grids,NULL);CHKERRQ(ierr);
1340*ce7c7f2fSMark Adams     i = (PetscInt)PCGAMG_LAYOUT_COMPACT;
1341*ce7c7f2fSMark Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_grid_layout_type","0: place reduced grids on processes in natural order; 1: distribute to whole machine for more memory bandwidth","PCGAMGSetCoarseGridLayoutType",i,&i,NULL);CHKERRQ(ierr);
1342*ce7c7f2fSMark Adams     pc_gamg->layout_type = (PCGAMGLayoutType)i;
134394ae4db5SBarry 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);
134494ae4db5SBarry 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);
1345a303c832SJed Brown     ierr = PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL);CHKERRQ(ierr);
1346c1eae691SMark Adams     n = PETSC_GAMG_MAXLEVELS;
1347c1eae691SMark Adams     ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr);
1348efd3c5ceSMark Adams     if (!flag || n < PETSC_GAMG_MAXLEVELS) {
1349efd3c5ceSMark Adams       if (!flag) n = 1;
1350c1eae691SMark Adams       i = n;
1351c1eae691SMark Adams       do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1352c1eae691SMark Adams     }
135394ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1354b7cbab4eSMark Adams 
1355b7cbab4eSMark Adams     /* set options for subtype */
1356e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1357676e1743SMark F. Adams   }
135814a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
135914a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
1360676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
13615b89ad90SMark F. Adams   PetscFunctionReturn(0);
13625b89ad90SMark F. Adams }
13635b89ad90SMark F. Adams 
13645b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
13655b89ad90SMark F. Adams /*MC
13661cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
13675b89ad90SMark F. Adams 
1368280d9858SJed Brown    Options Database Keys:
1369cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1370cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1371cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1372cab9ed1eSBarry 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
1373cab9ed1eSBarry 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>
1374cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1375cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
13766008e27bSRichard Tran Mills .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1377c1eae691SMark Adams -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1378cab9ed1eSBarry Smith 
1379cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1380cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1381cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1382cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1383cab9ed1eSBarry Smith 
1384db9745e2SBarry Smith    Multigrid options:
1385db9745e2SBarry Smith +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1386db9745e2SBarry Smith .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1387db9745e2SBarry Smith .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1388cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
13895b89ad90SMark F. Adams 
13901cc46a46SBarry Smith 
139195452b02SPatrick Sanan   Notes:
139295452b02SPatrick Sanan     In order to obtain good performance for PCGAMG for vector valued problems you must
1393db9745e2SBarry Smith        Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1394db9745e2SBarry Smith        Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1395db9745e2SBarry Smith        See the Users Manual Chapter 4 for more details
13961cc46a46SBarry Smith 
13975b89ad90SMark F. Adams   Level: intermediate
1398280d9858SJed Brown 
13991cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1400171cca9aSMark Adams            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation()
14015b89ad90SMark F. Adams M*/
1402b2573a8aSBarry Smith 
14038cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
14045b89ad90SMark F. Adams {
1405c1eae691SMark Adams   PetscErrorCode ierr,i;
14065b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
14075b89ad90SMark F. Adams   PC_MG          *mg;
14085b89ad90SMark F. Adams 
14095b89ad90SMark F. Adams   PetscFunctionBegin;
14101c1aac46SBarry Smith    /* register AMG type */
14111c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
14121c1aac46SBarry Smith 
14135b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14141c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
14155b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14165b89ad90SMark F. Adams 
14175b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1418b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
141969aca0b8SBarry Smith   ierr         = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
14205b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
14215b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
14225b89ad90SMark F. Adams 
1423b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
14241ab5ffc9SJed Brown 
14259d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14269d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14279d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14289d5b6da9SMark F. Adams   pc_gamg->data    = 0;
14295b89ad90SMark F. Adams 
14309d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
14315b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
14325b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
14335b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
14345b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
14355adeb434SBarry Smith   mg->view                = PCView_GAMG;
14365b89ad90SMark F. Adams 
143797d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
143897d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
1439bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1440bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1441cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr);
14421cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1443cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr);
1444171cca9aSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr);
1445*ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);CHKERRQ(ierr);
1446*ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);CHKERRQ(ierr);
1447bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1448c1eae691SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr);
1449bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1450c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1451bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
14529d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1453d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
14540c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1455171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1456*ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = PETSC_TRUE;
1457*ce7c7f2fSMark Adams   pc_gamg->layout_type      = PCGAMG_LAYOUT_COMPACT;
1458038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
145925a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1460c1eae691SMark Adams   for (i=0;i<PETSC_GAMG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1461c1eae691SMark Adams   pc_gamg->threshold_scale = 1.;
1462c1eae691SMark Adams   pc_gamg->Nlevels          = PETSC_GAMG_MAXLEVELS;
14639ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
1464c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
14659d5b6da9SMark F. Adams 
1466bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1467bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
14685b89ad90SMark F. Adams   PetscFunctionReturn(0);
14695b89ad90SMark F. Adams }
14703e3471ccSMark Adams 
14713e3471ccSMark Adams /*@C
14723e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
14738a690491SBarry Smith     from PCInitializePackage().
14743e3471ccSMark Adams 
14753e3471ccSMark Adams  Level: developer
14763e3471ccSMark Adams 
14773e3471ccSMark Adams  .seealso: PetscInitialize()
14783e3471ccSMark Adams @*/
14793e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
14803e3471ccSMark Adams {
14813e3471ccSMark Adams   PetscErrorCode ierr;
14823e3471ccSMark Adams 
14833e3471ccSMark Adams   PetscFunctionBegin;
14843e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
14853e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
14863e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
14873e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
14888e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
14893e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1490c1c463dbSMark Adams 
1491c1c463dbSMark Adams   /* general events */
1492fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1493fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1494fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1495fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1496c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1497c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1498fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1499c1c463dbSMark Adams 
15005b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG
15015b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
15025b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
15035b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
15045b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
15055b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
15065b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
15075b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
15085b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1509bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
15105b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
15115b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
15125b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
15135b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
15145b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
15155b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
15165b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
15175b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
15185b89ad90SMark F. Adams 
15195b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
15205b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
15215b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
15225b89ad90SMark F. Adams   /* create timer stages */
15235b89ad90SMark F. Adams #if defined GAMG_STAGES
15245b89ad90SMark F. Adams   {
15255b89ad90SMark F. Adams     char     str[32];
15265b89ad90SMark F. Adams     PetscInt lidx;
15275b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
15285b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
15295b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
15305b89ad90SMark F. Adams       sprintf(str,"MG Level %d",lidx);
15315b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
15325b89ad90SMark F. Adams     }
15335b89ad90SMark F. Adams   }
15345b89ad90SMark F. Adams #endif
15355b89ad90SMark F. Adams #endif
15363e3471ccSMark Adams   PetscFunctionReturn(0);
15373e3471ccSMark Adams }
15383e3471ccSMark Adams 
15393e3471ccSMark Adams /*@C
15401c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
15411c1aac46SBarry Smith     called from PetscFinalize() automatically.
15423e3471ccSMark Adams 
15433e3471ccSMark Adams  Level: developer
15443e3471ccSMark Adams 
15453e3471ccSMark Adams  .seealso: PetscFinalize()
15463e3471ccSMark Adams @*/
15473e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
15483e3471ccSMark Adams {
15493e3471ccSMark Adams   PetscErrorCode ierr;
15503e3471ccSMark Adams 
15513e3471ccSMark Adams   PetscFunctionBegin;
15523e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
15533e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
15543e3471ccSMark Adams   PetscFunctionReturn(0);
15553e3471ccSMark Adams }
1556a36cf38bSToby Isaac 
1557a36cf38bSToby Isaac /*@C
1558a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1559a36cf38bSToby Isaac 
1560a36cf38bSToby Isaac  Input Parameters:
1561a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1562a36cf38bSToby Isaac  - create - function for creating the gamg context.
1563a36cf38bSToby Isaac 
1564a36cf38bSToby Isaac   Level: advanced
1565a36cf38bSToby Isaac 
15661c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1567a36cf38bSToby Isaac @*/
1568a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1569a36cf38bSToby Isaac {
1570a36cf38bSToby Isaac   PetscErrorCode ierr;
1571a36cf38bSToby Isaac 
1572a36cf38bSToby Isaac   PetscFunctionBegin;
1573a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1574a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1575a36cf38bSToby Isaac   PetscFunctionReturn(0);
1576a36cf38bSToby Isaac }
1577a36cf38bSToby Isaac 
1578