xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 9b8ffb577ff6425795f77eba6c25eb9ef8e7c969)
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;
313e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
329d5b6da9SMark F. Adams 
33d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
34d3d6bff4SMark F. Adams #undef __FUNCT__
35d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
36d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
37d3d6bff4SMark F. Adams {
38d3d6bff4SMark F. Adams   PetscErrorCode ierr;
39d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
40d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
41d3d6bff4SMark F. Adams 
42d3d6bff4SMark F. Adams   PetscFunctionBegin;
43a2f3521dSMark F. Adams   if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */
44ce94432eSBarry Smith     PetscPrintf(PetscObjectComm((PetscObject)pc),"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__);
459d5b6da9SMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
4658471d46SMark F. Adams   }
470298fd71SBarry Smith   pc_gamg->data = NULL; pc_gamg->data_sz = 0;
48878e152fSMark F. Adams 
49878e152fSMark F. Adams   if (pc_gamg->orig_data) {
50878e152fSMark F. Adams     ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
51878e152fSMark F. Adams   }
52a2f3521dSMark F. Adams   PetscFunctionReturn(0);
53a2f3521dSMark F. Adams }
54a2f3521dSMark F. Adams 
55a2f3521dSMark F. Adams /* private 2x2 Mat Nest for Stokes */
56a2f3521dSMark F. Adams typedef struct {
57a2f3521dSMark F. Adams   Mat A11,A21,A12,Amat;
58a2f3521dSMark F. Adams   IS  prim_is,constr_is;
59a2f3521dSMark F. Adams } GAMGKKTMat;
60a2f3521dSMark F. Adams 
61a2f3521dSMark F. Adams #undef __FUNCT__
62a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatCreate"
63a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatCreate(Mat A, PetscBool iskkt, GAMGKKTMat *out)
64a2f3521dSMark F. Adams {
65a2f3521dSMark F. Adams   PetscFunctionBegin;
66a2f3521dSMark F. Adams   out->Amat = A;
67a2f3521dSMark F. Adams   if (iskkt) {
68a2f3521dSMark F. Adams     PetscErrorCode ierr;
69a2f3521dSMark F. Adams     IS             is_constraint, is_prime;
70a2f3521dSMark F. Adams     PetscInt       nmin,nmax;
71a2f3521dSMark F. Adams 
72a2f3521dSMark F. Adams     ierr = MatGetOwnershipRange(A, &nmin, &nmax);CHKERRQ(ierr);
73a2f3521dSMark F. Adams     ierr = MatFindZeroDiagonals(A, &is_constraint);CHKERRQ(ierr);
74a2f3521dSMark F. Adams     ierr = ISComplement(is_constraint, nmin, nmax, &is_prime);CHKERRQ(ierr);
752fa5cd67SKarl Rupp 
76a2f3521dSMark F. Adams     out->prim_is   = is_prime;
77a2f3521dSMark F. Adams     out->constr_is = is_constraint;
78a2f3521dSMark F. Adams 
79a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_prime, is_prime,      MAT_INITIAL_MATRIX, &out->A11);CHKERRQ(ierr);
80a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12);CHKERRQ(ierr);
81a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21);CHKERRQ(ierr);
82806fa848SBarry Smith   } else {
83a2f3521dSMark F. Adams     out->A11       = A;
840298fd71SBarry Smith     out->A21       = NULL;
850298fd71SBarry Smith     out->A12       = NULL;
860298fd71SBarry Smith     out->prim_is   = NULL;
870298fd71SBarry Smith     out->constr_is = NULL;
88a2f3521dSMark F. Adams   }
89a2f3521dSMark F. Adams   PetscFunctionReturn(0);
90a2f3521dSMark F. Adams }
91a2f3521dSMark F. Adams 
92a2f3521dSMark F. Adams #undef __FUNCT__
93a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatDestroy"
94a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatDestroy(GAMGKKTMat *mat)
95a2f3521dSMark F. Adams {
96a2f3521dSMark F. Adams   PetscErrorCode ierr;
97a2f3521dSMark F. Adams 
98a2f3521dSMark F. Adams   PetscFunctionBegin;
99a2f3521dSMark F. Adams   if (mat->A11 && mat->A11 != mat->Amat) {
100a2f3521dSMark F. Adams     ierr = MatDestroy(&mat->A11);CHKERRQ(ierr);
101a2f3521dSMark F. Adams   }
102a2f3521dSMark F. Adams   ierr = MatDestroy(&mat->A21);CHKERRQ(ierr);
103a2f3521dSMark F. Adams   ierr = MatDestroy(&mat->A12);CHKERRQ(ierr);
104a2f3521dSMark F. Adams 
105a2f3521dSMark F. Adams   ierr = ISDestroy(&mat->prim_is);CHKERRQ(ierr);
106a2f3521dSMark F. Adams   ierr = ISDestroy(&mat->constr_is);CHKERRQ(ierr);
107d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
108d3d6bff4SMark F. Adams }
109d3d6bff4SMark F. Adams 
1105b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1115b89ad90SMark F. Adams /*
112a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
113a147abb0SMark F. Adams      of active processors.
1145b89ad90SMark F. Adams 
1155b89ad90SMark F. Adams    Input Parameter:
116a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
117a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
1189d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
119c5bfad50SMark F. Adams    . cr_bs - coarse block size
1209d5b6da9SMark F. Adams    . isLast -
121a2f3521dSMark F. Adams    . stokes -
1223530afc2SMark F. Adams    In/Output Parameter:
123a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
124afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
12511e60469SMark F. Adams    Output Parameter:
1263530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1275b89ad90SMark F. Adams */
1285cb416c2SMark F. Adams 
1295b89ad90SMark F. Adams #undef __FUNCT__
1308263b398SMark F. Adams #define __FUNCT__ "createLevel"
1312fa5cd67SKarl 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)
1325b89ad90SMark F. Adams {
133a2f3521dSMark F. Adams   PetscErrorCode  ierr;
1349d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
135486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
1369d5b6da9SMark F. Adams   const PetscBool repart      = pc_gamg->repart;
1379d5b6da9SMark F. Adams   const PetscInt  min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
138a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
1393b4367a7SBarry Smith   MPI_Comm        comm;
140c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
141c5bfad50SMark F. Adams   PetscInt        ncrs_eq,ncrs_prim,f_bs;
1425b89ad90SMark F. Adams 
1435b89ad90SMark F. Adams   PetscFunctionBegin;
1443b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
1453b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1463b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
147c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
14811e60469SMark F. Adams   /* RAP */
1499d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
150038e3b61SMark F. Adams 
151a2f3521dSMark F. Adams   /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/
152a2f3521dSMark F. Adams   ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
15371959b99SBarry 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);
1540298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
155a2f3521dSMark F. Adams 
156c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
157a2f3521dSMark F. Adams   {
158472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
1590298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
160c5df96a5SBarry Smith     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
161c5df96a5SBarry Smith     if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1;
162c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
163c5df96a5SBarry Smith     if (isLast) new_size = 1;
164a2f3521dSMark F. Adams   }
165f852f58cSMark F. Adams 
1662fa5cd67SKarl Rupp   if (!repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1672fa5cd67SKarl Rupp   else {
168a2f3521dSMark 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;
169a2f3521dSMark F. Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
170a2f3521dSMark F. Adams     IS             is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
1715a9b9e01SMark F. Adams     VecScatter     vecscat;
17222063be5SMark F. Adams     PetscScalar    *array;
17322063be5SMark F. Adams     Vec            src_crd, dest_crd;
174e33ef3b1SMark F. Adams 
17571959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
17671959b99SBarry 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);
1770cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1780cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
179b4fbaa2aSMark F. Adams #endif
180a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
181c5df96a5SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt), &counts);CHKERRQ(ierr);
182a2f3521dSMark F. Adams     if (repart && !stokes) {
183a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1845a9b9e01SMark F. Adams       Mat adj;
1855a9b9e01SMark F. Adams 
186a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
1873b4367a7SBarry 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);
188a2f3521dSMark F. Adams         else {
189a2f3521dSMark F. Adams           PetscInt n;
1903b4367a7SBarry Smith           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
1913b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n);
192a2f3521dSMark F. Adams         }
193a2f3521dSMark F. Adams       }
1945a9b9e01SMark F. Adams 
195a2f3521dSMark F. Adams       /* get 'adj' */
196c5bfad50SMark F. Adams       if (cr_bs == 1) {
197038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
198806fa848SBarry Smith       } else {
199a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
200eb07cef2SMark F. Adams         Mat               tMat;
201a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
202b4fbaa2aSMark F. Adams         const PetscScalar *vals;
203b4fbaa2aSMark F. Adams         const PetscInt    *idx;
204a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
2059057884aSMark F. Adams         static PetscInt   llev = 0;
206b4fbaa2aSMark F. Adams 
207a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &d_nnz);CHKERRQ(ierr);
208a2f3521dSMark F. Adams         ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &o_nnz);CHKERRQ(ierr);
209a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
210a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
211c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
21258471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
213c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
214c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
21558471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
216a2f3521dSMark F. Adams           if (d_nnz[jj] > ncrs_prim) d_nnz[jj] = ncrs_prim;
217c5bfad50SMark F. Adams           if (o_nnz[jj] > (M/cr_bs-ncrs_prim)) o_nnz[jj] = M/cr_bs-ncrs_prim;
21858471d46SMark F. Adams         }
2196876a03eSMark F. Adams 
2203b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
221806fa848SBarry Smith         ierr = MatSetSizes(tMat, ncrs_prim, ncrs_prim,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
222a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
223a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
224a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
22558471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
2265f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
227eb07cef2SMark F. Adams 
228a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
229c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
23022063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
231eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
232c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
233eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
234eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
235eb07cef2SMark F. Adams           }
23622063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
237eb07cef2SMark F. Adams         }
238eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
239eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
240eb07cef2SMark F. Adams 
241b4fbaa2aSMark F. Adams         if (llev++ == -1) {
242b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2438caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
2443b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
245b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2463bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
247b4fbaa2aSMark F. Adams         }
248b4fbaa2aSMark F. Adams 
249eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
250eb07cef2SMark F. Adams 
251eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
252a2f3521dSMark F. Adams       } /* create 'adj' */
253f150b916SMark F. Adams 
254a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2555a9b9e01SMark F. Adams         char            prefix[256];
2565a9b9e01SMark F. Adams         const char      *pcpre;
257b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
258b4fbaa2aSMark F. Adams         MatPartitioning mpart;
259a4b7d37bSMark F. Adams         IS              proc_is;
260a2f3521dSMark F. Adams         PetscInt        targetPE;
2612f03bc48SMark F. Adams 
2623b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2635ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2649d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2658caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
26659a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
26711e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
268c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
269a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
27011e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2715a9b9e01SMark F. Adams 
2725ef31b24SMark F. Adams         /* collect IS info */
273a2f3521dSMark F. Adams         ierr     = PetscMalloc(ncrs_eq*sizeof(PetscInt), &newproc_idx);CHKERRQ(ierr);
274a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
275a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
276c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
277a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
278c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
279a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
280eb07cef2SMark F. Adams           }
2815ef31b24SMark F. Adams         }
282a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
283a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2845ef31b24SMark F. Adams       }
2855ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2865a9b9e01SMark F. Adams 
2873b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2888263b398SMark F. Adams       if (newproc_idx != 0) {
2898263b398SMark F. Adams         ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2905ef31b24SMark F. Adams       }
291806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
292a2f3521dSMark F. Adams 
293a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
2945a9b9e01SMark F. Adams       /* find factor */
295c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2965a9b9e01SMark F. Adams       else {
2975a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2985a9b9e01SMark F. Adams         jj = -1;
299c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
300c5df96a5SBarry Smith           if (size%kk==0) { /* a candidate */
301c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
3025a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
3035a9b9e01SMark F. Adams             if (fact > best_fact) {
3045a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
3055a9b9e01SMark F. Adams             }
3065a9b9e01SMark F. Adams           }
3075a9b9e01SMark F. Adams         }
3085a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
309a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
3105a9b9e01SMark F. Adams       }
311c5df96a5SBarry Smith       new_size = size/rfactor;
3125a9b9e01SMark F. Adams 
313c5df96a5SBarry Smith       if (new_size==nactive) {
314a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
3155a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
316a2f3521dSMark F. Adams         if (pc_gamg->verbose>0) {
3173b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq);
318a2f3521dSMark F. Adams         }
3195a9b9e01SMark F. Adams         PetscFunctionReturn(0);
3205a9b9e01SMark F. Adams       }
3215a9b9e01SMark F. Adams 
3223b4367a7SBarry Smith       if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq);
323c5df96a5SBarry Smith       targetPE = rank/rfactor;
3243b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
3255a9b9e01SMark F. Adams 
326a2f3521dSMark F. Adams       if (stokes) {
3273b4367a7SBarry Smith         ierr = ISCreateStride(comm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim);CHKERRQ(ierr);
3285a9b9e01SMark F. Adams       }
329a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
330e33ef3b1SMark F. Adams 
33111e60469SMark F. Adams     /*
332a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
33311e60469SMark F. Adams      */
334a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
335a2f3521dSMark F. Adams     if (stokes) {
336a2f3521dSMark F. Adams       ierr = ISPartitioningToNumbering(is_eq_newproc_prim, &is_eq_num_prim);CHKERRQ(ierr);
337806fa848SBarry Smith     } else is_eq_num_prim = is_eq_num;
33811e60469SMark F. Adams     /*
339a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
34011e60469SMark F. Adams      */
341c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
342c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
343a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
344a2f3521dSMark F. Adams     if (stokes) {
345c5df96a5SBarry Smith       ierr          = ISPartitioningCount(is_eq_newproc_prim, size, counts);CHKERRQ(ierr);
346a2f3521dSMark F. Adams       ierr          = ISDestroy(&is_eq_newproc_prim);CHKERRQ(ierr);
347c5df96a5SBarry Smith       ncrs_prim_new = counts[rank]/cr_bs; /* nodes */
348806fa848SBarry Smith     } else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */
349a2f3521dSMark F. Adams 
350a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
3510cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3520cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
353b4fbaa2aSMark F. Adams #endif
354a2f3521dSMark F. Adams 
355a2f3521dSMark F. Adams     /* move data (for primal equations only) */
35622063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
3573b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
358a2f3521dSMark F. Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE);CHKERRQ(ierr);
359470e26f8SMark F. Adams     ierr = VecSetFromOptions(dest_crd);CHKERRQ(ierr); /* this is needed! */
36011e60469SMark F. Adams     /*
3619d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
362c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
36311e60469SMark F. Adams      */
364a2f3521dSMark F. Adams     ierr = PetscMalloc((ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx);CHKERRQ(ierr);
365a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
366a2f3521dSMark F. Adams     for (ii=0,jj=0; ii<ncrs_prim; ii++) {
367c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
368a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
36911e60469SMark F. Adams     }
370a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3713b4367a7SBarry Smith     ierr = ISCreateGeneral(comm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
37292a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
37311e60469SMark F. Adams     /*
37411e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
37511e60469SMark F. Adams      */
376a2f3521dSMark F. Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd);CHKERRQ(ierr);
3779d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
378a2f3521dSMark F. Adams       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
379a2f3521dSMark F. Adams       for (ii=0; ii<ncrs_prim; ii++) {
3809d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
381a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
382c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
383676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
384d3d6bff4SMark F. Adams         }
385038e3b61SMark F. Adams       }
386eb07cef2SMark F. Adams     }
387eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
388eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
38911e60469SMark F. Adams     /*
39011e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
39111e60469SMark F. Adams       to the correct processor
39211e60469SMark F. Adams     */
3930298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
39411e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
39511e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39611e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39711e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
39811e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
39911e60469SMark F. Adams     /*
40011e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
40111e60469SMark F. Adams     */
402c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
403a2f3521dSMark F. Adams     ierr = PetscMalloc(node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
4042fa5cd67SKarl Rupp 
405a2f3521dSMark F. Adams     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
406a2f3521dSMark F. Adams     strideNew        = ncrs_prim_new*ndata_rows;
4072fa5cd67SKarl Rupp 
40811e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
4099d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
410a2f3521dSMark F. Adams       for (ii=0; ii<ncrs_prim_new; ii++) {
4119d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
412a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
413c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
414d3d6bff4SMark F. Adams         }
415038e3b61SMark F. Adams       }
416038e3b61SMark F. Adams     }
41711e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
41811e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
419a2f3521dSMark F. Adams 
420a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
4210cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4220cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
423ed3f9983SMark F. Adams #endif
424a2f3521dSMark F. Adams 
42511e60469SMark F. Adams     /*
42611e60469SMark F. Adams       Invert for MatGetSubMatrix
42711e60469SMark F. Adams     */
428a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
429a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
430c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
431a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
432a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
433a2f3521dSMark F. Adams     }
434a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
4350cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4360cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
4370cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
438ed3f9983SMark F. Adams #endif
439a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
440a2f3521dSMark F. Adams     {
441a2f3521dSMark F. Adams       Mat mat;
442806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
443a2f3521dSMark F. Adams       *a_Amat_crs = mat;
444c5bfad50SMark F. Adams 
445c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
446c5bfad50SMark F. Adams         PetscInt cbs, rbs;
447c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
448c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
449c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
450c5df96a5SBarry 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);
451c5bfad50SMark F. Adams       }
452a2f3521dSMark F. Adams     }
453038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
454a2f3521dSMark F. Adams 
4550cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4560cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
457ed3f9983SMark F. Adams #endif
45811e60469SMark F. Adams     /* prolongator */
45911e60469SMark F. Adams     {
46011e60469SMark F. Adams       IS       findices;
461a2f3521dSMark F. Adams       PetscInt Istart,Iend;
462a2f3521dSMark F. Adams       Mat      Pnew;
463a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4640cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4650cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
466ed3f9983SMark F. Adams #endif
4673b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
468c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
469806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
47011e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
471c5bfad50SMark F. Adams 
472c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
473c5bfad50SMark F. Adams         PetscInt cbs, rbs;
474c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
475c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
476c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
477c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
478c5bfad50SMark F. Adams       }
4790cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4800cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
481ed3f9983SMark F. Adams #endif
4823530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
483a2f3521dSMark F. Adams 
484a2f3521dSMark F. Adams       /* output - repartitioned */
485a2f3521dSMark F. Adams       *a_P_inout = Pnew;
486e33ef3b1SMark F. Adams     }
487a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4885b89ad90SMark F. Adams 
489c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
490a2f3521dSMark F. Adams   }
4915a9b9e01SMark F. Adams 
492a2f3521dSMark F. Adams   /* outout matrix data */
493c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
494c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
495c8b0795cSMark F. Adams     if (llev==0) {
496c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
4973b4367a7SBarry Smith       PetscViewerASCIIOpen(comm,fname,&viewer);
498c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
499c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
500c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
501c8b0795cSMark F. Adams     }
502c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
5033b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
504c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
505c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
506c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
507c8b0795cSMark F. Adams   }
5085b89ad90SMark F. Adams   PetscFunctionReturn(0);
5095b89ad90SMark F. Adams }
5105b89ad90SMark F. Adams 
5115b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5125b89ad90SMark F. Adams /*
5135b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5145b89ad90SMark F. Adams                     by setting data structures and options.
5155b89ad90SMark F. Adams 
5165b89ad90SMark F. Adams    Input Parameter:
5175b89ad90SMark F. Adams .  pc - the preconditioner context
5185b89ad90SMark F. Adams 
5195b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
5205b89ad90SMark F. Adams 
5215b89ad90SMark F. Adams    Notes:
5225b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
5235b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
5245b89ad90SMark F. Adams */
5255b89ad90SMark F. Adams #undef __FUNCT__
5265b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
5279d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5285b89ad90SMark F. Adams {
5295b89ad90SMark F. Adams   PetscErrorCode ierr;
5309d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
5315b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
5322adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
533a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
5343b4367a7SBarry Smith   MPI_Comm       comm;
535c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
536587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
537c8b0795cSMark F. Adams   PetscReal      emaxs[GAMG_MAXLEVELS];
538e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
539a2f3521dSMark F. Adams   GAMGKKTMat     kktMatsArr[GAMG_MAXLEVELS];
540a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
541737a81a9SMark F. Adams   MatInfo        info;
5425169fedaSMark F. Adams   PetscBool      stokes = PETSC_FALSE, redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);
5435ef31b24SMark F. Adams 
5445b89ad90SMark F. Adams   PetscFunctionBegin;
5453b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
5463b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
5473b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
548dfd5c07aSMark F. Adams 
5493b4367a7SBarry 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);
550dfd5c07aSMark F. Adams 
55184d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
552878e152fSMark F. Adams     if (redo_mesh_setup) {
553878e152fSMark F. Adams       /* reset everything */
554878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
555878e152fSMark F. Adams       pc->setupcalled = 0;
556806fa848SBarry Smith     } else {
55784d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
55803a628feSMark F. Adams       /* just do Galerkin grids */
55958471d46SMark F. Adams       Mat          B,dA,dB;
56058471d46SMark F. Adams 
56171959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
5629d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
56358471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5640298fd71SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,NULL);CHKERRQ(ierr);
56558471d46SMark F. Adams         /* (re)set to get dirty flag */
5669d5b6da9SMark F. Adams         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
56758471d46SMark F. Adams 
5682fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
56903a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
5702fb0b348SMark F. Adams           if (pc_gamg->setup_count==2 && (pc_gamg->repart || level==0)) {
57103a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
572084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
5732fa5cd67SKarl Rupp 
57403a628feSMark F. Adams             mglevels[level]->A = B;
575806fa848SBarry Smith           } else {
5760298fd71SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B,NULL);CHKERRQ(ierr);
57758471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
57803a628feSMark F. Adams           }
57958471d46SMark F. Adams           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
58058471d46SMark F. Adams           dB   = B;
58158471d46SMark F. Adams         }
5825f8cf99dSMark F. Adams       }
583d5280255SMark F. Adams 
584d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
585d5280255SMark F. Adams 
586d5280255SMark F. Adams       /* PCSetUp_MG seems to insists on setting this to GMRES */
5871ac9965eSMark F. Adams       ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr);
58858471d46SMark F. Adams       PetscFunctionReturn(0);
589eb07cef2SMark F. Adams     }
590878e152fSMark F. Adams   }
591f6536408SMark F. Adams 
5920298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
593eb07cef2SMark F. Adams 
5949d1b15c3SMark F. Adams   ierr = GAMGKKTMatCreate(Pmat, stokes, &kktMatsArr[0]);CHKERRQ(ierr);
5959d1b15c3SMark F. Adams 
596878e152fSMark F. Adams   if (!pc_gamg->data) {
597878e152fSMark F. Adams     if (pc_gamg->orig_data) {
598878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5990298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
6002fa5cd67SKarl Rupp 
601878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
602878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
603878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
6042fa5cd67SKarl Rupp 
605878e152fSMark F. Adams       ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
606878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
607806fa848SBarry Smith     } else {
6081ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
6093b4367a7SBarry Smith       if (stokes) SETERRQ(comm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
6101ab5ffc9SJed Brown       ierr = pc_gamg->ops->createdefaultdata(pc, kktMatsArr[0].A11);CHKERRQ(ierr);
6119d5b6da9SMark F. Adams     }
612878e152fSMark F. Adams   }
613878e152fSMark F. Adams 
614878e152fSMark F. Adams   /* cache original data for reuse */
615878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
616878e152fSMark F. Adams     ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr);
617878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
618878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
619878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
620878e152fSMark F. Adams   }
621038e3b61SMark F. Adams 
622302f38e8SMark F. Adams   /* get basic dims */
6232fa5cd67SKarl Rupp   if (stokes) bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
6242fa5cd67SKarl Rupp   else {
625302f38e8SMark F. Adams     ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
626302f38e8SMark F. Adams   }
627a2f3521dSMark F. Adams 
628a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
629c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
63084f9421dSMark F. Adams     PetscInt NN = M;
63184f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
63284f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
6333bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
634806fa848SBarry Smith     } else {
635806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
63684f9421dSMark F. Adams     }
637b2a4f308SMark F. Adams     nnz0   = info.nz_used;
638b2a4f308SMark F. Adams     nnztot = info.nz_used;
6393b4367a7SBarry 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",
640c5df96a5SBarry Smith                          rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
641c5df96a5SBarry Smith                          (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr);
642c8b0795cSMark F. Adams   }
64384d3f75bSMark F. Adams 
644a2f3521dSMark F. Adams   /* Get A_i and R_i */
645c5df96a5SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
646c5df96a5SBarry Smith        level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit);  /* && (size==1 || nactivepe>1); */
6470205a208SMark F. Adams        level++) {
6485b89ad90SMark F. Adams     level1 = level + 1;
6490cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6500cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
651a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
652a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
653b4fbaa2aSMark F. Adams #endif
654a2f3521dSMark F. Adams #endif
655a2f3521dSMark F. Adams     /* deal with Stokes, get sub matrices */
6569d1b15c3SMark F. Adams     if (level > 0) {
657a2f3521dSMark F. Adams       ierr = GAMGKKTMatCreate(Aarr[level], stokes, &kktMatsArr[level]);CHKERRQ(ierr);
6589d1b15c3SMark F. Adams     }
659c8b0795cSMark F. Adams     { /* construct prolongator */
660725b86d8SJed Brown       Mat              Gmat;
6610cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
662a2f3521dSMark F. Adams       Mat              Prol11,Prol22;
663c8b0795cSMark F. Adams 
6641ab5ffc9SJed Brown       ierr = pc_gamg->ops->graph(pc,kktMatsArr[level].A11, &Gmat);CHKERRQ(ierr);
6651ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
6661ab5ffc9SJed Brown       ierr = pc_gamg->ops->prolongator(pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11);CHKERRQ(ierr);
667c8b0795cSMark F. Adams 
668a2f3521dSMark F. Adams       /* could have failed to create new level */
669a2f3521dSMark F. Adams       if (Prol11) {
6709d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6710298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
672a2f3521dSMark F. Adams 
673a2f3521dSMark F. Adams         if (stokes) {
6741ab5ffc9SJed Brown           if (!pc_gamg->ops->formkktprol) SETERRQ(comm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
675a2f3521dSMark F. Adams           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
6761ab5ffc9SJed Brown           ierr = pc_gamg->ops->formkktprol(pc, Prol11, kktMatsArr[level].A21, &Prol22);CHKERRQ(ierr);
677c8b0795cSMark F. Adams         }
678c8b0795cSMark F. Adams 
6791ab5ffc9SJed Brown         if (pc_gamg->ops->optprol) {
680c8b0795cSMark F. Adams           /* smooth */
6811ab5ffc9SJed Brown           ierr = pc_gamg->ops->optprol(pc, kktMatsArr[level].A11, &Prol11);CHKERRQ(ierr);
682c8b0795cSMark F. Adams         }
683c8b0795cSMark F. Adams 
684a2f3521dSMark F. Adams         if (stokes) {
685da80777bSKarl Rupp           IS  is_row[2];
686da80777bSKarl Rupp           Mat a[4];
687da80777bSKarl Rupp 
688da80777bSKarl Rupp           is_row[0] = kktMatsArr[level].prim_is; is_row[1] = kktMatsArr[level].constr_is;
6890298fd71SBarry Smith           a[0]      = Prol11; a[1] = NULL; a[2] = NULL; a[3] = Prol22;
6903b4367a7SBarry Smith           ierr      = MatCreateNest(comm,2,is_row, 2, is_row, a, &Parr[level1]);CHKERRQ(ierr);
6912fa5cd67SKarl Rupp         } else Parr[level1] = Prol11;
6920298fd71SBarry Smith       } else Parr[level1] = NULL;
693ffc955d6SMark F. Adams 
694ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
695806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
696ffc955d6SMark F. Adams       }
697ffc955d6SMark F. Adams 
698a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
69941b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
700a2f3521dSMark F. Adams     } /* construct prolongator scope */
7010cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7020cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
703c8b0795cSMark F. Adams #endif
7049d5b6da9SMark F. Adams     /* cache eigen estimate */
7059d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
7069d5b6da9SMark F. Adams       PetscBool flag;
707806fa848SBarry Smith       ierr = PetscObjectComposedDataGetReal((PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
7089d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
709806fa848SBarry Smith     } else emaxs[level] = -1.;
7102adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
711c8b0795cSMark F. Adams     if (!Parr[level1]) {
712806fa848SBarry Smith       if (pc_gamg->verbose) {
7133b4367a7SBarry Smith         ierr =  PetscPrintf(comm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr);
714806fa848SBarry Smith       }
715c8b0795cSMark F. Adams       break;
716c8b0795cSMark F. Adams     }
7170cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7180cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
719b4fbaa2aSMark F. Adams #endif
720a2f3521dSMark F. Adams 
721a2f3521dSMark F. Adams     ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
722806fa848SBarry Smith                        stokes, &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr);
723a2f3521dSMark F. Adams 
7240cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7250cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
726b4fbaa2aSMark F. Adams #endif
727a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
728a2f3521dSMark F. Adams 
729a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0) {
7300cbbd2e1SMark F. Adams       PetscInt NN = M;
7310cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
732a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
7333bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
734806fa848SBarry Smith       } else {
735806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
7360cbbd2e1SMark F. Adams       }
737a2f3521dSMark F. Adams 
7380cbbd2e1SMark F. Adams       nnztot += info.nz_used;
7393b4367a7SBarry Smith       ierr    = PetscPrintf(comm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
740c5df96a5SBarry Smith                             rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
741806fa848SBarry Smith                             (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
742c8b0795cSMark F. Adams     }
743a2f3521dSMark F. Adams 
744a2f3521dSMark F. Adams     /* stop if one node -- could pull back for singular problems */
745c8b0795cSMark F. Adams     if (M/pc_gamg->data_cell_cols < 2) {
74681708dccSMark F. Adams       level++;
74781708dccSMark F. Adams       break;
74881708dccSMark F. Adams     }
7490cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
750b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
751b4fbaa2aSMark F. Adams #endif
752c8b0795cSMark F. Adams   } /* levels */
753c8b0795cSMark F. Adams 
754c8b0795cSMark F. Adams   if (pc_gamg->data) {
755c8b0795cSMark F. Adams     ierr          = PetscFree(pc_gamg->data);CHKERRQ(ierr);
7560298fd71SBarry Smith     pc_gamg->data = NULL;
7575b89ad90SMark F. Adams   }
758c8b0795cSMark F. Adams 
7593b4367a7SBarry Smith   if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
7609d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7615b89ad90SMark F. Adams   fine_level       = level;
7620298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
7635b89ad90SMark F. Adams 
76484d3f75bSMark F. Adams   /* simple setup */
76584d3f75bSMark F. Adams   if (!PETSC_TRUE) {
76684d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
76784d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
76884d3f75bSMark F. Adams          lidx<fine_level;
76984d3f75bSMark F. Adams          lidx++, level--) {
77084d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
77184d3f75bSMark F. Adams       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77284d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
77384d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
77484d3f75bSMark F. Adams     }
77584d3f75bSMark F. Adams     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77684d3f75bSMark F. Adams 
77784d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
778806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
779d5280255SMark F. Adams     /* set default smoothers & set operators */
7809d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
781587fa25dSMark F. Adams          lidx <= fine_level;
782587fa25dSMark F. Adams          lidx++, level--) {
783ffc955d6SMark F. Adams       KSP smoother;
784ffc955d6SMark F. Adams       PC  subpc;
785a2f3521dSMark F. Adams 
7869d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
787f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
788ffc955d6SMark F. Adams 
789a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
790a2f3521dSMark F. Adams       /* set ops */
791a2f3521dSMark F. Adams       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
792a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
793a2f3521dSMark F. Adams 
794a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
795a2f3521dSMark F. Adams       if (stokes) {
796a2f3521dSMark F. Adams         KSP      *ksps;
797a2f3521dSMark F. Adams         PetscInt nn;
798a2f3521dSMark F. Adams         ierr     = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);CHKERRQ(ierr);
799a2f3521dSMark F. Adams         ierr     = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is);CHKERRQ(ierr);
800a2f3521dSMark F. Adams         ierr     = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
801a2f3521dSMark F. Adams         smoother = ksps[0];
802a2f3521dSMark F. Adams         ierr     = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
803a2f3521dSMark F. Adams         ierr     = PetscFree(ksps);CHKERRQ(ierr);
804a2f3521dSMark F. Adams       }
805a2f3521dSMark F. Adams       ierr = GAMGKKTMatDestroy(&kktMatsArr[level]);CHKERRQ(ierr);
806a2f3521dSMark F. Adams 
807a2f3521dSMark F. Adams       /* set defaults */
8086c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
809a2f3521dSMark F. Adams 
810ffc955d6SMark F. Adams       /* override defaults and command line args (!) */
811ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
8122d3561bbSSatish Balay         PetscInt sz;
8132d3561bbSSatish Balay         IS       *is;
814a2f3521dSMark F. Adams 
8152d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
8162d3561bbSSatish Balay         is   = ASMLocalIDsArr[level];
817ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
818ffc955d6SMark F. Adams         if (sz==0) {
819ffc955d6SMark F. Adams           IS       is;
820ffc955d6SMark F. Adams           PetscInt my0,kk;
821ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
822ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
8230298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr);
824a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
825806fa848SBarry Smith         } else {
826a94c3b12SMark F. Adams           PetscInt kk;
8270298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr);
828a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
829a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
830a94c3b12SMark F. Adams           }
831ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
832ffc955d6SMark F. Adams         }
833ffc955d6SMark F. Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
834ffc955d6SMark F. Adams 
8350298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
836ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
837ffc955d6SMark F. Adams         ierr                  = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
838806fa848SBarry Smith       } else {
8399d5b6da9SMark F. Adams         ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr);
840ffc955d6SMark F. Adams       }
841d5280255SMark F. Adams     }
842d5280255SMark F. Adams     {
843d5280255SMark F. Adams       /* coarse grid */
844d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
845d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
846d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
847d5280255SMark F. Adams       ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
848d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
849d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
850d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
851d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
85271959b99SBarry Smith       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
85371959b99SBarry Smith       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
854d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
855d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
8569dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
8572fb0b348SMark F. Adams       ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
858d5280255SMark F. Adams     }
859d5280255SMark F. Adams 
860d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
861d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
862d5280255SMark F. Adams     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr);
863d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
8643b4367a7SBarry Smith     if (mg->galerkin != 2) SETERRQ(comm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
865d5280255SMark F. Adams 
866d5280255SMark F. Adams     /* create cheby smoothers */
867d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
868d5280255SMark F. Adams          lidx <= fine_level;
869d5280255SMark F. Adams          lidx++, level--) {
870d5280255SMark F. Adams       KSP       smoother;
871ffc955d6SMark F. Adams       PetscBool flag;
872d5280255SMark F. Adams       PC        subpc;
873d5280255SMark F. Adams 
874ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
875a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
876a2f3521dSMark F. Adams 
877a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
878a2f3521dSMark F. Adams       if (stokes) {
879a2f3521dSMark F. Adams         KSP      *ksps;
880a2f3521dSMark F. Adams         PetscInt nn;
881a2f3521dSMark F. Adams         ierr     = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
882a2f3521dSMark F. Adams         smoother = ksps[0];
883a2f3521dSMark F. Adams         ierr     = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
884a2f3521dSMark F. Adams         ierr     = PetscFree(ksps);CHKERRQ(ierr);
885a2f3521dSMark F. Adams       }
886ffc955d6SMark F. Adams 
887ffc955d6SMark F. Adams       /* do my own cheby */
8886c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
889ffc955d6SMark F. Adams       if (flag) {
890ffc955d6SMark F. Adams         PetscReal emax, emin;
891251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
892ffc955d6SMark F. Adams         if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
893587fa25dSMark F. Adams         else { /* eigen estimate 'emax' */
894e696ed0bSMark F. Adams           KSP eksp;
895e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
896b2a4f308SMark F. Adams           Vec bb, xx;
897038e3b61SMark F. Adams 
8985745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr);
8995745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr);
900fc4362bfSMark F. Adams           {
901fc4362bfSMark F. Adams             PetscRandom rctx;
9023b4367a7SBarry Smith             ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
903fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
904fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
905fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
9065b89ad90SMark F. Adams           }
907a2f3521dSMark F. Adams 
908e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
909e696ed0bSMark F. Adams           {
910e696ed0bSMark F. Adams             PetscInt    Istart,Iend,ncols,jj,Ii;
911e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
912e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
913e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
914e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
915e696ed0bSMark F. Adams               if (ncols <= 1) {
916e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
917a94c3b12SMark F. Adams               }
918e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
919a94c3b12SMark F. Adams             }
920a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
921a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
922a94c3b12SMark F. Adams           }
923e696ed0bSMark F. Adams 
9243b4367a7SBarry Smith           ierr = KSPCreate(comm, &eksp);CHKERRQ(ierr);
925806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
926fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
9271a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
9281a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
929f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
930f6536408SMark F. Adams 
931f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
932f6536408SMark F. Adams           ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
933fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
9345a9b9e01SMark F. Adams 
935d3d0db20SJed Brown           /* set PC type to be same as smoother */
936ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
937b2a4f308SMark F. Adams 
9385a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
9395a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
9405a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
941fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
9425a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
9435a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
9445a9b9e01SMark F. Adams 
945fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
9465a9b9e01SMark F. Adams 
947fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
948fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
949fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
950f6536408SMark F. Adams 
951ffc955d6SMark F. Adams           if (pc_gamg->verbose > 0) {
952a94c3b12SMark F. Adams             PetscInt N1, tt;
953a94c3b12SMark F. Adams             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
9543b4367a7SBarry 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);
955f6536408SMark F. Adams           }
956fc4362bfSMark F. Adams         }
957038e3b61SMark F. Adams         {
958c5bfad50SMark F. Adams           PetscInt N1, N0;
9590298fd71SBarry Smith           ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr);
9600298fd71SBarry Smith           ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr);
961f6536408SMark F. Adams           /* heuristic - is this crap? */
962b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
9635e7c91beSJed Brown           emin  = emax * pc_gamg->eigtarget[0];
9645e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
965038e3b61SMark F. Adams         }
9666c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
967ffc955d6SMark F. Adams       } /* setup checby flag */
968ffc955d6SMark F. Adams     } /* non-coarse levels */
969737a81a9SMark F. Adams 
970d5280255SMark F. Adams     /* clean up */
971d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
972587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
973587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
9745b89ad90SMark F. Adams     }
9750cbbd2e1SMark F. Adams 
9760cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9770cbbd2e1SMark F. Adams 
9781ac9965eSMark F. Adams     if (PETSC_TRUE) {
97958471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
9809d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
98158471d46SMark F. Adams       ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
98258471d46SMark F. Adams     }
983806fa848SBarry Smith   } else {
9845f8cf99dSMark F. Adams     KSP smoother;
9853b4367a7SBarry 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__);
9869d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
9875f8cf99dSMark F. Adams     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9885f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
9899d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9905f8cf99dSMark F. Adams   }
9915b89ad90SMark F. Adams   PetscFunctionReturn(0);
9925b89ad90SMark F. Adams }
9935b89ad90SMark F. Adams 
994eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
9955b89ad90SMark F. Adams /*
9965b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
9975b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
9985b89ad90SMark F. Adams 
9995b89ad90SMark F. Adams    Input Parameter:
10005b89ad90SMark F. Adams .  pc - the preconditioner context
10015b89ad90SMark F. Adams 
10025b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
10035b89ad90SMark F. Adams */
10045b89ad90SMark F. Adams #undef __FUNCT__
10055b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
10065b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
10075b89ad90SMark F. Adams {
10085b89ad90SMark F. Adams   PetscErrorCode ierr;
10095b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
10105b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
10115b89ad90SMark F. Adams 
10125b89ad90SMark F. Adams   PetscFunctionBegin;
10135b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
1014*9b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
1015*9b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
1016*9b8ffb57SJed Brown   }
10171ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
10181ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
10195b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
10205b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
10215b89ad90SMark F. Adams   PetscFunctionReturn(0);
10225b89ad90SMark F. Adams }
10235b89ad90SMark F. Adams 
1024676e1743SMark F. Adams 
1025676e1743SMark F. Adams #undef __FUNCT__
1026676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
1027676e1743SMark F. Adams /*@
1028676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
1029676e1743SMark F. Adams    processor reduction.
1030676e1743SMark F. Adams 
1031676e1743SMark F. Adams    Not Collective on PC
1032676e1743SMark F. Adams 
1033676e1743SMark F. Adams    Input Parameters:
1034676e1743SMark F. Adams .  pc - the preconditioner context
1035676e1743SMark F. Adams 
1036676e1743SMark F. Adams 
1037676e1743SMark F. Adams    Options Database Key:
1038676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
1039676e1743SMark F. Adams 
1040676e1743SMark F. Adams    Level: intermediate
1041676e1743SMark F. Adams 
1042676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1043676e1743SMark F. Adams 
1044676e1743SMark F. Adams .seealso: ()
1045676e1743SMark F. Adams @*/
1046676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
1047676e1743SMark F. Adams {
1048676e1743SMark F. Adams   PetscErrorCode ierr;
1049676e1743SMark F. Adams 
1050676e1743SMark F. Adams   PetscFunctionBegin;
1051676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1052676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1053676e1743SMark F. Adams   PetscFunctionReturn(0);
1054676e1743SMark F. Adams }
1055676e1743SMark F. Adams 
1056676e1743SMark F. Adams #undef __FUNCT__
1057676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
10581e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1059676e1743SMark F. Adams {
1060c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1061c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1062676e1743SMark F. Adams 
1063676e1743SMark F. Adams   PetscFunctionBegin;
10649d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
1065676e1743SMark F. Adams   PetscFunctionReturn(0);
1066676e1743SMark F. Adams }
1067676e1743SMark F. Adams 
1068676e1743SMark F. Adams #undef __FUNCT__
1069389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
1070389730f3SMark F. Adams /*@
1071389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
1072389730f3SMark F. Adams 
1073389730f3SMark F. Adams  Collective on PC
1074389730f3SMark F. Adams 
1075389730f3SMark F. Adams    Input Parameters:
1076389730f3SMark F. Adams .  pc - the preconditioner context
1077389730f3SMark F. Adams 
1078389730f3SMark F. Adams 
1079389730f3SMark F. Adams    Options Database Key:
1080389730f3SMark F. Adams .  -pc_gamg_coarse_eq_limit
1081389730f3SMark F. Adams 
1082389730f3SMark F. Adams    Level: intermediate
1083389730f3SMark F. Adams 
1084389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1085389730f3SMark F. Adams 
1086389730f3SMark F. Adams .seealso: ()
1087389730f3SMark F. Adams  @*/
1088389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
1089389730f3SMark F. Adams {
1090389730f3SMark F. Adams   PetscErrorCode ierr;
1091389730f3SMark F. Adams 
1092389730f3SMark F. Adams   PetscFunctionBegin;
1093389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1094389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1095389730f3SMark F. Adams   PetscFunctionReturn(0);
1096389730f3SMark F. Adams }
1097389730f3SMark F. Adams 
1098389730f3SMark F. Adams #undef __FUNCT__
1099389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
11001e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1101389730f3SMark F. Adams {
1102389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1103389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1104389730f3SMark F. Adams 
1105389730f3SMark F. Adams   PetscFunctionBegin;
11069d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
1107389730f3SMark F. Adams   PetscFunctionReturn(0);
1108389730f3SMark F. Adams }
1109389730f3SMark F. Adams 
1110389730f3SMark F. Adams #undef __FUNCT__
11118263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1112676e1743SMark F. Adams /*@
11138263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1114676e1743SMark F. Adams 
1115676e1743SMark F. Adams    Collective on PC
1116676e1743SMark F. Adams 
1117676e1743SMark F. Adams    Input Parameters:
1118676e1743SMark F. Adams .  pc - the preconditioner context
1119676e1743SMark F. Adams 
1120676e1743SMark F. Adams 
1121676e1743SMark F. Adams    Options Database Key:
11228263b398SMark F. Adams .  -pc_gamg_repartition
1123676e1743SMark F. Adams 
1124676e1743SMark F. Adams    Level: intermediate
1125676e1743SMark F. Adams 
1126676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1127676e1743SMark F. Adams 
1128676e1743SMark F. Adams .seealso: ()
1129676e1743SMark F. Adams @*/
11308263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1131676e1743SMark F. Adams {
1132676e1743SMark F. Adams   PetscErrorCode ierr;
1133676e1743SMark F. Adams 
1134676e1743SMark F. Adams   PetscFunctionBegin;
1135676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11368263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1137676e1743SMark F. Adams   PetscFunctionReturn(0);
1138676e1743SMark F. Adams }
1139676e1743SMark F. Adams 
1140676e1743SMark F. Adams #undef __FUNCT__
11418263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
11421e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1143676e1743SMark F. Adams {
1144c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1145c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1146676e1743SMark F. Adams 
1147676e1743SMark F. Adams   PetscFunctionBegin;
11489d5b6da9SMark F. Adams   pc_gamg->repart = n;
1149676e1743SMark F. Adams   PetscFunctionReturn(0);
1150676e1743SMark F. Adams }
1151676e1743SMark F. Adams 
1152676e1743SMark F. Adams #undef __FUNCT__
1153dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl"
1154dfd5c07aSMark F. Adams /*@
1155dfd5c07aSMark F. Adams    PCGAMGSetReuseProl - Reuse prlongation
1156dfd5c07aSMark F. Adams 
1157dfd5c07aSMark F. Adams    Collective on PC
1158dfd5c07aSMark F. Adams 
1159dfd5c07aSMark F. Adams    Input Parameters:
1160dfd5c07aSMark F. Adams .  pc - the preconditioner context
1161dfd5c07aSMark F. Adams 
1162dfd5c07aSMark F. Adams 
1163dfd5c07aSMark F. Adams    Options Database Key:
1164dfd5c07aSMark F. Adams .  -pc_gamg_reuse_interpolation
1165dfd5c07aSMark F. Adams 
1166dfd5c07aSMark F. Adams    Level: intermediate
1167dfd5c07aSMark F. Adams 
1168dfd5c07aSMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1169dfd5c07aSMark F. Adams 
1170dfd5c07aSMark F. Adams .seealso: ()
1171dfd5c07aSMark F. Adams @*/
1172dfd5c07aSMark F. Adams PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n)
1173dfd5c07aSMark F. Adams {
1174dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1175dfd5c07aSMark F. Adams 
1176dfd5c07aSMark F. Adams   PetscFunctionBegin;
1177dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1178dfd5c07aSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1179dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1180dfd5c07aSMark F. Adams }
1181dfd5c07aSMark F. Adams 
1182dfd5c07aSMark F. Adams #undef __FUNCT__
1183dfd5c07aSMark F. Adams #define __FUNCT__ "PCGAMGSetReuseProl_GAMG"
11841e6b0712SBarry Smith static PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n)
1185dfd5c07aSMark F. Adams {
1186dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1187dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1188dfd5c07aSMark F. Adams 
1189dfd5c07aSMark F. Adams   PetscFunctionBegin;
1190dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1191dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1192dfd5c07aSMark F. Adams }
1193dfd5c07aSMark F. Adams 
1194dfd5c07aSMark F. Adams #undef __FUNCT__
1195ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1196ffc955d6SMark F. Adams /*@
1197ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1198ffc955d6SMark F. Adams 
1199ffc955d6SMark F. Adams    Collective on PC
1200ffc955d6SMark F. Adams 
1201ffc955d6SMark F. Adams    Input Parameters:
1202ffc955d6SMark F. Adams .  pc - the preconditioner context
1203ffc955d6SMark F. Adams 
1204ffc955d6SMark F. Adams 
1205ffc955d6SMark F. Adams    Options Database Key:
1206ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1207ffc955d6SMark F. Adams 
1208ffc955d6SMark F. Adams    Level: intermediate
1209ffc955d6SMark F. Adams 
1210ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1211ffc955d6SMark F. Adams 
1212ffc955d6SMark F. Adams .seealso: ()
1213ffc955d6SMark F. Adams @*/
1214ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1215ffc955d6SMark F. Adams {
1216ffc955d6SMark F. Adams   PetscErrorCode ierr;
1217ffc955d6SMark F. Adams 
1218ffc955d6SMark F. Adams   PetscFunctionBegin;
1219ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1220ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1221ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1222ffc955d6SMark F. Adams }
1223ffc955d6SMark F. Adams 
1224ffc955d6SMark F. Adams #undef __FUNCT__
1225ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
12261e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1227ffc955d6SMark F. Adams {
1228ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1229ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1230ffc955d6SMark F. Adams 
1231ffc955d6SMark F. Adams   PetscFunctionBegin;
1232ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1233ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1234ffc955d6SMark F. Adams }
1235ffc955d6SMark F. Adams 
1236ffc955d6SMark F. Adams #undef __FUNCT__
12374ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
12384ef23d27SMark F. Adams /*@
12394ef23d27SMark F. Adams    PCGAMGSetNlevels -
12404ef23d27SMark F. Adams 
12414ef23d27SMark F. Adams    Not collective on PC
12424ef23d27SMark F. Adams 
12434ef23d27SMark F. Adams    Input Parameters:
12444ef23d27SMark F. Adams .  pc - the preconditioner context
12454ef23d27SMark F. Adams 
12464ef23d27SMark F. Adams 
12474ef23d27SMark F. Adams    Options Database Key:
12484ef23d27SMark F. Adams .  -pc_mg_levels
12494ef23d27SMark F. Adams 
12504ef23d27SMark F. Adams    Level: intermediate
12514ef23d27SMark F. Adams 
12524ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12534ef23d27SMark F. Adams 
12544ef23d27SMark F. Adams .seealso: ()
12554ef23d27SMark F. Adams @*/
12564ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12574ef23d27SMark F. Adams {
12584ef23d27SMark F. Adams   PetscErrorCode ierr;
12594ef23d27SMark F. Adams 
12604ef23d27SMark F. Adams   PetscFunctionBegin;
12614ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12624ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
12634ef23d27SMark F. Adams   PetscFunctionReturn(0);
12644ef23d27SMark F. Adams }
12654ef23d27SMark F. Adams 
12664ef23d27SMark F. Adams #undef __FUNCT__
12674ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
12681e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12694ef23d27SMark F. Adams {
12704ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
12714ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12724ef23d27SMark F. Adams 
12734ef23d27SMark F. Adams   PetscFunctionBegin;
12749d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12754ef23d27SMark F. Adams   PetscFunctionReturn(0);
12764ef23d27SMark F. Adams }
12774ef23d27SMark F. Adams 
12784ef23d27SMark F. Adams #undef __FUNCT__
12793542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
12803542efc5SMark F. Adams /*@
12813542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12823542efc5SMark F. Adams 
12833542efc5SMark F. Adams    Not collective on PC
12843542efc5SMark F. Adams 
12853542efc5SMark F. Adams    Input Parameters:
12863542efc5SMark F. Adams .  pc - the preconditioner context
12873542efc5SMark F. Adams 
12883542efc5SMark F. Adams 
12893542efc5SMark F. Adams    Options Database Key:
12903542efc5SMark F. Adams .  -pc_gamg_threshold
12913542efc5SMark F. Adams 
12923542efc5SMark F. Adams    Level: intermediate
12933542efc5SMark F. Adams 
12943542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12953542efc5SMark F. Adams 
12963542efc5SMark F. Adams .seealso: ()
12973542efc5SMark F. Adams @*/
12983542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
12993542efc5SMark F. Adams {
13003542efc5SMark F. Adams   PetscErrorCode ierr;
13013542efc5SMark F. Adams 
13023542efc5SMark F. Adams   PetscFunctionBegin;
13033542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13043542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
13053542efc5SMark F. Adams   PetscFunctionReturn(0);
13063542efc5SMark F. Adams }
13073542efc5SMark F. Adams 
13083542efc5SMark F. Adams #undef __FUNCT__
13093542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
13101e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
13113542efc5SMark F. Adams {
1312c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1313c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
13143542efc5SMark F. Adams 
13153542efc5SMark F. Adams   PetscFunctionBegin;
13169d5b6da9SMark F. Adams   pc_gamg->threshold = n;
13173542efc5SMark F. Adams   PetscFunctionReturn(0);
13183542efc5SMark F. Adams }
13193542efc5SMark F. Adams 
13203542efc5SMark F. Adams #undef __FUNCT__
13219d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1322676e1743SMark F. Adams /*@
13239d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1324676e1743SMark F. Adams 
1325676e1743SMark F. Adams    Collective on PC
1326676e1743SMark F. Adams 
1327676e1743SMark F. Adams    Input Parameters:
1328676e1743SMark F. Adams .  pc - the preconditioner context
1329676e1743SMark F. Adams 
1330676e1743SMark F. Adams 
1331676e1743SMark F. Adams    Options Database Key:
13323542efc5SMark F. Adams .  -pc_gamg_type
1333676e1743SMark F. Adams 
1334676e1743SMark F. Adams    Level: intermediate
1335676e1743SMark F. Adams 
1336676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1337676e1743SMark F. Adams 
1338676e1743SMark F. Adams .seealso: ()
1339676e1743SMark F. Adams @*/
134019fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1341676e1743SMark F. Adams {
1342676e1743SMark F. Adams   PetscErrorCode ierr;
1343676e1743SMark F. Adams 
1344676e1743SMark F. Adams   PetscFunctionBegin;
1345676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1346806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1347676e1743SMark F. Adams   PetscFunctionReturn(0);
1348676e1743SMark F. Adams }
1349676e1743SMark F. Adams 
1350676e1743SMark F. Adams #undef __FUNCT__
13519d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
13521e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1353676e1743SMark F. Adams {
13549d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
13551ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
13561ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1357676e1743SMark F. Adams 
1358676e1743SMark F. Adams   PetscFunctionBegin;
13591c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
13609d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
13611ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
13621ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
13631ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
13641ab5ffc9SJed Brown   }
13651ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
13661ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
13679d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1368676e1743SMark F. Adams   PetscFunctionReturn(0);
1369676e1743SMark F. Adams }
1370676e1743SMark F. Adams 
13715b89ad90SMark F. Adams #undef __FUNCT__
13725b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
13735b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
13745b89ad90SMark F. Adams {
1375676e1743SMark F. Adams   PetscErrorCode ierr;
1376676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1377676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1378676e1743SMark F. Adams   PetscBool      flag;
13795e7c91beSJed Brown   PetscInt       two   = 2;
13803b4367a7SBarry Smith   MPI_Comm       comm;
13815b89ad90SMark F. Adams 
13825b89ad90SMark F. Adams   PetscFunctionBegin;
13833b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1384676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
1385676e1743SMark F. Adams   {
1386b7cbab4eSMark Adams     /* -pc_gamg_type */
1387b7cbab4eSMark Adams     {
1388b7cbab4eSMark Adams       char tname[256] = PCGAMGAGG;
13891ab5ffc9SJed Brown       const char *deftype = pc_gamg->gamg_type_name ? pc_gamg->gamg_type_name : tname;
13901ab5ffc9SJed Brown       ierr = PetscOptionsList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1391b7cbab4eSMark Adams       /* call PCCreateGAMG_XYZ */
13921ab5ffc9SJed Brown       if (flag || !pc_gamg->gamg_type_name) {
13931ab5ffc9SJed Brown         ierr = PCGAMGSetType(pc, flag ? tname : deftype);CHKERRQ(ierr);
13941ab5ffc9SJed Brown       }
1395b7cbab4eSMark Adams     }
139675b74bdaSMark F. Adams     /* -pc_gamg_verbose */
13979d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
13989d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
13990298fd71SBarry Smith                            &pc_gamg->verbose, NULL);CHKERRQ(ierr);
14008263b398SMark F. Adams     /* -pc_gamg_repartition */
14018263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
14028263b398SMark F. Adams                             "Repartion coarse grids (false)",
14038263b398SMark F. Adams                             "PCGAMGRepartitioning",
14049d5b6da9SMark F. Adams                             pc_gamg->repart,
14059d5b6da9SMark F. Adams                             &pc_gamg->repart,
1406806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1407dfd5c07aSMark F. Adams     /* -pc_gamg_reuse_interpolation */
1408dfd5c07aSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation",
1409dfd5c07aSMark F. Adams                             "Reuse prolongation operator (true)",
1410dfd5c07aSMark F. Adams                             "PCGAMGReuseProl",
1411dfd5c07aSMark F. Adams                             pc_gamg->reuse_prol,
1412dfd5c07aSMark F. Adams                             &pc_gamg->reuse_prol,
1413dfd5c07aSMark F. Adams                             &flag);CHKERRQ(ierr);
1414ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
1415ffc955d6SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1416ffc955d6SMark F. Adams                             "Use aggregation agragates for GASM smoother (false)",
1417ffc955d6SMark F. Adams                             "PCGAMGUseASMAggs",
1418ffc955d6SMark F. Adams                             pc_gamg->use_aggs_in_gasm,
1419ffc955d6SMark F. Adams                             &pc_gamg->use_aggs_in_gasm,
1420806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1421c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1422676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1423676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1424676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
14259d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
14269d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1427806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1428389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1429389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1430389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1431389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
14329d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
14339d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1434806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1435c20e4228SMark F. Adams     /* -pc_gamg_threshold */
14363542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
14373542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
14383542efc5SMark F. Adams                             "PCGAMGSetThreshold",
14399d5b6da9SMark F. Adams                             pc_gamg->threshold,
14409d5b6da9SMark F. Adams                             &pc_gamg->threshold,
1441806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1442806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
14433b4367a7SBarry Smith       ierr = PetscPrintf(comm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1444806fa848SBarry Smith     }
1445b7cbab4eSMark Adams     /* -pc_gamg_eigtarget */
14460298fd71SBarry Smith     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr);
14474ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
14484ef23d27SMark F. Adams                            "Set number of MG levels",
14494ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
14509d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
14519d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
1452806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1453b7cbab4eSMark Adams 
1454b7cbab4eSMark Adams     /* set options for subtype */
14551ab5ffc9SJed Brown     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(pc);CHKERRQ(ierr);}
1456676e1743SMark F. Adams   }
1457676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
14585b89ad90SMark F. Adams   PetscFunctionReturn(0);
14595b89ad90SMark F. Adams }
14605b89ad90SMark F. Adams 
14615b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
14625b89ad90SMark F. Adams /*MC
1463856830bfSJed Brown      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
14649d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
14655b89ad90SMark F. Adams 
1466280d9858SJed Brown    Options Database Keys:
14675b89ad90SMark F. Adams    Multigrid options(inherited)
1468280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1469280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1470280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
14718c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
14725b89ad90SMark F. Adams 
14735b89ad90SMark F. Adams   Level: intermediate
1474280d9858SJed Brown 
14755b89ad90SMark F. Adams   Concepts: multigrid
14765b89ad90SMark F. Adams 
14775b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1478280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
14795b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
14805b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
14815b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
14825b89ad90SMark F. Adams M*/
1483b2573a8aSBarry Smith 
14845b89ad90SMark F. Adams #undef __FUNCT__
14855b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
14868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
14875b89ad90SMark F. Adams {
14885b89ad90SMark F. Adams   PetscErrorCode ierr;
14895b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
14905b89ad90SMark F. Adams   PC_MG          *mg;
14910cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14929d5b6da9SMark F. Adams   static long count = 0;
14935ee9c036SSatish Balay #endif
14945b89ad90SMark F. Adams 
14955b89ad90SMark F. Adams   PetscFunctionBegin;
14965b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14975b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
14985b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14995b89ad90SMark F. Adams 
15005b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
15015b89ad90SMark F. Adams   ierr         = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr);
15025b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
1503ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
15045b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
15055b89ad90SMark F. Adams 
15061ab5ffc9SJed Brown   ierr = PetscNewLog(pc,struct _PCGAMGOps,&pc_gamg->ops);CHKERRQ(ierr);
15071ab5ffc9SJed Brown 
15089d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
15099d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
15109d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
15119d5b6da9SMark F. Adams   pc_gamg->data    = 0;
15125b89ad90SMark F. Adams 
15139d5b6da9SMark F. Adams   /* register AMG type */
15143e3471ccSMark Adams #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
15153e3471ccSMark Adams   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
15163e3471ccSMark Adams #endif
15179d5b6da9SMark F. Adams 
15189d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
15195b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
15205b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
15215b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
15225b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
15235b89ad90SMark F. Adams 
1524bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1525bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1526bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
1527bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseProl_C",PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr);
1528bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1529bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1530bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1531bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
15329d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1533dfd5c07aSMark F. Adams   pc_gamg->reuse_prol       = PETSC_TRUE;
1534ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1535038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
1536c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit  = 800;
1537a2f3521dSMark F. Adams   pc_gamg->threshold        = 0.001;
15389d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
15399d5b6da9SMark F. Adams   pc_gamg->verbose          = 0;
15409d5b6da9SMark F. Adams   pc_gamg->emax_id          = -1;
15415e7c91beSJed Brown   pc_gamg->eigtarget[0]     = 0.05;
15425e7c91beSJed Brown   pc_gamg->eigtarget[1]     = 1.05;
15439d5b6da9SMark F. Adams 
15440cbbd2e1SMark F. Adams   /* private events */
15450cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1546785cba28SMark F. Adams   if (count++ == 0) {
1547806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1548806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
15490cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
15500cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
15510cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1552806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1553806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1554806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1555806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1556806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1557806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1558806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1559806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1560806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1561806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1562806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1563806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1564f852f58cSMark F. Adams 
15650cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
15660cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
15670cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1568b4fbaa2aSMark F. Adams     /* create timer stages */
1569b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1570b4fbaa2aSMark F. Adams     {
1571b4fbaa2aSMark F. Adams       char     str[32];
1572b4fbaa2aSMark F. Adams       PetscInt lidx;
1573806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1574806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1575b4fbaa2aSMark F. Adams       for (lidx=1; lidx<9; lidx++) {
1576b4fbaa2aSMark F. Adams         sprintf(str,"MG Level %d",lidx);
1577806fa848SBarry Smith         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1578b4fbaa2aSMark F. Adams       }
1579b4fbaa2aSMark F. Adams     }
1580b4fbaa2aSMark F. Adams #endif
1581b4fbaa2aSMark F. Adams   }
1582b4fbaa2aSMark F. Adams #endif
15830cbbd2e1SMark F. Adams   /* general events */
15840cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
1585806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1586806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1587806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1588806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1589806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1590806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1591806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1592806fa848SBarry Smith   ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr);
15930cbbd2e1SMark F. Adams #endif
15940cbbd2e1SMark F. Adams 
15955b89ad90SMark F. Adams   PetscFunctionReturn(0);
15965b89ad90SMark F. Adams }
15973e3471ccSMark Adams 
15983e3471ccSMark Adams #undef __FUNCT__
15993e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage"
16003e3471ccSMark Adams /*@C
16013e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
16023e3471ccSMark Adams  from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
16033e3471ccSMark Adams  when using static libraries.
16043e3471ccSMark Adams 
16053e3471ccSMark Adams  Level: developer
16063e3471ccSMark Adams 
16073e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
16083e3471ccSMark Adams  .seealso: PetscInitialize()
16093e3471ccSMark Adams @*/
16103e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
16113e3471ccSMark Adams {
16123e3471ccSMark Adams   PetscErrorCode ierr;
16133e3471ccSMark Adams 
16143e3471ccSMark Adams   PetscFunctionBegin;
16153e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
16163e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
16173e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
16183e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
16193e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
16203e3471ccSMark Adams   PetscFunctionReturn(0);
16213e3471ccSMark Adams }
16223e3471ccSMark Adams 
16233e3471ccSMark Adams #undef __FUNCT__
16243e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage"
16253e3471ccSMark Adams /*@C
16263e3471ccSMark Adams  PCGAMGFinalizePackage - This function destroys everything in the PCGAMG package. It is
16273e3471ccSMark Adams  called from PetscFinalize().
16283e3471ccSMark Adams 
16293e3471ccSMark Adams  Level: developer
16303e3471ccSMark Adams 
16313e3471ccSMark Adams  .keywords: Petsc, destroy, package
16323e3471ccSMark Adams  .seealso: PetscFinalize()
16333e3471ccSMark Adams @*/
16343e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
16353e3471ccSMark Adams {
16363e3471ccSMark Adams   PetscErrorCode ierr;
16373e3471ccSMark Adams 
16383e3471ccSMark Adams   PetscFunctionBegin;
16393e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
16403e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
16413e3471ccSMark Adams   PetscFunctionReturn(0);
16423e3471ccSMark Adams }
1643