xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 2fa5cd679192b9b390e47ae2d0650965e6b1d9fa)
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);
74*2fa5cd67SKarl Rupp 
75a2f3521dSMark F. Adams     out->prim_is   = is_prime;
76a2f3521dSMark F. Adams     out->constr_is = is_constraint;
77a2f3521dSMark F. Adams 
78a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_prime, is_prime,      MAT_INITIAL_MATRIX, &out->A11);CHKERRQ(ierr);
79a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12);CHKERRQ(ierr);
80a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21);CHKERRQ(ierr);
81806fa848SBarry Smith   } else {
82a2f3521dSMark F. Adams     out->A11       = A;
83a2f3521dSMark F. Adams     out->A21       = PETSC_NULL;
84a2f3521dSMark F. Adams     out->A12       = PETSC_NULL;
85a2f3521dSMark F. Adams     out->prim_is   = PETSC_NULL;
86a2f3521dSMark F. Adams     out->constr_is = PETSC_NULL;
87a2f3521dSMark F. Adams   }
88a2f3521dSMark F. Adams   PetscFunctionReturn(0);
89a2f3521dSMark F. Adams }
90a2f3521dSMark F. Adams 
91a2f3521dSMark F. Adams #undef __FUNCT__
92a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatDestroy"
93a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatDestroy(GAMGKKTMat *mat)
94a2f3521dSMark F. Adams {
95a2f3521dSMark F. Adams   PetscErrorCode ierr;
96a2f3521dSMark F. Adams 
97a2f3521dSMark F. Adams   PetscFunctionBegin;
98a2f3521dSMark F. Adams   if (mat->A11 && mat->A11 != mat->Amat) {
99a2f3521dSMark F. Adams     ierr = MatDestroy(&mat->A11);CHKERRQ(ierr);
100a2f3521dSMark F. Adams   }
101a2f3521dSMark F. Adams   ierr = MatDestroy(&mat->A21);CHKERRQ(ierr);
102a2f3521dSMark F. Adams   ierr = MatDestroy(&mat->A12);CHKERRQ(ierr);
103a2f3521dSMark F. Adams 
104a2f3521dSMark F. Adams   ierr = ISDestroy(&mat->prim_is);CHKERRQ(ierr);
105a2f3521dSMark F. Adams   ierr = ISDestroy(&mat->constr_is);CHKERRQ(ierr);
106d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
107d3d6bff4SMark F. Adams }
108d3d6bff4SMark F. Adams 
1095b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1105b89ad90SMark F. Adams /*
111a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
112a147abb0SMark F. Adams      of active processors.
1135b89ad90SMark F. Adams 
1145b89ad90SMark F. Adams    Input Parameter:
115a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
116a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
1179d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
118c5bfad50SMark F. Adams    . cr_bs - coarse block size
1199d5b6da9SMark F. Adams    . isLast -
120a2f3521dSMark F. Adams    . stokes -
1213530afc2SMark F. Adams    In/Output Parameter:
122a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
123afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
12411e60469SMark F. Adams    Output Parameter:
1253530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1265b89ad90SMark F. Adams */
1275cb416c2SMark F. Adams 
1285b89ad90SMark F. Adams #undef __FUNCT__
1298263b398SMark F. Adams #define __FUNCT__ "createLevel"
130*2fa5cd67SKarl Rupp static PetscErrorCode createLevel(const PC pc,const Mat Amat_fine,const PetscInt cr_bs,const PetscBool isLast,const PetscBool stokes,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc)
1315b89ad90SMark F. Adams {
132a2f3521dSMark F. Adams   PetscErrorCode  ierr;
1339d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
134486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
1359d5b6da9SMark F. Adams   const PetscBool repart      = pc_gamg->repart;
1369d5b6da9SMark F. Adams   const PetscInt  min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
137a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
1389d5b6da9SMark F. Adams   MPI_Comm        wcomm = ((PetscObject)Amat_fine)->comm;
139c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
140c5bfad50SMark F. Adams   PetscInt        ncrs_eq,ncrs_prim,f_bs;
1415b89ad90SMark F. Adams 
1425b89ad90SMark F. Adams   PetscFunctionBegin;
143c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
144c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
145c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
14611e60469SMark F. Adams   /* RAP */
1479d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
148038e3b61SMark F. Adams 
149a2f3521dSMark F. Adams   /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/
150a2f3521dSMark F. Adams   ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
151a2f3521dSMark F. Adams   assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)==0);
152a2f3521dSMark F. Adams   ierr = MatGetLocalSize(Cmat, &ncrs_eq, PETSC_NULL);CHKERRQ(ierr);
153a2f3521dSMark F. Adams 
154c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
155a2f3521dSMark F. Adams   {
156472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
157a2f3521dSMark F. Adams     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, PETSC_NULL);CHKERRQ(ierr);
158c5df96a5SBarry Smith     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
159c5df96a5SBarry Smith     if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1;
160c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
161c5df96a5SBarry Smith     if (isLast) new_size = 1;
162a2f3521dSMark F. Adams   }
163f852f58cSMark F. Adams 
164*2fa5cd67SKarl Rupp   if (!repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
165*2fa5cd67SKarl Rupp   else {
166a2f3521dSMark 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;
167a2f3521dSMark F. Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
168a2f3521dSMark F. Adams     IS             is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
1695a9b9e01SMark F. Adams     VecScatter     vecscat;
17022063be5SMark F. Adams     PetscScalar    *array;
17122063be5SMark F. Adams     Vec            src_crd, dest_crd;
172e33ef3b1SMark F. Adams 
173c5bfad50SMark F. Adams     nloc_old = ncrs_eq/cr_bs; assert(ncrs_eq%cr_bs==0);
1740cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1750cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
176b4fbaa2aSMark F. Adams #endif
177a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
178c5df96a5SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt), &counts);CHKERRQ(ierr);
179a2f3521dSMark F. Adams     if (repart && !stokes) {
180a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1815a9b9e01SMark F. Adams       Mat adj;
1825a9b9e01SMark F. Adams 
183a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
184c5df96a5SBarry 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);
185a2f3521dSMark F. Adams         else {
186a2f3521dSMark F. Adams           PetscInt n;
187a2f3521dSMark F. Adams           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm);CHKERRQ(ierr);
188c5df96a5SBarry Smith           PetscPrintf(wcomm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n);
189a2f3521dSMark F. Adams         }
190a2f3521dSMark F. Adams       }
1915a9b9e01SMark F. Adams 
192a2f3521dSMark F. Adams       /* get 'adj' */
193c5bfad50SMark F. Adams       if (cr_bs == 1) {
194038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
195806fa848SBarry Smith       } else {
196a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
197eb07cef2SMark F. Adams         Mat               tMat;
198a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
199b4fbaa2aSMark F. Adams         const PetscScalar *vals;
200b4fbaa2aSMark F. Adams         const PetscInt    *idx;
201a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
2029057884aSMark F. Adams         static PetscInt   llev = 0;
203b4fbaa2aSMark F. Adams 
204a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &d_nnz);CHKERRQ(ierr);
205a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &o_nnz);CHKERRQ(ierr);
206a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
207a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
208c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
20958471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
210c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
211c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
21258471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
213a2f3521dSMark F. Adams           if (d_nnz[jj] > ncrs_prim) d_nnz[jj] = ncrs_prim;
214c5bfad50SMark F. Adams           if (o_nnz[jj] > (M/cr_bs-ncrs_prim)) o_nnz[jj] = M/cr_bs-ncrs_prim;
21558471d46SMark F. Adams         }
2166876a03eSMark F. Adams 
217a2f3521dSMark F. Adams         ierr = MatCreate(wcomm, &tMat);CHKERRQ(ierr);
218806fa848SBarry Smith         ierr = MatSetSizes(tMat, ncrs_prim, ncrs_prim,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
219a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
220a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
221a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
22258471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
2235f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
224eb07cef2SMark F. Adams 
225a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
226c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
22722063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
228eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
229c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
230eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
231eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
232eb07cef2SMark F. Adams           }
23322063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
234eb07cef2SMark F. Adams         }
235eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
236eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
237eb07cef2SMark F. Adams 
238b4fbaa2aSMark F. Adams         if (llev++ == -1) {
239b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2408caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
241b4fbaa2aSMark F. Adams           PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
242b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2433bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
244b4fbaa2aSMark F. Adams         }
245b4fbaa2aSMark F. Adams 
246eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
247eb07cef2SMark F. Adams 
248eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
249a2f3521dSMark F. Adams       } /* create 'adj' */
250f150b916SMark F. Adams 
251a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2525a9b9e01SMark F. Adams         char            prefix[256];
2535a9b9e01SMark F. Adams         const char      *pcpre;
254b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
255b4fbaa2aSMark F. Adams         MatPartitioning mpart;
256a4b7d37bSMark F. Adams         IS              proc_is;
257a2f3521dSMark F. Adams         PetscInt        targetPE;
2582f03bc48SMark F. Adams 
2595a9b9e01SMark F. Adams         ierr = MatPartitioningCreate(wcomm, &mpart);CHKERRQ(ierr);
2605ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2619d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2628caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
26359a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
26411e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
265c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
266a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
26711e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2685a9b9e01SMark F. Adams 
2695ef31b24SMark F. Adams         /* collect IS info */
270a2f3521dSMark F. Adams         ierr     = PetscMalloc(ncrs_eq*sizeof(PetscInt), &newproc_idx);CHKERRQ(ierr);
271a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
272a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
273c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
274a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
275c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
276a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
277eb07cef2SMark F. Adams           }
2785ef31b24SMark F. Adams         }
279a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
280a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2815ef31b24SMark F. Adams       }
2825ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2835a9b9e01SMark F. Adams 
284806fa848SBarry Smith       ierr = ISCreateGeneral(wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2858263b398SMark F. Adams       if (newproc_idx != 0) {
2868263b398SMark F. Adams         ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2875ef31b24SMark F. Adams       }
288806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
289a2f3521dSMark F. Adams 
290a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
2915a9b9e01SMark F. Adams       /* find factor */
292c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2935a9b9e01SMark F. Adams       else {
2945a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2955a9b9e01SMark F. Adams         jj = -1;
296c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
297c5df96a5SBarry Smith           if (size%kk==0) { /* a candidate */
298c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2995a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
3005a9b9e01SMark F. Adams             if (fact > best_fact) {
3015a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
3025a9b9e01SMark F. Adams             }
3035a9b9e01SMark F. Adams           }
3045a9b9e01SMark F. Adams         }
3055a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
306a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
3075a9b9e01SMark F. Adams       }
308c5df96a5SBarry Smith       new_size = size/rfactor;
3095a9b9e01SMark F. Adams 
310c5df96a5SBarry Smith       if (new_size==nactive) {
311a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
3125a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
313a2f3521dSMark F. Adams         if (pc_gamg->verbose>0) {
314c5df96a5SBarry Smith           PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq);
315a2f3521dSMark F. Adams         }
3165a9b9e01SMark F. Adams         PetscFunctionReturn(0);
3175a9b9e01SMark F. Adams       }
3185a9b9e01SMark F. Adams 
319c5df96a5SBarry Smith       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq);
320c5df96a5SBarry Smith       targetPE = rank/rfactor;
321a2f3521dSMark F. Adams       ierr     = ISCreateStride(wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
3225a9b9e01SMark F. Adams 
323a2f3521dSMark F. Adams       if (stokes) {
324c5bfad50SMark F. Adams         ierr = ISCreateStride(wcomm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim);CHKERRQ(ierr);
3255a9b9e01SMark F. Adams       }
326a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
327e33ef3b1SMark F. Adams 
32811e60469SMark F. Adams     /*
329a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
33011e60469SMark F. Adams      */
331a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
332a2f3521dSMark F. Adams     if (stokes) {
333a2f3521dSMark F. Adams       ierr = ISPartitioningToNumbering(is_eq_newproc_prim, &is_eq_num_prim);CHKERRQ(ierr);
334806fa848SBarry Smith     } else is_eq_num_prim = is_eq_num;
33511e60469SMark F. Adams     /*
336a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
33711e60469SMark F. Adams      */
338c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
339c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
340a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
341a2f3521dSMark F. Adams     if (stokes) {
342c5df96a5SBarry Smith       ierr          = ISPartitioningCount(is_eq_newproc_prim, size, counts);CHKERRQ(ierr);
343a2f3521dSMark F. Adams       ierr          = ISDestroy(&is_eq_newproc_prim);CHKERRQ(ierr);
344c5df96a5SBarry Smith       ncrs_prim_new = counts[rank]/cr_bs; /* nodes */
345806fa848SBarry Smith     } else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */
346a2f3521dSMark F. Adams 
347a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
3480cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3490cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
350b4fbaa2aSMark F. Adams #endif
351a2f3521dSMark F. Adams 
352a2f3521dSMark F. Adams     /* move data (for primal equations only) */
35322063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
3543bf036e2SBarry Smith     ierr = VecCreate(wcomm, &dest_crd);CHKERRQ(ierr);
355a2f3521dSMark F. Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE);CHKERRQ(ierr);
356470e26f8SMark F. Adams     ierr = VecSetFromOptions(dest_crd);CHKERRQ(ierr); /* this is needed! */
35711e60469SMark F. Adams     /*
3589d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
359c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
36011e60469SMark F. Adams      */
361a2f3521dSMark F. Adams     ierr = PetscMalloc((ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx);CHKERRQ(ierr);
362a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
363a2f3521dSMark F. Adams     for (ii=0,jj=0; ii<ncrs_prim; ii++) {
364c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
365a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
36611e60469SMark F. Adams     }
367a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
368806fa848SBarry Smith     ierr = ISCreateGeneral(wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
36992a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
37011e60469SMark F. Adams     /*
37111e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
37211e60469SMark F. Adams      */
373a2f3521dSMark F. Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd);CHKERRQ(ierr);
3749d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
375a2f3521dSMark F. Adams       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
376a2f3521dSMark F. Adams       for (ii=0; ii<ncrs_prim; ii++) {
3779d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
378a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
379c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
380676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
381d3d6bff4SMark F. Adams         }
382038e3b61SMark F. Adams       }
383eb07cef2SMark F. Adams     }
384eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
385eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
38611e60469SMark F. Adams     /*
38711e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
38811e60469SMark F. Adams       to the correct processor
38911e60469SMark F. Adams     */
390806fa848SBarry Smith     ierr = VecScatterCreate(src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
39111e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
39211e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39311e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39411e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
39511e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
39611e60469SMark F. Adams     /*
39711e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
39811e60469SMark F. Adams     */
399c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
400a2f3521dSMark F. Adams     ierr = PetscMalloc(node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
401*2fa5cd67SKarl Rupp 
402a2f3521dSMark F. Adams     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
403a2f3521dSMark F. Adams     strideNew        = ncrs_prim_new*ndata_rows;
404*2fa5cd67SKarl Rupp 
40511e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
4069d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
407a2f3521dSMark F. Adams       for (ii=0; ii<ncrs_prim_new; ii++) {
4089d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
409a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
410c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
411d3d6bff4SMark F. Adams         }
412038e3b61SMark F. Adams       }
413038e3b61SMark F. Adams     }
41411e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
41511e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
416a2f3521dSMark F. Adams 
417a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
4180cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4190cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
420ed3f9983SMark F. Adams #endif
421a2f3521dSMark F. Adams 
42211e60469SMark F. Adams     /*
42311e60469SMark F. Adams       Invert for MatGetSubMatrix
42411e60469SMark F. Adams     */
425a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
426a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
427c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
428a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
429a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
430a2f3521dSMark F. Adams     }
431a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
4320cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4330cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
4340cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
435ed3f9983SMark F. Adams #endif
436a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
437a2f3521dSMark F. Adams     {
438a2f3521dSMark F. Adams       Mat mat;
439806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
440a2f3521dSMark F. Adams       *a_Amat_crs = mat;
441c5bfad50SMark F. Adams 
442c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
443c5bfad50SMark F. Adams         PetscInt cbs, rbs;
444c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
445c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
446c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
447c5df96a5SBarry 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);
448c5bfad50SMark F. Adams       }
449a2f3521dSMark F. Adams     }
450038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
451a2f3521dSMark F. Adams 
4520cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4530cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
454ed3f9983SMark F. Adams #endif
45511e60469SMark F. Adams     /* prolongator */
45611e60469SMark F. Adams     {
45711e60469SMark F. Adams       IS       findices;
458a2f3521dSMark F. Adams       PetscInt Istart,Iend;
459a2f3521dSMark F. Adams       Mat      Pnew;
460a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4610cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4620cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
463ed3f9983SMark F. Adams #endif
46411e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
465c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
466806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
46711e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
468c5bfad50SMark F. Adams 
469c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
470c5bfad50SMark F. Adams         PetscInt cbs, rbs;
471c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
472c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
473c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
474c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
475c5bfad50SMark F. Adams       }
4760cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4770cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
478ed3f9983SMark F. Adams #endif
4793530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
480a2f3521dSMark F. Adams 
481a2f3521dSMark F. Adams       /* output - repartitioned */
482a2f3521dSMark F. Adams       *a_P_inout = Pnew;
483e33ef3b1SMark F. Adams     }
484a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4855b89ad90SMark F. Adams 
486c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
487a2f3521dSMark F. Adams   }
4885a9b9e01SMark F. Adams 
489a2f3521dSMark F. Adams   /* outout matrix data */
490c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
491c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
492c8b0795cSMark F. Adams     if (llev==0) {
493c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
494c8b0795cSMark F. Adams       PetscViewerASCIIOpen(wcomm,fname,&viewer);
495c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
496c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
497c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
498c8b0795cSMark F. Adams     }
499c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
500c8b0795cSMark F. Adams     PetscViewerASCIIOpen(wcomm,fname,&viewer);
501c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
502c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
503c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
504c8b0795cSMark F. Adams   }
5055b89ad90SMark F. Adams   PetscFunctionReturn(0);
5065b89ad90SMark F. Adams }
5075b89ad90SMark F. Adams 
5085b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5095b89ad90SMark F. Adams /*
5105b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5115b89ad90SMark F. Adams                     by setting data structures and options.
5125b89ad90SMark F. Adams 
5135b89ad90SMark F. Adams    Input Parameter:
5145b89ad90SMark F. Adams .  pc - the preconditioner context
5155b89ad90SMark F. Adams 
5165b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
5175b89ad90SMark F. Adams 
5185b89ad90SMark F. Adams    Notes:
5195b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
5205b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
5215b89ad90SMark F. Adams */
5225b89ad90SMark F. Adams #undef __FUNCT__
5235b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
5249d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5255b89ad90SMark F. Adams {
5265b89ad90SMark F. Adams   PetscErrorCode ierr;
5279d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
5285b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
5292adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
530a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
5319d5b6da9SMark F. Adams   MPI_Comm       wcomm = ((PetscObject)pc)->comm;
532c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
533587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
534c8b0795cSMark F. Adams   PetscReal      emaxs[GAMG_MAXLEVELS];
535e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
536a2f3521dSMark F. Adams   GAMGKKTMat     kktMatsArr[GAMG_MAXLEVELS];
537a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
538737a81a9SMark F. Adams   MatInfo        info;
5395169fedaSMark F. Adams   PetscBool      stokes = PETSC_FALSE, redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);
5405ef31b24SMark F. Adams 
5415b89ad90SMark F. Adams   PetscFunctionBegin;
542c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm,&rank);CHKERRQ(ierr);
543c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm,&size);CHKERRQ(ierr);
544dfd5c07aSMark F. Adams 
545c5df96a5SBarry 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);
546dfd5c07aSMark F. Adams 
54784d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
548878e152fSMark F. Adams     if (redo_mesh_setup) {
549878e152fSMark F. Adams       /* reset everything */
550878e152fSMark F. Adams       ierr            = PCReset_MG(pc);CHKERRQ(ierr);
551878e152fSMark F. Adams       pc->setupcalled = 0;
552806fa848SBarry Smith     } else {
55384d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
55403a628feSMark F. Adams       /* just do Galerkin grids */
55558471d46SMark F. Adams       Mat B,dA,dB;
556d5280255SMark F. Adams       assert(pc->setupcalled);
55758471d46SMark F. Adams 
5589d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
55958471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5609d5b6da9SMark F. Adams         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
56158471d46SMark F. Adams         /* (re)set to get dirty flag */
5629d5b6da9SMark F. Adams         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
56358471d46SMark F. Adams 
5649d5b6da9SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>-1; level--) {
56503a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
56684d3f75bSMark F. Adams           if (pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
56703a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
568084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
569*2fa5cd67SKarl Rupp 
57003a628feSMark F. Adams             mglevels[level]->A = B;
571806fa848SBarry Smith           } else {
572084a8fe3SJed Brown             ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
57358471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
57403a628feSMark F. Adams           }
57558471d46SMark F. Adams           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
57658471d46SMark F. Adams           dB   = B;
57758471d46SMark F. Adams         }
5785f8cf99dSMark F. Adams       }
579d5280255SMark F. Adams 
580d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
581d5280255SMark F. Adams 
582d5280255SMark F. Adams       /* PCSetUp_MG seems to insists on setting this to GMRES */
5831ac9965eSMark F. Adams       ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr);
58458471d46SMark F. Adams       PetscFunctionReturn(0);
585eb07cef2SMark F. Adams     }
586878e152fSMark F. Adams   }
587f6536408SMark F. Adams 
588302f38e8SMark F. Adams   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
589eb07cef2SMark F. Adams 
5909d1b15c3SMark F. Adams   ierr = GAMGKKTMatCreate(Pmat, stokes, &kktMatsArr[0]);CHKERRQ(ierr);
5919d1b15c3SMark F. Adams 
592878e152fSMark F. Adams   if (!pc_gamg->data) {
593878e152fSMark F. Adams     if (pc_gamg->orig_data) {
594878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
595878e152fSMark F. Adams       ierr = MatGetLocalSize(Pmat, &qq, PETSC_NULL);CHKERRQ(ierr);
596*2fa5cd67SKarl Rupp 
597878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
598878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
599878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
600*2fa5cd67SKarl Rupp 
601878e152fSMark F. Adams       ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
602878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
603806fa848SBarry Smith     } else {
604f23aa3ddSBarry Smith       if (!pc_gamg->createdefaultdata) SETERRQ(wcomm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
605f23aa3ddSBarry Smith       if (stokes) SETERRQ(wcomm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
6069d1b15c3SMark F. Adams       ierr = pc_gamg->createdefaultdata(pc, kktMatsArr[0].A11);CHKERRQ(ierr);
6079d5b6da9SMark F. Adams     }
608878e152fSMark F. Adams   }
609878e152fSMark F. Adams 
610878e152fSMark F. Adams   /* cache original data for reuse */
611878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
612878e152fSMark F. Adams     ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr);
613878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
614878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
615878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
616878e152fSMark F. Adams   }
617038e3b61SMark F. Adams 
618302f38e8SMark F. Adams   /* get basic dims */
619*2fa5cd67SKarl Rupp   if (stokes) bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
620*2fa5cd67SKarl Rupp   else {
621302f38e8SMark F. Adams     ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
622302f38e8SMark F. Adams   }
623a2f3521dSMark F. Adams 
624a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
625c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
62684f9421dSMark F. Adams     PetscInt NN = M;
62784f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
62884f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
6293bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
630806fa848SBarry Smith     } else {
631806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
63284f9421dSMark F. Adams     }
633b2a4f308SMark F. Adams     nnz0   = info.nz_used;
634b2a4f308SMark F. Adams     nnztot = info.nz_used;
635806fa848SBarry 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",
636c5df96a5SBarry Smith                          rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
637c5df96a5SBarry Smith                          (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr);
638c8b0795cSMark F. Adams   }
63984d3f75bSMark F. Adams 
640a2f3521dSMark F. Adams   /* Get A_i and R_i */
641c5df96a5SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
642c5df96a5SBarry Smith        level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit);  /* && (size==1 || nactivepe>1); */
6430205a208SMark F. Adams        level++) {
6445b89ad90SMark F. Adams     level1 = level + 1;
6450cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6460cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
647a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
648a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
649b4fbaa2aSMark F. Adams #endif
650a2f3521dSMark F. Adams #endif
651a2f3521dSMark F. Adams     /* deal with Stokes, get sub matrices */
6529d1b15c3SMark F. Adams     if (level > 0) {
653a2f3521dSMark F. Adams       ierr = GAMGKKTMatCreate(Aarr[level], stokes, &kktMatsArr[level]);CHKERRQ(ierr);
6549d1b15c3SMark F. Adams     }
655c8b0795cSMark F. Adams     { /* construct prolongator */
656725b86d8SJed Brown       Mat              Gmat;
6570cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
658a2f3521dSMark F. Adams       Mat              Prol11,Prol22;
659c8b0795cSMark F. Adams 
660a2f3521dSMark F. Adams       ierr = pc_gamg->graph(pc,kktMatsArr[level].A11, &Gmat);CHKERRQ(ierr);
6610cbbd2e1SMark F. Adams       ierr = pc_gamg->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
662a2f3521dSMark F. Adams       ierr = pc_gamg->prolongator(pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11);CHKERRQ(ierr);
663c8b0795cSMark F. Adams 
664a2f3521dSMark F. Adams       /* could have failed to create new level */
665a2f3521dSMark F. Adams       if (Prol11) {
6669d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6679d1b15c3SMark F. Adams         ierr = MatGetBlockSizes(Prol11, PETSC_NULL, &bs);CHKERRQ(ierr);
668a2f3521dSMark F. Adams 
669a2f3521dSMark F. Adams         if (stokes) {
670a2f3521dSMark F. Adams           if (!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
671a2f3521dSMark F. Adams           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
672a2f3521dSMark F. Adams           ierr = pc_gamg->formkktprol(pc, Prol11, kktMatsArr[level].A21, &Prol22);CHKERRQ(ierr);
673c8b0795cSMark F. Adams         }
674c8b0795cSMark F. Adams 
675a2f3521dSMark F. Adams         if (pc_gamg->optprol) {
676c8b0795cSMark F. Adams           /* smooth */
677a2f3521dSMark F. Adams           ierr = pc_gamg->optprol(pc, kktMatsArr[level].A11, &Prol11);CHKERRQ(ierr);
678c8b0795cSMark F. Adams         }
679c8b0795cSMark F. Adams 
680a2f3521dSMark F. Adams         if (stokes) {
681da80777bSKarl Rupp           IS  is_row[2];
682da80777bSKarl Rupp           Mat a[4];
683da80777bSKarl Rupp 
684da80777bSKarl Rupp           is_row[0] = kktMatsArr[level].prim_is; is_row[1] = kktMatsArr[level].constr_is;
685da80777bSKarl Rupp           a[0]      = Prol11; a[1] = PETSC_NULL; a[2] = PETSC_NULL; a[3] = Prol22;
686a2f3521dSMark F. Adams           ierr      = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1]);CHKERRQ(ierr);
687*2fa5cd67SKarl Rupp         } else Parr[level1] = Prol11;
688806fa848SBarry Smith       } else Parr[level1] = PETSC_NULL;
689ffc955d6SMark F. Adams 
690ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
691806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
692ffc955d6SMark F. Adams       }
693ffc955d6SMark F. Adams 
694a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
69541b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
696a2f3521dSMark F. Adams     } /* construct prolongator scope */
6970cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6980cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
699c8b0795cSMark F. Adams #endif
7009d5b6da9SMark F. Adams     /* cache eigen estimate */
7019d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
7029d5b6da9SMark F. Adams       PetscBool flag;
703806fa848SBarry Smith       ierr = PetscObjectComposedDataGetReal((PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
7049d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
705806fa848SBarry Smith     } else emaxs[level] = -1.;
7062adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
707c8b0795cSMark F. Adams     if (!Parr[level1]) {
708806fa848SBarry Smith       if (pc_gamg->verbose) {
709c5df96a5SBarry Smith         ierr =  PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr);
710806fa848SBarry Smith       }
711c8b0795cSMark F. Adams       break;
712c8b0795cSMark F. Adams     }
7130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7140cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
715b4fbaa2aSMark F. Adams #endif
716a2f3521dSMark F. Adams 
717a2f3521dSMark F. Adams     ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
718806fa848SBarry Smith                        stokes, &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr);
719a2f3521dSMark F. Adams 
7200cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7210cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
722b4fbaa2aSMark F. Adams #endif
723a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
724a2f3521dSMark F. Adams 
725a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0) {
7260cbbd2e1SMark F. Adams       PetscInt NN = M;
7270cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
728a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
7293bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
730806fa848SBarry Smith       } else {
731806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
7320cbbd2e1SMark F. Adams       }
733a2f3521dSMark F. Adams 
7340cbbd2e1SMark F. Adams       nnztot += info.nz_used;
735806fa848SBarry Smith       ierr    = PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
736c5df96a5SBarry Smith                             rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
737806fa848SBarry Smith                             (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
738c8b0795cSMark F. Adams     }
739a2f3521dSMark F. Adams 
740a2f3521dSMark F. Adams     /* stop if one node -- could pull back for singular problems */
741c8b0795cSMark F. Adams     if (M/pc_gamg->data_cell_cols < 2) {
74281708dccSMark F. Adams       level++;
74381708dccSMark F. Adams       break;
74481708dccSMark F. Adams     }
7450cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
746b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
747b4fbaa2aSMark F. Adams #endif
748c8b0795cSMark F. Adams   } /* levels */
749c8b0795cSMark F. Adams 
750c8b0795cSMark F. Adams   if (pc_gamg->data) {
751c8b0795cSMark F. Adams     ierr          = PetscFree(pc_gamg->data);CHKERRQ(ierr);
752a2f3521dSMark F. Adams     pc_gamg->data = PETSC_NULL;
7535b89ad90SMark F. Adams   }
754c8b0795cSMark F. Adams 
7550cbbd2e1SMark F. Adams   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
7569d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7575b89ad90SMark F. Adams   fine_level       = level;
7589d5b6da9SMark F. Adams   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
7595b89ad90SMark F. Adams 
76084d3f75bSMark F. Adams   /* simple setup */
76184d3f75bSMark F. Adams   if (!PETSC_TRUE) {
76284d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
76384d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
76484d3f75bSMark F. Adams          lidx<fine_level;
76584d3f75bSMark F. Adams          lidx++, level--) {
76684d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
76784d3f75bSMark F. Adams       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
76884d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
76984d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
77084d3f75bSMark F. Adams     }
77184d3f75bSMark F. Adams     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77284d3f75bSMark F. Adams 
77384d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
774806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
775d5280255SMark F. Adams     /* set default smoothers & set operators */
7769d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
777587fa25dSMark F. Adams          lidx <= fine_level;
778587fa25dSMark F. Adams          lidx++, level--) {
779ffc955d6SMark F. Adams       KSP smoother;
780ffc955d6SMark F. Adams       PC  subpc;
781a2f3521dSMark F. Adams 
7829d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
783f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
784ffc955d6SMark F. Adams 
785a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
786a2f3521dSMark F. Adams       /* set ops */
787a2f3521dSMark F. Adams       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
788a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
789a2f3521dSMark F. Adams 
790a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
791a2f3521dSMark F. Adams       if (stokes) {
792a2f3521dSMark F. Adams         KSP      *ksps;
793a2f3521dSMark F. Adams         PetscInt nn;
794a2f3521dSMark F. Adams         ierr     = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);CHKERRQ(ierr);
795a2f3521dSMark F. Adams         ierr     = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is);CHKERRQ(ierr);
796a2f3521dSMark F. Adams         ierr     = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
797a2f3521dSMark F. Adams         smoother = ksps[0];
798a2f3521dSMark F. Adams         ierr     = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
799a2f3521dSMark F. Adams         ierr     = PetscFree(ksps);CHKERRQ(ierr);
800a2f3521dSMark F. Adams       }
801a2f3521dSMark F. Adams       ierr = GAMGKKTMatDestroy(&kktMatsArr[level]);CHKERRQ(ierr);
802a2f3521dSMark F. Adams 
803a2f3521dSMark F. Adams       /* set defaults */
8046c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
805a2f3521dSMark F. Adams 
806ffc955d6SMark F. Adams       /* override defaults and command line args (!) */
807ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
8082d3561bbSSatish Balay         PetscInt sz;
8092d3561bbSSatish Balay         IS       *is;
810a2f3521dSMark F. Adams 
8112d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
8122d3561bbSSatish Balay         is   = ASMLocalIDsArr[level];
813ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
814ffc955d6SMark F. Adams         if (sz==0) {
815ffc955d6SMark F. Adams           IS       is;
816ffc955d6SMark F. Adams           PetscInt my0,kk;
817ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
818ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
81906b43e7eSDmitry Karpeev           ierr = PCGASMSetSubdomains(subpc, 1, &is, PETSC_NULL);CHKERRQ(ierr);
820a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
821806fa848SBarry Smith         } else {
822a94c3b12SMark F. Adams           PetscInt kk;
82306b43e7eSDmitry Karpeev           ierr = PCGASMSetSubdomains(subpc, sz, is, PETSC_NULL);CHKERRQ(ierr);
824a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
825a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
826a94c3b12SMark F. Adams           }
827ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
828ffc955d6SMark F. Adams         }
829ffc955d6SMark F. Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
830ffc955d6SMark F. Adams 
831ffc955d6SMark F. Adams         ASMLocalIDsArr[level] = PETSC_NULL;
832ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
833ffc955d6SMark F. Adams         ierr                  = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
834806fa848SBarry Smith       } else {
8359d5b6da9SMark F. Adams         ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr);
836ffc955d6SMark F. Adams       }
837d5280255SMark F. Adams     }
838d5280255SMark F. Adams     {
839d5280255SMark F. Adams       /* coarse grid */
840d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
841d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
842d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
843d5280255SMark F. Adams       ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
844d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
845d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
846d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
847d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
848d5280255SMark F. Adams       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
849d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
850d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
8519dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
852d5280255SMark F. Adams     }
853d5280255SMark F. Adams 
854d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
855d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
856d5280255SMark F. Adams     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr);
857d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
858d5280255SMark 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.");
859d5280255SMark F. Adams 
860d5280255SMark F. Adams     /* create cheby smoothers */
861d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
862d5280255SMark F. Adams          lidx <= fine_level;
863d5280255SMark F. Adams          lidx++, level--) {
864d5280255SMark F. Adams       KSP       smoother;
865ffc955d6SMark F. Adams       PetscBool flag;
866d5280255SMark F. Adams       PC        subpc;
867d5280255SMark F. Adams 
868ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
869a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
870a2f3521dSMark F. Adams 
871a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
872a2f3521dSMark F. Adams       if (stokes) {
873a2f3521dSMark F. Adams         KSP      *ksps;
874a2f3521dSMark F. Adams         PetscInt nn;
875a2f3521dSMark F. Adams         ierr     = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
876a2f3521dSMark F. Adams         smoother = ksps[0];
877a2f3521dSMark F. Adams         ierr     = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
878a2f3521dSMark F. Adams         ierr     = PetscFree(ksps);CHKERRQ(ierr);
879a2f3521dSMark F. Adams       }
880ffc955d6SMark F. Adams 
881ffc955d6SMark F. Adams       /* do my own cheby */
8826c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
883ffc955d6SMark F. Adams       if (flag) {
884ffc955d6SMark F. Adams         PetscReal emax, emin;
885251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
886ffc955d6SMark F. Adams         if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
887587fa25dSMark F. Adams         else { /* eigen estimate 'emax' */
888e696ed0bSMark F. Adams           KSP eksp;
889e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
890b2a4f308SMark F. Adams           Vec bb, xx;
891038e3b61SMark F. Adams 
8925745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr);
8935745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr);
894fc4362bfSMark F. Adams           {
895fc4362bfSMark F. Adams             PetscRandom rctx;
896fc4362bfSMark F. Adams             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
897fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
898fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
899fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
9005b89ad90SMark F. Adams           }
901a2f3521dSMark F. Adams 
902e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
903e696ed0bSMark F. Adams           {
904e696ed0bSMark F. Adams             PetscInt    Istart,Iend,ncols,jj,Ii;
905e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
906e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
907e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
908e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
909e696ed0bSMark F. Adams               if (ncols <= 1) {
910e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
911a94c3b12SMark F. Adams               }
912e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
913a94c3b12SMark F. Adams             }
914a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
915a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
916a94c3b12SMark F. Adams           }
917e696ed0bSMark F. Adams 
918fc4362bfSMark F. Adams           ierr = KSPCreate(wcomm, &eksp);CHKERRQ(ierr);
919806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
920fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
9211a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
9221a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
923f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
924f6536408SMark F. Adams 
925f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
926f6536408SMark F. Adams           ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
927fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
9285a9b9e01SMark F. Adams 
929d3d0db20SJed Brown           /* set PC type to be same as smoother */
930ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
931b2a4f308SMark F. Adams 
9325a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
9335a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
9345a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
935fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
9365a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
9375a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
9385a9b9e01SMark F. Adams 
939fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
9405a9b9e01SMark F. Adams 
941fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
942fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
943fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
944f6536408SMark F. Adams 
945ffc955d6SMark F. Adams           if (pc_gamg->verbose > 0) {
946a94c3b12SMark F. Adams             PetscInt N1, tt;
947a94c3b12SMark F. Adams             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
948a94c3b12SMark 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);
949f6536408SMark F. Adams           }
950fc4362bfSMark F. Adams         }
951038e3b61SMark F. Adams         {
952c5bfad50SMark F. Adams           PetscInt N1, N0;
953c5bfad50SMark F. Adams           ierr = MatGetSize(Aarr[level], &N1, PETSC_NULL);CHKERRQ(ierr);
954c5bfad50SMark F. Adams           ierr = MatGetSize(Aarr[level+1], &N0, PETSC_NULL);CHKERRQ(ierr);
955f6536408SMark F. Adams           /* heuristic - is this crap? */
956b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
9575e7c91beSJed Brown           emin  = emax * pc_gamg->eigtarget[0];
9585e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
959038e3b61SMark F. Adams         }
9606c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
961ffc955d6SMark F. Adams       } /* setup checby flag */
962ffc955d6SMark F. Adams     } /* non-coarse levels */
963737a81a9SMark F. Adams 
964d5280255SMark F. Adams     /* clean up */
965d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
966587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
967587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
9685b89ad90SMark F. Adams     }
9690cbbd2e1SMark F. Adams 
9700cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9710cbbd2e1SMark F. Adams 
9721ac9965eSMark F. Adams     if (PETSC_TRUE) {
97358471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
9749d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
97558471d46SMark F. Adams       ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
97658471d46SMark F. Adams     }
977806fa848SBarry Smith   } else {
9785f8cf99dSMark F. Adams     KSP smoother;
979c5df96a5SBarry 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__);
9809d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
9815f8cf99dSMark F. Adams     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9825f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
9839d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9845f8cf99dSMark F. Adams   }
9855b89ad90SMark F. Adams   PetscFunctionReturn(0);
9865b89ad90SMark F. Adams }
9875b89ad90SMark F. Adams 
988eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
9895b89ad90SMark F. Adams /*
9905b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
9915b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
9925b89ad90SMark F. Adams 
9935b89ad90SMark F. Adams    Input Parameter:
9945b89ad90SMark F. Adams .  pc - the preconditioner context
9955b89ad90SMark F. Adams 
9965b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
9975b89ad90SMark F. Adams */
9985b89ad90SMark F. Adams #undef __FUNCT__
9995b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
10005b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
10015b89ad90SMark F. Adams {
10025b89ad90SMark F. Adams   PetscErrorCode ierr;
10035b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
10045b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
10055b89ad90SMark F. Adams 
10065b89ad90SMark F. Adams   PetscFunctionBegin;
10075b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
10085b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
10095b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
10105b89ad90SMark F. Adams   PetscFunctionReturn(0);
10115b89ad90SMark F. Adams }
10125b89ad90SMark F. Adams 
1013676e1743SMark F. Adams 
1014676e1743SMark F. Adams #undef __FUNCT__
1015676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
1016676e1743SMark F. Adams /*@
1017676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
1018676e1743SMark F. Adams    processor reduction.
1019676e1743SMark F. Adams 
1020676e1743SMark F. Adams    Not Collective on PC
1021676e1743SMark F. Adams 
1022676e1743SMark F. Adams    Input Parameters:
1023676e1743SMark F. Adams .  pc - the preconditioner context
1024676e1743SMark F. Adams 
1025676e1743SMark F. Adams 
1026676e1743SMark F. Adams    Options Database Key:
1027676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
1028676e1743SMark F. Adams 
1029676e1743SMark F. Adams    Level: intermediate
1030676e1743SMark F. Adams 
1031676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1032676e1743SMark F. Adams 
1033676e1743SMark F. Adams .seealso: ()
1034676e1743SMark F. Adams @*/
1035676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
1036676e1743SMark F. Adams {
1037676e1743SMark F. Adams   PetscErrorCode ierr;
1038676e1743SMark F. Adams 
1039676e1743SMark F. Adams   PetscFunctionBegin;
1040676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1041676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1042676e1743SMark F. Adams   PetscFunctionReturn(0);
1043676e1743SMark F. Adams }
1044676e1743SMark F. Adams 
1045676e1743SMark F. Adams EXTERN_C_BEGIN
1046676e1743SMark F. Adams #undef __FUNCT__
1047676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
1048676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1049676e1743SMark F. Adams {
1050c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1051c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1052676e1743SMark F. Adams 
1053676e1743SMark F. Adams   PetscFunctionBegin;
10549d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
1055676e1743SMark F. Adams   PetscFunctionReturn(0);
1056676e1743SMark F. Adams }
1057676e1743SMark F. Adams EXTERN_C_END
1058676e1743SMark F. Adams 
1059676e1743SMark F. Adams #undef __FUNCT__
1060389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
1061389730f3SMark F. Adams /*@
1062389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
1063389730f3SMark F. Adams 
1064389730f3SMark F. Adams  Collective on PC
1065389730f3SMark F. Adams 
1066389730f3SMark F. Adams    Input Parameters:
1067389730f3SMark F. Adams .  pc - the preconditioner context
1068389730f3SMark F. Adams 
1069389730f3SMark F. Adams 
1070389730f3SMark F. Adams    Options Database Key:
1071389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
1072389730f3SMark F. Adams 
1073389730f3SMark F. Adams    Level: intermediate
1074389730f3SMark F. Adams 
1075389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1076389730f3SMark F. Adams 
1077389730f3SMark F. Adams .seealso: ()
1078389730f3SMark F. Adams  @*/
1079389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
1080389730f3SMark F. Adams {
1081389730f3SMark F. Adams   PetscErrorCode ierr;
1082389730f3SMark F. Adams 
1083389730f3SMark F. Adams   PetscFunctionBegin;
1084389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1085389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1086389730f3SMark F. Adams   PetscFunctionReturn(0);
1087389730f3SMark F. Adams }
1088389730f3SMark F. Adams 
1089389730f3SMark F. Adams EXTERN_C_BEGIN
1090389730f3SMark F. Adams #undef __FUNCT__
1091389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
1092389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1093389730f3SMark F. Adams {
1094389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1095389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1096389730f3SMark F. Adams 
1097389730f3SMark F. Adams   PetscFunctionBegin;
10989d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
1099389730f3SMark F. Adams   PetscFunctionReturn(0);
1100389730f3SMark F. Adams }
1101389730f3SMark F. Adams EXTERN_C_END
1102389730f3SMark F. Adams 
1103389730f3SMark F. Adams #undef __FUNCT__
11048263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1105676e1743SMark F. Adams /*@
11068263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1107676e1743SMark F. Adams 
1108676e1743SMark F. Adams    Collective on PC
1109676e1743SMark F. Adams 
1110676e1743SMark F. Adams    Input Parameters:
1111676e1743SMark F. Adams .  pc - the preconditioner context
1112676e1743SMark F. Adams 
1113676e1743SMark F. Adams 
1114676e1743SMark F. Adams    Options Database Key:
11158263b398SMark F. Adams .  -pc_gamg_repartition
1116676e1743SMark F. Adams 
1117676e1743SMark F. Adams    Level: intermediate
1118676e1743SMark F. Adams 
1119676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1120676e1743SMark F. Adams 
1121676e1743SMark F. Adams .seealso: ()
1122676e1743SMark F. Adams @*/
11238263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1124676e1743SMark F. Adams {
1125676e1743SMark F. Adams   PetscErrorCode ierr;
1126676e1743SMark F. Adams 
1127676e1743SMark F. Adams   PetscFunctionBegin;
1128676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11298263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1130676e1743SMark F. Adams   PetscFunctionReturn(0);
1131676e1743SMark F. Adams }
1132676e1743SMark F. Adams 
1133676e1743SMark F. Adams EXTERN_C_BEGIN
1134676e1743SMark F. Adams #undef __FUNCT__
11358263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
11368263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1137676e1743SMark F. Adams {
1138c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1139c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1140676e1743SMark F. Adams 
1141676e1743SMark F. Adams   PetscFunctionBegin;
11429d5b6da9SMark F. Adams   pc_gamg->repart = n;
1143676e1743SMark F. Adams   PetscFunctionReturn(0);
1144676e1743SMark F. Adams }
1145676e1743SMark F. Adams EXTERN_C_END
1146676e1743SMark F. Adams 
1147676e1743SMark F. Adams #undef __FUNCT__
1148dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl"
1149dfd5c07aSMark F. Adams /*@
1150dfd5c07aSMark F. Adams    PCGAMGSetReuseProl - Reuse prlongation
1151dfd5c07aSMark F. Adams 
1152dfd5c07aSMark F. Adams    Collective on PC
1153dfd5c07aSMark F. Adams 
1154dfd5c07aSMark F. Adams    Input Parameters:
1155dfd5c07aSMark F. Adams .  pc - the preconditioner context
1156dfd5c07aSMark F. Adams 
1157dfd5c07aSMark F. Adams 
1158dfd5c07aSMark F. Adams    Options Database Key:
1159dfd5c07aSMark F. Adams .  -pc_gamg_reuse_interpolation
1160dfd5c07aSMark F. Adams 
1161dfd5c07aSMark F. Adams    Level: intermediate
1162dfd5c07aSMark F. Adams 
1163dfd5c07aSMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1164dfd5c07aSMark F. Adams 
1165dfd5c07aSMark F. Adams .seealso: ()
1166dfd5c07aSMark F. Adams @*/
1167dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n)
1168dfd5c07aSMark F. Adams {
1169dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1170dfd5c07aSMark F. Adams 
1171dfd5c07aSMark F. Adams   PetscFunctionBegin;
1172dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1173dfd5c07aSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1174dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1175dfd5c07aSMark F. Adams }
1176dfd5c07aSMark F. Adams 
1177dfd5c07aSMark F. Adams EXTERN_C_BEGIN
1178dfd5c07aSMark F. Adams #undef __FUNCT__
1179dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl_GAMG"
1180dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n)
1181dfd5c07aSMark F. Adams {
1182dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1183dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1184dfd5c07aSMark F. Adams 
1185dfd5c07aSMark F. Adams   PetscFunctionBegin;
1186dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1187dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1188dfd5c07aSMark F. Adams }
1189dfd5c07aSMark F. Adams EXTERN_C_END
1190dfd5c07aSMark F. Adams 
1191dfd5c07aSMark F. Adams #undef __FUNCT__
1192ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1193ffc955d6SMark F. Adams /*@
1194ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1195ffc955d6SMark F. Adams 
1196ffc955d6SMark F. Adams    Collective on PC
1197ffc955d6SMark F. Adams 
1198ffc955d6SMark F. Adams    Input Parameters:
1199ffc955d6SMark F. Adams .  pc - the preconditioner context
1200ffc955d6SMark F. Adams 
1201ffc955d6SMark F. Adams 
1202ffc955d6SMark F. Adams    Options Database Key:
1203ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1204ffc955d6SMark F. Adams 
1205ffc955d6SMark F. Adams    Level: intermediate
1206ffc955d6SMark F. Adams 
1207ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1208ffc955d6SMark F. Adams 
1209ffc955d6SMark F. Adams .seealso: ()
1210ffc955d6SMark F. Adams @*/
1211ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1212ffc955d6SMark F. Adams {
1213ffc955d6SMark F. Adams   PetscErrorCode ierr;
1214ffc955d6SMark F. Adams 
1215ffc955d6SMark F. Adams   PetscFunctionBegin;
1216ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1217ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1218ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1219ffc955d6SMark F. Adams }
1220ffc955d6SMark F. Adams 
1221ffc955d6SMark F. Adams EXTERN_C_BEGIN
1222ffc955d6SMark F. Adams #undef __FUNCT__
1223ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
1224ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1225ffc955d6SMark F. Adams {
1226ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1227ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1228ffc955d6SMark F. Adams 
1229ffc955d6SMark F. Adams   PetscFunctionBegin;
1230ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1231ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1232ffc955d6SMark F. Adams }
1233ffc955d6SMark F. Adams EXTERN_C_END
1234ffc955d6SMark F. Adams 
1235ffc955d6SMark F. Adams #undef __FUNCT__
12364ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
12374ef23d27SMark F. Adams /*@
12384ef23d27SMark F. Adams    PCGAMGSetNlevels -
12394ef23d27SMark F. Adams 
12404ef23d27SMark F. Adams    Not collective on PC
12414ef23d27SMark F. Adams 
12424ef23d27SMark F. Adams    Input Parameters:
12434ef23d27SMark F. Adams .  pc - the preconditioner context
12444ef23d27SMark F. Adams 
12454ef23d27SMark F. Adams 
12464ef23d27SMark F. Adams    Options Database Key:
12474ef23d27SMark F. Adams .  -pc_mg_levels
12484ef23d27SMark F. Adams 
12494ef23d27SMark F. Adams    Level: intermediate
12504ef23d27SMark F. Adams 
12514ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12524ef23d27SMark F. Adams 
12534ef23d27SMark F. Adams .seealso: ()
12544ef23d27SMark F. Adams @*/
12554ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12564ef23d27SMark F. Adams {
12574ef23d27SMark F. Adams   PetscErrorCode ierr;
12584ef23d27SMark F. Adams 
12594ef23d27SMark F. Adams   PetscFunctionBegin;
12604ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12614ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
12624ef23d27SMark F. Adams   PetscFunctionReturn(0);
12634ef23d27SMark F. Adams }
12644ef23d27SMark F. Adams 
12654ef23d27SMark F. Adams EXTERN_C_BEGIN
12664ef23d27SMark F. Adams #undef __FUNCT__
12674ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
12684ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12694ef23d27SMark F. Adams {
12704ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
12714ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12724ef23d27SMark F. Adams 
12734ef23d27SMark F. Adams   PetscFunctionBegin;
12749d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12754ef23d27SMark F. Adams   PetscFunctionReturn(0);
12764ef23d27SMark F. Adams }
12774ef23d27SMark F. Adams EXTERN_C_END
12784ef23d27SMark F. Adams 
12794ef23d27SMark F. Adams #undef __FUNCT__
12803542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
12813542efc5SMark F. Adams /*@
12823542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12833542efc5SMark F. Adams 
12843542efc5SMark F. Adams    Not collective on PC
12853542efc5SMark F. Adams 
12863542efc5SMark F. Adams    Input Parameters:
12873542efc5SMark F. Adams .  pc - the preconditioner context
12883542efc5SMark F. Adams 
12893542efc5SMark F. Adams 
12903542efc5SMark F. Adams    Options Database Key:
12913542efc5SMark F. Adams .  -pc_gamg_threshold
12923542efc5SMark F. Adams 
12933542efc5SMark F. Adams    Level: intermediate
12943542efc5SMark F. Adams 
12953542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12963542efc5SMark F. Adams 
12973542efc5SMark F. Adams .seealso: ()
12983542efc5SMark F. Adams @*/
12993542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
13003542efc5SMark F. Adams {
13013542efc5SMark F. Adams   PetscErrorCode ierr;
13023542efc5SMark F. Adams 
13033542efc5SMark F. Adams   PetscFunctionBegin;
13043542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13053542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
13063542efc5SMark F. Adams   PetscFunctionReturn(0);
13073542efc5SMark F. Adams }
13083542efc5SMark F. Adams 
13093542efc5SMark F. Adams EXTERN_C_BEGIN
13103542efc5SMark F. Adams #undef __FUNCT__
13113542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
13123542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
13133542efc5SMark F. Adams {
1314c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1315c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
13163542efc5SMark F. Adams 
13173542efc5SMark F. Adams   PetscFunctionBegin;
13189d5b6da9SMark F. Adams   pc_gamg->threshold = n;
13193542efc5SMark F. Adams   PetscFunctionReturn(0);
13203542efc5SMark F. Adams }
13213542efc5SMark F. Adams EXTERN_C_END
13223542efc5SMark F. Adams 
13233542efc5SMark F. Adams #undef __FUNCT__
13249d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1325676e1743SMark F. Adams /*@
13269d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1327676e1743SMark F. Adams 
1328676e1743SMark F. Adams    Collective on PC
1329676e1743SMark F. Adams 
1330676e1743SMark F. Adams    Input Parameters:
1331676e1743SMark F. Adams .  pc - the preconditioner context
1332676e1743SMark F. Adams 
1333676e1743SMark F. Adams 
1334676e1743SMark F. Adams    Options Database Key:
13353542efc5SMark F. Adams .  -pc_gamg_type
1336676e1743SMark F. Adams 
1337676e1743SMark F. Adams    Level: intermediate
1338676e1743SMark F. Adams 
1339676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1340676e1743SMark F. Adams 
1341676e1743SMark F. Adams .seealso: ()
1342676e1743SMark F. Adams @*/
134319fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1344676e1743SMark F. Adams {
1345676e1743SMark F. Adams   PetscErrorCode ierr;
1346676e1743SMark F. Adams 
1347676e1743SMark F. Adams   PetscFunctionBegin;
1348676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1349806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1350676e1743SMark F. Adams   PetscFunctionReturn(0);
1351676e1743SMark F. Adams }
1352676e1743SMark F. Adams 
1353676e1743SMark F. Adams EXTERN_C_BEGIN
1354676e1743SMark F. Adams #undef __FUNCT__
13559d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
135619fd82e9SBarry Smith PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1357676e1743SMark F. Adams {
13589d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
1359676e1743SMark F. Adams 
1360676e1743SMark F. Adams   PetscFunctionBegin;
1361140e18c1SBarry Smith   ierr = PetscFunctionListFind(((PetscObject)pc)->comm,GAMGList,type,PETSC_FALSE,(PetscVoidStarFunction)&r);CHKERRQ(ierr);
13629d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
13639d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1364676e1743SMark F. Adams   PetscFunctionReturn(0);
1365676e1743SMark F. Adams }
1366676e1743SMark F. Adams EXTERN_C_END
1367676e1743SMark F. Adams 
13685b89ad90SMark F. Adams #undef __FUNCT__
13695b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
13705b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
13715b89ad90SMark F. Adams {
1372676e1743SMark F. Adams   PetscErrorCode ierr;
1373676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1374676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1375676e1743SMark F. Adams   PetscBool      flag;
13765e7c91beSJed Brown   PetscInt       two   = 2;
1377b43b03e9SMark F. Adams   MPI_Comm       wcomm = ((PetscObject)pc)->comm;
13785b89ad90SMark F. Adams 
13795b89ad90SMark F. Adams   PetscFunctionBegin;
1380676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
1381676e1743SMark F. Adams   {
138275b74bdaSMark F. Adams     /* -pc_gamg_verbose */
13839d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
13849d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
1385806fa848SBarry Smith                            &pc_gamg->verbose, PETSC_NULL);CHKERRQ(ierr);
13868263b398SMark F. Adams     /* -pc_gamg_repartition */
13878263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
13888263b398SMark F. Adams                             "Repartion coarse grids (false)",
13898263b398SMark F. Adams                             "PCGAMGRepartitioning",
13909d5b6da9SMark F. Adams                             pc_gamg->repart,
13919d5b6da9SMark F. Adams                             &pc_gamg->repart,
1392806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1393dfd5c07aSMark F. Adams     /* -pc_gamg_reuse_interpolation */
1394dfd5c07aSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation",
1395dfd5c07aSMark F. Adams                             "Reuse prolongation operator (true)",
1396dfd5c07aSMark F. Adams                             "PCGAMGReuseProl",
1397dfd5c07aSMark F. Adams                             pc_gamg->reuse_prol,
1398dfd5c07aSMark F. Adams                             &pc_gamg->reuse_prol,
1399dfd5c07aSMark F. Adams                             &flag);CHKERRQ(ierr);
1400ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
1401ffc955d6SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1402ffc955d6SMark F. Adams                             "Use aggregation agragates for GASM smoother (false)",
1403ffc955d6SMark F. Adams                             "PCGAMGUseASMAggs",
1404ffc955d6SMark F. Adams                             pc_gamg->use_aggs_in_gasm,
1405ffc955d6SMark F. Adams                             &pc_gamg->use_aggs_in_gasm,
1406806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1407c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1408676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1409676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1410676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
14119d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
14129d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1413806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1414389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1415389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1416389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1417389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
14189d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
14199d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1420806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1421c20e4228SMark F. Adams     /* -pc_gamg_threshold */
14223542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
14233542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
14243542efc5SMark F. Adams                             "PCGAMGSetThreshold",
14259d5b6da9SMark F. Adams                             pc_gamg->threshold,
14269d5b6da9SMark F. Adams                             &pc_gamg->threshold,
1427806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1428806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
1429806fa848SBarry Smith       ierr = PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1430806fa848SBarry Smith     }
14314ef23d27SMark F. Adams 
14325e7c91beSJed Brown     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,PETSC_NULL);CHKERRQ(ierr);
14334ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
14344ef23d27SMark F. Adams                            "Set number of MG levels",
14354ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
14369d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
14379d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
1438806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1439676e1743SMark F. Adams   }
1440676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
14415b89ad90SMark F. Adams   PetscFunctionReturn(0);
14425b89ad90SMark F. Adams }
14435b89ad90SMark F. Adams 
14445b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
14455b89ad90SMark F. Adams /*MC
1446856830bfSJed Brown      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
14479d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
14485b89ad90SMark F. Adams 
1449280d9858SJed Brown    Options Database Keys:
14505b89ad90SMark F. Adams    Multigrid options(inherited)
1451280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1452280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1453280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
14548c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
14555b89ad90SMark F. Adams 
14565b89ad90SMark F. Adams   Level: intermediate
1457280d9858SJed Brown 
14585b89ad90SMark F. Adams   Concepts: multigrid
14595b89ad90SMark F. Adams 
14605b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1461280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
14625b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
14635b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
14645b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
14655b89ad90SMark F. Adams M*/
14665b89ad90SMark F. Adams EXTERN_C_BEGIN
14675b89ad90SMark F. Adams #undef __FUNCT__
14685b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
14695b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
14705b89ad90SMark F. Adams {
14715b89ad90SMark F. Adams   PetscErrorCode ierr;
14725b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
14735b89ad90SMark F. Adams   PC_MG          *mg;
14740cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14759d5b6da9SMark F. Adams   static long count = 0;
14765ee9c036SSatish Balay #endif
14775b89ad90SMark F. Adams 
14785b89ad90SMark F. Adams   PetscFunctionBegin;
14795b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14805b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
14815b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14825b89ad90SMark F. Adams 
14835b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
14845b89ad90SMark F. Adams   ierr         = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr);
14855b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
1486ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
14875b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
14885b89ad90SMark F. Adams 
14899d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14909d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14919d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14929d5b6da9SMark F. Adams   pc_gamg->data    = 0;
14935b89ad90SMark F. Adams 
14949d5b6da9SMark F. Adams   /* register AMG type */
14959d5b6da9SMark F. Adams   if (!GAMGList) {
1496140e18c1SBarry Smith     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void (*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
1497140e18c1SBarry Smith     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void (*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
14989d5b6da9SMark F. Adams   }
14999d5b6da9SMark F. Adams 
15009d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
15015b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
15025b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
15035b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
15045b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
15055b89ad90SMark F. Adams 
15065b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1507676e1743SMark F. Adams                                            "PCGAMGSetProcEqLim_C",
1508676e1743SMark F. Adams                                            "PCGAMGSetProcEqLim_GAMG",
1509806fa848SBarry Smith                                            PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1510676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1511389730f3SMark F. Adams                                            "PCGAMGSetCoarseEqLim_C",
1512389730f3SMark F. Adams                                            "PCGAMGSetCoarseEqLim_GAMG",
1513806fa848SBarry Smith                                            PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1514389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15158263b398SMark F. Adams                                            "PCGAMGSetRepartitioning_C",
15168263b398SMark F. Adams                                            "PCGAMGSetRepartitioning_GAMG",
1517806fa848SBarry Smith                                            PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
1518676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1519dfd5c07aSMark F. Adams                                            "PCGAMGSetReuseProl_C",
1520dfd5c07aSMark F. Adams                                            "PCGAMGSetReuseProl_GAMG",
1521dfd5c07aSMark F. Adams                                            PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr);
1522dfd5c07aSMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1523ffc955d6SMark F. Adams                                            "PCGAMGSetUseASMAggs_C",
1524ffc955d6SMark F. Adams                                            "PCGAMGSetUseASMAggs_GAMG",
1525806fa848SBarry Smith                                            PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1526ffc955d6SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15273542efc5SMark F. Adams                                            "PCGAMGSetThreshold_C",
15283542efc5SMark F. Adams                                            "PCGAMGSetThreshold_GAMG",
1529806fa848SBarry Smith                                            PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
15303542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15319d5b6da9SMark F. Adams                                            "PCGAMGSetType_C",
15329d5b6da9SMark F. Adams                                            "PCGAMGSetType_GAMG",
1533806fa848SBarry Smith                                            PCGAMGSetType_GAMG);CHKERRQ(ierr);
15349d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1535dfd5c07aSMark F. Adams   pc_gamg->reuse_prol       = PETSC_TRUE;
1536ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1537038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
1538c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit  = 800;
1539a2f3521dSMark F. Adams   pc_gamg->threshold        = 0.001;
15409d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
15419d5b6da9SMark F. Adams   pc_gamg->verbose          = 0;
15429d5b6da9SMark F. Adams   pc_gamg->emax_id          = -1;
15435e7c91beSJed Brown   pc_gamg->eigtarget[0]     = 0.05;
15445e7c91beSJed Brown   pc_gamg->eigtarget[1]     = 1.05;
15459d5b6da9SMark F. Adams 
15460cbbd2e1SMark F. Adams   /* private events */
15470cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1548785cba28SMark F. Adams   if (count++ == 0) {
1549806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1550806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
15510cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
15520cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
15530cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1554806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1555806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1556806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1557806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1558806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1559806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1560806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1561806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1562806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1563806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1564806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1565806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1566f852f58cSMark F. Adams 
15670cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
15680cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
15690cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1570b4fbaa2aSMark F. Adams     /* create timer stages */
1571b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1572b4fbaa2aSMark F. Adams     {
1573b4fbaa2aSMark F. Adams       char     str[32];
1574b4fbaa2aSMark F. Adams       PetscInt lidx;
1575806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1576806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1577b4fbaa2aSMark F. Adams       for (lidx=1; lidx<9; lidx++) {
1578b4fbaa2aSMark F. Adams         sprintf(str,"MG Level %d",lidx);
1579806fa848SBarry Smith         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1580b4fbaa2aSMark F. Adams       }
1581b4fbaa2aSMark F. Adams     }
1582b4fbaa2aSMark F. Adams #endif
1583b4fbaa2aSMark F. Adams   }
1584b4fbaa2aSMark F. Adams #endif
15850cbbd2e1SMark F. Adams   /* general events */
15860cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
1587806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1588806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1589806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1590806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1591806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1592806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1593806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1594806fa848SBarry Smith   ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr);
15950cbbd2e1SMark F. Adams #endif
15960cbbd2e1SMark F. Adams 
159760a6d8f8SMark F. Adams   /* instantiate derived type */
159860a6d8f8SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
159960a6d8f8SMark F. Adams   {
160060a6d8f8SMark F. Adams     char tname[256] = GAMGAGG;
1601806fa848SBarry Smith     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), PETSC_NULL);CHKERRQ(ierr);
160260a6d8f8SMark F. Adams     ierr = PCGAMGSetType(pc, tname);CHKERRQ(ierr);
160360a6d8f8SMark F. Adams   }
160460a6d8f8SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
16055b89ad90SMark F. Adams   PetscFunctionReturn(0);
16065b89ad90SMark F. Adams }
16075b89ad90SMark F. Adams EXTERN_C_END
1608