xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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 */
43*ce94432eSBarry Smith     PetscPrintf(PetscObjectComm((PetscObject)pc),"***[%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   }
460298fd71SBarry Smith   pc_gamg->data = NULL; pc_gamg->data_sz = 0;
47878e152fSMark F. Adams 
48878e152fSMark F. Adams   if (pc_gamg->orig_data) {
49878e152fSMark F. Adams     ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
50878e152fSMark F. Adams   }
51a2f3521dSMark F. Adams   PetscFunctionReturn(0);
52a2f3521dSMark F. Adams }
53a2f3521dSMark F. Adams 
54a2f3521dSMark F. Adams /* private 2x2 Mat Nest for Stokes */
55a2f3521dSMark F. Adams typedef struct {
56a2f3521dSMark F. Adams   Mat A11,A21,A12,Amat;
57a2f3521dSMark F. Adams   IS  prim_is,constr_is;
58a2f3521dSMark F. Adams } GAMGKKTMat;
59a2f3521dSMark F. Adams 
60a2f3521dSMark F. Adams #undef __FUNCT__
61a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatCreate"
62a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatCreate(Mat A, PetscBool iskkt, GAMGKKTMat *out)
63a2f3521dSMark F. Adams {
64a2f3521dSMark F. Adams   PetscFunctionBegin;
65a2f3521dSMark F. Adams   out->Amat = A;
66a2f3521dSMark F. Adams   if (iskkt) {
67a2f3521dSMark F. Adams     PetscErrorCode ierr;
68a2f3521dSMark F. Adams     IS             is_constraint, is_prime;
69a2f3521dSMark F. Adams     PetscInt       nmin,nmax;
70a2f3521dSMark F. Adams 
71a2f3521dSMark F. Adams     ierr = MatGetOwnershipRange(A, &nmin, &nmax);CHKERRQ(ierr);
72a2f3521dSMark F. Adams     ierr = MatFindZeroDiagonals(A, &is_constraint);CHKERRQ(ierr);
73a2f3521dSMark F. Adams     ierr = ISComplement(is_constraint, nmin, nmax, &is_prime);CHKERRQ(ierr);
742fa5cd67SKarl Rupp 
75a2f3521dSMark F. Adams     out->prim_is   = is_prime;
76a2f3521dSMark F. Adams     out->constr_is = is_constraint;
77a2f3521dSMark F. Adams 
78a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_prime, is_prime,      MAT_INITIAL_MATRIX, &out->A11);CHKERRQ(ierr);
79a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12);CHKERRQ(ierr);
80a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21);CHKERRQ(ierr);
81806fa848SBarry Smith   } else {
82a2f3521dSMark F. Adams     out->A11       = A;
830298fd71SBarry Smith     out->A21       = NULL;
840298fd71SBarry Smith     out->A12       = NULL;
850298fd71SBarry Smith     out->prim_is   = NULL;
860298fd71SBarry Smith     out->constr_is = NULL;
87a2f3521dSMark F. Adams   }
88a2f3521dSMark F. Adams   PetscFunctionReturn(0);
89a2f3521dSMark F. Adams }
90a2f3521dSMark F. Adams 
91a2f3521dSMark F. Adams #undef __FUNCT__
92a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatDestroy"
93a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatDestroy(GAMGKKTMat *mat)
94a2f3521dSMark F. Adams {
95a2f3521dSMark F. Adams   PetscErrorCode ierr;
96a2f3521dSMark F. Adams 
97a2f3521dSMark F. Adams   PetscFunctionBegin;
98a2f3521dSMark F. Adams   if (mat->A11 && mat->A11 != mat->Amat) {
99a2f3521dSMark F. Adams     ierr = MatDestroy(&mat->A11);CHKERRQ(ierr);
100a2f3521dSMark F. Adams   }
101a2f3521dSMark F. Adams   ierr = MatDestroy(&mat->A21);CHKERRQ(ierr);
102a2f3521dSMark F. Adams   ierr = MatDestroy(&mat->A12);CHKERRQ(ierr);
103a2f3521dSMark F. Adams 
104a2f3521dSMark F. Adams   ierr = ISDestroy(&mat->prim_is);CHKERRQ(ierr);
105a2f3521dSMark F. Adams   ierr = ISDestroy(&mat->constr_is);CHKERRQ(ierr);
106d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
107d3d6bff4SMark F. Adams }
108d3d6bff4SMark F. Adams 
1095b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1105b89ad90SMark F. Adams /*
111a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
112a147abb0SMark F. Adams      of active processors.
1135b89ad90SMark F. Adams 
1145b89ad90SMark F. Adams    Input Parameter:
115a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
116a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
1179d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
118c5bfad50SMark F. Adams    . cr_bs - coarse block size
1199d5b6da9SMark F. Adams    . isLast -
120a2f3521dSMark F. Adams    . stokes -
1213530afc2SMark F. Adams    In/Output Parameter:
122a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
123afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
12411e60469SMark F. Adams    Output Parameter:
1253530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1265b89ad90SMark F. Adams */
1275cb416c2SMark F. Adams 
1285b89ad90SMark F. Adams #undef __FUNCT__
1298263b398SMark F. Adams #define __FUNCT__ "createLevel"
1302fa5cd67SKarl Rupp static PetscErrorCode createLevel(const PC pc,const Mat Amat_fine,const PetscInt cr_bs,const PetscBool isLast,const PetscBool stokes,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc)
1315b89ad90SMark F. Adams {
132a2f3521dSMark F. Adams   PetscErrorCode  ierr;
1339d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
134486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
1359d5b6da9SMark F. Adams   const PetscBool repart      = pc_gamg->repart;
1369d5b6da9SMark F. Adams   const PetscInt  min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
137a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
1389d5b6da9SMark F. Adams   MPI_Comm        wcomm = ((PetscObject)Amat_fine)->comm;
139c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
140c5bfad50SMark F. Adams   PetscInt        ncrs_eq,ncrs_prim,f_bs;
1415b89ad90SMark F. Adams 
1425b89ad90SMark F. Adams   PetscFunctionBegin;
143c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
144c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
145c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
14611e60469SMark F. Adams   /* RAP */
1479d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
148038e3b61SMark F. Adams 
149a2f3521dSMark F. Adams   /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/
150a2f3521dSMark F. Adams   ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
15171959b99SBarry Smith   if (pc_gamg->data_sz % (pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data_sz %D not divisible by (pc_gamg->data_cell_cols %D *pc_gamg->data_cell_rows %D)",pc_gamg->data_sz,pc_gamg->data_cell_cols,pc_gamg->data_cell_rows);
1520298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
153a2f3521dSMark F. Adams 
154c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
155a2f3521dSMark F. Adams   {
156472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
1570298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
158c5df96a5SBarry Smith     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
159c5df96a5SBarry Smith     if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1;
160c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
161c5df96a5SBarry Smith     if (isLast) new_size = 1;
162a2f3521dSMark F. Adams   }
163f852f58cSMark F. Adams 
1642fa5cd67SKarl Rupp   if (!repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1652fa5cd67SKarl Rupp   else {
166a2f3521dSMark F. Adams     const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
167a2f3521dSMark F. Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
168a2f3521dSMark F. Adams     IS             is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
1695a9b9e01SMark F. Adams     VecScatter     vecscat;
17022063be5SMark F. Adams     PetscScalar    *array;
17122063be5SMark F. Adams     Vec            src_crd, dest_crd;
172e33ef3b1SMark F. Adams 
17371959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
17471959b99SBarry Smith     if (ncrs_eq % cr_bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs);
1750cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1760cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
177b4fbaa2aSMark F. Adams #endif
178a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
179c5df96a5SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt), &counts);CHKERRQ(ierr);
180a2f3521dSMark F. Adams     if (repart && !stokes) {
181a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1825a9b9e01SMark F. Adams       Mat adj;
1835a9b9e01SMark F. Adams 
184a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
185c5df96a5SBarry Smith         if (pc_gamg->verbose==1) PetscPrintf(wcomm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,ncrs_eq);
186a2f3521dSMark F. Adams         else {
187a2f3521dSMark F. Adams           PetscInt n;
188a2f3521dSMark F. Adams           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm);CHKERRQ(ierr);
189c5df96a5SBarry Smith           PetscPrintf(wcomm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n);
190a2f3521dSMark F. Adams         }
191a2f3521dSMark F. Adams       }
1925a9b9e01SMark F. Adams 
193a2f3521dSMark F. Adams       /* get 'adj' */
194c5bfad50SMark F. Adams       if (cr_bs == 1) {
195038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
196806fa848SBarry Smith       } else {
197a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
198eb07cef2SMark F. Adams         Mat               tMat;
199a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
200b4fbaa2aSMark F. Adams         const PetscScalar *vals;
201b4fbaa2aSMark F. Adams         const PetscInt    *idx;
202a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
2039057884aSMark F. Adams         static PetscInt   llev = 0;
204b4fbaa2aSMark F. Adams 
205a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &d_nnz);CHKERRQ(ierr);
206a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &o_nnz);CHKERRQ(ierr);
207a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
208a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
209c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
21058471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
211c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
212c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
21358471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
214a2f3521dSMark F. Adams           if (d_nnz[jj] > ncrs_prim) d_nnz[jj] = ncrs_prim;
215c5bfad50SMark F. Adams           if (o_nnz[jj] > (M/cr_bs-ncrs_prim)) o_nnz[jj] = M/cr_bs-ncrs_prim;
21658471d46SMark F. Adams         }
2176876a03eSMark F. Adams 
218a2f3521dSMark F. Adams         ierr = MatCreate(wcomm, &tMat);CHKERRQ(ierr);
219806fa848SBarry Smith         ierr = MatSetSizes(tMat, ncrs_prim, ncrs_prim,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
220a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
221a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
222a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
22358471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
2245f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
225eb07cef2SMark F. Adams 
226a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
227c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
22822063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
229eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
230c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
231eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
232eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
233eb07cef2SMark F. Adams           }
23422063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
235eb07cef2SMark F. Adams         }
236eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
237eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
238eb07cef2SMark F. Adams 
239b4fbaa2aSMark F. Adams         if (llev++ == -1) {
240b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2418caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
242b4fbaa2aSMark F. Adams           PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
243b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2443bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
245b4fbaa2aSMark F. Adams         }
246b4fbaa2aSMark F. Adams 
247eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
248eb07cef2SMark F. Adams 
249eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
250a2f3521dSMark F. Adams       } /* create 'adj' */
251f150b916SMark F. Adams 
252a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2535a9b9e01SMark F. Adams         char            prefix[256];
2545a9b9e01SMark F. Adams         const char      *pcpre;
255b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
256b4fbaa2aSMark F. Adams         MatPartitioning mpart;
257a4b7d37bSMark F. Adams         IS              proc_is;
258a2f3521dSMark F. Adams         PetscInt        targetPE;
2592f03bc48SMark F. Adams 
2605a9b9e01SMark F. Adams         ierr = MatPartitioningCreate(wcomm, &mpart);CHKERRQ(ierr);
2615ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2629d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2638caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
26459a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
26511e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
266c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
267a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
26811e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2695a9b9e01SMark F. Adams 
2705ef31b24SMark F. Adams         /* collect IS info */
271a2f3521dSMark F. Adams         ierr     = PetscMalloc(ncrs_eq*sizeof(PetscInt), &newproc_idx);CHKERRQ(ierr);
272a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
273a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
274c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
275a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
276c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
277a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
278eb07cef2SMark F. Adams           }
2795ef31b24SMark F. Adams         }
280a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
281a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2825ef31b24SMark F. Adams       }
2835ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2845a9b9e01SMark F. Adams 
285806fa848SBarry Smith       ierr = ISCreateGeneral(wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2868263b398SMark F. Adams       if (newproc_idx != 0) {
2878263b398SMark F. Adams         ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2885ef31b24SMark F. Adams       }
289806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
290a2f3521dSMark F. Adams 
291a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
2925a9b9e01SMark F. Adams       /* find factor */
293c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2945a9b9e01SMark F. Adams       else {
2955a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2965a9b9e01SMark F. Adams         jj = -1;
297c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
298c5df96a5SBarry Smith           if (size%kk==0) { /* a candidate */
299c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
3005a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
3015a9b9e01SMark F. Adams             if (fact > best_fact) {
3025a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
3035a9b9e01SMark F. Adams             }
3045a9b9e01SMark F. Adams           }
3055a9b9e01SMark F. Adams         }
3065a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
307a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
3085a9b9e01SMark F. Adams       }
309c5df96a5SBarry Smith       new_size = size/rfactor;
3105a9b9e01SMark F. Adams 
311c5df96a5SBarry Smith       if (new_size==nactive) {
312a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
3135a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
314a2f3521dSMark F. Adams         if (pc_gamg->verbose>0) {
315c5df96a5SBarry Smith           PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq);
316a2f3521dSMark F. Adams         }
3175a9b9e01SMark F. Adams         PetscFunctionReturn(0);
3185a9b9e01SMark F. Adams       }
3195a9b9e01SMark F. Adams 
320c5df96a5SBarry Smith       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq);
321c5df96a5SBarry Smith       targetPE = rank/rfactor;
322a2f3521dSMark F. Adams       ierr     = ISCreateStride(wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
3235a9b9e01SMark F. Adams 
324a2f3521dSMark F. Adams       if (stokes) {
325c5bfad50SMark F. Adams         ierr = ISCreateStride(wcomm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim);CHKERRQ(ierr);
3265a9b9e01SMark F. Adams       }
327a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
328e33ef3b1SMark F. Adams 
32911e60469SMark F. Adams     /*
330a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
33111e60469SMark F. Adams      */
332a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
333a2f3521dSMark F. Adams     if (stokes) {
334a2f3521dSMark F. Adams       ierr = ISPartitioningToNumbering(is_eq_newproc_prim, &is_eq_num_prim);CHKERRQ(ierr);
335806fa848SBarry Smith     } else is_eq_num_prim = is_eq_num;
33611e60469SMark F. Adams     /*
337a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
33811e60469SMark F. Adams      */
339c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
340c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
341a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
342a2f3521dSMark F. Adams     if (stokes) {
343c5df96a5SBarry Smith       ierr          = ISPartitioningCount(is_eq_newproc_prim, size, counts);CHKERRQ(ierr);
344a2f3521dSMark F. Adams       ierr          = ISDestroy(&is_eq_newproc_prim);CHKERRQ(ierr);
345c5df96a5SBarry Smith       ncrs_prim_new = counts[rank]/cr_bs; /* nodes */
346806fa848SBarry Smith     } else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */
347a2f3521dSMark F. Adams 
348a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
3490cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3500cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
351b4fbaa2aSMark F. Adams #endif
352a2f3521dSMark F. Adams 
353a2f3521dSMark F. Adams     /* move data (for primal equations only) */
35422063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
3553bf036e2SBarry Smith     ierr = VecCreate(wcomm, &dest_crd);CHKERRQ(ierr);
356a2f3521dSMark F. Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE);CHKERRQ(ierr);
357470e26f8SMark F. Adams     ierr = VecSetFromOptions(dest_crd);CHKERRQ(ierr); /* this is needed! */
35811e60469SMark F. Adams     /*
3599d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
360c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
36111e60469SMark F. Adams      */
362a2f3521dSMark F. Adams     ierr = PetscMalloc((ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx);CHKERRQ(ierr);
363a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
364a2f3521dSMark F. Adams     for (ii=0,jj=0; ii<ncrs_prim; ii++) {
365c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
366a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
36711e60469SMark F. Adams     }
368a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
369806fa848SBarry Smith     ierr = ISCreateGeneral(wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
37092a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
37111e60469SMark F. Adams     /*
37211e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
37311e60469SMark F. Adams      */
374a2f3521dSMark F. Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd);CHKERRQ(ierr);
3759d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
376a2f3521dSMark F. Adams       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
377a2f3521dSMark F. Adams       for (ii=0; ii<ncrs_prim; ii++) {
3789d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
379a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
380c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
381676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
382d3d6bff4SMark F. Adams         }
383038e3b61SMark F. Adams       }
384eb07cef2SMark F. Adams     }
385eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
386eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
38711e60469SMark F. Adams     /*
38811e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
38911e60469SMark F. Adams       to the correct processor
39011e60469SMark F. Adams     */
3910298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
39211e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
39311e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39411e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39511e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
39611e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
39711e60469SMark F. Adams     /*
39811e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
39911e60469SMark F. Adams     */
400c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
401a2f3521dSMark F. Adams     ierr = PetscMalloc(node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
4022fa5cd67SKarl Rupp 
403a2f3521dSMark F. Adams     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
404a2f3521dSMark F. Adams     strideNew        = ncrs_prim_new*ndata_rows;
4052fa5cd67SKarl Rupp 
40611e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
4079d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
408a2f3521dSMark F. Adams       for (ii=0; ii<ncrs_prim_new; ii++) {
4099d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
410a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
411c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
412d3d6bff4SMark F. Adams         }
413038e3b61SMark F. Adams       }
414038e3b61SMark F. Adams     }
41511e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
41611e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
417a2f3521dSMark F. Adams 
418a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
4190cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4200cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
421ed3f9983SMark F. Adams #endif
422a2f3521dSMark F. Adams 
42311e60469SMark F. Adams     /*
42411e60469SMark F. Adams       Invert for MatGetSubMatrix
42511e60469SMark F. Adams     */
426a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
427a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
428c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
429a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
430a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
431a2f3521dSMark F. Adams     }
432a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
4330cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4340cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
4350cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
436ed3f9983SMark F. Adams #endif
437a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
438a2f3521dSMark F. Adams     {
439a2f3521dSMark F. Adams       Mat mat;
440806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
441a2f3521dSMark F. Adams       *a_Amat_crs = mat;
442c5bfad50SMark F. Adams 
443c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
444c5bfad50SMark F. Adams         PetscInt cbs, rbs;
445c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
446c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
447c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
448c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",rank,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr);
449c5bfad50SMark F. Adams       }
450a2f3521dSMark F. Adams     }
451038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
452a2f3521dSMark F. Adams 
4530cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4540cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
455ed3f9983SMark F. Adams #endif
45611e60469SMark F. Adams     /* prolongator */
45711e60469SMark F. Adams     {
45811e60469SMark F. Adams       IS       findices;
459a2f3521dSMark F. Adams       PetscInt Istart,Iend;
460a2f3521dSMark F. Adams       Mat      Pnew;
461a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4620cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4630cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
464ed3f9983SMark F. Adams #endif
46511e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
466c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
467806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
46811e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
469c5bfad50SMark F. Adams 
470c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
471c5bfad50SMark F. Adams         PetscInt cbs, rbs;
472c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
473c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
474c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
475c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
476c5bfad50SMark F. Adams       }
4770cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4780cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
479ed3f9983SMark F. Adams #endif
4803530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
481a2f3521dSMark F. Adams 
482a2f3521dSMark F. Adams       /* output - repartitioned */
483a2f3521dSMark F. Adams       *a_P_inout = Pnew;
484e33ef3b1SMark F. Adams     }
485a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4865b89ad90SMark F. Adams 
487c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
488a2f3521dSMark F. Adams   }
4895a9b9e01SMark F. Adams 
490a2f3521dSMark F. Adams   /* outout matrix data */
491c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
492c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
493c8b0795cSMark F. Adams     if (llev==0) {
494c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
495c8b0795cSMark F. Adams       PetscViewerASCIIOpen(wcomm,fname,&viewer);
496c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
497c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
498c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
499c8b0795cSMark F. Adams     }
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(Cmat, viewer);CHKERRQ(ierr);
504c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
505c8b0795cSMark F. Adams   }
5065b89ad90SMark F. Adams   PetscFunctionReturn(0);
5075b89ad90SMark F. Adams }
5085b89ad90SMark F. Adams 
5095b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5105b89ad90SMark F. Adams /*
5115b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5125b89ad90SMark F. Adams                     by setting data structures and options.
5135b89ad90SMark F. Adams 
5145b89ad90SMark F. Adams    Input Parameter:
5155b89ad90SMark F. Adams .  pc - the preconditioner context
5165b89ad90SMark F. Adams 
5175b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
5185b89ad90SMark F. Adams 
5195b89ad90SMark F. Adams    Notes:
5205b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
5215b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
5225b89ad90SMark F. Adams */
5235b89ad90SMark F. Adams #undef __FUNCT__
5245b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
5259d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5265b89ad90SMark F. Adams {
5275b89ad90SMark F. Adams   PetscErrorCode ierr;
5289d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
5295b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
5302adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
531a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
532*ce94432eSBarry Smith   MPI_Comm       wcomm;
533c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
534587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
535c8b0795cSMark F. Adams   PetscReal      emaxs[GAMG_MAXLEVELS];
536e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
537a2f3521dSMark F. Adams   GAMGKKTMat     kktMatsArr[GAMG_MAXLEVELS];
538a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
539737a81a9SMark F. Adams   MatInfo        info;
5405169fedaSMark F. Adams   PetscBool      stokes = PETSC_FALSE, redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);
5415ef31b24SMark F. Adams 
5425b89ad90SMark F. Adams   PetscFunctionBegin;
543*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&wcomm);CHKERRQ(ierr);
544c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm,&rank);CHKERRQ(ierr);
545c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm,&size);CHKERRQ(ierr);
546dfd5c07aSMark F. Adams 
547c5df96a5SBarry Smith   if (pc_gamg->verbose>2) PetscPrintf(wcomm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",rank,__FUNCT__,pc_gamg->setup_count,pc->setupcalled);
548dfd5c07aSMark F. Adams 
54984d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
550878e152fSMark F. Adams     if (redo_mesh_setup) {
551878e152fSMark F. Adams       /* reset everything */
552878e152fSMark F. Adams       ierr            = PCReset_MG(pc);CHKERRQ(ierr);
553878e152fSMark F. Adams       pc->setupcalled = 0;
554806fa848SBarry Smith     } else {
55584d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
55603a628feSMark F. Adams       /* just do Galerkin grids */
55758471d46SMark F. Adams       Mat          B,dA,dB;
55858471d46SMark F. Adams 
55971959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
5609d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
56158471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5620298fd71SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,NULL);CHKERRQ(ierr);
56358471d46SMark F. Adams         /* (re)set to get dirty flag */
5649d5b6da9SMark F. Adams         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
56558471d46SMark F. Adams 
5669d5b6da9SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>-1; level--) {
56703a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
56884d3f75bSMark F. Adams           if (pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
56903a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
570084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
5712fa5cd67SKarl Rupp 
57203a628feSMark F. Adams             mglevels[level]->A = B;
573806fa848SBarry Smith           } else {
5740298fd71SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B,NULL);CHKERRQ(ierr);
57558471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
57603a628feSMark F. Adams           }
57758471d46SMark F. Adams           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
57858471d46SMark F. Adams           dB   = B;
57958471d46SMark F. Adams         }
5805f8cf99dSMark F. Adams       }
581d5280255SMark F. Adams 
582d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
583d5280255SMark F. Adams 
584d5280255SMark F. Adams       /* PCSetUp_MG seems to insists on setting this to GMRES */
5851ac9965eSMark F. Adams       ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr);
58658471d46SMark F. Adams       PetscFunctionReturn(0);
587eb07cef2SMark F. Adams     }
588878e152fSMark F. Adams   }
589f6536408SMark F. Adams 
5900298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
591eb07cef2SMark F. Adams 
5929d1b15c3SMark F. Adams   ierr = GAMGKKTMatCreate(Pmat, stokes, &kktMatsArr[0]);CHKERRQ(ierr);
5939d1b15c3SMark F. Adams 
594878e152fSMark F. Adams   if (!pc_gamg->data) {
595878e152fSMark F. Adams     if (pc_gamg->orig_data) {
596878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5970298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5982fa5cd67SKarl Rupp 
599878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
600878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
601878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
6022fa5cd67SKarl Rupp 
603878e152fSMark F. Adams       ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
604878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
605806fa848SBarry Smith     } else {
606f23aa3ddSBarry Smith       if (!pc_gamg->createdefaultdata) SETERRQ(wcomm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
607f23aa3ddSBarry Smith       if (stokes) SETERRQ(wcomm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
6089d1b15c3SMark F. Adams       ierr = pc_gamg->createdefaultdata(pc, kktMatsArr[0].A11);CHKERRQ(ierr);
6099d5b6da9SMark F. Adams     }
610878e152fSMark F. Adams   }
611878e152fSMark F. Adams 
612878e152fSMark F. Adams   /* cache original data for reuse */
613878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
614878e152fSMark F. Adams     ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr);
615878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
616878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
617878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
618878e152fSMark F. Adams   }
619038e3b61SMark F. Adams 
620302f38e8SMark F. Adams   /* get basic dims */
6212fa5cd67SKarl Rupp   if (stokes) bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
6222fa5cd67SKarl Rupp   else {
623302f38e8SMark F. Adams     ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
624302f38e8SMark F. Adams   }
625a2f3521dSMark F. Adams 
626a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
627c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
62884f9421dSMark F. Adams     PetscInt NN = M;
62984f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
63084f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
6313bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
632806fa848SBarry Smith     } else {
633806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
63484f9421dSMark F. Adams     }
635b2a4f308SMark F. Adams     nnz0   = info.nz_used;
636b2a4f308SMark F. Adams     nnztot = info.nz_used;
637806fa848SBarry 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",
638c5df96a5SBarry Smith                          rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
639c5df96a5SBarry Smith                          (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr);
640c8b0795cSMark F. Adams   }
64184d3f75bSMark F. Adams 
642a2f3521dSMark F. Adams   /* Get A_i and R_i */
643c5df96a5SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
644c5df96a5SBarry Smith        level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit);  /* && (size==1 || nactivepe>1); */
6450205a208SMark F. Adams        level++) {
6465b89ad90SMark F. Adams     level1 = level + 1;
6470cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6480cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
649a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
650a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
651b4fbaa2aSMark F. Adams #endif
652a2f3521dSMark F. Adams #endif
653a2f3521dSMark F. Adams     /* deal with Stokes, get sub matrices */
6549d1b15c3SMark F. Adams     if (level > 0) {
655a2f3521dSMark F. Adams       ierr = GAMGKKTMatCreate(Aarr[level], stokes, &kktMatsArr[level]);CHKERRQ(ierr);
6569d1b15c3SMark F. Adams     }
657c8b0795cSMark F. Adams     { /* construct prolongator */
658725b86d8SJed Brown       Mat              Gmat;
6590cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
660a2f3521dSMark F. Adams       Mat              Prol11,Prol22;
661c8b0795cSMark F. Adams 
662a2f3521dSMark F. Adams       ierr = pc_gamg->graph(pc,kktMatsArr[level].A11, &Gmat);CHKERRQ(ierr);
6630cbbd2e1SMark F. Adams       ierr = pc_gamg->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
664a2f3521dSMark F. Adams       ierr = pc_gamg->prolongator(pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11);CHKERRQ(ierr);
665c8b0795cSMark F. Adams 
666a2f3521dSMark F. Adams       /* could have failed to create new level */
667a2f3521dSMark F. Adams       if (Prol11) {
6689d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6690298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
670a2f3521dSMark F. Adams 
671a2f3521dSMark F. Adams         if (stokes) {
672a2f3521dSMark F. Adams           if (!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
673a2f3521dSMark F. Adams           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
674a2f3521dSMark F. Adams           ierr = pc_gamg->formkktprol(pc, Prol11, kktMatsArr[level].A21, &Prol22);CHKERRQ(ierr);
675c8b0795cSMark F. Adams         }
676c8b0795cSMark F. Adams 
677a2f3521dSMark F. Adams         if (pc_gamg->optprol) {
678c8b0795cSMark F. Adams           /* smooth */
679a2f3521dSMark F. Adams           ierr = pc_gamg->optprol(pc, kktMatsArr[level].A11, &Prol11);CHKERRQ(ierr);
680c8b0795cSMark F. Adams         }
681c8b0795cSMark F. Adams 
682a2f3521dSMark F. Adams         if (stokes) {
683da80777bSKarl Rupp           IS  is_row[2];
684da80777bSKarl Rupp           Mat a[4];
685da80777bSKarl Rupp 
686da80777bSKarl Rupp           is_row[0] = kktMatsArr[level].prim_is; is_row[1] = kktMatsArr[level].constr_is;
6870298fd71SBarry Smith           a[0]      = Prol11; a[1] = NULL; a[2] = NULL; a[3] = Prol22;
688a2f3521dSMark F. Adams           ierr      = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1]);CHKERRQ(ierr);
6892fa5cd67SKarl Rupp         } else Parr[level1] = Prol11;
6900298fd71SBarry Smith       } else Parr[level1] = NULL;
691ffc955d6SMark F. Adams 
692ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
693806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
694ffc955d6SMark F. Adams       }
695ffc955d6SMark F. Adams 
696a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
69741b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
698a2f3521dSMark F. Adams     } /* construct prolongator scope */
6990cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7000cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
701c8b0795cSMark F. Adams #endif
7029d5b6da9SMark F. Adams     /* cache eigen estimate */
7039d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
7049d5b6da9SMark F. Adams       PetscBool flag;
705806fa848SBarry Smith       ierr = PetscObjectComposedDataGetReal((PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
7069d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
707806fa848SBarry Smith     } else emaxs[level] = -1.;
7082adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
709c8b0795cSMark F. Adams     if (!Parr[level1]) {
710806fa848SBarry Smith       if (pc_gamg->verbose) {
711c5df96a5SBarry Smith         ierr =  PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr);
712806fa848SBarry Smith       }
713c8b0795cSMark F. Adams       break;
714c8b0795cSMark F. Adams     }
7150cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7160cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
717b4fbaa2aSMark F. Adams #endif
718a2f3521dSMark F. Adams 
719a2f3521dSMark F. Adams     ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
720806fa848SBarry Smith                        stokes, &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr);
721a2f3521dSMark F. Adams 
7220cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7230cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
724b4fbaa2aSMark F. Adams #endif
725a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
726a2f3521dSMark F. Adams 
727a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0) {
7280cbbd2e1SMark F. Adams       PetscInt NN = M;
7290cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
730a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
7313bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
732806fa848SBarry Smith       } else {
733806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
7340cbbd2e1SMark F. Adams       }
735a2f3521dSMark F. Adams 
7360cbbd2e1SMark F. Adams       nnztot += info.nz_used;
737806fa848SBarry Smith       ierr    = PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
738c5df96a5SBarry Smith                             rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
739806fa848SBarry Smith                             (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
740c8b0795cSMark F. Adams     }
741a2f3521dSMark F. Adams 
742a2f3521dSMark F. Adams     /* stop if one node -- could pull back for singular problems */
743c8b0795cSMark F. Adams     if (M/pc_gamg->data_cell_cols < 2) {
74481708dccSMark F. Adams       level++;
74581708dccSMark F. Adams       break;
74681708dccSMark F. Adams     }
7470cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
748b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
749b4fbaa2aSMark F. Adams #endif
750c8b0795cSMark F. Adams   } /* levels */
751c8b0795cSMark F. Adams 
752c8b0795cSMark F. Adams   if (pc_gamg->data) {
753c8b0795cSMark F. Adams     ierr          = PetscFree(pc_gamg->data);CHKERRQ(ierr);
7540298fd71SBarry Smith     pc_gamg->data = NULL;
7555b89ad90SMark F. Adams   }
756c8b0795cSMark F. Adams 
7570cbbd2e1SMark F. Adams   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
7589d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7595b89ad90SMark F. Adams   fine_level       = level;
7600298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
7615b89ad90SMark F. Adams 
76284d3f75bSMark F. Adams   /* simple setup */
76384d3f75bSMark F. Adams   if (!PETSC_TRUE) {
76484d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
76584d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
76684d3f75bSMark F. Adams          lidx<fine_level;
76784d3f75bSMark F. Adams          lidx++, level--) {
76884d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
76984d3f75bSMark F. Adams       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77084d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
77184d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
77284d3f75bSMark F. Adams     }
77384d3f75bSMark F. Adams     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77484d3f75bSMark F. Adams 
77584d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
776806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
777d5280255SMark F. Adams     /* set default smoothers & set operators */
7789d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
779587fa25dSMark F. Adams          lidx <= fine_level;
780587fa25dSMark F. Adams          lidx++, level--) {
781ffc955d6SMark F. Adams       KSP smoother;
782ffc955d6SMark F. Adams       PC  subpc;
783a2f3521dSMark F. Adams 
7849d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
785f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
786ffc955d6SMark F. Adams 
787a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
788a2f3521dSMark F. Adams       /* set ops */
789a2f3521dSMark F. Adams       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
790a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
791a2f3521dSMark F. Adams 
792a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
793a2f3521dSMark F. Adams       if (stokes) {
794a2f3521dSMark F. Adams         KSP      *ksps;
795a2f3521dSMark F. Adams         PetscInt nn;
796a2f3521dSMark F. Adams         ierr     = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);CHKERRQ(ierr);
797a2f3521dSMark F. Adams         ierr     = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is);CHKERRQ(ierr);
798a2f3521dSMark F. Adams         ierr     = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
799a2f3521dSMark F. Adams         smoother = ksps[0];
800a2f3521dSMark F. Adams         ierr     = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
801a2f3521dSMark F. Adams         ierr     = PetscFree(ksps);CHKERRQ(ierr);
802a2f3521dSMark F. Adams       }
803a2f3521dSMark F. Adams       ierr = GAMGKKTMatDestroy(&kktMatsArr[level]);CHKERRQ(ierr);
804a2f3521dSMark F. Adams 
805a2f3521dSMark F. Adams       /* set defaults */
8066c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
807a2f3521dSMark F. Adams 
808ffc955d6SMark F. Adams       /* override defaults and command line args (!) */
809ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
8102d3561bbSSatish Balay         PetscInt sz;
8112d3561bbSSatish Balay         IS       *is;
812a2f3521dSMark F. Adams 
8132d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
8142d3561bbSSatish Balay         is   = ASMLocalIDsArr[level];
815ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
816ffc955d6SMark F. Adams         if (sz==0) {
817ffc955d6SMark F. Adams           IS       is;
818ffc955d6SMark F. Adams           PetscInt my0,kk;
819ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
820ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
8210298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr);
822a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
823806fa848SBarry Smith         } else {
824a94c3b12SMark F. Adams           PetscInt kk;
8250298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr);
826a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
827a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
828a94c3b12SMark F. Adams           }
829ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
830ffc955d6SMark F. Adams         }
831ffc955d6SMark F. Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
832ffc955d6SMark F. Adams 
8330298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
834ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
835ffc955d6SMark F. Adams         ierr                  = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
836806fa848SBarry Smith       } else {
8379d5b6da9SMark F. Adams         ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr);
838ffc955d6SMark F. Adams       }
839d5280255SMark F. Adams     }
840d5280255SMark F. Adams     {
841d5280255SMark F. Adams       /* coarse grid */
842d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
843d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
844d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
845d5280255SMark F. Adams       ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
846d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
847d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
848d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
849d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
85071959b99SBarry Smith       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
85171959b99SBarry Smith       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
852d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
853d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
8549dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
855d5280255SMark F. Adams     }
856d5280255SMark F. Adams 
857d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
858d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
859d5280255SMark F. Adams     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr);
860d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
861d5280255SMark 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.");
862d5280255SMark F. Adams 
863d5280255SMark F. Adams     /* create cheby smoothers */
864d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
865d5280255SMark F. Adams          lidx <= fine_level;
866d5280255SMark F. Adams          lidx++, level--) {
867d5280255SMark F. Adams       KSP       smoother;
868ffc955d6SMark F. Adams       PetscBool flag;
869d5280255SMark F. Adams       PC        subpc;
870d5280255SMark F. Adams 
871ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
872a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
873a2f3521dSMark F. Adams 
874a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
875a2f3521dSMark F. Adams       if (stokes) {
876a2f3521dSMark F. Adams         KSP      *ksps;
877a2f3521dSMark F. Adams         PetscInt nn;
878a2f3521dSMark F. Adams         ierr     = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
879a2f3521dSMark F. Adams         smoother = ksps[0];
880a2f3521dSMark F. Adams         ierr     = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
881a2f3521dSMark F. Adams         ierr     = PetscFree(ksps);CHKERRQ(ierr);
882a2f3521dSMark F. Adams       }
883ffc955d6SMark F. Adams 
884ffc955d6SMark F. Adams       /* do my own cheby */
8856c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
886ffc955d6SMark F. Adams       if (flag) {
887ffc955d6SMark F. Adams         PetscReal emax, emin;
888251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
889ffc955d6SMark F. Adams         if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
890587fa25dSMark F. Adams         else { /* eigen estimate 'emax' */
891e696ed0bSMark F. Adams           KSP eksp;
892e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
893b2a4f308SMark F. Adams           Vec bb, xx;
894038e3b61SMark F. Adams 
8955745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr);
8965745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr);
897fc4362bfSMark F. Adams           {
898fc4362bfSMark F. Adams             PetscRandom rctx;
899fc4362bfSMark F. Adams             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
900fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
901fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
902fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
9035b89ad90SMark F. Adams           }
904a2f3521dSMark F. Adams 
905e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
906e696ed0bSMark F. Adams           {
907e696ed0bSMark F. Adams             PetscInt    Istart,Iend,ncols,jj,Ii;
908e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
909e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
910e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
911e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
912e696ed0bSMark F. Adams               if (ncols <= 1) {
913e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
914a94c3b12SMark F. Adams               }
915e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
916a94c3b12SMark F. Adams             }
917a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
918a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
919a94c3b12SMark F. Adams           }
920e696ed0bSMark F. Adams 
921fc4362bfSMark F. Adams           ierr = KSPCreate(wcomm, &eksp);CHKERRQ(ierr);
922806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
923fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
9241a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
9251a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
926f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
927f6536408SMark F. Adams 
928f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
929f6536408SMark F. Adams           ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
930fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
9315a9b9e01SMark F. Adams 
932d3d0db20SJed Brown           /* set PC type to be same as smoother */
933ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
934b2a4f308SMark F. Adams 
9355a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
9365a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
9375a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
938fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
9395a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
9405a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
9415a9b9e01SMark F. Adams 
942fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
9435a9b9e01SMark F. Adams 
944fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
945fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
946fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
947f6536408SMark F. Adams 
948ffc955d6SMark F. Adams           if (pc_gamg->verbose > 0) {
949a94c3b12SMark F. Adams             PetscInt N1, tt;
950a94c3b12SMark F. Adams             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
951a94c3b12SMark 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);
952f6536408SMark F. Adams           }
953fc4362bfSMark F. Adams         }
954038e3b61SMark F. Adams         {
955c5bfad50SMark F. Adams           PetscInt N1, N0;
9560298fd71SBarry Smith           ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr);
9570298fd71SBarry Smith           ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr);
958f6536408SMark F. Adams           /* heuristic - is this crap? */
959b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
9605e7c91beSJed Brown           emin  = emax * pc_gamg->eigtarget[0];
9615e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
962038e3b61SMark F. Adams         }
9636c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
964ffc955d6SMark F. Adams       } /* setup checby flag */
965ffc955d6SMark F. Adams     } /* non-coarse levels */
966737a81a9SMark F. Adams 
967d5280255SMark F. Adams     /* clean up */
968d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
969587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
970587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
9715b89ad90SMark F. Adams     }
9720cbbd2e1SMark F. Adams 
9730cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9740cbbd2e1SMark F. Adams 
9751ac9965eSMark F. Adams     if (PETSC_TRUE) {
97658471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
9779d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
97858471d46SMark F. Adams       ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
97958471d46SMark F. Adams     }
980806fa848SBarry Smith   } else {
9815f8cf99dSMark F. Adams     KSP smoother;
982c5df96a5SBarry Smith     if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",rank,__FUNCT__);
9839d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
9845f8cf99dSMark F. Adams     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9855f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
9869d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9875f8cf99dSMark F. Adams   }
9885b89ad90SMark F. Adams   PetscFunctionReturn(0);
9895b89ad90SMark F. Adams }
9905b89ad90SMark F. Adams 
991eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
9925b89ad90SMark F. Adams /*
9935b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
9945b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
9955b89ad90SMark F. Adams 
9965b89ad90SMark F. Adams    Input Parameter:
9975b89ad90SMark F. Adams .  pc - the preconditioner context
9985b89ad90SMark F. Adams 
9995b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
10005b89ad90SMark F. Adams */
10015b89ad90SMark F. Adams #undef __FUNCT__
10025b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
10035b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
10045b89ad90SMark F. Adams {
10055b89ad90SMark F. Adams   PetscErrorCode ierr;
10065b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
10075b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
10085b89ad90SMark F. Adams 
10095b89ad90SMark F. Adams   PetscFunctionBegin;
10105b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
10115b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
10125b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
10135b89ad90SMark F. Adams   PetscFunctionReturn(0);
10145b89ad90SMark F. Adams }
10155b89ad90SMark F. Adams 
1016676e1743SMark F. Adams 
1017676e1743SMark F. Adams #undef __FUNCT__
1018676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
1019676e1743SMark F. Adams /*@
1020676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
1021676e1743SMark F. Adams    processor reduction.
1022676e1743SMark F. Adams 
1023676e1743SMark F. Adams    Not Collective on PC
1024676e1743SMark F. Adams 
1025676e1743SMark F. Adams    Input Parameters:
1026676e1743SMark F. Adams .  pc - the preconditioner context
1027676e1743SMark F. Adams 
1028676e1743SMark F. Adams 
1029676e1743SMark F. Adams    Options Database Key:
1030676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
1031676e1743SMark F. Adams 
1032676e1743SMark F. Adams    Level: intermediate
1033676e1743SMark F. Adams 
1034676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1035676e1743SMark F. Adams 
1036676e1743SMark F. Adams .seealso: ()
1037676e1743SMark F. Adams @*/
1038676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
1039676e1743SMark F. Adams {
1040676e1743SMark F. Adams   PetscErrorCode ierr;
1041676e1743SMark F. Adams 
1042676e1743SMark F. Adams   PetscFunctionBegin;
1043676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1044676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1045676e1743SMark F. Adams   PetscFunctionReturn(0);
1046676e1743SMark F. Adams }
1047676e1743SMark F. Adams 
1048676e1743SMark F. Adams EXTERN_C_BEGIN
1049676e1743SMark F. Adams #undef __FUNCT__
1050676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
1051676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1052676e1743SMark F. Adams {
1053c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1054c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1055676e1743SMark F. Adams 
1056676e1743SMark F. Adams   PetscFunctionBegin;
10579d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
1058676e1743SMark F. Adams   PetscFunctionReturn(0);
1059676e1743SMark F. Adams }
1060676e1743SMark F. Adams EXTERN_C_END
1061676e1743SMark F. Adams 
1062676e1743SMark F. Adams #undef __FUNCT__
1063389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
1064389730f3SMark F. Adams /*@
1065389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
1066389730f3SMark F. Adams 
1067389730f3SMark F. Adams  Collective on PC
1068389730f3SMark F. Adams 
1069389730f3SMark F. Adams    Input Parameters:
1070389730f3SMark F. Adams .  pc - the preconditioner context
1071389730f3SMark F. Adams 
1072389730f3SMark F. Adams 
1073389730f3SMark F. Adams    Options Database Key:
1074389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
1075389730f3SMark F. Adams 
1076389730f3SMark F. Adams    Level: intermediate
1077389730f3SMark F. Adams 
1078389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1079389730f3SMark F. Adams 
1080389730f3SMark F. Adams .seealso: ()
1081389730f3SMark F. Adams  @*/
1082389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
1083389730f3SMark F. Adams {
1084389730f3SMark F. Adams   PetscErrorCode ierr;
1085389730f3SMark F. Adams 
1086389730f3SMark F. Adams   PetscFunctionBegin;
1087389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1088389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1089389730f3SMark F. Adams   PetscFunctionReturn(0);
1090389730f3SMark F. Adams }
1091389730f3SMark F. Adams 
1092389730f3SMark F. Adams EXTERN_C_BEGIN
1093389730f3SMark F. Adams #undef __FUNCT__
1094389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
1095389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1096389730f3SMark F. Adams {
1097389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1098389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1099389730f3SMark F. Adams 
1100389730f3SMark F. Adams   PetscFunctionBegin;
11019d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
1102389730f3SMark F. Adams   PetscFunctionReturn(0);
1103389730f3SMark F. Adams }
1104389730f3SMark F. Adams EXTERN_C_END
1105389730f3SMark F. Adams 
1106389730f3SMark F. Adams #undef __FUNCT__
11078263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1108676e1743SMark F. Adams /*@
11098263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1110676e1743SMark F. Adams 
1111676e1743SMark F. Adams    Collective on PC
1112676e1743SMark F. Adams 
1113676e1743SMark F. Adams    Input Parameters:
1114676e1743SMark F. Adams .  pc - the preconditioner context
1115676e1743SMark F. Adams 
1116676e1743SMark F. Adams 
1117676e1743SMark F. Adams    Options Database Key:
11188263b398SMark F. Adams .  -pc_gamg_repartition
1119676e1743SMark F. Adams 
1120676e1743SMark F. Adams    Level: intermediate
1121676e1743SMark F. Adams 
1122676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1123676e1743SMark F. Adams 
1124676e1743SMark F. Adams .seealso: ()
1125676e1743SMark F. Adams @*/
11268263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1127676e1743SMark F. Adams {
1128676e1743SMark F. Adams   PetscErrorCode ierr;
1129676e1743SMark F. Adams 
1130676e1743SMark F. Adams   PetscFunctionBegin;
1131676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11328263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1133676e1743SMark F. Adams   PetscFunctionReturn(0);
1134676e1743SMark F. Adams }
1135676e1743SMark F. Adams 
1136676e1743SMark F. Adams EXTERN_C_BEGIN
1137676e1743SMark F. Adams #undef __FUNCT__
11388263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
11398263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1140676e1743SMark F. Adams {
1141c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1142c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1143676e1743SMark F. Adams 
1144676e1743SMark F. Adams   PetscFunctionBegin;
11459d5b6da9SMark F. Adams   pc_gamg->repart = n;
1146676e1743SMark F. Adams   PetscFunctionReturn(0);
1147676e1743SMark F. Adams }
1148676e1743SMark F. Adams EXTERN_C_END
1149676e1743SMark F. Adams 
1150676e1743SMark F. Adams #undef __FUNCT__
1151dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl"
1152dfd5c07aSMark F. Adams /*@
1153dfd5c07aSMark F. Adams    PCGAMGSetReuseProl - Reuse prlongation
1154dfd5c07aSMark F. Adams 
1155dfd5c07aSMark F. Adams    Collective on PC
1156dfd5c07aSMark F. Adams 
1157dfd5c07aSMark F. Adams    Input Parameters:
1158dfd5c07aSMark F. Adams .  pc - the preconditioner context
1159dfd5c07aSMark F. Adams 
1160dfd5c07aSMark F. Adams 
1161dfd5c07aSMark F. Adams    Options Database Key:
1162dfd5c07aSMark F. Adams .  -pc_gamg_reuse_interpolation
1163dfd5c07aSMark F. Adams 
1164dfd5c07aSMark F. Adams    Level: intermediate
1165dfd5c07aSMark F. Adams 
1166dfd5c07aSMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1167dfd5c07aSMark F. Adams 
1168dfd5c07aSMark F. Adams .seealso: ()
1169dfd5c07aSMark F. Adams @*/
1170dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n)
1171dfd5c07aSMark F. Adams {
1172dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1173dfd5c07aSMark F. Adams 
1174dfd5c07aSMark F. Adams   PetscFunctionBegin;
1175dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1176dfd5c07aSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1177dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1178dfd5c07aSMark F. Adams }
1179dfd5c07aSMark F. Adams 
1180dfd5c07aSMark F. Adams EXTERN_C_BEGIN
1181dfd5c07aSMark F. Adams #undef __FUNCT__
1182dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl_GAMG"
1183dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n)
1184dfd5c07aSMark F. Adams {
1185dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1186dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1187dfd5c07aSMark F. Adams 
1188dfd5c07aSMark F. Adams   PetscFunctionBegin;
1189dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1190dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1191dfd5c07aSMark F. Adams }
1192dfd5c07aSMark F. Adams EXTERN_C_END
1193dfd5c07aSMark F. Adams 
1194dfd5c07aSMark F. Adams #undef __FUNCT__
1195ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1196ffc955d6SMark F. Adams /*@
1197ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1198ffc955d6SMark F. Adams 
1199ffc955d6SMark F. Adams    Collective on PC
1200ffc955d6SMark F. Adams 
1201ffc955d6SMark F. Adams    Input Parameters:
1202ffc955d6SMark F. Adams .  pc - the preconditioner context
1203ffc955d6SMark F. Adams 
1204ffc955d6SMark F. Adams 
1205ffc955d6SMark F. Adams    Options Database Key:
1206ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1207ffc955d6SMark F. Adams 
1208ffc955d6SMark F. Adams    Level: intermediate
1209ffc955d6SMark F. Adams 
1210ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1211ffc955d6SMark F. Adams 
1212ffc955d6SMark F. Adams .seealso: ()
1213ffc955d6SMark F. Adams @*/
1214ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1215ffc955d6SMark F. Adams {
1216ffc955d6SMark F. Adams   PetscErrorCode ierr;
1217ffc955d6SMark F. Adams 
1218ffc955d6SMark F. Adams   PetscFunctionBegin;
1219ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1220ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1221ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1222ffc955d6SMark F. Adams }
1223ffc955d6SMark F. Adams 
1224ffc955d6SMark F. Adams EXTERN_C_BEGIN
1225ffc955d6SMark F. Adams #undef __FUNCT__
1226ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
1227ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1228ffc955d6SMark F. Adams {
1229ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1230ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1231ffc955d6SMark F. Adams 
1232ffc955d6SMark F. Adams   PetscFunctionBegin;
1233ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1234ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1235ffc955d6SMark F. Adams }
1236ffc955d6SMark F. Adams EXTERN_C_END
1237ffc955d6SMark F. Adams 
1238ffc955d6SMark F. Adams #undef __FUNCT__
12394ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
12404ef23d27SMark F. Adams /*@
12414ef23d27SMark F. Adams    PCGAMGSetNlevels -
12424ef23d27SMark F. Adams 
12434ef23d27SMark F. Adams    Not collective on PC
12444ef23d27SMark F. Adams 
12454ef23d27SMark F. Adams    Input Parameters:
12464ef23d27SMark F. Adams .  pc - the preconditioner context
12474ef23d27SMark F. Adams 
12484ef23d27SMark F. Adams 
12494ef23d27SMark F. Adams    Options Database Key:
12504ef23d27SMark F. Adams .  -pc_mg_levels
12514ef23d27SMark F. Adams 
12524ef23d27SMark F. Adams    Level: intermediate
12534ef23d27SMark F. Adams 
12544ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12554ef23d27SMark F. Adams 
12564ef23d27SMark F. Adams .seealso: ()
12574ef23d27SMark F. Adams @*/
12584ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12594ef23d27SMark F. Adams {
12604ef23d27SMark F. Adams   PetscErrorCode ierr;
12614ef23d27SMark F. Adams 
12624ef23d27SMark F. Adams   PetscFunctionBegin;
12634ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12644ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
12654ef23d27SMark F. Adams   PetscFunctionReturn(0);
12664ef23d27SMark F. Adams }
12674ef23d27SMark F. Adams 
12684ef23d27SMark F. Adams EXTERN_C_BEGIN
12694ef23d27SMark F. Adams #undef __FUNCT__
12704ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
12714ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12724ef23d27SMark F. Adams {
12734ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
12744ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12754ef23d27SMark F. Adams 
12764ef23d27SMark F. Adams   PetscFunctionBegin;
12779d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12784ef23d27SMark F. Adams   PetscFunctionReturn(0);
12794ef23d27SMark F. Adams }
12804ef23d27SMark F. Adams EXTERN_C_END
12814ef23d27SMark F. Adams 
12824ef23d27SMark F. Adams #undef __FUNCT__
12833542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
12843542efc5SMark F. Adams /*@
12853542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12863542efc5SMark F. Adams 
12873542efc5SMark F. Adams    Not collective on PC
12883542efc5SMark F. Adams 
12893542efc5SMark F. Adams    Input Parameters:
12903542efc5SMark F. Adams .  pc - the preconditioner context
12913542efc5SMark F. Adams 
12923542efc5SMark F. Adams 
12933542efc5SMark F. Adams    Options Database Key:
12943542efc5SMark F. Adams .  -pc_gamg_threshold
12953542efc5SMark F. Adams 
12963542efc5SMark F. Adams    Level: intermediate
12973542efc5SMark F. Adams 
12983542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12993542efc5SMark F. Adams 
13003542efc5SMark F. Adams .seealso: ()
13013542efc5SMark F. Adams @*/
13023542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
13033542efc5SMark F. Adams {
13043542efc5SMark F. Adams   PetscErrorCode ierr;
13053542efc5SMark F. Adams 
13063542efc5SMark F. Adams   PetscFunctionBegin;
13073542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13083542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
13093542efc5SMark F. Adams   PetscFunctionReturn(0);
13103542efc5SMark F. Adams }
13113542efc5SMark F. Adams 
13123542efc5SMark F. Adams EXTERN_C_BEGIN
13133542efc5SMark F. Adams #undef __FUNCT__
13143542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
13153542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
13163542efc5SMark F. Adams {
1317c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1318c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
13193542efc5SMark F. Adams 
13203542efc5SMark F. Adams   PetscFunctionBegin;
13219d5b6da9SMark F. Adams   pc_gamg->threshold = n;
13223542efc5SMark F. Adams   PetscFunctionReturn(0);
13233542efc5SMark F. Adams }
13243542efc5SMark F. Adams EXTERN_C_END
13253542efc5SMark F. Adams 
13263542efc5SMark F. Adams #undef __FUNCT__
13279d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1328676e1743SMark F. Adams /*@
13299d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1330676e1743SMark F. Adams 
1331676e1743SMark F. Adams    Collective on PC
1332676e1743SMark F. Adams 
1333676e1743SMark F. Adams    Input Parameters:
1334676e1743SMark F. Adams .  pc - the preconditioner context
1335676e1743SMark F. Adams 
1336676e1743SMark F. Adams 
1337676e1743SMark F. Adams    Options Database Key:
13383542efc5SMark F. Adams .  -pc_gamg_type
1339676e1743SMark F. Adams 
1340676e1743SMark F. Adams    Level: intermediate
1341676e1743SMark F. Adams 
1342676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1343676e1743SMark F. Adams 
1344676e1743SMark F. Adams .seealso: ()
1345676e1743SMark F. Adams @*/
134619fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1347676e1743SMark F. Adams {
1348676e1743SMark F. Adams   PetscErrorCode ierr;
1349676e1743SMark F. Adams 
1350676e1743SMark F. Adams   PetscFunctionBegin;
1351676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1352806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1353676e1743SMark F. Adams   PetscFunctionReturn(0);
1354676e1743SMark F. Adams }
1355676e1743SMark F. Adams 
1356676e1743SMark F. Adams EXTERN_C_BEGIN
1357676e1743SMark F. Adams #undef __FUNCT__
13589d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
135919fd82e9SBarry Smith PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1360676e1743SMark F. Adams {
13619d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
1362676e1743SMark F. Adams 
1363676e1743SMark F. Adams   PetscFunctionBegin;
1364*ce94432eSBarry Smith   ierr = PetscFunctionListFind(PetscObjectComm((PetscObject)pc),GAMGList,type,PETSC_FALSE,(PetscVoidStarFunction)&r);CHKERRQ(ierr);
13659d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
13669d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1367676e1743SMark F. Adams   PetscFunctionReturn(0);
1368676e1743SMark F. Adams }
1369676e1743SMark F. Adams EXTERN_C_END
1370676e1743SMark F. Adams 
13715b89ad90SMark F. Adams #undef __FUNCT__
13725b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
13735b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
13745b89ad90SMark F. Adams {
1375676e1743SMark F. Adams   PetscErrorCode ierr;
1376676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1377676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1378676e1743SMark F. Adams   PetscBool      flag;
13795e7c91beSJed Brown   PetscInt       two   = 2;
1380*ce94432eSBarry Smith   MPI_Comm       wcomm;
13815b89ad90SMark F. Adams 
13825b89ad90SMark F. Adams   PetscFunctionBegin;
1383*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&wcomm);CHKERRQ(ierr);
1384676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
1385676e1743SMark F. Adams   {
138675b74bdaSMark F. Adams     /* -pc_gamg_verbose */
13879d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
13889d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
13890298fd71SBarry Smith                            &pc_gamg->verbose, NULL);CHKERRQ(ierr);
13908263b398SMark F. Adams     /* -pc_gamg_repartition */
13918263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
13928263b398SMark F. Adams                             "Repartion coarse grids (false)",
13938263b398SMark F. Adams                             "PCGAMGRepartitioning",
13949d5b6da9SMark F. Adams                             pc_gamg->repart,
13959d5b6da9SMark F. Adams                             &pc_gamg->repart,
1396806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1397dfd5c07aSMark F. Adams     /* -pc_gamg_reuse_interpolation */
1398dfd5c07aSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation",
1399dfd5c07aSMark F. Adams                             "Reuse prolongation operator (true)",
1400dfd5c07aSMark F. Adams                             "PCGAMGReuseProl",
1401dfd5c07aSMark F. Adams                             pc_gamg->reuse_prol,
1402dfd5c07aSMark F. Adams                             &pc_gamg->reuse_prol,
1403dfd5c07aSMark F. Adams                             &flag);CHKERRQ(ierr);
1404ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
1405ffc955d6SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1406ffc955d6SMark F. Adams                             "Use aggregation agragates for GASM smoother (false)",
1407ffc955d6SMark F. Adams                             "PCGAMGUseASMAggs",
1408ffc955d6SMark F. Adams                             pc_gamg->use_aggs_in_gasm,
1409ffc955d6SMark F. Adams                             &pc_gamg->use_aggs_in_gasm,
1410806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1411c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1412676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1413676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1414676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
14159d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
14169d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1417806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1418389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1419389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1420389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1421389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
14229d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
14239d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1424806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1425c20e4228SMark F. Adams     /* -pc_gamg_threshold */
14263542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
14273542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
14283542efc5SMark F. Adams                             "PCGAMGSetThreshold",
14299d5b6da9SMark F. Adams                             pc_gamg->threshold,
14309d5b6da9SMark F. Adams                             &pc_gamg->threshold,
1431806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1432806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
1433806fa848SBarry Smith       ierr = PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1434806fa848SBarry Smith     }
14354ef23d27SMark F. Adams 
14360298fd71SBarry Smith     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr);
14374ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
14384ef23d27SMark F. Adams                            "Set number of MG levels",
14394ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
14409d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
14419d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
1442806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1443676e1743SMark F. Adams   }
1444676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
14455b89ad90SMark F. Adams   PetscFunctionReturn(0);
14465b89ad90SMark F. Adams }
14475b89ad90SMark F. Adams 
14485b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
14495b89ad90SMark F. Adams /*MC
1450856830bfSJed Brown      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
14519d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
14525b89ad90SMark F. Adams 
1453280d9858SJed Brown    Options Database Keys:
14545b89ad90SMark F. Adams    Multigrid options(inherited)
1455280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1456280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1457280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
14588c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
14595b89ad90SMark F. Adams 
14605b89ad90SMark F. Adams   Level: intermediate
1461280d9858SJed Brown 
14625b89ad90SMark F. Adams   Concepts: multigrid
14635b89ad90SMark F. Adams 
14645b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1465280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
14665b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
14675b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
14685b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
14695b89ad90SMark F. Adams M*/
14705b89ad90SMark F. Adams EXTERN_C_BEGIN
14715b89ad90SMark F. Adams #undef __FUNCT__
14725b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
14735b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
14745b89ad90SMark F. Adams {
14755b89ad90SMark F. Adams   PetscErrorCode ierr;
14765b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
14775b89ad90SMark F. Adams   PC_MG          *mg;
14780cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14799d5b6da9SMark F. Adams   static long count = 0;
14805ee9c036SSatish Balay #endif
14815b89ad90SMark F. Adams 
14825b89ad90SMark F. Adams   PetscFunctionBegin;
14835b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14845b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
14855b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14865b89ad90SMark F. Adams 
14875b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
14885b89ad90SMark F. Adams   ierr         = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr);
14895b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
1490ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
14915b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
14925b89ad90SMark F. Adams 
14939d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14949d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14959d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14969d5b6da9SMark F. Adams   pc_gamg->data    = 0;
14975b89ad90SMark F. Adams 
14989d5b6da9SMark F. Adams   /* register AMG type */
14999d5b6da9SMark F. Adams   if (!GAMGList) {
1500140e18c1SBarry Smith     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void (*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
1501140e18c1SBarry Smith     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void (*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
15029d5b6da9SMark F. Adams   }
15039d5b6da9SMark F. Adams 
15049d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
15055b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
15065b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
15075b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
15085b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
15095b89ad90SMark F. Adams 
15105b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1511676e1743SMark F. Adams                                            "PCGAMGSetProcEqLim_C",
1512676e1743SMark F. Adams                                            "PCGAMGSetProcEqLim_GAMG",
1513806fa848SBarry Smith                                            PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1514676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1515389730f3SMark F. Adams                                            "PCGAMGSetCoarseEqLim_C",
1516389730f3SMark F. Adams                                            "PCGAMGSetCoarseEqLim_GAMG",
1517806fa848SBarry Smith                                            PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1518389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15198263b398SMark F. Adams                                            "PCGAMGSetRepartitioning_C",
15208263b398SMark F. Adams                                            "PCGAMGSetRepartitioning_GAMG",
1521806fa848SBarry Smith                                            PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
1522676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1523dfd5c07aSMark F. Adams                                            "PCGAMGSetReuseProl_C",
1524dfd5c07aSMark F. Adams                                            "PCGAMGSetReuseProl_GAMG",
1525dfd5c07aSMark F. Adams                                            PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr);
1526dfd5c07aSMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1527ffc955d6SMark F. Adams                                            "PCGAMGSetUseASMAggs_C",
1528ffc955d6SMark F. Adams                                            "PCGAMGSetUseASMAggs_GAMG",
1529806fa848SBarry Smith                                            PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1530ffc955d6SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15313542efc5SMark F. Adams                                            "PCGAMGSetThreshold_C",
15323542efc5SMark F. Adams                                            "PCGAMGSetThreshold_GAMG",
1533806fa848SBarry Smith                                            PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
15343542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
15359d5b6da9SMark F. Adams                                            "PCGAMGSetType_C",
15369d5b6da9SMark F. Adams                                            "PCGAMGSetType_GAMG",
1537806fa848SBarry Smith                                            PCGAMGSetType_GAMG);CHKERRQ(ierr);
15389d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1539dfd5c07aSMark F. Adams   pc_gamg->reuse_prol       = PETSC_TRUE;
1540ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1541038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
1542c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit  = 800;
1543a2f3521dSMark F. Adams   pc_gamg->threshold        = 0.001;
15449d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
15459d5b6da9SMark F. Adams   pc_gamg->verbose          = 0;
15469d5b6da9SMark F. Adams   pc_gamg->emax_id          = -1;
15475e7c91beSJed Brown   pc_gamg->eigtarget[0]     = 0.05;
15485e7c91beSJed Brown   pc_gamg->eigtarget[1]     = 1.05;
15499d5b6da9SMark F. Adams 
15500cbbd2e1SMark F. Adams   /* private events */
15510cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1552785cba28SMark F. Adams   if (count++ == 0) {
1553806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1554806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
15550cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
15560cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
15570cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1558806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1559806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1560806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1561806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1562806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1563806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1564806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1565806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1566806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1567806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1568806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1569806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1570f852f58cSMark F. Adams 
15710cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
15720cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
15730cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1574b4fbaa2aSMark F. Adams     /* create timer stages */
1575b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1576b4fbaa2aSMark F. Adams     {
1577b4fbaa2aSMark F. Adams       char     str[32];
1578b4fbaa2aSMark F. Adams       PetscInt lidx;
1579806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1580806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1581b4fbaa2aSMark F. Adams       for (lidx=1; lidx<9; lidx++) {
1582b4fbaa2aSMark F. Adams         sprintf(str,"MG Level %d",lidx);
1583806fa848SBarry Smith         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1584b4fbaa2aSMark F. Adams       }
1585b4fbaa2aSMark F. Adams     }
1586b4fbaa2aSMark F. Adams #endif
1587b4fbaa2aSMark F. Adams   }
1588b4fbaa2aSMark F. Adams #endif
15890cbbd2e1SMark F. Adams   /* general events */
15900cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
1591806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1592806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1593806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1594806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1595806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1596806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1597806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1598806fa848SBarry Smith   ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr);
15990cbbd2e1SMark F. Adams #endif
16000cbbd2e1SMark F. Adams 
160160a6d8f8SMark F. Adams   /* instantiate derived type */
160260a6d8f8SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
160360a6d8f8SMark F. Adams   {
160460a6d8f8SMark F. Adams     char tname[256] = GAMGAGG;
16050298fd71SBarry Smith     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), NULL);CHKERRQ(ierr);
160660a6d8f8SMark F. Adams     ierr = PCGAMGSetType(pc, tname);CHKERRQ(ierr);
160760a6d8f8SMark F. Adams   }
160860a6d8f8SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
16095b89ad90SMark F. Adams   PetscFunctionReturn(0);
16105b89ad90SMark F. Adams }
16115b89ad90SMark F. Adams EXTERN_C_END
1612