xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 1147fc2a6bf7922defbbca15fc6c33f234803ff0)
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   }
51878e152fSMark F. Adams 
52a2f3521dSMark F. Adams   PetscFunctionReturn(0);
53a2f3521dSMark F. Adams }
54a2f3521dSMark F. Adams 
55a2f3521dSMark F. Adams /* private 2x2 Mat Nest for Stokes */
56a2f3521dSMark F. Adams typedef struct{
57a2f3521dSMark F. Adams   Mat A11,A21,A12,Amat;
58a2f3521dSMark F. Adams   IS  prim_is,constr_is;
59a2f3521dSMark F. Adams }GAMGKKTMat;
60a2f3521dSMark F. Adams 
61a2f3521dSMark F. Adams #undef __FUNCT__
62a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatCreate"
63a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatCreate(Mat A, PetscBool iskkt, GAMGKKTMat *out)
64a2f3521dSMark F. Adams {
65a2f3521dSMark F. Adams   PetscFunctionBegin;
66a2f3521dSMark F. Adams   out->Amat = A;
67a2f3521dSMark F. Adams   if (iskkt) {
68a2f3521dSMark F. Adams     PetscErrorCode   ierr;
69a2f3521dSMark F. Adams     IS       is_constraint, is_prime;
70a2f3521dSMark F. Adams     PetscInt nmin,nmax;
71a2f3521dSMark F. Adams 
72a2f3521dSMark F. Adams     ierr = MatGetOwnershipRange(A, &nmin, &nmax);CHKERRQ(ierr);
73a2f3521dSMark F. Adams     ierr = MatFindZeroDiagonals(A, &is_constraint);CHKERRQ(ierr);
74a2f3521dSMark F. Adams     ierr = ISComplement(is_constraint, nmin, nmax, &is_prime);CHKERRQ(ierr);
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);
106a2f3521dSMark F. Adams 
107d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
108d3d6bff4SMark F. Adams }
109d3d6bff4SMark F. Adams 
1105b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1115b89ad90SMark F. Adams /*
112a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
113a147abb0SMark F. Adams      of active processors.
1145b89ad90SMark F. Adams 
1155b89ad90SMark F. Adams    Input Parameter:
116a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
117a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
1189d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
119c5bfad50SMark F. Adams    . cr_bs - coarse block size
1209d5b6da9SMark F. Adams    . isLast -
121a2f3521dSMark F. Adams    . stokes -
1223530afc2SMark F. Adams    In/Output Parameter:
123a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
124afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
12511e60469SMark F. Adams    Output Parameter:
1263530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1275b89ad90SMark F. Adams */
1285cb416c2SMark F. Adams 
1295b89ad90SMark F. Adams #undef __FUNCT__
1308263b398SMark F. Adams #define __FUNCT__ "createLevel"
1310cbbd2e1SMark F. Adams static PetscErrorCode createLevel(const PC pc,
1329d5b6da9SMark F. Adams                                    const Mat Amat_fine,
133c5bfad50SMark F. Adams                                    const PetscInt cr_bs,
1349d5b6da9SMark F. Adams                                    const PetscBool isLast,
135a2f3521dSMark F. Adams                                    const PetscBool stokes,
1363530afc2SMark F. Adams                                    Mat *a_P_inout,
137a2f3521dSMark F. Adams                                    Mat *a_Amat_crs,
138*1147fc2aSKarl Rupp                                    PetscMPIInt *a_nactive_proc)
1395b89ad90SMark F. Adams {
140a2f3521dSMark F. Adams   PetscErrorCode   ierr;
1419d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
142486a8d0bSJed Brown   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1439d5b6da9SMark F. Adams   const PetscBool  repart = pc_gamg->repart;
1449d5b6da9SMark F. Adams   const PetscInt   min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
145a2f3521dSMark F. Adams   Mat              Cmat,Pold=*a_P_inout;
1469d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)Amat_fine)->comm;
1475a9b9e01SMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive=*a_nactive_proc;
148c5bfad50SMark F. Adams   PetscInt         ncrs_eq,ncrs_prim,f_bs;
1495b89ad90SMark F. Adams 
1505b89ad90SMark F. Adams   PetscFunctionBegin;
15111e60469SMark F. Adams   ierr = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr);
15211e60469SMark F. Adams   ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr);
153c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
15411e60469SMark F. Adams   /* RAP */
1559d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
156038e3b61SMark F. Adams 
157a2f3521dSMark F. Adams   /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/
158a2f3521dSMark F. Adams   ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
159a2f3521dSMark F. Adams   assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)==0);
160a2f3521dSMark F. Adams   ierr = MatGetLocalSize(Cmat, &ncrs_eq, PETSC_NULL);CHKERRQ(ierr);
161a2f3521dSMark F. Adams 
162a2f3521dSMark F. Adams   /* get number of PEs to make active 'new_npe', reduce, can be any integer 1-P */
163a2f3521dSMark F. Adams   {
164472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
165a2f3521dSMark F. Adams     ierr = MatGetSize(Cmat, &ncrs_eq_glob, PETSC_NULL);CHKERRQ(ierr);
166472110cdSMark F. Adams     new_npe = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
167472110cdSMark F. Adams     if (new_npe == 0 || ncrs_eq_glob < coarse_max) new_npe = 1;
1685a9b9e01SMark F. Adams     else if (new_npe >= nactive) new_npe = nactive; /* no change, rare */
1699d5b6da9SMark F. Adams     if (isLast) new_npe = 1;
170a2f3521dSMark F. Adams   }
171f852f58cSMark F. Adams 
1725a9b9e01SMark F. Adams   if (!repart && new_npe==nactive) {
173a2f3521dSMark F. Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
174806fa848SBarry Smith   } else {
175a2f3521dSMark 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;
176a2f3521dSMark F. Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
177a2f3521dSMark F. Adams     IS              is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
1785a9b9e01SMark F. Adams     VecScatter      vecscat;
17922063be5SMark F. Adams     PetscScalar    *array;
18022063be5SMark F. Adams     Vec             src_crd, dest_crd;
181e33ef3b1SMark F. Adams 
182c5bfad50SMark F. Adams     nloc_old = ncrs_eq/cr_bs; assert(ncrs_eq%cr_bs==0);
1830cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1840cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
185b4fbaa2aSMark F. Adams #endif
186a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
187a2f3521dSMark F. Adams     ierr = PetscMalloc(npe*sizeof(PetscInt), &counts);CHKERRQ(ierr);
188a2f3521dSMark F. Adams     if (repart && !stokes) {
189a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1905a9b9e01SMark F. Adams       Mat             adj;
1915a9b9e01SMark F. Adams 
192a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
193a2f3521dSMark F. Adams         if (pc_gamg->verbose==1) PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,ncrs_eq);
194a2f3521dSMark F. Adams         else {
195a2f3521dSMark F. Adams           PetscInt n;
196a2f3521dSMark F. Adams           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm);CHKERRQ(ierr);
197a2f3521dSMark F. Adams           PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,n);
198a2f3521dSMark F. Adams         }
199a2f3521dSMark F. Adams       }
2005a9b9e01SMark F. Adams 
201a2f3521dSMark F. Adams       /* get 'adj' */
202c5bfad50SMark F. Adams       if (cr_bs == 1) {
203038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
204806fa848SBarry Smith       } else {
205a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
206eb07cef2SMark F. Adams         Mat tMat;
207a2f3521dSMark F. Adams         PetscInt Istart_crs,Iend_crs,ncols,jj,Ii;
208b4fbaa2aSMark F. Adams         const PetscScalar *vals;
209b4fbaa2aSMark F. Adams         const PetscInt *idx;
210a2f3521dSMark F. Adams         PetscInt *d_nnz, *o_nnz, M, N;
2119057884aSMark F. Adams         static PetscInt llev = 0;
212b4fbaa2aSMark F. Adams 
213a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &d_nnz);CHKERRQ(ierr);
214a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &o_nnz);CHKERRQ(ierr);
215a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
216a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
217c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0 ; Ii < Iend_crs ; Ii += cr_bs, jj++) {
21858471d46SMark F. Adams           ierr = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
219c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
220c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
22158471d46SMark F. Adams           ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
222a2f3521dSMark F. Adams           if (d_nnz[jj] > ncrs_prim) d_nnz[jj] = ncrs_prim;
223c5bfad50SMark F. Adams           if (o_nnz[jj] > (M/cr_bs-ncrs_prim)) o_nnz[jj] = M/cr_bs-ncrs_prim;
22458471d46SMark F. Adams         }
2256876a03eSMark F. Adams 
226a2f3521dSMark F. Adams         ierr = MatCreate(wcomm, &tMat);CHKERRQ(ierr);
227806fa848SBarry Smith         ierr = MatSetSizes(tMat, ncrs_prim, ncrs_prim,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
228a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
229a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
230a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
23158471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
2325f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
233eb07cef2SMark F. Adams 
234a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
235c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
23622063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
237eb07cef2SMark F. Adams           for (jj = 0 ; jj < ncols ; jj++) {
238c5bfad50SMark F. Adams             PetscInt dest_col = idx[jj]/cr_bs;
239eb07cef2SMark F. Adams             PetscScalar v = 1.0;
240eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
241eb07cef2SMark F. Adams           }
24222063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
243eb07cef2SMark F. Adams         }
244eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
245eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
246eb07cef2SMark F. Adams 
247b4fbaa2aSMark F. Adams         if (llev++ == -1) {
248b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2498caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
250b4fbaa2aSMark F. Adams           PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
251b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2523bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
253b4fbaa2aSMark F. Adams         }
254b4fbaa2aSMark F. Adams 
255eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
256eb07cef2SMark F. Adams 
257eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
258a2f3521dSMark F. Adams       } /* create 'adj' */
259f150b916SMark F. Adams 
260a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2615a9b9e01SMark F. Adams         char prefix[256];
2625a9b9e01SMark F. Adams         const char *pcpre;
263b4fbaa2aSMark F. Adams         const PetscInt *is_idx;
264b4fbaa2aSMark F. Adams         MatPartitioning  mpart;
265a4b7d37bSMark F. Adams         IS proc_is;
266a2f3521dSMark F. Adams         PetscInt targetPE;
2672f03bc48SMark F. Adams 
2685a9b9e01SMark F. Adams         ierr = MatPartitioningCreate(wcomm, &mpart);CHKERRQ(ierr);
2695ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2709d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2718caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
27259a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
27311e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
2745ef31b24SMark F. Adams         ierr = MatPartitioningSetNParts(mpart, new_npe);CHKERRQ(ierr);
275a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
27611e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2775a9b9e01SMark F. Adams 
2785ef31b24SMark F. Adams         /* collect IS info */
279a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_eq*sizeof(PetscInt), &newproc_idx);CHKERRQ(ierr);
280a4b7d37bSMark F. Adams         ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
281a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
282a2f3521dSMark F. Adams         /*targetPE = npe/new_npe;*/ /* spread partitioning across machine */
283a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
284c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
285a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
286eb07cef2SMark F. Adams           }
2875ef31b24SMark F. Adams         }
288a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
289a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2905ef31b24SMark F. Adams       }
2915ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2925a9b9e01SMark F. Adams 
293806fa848SBarry Smith       ierr = ISCreateGeneral(wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2948263b398SMark F. Adams       if (newproc_idx != 0) {
2958263b398SMark F. Adams         ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2965ef31b24SMark F. Adams       }
297806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
298a2f3521dSMark F. Adams 
299a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
3005a9b9e01SMark F. Adams       /* find factor */
3015a9b9e01SMark F. Adams       if (new_npe == 1) rfactor = npe; /* easy */
3025a9b9e01SMark F. Adams       else {
3035a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
3045a9b9e01SMark F. Adams         jj = -1;
3055a9b9e01SMark F. Adams         for (kk = 1 ; kk <= npe ; kk++) {
3065a9b9e01SMark F. Adams           if (npe%kk==0) { /* a candidate */
3075a9b9e01SMark F. Adams             PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe;
3085a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
3095a9b9e01SMark F. Adams             if (fact > best_fact) {
3105a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
3115a9b9e01SMark F. Adams             }
3125a9b9e01SMark F. Adams           }
3135a9b9e01SMark F. Adams         }
3145a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
315a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
3165a9b9e01SMark F. Adams       }
3175a9b9e01SMark F. Adams       new_npe = npe/rfactor;
3185a9b9e01SMark F. Adams 
3195a9b9e01SMark F. Adams       if (new_npe==nactive) {
320a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
3215a9b9e01SMark F. Adams         ierr = PetscFree(counts);CHKERRQ(ierr);
322a2f3521dSMark F. Adams         if (pc_gamg->verbose>0) {
323a2f3521dSMark F. Adams           PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq(loc)=%d\n",mype,__FUNCT__,new_npe,ncrs_eq);
324a2f3521dSMark F. Adams         }
3255a9b9e01SMark F. Adams         PetscFunctionReturn(0);
3265a9b9e01SMark F. Adams       }
3275a9b9e01SMark F. Adams 
328a2f3521dSMark F. Adams       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",mype,__FUNCT__,ncrs_eq);
3295a9b9e01SMark F. Adams       targetPE = mype/rfactor;
330a2f3521dSMark F. Adams       ierr = ISCreateStride(wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
3315a9b9e01SMark F. Adams 
332a2f3521dSMark F. Adams       if (stokes) {
333c5bfad50SMark F. Adams         ierr = ISCreateStride(wcomm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim);CHKERRQ(ierr);
3345a9b9e01SMark F. Adams       }
335a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
336e33ef3b1SMark F. Adams 
33711e60469SMark F. Adams     /*
338a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
33911e60469SMark F. Adams      */
340a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
341a2f3521dSMark F. Adams     if (stokes) {
342a2f3521dSMark F. Adams       ierr = ISPartitioningToNumbering(is_eq_newproc_prim, &is_eq_num_prim);CHKERRQ(ierr);
343806fa848SBarry Smith     } else is_eq_num_prim = is_eq_num;
34411e60469SMark F. Adams     /*
345a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
34611e60469SMark F. Adams      */
347a2f3521dSMark F. Adams     ierr = ISPartitioningCount(is_eq_newproc, npe, counts);CHKERRQ(ierr);
348a2f3521dSMark F. Adams     ncrs_eq_new = counts[mype];
349a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
350a2f3521dSMark F. Adams     if (stokes) {
351a2f3521dSMark F. Adams       ierr = ISPartitioningCount(is_eq_newproc_prim, npe, counts);CHKERRQ(ierr);
352a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_newproc_prim);CHKERRQ(ierr);
353c5bfad50SMark F. Adams       ncrs_prim_new = counts[mype]/cr_bs; /* nodes */
354806fa848SBarry Smith     } else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */
355a2f3521dSMark F. Adams 
356a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
3570cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3580cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
359b4fbaa2aSMark F. Adams #endif
360a2f3521dSMark F. Adams 
361a2f3521dSMark F. Adams     /* move data (for primal equations only) */
36222063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
3633bf036e2SBarry Smith     ierr = VecCreate(wcomm, &dest_crd);CHKERRQ(ierr);
364a2f3521dSMark F. Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE);CHKERRQ(ierr);
365470e26f8SMark F. Adams     ierr = VecSetFromOptions(dest_crd);CHKERRQ(ierr); /* this is needed! */
36611e60469SMark F. Adams     /*
3679d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
368c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
36911e60469SMark F. Adams      */
370a2f3521dSMark F. Adams     ierr = PetscMalloc((ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx);CHKERRQ(ierr);
371a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
372a2f3521dSMark F. Adams     for (ii=0,jj=0; ii<ncrs_prim ; ii++) {
373c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
374a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz ; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
37511e60469SMark F. Adams     }
376a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
377806fa848SBarry Smith     ierr = ISCreateGeneral(wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
37892a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
37911e60469SMark F. Adams     /*
38011e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
38111e60469SMark F. Adams      */
382a2f3521dSMark F. Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd);CHKERRQ(ierr);
3839d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols ; jj++) {
384a2f3521dSMark F. Adams       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
385a2f3521dSMark F. Adams       for (ii=0 ; ii<ncrs_prim ; ii++) {
3869d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows ; kk++) {
387a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
388c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
389676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
390d3d6bff4SMark F. Adams         }
391038e3b61SMark F. Adams       }
392eb07cef2SMark F. Adams     }
393eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
394eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
39511e60469SMark F. Adams     /*
39611e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
39711e60469SMark F. Adams       to the correct processor
39811e60469SMark F. Adams     */
399806fa848SBarry Smith     ierr = VecScatterCreate(src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
40011e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
40111e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40211e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40311e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
40411e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
40511e60469SMark F. Adams     /*
40611e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
40711e60469SMark F. Adams     */
408c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
409a2f3521dSMark F. Adams     ierr = PetscMalloc(node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
410a2f3521dSMark F. Adams     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
411a2f3521dSMark F. Adams     strideNew = ncrs_prim_new*ndata_rows;
41211e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
4139d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols ; jj++) {
414a2f3521dSMark F. Adams       for (ii=0 ; ii<ncrs_prim_new ; ii++) {
4159d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows ; kk++) {
416a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
417c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
418d3d6bff4SMark F. Adams         }
419038e3b61SMark F. Adams       }
420038e3b61SMark F. Adams     }
42111e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
42211e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
423a2f3521dSMark F. Adams 
424a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
4250cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4260cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
427ed3f9983SMark F. Adams #endif
428a2f3521dSMark F. Adams 
42911e60469SMark F. Adams     /*
43011e60469SMark F. Adams       Invert for MatGetSubMatrix
43111e60469SMark F. Adams     */
432a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
433a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
434c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
435a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
436a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
437a2f3521dSMark F. Adams     }
438a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
4390cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4400cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
4410cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
442ed3f9983SMark F. Adams #endif
443a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
444a2f3521dSMark F. Adams     {
445a2f3521dSMark F. Adams       Mat mat;
446806fa848SBarry Smith       ierr = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
447a2f3521dSMark F. Adams       *a_Amat_crs = mat;
448c5bfad50SMark F. Adams 
449c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
450c5bfad50SMark F. Adams         PetscInt cbs, rbs;
451c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
452806fa848SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
453c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
454806fa848SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",mype,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr);
455c5bfad50SMark F. Adams       }
456a2f3521dSMark F. Adams     }
457038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
458a2f3521dSMark F. Adams 
4590cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4600cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
461ed3f9983SMark F. Adams #endif
46211e60469SMark F. Adams     /* prolongator */
46311e60469SMark F. Adams     {
46411e60469SMark F. Adams       IS findices;
465a2f3521dSMark F. Adams       PetscInt Istart,Iend;
466a2f3521dSMark F. Adams       Mat Pnew;
467a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4680cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4690cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
470ed3f9983SMark F. Adams #endif
47111e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
472c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
473806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
47411e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
475c5bfad50SMark F. Adams 
476c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
477c5bfad50SMark F. Adams         PetscInt cbs, rbs;
478c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
479806fa848SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
480c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
481806fa848SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
482c5bfad50SMark F. Adams       }
4830cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4840cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
485ed3f9983SMark F. Adams #endif
4863530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
487a2f3521dSMark F. Adams 
488a2f3521dSMark F. Adams       /* output - repartitioned */
489a2f3521dSMark F. Adams       *a_P_inout = Pnew;
490e33ef3b1SMark F. Adams     }
491a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4925b89ad90SMark F. Adams 
4935a9b9e01SMark F. Adams     *a_nactive_proc = new_npe; /* output */
494a2f3521dSMark F. Adams   }
4955a9b9e01SMark F. Adams 
496a2f3521dSMark F. Adams   /* outout matrix data */
497c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
498c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
499c8b0795cSMark F. Adams     if (llev==0) {
500c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
501c8b0795cSMark F. Adams       PetscViewerASCIIOpen(wcomm,fname,&viewer);
502c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
503c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
504c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
505c8b0795cSMark F. Adams     }
506c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
507c8b0795cSMark F. Adams     PetscViewerASCIIOpen(wcomm,fname,&viewer);
508c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
509c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
510c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
511c8b0795cSMark F. Adams   }
512c8b0795cSMark F. Adams 
5135b89ad90SMark F. Adams   PetscFunctionReturn(0);
5145b89ad90SMark F. Adams }
5155b89ad90SMark F. Adams 
5165b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5175b89ad90SMark F. Adams /*
5185b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5195b89ad90SMark F. Adams                     by setting data structures and options.
5205b89ad90SMark F. Adams 
5215b89ad90SMark F. Adams    Input Parameter:
5225b89ad90SMark F. Adams .  pc - the preconditioner context
5235b89ad90SMark F. Adams 
5245b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
5255b89ad90SMark F. Adams 
5265b89ad90SMark F. Adams    Notes:
5275b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
5285b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
5295b89ad90SMark F. Adams */
5305b89ad90SMark F. Adams #undef __FUNCT__
5315b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
5329d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5335b89ad90SMark F. Adams {
5345b89ad90SMark F. Adams   PetscErrorCode  ierr;
5359d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
5365b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
5372adcac29SMark F. Adams   Mat              Pmat = pc->pmat;
538a2f3521dSMark F. Adams   PetscInt         fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
5399d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
5403530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
541587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
542c8b0795cSMark F. Adams   PetscReal        emaxs[GAMG_MAXLEVELS];
543e696ed0bSMark F. Adams   IS              *ASMLocalIDsArr[GAMG_MAXLEVELS];
544a2f3521dSMark F. Adams   GAMGKKTMat       kktMatsArr[GAMG_MAXLEVELS];
545a2f3521dSMark F. Adams   PetscLogDouble   nnz0=0.,nnztot=0.;
546737a81a9SMark F. Adams   MatInfo          info;
5475169fedaSMark F. Adams   PetscBool        stokes = PETSC_FALSE, redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);
5485ef31b24SMark F. Adams 
5495b89ad90SMark F. Adams   PetscFunctionBegin;
55084d3f75bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
55184d3f75bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
552dfd5c07aSMark F. Adams 
553a2f3521dSMark F. Adams   if (pc_gamg->verbose>2) PetscPrintf(wcomm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",mype,__FUNCT__,pc_gamg->setup_count,pc->setupcalled);
554dfd5c07aSMark F. Adams 
55584d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
556878e152fSMark F. Adams     if (redo_mesh_setup) {
557878e152fSMark F. Adams       /* reset everything */
558878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
559878e152fSMark F. Adams       pc->setupcalled = 0;
560806fa848SBarry Smith     } else {
56184d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
56203a628feSMark F. Adams       /* just do Galerkin grids */
56358471d46SMark F. Adams       Mat B,dA,dB;
564d5280255SMark F. Adams       assert(pc->setupcalled);
56558471d46SMark F. Adams 
5669d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
56758471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5689d5b6da9SMark F. Adams         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
56958471d46SMark F. Adams         /* (re)set to get dirty flag */
5709d5b6da9SMark F. Adams         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
57158471d46SMark F. Adams 
5729d5b6da9SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>-1; level--) {
57303a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
57484d3f75bSMark F. Adams           if (pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
57503a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
576084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
57703a628feSMark F. Adams             mglevels[level]->A = B;
578806fa848SBarry Smith           } else {
579084a8fe3SJed Brown             ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
58058471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
58103a628feSMark F. Adams           }
58258471d46SMark F. Adams           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
58358471d46SMark F. Adams           dB = B;
58458471d46SMark F. Adams         }
5855f8cf99dSMark F. Adams       }
586d5280255SMark F. Adams 
587d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
588d5280255SMark F. Adams 
589d5280255SMark F. Adams       /* PCSetUp_MG seems to insists on setting this to GMRES */
5901ac9965eSMark F. Adams       ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr);
591d5280255SMark F. Adams 
59258471d46SMark F. Adams       PetscFunctionReturn(0);
593eb07cef2SMark F. Adams     }
594878e152fSMark F. Adams   }
595f6536408SMark F. Adams 
596302f38e8SMark F. Adams   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
597eb07cef2SMark F. Adams 
5989d1b15c3SMark F. Adams   ierr = GAMGKKTMatCreate(Pmat, stokes, &kktMatsArr[0]);CHKERRQ(ierr);
5999d1b15c3SMark F. Adams 
600878e152fSMark F. Adams   if (!pc_gamg->data) {
601878e152fSMark F. Adams     if (pc_gamg->orig_data) {
602878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
603878e152fSMark F. Adams       ierr = MatGetLocalSize(Pmat, &qq, PETSC_NULL);CHKERRQ(ierr);
604878e152fSMark F. Adams       pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
605878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
606878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
607878e152fSMark F. Adams       ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
608878e152fSMark F. Adams       for (qq=0;qq<pc_gamg->data_sz;qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
609806fa848SBarry Smith     } else {
6109d5b6da9SMark F. Adams       if (!pc_gamg->createdefaultdata) {
611c666adbfSMark F. Adams         SETERRQ(wcomm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
612302f38e8SMark F. Adams       }
613302f38e8SMark F. Adams       if (stokes) {
614c666adbfSMark F. Adams         SETERRQ(wcomm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
615eb07cef2SMark F. Adams       }
6169d1b15c3SMark F. Adams       ierr = pc_gamg->createdefaultdata(pc, kktMatsArr[0].A11);CHKERRQ(ierr);
6179d5b6da9SMark F. Adams     }
618878e152fSMark F. Adams   }
619878e152fSMark F. Adams 
620878e152fSMark F. Adams   /* cache original data for reuse */
621878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
622878e152fSMark F. Adams     ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr);
623878e152fSMark F. Adams     for (qq=0;qq<pc_gamg->data_sz;qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
624878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
625878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
626878e152fSMark F. Adams   }
627038e3b61SMark F. Adams 
628302f38e8SMark F. Adams   /* get basic dims */
629302f38e8SMark F. Adams   if (stokes) {
630a2f3521dSMark F. Adams     bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
631806fa848SBarry Smith   } else {
632302f38e8SMark F. Adams     ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
633302f38e8SMark F. Adams   }
634a2f3521dSMark F. Adams 
635a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
636c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
63784f9421dSMark F. Adams     PetscInt NN = M;
63884f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
63984f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
6403bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
641806fa848SBarry Smith     } else {
642806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
64384f9421dSMark F. Adams     }
644b2a4f308SMark F. Adams     nnz0 = info.nz_used;
645b2a4f308SMark F. Adams     nnztot = info.nz_used;
646806fa848SBarry 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",
647a2f3521dSMark F. Adams                 mype,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
648806fa848SBarry Smith                 (int)(nnz0/(PetscReal)NN),npe);CHKERRQ(ierr);
649c8b0795cSMark F. Adams   }
65084d3f75bSMark F. Adams 
651a2f3521dSMark F. Adams   /* Get A_i and R_i */
6528f4b7eb5SMark F. Adams   for (level=0, Aarr[0]=Pmat, nactivepe = npe; /* hard wired stopping logic */
6539d5b6da9SMark F. Adams         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
6540205a208SMark F. Adams         level++) {
6555b89ad90SMark F. Adams     level1 = level + 1;
6560cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6570cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
658a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
659a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
660b4fbaa2aSMark F. Adams #endif
661a2f3521dSMark F. Adams #endif
662a2f3521dSMark F. Adams     /* deal with Stokes, get sub matrices */
6639d1b15c3SMark F. Adams     if (level > 0) {
664a2f3521dSMark F. Adams       ierr = GAMGKKTMatCreate(Aarr[level], stokes, &kktMatsArr[level]);CHKERRQ(ierr);
6659d1b15c3SMark F. Adams     }
666c8b0795cSMark F. Adams     { /* construct prolongator */
667725b86d8SJed Brown       Mat Gmat;
6680cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
669a2f3521dSMark F. Adams       Mat Prol11,Prol22;
670c8b0795cSMark F. Adams 
671a2f3521dSMark F. Adams       ierr = pc_gamg->graph(pc,kktMatsArr[level].A11, &Gmat);CHKERRQ(ierr);
6720cbbd2e1SMark F. Adams       ierr = pc_gamg->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
673a2f3521dSMark F. Adams       ierr = pc_gamg->prolongator(pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11);CHKERRQ(ierr);
674c8b0795cSMark F. Adams 
675a2f3521dSMark F. Adams       /* could have failed to create new level */
676a2f3521dSMark F. Adams       if (Prol11) {
6779d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6789d1b15c3SMark F. Adams         ierr = MatGetBlockSizes(Prol11, PETSC_NULL, &bs);CHKERRQ(ierr);
679a2f3521dSMark F. Adams 
680a2f3521dSMark F. Adams         if (stokes) {
681a2f3521dSMark F. Adams           if (!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
682a2f3521dSMark F. Adams           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
683a2f3521dSMark F. Adams           ierr = pc_gamg->formkktprol(pc, Prol11, kktMatsArr[level].A21, &Prol22);CHKERRQ(ierr);
684c8b0795cSMark F. Adams         }
685c8b0795cSMark F. Adams 
686a2f3521dSMark F. Adams         if (pc_gamg->optprol) {
687c8b0795cSMark F. Adams           /* smooth */
688a2f3521dSMark F. Adams           ierr = pc_gamg->optprol(pc, kktMatsArr[level].A11, &Prol11);CHKERRQ(ierr);
689c8b0795cSMark F. Adams         }
690c8b0795cSMark F. Adams 
691a2f3521dSMark F. Adams         if (stokes) {
692a2f3521dSMark F. Adams           IS is_row[2] = {kktMatsArr[level].prim_is,kktMatsArr[level].constr_is};
693a2f3521dSMark F. Adams           Mat a[4] = {Prol11, PETSC_NULL, PETSC_NULL, Prol22 };
694a2f3521dSMark F. Adams           ierr = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1]);CHKERRQ(ierr);
695806fa848SBarry Smith         } else {
696a2f3521dSMark F. Adams           Parr[level1] = Prol11;
697a2f3521dSMark F. Adams         }
698806fa848SBarry Smith       } else Parr[level1] = PETSC_NULL;
699ffc955d6SMark F. Adams 
700ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
701806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
702ffc955d6SMark F. Adams       }
703ffc955d6SMark F. Adams 
704a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
70541b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
706a2f3521dSMark F. Adams     } /* construct prolongator scope */
7070cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7080cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
709c8b0795cSMark F. Adams #endif
7109d5b6da9SMark F. Adams     /* cache eigen estimate */
7119d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
7129d5b6da9SMark F. Adams       PetscBool flag;
713806fa848SBarry Smith       ierr = PetscObjectComposedDataGetReal((PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
7149d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
715806fa848SBarry Smith     } else emaxs[level] = -1.;
7162adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
717c8b0795cSMark F. Adams     if (!Parr[level1]) {
718806fa848SBarry Smith       if (pc_gamg->verbose) {
7193bf036e2SBarry Smith         ierr =  PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);CHKERRQ(ierr);
720806fa848SBarry Smith       }
721c8b0795cSMark F. Adams       break;
722c8b0795cSMark F. Adams     }
7230cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7240cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
725b4fbaa2aSMark F. Adams #endif
726a2f3521dSMark F. Adams 
727a2f3521dSMark F. Adams     ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
728806fa848SBarry Smith                         stokes, &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr);
729a2f3521dSMark F. Adams 
7300cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7310cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
732b4fbaa2aSMark F. Adams #endif
733a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
734a2f3521dSMark F. Adams 
735a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0) {
7360cbbd2e1SMark F. Adams       PetscInt NN = M;
7370cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
738a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
7393bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
740806fa848SBarry Smith       } else {
741806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
7420cbbd2e1SMark F. Adams       }
743a2f3521dSMark F. Adams 
7440cbbd2e1SMark F. Adams       nnztot += info.nz_used;
745806fa848SBarry Smith       ierr = PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
746a2f3521dSMark F. Adams                   mype,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
747806fa848SBarry Smith                   (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
748c8b0795cSMark F. Adams     }
749a2f3521dSMark F. Adams 
750a2f3521dSMark F. Adams     /* stop if one node -- could pull back for singular problems */
751c8b0795cSMark F. Adams     if (M/pc_gamg->data_cell_cols < 2) {
75281708dccSMark F. Adams       level++;
75381708dccSMark F. Adams       break;
75481708dccSMark F. Adams     }
7550cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
756b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
757b4fbaa2aSMark F. Adams #endif
758c8b0795cSMark F. Adams   } /* levels */
759c8b0795cSMark F. Adams 
760c8b0795cSMark F. Adams   if (pc_gamg->data) {
761c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
762a2f3521dSMark F. Adams     pc_gamg->data = PETSC_NULL;
7635b89ad90SMark F. Adams   }
764c8b0795cSMark F. Adams 
7650cbbd2e1SMark F. Adams   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
7669d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7675b89ad90SMark F. Adams   fine_level = level;
7689d5b6da9SMark F. Adams   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
7695b89ad90SMark F. Adams 
77084d3f75bSMark F. Adams   /* simple setup */
77184d3f75bSMark F. Adams   if (!PETSC_TRUE) {
77284d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
77384d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
77484d3f75bSMark F. Adams          lidx<fine_level;
77584d3f75bSMark F. Adams          lidx++, level--) {
77684d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
77784d3f75bSMark F. Adams       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77884d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
77984d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
78084d3f75bSMark F. Adams     }
78184d3f75bSMark F. Adams     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
78284d3f75bSMark F. Adams 
78384d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
784806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
785d5280255SMark F. Adams     /* set default smoothers & set operators */
7869d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
787587fa25dSMark F. Adams           lidx <= fine_level;
788587fa25dSMark F. Adams           lidx++, level--) {
789ffc955d6SMark F. Adams       KSP smoother;
790ffc955d6SMark F. Adams       PC subpc;
791a2f3521dSMark F. Adams 
7929d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
793f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
794ffc955d6SMark F. Adams 
795a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
796a2f3521dSMark F. Adams       /* set ops */
797a2f3521dSMark F. Adams       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
798a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
799a2f3521dSMark F. Adams 
800a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
801a2f3521dSMark F. Adams       if (stokes) {
802a2f3521dSMark F. Adams         KSP *ksps;
803a2f3521dSMark F. Adams         PetscInt nn;
804a2f3521dSMark F. Adams         ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);CHKERRQ(ierr);
805a2f3521dSMark F. Adams         ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is);CHKERRQ(ierr);
806a2f3521dSMark F. Adams         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
807a2f3521dSMark F. Adams         smoother = ksps[0];
808a2f3521dSMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
809a2f3521dSMark F. Adams         ierr = PetscFree(ksps);CHKERRQ(ierr);
810a2f3521dSMark F. Adams       }
811a2f3521dSMark F. Adams       ierr = GAMGKKTMatDestroy(&kktMatsArr[level]);CHKERRQ(ierr);
812a2f3521dSMark F. Adams 
813a2f3521dSMark F. Adams       /* set defaults */
8146c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
815a2f3521dSMark F. Adams 
816ffc955d6SMark F. Adams       /* override defaults and command line args (!) */
817ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
8182d3561bbSSatish Balay         PetscInt sz;
8192d3561bbSSatish Balay         IS *is;
820a2f3521dSMark F. Adams 
8212d3561bbSSatish Balay         sz = nASMBlocksArr[level];
8222d3561bbSSatish Balay         is = ASMLocalIDsArr[level];
823ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
824ffc955d6SMark F. Adams         if (sz==0) {
825ffc955d6SMark F. Adams           IS is;
826ffc955d6SMark F. Adams           PetscInt my0,kk;
827ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
828ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
82906b43e7eSDmitry Karpeev           ierr = PCGASMSetSubdomains(subpc, 1, &is, PETSC_NULL);CHKERRQ(ierr);
830a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
831806fa848SBarry Smith         } else {
832a94c3b12SMark F. Adams           PetscInt kk;
83306b43e7eSDmitry Karpeev           ierr = PCGASMSetSubdomains(subpc, sz, is, PETSC_NULL);CHKERRQ(ierr);
834a94c3b12SMark F. Adams           for (kk=0;kk<sz;kk++) {
835a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
836a94c3b12SMark F. Adams           }
837ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
838ffc955d6SMark F. Adams         }
839ffc955d6SMark F. Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
840ffc955d6SMark F. Adams 
841ffc955d6SMark F. Adams         ASMLocalIDsArr[level] = PETSC_NULL;
842ffc955d6SMark F. Adams         nASMBlocksArr[level] = 0;
843ffc955d6SMark F. Adams         ierr = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
844806fa848SBarry Smith       } else {
8459d5b6da9SMark F. Adams         ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr);
846ffc955d6SMark F. Adams       }
847d5280255SMark F. Adams     }
848d5280255SMark F. Adams     {
849d5280255SMark F. Adams       /* coarse grid */
850d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
851d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
852d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
853d5280255SMark F. Adams       ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
854d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
855d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
856d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
857d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
858d5280255SMark F. Adams       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
859d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
860d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
861d5280255SMark F. Adams     }
862d5280255SMark F. Adams 
863d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
864d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
865d5280255SMark F. Adams     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr);
866d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
867d5280255SMark 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.");
868d5280255SMark F. Adams 
869d5280255SMark F. Adams     /* create cheby smoothers */
870d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
871d5280255SMark F. Adams           lidx <= fine_level;
872d5280255SMark F. Adams           lidx++, level--) {
873d5280255SMark F. Adams       KSP smoother;
874ffc955d6SMark F. Adams       PetscBool flag;
875d5280255SMark F. Adams       PC subpc;
876d5280255SMark F. Adams 
877ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
878a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
879a2f3521dSMark F. Adams 
880a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
881a2f3521dSMark F. Adams       if (stokes) {
882a2f3521dSMark F. Adams         KSP *ksps;
883a2f3521dSMark F. Adams         PetscInt nn;
884a2f3521dSMark F. Adams         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
885a2f3521dSMark F. Adams         smoother = ksps[0];
886a2f3521dSMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
887a2f3521dSMark F. Adams         ierr = PetscFree(ksps);CHKERRQ(ierr);
888a2f3521dSMark F. Adams       }
889ffc955d6SMark F. Adams 
890ffc955d6SMark F. Adams       /* do my own cheby */
8916c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
892ffc955d6SMark F. Adams       if (flag) {
893ffc955d6SMark F. Adams         PetscReal emax, emin;
894251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
895ffc955d6SMark F. Adams         if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
896587fa25dSMark F. Adams         else { /* eigen estimate 'emax' */
897e696ed0bSMark F. Adams           KSP eksp;
898e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
899b2a4f308SMark F. Adams           Vec bb, xx;
900038e3b61SMark F. Adams 
9015745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr);
9025745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr);
903fc4362bfSMark F. Adams           {
904fc4362bfSMark F. Adams             PetscRandom    rctx;
905fc4362bfSMark F. Adams             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
906fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
907fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
908fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
9095b89ad90SMark F. Adams           }
910a2f3521dSMark F. Adams 
911e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
912e696ed0bSMark F. Adams           {
913e696ed0bSMark F. Adams             PetscInt Istart,Iend,ncols,jj,Ii;
914e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
915e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
916e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0 ; Ii < Iend ; Ii++, jj++) {
917e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
918e696ed0bSMark F. Adams               if (ncols <= 1) {
919e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
920a94c3b12SMark F. Adams               }
921e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
922a94c3b12SMark F. Adams             }
923a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
924a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
925a94c3b12SMark F. Adams           }
926e696ed0bSMark F. Adams 
927fc4362bfSMark F. Adams           ierr = KSPCreate(wcomm, &eksp);CHKERRQ(ierr);
928806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
929fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
9301a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
9311a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
932f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
933f6536408SMark F. Adams 
934f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
935f6536408SMark F. Adams           ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
936fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
9375a9b9e01SMark F. Adams 
938d3d0db20SJed Brown           /* set PC type to be same as smoother */
939ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
940b2a4f308SMark F. Adams 
9415a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
9425a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
9435a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
944fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
9455a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
9465a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
9475a9b9e01SMark F. Adams 
948fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
9495a9b9e01SMark F. Adams 
950fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
951fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
952fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
953f6536408SMark F. Adams 
954ffc955d6SMark F. Adams           if (pc_gamg->verbose > 0) {
955a94c3b12SMark F. Adams             PetscInt N1, tt;
956a94c3b12SMark F. Adams             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
957a94c3b12SMark 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);
958f6536408SMark F. Adams           }
959fc4362bfSMark F. Adams         }
960038e3b61SMark F. Adams         {
961c5bfad50SMark F. Adams           PetscInt N1, N0;
962c5bfad50SMark F. Adams           ierr = MatGetSize(Aarr[level], &N1, PETSC_NULL);CHKERRQ(ierr);
963c5bfad50SMark F. Adams           ierr = MatGetSize(Aarr[level+1], &N0, PETSC_NULL);CHKERRQ(ierr);
964f6536408SMark F. Adams           /* heuristic - is this crap? */
965b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
9665e7c91beSJed Brown           emin = emax * pc_gamg->eigtarget[0];
9675e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
968038e3b61SMark F. Adams         }
9696c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
970ffc955d6SMark F. Adams       } /* setup checby flag */
971ffc955d6SMark F. Adams     } /* non-coarse levels */
972737a81a9SMark F. Adams 
973d5280255SMark F. Adams     /* clean up */
974d5280255SMark F. Adams     for (level=1;level<pc_gamg->Nlevels;level++) {
975587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
976587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
9775b89ad90SMark F. Adams     }
9780cbbd2e1SMark F. Adams 
9790cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9800cbbd2e1SMark F. Adams 
9811ac9965eSMark F. Adams     if (PETSC_TRUE) {
98258471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
9839d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
98458471d46SMark F. Adams       ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
98558471d46SMark F. Adams     }
986806fa848SBarry Smith   } else {
9875f8cf99dSMark F. Adams     KSP smoother;
988b43b03e9SMark F. Adams     if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__);
9899d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
9905f8cf99dSMark F. Adams     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9915f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
9929d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9935f8cf99dSMark F. Adams   }
99484d3f75bSMark F. Adams 
9955b89ad90SMark F. Adams   PetscFunctionReturn(0);
9965b89ad90SMark F. Adams }
9975b89ad90SMark F. Adams 
998eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
9995b89ad90SMark F. Adams /*
10005b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
10015b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
10025b89ad90SMark F. Adams 
10035b89ad90SMark F. Adams    Input Parameter:
10045b89ad90SMark F. Adams .  pc - the preconditioner context
10055b89ad90SMark F. Adams 
10065b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
10075b89ad90SMark F. Adams */
10085b89ad90SMark F. Adams #undef __FUNCT__
10095b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
10105b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
10115b89ad90SMark F. Adams {
10125b89ad90SMark F. Adams   PetscErrorCode  ierr;
10135b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
10145b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
10155b89ad90SMark F. Adams 
10165b89ad90SMark F. Adams   PetscFunctionBegin;
10175b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
10185b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
10195b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
10205b89ad90SMark F. Adams   PetscFunctionReturn(0);
10215b89ad90SMark F. Adams }
10225b89ad90SMark F. Adams 
1023676e1743SMark F. Adams 
1024676e1743SMark F. Adams #undef __FUNCT__
1025676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
1026676e1743SMark F. Adams /*@
1027676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
1028676e1743SMark F. Adams    processor reduction.
1029676e1743SMark F. Adams 
1030676e1743SMark F. Adams    Not Collective on PC
1031676e1743SMark F. Adams 
1032676e1743SMark F. Adams    Input Parameters:
1033676e1743SMark F. Adams .  pc - the preconditioner context
1034676e1743SMark F. Adams 
1035676e1743SMark F. Adams 
1036676e1743SMark F. Adams    Options Database Key:
1037676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
1038676e1743SMark F. Adams 
1039676e1743SMark F. Adams    Level: intermediate
1040676e1743SMark F. Adams 
1041676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1042676e1743SMark F. Adams 
1043676e1743SMark F. Adams .seealso: ()
1044676e1743SMark F. Adams @*/
1045676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
1046676e1743SMark F. Adams {
1047676e1743SMark F. Adams   PetscErrorCode ierr;
1048676e1743SMark F. Adams 
1049676e1743SMark F. Adams   PetscFunctionBegin;
1050676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1051676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1052676e1743SMark F. Adams   PetscFunctionReturn(0);
1053676e1743SMark F. Adams }
1054676e1743SMark F. Adams 
1055676e1743SMark F. Adams EXTERN_C_BEGIN
1056676e1743SMark F. Adams #undef __FUNCT__
1057676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
1058676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1059676e1743SMark F. Adams {
1060c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1061c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1062676e1743SMark F. Adams 
1063676e1743SMark F. Adams   PetscFunctionBegin;
10649d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
1065676e1743SMark F. Adams   PetscFunctionReturn(0);
1066676e1743SMark F. Adams }
1067676e1743SMark F. Adams EXTERN_C_END
1068676e1743SMark F. Adams 
1069676e1743SMark F. Adams #undef __FUNCT__
1070389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
1071389730f3SMark F. Adams /*@
1072389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
1073389730f3SMark F. Adams 
1074389730f3SMark F. Adams  Collective on PC
1075389730f3SMark F. Adams 
1076389730f3SMark F. Adams    Input Parameters:
1077389730f3SMark F. Adams .  pc - the preconditioner context
1078389730f3SMark F. Adams 
1079389730f3SMark F. Adams 
1080389730f3SMark F. Adams    Options Database Key:
1081389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
1082389730f3SMark F. Adams 
1083389730f3SMark F. Adams    Level: intermediate
1084389730f3SMark F. Adams 
1085389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1086389730f3SMark F. Adams 
1087389730f3SMark F. Adams .seealso: ()
1088389730f3SMark F. Adams  @*/
1089389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
1090389730f3SMark F. Adams {
1091389730f3SMark F. Adams   PetscErrorCode ierr;
1092389730f3SMark F. Adams 
1093389730f3SMark F. Adams   PetscFunctionBegin;
1094389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1095389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1096389730f3SMark F. Adams   PetscFunctionReturn(0);
1097389730f3SMark F. Adams }
1098389730f3SMark F. Adams 
1099389730f3SMark F. Adams EXTERN_C_BEGIN
1100389730f3SMark F. Adams #undef __FUNCT__
1101389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
1102389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1103389730f3SMark F. Adams {
1104389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1105389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1106389730f3SMark F. Adams 
1107389730f3SMark F. Adams   PetscFunctionBegin;
11089d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
1109389730f3SMark F. Adams   PetscFunctionReturn(0);
1110389730f3SMark F. Adams }
1111389730f3SMark F. Adams EXTERN_C_END
1112389730f3SMark F. Adams 
1113389730f3SMark F. Adams #undef __FUNCT__
11148263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1115676e1743SMark F. Adams /*@
11168263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1117676e1743SMark F. Adams 
1118676e1743SMark F. Adams    Collective on PC
1119676e1743SMark F. Adams 
1120676e1743SMark F. Adams    Input Parameters:
1121676e1743SMark F. Adams .  pc - the preconditioner context
1122676e1743SMark F. Adams 
1123676e1743SMark F. Adams 
1124676e1743SMark F. Adams    Options Database Key:
11258263b398SMark F. Adams .  -pc_gamg_repartition
1126676e1743SMark F. Adams 
1127676e1743SMark F. Adams    Level: intermediate
1128676e1743SMark F. Adams 
1129676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1130676e1743SMark F. Adams 
1131676e1743SMark F. Adams .seealso: ()
1132676e1743SMark F. Adams @*/
11338263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1134676e1743SMark F. Adams {
1135676e1743SMark F. Adams   PetscErrorCode ierr;
1136676e1743SMark F. Adams 
1137676e1743SMark F. Adams   PetscFunctionBegin;
1138676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11398263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1140676e1743SMark F. Adams   PetscFunctionReturn(0);
1141676e1743SMark F. Adams }
1142676e1743SMark F. Adams 
1143676e1743SMark F. Adams EXTERN_C_BEGIN
1144676e1743SMark F. Adams #undef __FUNCT__
11458263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
11468263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1147676e1743SMark F. Adams {
1148c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1149c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1150676e1743SMark F. Adams 
1151676e1743SMark F. Adams   PetscFunctionBegin;
11529d5b6da9SMark F. Adams   pc_gamg->repart = n;
1153676e1743SMark F. Adams   PetscFunctionReturn(0);
1154676e1743SMark F. Adams }
1155676e1743SMark F. Adams EXTERN_C_END
1156676e1743SMark F. Adams 
1157676e1743SMark F. Adams #undef __FUNCT__
1158dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl"
1159dfd5c07aSMark F. Adams /*@
1160dfd5c07aSMark F. Adams    PCGAMGSetReuseProl - Reuse prlongation
1161dfd5c07aSMark F. Adams 
1162dfd5c07aSMark F. Adams    Collective on PC
1163dfd5c07aSMark F. Adams 
1164dfd5c07aSMark F. Adams    Input Parameters:
1165dfd5c07aSMark F. Adams .  pc - the preconditioner context
1166dfd5c07aSMark F. Adams 
1167dfd5c07aSMark F. Adams 
1168dfd5c07aSMark F. Adams    Options Database Key:
1169dfd5c07aSMark F. Adams .  -pc_gamg_reuse_interpolation
1170dfd5c07aSMark F. Adams 
1171dfd5c07aSMark F. Adams    Level: intermediate
1172dfd5c07aSMark F. Adams 
1173dfd5c07aSMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1174dfd5c07aSMark F. Adams 
1175dfd5c07aSMark F. Adams .seealso: ()
1176dfd5c07aSMark F. Adams @*/
1177dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n)
1178dfd5c07aSMark F. Adams {
1179dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1180dfd5c07aSMark F. Adams 
1181dfd5c07aSMark F. Adams   PetscFunctionBegin;
1182dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1183dfd5c07aSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1184dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1185dfd5c07aSMark F. Adams }
1186dfd5c07aSMark F. Adams 
1187dfd5c07aSMark F. Adams EXTERN_C_BEGIN
1188dfd5c07aSMark F. Adams #undef __FUNCT__
1189dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl_GAMG"
1190dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n)
1191dfd5c07aSMark F. Adams {
1192dfd5c07aSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1193dfd5c07aSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1194dfd5c07aSMark F. Adams 
1195dfd5c07aSMark F. Adams   PetscFunctionBegin;
1196dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1197dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1198dfd5c07aSMark F. Adams }
1199dfd5c07aSMark F. Adams EXTERN_C_END
1200dfd5c07aSMark F. Adams 
1201dfd5c07aSMark F. Adams #undef __FUNCT__
1202ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1203ffc955d6SMark F. Adams /*@
1204ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1205ffc955d6SMark F. Adams 
1206ffc955d6SMark F. Adams    Collective on PC
1207ffc955d6SMark F. Adams 
1208ffc955d6SMark F. Adams    Input Parameters:
1209ffc955d6SMark F. Adams .  pc - the preconditioner context
1210ffc955d6SMark F. Adams 
1211ffc955d6SMark F. Adams 
1212ffc955d6SMark F. Adams    Options Database Key:
1213ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1214ffc955d6SMark F. Adams 
1215ffc955d6SMark F. Adams    Level: intermediate
1216ffc955d6SMark F. Adams 
1217ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1218ffc955d6SMark F. Adams 
1219ffc955d6SMark F. Adams .seealso: ()
1220ffc955d6SMark F. Adams @*/
1221ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1222ffc955d6SMark F. Adams {
1223ffc955d6SMark F. Adams   PetscErrorCode ierr;
1224ffc955d6SMark F. Adams 
1225ffc955d6SMark F. Adams   PetscFunctionBegin;
1226ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1227ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1228ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1229ffc955d6SMark F. Adams }
1230ffc955d6SMark F. Adams 
1231ffc955d6SMark F. Adams EXTERN_C_BEGIN
1232ffc955d6SMark F. Adams #undef __FUNCT__
1233ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
1234ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1235ffc955d6SMark F. Adams {
1236ffc955d6SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1237ffc955d6SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1238ffc955d6SMark F. Adams 
1239ffc955d6SMark F. Adams   PetscFunctionBegin;
1240ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1241ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1242ffc955d6SMark F. Adams }
1243ffc955d6SMark F. Adams EXTERN_C_END
1244ffc955d6SMark F. Adams 
1245ffc955d6SMark F. Adams #undef __FUNCT__
12464ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
12474ef23d27SMark F. Adams /*@
12484ef23d27SMark F. Adams    PCGAMGSetNlevels -
12494ef23d27SMark F. Adams 
12504ef23d27SMark F. Adams    Not collective on PC
12514ef23d27SMark F. Adams 
12524ef23d27SMark F. Adams    Input Parameters:
12534ef23d27SMark F. Adams .  pc - the preconditioner context
12544ef23d27SMark F. Adams 
12554ef23d27SMark F. Adams 
12564ef23d27SMark F. Adams    Options Database Key:
12574ef23d27SMark F. Adams .  -pc_mg_levels
12584ef23d27SMark F. Adams 
12594ef23d27SMark F. Adams    Level: intermediate
12604ef23d27SMark F. Adams 
12614ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12624ef23d27SMark F. Adams 
12634ef23d27SMark F. Adams .seealso: ()
12644ef23d27SMark F. Adams @*/
12654ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12664ef23d27SMark F. Adams {
12674ef23d27SMark F. Adams   PetscErrorCode ierr;
12684ef23d27SMark F. Adams 
12694ef23d27SMark F. Adams   PetscFunctionBegin;
12704ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12714ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
12724ef23d27SMark F. Adams   PetscFunctionReturn(0);
12734ef23d27SMark F. Adams }
12744ef23d27SMark F. Adams 
12754ef23d27SMark F. Adams EXTERN_C_BEGIN
12764ef23d27SMark F. Adams #undef __FUNCT__
12774ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
12784ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12794ef23d27SMark F. Adams {
12804ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
12814ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
12824ef23d27SMark F. Adams 
12834ef23d27SMark F. Adams   PetscFunctionBegin;
12849d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12854ef23d27SMark F. Adams   PetscFunctionReturn(0);
12864ef23d27SMark F. Adams }
12874ef23d27SMark F. Adams EXTERN_C_END
12884ef23d27SMark F. Adams 
12894ef23d27SMark F. Adams #undef __FUNCT__
12903542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
12913542efc5SMark F. Adams /*@
12923542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12933542efc5SMark F. Adams 
12943542efc5SMark F. Adams    Not collective on PC
12953542efc5SMark F. Adams 
12963542efc5SMark F. Adams    Input Parameters:
12973542efc5SMark F. Adams .  pc - the preconditioner context
12983542efc5SMark F. Adams 
12993542efc5SMark F. Adams 
13003542efc5SMark F. Adams    Options Database Key:
13013542efc5SMark F. Adams .  -pc_gamg_threshold
13023542efc5SMark F. Adams 
13033542efc5SMark F. Adams    Level: intermediate
13043542efc5SMark F. Adams 
13053542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
13063542efc5SMark F. Adams 
13073542efc5SMark F. Adams .seealso: ()
13083542efc5SMark F. Adams @*/
13093542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
13103542efc5SMark F. Adams {
13113542efc5SMark F. Adams   PetscErrorCode ierr;
13123542efc5SMark F. Adams 
13133542efc5SMark F. Adams   PetscFunctionBegin;
13143542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13153542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
13163542efc5SMark F. Adams   PetscFunctionReturn(0);
13173542efc5SMark F. Adams }
13183542efc5SMark F. Adams 
13193542efc5SMark F. Adams EXTERN_C_BEGIN
13203542efc5SMark F. Adams #undef __FUNCT__
13213542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
13223542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
13233542efc5SMark F. Adams {
1324c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1325c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
13263542efc5SMark F. Adams 
13273542efc5SMark F. Adams   PetscFunctionBegin;
13289d5b6da9SMark F. Adams   pc_gamg->threshold = n;
13293542efc5SMark F. Adams   PetscFunctionReturn(0);
13303542efc5SMark F. Adams }
13313542efc5SMark F. Adams EXTERN_C_END
13323542efc5SMark F. Adams 
13333542efc5SMark F. Adams #undef __FUNCT__
13349d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1335676e1743SMark F. Adams /*@
13369d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1337676e1743SMark F. Adams 
1338676e1743SMark F. Adams    Collective on PC
1339676e1743SMark F. Adams 
1340676e1743SMark F. Adams    Input Parameters:
1341676e1743SMark F. Adams .  pc - the preconditioner context
1342676e1743SMark F. Adams 
1343676e1743SMark F. Adams 
1344676e1743SMark F. Adams    Options Database Key:
13453542efc5SMark F. Adams .  -pc_gamg_type
1346676e1743SMark F. Adams 
1347676e1743SMark F. Adams    Level: intermediate
1348676e1743SMark F. Adams 
1349676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1350676e1743SMark F. Adams 
1351676e1743SMark F. Adams .seealso: ()
1352676e1743SMark F. Adams @*/
135319fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1354676e1743SMark F. Adams {
1355676e1743SMark F. Adams   PetscErrorCode ierr;
1356676e1743SMark F. Adams 
1357676e1743SMark F. Adams   PetscFunctionBegin;
1358676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1359806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1360676e1743SMark F. Adams   PetscFunctionReturn(0);
1361676e1743SMark F. Adams }
1362676e1743SMark F. Adams 
1363676e1743SMark F. Adams EXTERN_C_BEGIN
1364676e1743SMark F. Adams #undef __FUNCT__
13659d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
136619fd82e9SBarry Smith PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1367676e1743SMark F. Adams {
13689d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
1369676e1743SMark F. Adams 
1370676e1743SMark F. Adams   PetscFunctionBegin;
1371140e18c1SBarry Smith   ierr = PetscFunctionListFind(((PetscObject)pc)->comm,GAMGList,type,PETSC_FALSE,(PetscVoidStarFunction)&r);CHKERRQ(ierr);
13729d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
13739d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1374676e1743SMark F. Adams   PetscFunctionReturn(0);
1375676e1743SMark F. Adams }
1376676e1743SMark F. Adams EXTERN_C_END
1377676e1743SMark F. Adams 
13785b89ad90SMark F. Adams #undef __FUNCT__
13795b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
13805b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
13815b89ad90SMark F. Adams {
1382676e1743SMark F. Adams   PetscErrorCode  ierr;
1383676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1384676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1385676e1743SMark F. Adams   PetscBool        flag;
13865e7c91beSJed Brown   PetscInt         two = 2;
1387b43b03e9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
13885b89ad90SMark F. Adams 
13895b89ad90SMark F. Adams   PetscFunctionBegin;
1390676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
1391676e1743SMark F. Adams   {
139275b74bdaSMark F. Adams     /* -pc_gamg_verbose */
13939d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
13949d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
1395806fa848SBarry Smith                            &pc_gamg->verbose, PETSC_NULL);CHKERRQ(ierr);
13968263b398SMark F. Adams     /* -pc_gamg_repartition */
13978263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
13988263b398SMark F. Adams                             "Repartion coarse grids (false)",
13998263b398SMark F. Adams                             "PCGAMGRepartitioning",
14009d5b6da9SMark F. Adams                             pc_gamg->repart,
14019d5b6da9SMark F. Adams                             &pc_gamg->repart,
1402806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1403dfd5c07aSMark F. Adams     /* -pc_gamg_reuse_interpolation */
1404dfd5c07aSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation",
1405dfd5c07aSMark F. Adams                             "Reuse prolongation operator (true)",
1406dfd5c07aSMark F. Adams                             "PCGAMGReuseProl",
1407dfd5c07aSMark F. Adams                             pc_gamg->reuse_prol,
1408dfd5c07aSMark F. Adams                             &pc_gamg->reuse_prol,
1409dfd5c07aSMark F. Adams                             &flag);CHKERRQ(ierr);
1410ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
1411ffc955d6SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1412ffc955d6SMark F. Adams                             "Use aggregation agragates for GASM smoother (false)",
1413ffc955d6SMark F. Adams                             "PCGAMGUseASMAggs",
1414ffc955d6SMark F. Adams                             pc_gamg->use_aggs_in_gasm,
1415ffc955d6SMark F. Adams                             &pc_gamg->use_aggs_in_gasm,
1416806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1417c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1418676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1419676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1420676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
14219d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
14229d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1423806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1424389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1425389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1426389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1427389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
14289d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
14299d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1430806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1431c20e4228SMark F. Adams     /* -pc_gamg_threshold */
14323542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
14333542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
14343542efc5SMark F. Adams                             "PCGAMGSetThreshold",
14359d5b6da9SMark F. Adams                             pc_gamg->threshold,
14369d5b6da9SMark F. Adams                             &pc_gamg->threshold,
1437806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1438806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
1439806fa848SBarry Smith       ierr = PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1440806fa848SBarry Smith     }
14414ef23d27SMark F. Adams 
14425e7c91beSJed Brown     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,PETSC_NULL);CHKERRQ(ierr);
14434ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
14444ef23d27SMark F. Adams                            "Set number of MG levels",
14454ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
14469d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
14479d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
1448806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1449676e1743SMark F. Adams   }
1450676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1451676e1743SMark F. Adams 
14525b89ad90SMark F. Adams   PetscFunctionReturn(0);
14535b89ad90SMark F. Adams }
14545b89ad90SMark F. Adams 
14555b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
14565b89ad90SMark F. Adams /*MC
1457856830bfSJed Brown      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
14589d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
14595b89ad90SMark F. Adams 
1460280d9858SJed Brown    Options Database Keys:
14615b89ad90SMark F. Adams    Multigrid options(inherited)
1462280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1463280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1464280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
14658c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
14665b89ad90SMark F. Adams 
14675b89ad90SMark F. Adams   Level: intermediate
1468280d9858SJed Brown 
14695b89ad90SMark F. Adams   Concepts: multigrid
14705b89ad90SMark F. Adams 
14715b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1472280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
14735b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
14745b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
14755b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
14765b89ad90SMark F. Adams M*/
14775b89ad90SMark F. Adams EXTERN_C_BEGIN
14785b89ad90SMark F. Adams #undef __FUNCT__
14795b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
14805b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
14815b89ad90SMark F. Adams {
14825b89ad90SMark F. Adams   PetscErrorCode  ierr;
14835b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
14845b89ad90SMark F. Adams   PC_MG           *mg;
14850cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14869d5b6da9SMark F. Adams   static long count = 0;
14875ee9c036SSatish Balay #endif
14885b89ad90SMark F. Adams 
14895b89ad90SMark F. Adams   PetscFunctionBegin;
149060a6d8f8SMark F. Adams 
14915b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14925b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
14935b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14945b89ad90SMark F. Adams 
14955b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
14965b89ad90SMark F. Adams   ierr = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr);
14975b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
1498ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
14995b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
15005b89ad90SMark F. Adams 
15019d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
15029d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
15039d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
15049d5b6da9SMark F. Adams   pc_gamg->data = 0;
15055b89ad90SMark F. Adams 
15069d5b6da9SMark F. Adams   /* register AMG type */
15079d5b6da9SMark F. Adams   if (!GAMGList) {
1508140e18c1SBarry Smith     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
1509140e18c1SBarry Smith     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
15109d5b6da9SMark F. Adams   }
15119d5b6da9SMark F. Adams 
15129d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
15135b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
15145b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
15155b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
15165b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
15175b89ad90SMark F. Adams 
15185b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1519676e1743SMark F. Adams                                             "PCGAMGSetProcEqLim_C",
1520676e1743SMark F. Adams                                             "PCGAMGSetProcEqLim_GAMG",
1521806fa848SBarry Smith                                             PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1522676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1523389730f3SMark F. Adams                                             "PCGAMGSetCoarseEqLim_C",
1524389730f3SMark F. Adams                                             "PCGAMGSetCoarseEqLim_GAMG",
1525806fa848SBarry Smith                                             PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1526389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15278263b398SMark F. Adams                                             "PCGAMGSetRepartitioning_C",
15288263b398SMark F. Adams                                             "PCGAMGSetRepartitioning_GAMG",
1529806fa848SBarry Smith                                             PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
1530676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1531dfd5c07aSMark F. Adams                                             "PCGAMGSetReuseProl_C",
1532dfd5c07aSMark F. Adams                                             "PCGAMGSetReuseProl_GAMG",
1533dfd5c07aSMark F. Adams                                             PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr);
1534dfd5c07aSMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1535ffc955d6SMark F. Adams                                             "PCGAMGSetUseASMAggs_C",
1536ffc955d6SMark F. Adams                                             "PCGAMGSetUseASMAggs_GAMG",
1537806fa848SBarry Smith                                             PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1538ffc955d6SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15393542efc5SMark F. Adams                                             "PCGAMGSetThreshold_C",
15403542efc5SMark F. Adams                                             "PCGAMGSetThreshold_GAMG",
1541806fa848SBarry Smith                                             PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
15423542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15439d5b6da9SMark F. Adams                                             "PCGAMGSetType_C",
15449d5b6da9SMark F. Adams                                             "PCGAMGSetType_GAMG",
1545806fa848SBarry Smith                                             PCGAMGSetType_GAMG);CHKERRQ(ierr);
15469d5b6da9SMark F. Adams   pc_gamg->repart = PETSC_FALSE;
1547dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = PETSC_TRUE;
1548ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
15499d5b6da9SMark F. Adams   pc_gamg->min_eq_proc = 100;
1550c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit = 800;
1551a2f3521dSMark F. Adams   pc_gamg->threshold = 0.001;
15529d5b6da9SMark F. Adams   pc_gamg->Nlevels = GAMG_MAXLEVELS;
15539d5b6da9SMark F. Adams   pc_gamg->verbose = 0;
15549d5b6da9SMark F. Adams   pc_gamg->emax_id = -1;
15555e7c91beSJed Brown   pc_gamg->eigtarget[0] = 0.05;
15565e7c91beSJed Brown   pc_gamg->eigtarget[1] = 1.05;
15579d5b6da9SMark F. Adams 
15580cbbd2e1SMark F. Adams   /* private events */
15590cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1560785cba28SMark F. Adams   if (count++ == 0) {
1561806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1562806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
15630cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
15640cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
15650cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1566806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1567806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1568806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1569806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1570806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1571806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1572806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1573806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1574806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1575806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1576806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1577806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1578f852f58cSMark F. Adams 
15790cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
15800cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
15810cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1582b4fbaa2aSMark F. Adams     /* create timer stages */
1583b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1584b4fbaa2aSMark F. Adams     {
1585b4fbaa2aSMark F. Adams       char     str[32];
1586b4fbaa2aSMark F. Adams       PetscInt lidx;
1587806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1588806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1589b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++) {
1590b4fbaa2aSMark F. Adams         sprintf(str,"MG Level %d",lidx);
1591806fa848SBarry Smith         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1592b4fbaa2aSMark F. Adams       }
1593b4fbaa2aSMark F. Adams     }
1594b4fbaa2aSMark F. Adams #endif
1595b4fbaa2aSMark F. Adams   }
1596b4fbaa2aSMark F. Adams #endif
15970cbbd2e1SMark F. Adams   /* general events */
15980cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
1599806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1600806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1601806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1602806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1603806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1604806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1605806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1606806fa848SBarry Smith   ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr);
16070cbbd2e1SMark F. Adams #endif
16080cbbd2e1SMark F. Adams 
160960a6d8f8SMark F. Adams   /* instantiate derived type */
161060a6d8f8SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
161160a6d8f8SMark F. Adams   {
161260a6d8f8SMark F. Adams     char tname[256] = GAMGAGG;
1613806fa848SBarry Smith     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), PETSC_NULL);CHKERRQ(ierr);
161460a6d8f8SMark F. Adams     ierr = PCGAMGSetType(pc, tname);CHKERRQ(ierr);
161560a6d8f8SMark F. Adams   }
161660a6d8f8SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
16175b89ad90SMark F. Adams   PetscFunctionReturn(0);
16185b89ad90SMark F. Adams }
16195b89ad90SMark F. Adams EXTERN_C_END
1620