xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 2fb0b3480a5e84094535aa38efa0c85e1d8a6cbc)
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 */
43ce94432eSBarry 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;
1383b4367a7SBarry Smith   MPI_Comm        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;
1433b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
1443b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1453b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
146c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
14711e60469SMark F. Adams   /* RAP */
1489d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
149038e3b61SMark F. Adams 
150a2f3521dSMark F. Adams   /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/
151a2f3521dSMark F. Adams   ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
15271959b99SBarry 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);
1530298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
154a2f3521dSMark F. Adams 
155c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
156a2f3521dSMark F. Adams   {
157472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
1580298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
159c5df96a5SBarry Smith     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
160c5df96a5SBarry Smith     if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1;
161c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
162c5df96a5SBarry Smith     if (isLast) new_size = 1;
163a2f3521dSMark F. Adams   }
164f852f58cSMark F. Adams 
1652fa5cd67SKarl Rupp   if (!repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1662fa5cd67SKarl Rupp   else {
167a2f3521dSMark 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;
168a2f3521dSMark F. Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
169a2f3521dSMark F. Adams     IS             is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
1705a9b9e01SMark F. Adams     VecScatter     vecscat;
17122063be5SMark F. Adams     PetscScalar    *array;
17222063be5SMark F. Adams     Vec            src_crd, dest_crd;
173e33ef3b1SMark F. Adams 
17471959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
17571959b99SBarry 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);
1760cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1770cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
178b4fbaa2aSMark F. Adams #endif
179a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
180c5df96a5SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt), &counts);CHKERRQ(ierr);
181a2f3521dSMark F. Adams     if (repart && !stokes) {
182a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1835a9b9e01SMark F. Adams       Mat adj;
1845a9b9e01SMark F. Adams 
185a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
1863b4367a7SBarry Smith         if (pc_gamg->verbose==1) PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,ncrs_eq);
187a2f3521dSMark F. Adams         else {
188a2f3521dSMark F. Adams           PetscInt n;
1893b4367a7SBarry Smith           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
1903b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n);
191a2f3521dSMark F. Adams         }
192a2f3521dSMark F. Adams       }
1935a9b9e01SMark F. Adams 
194a2f3521dSMark F. Adams       /* get 'adj' */
195c5bfad50SMark F. Adams       if (cr_bs == 1) {
196038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
197806fa848SBarry Smith       } else {
198a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
199eb07cef2SMark F. Adams         Mat               tMat;
200a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
201b4fbaa2aSMark F. Adams         const PetscScalar *vals;
202b4fbaa2aSMark F. Adams         const PetscInt    *idx;
203a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
2049057884aSMark F. Adams         static PetscInt   llev = 0;
205b4fbaa2aSMark F. Adams 
206a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &d_nnz);CHKERRQ(ierr);
207a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &o_nnz);CHKERRQ(ierr);
208a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
209a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
210c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
21158471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
212c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
213c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
21458471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
215a2f3521dSMark F. Adams           if (d_nnz[jj] > ncrs_prim) d_nnz[jj] = ncrs_prim;
216c5bfad50SMark F. Adams           if (o_nnz[jj] > (M/cr_bs-ncrs_prim)) o_nnz[jj] = M/cr_bs-ncrs_prim;
21758471d46SMark F. Adams         }
2186876a03eSMark F. Adams 
2193b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
220806fa848SBarry Smith         ierr = MatSetSizes(tMat, ncrs_prim, ncrs_prim,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
221a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
222a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
223a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
22458471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
2255f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
226eb07cef2SMark F. Adams 
227a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
228c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
22922063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
230eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
231c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
232eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
233eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
234eb07cef2SMark F. Adams           }
23522063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
236eb07cef2SMark F. Adams         }
237eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
238eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
239eb07cef2SMark F. Adams 
240b4fbaa2aSMark F. Adams         if (llev++ == -1) {
241b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2428caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
2433b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
244b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2453bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
246b4fbaa2aSMark F. Adams         }
247b4fbaa2aSMark F. Adams 
248eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
249eb07cef2SMark F. Adams 
250eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
251a2f3521dSMark F. Adams       } /* create 'adj' */
252f150b916SMark F. Adams 
253a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2545a9b9e01SMark F. Adams         char            prefix[256];
2555a9b9e01SMark F. Adams         const char      *pcpre;
256b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
257b4fbaa2aSMark F. Adams         MatPartitioning mpart;
258a4b7d37bSMark F. Adams         IS              proc_is;
259a2f3521dSMark F. Adams         PetscInt        targetPE;
2602f03bc48SMark F. Adams 
2613b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2625ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2639d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2648caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
26559a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
26611e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
267c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
268a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
26911e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2705a9b9e01SMark F. Adams 
2715ef31b24SMark F. Adams         /* collect IS info */
272a2f3521dSMark F. Adams         ierr     = PetscMalloc(ncrs_eq*sizeof(PetscInt), &newproc_idx);CHKERRQ(ierr);
273a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
274a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
275c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
276a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
277c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
278a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
279eb07cef2SMark F. Adams           }
2805ef31b24SMark F. Adams         }
281a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
282a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2835ef31b24SMark F. Adams       }
2845ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2855a9b9e01SMark F. Adams 
2863b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2878263b398SMark F. Adams       if (newproc_idx != 0) {
2888263b398SMark F. Adams         ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2895ef31b24SMark F. Adams       }
290806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
291a2f3521dSMark F. Adams 
292a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
2935a9b9e01SMark F. Adams       /* find factor */
294c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2955a9b9e01SMark F. Adams       else {
2965a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2975a9b9e01SMark F. Adams         jj = -1;
298c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
299c5df96a5SBarry Smith           if (size%kk==0) { /* a candidate */
300c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
3015a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
3025a9b9e01SMark F. Adams             if (fact > best_fact) {
3035a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
3045a9b9e01SMark F. Adams             }
3055a9b9e01SMark F. Adams           }
3065a9b9e01SMark F. Adams         }
3075a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
308a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
3095a9b9e01SMark F. Adams       }
310c5df96a5SBarry Smith       new_size = size/rfactor;
3115a9b9e01SMark F. Adams 
312c5df96a5SBarry Smith       if (new_size==nactive) {
313a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
3145a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
315a2f3521dSMark F. Adams         if (pc_gamg->verbose>0) {
3163b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq);
317a2f3521dSMark F. Adams         }
3185a9b9e01SMark F. Adams         PetscFunctionReturn(0);
3195a9b9e01SMark F. Adams       }
3205a9b9e01SMark F. Adams 
3213b4367a7SBarry Smith       if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq);
322c5df96a5SBarry Smith       targetPE = rank/rfactor;
3233b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
3245a9b9e01SMark F. Adams 
325a2f3521dSMark F. Adams       if (stokes) {
3263b4367a7SBarry Smith         ierr = ISCreateStride(comm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim);CHKERRQ(ierr);
3275a9b9e01SMark F. Adams       }
328a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
329e33ef3b1SMark F. Adams 
33011e60469SMark F. Adams     /*
331a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
33211e60469SMark F. Adams      */
333a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
334a2f3521dSMark F. Adams     if (stokes) {
335a2f3521dSMark F. Adams       ierr = ISPartitioningToNumbering(is_eq_newproc_prim, &is_eq_num_prim);CHKERRQ(ierr);
336806fa848SBarry Smith     } else is_eq_num_prim = is_eq_num;
33711e60469SMark F. Adams     /*
338a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
33911e60469SMark F. Adams      */
340c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
341c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
342a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
343a2f3521dSMark F. Adams     if (stokes) {
344c5df96a5SBarry Smith       ierr          = ISPartitioningCount(is_eq_newproc_prim, size, counts);CHKERRQ(ierr);
345a2f3521dSMark F. Adams       ierr          = ISDestroy(&is_eq_newproc_prim);CHKERRQ(ierr);
346c5df96a5SBarry Smith       ncrs_prim_new = counts[rank]/cr_bs; /* nodes */
347806fa848SBarry Smith     } else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */
348a2f3521dSMark F. Adams 
349a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
3500cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3510cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
352b4fbaa2aSMark F. Adams #endif
353a2f3521dSMark F. Adams 
354a2f3521dSMark F. Adams     /* move data (for primal equations only) */
35522063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
3563b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
357a2f3521dSMark F. Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE);CHKERRQ(ierr);
358470e26f8SMark F. Adams     ierr = VecSetFromOptions(dest_crd);CHKERRQ(ierr); /* this is needed! */
35911e60469SMark F. Adams     /*
3609d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
361c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
36211e60469SMark F. Adams      */
363a2f3521dSMark F. Adams     ierr = PetscMalloc((ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx);CHKERRQ(ierr);
364a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
365a2f3521dSMark F. Adams     for (ii=0,jj=0; ii<ncrs_prim; ii++) {
366c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
367a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
36811e60469SMark F. Adams     }
369a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3703b4367a7SBarry Smith     ierr = ISCreateGeneral(comm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
37192a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
37211e60469SMark F. Adams     /*
37311e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
37411e60469SMark F. Adams      */
375a2f3521dSMark F. Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd);CHKERRQ(ierr);
3769d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
377a2f3521dSMark F. Adams       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
378a2f3521dSMark F. Adams       for (ii=0; ii<ncrs_prim; ii++) {
3799d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
380a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
381c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
382676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
383d3d6bff4SMark F. Adams         }
384038e3b61SMark F. Adams       }
385eb07cef2SMark F. Adams     }
386eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
387eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
38811e60469SMark F. Adams     /*
38911e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
39011e60469SMark F. Adams       to the correct processor
39111e60469SMark F. Adams     */
3920298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
39311e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
39411e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39511e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39611e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
39711e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
39811e60469SMark F. Adams     /*
39911e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
40011e60469SMark F. Adams     */
401c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
402a2f3521dSMark F. Adams     ierr = PetscMalloc(node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
4032fa5cd67SKarl Rupp 
404a2f3521dSMark F. Adams     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
405a2f3521dSMark F. Adams     strideNew        = ncrs_prim_new*ndata_rows;
4062fa5cd67SKarl Rupp 
40711e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
4089d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
409a2f3521dSMark F. Adams       for (ii=0; ii<ncrs_prim_new; ii++) {
4109d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
411a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
412c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
413d3d6bff4SMark F. Adams         }
414038e3b61SMark F. Adams       }
415038e3b61SMark F. Adams     }
41611e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
41711e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
418a2f3521dSMark F. Adams 
419a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
4200cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4210cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
422ed3f9983SMark F. Adams #endif
423a2f3521dSMark F. Adams 
42411e60469SMark F. Adams     /*
42511e60469SMark F. Adams       Invert for MatGetSubMatrix
42611e60469SMark F. Adams     */
427a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
428a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
429c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
430a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
431a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
432a2f3521dSMark F. Adams     }
433a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
4340cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4350cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
4360cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
437ed3f9983SMark F. Adams #endif
438a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
439a2f3521dSMark F. Adams     {
440a2f3521dSMark F. Adams       Mat mat;
441806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
442a2f3521dSMark F. Adams       *a_Amat_crs = mat;
443c5bfad50SMark F. Adams 
444c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
445c5bfad50SMark F. Adams         PetscInt cbs, rbs;
446c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
447c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
448c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
449c5df96a5SBarry 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);
450c5bfad50SMark F. Adams       }
451a2f3521dSMark F. Adams     }
452038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
453a2f3521dSMark F. Adams 
4540cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4550cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
456ed3f9983SMark F. Adams #endif
45711e60469SMark F. Adams     /* prolongator */
45811e60469SMark F. Adams     {
45911e60469SMark F. Adams       IS       findices;
460a2f3521dSMark F. Adams       PetscInt Istart,Iend;
461a2f3521dSMark F. Adams       Mat      Pnew;
462a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4630cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4640cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
465ed3f9983SMark F. Adams #endif
4663b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
467c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
468806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
46911e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
470c5bfad50SMark F. Adams 
471c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
472c5bfad50SMark F. Adams         PetscInt cbs, rbs;
473c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
474c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
475c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
476c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
477c5bfad50SMark F. Adams       }
4780cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4790cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
480ed3f9983SMark F. Adams #endif
4813530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
482a2f3521dSMark F. Adams 
483a2f3521dSMark F. Adams       /* output - repartitioned */
484a2f3521dSMark F. Adams       *a_P_inout = Pnew;
485e33ef3b1SMark F. Adams     }
486a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4875b89ad90SMark F. Adams 
488c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
489a2f3521dSMark F. Adams   }
4905a9b9e01SMark F. Adams 
491a2f3521dSMark F. Adams   /* outout matrix data */
492c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
493c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
494c8b0795cSMark F. Adams     if (llev==0) {
495c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
4963b4367a7SBarry Smith       PetscViewerASCIIOpen(comm,fname,&viewer);
497c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
498c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
499c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
500c8b0795cSMark F. Adams     }
501c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
5023b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
503c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
504c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
505c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
506c8b0795cSMark F. Adams   }
5075b89ad90SMark F. Adams   PetscFunctionReturn(0);
5085b89ad90SMark F. Adams }
5095b89ad90SMark F. Adams 
5105b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5115b89ad90SMark F. Adams /*
5125b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5135b89ad90SMark F. Adams                     by setting data structures and options.
5145b89ad90SMark F. Adams 
5155b89ad90SMark F. Adams    Input Parameter:
5165b89ad90SMark F. Adams .  pc - the preconditioner context
5175b89ad90SMark F. Adams 
5185b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
5195b89ad90SMark F. Adams 
5205b89ad90SMark F. Adams    Notes:
5215b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
5225b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
5235b89ad90SMark F. Adams */
5245b89ad90SMark F. Adams #undef __FUNCT__
5255b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
5269d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5275b89ad90SMark F. Adams {
5285b89ad90SMark F. Adams   PetscErrorCode ierr;
5299d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
5305b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
5312adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
532a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
5333b4367a7SBarry Smith   MPI_Comm       comm;
534c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
535587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
536c8b0795cSMark F. Adams   PetscReal      emaxs[GAMG_MAXLEVELS];
537e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
538a2f3521dSMark F. Adams   GAMGKKTMat     kktMatsArr[GAMG_MAXLEVELS];
539a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
540737a81a9SMark F. Adams   MatInfo        info;
5415169fedaSMark F. Adams   PetscBool      stokes = PETSC_FALSE, redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);
5425ef31b24SMark F. Adams 
5435b89ad90SMark F. Adams   PetscFunctionBegin;
5443b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
5453b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
5463b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
547dfd5c07aSMark F. Adams 
5483b4367a7SBarry Smith   if (pc_gamg->verbose>2) PetscPrintf(comm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",rank,__FUNCT__,pc_gamg->setup_count,pc->setupcalled);
549dfd5c07aSMark F. Adams 
55084d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
551878e152fSMark F. Adams     if (redo_mesh_setup) {
552878e152fSMark F. Adams       /* reset everything */
553878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
554878e152fSMark F. Adams       pc->setupcalled = 0;
555806fa848SBarry Smith     } else {
55684d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
55703a628feSMark F. Adams       /* just do Galerkin grids */
55858471d46SMark F. Adams       Mat          B,dA,dB;
55958471d46SMark F. Adams 
56071959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
5619d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
56258471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5630298fd71SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,NULL);CHKERRQ(ierr);
56458471d46SMark F. Adams         /* (re)set to get dirty flag */
5659d5b6da9SMark F. Adams         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
56658471d46SMark F. Adams 
567*2fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
56803a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
569*2fb0b348SMark F. Adams           if (pc_gamg->setup_count==2 && (pc_gamg->repart || level==0)) {
57003a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
571084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
5722fa5cd67SKarl Rupp 
57303a628feSMark F. Adams             mglevels[level]->A = B;
574806fa848SBarry Smith           } else {
5750298fd71SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B,NULL);CHKERRQ(ierr);
57658471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
57703a628feSMark F. Adams           }
57858471d46SMark F. Adams           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
57958471d46SMark F. Adams           dB   = B;
58058471d46SMark F. Adams         }
5815f8cf99dSMark F. Adams       }
582d5280255SMark F. Adams 
583d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
584d5280255SMark F. Adams 
585d5280255SMark F. Adams       /* PCSetUp_MG seems to insists on setting this to GMRES */
5861ac9965eSMark F. Adams       ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr);
58758471d46SMark F. Adams       PetscFunctionReturn(0);
588eb07cef2SMark F. Adams     }
589878e152fSMark F. Adams   }
590f6536408SMark F. Adams 
5910298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
592eb07cef2SMark F. Adams 
5939d1b15c3SMark F. Adams   ierr = GAMGKKTMatCreate(Pmat, stokes, &kktMatsArr[0]);CHKERRQ(ierr);
5949d1b15c3SMark F. Adams 
595878e152fSMark F. Adams   if (!pc_gamg->data) {
596878e152fSMark F. Adams     if (pc_gamg->orig_data) {
597878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5980298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5992fa5cd67SKarl Rupp 
600878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
601878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
602878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
6032fa5cd67SKarl Rupp 
604878e152fSMark F. Adams       ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
605878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
606806fa848SBarry Smith     } else {
6073b4367a7SBarry Smith       if (!pc_gamg->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
6083b4367a7SBarry Smith       if (stokes) SETERRQ(comm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
6099d1b15c3SMark F. Adams       ierr = pc_gamg->createdefaultdata(pc, kktMatsArr[0].A11);CHKERRQ(ierr);
6109d5b6da9SMark F. Adams     }
611878e152fSMark F. Adams   }
612878e152fSMark F. Adams 
613878e152fSMark F. Adams   /* cache original data for reuse */
614878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
615878e152fSMark F. Adams     ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr);
616878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
617878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
618878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
619878e152fSMark F. Adams   }
620038e3b61SMark F. Adams 
621302f38e8SMark F. Adams   /* get basic dims */
6222fa5cd67SKarl Rupp   if (stokes) bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
6232fa5cd67SKarl Rupp   else {
624302f38e8SMark F. Adams     ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
625302f38e8SMark F. Adams   }
626a2f3521dSMark F. Adams 
627a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
628c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
62984f9421dSMark F. Adams     PetscInt NN = M;
63084f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
63184f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
6323bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
633806fa848SBarry Smith     } else {
634806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
63584f9421dSMark F. Adams     }
636b2a4f308SMark F. Adams     nnz0   = info.nz_used;
637b2a4f308SMark F. Adams     nnztot = info.nz_used;
6383b4367a7SBarry Smith     ierr   = PetscPrintf(comm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
639c5df96a5SBarry Smith                          rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
640c5df96a5SBarry Smith                          (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr);
641c8b0795cSMark F. Adams   }
64284d3f75bSMark F. Adams 
643a2f3521dSMark F. Adams   /* Get A_i and R_i */
644c5df96a5SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
645c5df96a5SBarry Smith        level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit);  /* && (size==1 || nactivepe>1); */
6460205a208SMark F. Adams        level++) {
6475b89ad90SMark F. Adams     level1 = level + 1;
6480cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6490cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
650a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
651a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
652b4fbaa2aSMark F. Adams #endif
653a2f3521dSMark F. Adams #endif
654a2f3521dSMark F. Adams     /* deal with Stokes, get sub matrices */
6559d1b15c3SMark F. Adams     if (level > 0) {
656a2f3521dSMark F. Adams       ierr = GAMGKKTMatCreate(Aarr[level], stokes, &kktMatsArr[level]);CHKERRQ(ierr);
6579d1b15c3SMark F. Adams     }
658c8b0795cSMark F. Adams     { /* construct prolongator */
659725b86d8SJed Brown       Mat              Gmat;
6600cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
661a2f3521dSMark F. Adams       Mat              Prol11,Prol22;
662c8b0795cSMark F. Adams 
663a2f3521dSMark F. Adams       ierr = pc_gamg->graph(pc,kktMatsArr[level].A11, &Gmat);CHKERRQ(ierr);
6640cbbd2e1SMark F. Adams       ierr = pc_gamg->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
665a2f3521dSMark F. Adams       ierr = pc_gamg->prolongator(pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11);CHKERRQ(ierr);
666c8b0795cSMark F. Adams 
667a2f3521dSMark F. Adams       /* could have failed to create new level */
668a2f3521dSMark F. Adams       if (Prol11) {
6699d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6700298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
671a2f3521dSMark F. Adams 
672a2f3521dSMark F. Adams         if (stokes) {
6733b4367a7SBarry Smith           if (!pc_gamg->formkktprol) SETERRQ(comm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
674a2f3521dSMark F. Adams           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
675a2f3521dSMark F. Adams           ierr = pc_gamg->formkktprol(pc, Prol11, kktMatsArr[level].A21, &Prol22);CHKERRQ(ierr);
676c8b0795cSMark F. Adams         }
677c8b0795cSMark F. Adams 
678a2f3521dSMark F. Adams         if (pc_gamg->optprol) {
679c8b0795cSMark F. Adams           /* smooth */
680a2f3521dSMark F. Adams           ierr = pc_gamg->optprol(pc, kktMatsArr[level].A11, &Prol11);CHKERRQ(ierr);
681c8b0795cSMark F. Adams         }
682c8b0795cSMark F. Adams 
683a2f3521dSMark F. Adams         if (stokes) {
684da80777bSKarl Rupp           IS  is_row[2];
685da80777bSKarl Rupp           Mat a[4];
686da80777bSKarl Rupp 
687da80777bSKarl Rupp           is_row[0] = kktMatsArr[level].prim_is; is_row[1] = kktMatsArr[level].constr_is;
6880298fd71SBarry Smith           a[0]      = Prol11; a[1] = NULL; a[2] = NULL; a[3] = Prol22;
6893b4367a7SBarry Smith           ierr      = MatCreateNest(comm,2,is_row, 2, is_row, a, &Parr[level1]);CHKERRQ(ierr);
6902fa5cd67SKarl Rupp         } else Parr[level1] = Prol11;
6910298fd71SBarry Smith       } else Parr[level1] = NULL;
692ffc955d6SMark F. Adams 
693ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
694806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
695ffc955d6SMark F. Adams       }
696ffc955d6SMark F. Adams 
697a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
69841b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
699a2f3521dSMark F. Adams     } /* construct prolongator scope */
7000cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7010cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
702c8b0795cSMark F. Adams #endif
7039d5b6da9SMark F. Adams     /* cache eigen estimate */
7049d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
7059d5b6da9SMark F. Adams       PetscBool flag;
706806fa848SBarry Smith       ierr = PetscObjectComposedDataGetReal((PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
7079d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
708806fa848SBarry Smith     } else emaxs[level] = -1.;
7092adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
710c8b0795cSMark F. Adams     if (!Parr[level1]) {
711806fa848SBarry Smith       if (pc_gamg->verbose) {
7123b4367a7SBarry Smith         ierr =  PetscPrintf(comm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr);
713806fa848SBarry Smith       }
714c8b0795cSMark F. Adams       break;
715c8b0795cSMark F. Adams     }
7160cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7170cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
718b4fbaa2aSMark F. Adams #endif
719a2f3521dSMark F. Adams 
720a2f3521dSMark F. Adams     ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
721806fa848SBarry Smith                        stokes, &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr);
722a2f3521dSMark F. Adams 
7230cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7240cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
725b4fbaa2aSMark F. Adams #endif
726a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
727a2f3521dSMark F. Adams 
728a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0) {
7290cbbd2e1SMark F. Adams       PetscInt NN = M;
7300cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
731a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
7323bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
733806fa848SBarry Smith       } else {
734806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
7350cbbd2e1SMark F. Adams       }
736a2f3521dSMark F. Adams 
7370cbbd2e1SMark F. Adams       nnztot += info.nz_used;
7383b4367a7SBarry Smith       ierr    = PetscPrintf(comm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
739c5df96a5SBarry Smith                             rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
740806fa848SBarry Smith                             (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
741c8b0795cSMark F. Adams     }
742a2f3521dSMark F. Adams 
743a2f3521dSMark F. Adams     /* stop if one node -- could pull back for singular problems */
744c8b0795cSMark F. Adams     if (M/pc_gamg->data_cell_cols < 2) {
74581708dccSMark F. Adams       level++;
74681708dccSMark F. Adams       break;
74781708dccSMark F. Adams     }
7480cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
749b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
750b4fbaa2aSMark F. Adams #endif
751c8b0795cSMark F. Adams   } /* levels */
752c8b0795cSMark F. Adams 
753c8b0795cSMark F. Adams   if (pc_gamg->data) {
754c8b0795cSMark F. Adams     ierr          = PetscFree(pc_gamg->data);CHKERRQ(ierr);
7550298fd71SBarry Smith     pc_gamg->data = NULL;
7565b89ad90SMark F. Adams   }
757c8b0795cSMark F. Adams 
7583b4367a7SBarry Smith   if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
7599d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7605b89ad90SMark F. Adams   fine_level       = level;
7610298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
7625b89ad90SMark F. Adams 
76384d3f75bSMark F. Adams   /* simple setup */
76484d3f75bSMark F. Adams   if (!PETSC_TRUE) {
76584d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
76684d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
76784d3f75bSMark F. Adams          lidx<fine_level;
76884d3f75bSMark F. Adams          lidx++, level--) {
76984d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
77084d3f75bSMark F. Adams       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77184d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
77284d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
77384d3f75bSMark F. Adams     }
77484d3f75bSMark F. Adams     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77584d3f75bSMark F. Adams 
77684d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
777806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
778d5280255SMark F. Adams     /* set default smoothers & set operators */
7799d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
780587fa25dSMark F. Adams          lidx <= fine_level;
781587fa25dSMark F. Adams          lidx++, level--) {
782ffc955d6SMark F. Adams       KSP smoother;
783ffc955d6SMark F. Adams       PC  subpc;
784a2f3521dSMark F. Adams 
7859d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
786f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
787ffc955d6SMark F. Adams 
788a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
789a2f3521dSMark F. Adams       /* set ops */
790a2f3521dSMark F. Adams       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
791a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
792a2f3521dSMark F. Adams 
793a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
794a2f3521dSMark F. Adams       if (stokes) {
795a2f3521dSMark F. Adams         KSP      *ksps;
796a2f3521dSMark F. Adams         PetscInt nn;
797a2f3521dSMark F. Adams         ierr     = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);CHKERRQ(ierr);
798a2f3521dSMark F. Adams         ierr     = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is);CHKERRQ(ierr);
799a2f3521dSMark F. Adams         ierr     = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
800a2f3521dSMark F. Adams         smoother = ksps[0];
801a2f3521dSMark F. Adams         ierr     = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
802a2f3521dSMark F. Adams         ierr     = PetscFree(ksps);CHKERRQ(ierr);
803a2f3521dSMark F. Adams       }
804a2f3521dSMark F. Adams       ierr = GAMGKKTMatDestroy(&kktMatsArr[level]);CHKERRQ(ierr);
805a2f3521dSMark F. Adams 
806a2f3521dSMark F. Adams       /* set defaults */
8076c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
808a2f3521dSMark F. Adams 
809ffc955d6SMark F. Adams       /* override defaults and command line args (!) */
810ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
8112d3561bbSSatish Balay         PetscInt sz;
8122d3561bbSSatish Balay         IS       *is;
813a2f3521dSMark F. Adams 
8142d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
8152d3561bbSSatish Balay         is   = ASMLocalIDsArr[level];
816ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
817ffc955d6SMark F. Adams         if (sz==0) {
818ffc955d6SMark F. Adams           IS       is;
819ffc955d6SMark F. Adams           PetscInt my0,kk;
820ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
821ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
8220298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr);
823a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
824806fa848SBarry Smith         } else {
825a94c3b12SMark F. Adams           PetscInt kk;
8260298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr);
827a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
828a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
829a94c3b12SMark F. Adams           }
830ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
831ffc955d6SMark F. Adams         }
832ffc955d6SMark F. Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
833ffc955d6SMark F. Adams 
8340298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
835ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
836ffc955d6SMark F. Adams         ierr                  = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
837806fa848SBarry Smith       } else {
8389d5b6da9SMark F. Adams         ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr);
839ffc955d6SMark F. Adams       }
840d5280255SMark F. Adams     }
841d5280255SMark F. Adams     {
842d5280255SMark F. Adams       /* coarse grid */
843d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
844d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
845d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
846d5280255SMark F. Adams       ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
847d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
848d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
849d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
850d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
85171959b99SBarry Smith       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
85271959b99SBarry Smith       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
853d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
854d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
8559dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
856*2fb0b348SMark F. Adams       ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
857d5280255SMark F. Adams     }
858d5280255SMark F. Adams 
859d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
860d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
861d5280255SMark F. Adams     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr);
862d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
8633b4367a7SBarry Smith     if (mg->galerkin != 2) SETERRQ(comm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
864d5280255SMark F. Adams 
865d5280255SMark F. Adams     /* create cheby smoothers */
866d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
867d5280255SMark F. Adams          lidx <= fine_level;
868d5280255SMark F. Adams          lidx++, level--) {
869d5280255SMark F. Adams       KSP       smoother;
870ffc955d6SMark F. Adams       PetscBool flag;
871d5280255SMark F. Adams       PC        subpc;
872d5280255SMark F. Adams 
873ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
874a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
875a2f3521dSMark F. Adams 
876a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
877a2f3521dSMark F. Adams       if (stokes) {
878a2f3521dSMark F. Adams         KSP      *ksps;
879a2f3521dSMark F. Adams         PetscInt nn;
880a2f3521dSMark F. Adams         ierr     = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
881a2f3521dSMark F. Adams         smoother = ksps[0];
882a2f3521dSMark F. Adams         ierr     = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
883a2f3521dSMark F. Adams         ierr     = PetscFree(ksps);CHKERRQ(ierr);
884a2f3521dSMark F. Adams       }
885ffc955d6SMark F. Adams 
886ffc955d6SMark F. Adams       /* do my own cheby */
8876c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
888ffc955d6SMark F. Adams       if (flag) {
889ffc955d6SMark F. Adams         PetscReal emax, emin;
890251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
891ffc955d6SMark F. Adams         if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
892587fa25dSMark F. Adams         else { /* eigen estimate 'emax' */
893e696ed0bSMark F. Adams           KSP eksp;
894e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
895b2a4f308SMark F. Adams           Vec bb, xx;
896038e3b61SMark F. Adams 
8975745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr);
8985745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr);
899fc4362bfSMark F. Adams           {
900fc4362bfSMark F. Adams             PetscRandom rctx;
9013b4367a7SBarry Smith             ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
902fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
903fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
904fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
9055b89ad90SMark F. Adams           }
906a2f3521dSMark F. Adams 
907e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
908e696ed0bSMark F. Adams           {
909e696ed0bSMark F. Adams             PetscInt    Istart,Iend,ncols,jj,Ii;
910e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
911e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
912e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
913e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
914e696ed0bSMark F. Adams               if (ncols <= 1) {
915e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
916a94c3b12SMark F. Adams               }
917e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
918a94c3b12SMark F. Adams             }
919a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
920a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
921a94c3b12SMark F. Adams           }
922e696ed0bSMark F. Adams 
9233b4367a7SBarry Smith           ierr = KSPCreate(comm, &eksp);CHKERRQ(ierr);
924806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
925fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
9261a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
9271a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
928f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
929f6536408SMark F. Adams 
930f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
931f6536408SMark F. Adams           ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
932fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
9335a9b9e01SMark F. Adams 
934d3d0db20SJed Brown           /* set PC type to be same as smoother */
935ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
936b2a4f308SMark F. Adams 
9375a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
9385a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
9395a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
940fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
9415a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
9425a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
9435a9b9e01SMark F. Adams 
944fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
9455a9b9e01SMark F. Adams 
946fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
947fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
948fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
949f6536408SMark F. Adams 
950ffc955d6SMark F. Adams           if (pc_gamg->verbose > 0) {
951a94c3b12SMark F. Adams             PetscInt N1, tt;
952a94c3b12SMark F. Adams             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
9533b4367a7SBarry Smith             PetscPrintf(comm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1);
954f6536408SMark F. Adams           }
955fc4362bfSMark F. Adams         }
956038e3b61SMark F. Adams         {
957c5bfad50SMark F. Adams           PetscInt N1, N0;
9580298fd71SBarry Smith           ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr);
9590298fd71SBarry Smith           ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr);
960f6536408SMark F. Adams           /* heuristic - is this crap? */
961b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
9625e7c91beSJed Brown           emin  = emax * pc_gamg->eigtarget[0];
9635e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
964038e3b61SMark F. Adams         }
9656c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
966ffc955d6SMark F. Adams       } /* setup checby flag */
967ffc955d6SMark F. Adams     } /* non-coarse levels */
968737a81a9SMark F. Adams 
969d5280255SMark F. Adams     /* clean up */
970d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
971587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
972587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
9735b89ad90SMark F. Adams     }
9740cbbd2e1SMark F. Adams 
9750cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9760cbbd2e1SMark F. Adams 
9771ac9965eSMark F. Adams     if (PETSC_TRUE) {
97858471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
9799d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
98058471d46SMark F. Adams       ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
98158471d46SMark F. Adams     }
982806fa848SBarry Smith   } else {
9835f8cf99dSMark F. Adams     KSP smoother;
9843b4367a7SBarry Smith     if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",rank,__FUNCT__);
9859d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
9865f8cf99dSMark F. Adams     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9875f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
9889d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9895f8cf99dSMark F. Adams   }
9905b89ad90SMark F. Adams   PetscFunctionReturn(0);
9915b89ad90SMark F. Adams }
9925b89ad90SMark F. Adams 
993eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
9945b89ad90SMark F. Adams /*
9955b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
9965b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
9975b89ad90SMark F. Adams 
9985b89ad90SMark F. Adams    Input Parameter:
9995b89ad90SMark F. Adams .  pc - the preconditioner context
10005b89ad90SMark F. Adams 
10015b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
10025b89ad90SMark F. Adams */
10035b89ad90SMark F. Adams #undef __FUNCT__
10045b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
10055b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
10065b89ad90SMark F. Adams {
10075b89ad90SMark F. Adams   PetscErrorCode ierr;
10085b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
10095b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
10105b89ad90SMark F. Adams 
10115b89ad90SMark F. Adams   PetscFunctionBegin;
10125b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
10135b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
10145b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
10155b89ad90SMark F. Adams   PetscFunctionReturn(0);
10165b89ad90SMark F. Adams }
10175b89ad90SMark F. Adams 
1018676e1743SMark F. Adams 
1019676e1743SMark F. Adams #undef __FUNCT__
1020676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
1021676e1743SMark F. Adams /*@
1022676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
1023676e1743SMark F. Adams    processor reduction.
1024676e1743SMark F. Adams 
1025676e1743SMark F. Adams    Not Collective on PC
1026676e1743SMark F. Adams 
1027676e1743SMark F. Adams    Input Parameters:
1028676e1743SMark F. Adams .  pc - the preconditioner context
1029676e1743SMark F. Adams 
1030676e1743SMark F. Adams 
1031676e1743SMark F. Adams    Options Database Key:
1032676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
1033676e1743SMark F. Adams 
1034676e1743SMark F. Adams    Level: intermediate
1035676e1743SMark F. Adams 
1036676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1037676e1743SMark F. Adams 
1038676e1743SMark F. Adams .seealso: ()
1039676e1743SMark F. Adams @*/
1040676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
1041676e1743SMark F. Adams {
1042676e1743SMark F. Adams   PetscErrorCode ierr;
1043676e1743SMark F. Adams 
1044676e1743SMark F. Adams   PetscFunctionBegin;
1045676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1046676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1047676e1743SMark F. Adams   PetscFunctionReturn(0);
1048676e1743SMark F. Adams }
1049676e1743SMark F. Adams 
1050676e1743SMark F. Adams #undef __FUNCT__
1051676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
10521e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1053676e1743SMark F. Adams {
1054c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1055c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1056676e1743SMark F. Adams 
1057676e1743SMark F. Adams   PetscFunctionBegin;
10589d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
1059676e1743SMark F. Adams   PetscFunctionReturn(0);
1060676e1743SMark F. Adams }
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 #undef __FUNCT__
1093389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
10941e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1095389730f3SMark F. Adams {
1096389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1097389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1098389730f3SMark F. Adams 
1099389730f3SMark F. Adams   PetscFunctionBegin;
11009d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
1101389730f3SMark F. Adams   PetscFunctionReturn(0);
1102389730f3SMark F. Adams }
1103389730f3SMark F. Adams 
1104389730f3SMark F. Adams #undef __FUNCT__
11058263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1106676e1743SMark F. Adams /*@
11078263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1108676e1743SMark F. Adams 
1109676e1743SMark F. Adams    Collective on PC
1110676e1743SMark F. Adams 
1111676e1743SMark F. Adams    Input Parameters:
1112676e1743SMark F. Adams .  pc - the preconditioner context
1113676e1743SMark F. Adams 
1114676e1743SMark F. Adams 
1115676e1743SMark F. Adams    Options Database Key:
11168263b398SMark F. Adams .  -pc_gamg_repartition
1117676e1743SMark F. Adams 
1118676e1743SMark F. Adams    Level: intermediate
1119676e1743SMark F. Adams 
1120676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1121676e1743SMark F. Adams 
1122676e1743SMark F. Adams .seealso: ()
1123676e1743SMark F. Adams @*/
11248263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1125676e1743SMark F. Adams {
1126676e1743SMark F. Adams   PetscErrorCode ierr;
1127676e1743SMark F. Adams 
1128676e1743SMark F. Adams   PetscFunctionBegin;
1129676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11308263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1131676e1743SMark F. Adams   PetscFunctionReturn(0);
1132676e1743SMark F. Adams }
1133676e1743SMark F. Adams 
1134676e1743SMark F. Adams #undef __FUNCT__
11358263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
11361e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1137676e1743SMark F. Adams {
1138c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1139c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1140676e1743SMark F. Adams 
1141676e1743SMark F. Adams   PetscFunctionBegin;
11429d5b6da9SMark F. Adams   pc_gamg->repart = n;
1143676e1743SMark F. Adams   PetscFunctionReturn(0);
1144676e1743SMark F. Adams }
1145676e1743SMark F. Adams 
1146676e1743SMark F. Adams #undef __FUNCT__
1147dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl"
1148dfd5c07aSMark F. Adams /*@
1149dfd5c07aSMark F. Adams    PCGAMGSetReuseProl - Reuse prlongation
1150dfd5c07aSMark F. Adams 
1151dfd5c07aSMark F. Adams    Collective on PC
1152dfd5c07aSMark F. Adams 
1153dfd5c07aSMark F. Adams    Input Parameters:
1154dfd5c07aSMark F. Adams .  pc - the preconditioner context
1155dfd5c07aSMark F. Adams 
1156dfd5c07aSMark F. Adams 
1157dfd5c07aSMark F. Adams    Options Database Key:
1158dfd5c07aSMark F. Adams .  -pc_gamg_reuse_interpolation
1159dfd5c07aSMark F. Adams 
1160dfd5c07aSMark F. Adams    Level: intermediate
1161dfd5c07aSMark F. Adams 
1162dfd5c07aSMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1163dfd5c07aSMark F. Adams 
1164dfd5c07aSMark F. Adams .seealso: ()
1165dfd5c07aSMark F. Adams @*/
1166dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n)
1167dfd5c07aSMark F. Adams {
1168dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1169dfd5c07aSMark F. Adams 
1170dfd5c07aSMark F. Adams   PetscFunctionBegin;
1171dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1172dfd5c07aSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1173dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1174dfd5c07aSMark F. Adams }
1175dfd5c07aSMark F. Adams 
1176dfd5c07aSMark F. Adams #undef __FUNCT__
1177dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl_GAMG"
11781e6b0712SBarry Smith static PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n)
1179dfd5c07aSMark F. Adams {
1180dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1181dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1182dfd5c07aSMark F. Adams 
1183dfd5c07aSMark F. Adams   PetscFunctionBegin;
1184dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1185dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1186dfd5c07aSMark F. Adams }
1187dfd5c07aSMark F. Adams 
1188dfd5c07aSMark F. Adams #undef __FUNCT__
1189ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1190ffc955d6SMark F. Adams /*@
1191ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1192ffc955d6SMark F. Adams 
1193ffc955d6SMark F. Adams    Collective on PC
1194ffc955d6SMark F. Adams 
1195ffc955d6SMark F. Adams    Input Parameters:
1196ffc955d6SMark F. Adams .  pc - the preconditioner context
1197ffc955d6SMark F. Adams 
1198ffc955d6SMark F. Adams 
1199ffc955d6SMark F. Adams    Options Database Key:
1200ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1201ffc955d6SMark F. Adams 
1202ffc955d6SMark F. Adams    Level: intermediate
1203ffc955d6SMark F. Adams 
1204ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1205ffc955d6SMark F. Adams 
1206ffc955d6SMark F. Adams .seealso: ()
1207ffc955d6SMark F. Adams @*/
1208ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1209ffc955d6SMark F. Adams {
1210ffc955d6SMark F. Adams   PetscErrorCode ierr;
1211ffc955d6SMark F. Adams 
1212ffc955d6SMark F. Adams   PetscFunctionBegin;
1213ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1214ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1215ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1216ffc955d6SMark F. Adams }
1217ffc955d6SMark F. Adams 
1218ffc955d6SMark F. Adams #undef __FUNCT__
1219ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
12201e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1221ffc955d6SMark F. Adams {
1222ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1223ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1224ffc955d6SMark F. Adams 
1225ffc955d6SMark F. Adams   PetscFunctionBegin;
1226ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1227ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1228ffc955d6SMark F. Adams }
1229ffc955d6SMark F. Adams 
1230ffc955d6SMark F. Adams #undef __FUNCT__
12314ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
12324ef23d27SMark F. Adams /*@
12334ef23d27SMark F. Adams    PCGAMGSetNlevels -
12344ef23d27SMark F. Adams 
12354ef23d27SMark F. Adams    Not collective on PC
12364ef23d27SMark F. Adams 
12374ef23d27SMark F. Adams    Input Parameters:
12384ef23d27SMark F. Adams .  pc - the preconditioner context
12394ef23d27SMark F. Adams 
12404ef23d27SMark F. Adams 
12414ef23d27SMark F. Adams    Options Database Key:
12424ef23d27SMark F. Adams .  -pc_mg_levels
12434ef23d27SMark F. Adams 
12444ef23d27SMark F. Adams    Level: intermediate
12454ef23d27SMark F. Adams 
12464ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12474ef23d27SMark F. Adams 
12484ef23d27SMark F. Adams .seealso: ()
12494ef23d27SMark F. Adams @*/
12504ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12514ef23d27SMark F. Adams {
12524ef23d27SMark F. Adams   PetscErrorCode ierr;
12534ef23d27SMark F. Adams 
12544ef23d27SMark F. Adams   PetscFunctionBegin;
12554ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12564ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
12574ef23d27SMark F. Adams   PetscFunctionReturn(0);
12584ef23d27SMark F. Adams }
12594ef23d27SMark F. Adams 
12604ef23d27SMark F. Adams #undef __FUNCT__
12614ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
12621e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12634ef23d27SMark F. Adams {
12644ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
12654ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12664ef23d27SMark F. Adams 
12674ef23d27SMark F. Adams   PetscFunctionBegin;
12689d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12694ef23d27SMark F. Adams   PetscFunctionReturn(0);
12704ef23d27SMark F. Adams }
12714ef23d27SMark F. Adams 
12724ef23d27SMark F. Adams #undef __FUNCT__
12733542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
12743542efc5SMark F. Adams /*@
12753542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12763542efc5SMark F. Adams 
12773542efc5SMark F. Adams    Not collective on PC
12783542efc5SMark F. Adams 
12793542efc5SMark F. Adams    Input Parameters:
12803542efc5SMark F. Adams .  pc - the preconditioner context
12813542efc5SMark F. Adams 
12823542efc5SMark F. Adams 
12833542efc5SMark F. Adams    Options Database Key:
12843542efc5SMark F. Adams .  -pc_gamg_threshold
12853542efc5SMark F. Adams 
12863542efc5SMark F. Adams    Level: intermediate
12873542efc5SMark F. Adams 
12883542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12893542efc5SMark F. Adams 
12903542efc5SMark F. Adams .seealso: ()
12913542efc5SMark F. Adams @*/
12923542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
12933542efc5SMark F. Adams {
12943542efc5SMark F. Adams   PetscErrorCode ierr;
12953542efc5SMark F. Adams 
12963542efc5SMark F. Adams   PetscFunctionBegin;
12973542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12983542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
12993542efc5SMark F. Adams   PetscFunctionReturn(0);
13003542efc5SMark F. Adams }
13013542efc5SMark F. Adams 
13023542efc5SMark F. Adams #undef __FUNCT__
13033542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
13041e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
13053542efc5SMark F. Adams {
1306c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1307c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
13083542efc5SMark F. Adams 
13093542efc5SMark F. Adams   PetscFunctionBegin;
13109d5b6da9SMark F. Adams   pc_gamg->threshold = n;
13113542efc5SMark F. Adams   PetscFunctionReturn(0);
13123542efc5SMark F. Adams }
13133542efc5SMark F. Adams 
13143542efc5SMark F. Adams #undef __FUNCT__
13159d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1316676e1743SMark F. Adams /*@
13179d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1318676e1743SMark F. Adams 
1319676e1743SMark F. Adams    Collective on PC
1320676e1743SMark F. Adams 
1321676e1743SMark F. Adams    Input Parameters:
1322676e1743SMark F. Adams .  pc - the preconditioner context
1323676e1743SMark F. Adams 
1324676e1743SMark F. Adams 
1325676e1743SMark F. Adams    Options Database Key:
13263542efc5SMark F. Adams .  -pc_gamg_type
1327676e1743SMark F. Adams 
1328676e1743SMark F. Adams    Level: intermediate
1329676e1743SMark F. Adams 
1330676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1331676e1743SMark F. Adams 
1332676e1743SMark F. Adams .seealso: ()
1333676e1743SMark F. Adams @*/
133419fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1335676e1743SMark F. Adams {
1336676e1743SMark F. Adams   PetscErrorCode ierr;
1337676e1743SMark F. Adams 
1338676e1743SMark F. Adams   PetscFunctionBegin;
1339676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1340806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1341676e1743SMark F. Adams   PetscFunctionReturn(0);
1342676e1743SMark F. Adams }
1343676e1743SMark F. Adams 
1344676e1743SMark F. Adams #undef __FUNCT__
13459d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
13461e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1347676e1743SMark F. Adams {
13489d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
1349676e1743SMark F. Adams 
1350676e1743SMark F. Adams   PetscFunctionBegin;
1351ce94432eSBarry Smith   ierr = PetscFunctionListFind(PetscObjectComm((PetscObject)pc),GAMGList,type,PETSC_FALSE,(PetscVoidStarFunction)&r);CHKERRQ(ierr);
13529d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
13539d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1354676e1743SMark F. Adams   PetscFunctionReturn(0);
1355676e1743SMark F. Adams }
1356676e1743SMark F. Adams 
13575b89ad90SMark F. Adams #undef __FUNCT__
13585b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
13595b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
13605b89ad90SMark F. Adams {
1361676e1743SMark F. Adams   PetscErrorCode ierr;
1362676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1363676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1364676e1743SMark F. Adams   PetscBool      flag;
13655e7c91beSJed Brown   PetscInt       two   = 2;
13663b4367a7SBarry Smith   MPI_Comm       comm;
13675b89ad90SMark F. Adams 
13685b89ad90SMark F. Adams   PetscFunctionBegin;
13693b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1370676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
1371676e1743SMark F. Adams   {
137275b74bdaSMark F. Adams     /* -pc_gamg_verbose */
13739d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
13749d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
13750298fd71SBarry Smith                            &pc_gamg->verbose, NULL);CHKERRQ(ierr);
13768263b398SMark F. Adams     /* -pc_gamg_repartition */
13778263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
13788263b398SMark F. Adams                             "Repartion coarse grids (false)",
13798263b398SMark F. Adams                             "PCGAMGRepartitioning",
13809d5b6da9SMark F. Adams                             pc_gamg->repart,
13819d5b6da9SMark F. Adams                             &pc_gamg->repart,
1382806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1383dfd5c07aSMark F. Adams     /* -pc_gamg_reuse_interpolation */
1384dfd5c07aSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation",
1385dfd5c07aSMark F. Adams                             "Reuse prolongation operator (true)",
1386dfd5c07aSMark F. Adams                             "PCGAMGReuseProl",
1387dfd5c07aSMark F. Adams                             pc_gamg->reuse_prol,
1388dfd5c07aSMark F. Adams                             &pc_gamg->reuse_prol,
1389dfd5c07aSMark F. Adams                             &flag);CHKERRQ(ierr);
1390ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
1391ffc955d6SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1392ffc955d6SMark F. Adams                             "Use aggregation agragates for GASM smoother (false)",
1393ffc955d6SMark F. Adams                             "PCGAMGUseASMAggs",
1394ffc955d6SMark F. Adams                             pc_gamg->use_aggs_in_gasm,
1395ffc955d6SMark F. Adams                             &pc_gamg->use_aggs_in_gasm,
1396806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1397c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1398676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1399676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1400676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
14019d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
14029d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1403806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1404389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1405389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1406389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1407389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
14089d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
14099d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1410806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1411c20e4228SMark F. Adams     /* -pc_gamg_threshold */
14123542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
14133542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
14143542efc5SMark F. Adams                             "PCGAMGSetThreshold",
14159d5b6da9SMark F. Adams                             pc_gamg->threshold,
14169d5b6da9SMark F. Adams                             &pc_gamg->threshold,
1417806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1418806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
14193b4367a7SBarry Smith       ierr = PetscPrintf(comm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1420806fa848SBarry Smith     }
14214ef23d27SMark F. Adams 
14220298fd71SBarry Smith     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr);
14234ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
14244ef23d27SMark F. Adams                            "Set number of MG levels",
14254ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
14269d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
14279d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
1428806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1429676e1743SMark F. Adams   }
1430676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
14315b89ad90SMark F. Adams   PetscFunctionReturn(0);
14325b89ad90SMark F. Adams }
14335b89ad90SMark F. Adams 
14345b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
14355b89ad90SMark F. Adams /*MC
1436856830bfSJed Brown      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
14379d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
14385b89ad90SMark F. Adams 
1439280d9858SJed Brown    Options Database Keys:
14405b89ad90SMark F. Adams    Multigrid options(inherited)
1441280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1442280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1443280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
14448c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
14455b89ad90SMark F. Adams 
14465b89ad90SMark F. Adams   Level: intermediate
1447280d9858SJed Brown 
14485b89ad90SMark F. Adams   Concepts: multigrid
14495b89ad90SMark F. Adams 
14505b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1451280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
14525b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
14535b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
14545b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
14555b89ad90SMark F. Adams M*/
1456b2573a8aSBarry Smith 
14575b89ad90SMark F. Adams #undef __FUNCT__
14585b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
14598cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
14605b89ad90SMark F. Adams {
14615b89ad90SMark F. Adams   PetscErrorCode ierr;
14625b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
14635b89ad90SMark F. Adams   PC_MG          *mg;
14640cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14659d5b6da9SMark F. Adams   static long count = 0;
14665ee9c036SSatish Balay #endif
14675b89ad90SMark F. Adams 
14685b89ad90SMark F. Adams   PetscFunctionBegin;
14695b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14705b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
14715b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14725b89ad90SMark F. Adams 
14735b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
14745b89ad90SMark F. Adams   ierr         = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr);
14755b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
1476ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
14775b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
14785b89ad90SMark F. Adams 
14799d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14809d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14819d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14829d5b6da9SMark F. Adams   pc_gamg->data    = 0;
14835b89ad90SMark F. Adams 
14849d5b6da9SMark F. Adams   /* register AMG type */
14859d5b6da9SMark F. Adams   if (!GAMGList) {
1486140e18c1SBarry Smith     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void (*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
1487140e18c1SBarry Smith     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void (*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
14889d5b6da9SMark F. Adams   }
14899d5b6da9SMark F. Adams 
14909d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
14915b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
14925b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
14935b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
14945b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
14955b89ad90SMark F. Adams 
149600de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C","PCGAMGSetProcEqLim_GAMG",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
149700de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C","PCGAMGSetCoarseEqLim_GAMG",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
149800de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C","PCGAMGSetRepartitioning_GAMG",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
149900de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseProl_C","PCGAMGSetReuseProl_GAMG",PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr);
150000de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C","PCGAMGSetUseASMAggs_GAMG",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
150100de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C","PCGAMGSetThreshold_GAMG",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
150200de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C","PCGAMGSetType_GAMG",PCGAMGSetType_GAMG);CHKERRQ(ierr);
15035a576424SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C","PCGAMGSetNlevels_GAMG",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
15049d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1505dfd5c07aSMark F. Adams   pc_gamg->reuse_prol       = PETSC_TRUE;
1506ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1507038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
1508c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit  = 800;
1509a2f3521dSMark F. Adams   pc_gamg->threshold        = 0.001;
15109d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
15119d5b6da9SMark F. Adams   pc_gamg->verbose          = 0;
15129d5b6da9SMark F. Adams   pc_gamg->emax_id          = -1;
15135e7c91beSJed Brown   pc_gamg->eigtarget[0]     = 0.05;
15145e7c91beSJed Brown   pc_gamg->eigtarget[1]     = 1.05;
15159d5b6da9SMark F. Adams 
15160cbbd2e1SMark F. Adams   /* private events */
15170cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1518785cba28SMark F. Adams   if (count++ == 0) {
1519806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1520806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
15210cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
15220cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
15230cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1524806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1525806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1526806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1527806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1528806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1529806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1530806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1531806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1532806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1533806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1534806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1535806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1536f852f58cSMark F. Adams 
15370cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
15380cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
15390cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1540b4fbaa2aSMark F. Adams     /* create timer stages */
1541b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1542b4fbaa2aSMark F. Adams     {
1543b4fbaa2aSMark F. Adams       char     str[32];
1544b4fbaa2aSMark F. Adams       PetscInt lidx;
1545806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1546806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1547b4fbaa2aSMark F. Adams       for (lidx=1; lidx<9; lidx++) {
1548b4fbaa2aSMark F. Adams         sprintf(str,"MG Level %d",lidx);
1549806fa848SBarry Smith         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1550b4fbaa2aSMark F. Adams       }
1551b4fbaa2aSMark F. Adams     }
1552b4fbaa2aSMark F. Adams #endif
1553b4fbaa2aSMark F. Adams   }
1554b4fbaa2aSMark F. Adams #endif
15550cbbd2e1SMark F. Adams   /* general events */
15560cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
1557806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1558806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1559806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1560806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1561806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1562806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1563806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1564806fa848SBarry Smith   ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr);
15650cbbd2e1SMark F. Adams #endif
15660cbbd2e1SMark F. Adams 
156760a6d8f8SMark F. Adams   /* instantiate derived type */
156860a6d8f8SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
156960a6d8f8SMark F. Adams   {
157060a6d8f8SMark F. Adams     char tname[256] = GAMGAGG;
15710298fd71SBarry Smith     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), NULL);CHKERRQ(ierr);
157260a6d8f8SMark F. Adams     ierr = PCGAMGSetType(pc, tname);CHKERRQ(ierr);
157360a6d8f8SMark F. Adams   }
157460a6d8f8SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
15755b89ad90SMark F. Adams   PetscFunctionReturn(0);
15765b89ad90SMark F. Adams }
1577