xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 038f3aa429affcd7cc7e66f572dbd0cce54afe4d)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4b45d2f2cSJed Brown #include "petsc-private/matimpl.h"
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
6b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
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
130cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_AGG;
140cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_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;
190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGOptprol_AGG;
20a2f3521dSMark F. Adams PetscLogEvent PC_GAMGKKTProl_AGG;
210cbbd2e1SMark F. Adams #endif
220cbbd2e1SMark F. Adams 
23b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
24b4fbaa2aSMark F. Adams 
25b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
260cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
27b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
28b4fbaa2aSMark F. Adams #endif
29f96513f1SMatthew G Knepley 
30140e18c1SBarry Smith static PetscFunctionList GAMGList = 0;
319d5b6da9SMark F. Adams 
32d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
33d3d6bff4SMark F. Adams #undef __FUNCT__
34d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
35d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
36d3d6bff4SMark F. Adams {
37d3d6bff4SMark F. Adams   PetscErrorCode  ierr;
38d3d6bff4SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
39d3d6bff4SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
40d3d6bff4SMark F. Adams 
41d3d6bff4SMark F. Adams   PetscFunctionBegin;
42a2f3521dSMark F. Adams   if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */
43878e152fSMark F. Adams     PetscPrintf(((PetscObject)pc)->comm,"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__);
449d5b6da9SMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
4558471d46SMark F. Adams   }
46a2f3521dSMark F. Adams   pc_gamg->data = PETSC_NULL; pc_gamg->data_sz = 0;
47878e152fSMark F. Adams 
48878e152fSMark F. Adams   if (pc_gamg->orig_data) {
49878e152fSMark F. Adams     ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
50878e152fSMark F. Adams   }
51a2f3521dSMark F. Adams   PetscFunctionReturn(0);
52a2f3521dSMark F. Adams }
53a2f3521dSMark F. Adams 
54a2f3521dSMark F. Adams /* private 2x2 Mat Nest for Stokes */
55a2f3521dSMark F. Adams typedef struct{
56a2f3521dSMark F. Adams   Mat A11,A21,A12,Amat;
57a2f3521dSMark F. Adams   IS  prim_is,constr_is;
58a2f3521dSMark F. Adams }GAMGKKTMat;
59a2f3521dSMark F. Adams 
60a2f3521dSMark F. Adams #undef __FUNCT__
61a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatCreate"
62a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatCreate(Mat A, PetscBool iskkt, GAMGKKTMat *out)
63a2f3521dSMark F. Adams {
64a2f3521dSMark F. Adams   PetscFunctionBegin;
65a2f3521dSMark F. Adams   out->Amat = A;
66a2f3521dSMark F. Adams   if (iskkt) {
67a2f3521dSMark F. Adams     PetscErrorCode   ierr;
68a2f3521dSMark F. Adams     IS       is_constraint, is_prime;
69a2f3521dSMark F. Adams     PetscInt nmin,nmax;
70a2f3521dSMark F. Adams 
71a2f3521dSMark F. Adams     ierr = MatGetOwnershipRange(A, &nmin, &nmax);CHKERRQ(ierr);
72a2f3521dSMark F. Adams     ierr = MatFindZeroDiagonals(A, &is_constraint);CHKERRQ(ierr);
73a2f3521dSMark F. Adams     ierr = ISComplement(is_constraint, nmin, nmax, &is_prime);CHKERRQ(ierr);
74a2f3521dSMark F. Adams     out->prim_is = is_prime;
75a2f3521dSMark F. Adams     out->constr_is = is_constraint;
76a2f3521dSMark F. Adams 
77a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_prime, is_prime,      MAT_INITIAL_MATRIX, &out->A11);CHKERRQ(ierr);
78a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12);CHKERRQ(ierr);
79a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21);CHKERRQ(ierr);
80806fa848SBarry Smith   } else {
81a2f3521dSMark F. Adams     out->A11 = A;
82a2f3521dSMark F. Adams     out->A21 = PETSC_NULL;
83a2f3521dSMark F. Adams     out->A12 = PETSC_NULL;
84a2f3521dSMark F. Adams     out->prim_is = PETSC_NULL;
85a2f3521dSMark F. Adams     out->constr_is = PETSC_NULL;
86a2f3521dSMark F. Adams   }
87a2f3521dSMark F. Adams   PetscFunctionReturn(0);
88a2f3521dSMark F. Adams }
89a2f3521dSMark F. Adams 
90a2f3521dSMark F. Adams #undef __FUNCT__
91a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatDestroy"
92a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatDestroy(GAMGKKTMat *mat)
93a2f3521dSMark F. Adams {
94a2f3521dSMark F. Adams   PetscErrorCode   ierr;
95a2f3521dSMark F. Adams 
96a2f3521dSMark F. Adams   PetscFunctionBegin;
97a2f3521dSMark F. Adams   if (mat->A11 && mat->A11 != mat->Amat) {
98a2f3521dSMark F. Adams     ierr = MatDestroy(&mat->A11);CHKERRQ(ierr);
99a2f3521dSMark F. Adams   }
100a2f3521dSMark F. Adams   ierr = MatDestroy(&mat->A21);CHKERRQ(ierr);
101a2f3521dSMark F. Adams   ierr = MatDestroy(&mat->A12);CHKERRQ(ierr);
102a2f3521dSMark F. Adams 
103a2f3521dSMark F. Adams   ierr = ISDestroy(&mat->prim_is);CHKERRQ(ierr);
104a2f3521dSMark F. Adams   ierr = ISDestroy(&mat->constr_is);CHKERRQ(ierr);
105d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
106d3d6bff4SMark F. Adams }
107d3d6bff4SMark F. Adams 
1085b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1095b89ad90SMark F. Adams /*
110a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
111a147abb0SMark F. Adams      of active processors.
1125b89ad90SMark F. Adams 
1135b89ad90SMark F. Adams    Input Parameter:
114a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
115a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
1169d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
117c5bfad50SMark F. Adams    . cr_bs - coarse block size
1189d5b6da9SMark F. Adams    . isLast -
119a2f3521dSMark F. Adams    . stokes -
1203530afc2SMark F. Adams    In/Output Parameter:
121a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
122afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
12311e60469SMark F. Adams    Output Parameter:
1243530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1255b89ad90SMark F. Adams */
1265cb416c2SMark F. Adams 
1275b89ad90SMark F. Adams #undef __FUNCT__
1288263b398SMark F. Adams #define __FUNCT__ "createLevel"
1290cbbd2e1SMark F. Adams static PetscErrorCode createLevel(const PC pc,
1309d5b6da9SMark F. Adams                                    const Mat Amat_fine,
131c5bfad50SMark F. Adams                                    const PetscInt cr_bs,
1329d5b6da9SMark F. Adams                                    const PetscBool isLast,
133a2f3521dSMark F. Adams                                    const PetscBool stokes,
1343530afc2SMark F. Adams                                    Mat *a_P_inout,
135a2f3521dSMark F. Adams                                    Mat *a_Amat_crs,
1361147fc2aSKarl Rupp                                    PetscMPIInt *a_nactive_proc)
1375b89ad90SMark F. Adams {
138a2f3521dSMark F. Adams   PetscErrorCode   ierr;
1399d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
140486a8d0bSJed Brown   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1419d5b6da9SMark F. Adams   const PetscBool  repart = pc_gamg->repart;
1429d5b6da9SMark F. Adams   const PetscInt   min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
143a2f3521dSMark F. Adams   Mat              Cmat,Pold=*a_P_inout;
1449d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)Amat_fine)->comm;
145c5df96a5SBarry Smith   PetscMPIInt      rank,size,new_size,nactive=*a_nactive_proc;
146c5bfad50SMark F. Adams   PetscInt         ncrs_eq,ncrs_prim,f_bs;
1475b89ad90SMark F. Adams 
1485b89ad90SMark F. Adams   PetscFunctionBegin;
149c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
150c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
151c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
15211e60469SMark F. Adams   /* RAP */
1539d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
154038e3b61SMark F. Adams 
155a2f3521dSMark F. Adams   /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/
156a2f3521dSMark F. Adams   ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
157a2f3521dSMark F. Adams   assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)==0);
158a2f3521dSMark F. Adams   ierr = MatGetLocalSize(Cmat, &ncrs_eq, PETSC_NULL);CHKERRQ(ierr);
159a2f3521dSMark F. Adams 
160c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
161a2f3521dSMark F. Adams   {
162472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
163a2f3521dSMark F. Adams     ierr = MatGetSize(Cmat, &ncrs_eq_glob, PETSC_NULL);CHKERRQ(ierr);
164c5df96a5SBarry Smith     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
165c5df96a5SBarry Smith     if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1;
166c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
167c5df96a5SBarry Smith     if (isLast) new_size = 1;
168a2f3521dSMark F. Adams   }
169f852f58cSMark F. Adams 
170c5df96a5SBarry Smith   if (!repart && new_size==nactive) {
171a2f3521dSMark F. Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
172806fa848SBarry Smith   } else {
173a2f3521dSMark F. 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;
174a2f3521dSMark F. Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
175a2f3521dSMark F. Adams     IS              is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
1765a9b9e01SMark F. Adams     VecScatter      vecscat;
17722063be5SMark F. Adams     PetscScalar    *array;
17822063be5SMark F. Adams     Vec             src_crd, dest_crd;
179e33ef3b1SMark F. Adams 
180c5bfad50SMark F. Adams     nloc_old = ncrs_eq/cr_bs; assert(ncrs_eq%cr_bs==0);
1810cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1820cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
183b4fbaa2aSMark F. Adams #endif
184a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
185c5df96a5SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt), &counts);CHKERRQ(ierr);
186a2f3521dSMark F. Adams     if (repart && !stokes) {
187a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1885a9b9e01SMark F. Adams       Mat             adj;
1895a9b9e01SMark F. Adams 
190a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
191c5df96a5SBarry Smith         if (pc_gamg->verbose==1) PetscPrintf(wcomm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,ncrs_eq);
192a2f3521dSMark F. Adams         else {
193a2f3521dSMark F. Adams           PetscInt n;
194a2f3521dSMark F. Adams           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm);CHKERRQ(ierr);
195c5df96a5SBarry Smith           PetscPrintf(wcomm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n);
196a2f3521dSMark F. Adams         }
197a2f3521dSMark F. Adams       }
1985a9b9e01SMark F. Adams 
199a2f3521dSMark F. Adams       /* get 'adj' */
200c5bfad50SMark F. Adams       if (cr_bs == 1) {
201038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
202806fa848SBarry Smith       } else {
203a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
204eb07cef2SMark F. Adams         Mat tMat;
205a2f3521dSMark F. Adams         PetscInt Istart_crs,Iend_crs,ncols,jj,Ii;
206b4fbaa2aSMark F. Adams         const PetscScalar *vals;
207b4fbaa2aSMark F. Adams         const PetscInt *idx;
208a2f3521dSMark F. Adams         PetscInt *d_nnz, *o_nnz, M, N;
2099057884aSMark F. Adams         static PetscInt llev = 0;
210b4fbaa2aSMark F. Adams 
211a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &d_nnz);CHKERRQ(ierr);
212a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &o_nnz);CHKERRQ(ierr);
213a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
214a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
215c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0 ; Ii < Iend_crs ; Ii += cr_bs, jj++) {
21658471d46SMark F. Adams           ierr = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
217c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
218c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
21958471d46SMark F. Adams           ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
220a2f3521dSMark F. Adams           if (d_nnz[jj] > ncrs_prim) d_nnz[jj] = ncrs_prim;
221c5bfad50SMark F. Adams           if (o_nnz[jj] > (M/cr_bs-ncrs_prim)) o_nnz[jj] = M/cr_bs-ncrs_prim;
22258471d46SMark F. Adams         }
2236876a03eSMark F. Adams 
224a2f3521dSMark F. Adams         ierr = MatCreate(wcomm, &tMat);CHKERRQ(ierr);
225806fa848SBarry Smith         ierr = MatSetSizes(tMat, ncrs_prim, ncrs_prim,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
226a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
227a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
228a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
22958471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
2305f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
231eb07cef2SMark F. Adams 
232a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
233c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
23422063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
235eb07cef2SMark F. Adams           for (jj = 0 ; jj < ncols ; jj++) {
236c5bfad50SMark F. Adams             PetscInt dest_col = idx[jj]/cr_bs;
237eb07cef2SMark F. Adams             PetscScalar v = 1.0;
238eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
239eb07cef2SMark F. Adams           }
24022063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
241eb07cef2SMark F. Adams         }
242eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
243eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
244eb07cef2SMark F. Adams 
245b4fbaa2aSMark F. Adams         if (llev++ == -1) {
246b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2478caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
248b4fbaa2aSMark F. Adams           PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
249b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2503bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
251b4fbaa2aSMark F. Adams         }
252b4fbaa2aSMark F. Adams 
253eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
254eb07cef2SMark F. Adams 
255eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
256a2f3521dSMark F. Adams       } /* create 'adj' */
257f150b916SMark F. Adams 
258a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2595a9b9e01SMark F. Adams         char prefix[256];
2605a9b9e01SMark F. Adams         const char *pcpre;
261b4fbaa2aSMark F. Adams         const PetscInt *is_idx;
262b4fbaa2aSMark F. Adams         MatPartitioning  mpart;
263a4b7d37bSMark F. Adams         IS proc_is;
264a2f3521dSMark F. Adams         PetscInt targetPE;
2652f03bc48SMark F. Adams 
2665a9b9e01SMark F. Adams         ierr = MatPartitioningCreate(wcomm, &mpart);CHKERRQ(ierr);
2675ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2689d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2698caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
27059a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
27111e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
272c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
273a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
27411e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2755a9b9e01SMark F. Adams 
2765ef31b24SMark F. Adams         /* collect IS info */
277a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_eq*sizeof(PetscInt), &newproc_idx);CHKERRQ(ierr);
278a4b7d37bSMark F. Adams         ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
279a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
280c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
281a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
282c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
283a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
284eb07cef2SMark F. Adams           }
2855ef31b24SMark F. Adams         }
286a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
287a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2885ef31b24SMark F. Adams       }
2895ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2905a9b9e01SMark F. Adams 
291806fa848SBarry Smith       ierr = ISCreateGeneral(wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2928263b398SMark F. Adams       if (newproc_idx != 0) {
2938263b398SMark F. Adams         ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2945ef31b24SMark F. Adams       }
295806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
296a2f3521dSMark F. Adams 
297a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
2985a9b9e01SMark F. Adams       /* find factor */
299c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
3005a9b9e01SMark F. Adams       else {
3015a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
3025a9b9e01SMark F. Adams         jj = -1;
303c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
304c5df96a5SBarry Smith           if (size%kk==0) { /* a candidate */
305c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
3065a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
3075a9b9e01SMark F. Adams             if (fact > best_fact) {
3085a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
3095a9b9e01SMark F. Adams             }
3105a9b9e01SMark F. Adams           }
3115a9b9e01SMark F. Adams         }
3125a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
313a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
3145a9b9e01SMark F. Adams       }
315c5df96a5SBarry Smith       new_size = size/rfactor;
3165a9b9e01SMark F. Adams 
317c5df96a5SBarry Smith       if (new_size==nactive) {
318a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
3195a9b9e01SMark F. Adams         ierr = PetscFree(counts);CHKERRQ(ierr);
320a2f3521dSMark F. Adams         if (pc_gamg->verbose>0) {
321c5df96a5SBarry Smith           PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq);
322a2f3521dSMark F. Adams         }
3235a9b9e01SMark F. Adams         PetscFunctionReturn(0);
3245a9b9e01SMark F. Adams       }
3255a9b9e01SMark F. Adams 
326c5df96a5SBarry Smith       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq);
327c5df96a5SBarry Smith       targetPE = rank/rfactor;
328a2f3521dSMark F. Adams       ierr = ISCreateStride(wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
3295a9b9e01SMark F. Adams 
330a2f3521dSMark F. Adams       if (stokes) {
331c5bfad50SMark F. Adams         ierr = ISCreateStride(wcomm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim);CHKERRQ(ierr);
3325a9b9e01SMark F. Adams       }
333a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
334e33ef3b1SMark F. Adams 
33511e60469SMark F. Adams     /*
336a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
33711e60469SMark F. Adams      */
338a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
339a2f3521dSMark F. Adams     if (stokes) {
340a2f3521dSMark F. Adams       ierr = ISPartitioningToNumbering(is_eq_newproc_prim, &is_eq_num_prim);CHKERRQ(ierr);
341806fa848SBarry Smith     } else is_eq_num_prim = is_eq_num;
34211e60469SMark F. Adams     /*
343a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
34411e60469SMark F. Adams      */
345c5df96a5SBarry Smith     ierr = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
346c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
347a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
348a2f3521dSMark F. Adams     if (stokes) {
349c5df96a5SBarry Smith       ierr = ISPartitioningCount(is_eq_newproc_prim, size, counts);CHKERRQ(ierr);
350a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_newproc_prim);CHKERRQ(ierr);
351c5df96a5SBarry Smith       ncrs_prim_new = counts[rank]/cr_bs; /* nodes */
352806fa848SBarry Smith     } else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */
353a2f3521dSMark F. Adams 
354a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
3550cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3560cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
357b4fbaa2aSMark F. Adams #endif
358a2f3521dSMark F. Adams 
359a2f3521dSMark F. Adams     /* move data (for primal equations only) */
36022063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
3613bf036e2SBarry Smith     ierr = VecCreate(wcomm, &dest_crd);CHKERRQ(ierr);
362a2f3521dSMark F. Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE);CHKERRQ(ierr);
363470e26f8SMark F. Adams     ierr = VecSetFromOptions(dest_crd);CHKERRQ(ierr); /* this is needed! */
36411e60469SMark F. Adams     /*
3659d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
366c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
36711e60469SMark F. Adams      */
368a2f3521dSMark F. Adams     ierr = PetscMalloc((ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx);CHKERRQ(ierr);
369a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
370a2f3521dSMark F. Adams     for (ii=0,jj=0; ii<ncrs_prim ; ii++) {
371c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
372a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz ; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
37311e60469SMark F. Adams     }
374a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
375806fa848SBarry Smith     ierr = ISCreateGeneral(wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
37692a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
37711e60469SMark F. Adams     /*
37811e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
37911e60469SMark F. Adams      */
380a2f3521dSMark F. Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd);CHKERRQ(ierr);
3819d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols ; jj++) {
382a2f3521dSMark F. Adams       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
383a2f3521dSMark F. Adams       for (ii=0 ; ii<ncrs_prim ; ii++) {
3849d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows ; kk++) {
385a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
386c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
387676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
388d3d6bff4SMark F. Adams         }
389038e3b61SMark F. Adams       }
390eb07cef2SMark F. Adams     }
391eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
392eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
39311e60469SMark F. Adams     /*
39411e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
39511e60469SMark F. Adams       to the correct processor
39611e60469SMark F. Adams     */
397806fa848SBarry Smith     ierr = VecScatterCreate(src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
39811e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
39911e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40011e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40111e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
40211e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
40311e60469SMark F. Adams     /*
40411e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
40511e60469SMark F. Adams     */
406c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
407a2f3521dSMark F. Adams     ierr = PetscMalloc(node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
408a2f3521dSMark F. Adams     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
409a2f3521dSMark F. Adams     strideNew = ncrs_prim_new*ndata_rows;
41011e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
4119d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols ; jj++) {
412a2f3521dSMark F. Adams       for (ii=0 ; ii<ncrs_prim_new ; ii++) {
4139d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows ; kk++) {
414a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
415c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
416d3d6bff4SMark F. Adams         }
417038e3b61SMark F. Adams       }
418038e3b61SMark F. Adams     }
41911e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
42011e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
421a2f3521dSMark F. Adams 
422a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
4230cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4240cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
425ed3f9983SMark F. Adams #endif
426a2f3521dSMark F. Adams 
42711e60469SMark F. Adams     /*
42811e60469SMark F. Adams       Invert for MatGetSubMatrix
42911e60469SMark F. Adams     */
430a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
431a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
432c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
433a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
434a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
435a2f3521dSMark F. Adams     }
436a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
4370cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4380cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
4390cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
440ed3f9983SMark F. Adams #endif
441a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
442a2f3521dSMark F. Adams     {
443a2f3521dSMark F. Adams       Mat mat;
444806fa848SBarry Smith       ierr = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
445a2f3521dSMark F. Adams       *a_Amat_crs = mat;
446c5bfad50SMark F. Adams 
447c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
448c5bfad50SMark F. Adams         PetscInt cbs, rbs;
449c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
450c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
451c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
452c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",rank,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr);
453c5bfad50SMark F. Adams       }
454a2f3521dSMark F. Adams     }
455038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
456a2f3521dSMark F. Adams 
4570cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4580cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
459ed3f9983SMark F. Adams #endif
46011e60469SMark F. Adams     /* prolongator */
46111e60469SMark F. Adams     {
46211e60469SMark F. Adams       IS findices;
463a2f3521dSMark F. Adams       PetscInt Istart,Iend;
464a2f3521dSMark F. Adams       Mat Pnew;
465a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4660cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4670cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
468ed3f9983SMark F. Adams #endif
46911e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
470c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
471806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
47211e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
473c5bfad50SMark F. Adams 
474c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
475c5bfad50SMark F. Adams         PetscInt cbs, rbs;
476c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
477c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
478c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
479c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
480c5bfad50SMark F. Adams       }
4810cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4820cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
483ed3f9983SMark F. Adams #endif
4843530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
485a2f3521dSMark F. Adams 
486a2f3521dSMark F. Adams       /* output - repartitioned */
487a2f3521dSMark F. Adams       *a_P_inout = Pnew;
488e33ef3b1SMark F. Adams     }
489a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4905b89ad90SMark F. Adams 
491c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
492a2f3521dSMark F. Adams   }
4935a9b9e01SMark F. Adams 
494a2f3521dSMark F. Adams   /* outout matrix data */
495c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
496c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
497c8b0795cSMark F. Adams     if (llev==0) {
498c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
499c8b0795cSMark F. Adams       PetscViewerASCIIOpen(wcomm,fname,&viewer);
500c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
501c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
502c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
503c8b0795cSMark F. Adams     }
504c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
505c8b0795cSMark F. Adams     PetscViewerASCIIOpen(wcomm,fname,&viewer);
506c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
507c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
508c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
509c8b0795cSMark F. Adams   }
5105b89ad90SMark F. Adams   PetscFunctionReturn(0);
5115b89ad90SMark F. Adams }
5125b89ad90SMark F. Adams 
5135b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5145b89ad90SMark F. Adams /*
5155b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5165b89ad90SMark F. Adams                     by setting data structures and options.
5175b89ad90SMark F. Adams 
5185b89ad90SMark F. Adams    Input Parameter:
5195b89ad90SMark F. Adams .  pc - the preconditioner context
5205b89ad90SMark F. Adams 
5215b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
5225b89ad90SMark F. Adams 
5235b89ad90SMark F. Adams    Notes:
5245b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
5255b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
5265b89ad90SMark F. Adams */
5275b89ad90SMark F. Adams #undef __FUNCT__
5285b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
5299d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5305b89ad90SMark F. Adams {
5315b89ad90SMark F. Adams   PetscErrorCode  ierr;
5329d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
5335b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
5342adcac29SMark F. Adams   Mat              Pmat = pc->pmat;
535a2f3521dSMark F. Adams   PetscInt         fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
5369d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
537c5df96a5SBarry Smith   PetscMPIInt      rank,size,nactivepe;
538587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
539c8b0795cSMark F. Adams   PetscReal        emaxs[GAMG_MAXLEVELS];
540e696ed0bSMark F. Adams   IS              *ASMLocalIDsArr[GAMG_MAXLEVELS];
541a2f3521dSMark F. Adams   GAMGKKTMat       kktMatsArr[GAMG_MAXLEVELS];
542a2f3521dSMark F. Adams   PetscLogDouble   nnz0=0.,nnztot=0.;
543737a81a9SMark F. Adams   MatInfo          info;
5445169fedaSMark F. Adams   PetscBool        stokes = PETSC_FALSE, redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);
5455ef31b24SMark F. Adams 
5465b89ad90SMark F. Adams   PetscFunctionBegin;
547c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm,&rank);CHKERRQ(ierr);
548c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm,&size);CHKERRQ(ierr);
549dfd5c07aSMark F. Adams 
550c5df96a5SBarry Smith   if (pc_gamg->verbose>2) PetscPrintf(wcomm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",rank,__FUNCT__,pc_gamg->setup_count,pc->setupcalled);
551dfd5c07aSMark F. Adams 
55284d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
553878e152fSMark F. Adams     if (redo_mesh_setup) {
554878e152fSMark F. Adams       /* reset everything */
555878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
556878e152fSMark F. Adams       pc->setupcalled = 0;
557806fa848SBarry Smith     } else {
55884d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
55903a628feSMark F. Adams       /* just do Galerkin grids */
56058471d46SMark F. Adams       Mat B,dA,dB;
561d5280255SMark F. Adams       assert(pc->setupcalled);
56258471d46SMark F. Adams 
5639d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
56458471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5659d5b6da9SMark F. Adams         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
56658471d46SMark F. Adams         /* (re)set to get dirty flag */
5679d5b6da9SMark F. Adams         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
56858471d46SMark F. Adams 
5699d5b6da9SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>-1; level--) {
57003a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
57184d3f75bSMark F. Adams           if (pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
57203a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
573084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
57403a628feSMark F. Adams             mglevels[level]->A = B;
575806fa848SBarry Smith           } else {
576084a8fe3SJed Brown             ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
57758471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
57803a628feSMark F. Adams           }
57958471d46SMark F. Adams           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
58058471d46SMark F. Adams           dB = B;
58158471d46SMark F. Adams         }
5825f8cf99dSMark F. Adams       }
583d5280255SMark F. Adams 
584d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
585d5280255SMark F. Adams 
586d5280255SMark F. Adams       /* PCSetUp_MG seems to insists on setting this to GMRES */
5871ac9965eSMark F. Adams       ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr);
58858471d46SMark F. Adams       PetscFunctionReturn(0);
589eb07cef2SMark F. Adams     }
590878e152fSMark F. Adams   }
591f6536408SMark F. Adams 
592302f38e8SMark F. Adams   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
593eb07cef2SMark F. Adams 
5949d1b15c3SMark F. Adams   ierr = GAMGKKTMatCreate(Pmat, stokes, &kktMatsArr[0]);CHKERRQ(ierr);
5959d1b15c3SMark F. Adams 
596878e152fSMark F. Adams   if (!pc_gamg->data) {
597878e152fSMark F. Adams     if (pc_gamg->orig_data) {
598878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
599878e152fSMark F. Adams       ierr = MatGetLocalSize(Pmat, &qq, PETSC_NULL);CHKERRQ(ierr);
600878e152fSMark F. Adams       pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
601878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
602878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
603878e152fSMark F. Adams       ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
604878e152fSMark F. Adams       for (qq=0;qq<pc_gamg->data_sz;qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
605806fa848SBarry Smith     } else {
606f23aa3ddSBarry Smith       if (!pc_gamg->createdefaultdata) SETERRQ(wcomm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
607f23aa3ddSBarry Smith       if (stokes) SETERRQ(wcomm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
6089d1b15c3SMark F. Adams       ierr = pc_gamg->createdefaultdata(pc, kktMatsArr[0].A11);CHKERRQ(ierr);
6099d5b6da9SMark F. Adams     }
610878e152fSMark F. Adams   }
611878e152fSMark F. Adams 
612878e152fSMark F. Adams   /* cache original data for reuse */
613878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
614878e152fSMark F. Adams     ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr);
615878e152fSMark F. Adams     for (qq=0;qq<pc_gamg->data_sz;qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
616878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
617878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
618878e152fSMark F. Adams   }
619038e3b61SMark F. Adams 
620302f38e8SMark F. Adams   /* get basic dims */
621302f38e8SMark F. Adams   if (stokes) {
622a2f3521dSMark F. Adams     bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
623806fa848SBarry Smith   } else {
624302f38e8SMark F. Adams     ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
625302f38e8SMark F. Adams   }
626a2f3521dSMark F. Adams 
627a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
628c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
62984f9421dSMark F. Adams     PetscInt NN = M;
63084f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
63184f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
6323bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
633806fa848SBarry Smith     } else {
634806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
63584f9421dSMark F. Adams     }
636b2a4f308SMark F. Adams     nnz0 = info.nz_used;
637b2a4f308SMark F. Adams     nnztot = info.nz_used;
638806fa848SBarry Smith     ierr = PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
639c5df96a5SBarry Smith                 rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
640c5df96a5SBarry Smith                 (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr);
641c8b0795cSMark F. Adams   }
64284d3f75bSMark F. Adams 
643a2f3521dSMark F. Adams   /* Get A_i and R_i */
644c5df96a5SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
645c5df96a5SBarry Smith         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (size==1 || nactivepe>1); */
6460205a208SMark F. Adams         level++) {
6475b89ad90SMark F. Adams     level1 = level + 1;
6480cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6490cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
650a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
651a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
652b4fbaa2aSMark F. Adams #endif
653a2f3521dSMark F. Adams #endif
654a2f3521dSMark F. Adams     /* deal with Stokes, get sub matrices */
6559d1b15c3SMark F. Adams     if (level > 0) {
656a2f3521dSMark F. Adams       ierr = GAMGKKTMatCreate(Aarr[level], stokes, &kktMatsArr[level]);CHKERRQ(ierr);
6579d1b15c3SMark F. Adams     }
658c8b0795cSMark F. Adams     { /* construct prolongator */
659725b86d8SJed Brown       Mat Gmat;
6600cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
661a2f3521dSMark F. Adams       Mat Prol11,Prol22;
662c8b0795cSMark F. Adams 
663a2f3521dSMark F. Adams       ierr = pc_gamg->graph(pc,kktMatsArr[level].A11, &Gmat);CHKERRQ(ierr);
6640cbbd2e1SMark F. Adams       ierr = pc_gamg->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
665a2f3521dSMark F. Adams       ierr = pc_gamg->prolongator(pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11);CHKERRQ(ierr);
666c8b0795cSMark F. Adams 
667a2f3521dSMark F. Adams       /* could have failed to create new level */
668a2f3521dSMark F. Adams       if (Prol11) {
6699d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6709d1b15c3SMark F. Adams         ierr = MatGetBlockSizes(Prol11, PETSC_NULL, &bs);CHKERRQ(ierr);
671a2f3521dSMark F. Adams 
672a2f3521dSMark F. Adams         if (stokes) {
673a2f3521dSMark F. Adams           if (!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
674a2f3521dSMark F. Adams           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
675a2f3521dSMark F. Adams           ierr = pc_gamg->formkktprol(pc, Prol11, kktMatsArr[level].A21, &Prol22);CHKERRQ(ierr);
676c8b0795cSMark F. Adams         }
677c8b0795cSMark F. Adams 
678a2f3521dSMark F. Adams         if (pc_gamg->optprol) {
679c8b0795cSMark F. Adams           /* smooth */
680a2f3521dSMark F. Adams           ierr = pc_gamg->optprol(pc, kktMatsArr[level].A11, &Prol11);CHKERRQ(ierr);
681c8b0795cSMark F. Adams         }
682c8b0795cSMark F. Adams 
683a2f3521dSMark F. Adams         if (stokes) {
684da80777bSKarl Rupp           IS is_row[2];
685da80777bSKarl Rupp           Mat a[4];
686da80777bSKarl Rupp 
687da80777bSKarl Rupp           is_row[0] = kktMatsArr[level].prim_is; is_row[1] = kktMatsArr[level].constr_is;
688da80777bSKarl Rupp           a[0] = Prol11; a[1] = PETSC_NULL; a[2] = PETSC_NULL; a[3] = Prol22;
689a2f3521dSMark F. Adams           ierr = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1]);CHKERRQ(ierr);
690806fa848SBarry Smith         } else {
691a2f3521dSMark F. Adams           Parr[level1] = Prol11;
692a2f3521dSMark F. Adams         }
693806fa848SBarry Smith       } else Parr[level1] = PETSC_NULL;
694ffc955d6SMark F. Adams 
695ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
696806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
697ffc955d6SMark F. Adams       }
698ffc955d6SMark F. Adams 
699a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
70041b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
701a2f3521dSMark F. Adams     } /* construct prolongator scope */
7020cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7030cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
704c8b0795cSMark F. Adams #endif
7059d5b6da9SMark F. Adams     /* cache eigen estimate */
7069d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
7079d5b6da9SMark F. Adams       PetscBool flag;
708806fa848SBarry Smith       ierr = PetscObjectComposedDataGetReal((PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
7099d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
710806fa848SBarry Smith     } else emaxs[level] = -1.;
7112adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
712c8b0795cSMark F. Adams     if (!Parr[level1]) {
713806fa848SBarry Smith       if (pc_gamg->verbose) {
714c5df96a5SBarry Smith         ierr =  PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr);
715806fa848SBarry Smith       }
716c8b0795cSMark F. Adams       break;
717c8b0795cSMark F. Adams     }
7180cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7190cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
720b4fbaa2aSMark F. Adams #endif
721a2f3521dSMark F. Adams 
722a2f3521dSMark F. Adams     ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
723806fa848SBarry Smith                         stokes, &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr);
724a2f3521dSMark F. Adams 
7250cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7260cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
727b4fbaa2aSMark F. Adams #endif
728a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
729a2f3521dSMark F. Adams 
730a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0) {
7310cbbd2e1SMark F. Adams       PetscInt NN = M;
7320cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
733a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
7343bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
735806fa848SBarry Smith       } else {
736806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
7370cbbd2e1SMark F. Adams       }
738a2f3521dSMark F. Adams 
7390cbbd2e1SMark F. Adams       nnztot += info.nz_used;
740806fa848SBarry Smith       ierr = PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
741c5df96a5SBarry Smith                   rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
742806fa848SBarry Smith                   (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
743c8b0795cSMark F. Adams     }
744a2f3521dSMark F. Adams 
745a2f3521dSMark F. Adams     /* stop if one node -- could pull back for singular problems */
746c8b0795cSMark F. Adams     if (M/pc_gamg->data_cell_cols < 2) {
74781708dccSMark F. Adams       level++;
74881708dccSMark F. Adams       break;
74981708dccSMark F. Adams     }
7500cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
751b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
752b4fbaa2aSMark F. Adams #endif
753c8b0795cSMark F. Adams   } /* levels */
754c8b0795cSMark F. Adams 
755c8b0795cSMark F. Adams   if (pc_gamg->data) {
756c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
757a2f3521dSMark F. Adams     pc_gamg->data = PETSC_NULL;
7585b89ad90SMark F. Adams   }
759c8b0795cSMark F. Adams 
7600cbbd2e1SMark F. Adams   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
7619d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7625b89ad90SMark F. Adams   fine_level = level;
7639d5b6da9SMark F. Adams   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
7645b89ad90SMark F. Adams 
76584d3f75bSMark F. Adams   /* simple setup */
76684d3f75bSMark F. Adams   if (!PETSC_TRUE) {
76784d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
76884d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
76984d3f75bSMark F. Adams          lidx<fine_level;
77084d3f75bSMark F. Adams          lidx++, level--) {
77184d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
77284d3f75bSMark F. Adams       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77384d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
77484d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
77584d3f75bSMark F. Adams     }
77684d3f75bSMark F. Adams     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77784d3f75bSMark F. Adams 
77884d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
779806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
780d5280255SMark F. Adams     /* set default smoothers & set operators */
7819d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
782587fa25dSMark F. Adams           lidx <= fine_level;
783587fa25dSMark F. Adams           lidx++, level--) {
784ffc955d6SMark F. Adams       KSP smoother;
785ffc955d6SMark F. Adams       PC subpc;
786a2f3521dSMark F. Adams 
7879d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
788f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
789ffc955d6SMark F. Adams 
790a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
791a2f3521dSMark F. Adams       /* set ops */
792a2f3521dSMark F. Adams       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
793a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
794a2f3521dSMark F. Adams 
795a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
796a2f3521dSMark F. Adams       if (stokes) {
797a2f3521dSMark F. Adams         KSP *ksps;
798a2f3521dSMark F. Adams         PetscInt nn;
799a2f3521dSMark F. Adams         ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);CHKERRQ(ierr);
800a2f3521dSMark F. Adams         ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is);CHKERRQ(ierr);
801a2f3521dSMark F. Adams         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
802a2f3521dSMark F. Adams         smoother = ksps[0];
803a2f3521dSMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
804a2f3521dSMark F. Adams         ierr = PetscFree(ksps);CHKERRQ(ierr);
805a2f3521dSMark F. Adams       }
806a2f3521dSMark F. Adams       ierr = GAMGKKTMatDestroy(&kktMatsArr[level]);CHKERRQ(ierr);
807a2f3521dSMark F. Adams 
808a2f3521dSMark F. Adams       /* set defaults */
8096c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
810a2f3521dSMark F. Adams 
811ffc955d6SMark F. Adams       /* override defaults and command line args (!) */
812ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
8132d3561bbSSatish Balay         PetscInt sz;
8142d3561bbSSatish Balay         IS *is;
815a2f3521dSMark F. Adams 
8162d3561bbSSatish Balay         sz = nASMBlocksArr[level];
8172d3561bbSSatish Balay         is = ASMLocalIDsArr[level];
818ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
819ffc955d6SMark F. Adams         if (sz==0) {
820ffc955d6SMark F. Adams           IS is;
821ffc955d6SMark F. Adams           PetscInt my0,kk;
822ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
823ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
82406b43e7eSDmitry Karpeev           ierr = PCGASMSetSubdomains(subpc, 1, &is, PETSC_NULL);CHKERRQ(ierr);
825a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
826806fa848SBarry Smith         } else {
827a94c3b12SMark F. Adams           PetscInt kk;
82806b43e7eSDmitry Karpeev           ierr = PCGASMSetSubdomains(subpc, sz, is, PETSC_NULL);CHKERRQ(ierr);
829a94c3b12SMark F. Adams           for (kk=0;kk<sz;kk++) {
830a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
831a94c3b12SMark F. Adams           }
832ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
833ffc955d6SMark F. Adams         }
834ffc955d6SMark F. Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
835ffc955d6SMark F. Adams 
836ffc955d6SMark F. Adams         ASMLocalIDsArr[level] = PETSC_NULL;
837ffc955d6SMark F. Adams         nASMBlocksArr[level] = 0;
838ffc955d6SMark F. Adams         ierr = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
839806fa848SBarry Smith       } else {
8409d5b6da9SMark F. Adams         ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr);
841ffc955d6SMark F. Adams       }
842d5280255SMark F. Adams     }
843d5280255SMark F. Adams     {
844d5280255SMark F. Adams       /* coarse grid */
845d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
846d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
847d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
848d5280255SMark F. Adams       ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
849d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
850d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
851d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
852d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
853d5280255SMark F. Adams       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
854d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
855d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
8569dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
857d5280255SMark F. Adams     }
858d5280255SMark F. Adams 
859d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
860d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
861d5280255SMark F. Adams     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr);
862d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
863d5280255SMark F. Adams     if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
864d5280255SMark F. Adams 
865d5280255SMark F. Adams     /* create cheby smoothers */
866d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
867d5280255SMark F. Adams           lidx <= fine_level;
868d5280255SMark F. Adams           lidx++, level--) {
869d5280255SMark F. Adams       KSP smoother;
870ffc955d6SMark F. Adams       PetscBool flag;
871d5280255SMark F. Adams       PC subpc;
872d5280255SMark F. Adams 
873ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
874a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
875a2f3521dSMark F. Adams 
876a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
877a2f3521dSMark F. Adams       if (stokes) {
878a2f3521dSMark F. Adams         KSP *ksps;
879a2f3521dSMark F. Adams         PetscInt nn;
880a2f3521dSMark F. Adams         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
881a2f3521dSMark F. Adams         smoother = ksps[0];
882a2f3521dSMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
883a2f3521dSMark F. Adams         ierr = PetscFree(ksps);CHKERRQ(ierr);
884a2f3521dSMark F. Adams       }
885ffc955d6SMark F. Adams 
886ffc955d6SMark F. Adams       /* do my own cheby */
8876c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
888ffc955d6SMark F. Adams       if (flag) {
889ffc955d6SMark F. Adams         PetscReal emax, emin;
890251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
891ffc955d6SMark F. Adams         if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
892587fa25dSMark F. Adams         else { /* eigen estimate 'emax' */
893e696ed0bSMark F. Adams           KSP eksp;
894e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
895b2a4f308SMark F. Adams           Vec bb, xx;
896038e3b61SMark F. Adams 
8975745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr);
8985745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr);
899fc4362bfSMark F. Adams           {
900fc4362bfSMark F. Adams             PetscRandom    rctx;
901fc4362bfSMark F. Adams             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
902fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
903fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
904fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
9055b89ad90SMark F. Adams           }
906a2f3521dSMark F. Adams 
907e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
908e696ed0bSMark F. Adams           {
909e696ed0bSMark F. Adams             PetscInt Istart,Iend,ncols,jj,Ii;
910e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
911e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
912e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0 ; Ii < Iend ; Ii++, jj++) {
913e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
914e696ed0bSMark F. Adams               if (ncols <= 1) {
915e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
916a94c3b12SMark F. Adams               }
917e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
918a94c3b12SMark F. Adams             }
919a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
920a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
921a94c3b12SMark F. Adams           }
922e696ed0bSMark F. Adams 
923fc4362bfSMark F. Adams           ierr = KSPCreate(wcomm, &eksp);CHKERRQ(ierr);
924806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
925fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
9261a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
9271a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
928f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
929f6536408SMark F. Adams 
930f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
931f6536408SMark F. Adams           ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
932fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
9335a9b9e01SMark F. Adams 
934d3d0db20SJed Brown           /* set PC type to be same as smoother */
935ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
936b2a4f308SMark F. Adams 
9375a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
9385a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
9395a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
940fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
9415a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
9425a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
9435a9b9e01SMark F. Adams 
944fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
9455a9b9e01SMark F. Adams 
946fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
947fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
948fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
949f6536408SMark F. Adams 
950ffc955d6SMark F. Adams           if (pc_gamg->verbose > 0) {
951a94c3b12SMark F. Adams             PetscInt N1, tt;
952a94c3b12SMark F. Adams             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
953a94c3b12SMark F. Adams             PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1);
954f6536408SMark F. Adams           }
955fc4362bfSMark F. Adams         }
956038e3b61SMark F. Adams         {
957c5bfad50SMark F. Adams           PetscInt N1, N0;
958c5bfad50SMark F. Adams           ierr = MatGetSize(Aarr[level], &N1, PETSC_NULL);CHKERRQ(ierr);
959c5bfad50SMark F. Adams           ierr = MatGetSize(Aarr[level+1], &N0, PETSC_NULL);CHKERRQ(ierr);
960f6536408SMark F. Adams           /* heuristic - is this crap? */
961b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
9625e7c91beSJed Brown           emin = emax * pc_gamg->eigtarget[0];
9635e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
964038e3b61SMark F. Adams         }
9656c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
966ffc955d6SMark F. Adams       } /* setup checby flag */
967ffc955d6SMark F. Adams     } /* non-coarse levels */
968737a81a9SMark F. Adams 
969d5280255SMark F. Adams     /* clean up */
970d5280255SMark F. Adams     for (level=1;level<pc_gamg->Nlevels;level++) {
971587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
972587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
9735b89ad90SMark F. Adams     }
9740cbbd2e1SMark F. Adams 
9750cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9760cbbd2e1SMark F. Adams 
9771ac9965eSMark F. Adams     if (PETSC_TRUE) {
97858471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
9799d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
98058471d46SMark F. Adams       ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
98158471d46SMark F. Adams     }
982806fa848SBarry Smith   } else {
9835f8cf99dSMark F. Adams     KSP smoother;
984c5df96a5SBarry Smith     if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",rank,__FUNCT__);
9859d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
9865f8cf99dSMark F. Adams     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9875f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
9889d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9895f8cf99dSMark F. Adams   }
9905b89ad90SMark F. Adams   PetscFunctionReturn(0);
9915b89ad90SMark F. Adams }
9925b89ad90SMark F. Adams 
993eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
9945b89ad90SMark F. Adams /*
9955b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
9965b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
9975b89ad90SMark F. Adams 
9985b89ad90SMark F. Adams    Input Parameter:
9995b89ad90SMark F. Adams .  pc - the preconditioner context
10005b89ad90SMark F. Adams 
10015b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
10025b89ad90SMark F. Adams */
10035b89ad90SMark F. Adams #undef __FUNCT__
10045b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
10055b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
10065b89ad90SMark F. Adams {
10075b89ad90SMark F. Adams   PetscErrorCode  ierr;
10085b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
10095b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
10105b89ad90SMark F. Adams 
10115b89ad90SMark F. Adams   PetscFunctionBegin;
10125b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
10135b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
10145b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
10155b89ad90SMark F. Adams   PetscFunctionReturn(0);
10165b89ad90SMark F. Adams }
10175b89ad90SMark F. Adams 
1018676e1743SMark F. Adams 
1019676e1743SMark F. Adams #undef __FUNCT__
1020676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
1021676e1743SMark F. Adams /*@
1022676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
1023676e1743SMark F. Adams    processor reduction.
1024676e1743SMark F. Adams 
1025676e1743SMark F. Adams    Not Collective on PC
1026676e1743SMark F. Adams 
1027676e1743SMark F. Adams    Input Parameters:
1028676e1743SMark F. Adams .  pc - the preconditioner context
1029676e1743SMark F. Adams 
1030676e1743SMark F. Adams 
1031676e1743SMark F. Adams    Options Database Key:
1032676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
1033676e1743SMark F. Adams 
1034676e1743SMark F. Adams    Level: intermediate
1035676e1743SMark F. Adams 
1036676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1037676e1743SMark F. Adams 
1038676e1743SMark F. Adams .seealso: ()
1039676e1743SMark F. Adams @*/
1040676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
1041676e1743SMark F. Adams {
1042676e1743SMark F. Adams   PetscErrorCode ierr;
1043676e1743SMark F. Adams 
1044676e1743SMark F. Adams   PetscFunctionBegin;
1045676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1046676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1047676e1743SMark F. Adams   PetscFunctionReturn(0);
1048676e1743SMark F. Adams }
1049676e1743SMark F. Adams 
1050676e1743SMark F. Adams EXTERN_C_BEGIN
1051676e1743SMark F. Adams #undef __FUNCT__
1052676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
1053676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1054676e1743SMark F. Adams {
1055c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1056c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1057676e1743SMark F. Adams 
1058676e1743SMark F. Adams   PetscFunctionBegin;
10599d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
1060676e1743SMark F. Adams   PetscFunctionReturn(0);
1061676e1743SMark F. Adams }
1062676e1743SMark F. Adams EXTERN_C_END
1063676e1743SMark F. Adams 
1064676e1743SMark F. Adams #undef __FUNCT__
1065389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
1066389730f3SMark F. Adams /*@
1067389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
1068389730f3SMark F. Adams 
1069389730f3SMark F. Adams  Collective on PC
1070389730f3SMark F. Adams 
1071389730f3SMark F. Adams    Input Parameters:
1072389730f3SMark F. Adams .  pc - the preconditioner context
1073389730f3SMark F. Adams 
1074389730f3SMark F. Adams 
1075389730f3SMark F. Adams    Options Database Key:
1076389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
1077389730f3SMark F. Adams 
1078389730f3SMark F. Adams    Level: intermediate
1079389730f3SMark F. Adams 
1080389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1081389730f3SMark F. Adams 
1082389730f3SMark F. Adams .seealso: ()
1083389730f3SMark F. Adams  @*/
1084389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
1085389730f3SMark F. Adams {
1086389730f3SMark F. Adams   PetscErrorCode ierr;
1087389730f3SMark F. Adams 
1088389730f3SMark F. Adams   PetscFunctionBegin;
1089389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1090389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1091389730f3SMark F. Adams   PetscFunctionReturn(0);
1092389730f3SMark F. Adams }
1093389730f3SMark F. Adams 
1094389730f3SMark F. Adams EXTERN_C_BEGIN
1095389730f3SMark F. Adams #undef __FUNCT__
1096389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
1097389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1098389730f3SMark F. Adams {
1099389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1100389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1101389730f3SMark F. Adams 
1102389730f3SMark F. Adams   PetscFunctionBegin;
11039d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
1104389730f3SMark F. Adams   PetscFunctionReturn(0);
1105389730f3SMark F. Adams }
1106389730f3SMark F. Adams EXTERN_C_END
1107389730f3SMark F. Adams 
1108389730f3SMark F. Adams #undef __FUNCT__
11098263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1110676e1743SMark F. Adams /*@
11118263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1112676e1743SMark F. Adams 
1113676e1743SMark F. Adams    Collective on PC
1114676e1743SMark F. Adams 
1115676e1743SMark F. Adams    Input Parameters:
1116676e1743SMark F. Adams .  pc - the preconditioner context
1117676e1743SMark F. Adams 
1118676e1743SMark F. Adams 
1119676e1743SMark F. Adams    Options Database Key:
11208263b398SMark F. Adams .  -pc_gamg_repartition
1121676e1743SMark F. Adams 
1122676e1743SMark F. Adams    Level: intermediate
1123676e1743SMark F. Adams 
1124676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1125676e1743SMark F. Adams 
1126676e1743SMark F. Adams .seealso: ()
1127676e1743SMark F. Adams @*/
11288263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1129676e1743SMark F. Adams {
1130676e1743SMark F. Adams   PetscErrorCode ierr;
1131676e1743SMark F. Adams 
1132676e1743SMark F. Adams   PetscFunctionBegin;
1133676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11348263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1135676e1743SMark F. Adams   PetscFunctionReturn(0);
1136676e1743SMark F. Adams }
1137676e1743SMark F. Adams 
1138676e1743SMark F. Adams EXTERN_C_BEGIN
1139676e1743SMark F. Adams #undef __FUNCT__
11408263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
11418263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1142676e1743SMark F. Adams {
1143c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1144c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1145676e1743SMark F. Adams 
1146676e1743SMark F. Adams   PetscFunctionBegin;
11479d5b6da9SMark F. Adams   pc_gamg->repart = n;
1148676e1743SMark F. Adams   PetscFunctionReturn(0);
1149676e1743SMark F. Adams }
1150676e1743SMark F. Adams EXTERN_C_END
1151676e1743SMark F. Adams 
1152676e1743SMark F. Adams #undef __FUNCT__
1153dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl"
1154dfd5c07aSMark F. Adams /*@
1155dfd5c07aSMark F. Adams    PCGAMGSetReuseProl - Reuse prlongation
1156dfd5c07aSMark F. Adams 
1157dfd5c07aSMark F. Adams    Collective on PC
1158dfd5c07aSMark F. Adams 
1159dfd5c07aSMark F. Adams    Input Parameters:
1160dfd5c07aSMark F. Adams .  pc - the preconditioner context
1161dfd5c07aSMark F. Adams 
1162dfd5c07aSMark F. Adams 
1163dfd5c07aSMark F. Adams    Options Database Key:
1164dfd5c07aSMark F. Adams .  -pc_gamg_reuse_interpolation
1165dfd5c07aSMark F. Adams 
1166dfd5c07aSMark F. Adams    Level: intermediate
1167dfd5c07aSMark F. Adams 
1168dfd5c07aSMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1169dfd5c07aSMark F. Adams 
1170dfd5c07aSMark F. Adams .seealso: ()
1171dfd5c07aSMark F. Adams @*/
1172dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n)
1173dfd5c07aSMark F. Adams {
1174dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1175dfd5c07aSMark F. Adams 
1176dfd5c07aSMark F. Adams   PetscFunctionBegin;
1177dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1178dfd5c07aSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1179dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1180dfd5c07aSMark F. Adams }
1181dfd5c07aSMark F. Adams 
1182dfd5c07aSMark F. Adams EXTERN_C_BEGIN
1183dfd5c07aSMark F. Adams #undef __FUNCT__
1184dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl_GAMG"
1185dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n)
1186dfd5c07aSMark F. Adams {
1187dfd5c07aSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1188dfd5c07aSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1189dfd5c07aSMark F. Adams 
1190dfd5c07aSMark F. Adams   PetscFunctionBegin;
1191dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1192dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1193dfd5c07aSMark F. Adams }
1194dfd5c07aSMark F. Adams EXTERN_C_END
1195dfd5c07aSMark F. Adams 
1196dfd5c07aSMark F. Adams #undef __FUNCT__
1197ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1198ffc955d6SMark F. Adams /*@
1199ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1200ffc955d6SMark F. Adams 
1201ffc955d6SMark F. Adams    Collective on PC
1202ffc955d6SMark F. Adams 
1203ffc955d6SMark F. Adams    Input Parameters:
1204ffc955d6SMark F. Adams .  pc - the preconditioner context
1205ffc955d6SMark F. Adams 
1206ffc955d6SMark F. Adams 
1207ffc955d6SMark F. Adams    Options Database Key:
1208ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1209ffc955d6SMark F. Adams 
1210ffc955d6SMark F. Adams    Level: intermediate
1211ffc955d6SMark F. Adams 
1212ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1213ffc955d6SMark F. Adams 
1214ffc955d6SMark F. Adams .seealso: ()
1215ffc955d6SMark F. Adams @*/
1216ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1217ffc955d6SMark F. Adams {
1218ffc955d6SMark F. Adams   PetscErrorCode ierr;
1219ffc955d6SMark F. Adams 
1220ffc955d6SMark F. Adams   PetscFunctionBegin;
1221ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1222ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1223ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1224ffc955d6SMark F. Adams }
1225ffc955d6SMark F. Adams 
1226ffc955d6SMark F. Adams EXTERN_C_BEGIN
1227ffc955d6SMark F. Adams #undef __FUNCT__
1228ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
1229ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1230ffc955d6SMark F. Adams {
1231ffc955d6SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1232ffc955d6SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1233ffc955d6SMark F. Adams 
1234ffc955d6SMark F. Adams   PetscFunctionBegin;
1235ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1236ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1237ffc955d6SMark F. Adams }
1238ffc955d6SMark F. Adams EXTERN_C_END
1239ffc955d6SMark F. Adams 
1240ffc955d6SMark F. Adams #undef __FUNCT__
12414ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
12424ef23d27SMark F. Adams /*@
12434ef23d27SMark F. Adams    PCGAMGSetNlevels -
12444ef23d27SMark F. Adams 
12454ef23d27SMark F. Adams    Not collective on PC
12464ef23d27SMark F. Adams 
12474ef23d27SMark F. Adams    Input Parameters:
12484ef23d27SMark F. Adams .  pc - the preconditioner context
12494ef23d27SMark F. Adams 
12504ef23d27SMark F. Adams 
12514ef23d27SMark F. Adams    Options Database Key:
12524ef23d27SMark F. Adams .  -pc_mg_levels
12534ef23d27SMark F. Adams 
12544ef23d27SMark F. Adams    Level: intermediate
12554ef23d27SMark F. Adams 
12564ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12574ef23d27SMark F. Adams 
12584ef23d27SMark F. Adams .seealso: ()
12594ef23d27SMark F. Adams @*/
12604ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12614ef23d27SMark F. Adams {
12624ef23d27SMark F. Adams   PetscErrorCode ierr;
12634ef23d27SMark F. Adams 
12644ef23d27SMark F. Adams   PetscFunctionBegin;
12654ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12664ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
12674ef23d27SMark F. Adams   PetscFunctionReturn(0);
12684ef23d27SMark F. Adams }
12694ef23d27SMark F. Adams 
12704ef23d27SMark F. Adams EXTERN_C_BEGIN
12714ef23d27SMark F. Adams #undef __FUNCT__
12724ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
12734ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12744ef23d27SMark F. Adams {
12754ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
12764ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
12774ef23d27SMark F. Adams 
12784ef23d27SMark F. Adams   PetscFunctionBegin;
12799d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12804ef23d27SMark F. Adams   PetscFunctionReturn(0);
12814ef23d27SMark F. Adams }
12824ef23d27SMark F. Adams EXTERN_C_END
12834ef23d27SMark F. Adams 
12844ef23d27SMark F. Adams #undef __FUNCT__
12853542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
12863542efc5SMark F. Adams /*@
12873542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12883542efc5SMark F. Adams 
12893542efc5SMark F. Adams    Not collective on PC
12903542efc5SMark F. Adams 
12913542efc5SMark F. Adams    Input Parameters:
12923542efc5SMark F. Adams .  pc - the preconditioner context
12933542efc5SMark F. Adams 
12943542efc5SMark F. Adams 
12953542efc5SMark F. Adams    Options Database Key:
12963542efc5SMark F. Adams .  -pc_gamg_threshold
12973542efc5SMark F. Adams 
12983542efc5SMark F. Adams    Level: intermediate
12993542efc5SMark F. Adams 
13003542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
13013542efc5SMark F. Adams 
13023542efc5SMark F. Adams .seealso: ()
13033542efc5SMark F. Adams @*/
13043542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
13053542efc5SMark F. Adams {
13063542efc5SMark F. Adams   PetscErrorCode ierr;
13073542efc5SMark F. Adams 
13083542efc5SMark F. Adams   PetscFunctionBegin;
13093542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13103542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
13113542efc5SMark F. Adams   PetscFunctionReturn(0);
13123542efc5SMark F. Adams }
13133542efc5SMark F. Adams 
13143542efc5SMark F. Adams EXTERN_C_BEGIN
13153542efc5SMark F. Adams #undef __FUNCT__
13163542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
13173542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
13183542efc5SMark F. Adams {
1319c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1320c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
13213542efc5SMark F. Adams 
13223542efc5SMark F. Adams   PetscFunctionBegin;
13239d5b6da9SMark F. Adams   pc_gamg->threshold = n;
13243542efc5SMark F. Adams   PetscFunctionReturn(0);
13253542efc5SMark F. Adams }
13263542efc5SMark F. Adams EXTERN_C_END
13273542efc5SMark F. Adams 
13283542efc5SMark F. Adams #undef __FUNCT__
13299d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1330676e1743SMark F. Adams /*@
13319d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1332676e1743SMark F. Adams 
1333676e1743SMark F. Adams    Collective on PC
1334676e1743SMark F. Adams 
1335676e1743SMark F. Adams    Input Parameters:
1336676e1743SMark F. Adams .  pc - the preconditioner context
1337676e1743SMark F. Adams 
1338676e1743SMark F. Adams 
1339676e1743SMark F. Adams    Options Database Key:
13403542efc5SMark F. Adams .  -pc_gamg_type
1341676e1743SMark F. Adams 
1342676e1743SMark F. Adams    Level: intermediate
1343676e1743SMark F. Adams 
1344676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1345676e1743SMark F. Adams 
1346676e1743SMark F. Adams .seealso: ()
1347676e1743SMark F. Adams @*/
134819fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1349676e1743SMark F. Adams {
1350676e1743SMark F. Adams   PetscErrorCode ierr;
1351676e1743SMark F. Adams 
1352676e1743SMark F. Adams   PetscFunctionBegin;
1353676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1354806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1355676e1743SMark F. Adams   PetscFunctionReturn(0);
1356676e1743SMark F. Adams }
1357676e1743SMark F. Adams 
1358676e1743SMark F. Adams EXTERN_C_BEGIN
1359676e1743SMark F. Adams #undef __FUNCT__
13609d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
136119fd82e9SBarry Smith PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1362676e1743SMark F. Adams {
13639d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
1364676e1743SMark F. Adams 
1365676e1743SMark F. Adams   PetscFunctionBegin;
1366140e18c1SBarry Smith   ierr = PetscFunctionListFind(((PetscObject)pc)->comm,GAMGList,type,PETSC_FALSE,(PetscVoidStarFunction)&r);CHKERRQ(ierr);
13679d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
13689d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1369676e1743SMark F. Adams   PetscFunctionReturn(0);
1370676e1743SMark F. Adams }
1371676e1743SMark F. Adams EXTERN_C_END
1372676e1743SMark F. Adams 
13735b89ad90SMark F. Adams #undef __FUNCT__
13745b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
13755b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
13765b89ad90SMark F. Adams {
1377676e1743SMark F. Adams   PetscErrorCode  ierr;
1378676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1379676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1380676e1743SMark F. Adams   PetscBool        flag;
13815e7c91beSJed Brown   PetscInt         two = 2;
1382b43b03e9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
13835b89ad90SMark F. Adams 
13845b89ad90SMark F. Adams   PetscFunctionBegin;
1385676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
1386676e1743SMark F. Adams   {
138775b74bdaSMark F. Adams     /* -pc_gamg_verbose */
13889d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
13899d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
1390806fa848SBarry Smith                            &pc_gamg->verbose, PETSC_NULL);CHKERRQ(ierr);
13918263b398SMark F. Adams     /* -pc_gamg_repartition */
13928263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
13938263b398SMark F. Adams                             "Repartion coarse grids (false)",
13948263b398SMark F. Adams                             "PCGAMGRepartitioning",
13959d5b6da9SMark F. Adams                             pc_gamg->repart,
13969d5b6da9SMark F. Adams                             &pc_gamg->repart,
1397806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1398dfd5c07aSMark F. Adams     /* -pc_gamg_reuse_interpolation */
1399dfd5c07aSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation",
1400dfd5c07aSMark F. Adams                             "Reuse prolongation operator (true)",
1401dfd5c07aSMark F. Adams                             "PCGAMGReuseProl",
1402dfd5c07aSMark F. Adams                             pc_gamg->reuse_prol,
1403dfd5c07aSMark F. Adams                             &pc_gamg->reuse_prol,
1404dfd5c07aSMark F. Adams                             &flag);CHKERRQ(ierr);
1405ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
1406ffc955d6SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1407ffc955d6SMark F. Adams                             "Use aggregation agragates for GASM smoother (false)",
1408ffc955d6SMark F. Adams                             "PCGAMGUseASMAggs",
1409ffc955d6SMark F. Adams                             pc_gamg->use_aggs_in_gasm,
1410ffc955d6SMark F. Adams                             &pc_gamg->use_aggs_in_gasm,
1411806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1412c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1413676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1414676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1415676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
14169d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
14179d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1418806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1419389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1420389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1421389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1422389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
14239d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
14249d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1425806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1426c20e4228SMark F. Adams     /* -pc_gamg_threshold */
14273542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
14283542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
14293542efc5SMark F. Adams                             "PCGAMGSetThreshold",
14309d5b6da9SMark F. Adams                             pc_gamg->threshold,
14319d5b6da9SMark F. Adams                             &pc_gamg->threshold,
1432806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1433806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
1434806fa848SBarry Smith       ierr = PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1435806fa848SBarry Smith     }
14364ef23d27SMark F. Adams 
14375e7c91beSJed Brown     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,PETSC_NULL);CHKERRQ(ierr);
14384ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
14394ef23d27SMark F. Adams                            "Set number of MG levels",
14404ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
14419d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
14429d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
1443806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1444676e1743SMark F. Adams   }
1445676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
14465b89ad90SMark F. Adams   PetscFunctionReturn(0);
14475b89ad90SMark F. Adams }
14485b89ad90SMark F. Adams 
14495b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
14505b89ad90SMark F. Adams /*MC
1451856830bfSJed Brown      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
14529d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
14535b89ad90SMark F. Adams 
1454280d9858SJed Brown    Options Database Keys:
14555b89ad90SMark F. Adams    Multigrid options(inherited)
1456280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1457280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1458280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
14598c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
14605b89ad90SMark F. Adams 
14615b89ad90SMark F. Adams   Level: intermediate
1462280d9858SJed Brown 
14635b89ad90SMark F. Adams   Concepts: multigrid
14645b89ad90SMark F. Adams 
14655b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1466280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
14675b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
14685b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
14695b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
14705b89ad90SMark F. Adams M*/
14715b89ad90SMark F. Adams EXTERN_C_BEGIN
14725b89ad90SMark F. Adams #undef __FUNCT__
14735b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
14745b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
14755b89ad90SMark F. Adams {
14765b89ad90SMark F. Adams   PetscErrorCode  ierr;
14775b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
14785b89ad90SMark F. Adams   PC_MG           *mg;
14790cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14809d5b6da9SMark F. Adams   static long count = 0;
14815ee9c036SSatish Balay #endif
14825b89ad90SMark F. Adams 
14835b89ad90SMark F. Adams   PetscFunctionBegin;
14845b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14855b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
14865b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14875b89ad90SMark F. Adams 
14885b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
14895b89ad90SMark F. Adams   ierr = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr);
14905b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
1491ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
14925b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
14935b89ad90SMark F. Adams 
14949d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14959d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14969d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14979d5b6da9SMark F. Adams   pc_gamg->data = 0;
14985b89ad90SMark F. Adams 
14999d5b6da9SMark F. Adams   /* register AMG type */
15009d5b6da9SMark F. Adams   if (!GAMGList) {
1501140e18c1SBarry Smith     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
1502140e18c1SBarry Smith     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
15039d5b6da9SMark F. Adams   }
15049d5b6da9SMark F. Adams 
15059d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
15065b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
15075b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
15085b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
15095b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
15105b89ad90SMark F. Adams 
15115b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1512676e1743SMark F. Adams                                             "PCGAMGSetProcEqLim_C",
1513676e1743SMark F. Adams                                             "PCGAMGSetProcEqLim_GAMG",
1514806fa848SBarry Smith                                             PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1515676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1516389730f3SMark F. Adams                                             "PCGAMGSetCoarseEqLim_C",
1517389730f3SMark F. Adams                                             "PCGAMGSetCoarseEqLim_GAMG",
1518806fa848SBarry Smith                                             PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1519389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15208263b398SMark F. Adams                                             "PCGAMGSetRepartitioning_C",
15218263b398SMark F. Adams                                             "PCGAMGSetRepartitioning_GAMG",
1522806fa848SBarry Smith                                             PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
1523676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1524dfd5c07aSMark F. Adams                                             "PCGAMGSetReuseProl_C",
1525dfd5c07aSMark F. Adams                                             "PCGAMGSetReuseProl_GAMG",
1526dfd5c07aSMark F. Adams                                             PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr);
1527dfd5c07aSMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1528ffc955d6SMark F. Adams                                             "PCGAMGSetUseASMAggs_C",
1529ffc955d6SMark F. Adams                                             "PCGAMGSetUseASMAggs_GAMG",
1530806fa848SBarry Smith                                             PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1531ffc955d6SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15323542efc5SMark F. Adams                                             "PCGAMGSetThreshold_C",
15333542efc5SMark F. Adams                                             "PCGAMGSetThreshold_GAMG",
1534806fa848SBarry Smith                                             PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
15353542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15369d5b6da9SMark F. Adams                                             "PCGAMGSetType_C",
15379d5b6da9SMark F. Adams                                             "PCGAMGSetType_GAMG",
1538806fa848SBarry Smith                                             PCGAMGSetType_GAMG);CHKERRQ(ierr);
15399d5b6da9SMark F. Adams   pc_gamg->repart = PETSC_FALSE;
1540dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = PETSC_TRUE;
1541ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1542*038f3aa4SMark F. Adams   pc_gamg->min_eq_proc = 50;
1543c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit = 800;
1544a2f3521dSMark F. Adams   pc_gamg->threshold = 0.001;
15459d5b6da9SMark F. Adams   pc_gamg->Nlevels = GAMG_MAXLEVELS;
15469d5b6da9SMark F. Adams   pc_gamg->verbose = 0;
15479d5b6da9SMark F. Adams   pc_gamg->emax_id = -1;
15485e7c91beSJed Brown   pc_gamg->eigtarget[0] = 0.05;
15495e7c91beSJed Brown   pc_gamg->eigtarget[1] = 1.05;
15509d5b6da9SMark F. Adams 
15510cbbd2e1SMark F. Adams   /* private events */
15520cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1553785cba28SMark F. Adams   if (count++ == 0) {
1554806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1555806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
15560cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
15570cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
15580cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1559806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1560806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1561806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1562806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1563806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1564806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1565806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1566806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1567806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1568806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1569806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1570806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1571f852f58cSMark F. Adams 
15720cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
15730cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
15740cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1575b4fbaa2aSMark F. Adams     /* create timer stages */
1576b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1577b4fbaa2aSMark F. Adams     {
1578b4fbaa2aSMark F. Adams       char     str[32];
1579b4fbaa2aSMark F. Adams       PetscInt lidx;
1580806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1581806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1582b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++) {
1583b4fbaa2aSMark F. Adams         sprintf(str,"MG Level %d",lidx);
1584806fa848SBarry Smith         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1585b4fbaa2aSMark F. Adams       }
1586b4fbaa2aSMark F. Adams     }
1587b4fbaa2aSMark F. Adams #endif
1588b4fbaa2aSMark F. Adams   }
1589b4fbaa2aSMark F. Adams #endif
15900cbbd2e1SMark F. Adams   /* general events */
15910cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
1592806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1593806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1594806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1595806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1596806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1597806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1598806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1599806fa848SBarry Smith   ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr);
16000cbbd2e1SMark F. Adams #endif
16010cbbd2e1SMark F. Adams 
160260a6d8f8SMark F. Adams   /* instantiate derived type */
160360a6d8f8SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
160460a6d8f8SMark F. Adams   {
160560a6d8f8SMark F. Adams     char tname[256] = GAMGAGG;
1606806fa848SBarry Smith     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), PETSC_NULL);CHKERRQ(ierr);
160760a6d8f8SMark F. Adams     ierr = PCGAMGSetType(pc, tname);CHKERRQ(ierr);
160860a6d8f8SMark F. Adams   }
160960a6d8f8SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
16105b89ad90SMark F. Adams   PetscFunctionReturn(0);
16115b89ad90SMark F. Adams }
16125b89ad90SMark F. Adams EXTERN_C_END
1613