xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 3bf036e263fd78807e2931ff42d9ddcd8aae3fd4)
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 
309d5b6da9SMark F. Adams static PetscFList GAMGList = 0;
319d5b6da9SMark F. Adams 
32d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
33d3d6bff4SMark F. Adams #undef __FUNCT__
34d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
35d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
36d3d6bff4SMark F. Adams {
37d3d6bff4SMark F. Adams   PetscErrorCode  ierr;
38d3d6bff4SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
39d3d6bff4SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
40d3d6bff4SMark F. Adams 
41d3d6bff4SMark F. Adams   PetscFunctionBegin;
42a2f3521dSMark F. Adams   if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */
43878e152fSMark F. Adams     PetscPrintf(((PetscObject)pc)->comm,"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__);
449d5b6da9SMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
4558471d46SMark F. Adams   }
46a2f3521dSMark F. Adams   pc_gamg->data = PETSC_NULL; pc_gamg->data_sz = 0;
47878e152fSMark F. Adams 
48878e152fSMark F. Adams   if (pc_gamg->orig_data) {
49878e152fSMark F. Adams     ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
50878e152fSMark F. Adams   }
51878e152fSMark F. Adams 
52a2f3521dSMark F. Adams   PetscFunctionReturn(0);
53a2f3521dSMark F. Adams }
54a2f3521dSMark F. Adams 
55a2f3521dSMark F. Adams /* private 2x2 Mat Nest for Stokes */
56a2f3521dSMark F. Adams typedef struct{
57a2f3521dSMark F. Adams   Mat A11,A21,A12,Amat;
58a2f3521dSMark F. Adams   IS  prim_is,constr_is;
59a2f3521dSMark F. Adams }GAMGKKTMat;
60a2f3521dSMark F. Adams 
61a2f3521dSMark F. Adams #undef __FUNCT__
62a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatCreate"
63a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatCreate(Mat A, PetscBool iskkt, GAMGKKTMat *out)
64a2f3521dSMark F. Adams {
65a2f3521dSMark F. Adams   PetscFunctionBegin;
66a2f3521dSMark F. Adams   out->Amat = A;
67a2f3521dSMark F. Adams   if (iskkt){
68a2f3521dSMark F. Adams     PetscErrorCode   ierr;
69a2f3521dSMark F. Adams     IS       is_constraint, is_prime;
70a2f3521dSMark F. Adams     PetscInt nmin,nmax;
71a2f3521dSMark F. Adams 
72a2f3521dSMark F. Adams     ierr = MatGetOwnershipRange(A, &nmin, &nmax);CHKERRQ(ierr);
73a2f3521dSMark F. Adams     ierr = MatFindZeroDiagonals(A, &is_constraint);CHKERRQ(ierr);
74a2f3521dSMark F. Adams     ierr = ISComplement(is_constraint, nmin, nmax, &is_prime);CHKERRQ(ierr);
75a2f3521dSMark F. Adams     out->prim_is = is_prime;
76a2f3521dSMark F. Adams     out->constr_is = is_constraint;
77a2f3521dSMark F. Adams 
78a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_prime, is_prime,      MAT_INITIAL_MATRIX, &out->A11);CHKERRQ(ierr);
79a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12);CHKERRQ(ierr);
80a2f3521dSMark F. Adams     ierr = MatGetSubMatrix(A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21);CHKERRQ(ierr);
81806fa848SBarry Smith   } else {
82a2f3521dSMark F. Adams     out->A11 = A;
83a2f3521dSMark F. Adams     out->A21 = PETSC_NULL;
84a2f3521dSMark F. Adams     out->A12 = PETSC_NULL;
85a2f3521dSMark F. Adams     out->prim_is = PETSC_NULL;
86a2f3521dSMark F. Adams     out->constr_is = PETSC_NULL;
87a2f3521dSMark F. Adams   }
88a2f3521dSMark F. Adams   PetscFunctionReturn(0);
89a2f3521dSMark F. Adams }
90a2f3521dSMark F. Adams 
91a2f3521dSMark F. Adams #undef __FUNCT__
92a2f3521dSMark F. Adams #define __FUNCT__ "GAMGKKTMatDestroy"
93a2f3521dSMark F. Adams static PetscErrorCode GAMGKKTMatDestroy(GAMGKKTMat *mat)
94a2f3521dSMark F. Adams {
95a2f3521dSMark F. Adams   PetscErrorCode   ierr;
96a2f3521dSMark F. Adams 
97a2f3521dSMark F. Adams   PetscFunctionBegin;
98a2f3521dSMark F. Adams   if (mat->A11 && mat->A11 != mat->Amat) {
99a2f3521dSMark F. Adams     ierr = MatDestroy(&mat->A11);CHKERRQ(ierr);
100a2f3521dSMark F. Adams   }
101a2f3521dSMark F. Adams   ierr = MatDestroy(&mat->A21);CHKERRQ(ierr);
102a2f3521dSMark F. Adams   ierr = MatDestroy(&mat->A12);CHKERRQ(ierr);
103a2f3521dSMark F. Adams 
104a2f3521dSMark F. Adams   ierr = ISDestroy(&mat->prim_is);CHKERRQ(ierr);
105a2f3521dSMark F. Adams   ierr = ISDestroy(&mat->constr_is);CHKERRQ(ierr);
106a2f3521dSMark F. Adams 
107d3d6bff4SMark F. Adams   PetscFunctionReturn(0);
108d3d6bff4SMark F. Adams }
109d3d6bff4SMark F. Adams 
1105b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1115b89ad90SMark F. Adams /*
112a147abb0SMark F. Adams    createLevel: create coarse op with RAP.  repartition and/or reduce number
113a147abb0SMark F. Adams      of active processors.
1145b89ad90SMark F. Adams 
1155b89ad90SMark F. Adams    Input Parameter:
116a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
117a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
1189d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
119c5bfad50SMark F. Adams    . cr_bs - coarse block size
1209d5b6da9SMark F. Adams    . isLast -
121a2f3521dSMark F. Adams    . stokes -
1223530afc2SMark F. Adams    In/Output Parameter:
123a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
124afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
12511e60469SMark F. Adams    Output Parameter:
1263530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
1275b89ad90SMark F. Adams */
1285cb416c2SMark F. Adams 
1295b89ad90SMark F. Adams #undef __FUNCT__
1308263b398SMark F. Adams #define __FUNCT__ "createLevel"
1310cbbd2e1SMark F. Adams static PetscErrorCode createLevel(const PC pc,
1329d5b6da9SMark F. Adams                                    const Mat Amat_fine,
133c5bfad50SMark F. Adams                                    const PetscInt cr_bs,
1349d5b6da9SMark F. Adams                                    const PetscBool isLast,
135a2f3521dSMark F. Adams                                    const PetscBool stokes,
1363530afc2SMark F. Adams                                    Mat *a_P_inout,
137a2f3521dSMark F. Adams                                    Mat *a_Amat_crs,
138a2f3521dSMark F. Adams                                    PetscMPIInt *a_nactive_proc
13911e60469SMark F. Adams                                   )
1405b89ad90SMark F. Adams {
141a2f3521dSMark F. Adams   PetscErrorCode   ierr;
1429d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
143486a8d0bSJed Brown   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1449d5b6da9SMark F. Adams   const PetscBool  repart = pc_gamg->repart;
1459d5b6da9SMark F. Adams   const PetscInt   min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
146a2f3521dSMark F. Adams   Mat              Cmat,Pold=*a_P_inout;
1479d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)Amat_fine)->comm;
1485a9b9e01SMark F. Adams   PetscMPIInt      mype,npe,new_npe,nactive=*a_nactive_proc;
149c5bfad50SMark F. Adams   PetscInt         ncrs_eq,ncrs_prim,f_bs;
1505b89ad90SMark F. Adams 
1515b89ad90SMark F. Adams   PetscFunctionBegin;
15211e60469SMark F. Adams   ierr = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr);
15311e60469SMark F. Adams   ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr);
154c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
15511e60469SMark F. Adams   /* RAP */
1569d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
157038e3b61SMark F. Adams 
158a2f3521dSMark F. Adams   /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/
159a2f3521dSMark F. Adams   ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
160a2f3521dSMark F. Adams   assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)==0);
161a2f3521dSMark F. Adams   ierr = MatGetLocalSize(Cmat, &ncrs_eq, PETSC_NULL);CHKERRQ(ierr);
162a2f3521dSMark F. Adams 
163a2f3521dSMark F. Adams   /* get number of PEs to make active 'new_npe', reduce, can be any integer 1-P */
164a2f3521dSMark F. Adams   {
165472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
166a2f3521dSMark F. Adams     ierr = MatGetSize(Cmat, &ncrs_eq_glob, PETSC_NULL);CHKERRQ(ierr);
167472110cdSMark F. Adams     new_npe = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
168472110cdSMark F. Adams     if (new_npe == 0 || ncrs_eq_glob < coarse_max) new_npe = 1;
1695a9b9e01SMark F. Adams     else if (new_npe >= nactive) new_npe = nactive; /* no change, rare */
1709d5b6da9SMark F. Adams     if (isLast) new_npe = 1;
171a2f3521dSMark F. Adams   }
172f852f58cSMark F. Adams 
1735a9b9e01SMark F. Adams   if (!repart && new_npe==nactive) {
174a2f3521dSMark F. Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
175806fa848SBarry Smith   } else {
176a2f3521dSMark 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;
177a2f3521dSMark F. Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
178a2f3521dSMark F. Adams     IS              is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
1795a9b9e01SMark F. Adams     VecScatter      vecscat;
18022063be5SMark F. Adams     PetscScalar    *array;
18122063be5SMark F. Adams     Vec             src_crd, dest_crd;
182e33ef3b1SMark F. Adams 
183c5bfad50SMark F. Adams     nloc_old = ncrs_eq/cr_bs; assert(ncrs_eq%cr_bs==0);
1840cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1850cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
186b4fbaa2aSMark F. Adams #endif
187a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
188a2f3521dSMark F. Adams     ierr = PetscMalloc(npe*sizeof(PetscInt), &counts);CHKERRQ(ierr);
189a2f3521dSMark F. Adams     if (repart && !stokes) {
190a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1915a9b9e01SMark F. Adams       Mat             adj;
1925a9b9e01SMark F. Adams 
193a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
194a2f3521dSMark F. Adams         if (pc_gamg->verbose==1) PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,ncrs_eq);
195a2f3521dSMark F. Adams         else {
196a2f3521dSMark F. Adams           PetscInt n;
197a2f3521dSMark F. Adams           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm);CHKERRQ(ierr);
198a2f3521dSMark F. Adams           PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,n);
199a2f3521dSMark F. Adams         }
200a2f3521dSMark F. Adams       }
2015a9b9e01SMark F. Adams 
202a2f3521dSMark F. Adams       /* get 'adj' */
203c5bfad50SMark F. Adams       if (cr_bs == 1) {
204038e3b61SMark F. Adams 	ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
205806fa848SBarry Smith       } else{
206a2f3521dSMark F. Adams 	/* make a scalar matrix to partition (no Stokes here) */
207eb07cef2SMark F. Adams 	Mat tMat;
208a2f3521dSMark F. Adams 	PetscInt Istart_crs,Iend_crs,ncols,jj,Ii;
209b4fbaa2aSMark F. Adams 	const PetscScalar *vals;
210b4fbaa2aSMark F. Adams 	const PetscInt *idx;
211a2f3521dSMark F. Adams 	PetscInt *d_nnz, *o_nnz, M, N;
2129057884aSMark F. Adams 	static PetscInt llev = 0;
213b4fbaa2aSMark F. Adams 
214a2f3521dSMark F. Adams 	ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &d_nnz);CHKERRQ(ierr);
215a2f3521dSMark F. Adams 	ierr = PetscMalloc(ncrs_prim*sizeof(PetscInt), &o_nnz);CHKERRQ(ierr);
216a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
217a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
218c5bfad50SMark F. Adams 	for (Ii = Istart_crs, jj = 0 ; Ii < Iend_crs ; Ii += cr_bs, jj++) {
21958471d46SMark F. Adams 	  ierr = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
220c5bfad50SMark F. Adams 	  d_nnz[jj] = ncols/cr_bs;
221c5bfad50SMark F. Adams 	  o_nnz[jj] = ncols/cr_bs;
22258471d46SMark F. Adams 	  ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
223a2f3521dSMark F. Adams 	  if (d_nnz[jj] > ncrs_prim) d_nnz[jj] = ncrs_prim;
224c5bfad50SMark F. Adams 	  if (o_nnz[jj] > (M/cr_bs-ncrs_prim)) o_nnz[jj] = M/cr_bs-ncrs_prim;
22558471d46SMark F. Adams 	}
2266876a03eSMark F. Adams 
227a2f3521dSMark F. Adams 	ierr = MatCreate(wcomm, &tMat);CHKERRQ(ierr);
228806fa848SBarry Smith 	ierr = MatSetSizes(tMat, ncrs_prim, ncrs_prim,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
229a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
230a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
231a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
23258471d46SMark F. Adams 	ierr = PetscFree(d_nnz);CHKERRQ(ierr);
2335f8cf99dSMark F. Adams 	ierr = PetscFree(o_nnz);CHKERRQ(ierr);
234eb07cef2SMark F. Adams 
235a2f3521dSMark F. Adams 	for (ii = Istart_crs; ii < Iend_crs; ii++) {
236c5bfad50SMark F. Adams 	  PetscInt dest_row = ii/cr_bs;
23722063be5SMark F. Adams 	  ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
238eb07cef2SMark F. Adams 	  for (jj = 0 ; jj < ncols ; jj++){
239c5bfad50SMark F. Adams 	    PetscInt dest_col = idx[jj]/cr_bs;
240eb07cef2SMark F. Adams 	    PetscScalar v = 1.0;
241eb07cef2SMark F. Adams 	    ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
242eb07cef2SMark F. Adams 	  }
24322063be5SMark F. Adams 	  ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
244eb07cef2SMark F. Adams 	}
245eb07cef2SMark F. Adams 	ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
246eb07cef2SMark F. Adams 	ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
247eb07cef2SMark F. Adams 
248b4fbaa2aSMark F. Adams 	if (llev++ == -1) {
249b4fbaa2aSMark F. Adams 	  PetscViewer viewer; char fname[32];
2508caf3d72SBarry Smith 	  ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
251b4fbaa2aSMark F. Adams 	  PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
252b4fbaa2aSMark F. Adams 	  ierr = MatView(tMat, viewer);CHKERRQ(ierr);
253*3bf036e2SBarry Smith 	  ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
254b4fbaa2aSMark F. Adams 	}
255b4fbaa2aSMark F. Adams 
256eb07cef2SMark F. Adams 	ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
257eb07cef2SMark F. Adams 
258eb07cef2SMark F. Adams 	ierr = MatDestroy(&tMat);CHKERRQ(ierr);
259a2f3521dSMark F. Adams       } /* create 'adj' */
260f150b916SMark F. Adams 
261a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2625a9b9e01SMark F. Adams 	char prefix[256];
2635a9b9e01SMark F. Adams 	const char *pcpre;
264b4fbaa2aSMark F. Adams 	const PetscInt *is_idx;
265b4fbaa2aSMark F. Adams 	MatPartitioning  mpart;
266a4b7d37bSMark F. Adams 	IS proc_is;
267a2f3521dSMark F. Adams 	PetscInt targetPE;
2682f03bc48SMark F. Adams 
2695a9b9e01SMark F. Adams 	ierr = MatPartitioningCreate(wcomm, &mpart);CHKERRQ(ierr);
2705ef31b24SMark F. Adams 	ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2719d5b6da9SMark F. Adams 	ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2728caf3d72SBarry Smith 	ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
27359a0be82SJed Brown 	ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
27411e60469SMark F. Adams 	ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
2755ef31b24SMark F. Adams 	ierr = MatPartitioningSetNParts(mpart, new_npe);CHKERRQ(ierr);
276a4b7d37bSMark F. Adams 	ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
27711e60469SMark F. Adams 	ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2785a9b9e01SMark F. Adams 
2795ef31b24SMark F. Adams 	/* collect IS info */
280a2f3521dSMark F. Adams 	ierr = PetscMalloc(ncrs_eq*sizeof(PetscInt), &newproc_idx);CHKERRQ(ierr);
281a4b7d37bSMark F. Adams 	ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
282a2f3521dSMark F. Adams 	targetPE = 1; /* bring to "front" of machine */
283a2f3521dSMark F. Adams 	/*targetPE = npe/new_npe;*/ /* spread partitioning across machine */
284a2f3521dSMark F. Adams 	for (kk = jj = 0 ; kk < nloc_old ; kk++){
285c5bfad50SMark F. Adams 	  for (ii = 0 ; ii < cr_bs ; ii++, jj++){
286a2f3521dSMark F. Adams 	    newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
287eb07cef2SMark F. Adams 	  }
2885ef31b24SMark F. Adams 	}
289a4b7d37bSMark F. Adams 	ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
290a4b7d37bSMark F. Adams 	ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2915ef31b24SMark F. Adams       }
2925ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2935a9b9e01SMark F. Adams 
294806fa848SBarry Smith       ierr = ISCreateGeneral(wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2958263b398SMark F. Adams       if (newproc_idx != 0) {
2968263b398SMark F. Adams 	ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2975ef31b24SMark F. Adams       }
298806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
299a2f3521dSMark F. Adams 
300a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
3015a9b9e01SMark F. Adams       /* find factor */
3025a9b9e01SMark F. Adams       if (new_npe == 1) rfactor = npe; /* easy */
3035a9b9e01SMark F. Adams       else {
3045a9b9e01SMark F. Adams 	PetscReal best_fact = 0.;
3055a9b9e01SMark F. Adams 	jj = -1;
3065a9b9e01SMark F. Adams 	for (kk = 1 ; kk <= npe ; kk++){
3075a9b9e01SMark F. Adams 	  if (npe%kk==0) { /* a candidate */
3085a9b9e01SMark F. Adams 	    PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe;
3095a9b9e01SMark F. Adams 	    if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
3105a9b9e01SMark F. Adams 	    if (fact > best_fact) {
3115a9b9e01SMark F. Adams 	      best_fact = fact; jj = kk;
3125a9b9e01SMark F. Adams 	    }
3135a9b9e01SMark F. Adams 	  }
3145a9b9e01SMark F. Adams 	}
3155a9b9e01SMark F. Adams 	if (jj != -1) rfactor = jj;
316a2f3521dSMark F. Adams 	else rfactor = 1; /* does this happen .. a prime */
3175a9b9e01SMark F. Adams       }
3185a9b9e01SMark F. Adams       new_npe = npe/rfactor;
3195a9b9e01SMark F. Adams 
3205a9b9e01SMark F. Adams       if (new_npe==nactive) {
321a2f3521dSMark F. Adams 	*a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
3225a9b9e01SMark F. Adams 	ierr = PetscFree(counts);CHKERRQ(ierr);
323a2f3521dSMark F. Adams 	if (pc_gamg->verbose>0){
324a2f3521dSMark F. Adams           PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq(loc)=%d\n",mype,__FUNCT__,new_npe,ncrs_eq);
325a2f3521dSMark F. Adams         }
3265a9b9e01SMark F. Adams 	PetscFunctionReturn(0);
3275a9b9e01SMark F. Adams       }
3285a9b9e01SMark F. Adams 
329a2f3521dSMark F. Adams       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",mype,__FUNCT__,ncrs_eq);
3305a9b9e01SMark F. Adams       targetPE = mype/rfactor;
331a2f3521dSMark F. Adams       ierr = ISCreateStride(wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
3325a9b9e01SMark F. Adams 
333a2f3521dSMark F. Adams       if (stokes) {
334c5bfad50SMark F. Adams         ierr = ISCreateStride(wcomm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim);CHKERRQ(ierr);
3355a9b9e01SMark F. Adams       }
336a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
337e33ef3b1SMark F. Adams 
33811e60469SMark F. Adams     /*
339a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
34011e60469SMark F. Adams      */
341a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
342a2f3521dSMark F. Adams     if (stokes) {
343a2f3521dSMark F. Adams       ierr = ISPartitioningToNumbering(is_eq_newproc_prim, &is_eq_num_prim);CHKERRQ(ierr);
344806fa848SBarry Smith     } else is_eq_num_prim = is_eq_num;
34511e60469SMark F. Adams     /*
346a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
34711e60469SMark F. Adams      */
348a2f3521dSMark F. Adams     ierr = ISPartitioningCount(is_eq_newproc, npe, counts);CHKERRQ(ierr);
349a2f3521dSMark F. Adams     ncrs_eq_new = counts[mype];
350a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
351a2f3521dSMark F. Adams     if (stokes) {
352a2f3521dSMark F. Adams       ierr = ISPartitioningCount(is_eq_newproc_prim, npe, counts);CHKERRQ(ierr);
353a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_newproc_prim);CHKERRQ(ierr);
354c5bfad50SMark F. Adams       ncrs_prim_new = counts[mype]/cr_bs; /* nodes */
355806fa848SBarry Smith     } else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */
356a2f3521dSMark F. Adams 
357a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
3580cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3590cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
360b4fbaa2aSMark F. Adams #endif
361a2f3521dSMark F. Adams 
362a2f3521dSMark F. Adams     /* move data (for primal equations only) */
36322063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
364*3bf036e2SBarry Smith     ierr = VecCreate(wcomm, &dest_crd);CHKERRQ(ierr);
365a2f3521dSMark F. Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE);CHKERRQ(ierr);
366470e26f8SMark F. Adams     ierr = VecSetFromOptions(dest_crd);CHKERRQ(ierr); /* this is needed! */
36711e60469SMark F. Adams     /*
3689d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
369c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
37011e60469SMark F. Adams      */
371a2f3521dSMark F. Adams     ierr = PetscMalloc((ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx);CHKERRQ(ierr);
372a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
373a2f3521dSMark F. Adams     for (ii=0,jj=0; ii<ncrs_prim ; ii++) {
374c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
375a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz ; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
37611e60469SMark F. Adams     }
377a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
378806fa848SBarry Smith     ierr = ISCreateGeneral(wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
37992a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
38011e60469SMark F. Adams     /*
38111e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
38211e60469SMark F. Adams      */
383a2f3521dSMark F. Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd);CHKERRQ(ierr);
3849d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols ; jj++) {
385a2f3521dSMark F. Adams       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
386a2f3521dSMark F. Adams       for (ii=0 ; ii<ncrs_prim ; ii++) {
3879d5b6da9SMark F. Adams 	for (kk=0; kk<ndata_rows ; kk++) {
388a2f3521dSMark F. Adams 	  PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
389c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
390676e1743SMark F. Adams 	  ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
391d3d6bff4SMark F. Adams 	}
392038e3b61SMark F. Adams       }
393eb07cef2SMark F. Adams     }
394eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
395eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
39611e60469SMark F. Adams     /*
39711e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
39811e60469SMark F. Adams       to the correct processor
39911e60469SMark F. Adams     */
400806fa848SBarry Smith     ierr = VecScatterCreate(src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
40111e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
40211e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40311e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40411e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
40511e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
40611e60469SMark F. Adams     /*
40711e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
40811e60469SMark F. Adams     */
409c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
410a2f3521dSMark F. Adams     ierr = PetscMalloc(node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
411a2f3521dSMark F. Adams     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
412a2f3521dSMark F. Adams     strideNew = ncrs_prim_new*ndata_rows;
41311e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
4149d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols ; jj++) {
415a2f3521dSMark F. Adams       for (ii=0 ; ii<ncrs_prim_new ; ii++) {
4169d5b6da9SMark F. Adams 	for (kk=0; kk<ndata_rows ; kk++) {
417a2f3521dSMark F. Adams 	  PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
418c8b0795cSMark F. Adams 	  pc_gamg->data[ix] = PetscRealPart(array[jx]);
419d3d6bff4SMark F. Adams 	}
420038e3b61SMark F. Adams       }
421038e3b61SMark F. Adams     }
42211e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
42311e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
424a2f3521dSMark F. Adams 
425a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
4260cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4270cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
428ed3f9983SMark F. Adams #endif
429a2f3521dSMark F. Adams 
43011e60469SMark F. Adams     /*
43111e60469SMark F. Adams       Invert for MatGetSubMatrix
43211e60469SMark F. Adams     */
433a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
434a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
435c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
436a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
437a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
438a2f3521dSMark F. Adams     }
439a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
4400cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4410cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
4420cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
443ed3f9983SMark F. Adams #endif
444a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
445a2f3521dSMark F. Adams     {
446a2f3521dSMark F. Adams       Mat mat;
447806fa848SBarry Smith       ierr = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
448a2f3521dSMark F. Adams       *a_Amat_crs = mat;
449c5bfad50SMark F. Adams 
450c5bfad50SMark F. Adams       if (!PETSC_TRUE){
451c5bfad50SMark F. Adams         PetscInt cbs, rbs;
452c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
453806fa848SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
454c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
455806fa848SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",mype,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr);
456c5bfad50SMark F. Adams       }
457a2f3521dSMark F. Adams     }
458038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
459a2f3521dSMark F. Adams 
4600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4610cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
462ed3f9983SMark F. Adams #endif
46311e60469SMark F. Adams     /* prolongator */
46411e60469SMark F. Adams     {
46511e60469SMark F. Adams       IS findices;
466a2f3521dSMark F. Adams       PetscInt Istart,Iend;
467a2f3521dSMark F. Adams       Mat Pnew;
468a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4690cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4700cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
471ed3f9983SMark F. Adams #endif
47211e60469SMark F. Adams       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
473c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
474806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
47511e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
476c5bfad50SMark F. Adams 
477c5bfad50SMark F. Adams       if (!PETSC_TRUE){
478c5bfad50SMark F. Adams         PetscInt cbs, rbs;
479c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
480806fa848SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
481c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
482806fa848SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
483c5bfad50SMark F. Adams       }
4840cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4850cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
486ed3f9983SMark F. Adams #endif
4873530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
488a2f3521dSMark F. Adams 
489a2f3521dSMark F. Adams       /* output - repartitioned */
490a2f3521dSMark F. Adams       *a_P_inout = Pnew;
491e33ef3b1SMark F. Adams     }
492a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4935b89ad90SMark F. Adams 
4945a9b9e01SMark F. Adams     *a_nactive_proc = new_npe; /* output */
495a2f3521dSMark F. Adams   }
4965a9b9e01SMark F. Adams 
497a2f3521dSMark F. Adams   /* outout matrix data */
498c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
499c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
500c8b0795cSMark F. Adams     if (llev==0) {
501c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
502c8b0795cSMark F. Adams       PetscViewerASCIIOpen(wcomm,fname,&viewer);
503c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
504c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
505c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
506c8b0795cSMark F. Adams     }
507c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
508c8b0795cSMark F. Adams     PetscViewerASCIIOpen(wcomm,fname,&viewer);
509c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
510c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
511c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
512c8b0795cSMark F. Adams   }
513c8b0795cSMark F. Adams 
5145b89ad90SMark F. Adams   PetscFunctionReturn(0);
5155b89ad90SMark F. Adams }
5165b89ad90SMark F. Adams 
5175b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5185b89ad90SMark F. Adams /*
5195b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5205b89ad90SMark F. Adams                     by setting data structures and options.
5215b89ad90SMark F. Adams 
5225b89ad90SMark F. Adams    Input Parameter:
5235b89ad90SMark F. Adams .  pc - the preconditioner context
5245b89ad90SMark F. Adams 
5255b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
5265b89ad90SMark F. Adams 
5275b89ad90SMark F. Adams    Notes:
5285b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
5295b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
5305b89ad90SMark F. Adams */
5315b89ad90SMark F. Adams #undef __FUNCT__
5325b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
5339d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5345b89ad90SMark F. Adams {
5355b89ad90SMark F. Adams   PetscErrorCode  ierr;
5369d5b6da9SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
5375b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
5382adcac29SMark F. Adams   Mat              Pmat = pc->pmat;
539a2f3521dSMark F. Adams   PetscInt         fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
5409d5b6da9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
5413530afc2SMark F. Adams   PetscMPIInt      mype,npe,nactivepe;
542587fa25dSMark F. Adams   Mat              Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
543c8b0795cSMark F. Adams   PetscReal        emaxs[GAMG_MAXLEVELS];
544e696ed0bSMark F. Adams   IS              *ASMLocalIDsArr[GAMG_MAXLEVELS];
545a2f3521dSMark F. Adams   GAMGKKTMat       kktMatsArr[GAMG_MAXLEVELS];
546a2f3521dSMark F. Adams   PetscLogDouble   nnz0=0.,nnztot=0.;
547737a81a9SMark F. Adams   MatInfo          info;
548878e152fSMark F. Adams   PetscBool        stokes = PETSC_FALSE, redo_mesh_setup = PETSC_FALSE;
5495ef31b24SMark F. Adams 
5505b89ad90SMark F. Adams   PetscFunctionBegin;
55184d3f75bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
55284d3f75bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
553a2f3521dSMark F. Adams   if (pc_gamg->verbose>2) PetscPrintf(wcomm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",mype,__FUNCT__,pc_gamg->setup_count,pc->setupcalled);
55484d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
555878e152fSMark F. Adams     if (redo_mesh_setup) {
556878e152fSMark F. Adams       /* reset everything */
557878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
558878e152fSMark F. Adams       pc->setupcalled = 0;
559806fa848SBarry Smith     } else {
56084d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
56103a628feSMark F. Adams       /* just do Galerkin grids */
56258471d46SMark F. Adams       Mat B,dA,dB;
563d5280255SMark F. Adams       assert(pc->setupcalled);
56458471d46SMark F. Adams 
5659d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
56658471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5679d5b6da9SMark F. Adams         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
56858471d46SMark F. Adams         /* (re)set to get dirty flag */
5699d5b6da9SMark F. Adams         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
57058471d46SMark F. Adams 
5719d5b6da9SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>-1; level--) {
57203a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
57384d3f75bSMark F. Adams           if (pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
57403a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
575084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
57603a628feSMark F. Adams             mglevels[level]->A = B;
577806fa848SBarry Smith           } else {
578084a8fe3SJed Brown             ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
57958471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
58003a628feSMark F. Adams           }
58158471d46SMark F. Adams         ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
58258471d46SMark F. Adams         dB = B;
58358471d46SMark F. Adams         }
5845f8cf99dSMark F. Adams       }
585d5280255SMark F. Adams 
586d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
587d5280255SMark F. Adams 
588d5280255SMark F. Adams       /* PCSetUp_MG seems to insists on setting this to GMRES */
589d5280255SMark F. Adams       ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr);
590d5280255SMark F. Adams 
59158471d46SMark F. Adams       PetscFunctionReturn(0);
592eb07cef2SMark F. Adams     }
593878e152fSMark F. Adams   }
594f6536408SMark F. Adams 
595302f38e8SMark F. Adams   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
596eb07cef2SMark F. Adams 
5979d1b15c3SMark F. Adams   ierr = GAMGKKTMatCreate(Pmat, stokes, &kktMatsArr[0]);CHKERRQ(ierr);
5989d1b15c3SMark F. Adams 
599878e152fSMark F. Adams   if (!pc_gamg->data) {
600878e152fSMark F. Adams     if (pc_gamg->orig_data) {
601878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
602878e152fSMark F. Adams       ierr = MatGetLocalSize(Pmat, &qq, PETSC_NULL);CHKERRQ(ierr);
603878e152fSMark F. Adams       pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
604878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
605878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
606878e152fSMark F. Adams       ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
607878e152fSMark F. Adams       for (qq=0;qq<pc_gamg->data_sz;qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
608806fa848SBarry Smith     } else {
6099d5b6da9SMark F. Adams       if (!pc_gamg->createdefaultdata){
610c666adbfSMark F. Adams         SETERRQ(wcomm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
611302f38e8SMark F. Adams       }
612302f38e8SMark F. Adams       if (stokes) {
613c666adbfSMark F. Adams         SETERRQ(wcomm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
614eb07cef2SMark F. Adams       }
6159d1b15c3SMark F. Adams       ierr = pc_gamg->createdefaultdata(pc, kktMatsArr[0].A11);CHKERRQ(ierr);
6169d5b6da9SMark F. Adams     }
617878e152fSMark F. Adams   }
618878e152fSMark F. Adams 
619878e152fSMark F. Adams   /* cache original data for reuse */
620878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
621878e152fSMark F. Adams     ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr);
622878e152fSMark F. Adams     for (qq=0;qq<pc_gamg->data_sz;qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
623878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
624878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
625878e152fSMark F. Adams   }
626038e3b61SMark F. Adams 
627302f38e8SMark F. Adams   /* get basic dims */
628302f38e8SMark F. Adams   if (stokes) {
629a2f3521dSMark F. Adams     bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
630806fa848SBarry Smith   } else {
631302f38e8SMark F. Adams     ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
632302f38e8SMark F. Adams   }
633a2f3521dSMark F. Adams 
634a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
635c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
63684f9421dSMark F. Adams     PetscInt NN = M;
63784f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
63884f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
639*3bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
640806fa848SBarry Smith     } else {
641806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
64284f9421dSMark F. Adams     }
643b2a4f308SMark F. Adams     nnz0 = info.nz_used;
644b2a4f308SMark F. Adams     nnztot = info.nz_used;
645806fa848SBarry Smith     ierr = PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
646a2f3521dSMark F. Adams                 mype,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
647806fa848SBarry Smith                 (int)(nnz0/(PetscReal)NN),npe);CHKERRQ(ierr);
648c8b0795cSMark F. Adams   }
64984d3f75bSMark F. Adams 
650a2f3521dSMark F. Adams   /* Get A_i and R_i */
6518f4b7eb5SMark F. Adams   for (level=0, Aarr[0]=Pmat, nactivepe = npe; /* hard wired stopping logic */
6529d5b6da9SMark F. Adams         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
6530205a208SMark F. Adams         level++){
6545b89ad90SMark F. Adams     level1 = level + 1;
6550cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6560cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
657a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
658a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
659b4fbaa2aSMark F. Adams #endif
660a2f3521dSMark F. Adams #endif
661a2f3521dSMark F. Adams     /* deal with Stokes, get sub matrices */
6629d1b15c3SMark F. Adams     if (level > 0) {
663a2f3521dSMark F. Adams       ierr = GAMGKKTMatCreate(Aarr[level], stokes, &kktMatsArr[level]);CHKERRQ(ierr);
6649d1b15c3SMark F. Adams     }
665c8b0795cSMark F. Adams     { /* construct prolongator */
666725b86d8SJed Brown       Mat Gmat;
6670cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
668a2f3521dSMark F. Adams       Mat Prol11,Prol22;
669c8b0795cSMark F. Adams 
670a2f3521dSMark F. Adams       ierr = pc_gamg->graph(pc,kktMatsArr[level].A11, &Gmat);CHKERRQ(ierr);
6710cbbd2e1SMark F. Adams       ierr = pc_gamg->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
672a2f3521dSMark F. Adams       ierr = pc_gamg->prolongator(pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11);CHKERRQ(ierr);
673c8b0795cSMark F. Adams 
674a2f3521dSMark F. Adams       /* could have failed to create new level */
675a2f3521dSMark F. Adams       if (Prol11){
6769d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6779d1b15c3SMark F. Adams         ierr = MatGetBlockSizes(Prol11, PETSC_NULL, &bs);CHKERRQ(ierr);
678a2f3521dSMark F. Adams 
679a2f3521dSMark F. Adams         if (stokes) {
680a2f3521dSMark F. Adams           if (!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
681a2f3521dSMark F. Adams           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
682a2f3521dSMark F. Adams           ierr = pc_gamg->formkktprol(pc, Prol11, kktMatsArr[level].A21, &Prol22);CHKERRQ(ierr);
683c8b0795cSMark F. Adams         }
684c8b0795cSMark F. Adams 
685a2f3521dSMark F. Adams         if (pc_gamg->optprol){
686c8b0795cSMark F. Adams           /* smooth */
687a2f3521dSMark F. Adams           ierr = pc_gamg->optprol(pc, kktMatsArr[level].A11, &Prol11);CHKERRQ(ierr);
688c8b0795cSMark F. Adams         }
689c8b0795cSMark F. Adams 
690a2f3521dSMark F. Adams         if (stokes) {
691a2f3521dSMark F. Adams           IS is_row[2] = {kktMatsArr[level].prim_is,kktMatsArr[level].constr_is};
692a2f3521dSMark F. Adams           Mat a[4] = {Prol11, PETSC_NULL, PETSC_NULL, Prol22 };
693a2f3521dSMark F. Adams           ierr = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1]);CHKERRQ(ierr);
694806fa848SBarry Smith         } else {
695a2f3521dSMark F. Adams           Parr[level1] = Prol11;
696a2f3521dSMark F. Adams         }
697806fa848SBarry Smith       } else Parr[level1] = PETSC_NULL;
698ffc955d6SMark F. Adams 
699ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
700806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
701ffc955d6SMark F. Adams       }
702ffc955d6SMark F. Adams 
703a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
70441b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
705a2f3521dSMark F. Adams     } /* construct prolongator scope */
7060cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7070cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
708c8b0795cSMark F. Adams #endif
7099d5b6da9SMark F. Adams     /* cache eigen estimate */
7109d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1){
7119d5b6da9SMark F. Adams       PetscBool flag;
712806fa848SBarry Smith       ierr = PetscObjectComposedDataGetReal((PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
7139d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
714806fa848SBarry Smith     } else emaxs[level] = -1.;
7152adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
716c8b0795cSMark F. Adams     if (!Parr[level1]) {
717806fa848SBarry Smith       if (pc_gamg->verbose) {
718*3bf036e2SBarry Smith         ierr =  PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);CHKERRQ(ierr);
719806fa848SBarry Smith       }
720c8b0795cSMark F. Adams       break;
721c8b0795cSMark F. Adams     }
7220cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7230cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
724b4fbaa2aSMark F. Adams #endif
725a2f3521dSMark F. Adams 
726a2f3521dSMark F. Adams     ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
727806fa848SBarry Smith                         stokes, &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr);
728a2f3521dSMark F. Adams 
7290cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7300cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
731b4fbaa2aSMark F. Adams #endif
732a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
733a2f3521dSMark F. Adams 
734a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0){
7350cbbd2e1SMark F. Adams       PetscInt NN = M;
7360cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
737a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
738*3bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
739806fa848SBarry Smith       } else {
740806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
7410cbbd2e1SMark F. Adams       }
742a2f3521dSMark F. Adams 
7430cbbd2e1SMark F. Adams       nnztot += info.nz_used;
744806fa848SBarry Smith       ierr = PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
745a2f3521dSMark F. Adams                   mype,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
746806fa848SBarry Smith                   (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
747c8b0795cSMark F. Adams     }
748a2f3521dSMark F. Adams 
749a2f3521dSMark F. Adams     /* stop if one node -- could pull back for singular problems */
750c8b0795cSMark F. Adams     if (M/pc_gamg->data_cell_cols < 2) {
75181708dccSMark F. Adams       level++;
75281708dccSMark F. Adams       break;
75381708dccSMark F. Adams     }
7540cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
755b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
756b4fbaa2aSMark F. Adams #endif
757c8b0795cSMark F. Adams   } /* levels */
758c8b0795cSMark F. Adams 
759c8b0795cSMark F. Adams   if (pc_gamg->data) {
760c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
761a2f3521dSMark F. Adams     pc_gamg->data = PETSC_NULL;
7625b89ad90SMark F. Adams   }
763c8b0795cSMark F. Adams 
7640cbbd2e1SMark F. Adams   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
7659d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7665b89ad90SMark F. Adams   fine_level = level;
7679d5b6da9SMark F. Adams   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
7685b89ad90SMark F. Adams 
76984d3f75bSMark F. Adams   /* simple setup */
77084d3f75bSMark F. Adams   if (!PETSC_TRUE){
77184d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
77284d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
77384d3f75bSMark F. Adams          lidx<fine_level;
77484d3f75bSMark F. Adams          lidx++, level--){
77584d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
77684d3f75bSMark F. Adams       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77784d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
77884d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
77984d3f75bSMark F. Adams     }
78084d3f75bSMark F. Adams     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
78184d3f75bSMark F. Adams 
78284d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
783806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
784d5280255SMark F. Adams     /* set default smoothers & set operators */
7859d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
786587fa25dSMark F. Adams           lidx <= fine_level;
787587fa25dSMark F. Adams           lidx++, level--) {
788ffc955d6SMark F. Adams       KSP smoother;
789ffc955d6SMark F. Adams       PC subpc;
790a2f3521dSMark F. Adams 
7919d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
792f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
793ffc955d6SMark F. Adams 
794a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
795a2f3521dSMark F. Adams       /* set ops */
796a2f3521dSMark F. Adams       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
797a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
798a2f3521dSMark F. Adams 
799a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
800a2f3521dSMark F. Adams       if (stokes) {
801a2f3521dSMark F. Adams         KSP *ksps;
802a2f3521dSMark F. Adams         PetscInt nn;
803a2f3521dSMark F. Adams         ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);CHKERRQ(ierr);
804a2f3521dSMark F. Adams         ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is);CHKERRQ(ierr);
805a2f3521dSMark F. Adams         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
806a2f3521dSMark F. Adams         smoother = ksps[0];
807a2f3521dSMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
808a2f3521dSMark F. Adams         ierr = PetscFree(ksps);CHKERRQ(ierr);
809a2f3521dSMark F. Adams       }
810a2f3521dSMark F. Adams       ierr = GAMGKKTMatDestroy(&kktMatsArr[level]);CHKERRQ(ierr);
811a2f3521dSMark F. Adams 
812a2f3521dSMark F. Adams       /* set defaults */
8136c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
814a2f3521dSMark F. Adams 
815ffc955d6SMark F. Adams       /* override defaults and command line args (!) */
816ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
8172d3561bbSSatish Balay         PetscInt sz;
8182d3561bbSSatish Balay         IS *is;
819a2f3521dSMark F. Adams 
8202d3561bbSSatish Balay         sz = nASMBlocksArr[level];
8212d3561bbSSatish Balay         is = ASMLocalIDsArr[level];
822ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
823ffc955d6SMark F. Adams         if (sz==0){
824ffc955d6SMark F. Adams           IS is;
825ffc955d6SMark F. Adams           PetscInt my0,kk;
826ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
827ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
82806b43e7eSDmitry Karpeev           ierr = PCGASMSetSubdomains(subpc, 1, &is, PETSC_NULL);CHKERRQ(ierr);
829a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
830806fa848SBarry Smith         } else {
831a94c3b12SMark F. Adams           PetscInt kk;
83206b43e7eSDmitry Karpeev           ierr = PCGASMSetSubdomains(subpc, sz, is, PETSC_NULL);CHKERRQ(ierr);
833a94c3b12SMark F. Adams           for (kk=0;kk<sz;kk++){
834a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
835a94c3b12SMark F. Adams           }
836ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
837ffc955d6SMark F. Adams         }
838ffc955d6SMark F. Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
839ffc955d6SMark F. Adams 
840ffc955d6SMark F. Adams         ASMLocalIDsArr[level] = PETSC_NULL;
841ffc955d6SMark F. Adams         nASMBlocksArr[level] = 0;
842ffc955d6SMark F. Adams         ierr = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
843806fa848SBarry Smith       } else {
8449d5b6da9SMark F. Adams         ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr);
845ffc955d6SMark F. Adams       }
846d5280255SMark F. Adams     }
847d5280255SMark F. Adams     {
848d5280255SMark F. Adams       /* coarse grid */
849d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
850d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
851d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
852d5280255SMark F. Adams       ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
853d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
854d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
855d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
856d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
857d5280255SMark F. Adams       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
858d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
859d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
860d5280255SMark F. Adams     }
861d5280255SMark F. Adams 
862d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
863d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
864d5280255SMark F. Adams     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr);
865d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
866d5280255SMark F. Adams     if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
867d5280255SMark F. Adams 
868d5280255SMark F. Adams     /* create cheby smoothers */
869d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
870d5280255SMark F. Adams           lidx <= fine_level;
871d5280255SMark F. Adams           lidx++, level--) {
872d5280255SMark F. Adams       KSP smoother;
873ffc955d6SMark F. Adams       PetscBool flag;
874d5280255SMark F. Adams       PC subpc;
875d5280255SMark F. Adams 
876ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
877a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
878a2f3521dSMark F. Adams 
879a2f3521dSMark F. Adams       /* create field split PC, get subsmoother */
880a2f3521dSMark F. Adams       if (stokes) {
881a2f3521dSMark F. Adams         KSP *ksps;
882a2f3521dSMark F. Adams         PetscInt nn;
883a2f3521dSMark F. Adams         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
884a2f3521dSMark F. Adams         smoother = ksps[0];
885a2f3521dSMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
886a2f3521dSMark F. Adams         ierr = PetscFree(ksps);CHKERRQ(ierr);
887a2f3521dSMark F. Adams       }
888ffc955d6SMark F. Adams 
889ffc955d6SMark F. Adams       /* do my own cheby */
8906c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
891ffc955d6SMark F. Adams       if (flag) {
892ffc955d6SMark F. Adams         PetscReal emax, emin;
893251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
894ffc955d6SMark F. Adams         if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
895587fa25dSMark F. Adams         else{ /* eigen estimate 'emax' */
896e696ed0bSMark F. Adams           KSP eksp;
897e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
898b2a4f308SMark F. Adams           Vec bb, xx;
899038e3b61SMark F. Adams 
9005745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr);
9015745e0f5SMark F. Adams           ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr);
902fc4362bfSMark F. Adams           {
903fc4362bfSMark F. Adams             PetscRandom    rctx;
904fc4362bfSMark F. Adams             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
905fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
906fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
907fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
9085b89ad90SMark F. Adams           }
909a2f3521dSMark F. Adams 
910e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
911e696ed0bSMark F. Adams           {
912e696ed0bSMark F. Adams             PetscInt Istart,Iend,ncols,jj,Ii;
913e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
914e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
915e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0 ; Ii < Iend ; Ii++, jj++) {
916e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
917e696ed0bSMark F. Adams               if(ncols <= 1) {
918e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
919a94c3b12SMark F. Adams               }
920e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
921a94c3b12SMark F. Adams             }
922a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
923a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
924a94c3b12SMark F. Adams           }
925e696ed0bSMark F. Adams 
926fc4362bfSMark F. Adams           ierr = KSPCreate(wcomm, &eksp);CHKERRQ(ierr);
927806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
928fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
9291a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
9301a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
931f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
932f6536408SMark F. Adams 
933f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
934f6536408SMark F. Adams           ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
935fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
9365a9b9e01SMark F. Adams 
937d3d0db20SJed Brown           /* set PC type to be same as smoother */
938ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
939b2a4f308SMark F. Adams 
9405a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
9415a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
9425a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
943fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
9445a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
9455a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
9465a9b9e01SMark F. Adams 
947fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
9485a9b9e01SMark F. Adams 
949fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
950fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
951fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
952f6536408SMark F. Adams 
953ffc955d6SMark F. Adams           if (pc_gamg->verbose > 0) {
954a94c3b12SMark F. Adams             PetscInt N1, tt;
955a94c3b12SMark F. Adams             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
956a94c3b12SMark F. Adams             PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1);
957f6536408SMark F. Adams           }
958fc4362bfSMark F. Adams         }
959038e3b61SMark F. Adams         {
960c5bfad50SMark F. Adams           PetscInt N1, N0;
961c5bfad50SMark F. Adams           ierr = MatGetSize(Aarr[level], &N1, PETSC_NULL);CHKERRQ(ierr);
962c5bfad50SMark F. Adams           ierr = MatGetSize(Aarr[level+1], &N0, PETSC_NULL);CHKERRQ(ierr);
963f6536408SMark F. Adams           /* heuristic - is this crap? */
964b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
9655e7c91beSJed Brown 	  emin = emax * pc_gamg->eigtarget[0];
9665e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
967038e3b61SMark F. Adams         }
9686c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
969ffc955d6SMark F. Adams       } /* setup checby flag */
970ffc955d6SMark F. Adams     } /* non-coarse levels */
971737a81a9SMark F. Adams 
972d5280255SMark F. Adams     /* clean up */
973d5280255SMark F. Adams     for (level=1;level<pc_gamg->Nlevels;level++){
974587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
975587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
9765b89ad90SMark F. Adams     }
9770cbbd2e1SMark F. Adams 
9780cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9790cbbd2e1SMark F. Adams 
980a94c3b12SMark F. Adams     if (PETSC_FALSE){
98158471d46SMark F. Adams       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
9829d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
98358471d46SMark F. Adams       ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
98458471d46SMark F. Adams     }
985806fa848SBarry Smith   } else {
9865f8cf99dSMark F. Adams     KSP smoother;
987b43b03e9SMark F. Adams     if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__);
9889d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
9895f8cf99dSMark F. Adams     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9905f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
9919d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
9925f8cf99dSMark F. Adams   }
99384d3f75bSMark F. Adams 
9945b89ad90SMark F. Adams   PetscFunctionReturn(0);
9955b89ad90SMark F. Adams }
9965b89ad90SMark F. Adams 
997eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
9985b89ad90SMark F. Adams /*
9995b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
10005b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
10015b89ad90SMark F. Adams 
10025b89ad90SMark F. Adams    Input Parameter:
10035b89ad90SMark F. Adams .  pc - the preconditioner context
10045b89ad90SMark F. Adams 
10055b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
10065b89ad90SMark F. Adams */
10075b89ad90SMark F. Adams #undef __FUNCT__
10085b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
10095b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
10105b89ad90SMark F. Adams {
10115b89ad90SMark F. Adams   PetscErrorCode  ierr;
10125b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
10135b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
10145b89ad90SMark F. Adams 
10155b89ad90SMark F. Adams   PetscFunctionBegin;
10165b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
10175b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
10185b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
10195b89ad90SMark F. Adams   PetscFunctionReturn(0);
10205b89ad90SMark F. Adams }
10215b89ad90SMark F. Adams 
1022676e1743SMark F. Adams 
1023676e1743SMark F. Adams #undef __FUNCT__
1024676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
1025676e1743SMark F. Adams /*@
1026676e1743SMark F. Adams    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
1027676e1743SMark F. Adams    processor reduction.
1028676e1743SMark F. Adams 
1029676e1743SMark F. Adams    Not Collective on PC
1030676e1743SMark F. Adams 
1031676e1743SMark F. Adams    Input Parameters:
1032676e1743SMark F. Adams .  pc - the preconditioner context
1033676e1743SMark F. Adams 
1034676e1743SMark F. Adams 
1035676e1743SMark F. Adams    Options Database Key:
1036676e1743SMark F. Adams .  -pc_gamg_process_eq_limit
1037676e1743SMark F. Adams 
1038676e1743SMark F. Adams    Level: intermediate
1039676e1743SMark F. Adams 
1040676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1041676e1743SMark F. Adams 
1042676e1743SMark F. Adams .seealso: ()
1043676e1743SMark F. Adams @*/
1044676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
1045676e1743SMark F. Adams {
1046676e1743SMark F. Adams   PetscErrorCode ierr;
1047676e1743SMark F. Adams 
1048676e1743SMark F. Adams   PetscFunctionBegin;
1049676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1050676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1051676e1743SMark F. Adams   PetscFunctionReturn(0);
1052676e1743SMark F. Adams }
1053676e1743SMark F. Adams 
1054676e1743SMark F. Adams EXTERN_C_BEGIN
1055676e1743SMark F. Adams #undef __FUNCT__
1056676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
1057676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1058676e1743SMark F. Adams {
1059c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1060c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1061676e1743SMark F. Adams 
1062676e1743SMark F. Adams   PetscFunctionBegin;
10639d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
1064676e1743SMark F. Adams   PetscFunctionReturn(0);
1065676e1743SMark F. Adams }
1066676e1743SMark F. Adams EXTERN_C_END
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 EXTERN_C_BEGIN
1099389730f3SMark F. Adams #undef __FUNCT__
1100389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
1101389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1102389730f3SMark F. Adams {
1103389730f3SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1104389730f3SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1105389730f3SMark F. Adams 
1106389730f3SMark F. Adams   PetscFunctionBegin;
11079d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
1108389730f3SMark F. Adams   PetscFunctionReturn(0);
1109389730f3SMark F. Adams }
1110389730f3SMark F. Adams EXTERN_C_END
1111389730f3SMark F. Adams 
1112389730f3SMark F. Adams #undef __FUNCT__
11138263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1114676e1743SMark F. Adams /*@
11158263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1116676e1743SMark F. Adams 
1117676e1743SMark F. Adams    Collective on PC
1118676e1743SMark F. Adams 
1119676e1743SMark F. Adams    Input Parameters:
1120676e1743SMark F. Adams .  pc - the preconditioner context
1121676e1743SMark F. Adams 
1122676e1743SMark F. Adams 
1123676e1743SMark F. Adams    Options Database Key:
11248263b398SMark F. Adams .  -pc_gamg_repartition
1125676e1743SMark F. Adams 
1126676e1743SMark F. Adams    Level: intermediate
1127676e1743SMark F. Adams 
1128676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1129676e1743SMark F. Adams 
1130676e1743SMark F. Adams .seealso: ()
1131676e1743SMark F. Adams @*/
11328263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1133676e1743SMark F. Adams {
1134676e1743SMark F. Adams   PetscErrorCode ierr;
1135676e1743SMark F. Adams 
1136676e1743SMark F. Adams   PetscFunctionBegin;
1137676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11388263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1139676e1743SMark F. Adams   PetscFunctionReturn(0);
1140676e1743SMark F. Adams }
1141676e1743SMark F. Adams 
1142676e1743SMark F. Adams EXTERN_C_BEGIN
1143676e1743SMark F. Adams #undef __FUNCT__
11448263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
11458263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1146676e1743SMark F. Adams {
1147c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1148c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1149676e1743SMark F. Adams 
1150676e1743SMark F. Adams   PetscFunctionBegin;
11519d5b6da9SMark F. Adams   pc_gamg->repart = n;
1152676e1743SMark F. Adams   PetscFunctionReturn(0);
1153676e1743SMark F. Adams }
1154676e1743SMark F. Adams EXTERN_C_END
1155676e1743SMark F. Adams 
1156676e1743SMark F. Adams #undef __FUNCT__
1157ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1158ffc955d6SMark F. Adams /*@
1159ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1160ffc955d6SMark F. Adams 
1161ffc955d6SMark F. Adams    Collective on PC
1162ffc955d6SMark F. Adams 
1163ffc955d6SMark F. Adams    Input Parameters:
1164ffc955d6SMark F. Adams .  pc - the preconditioner context
1165ffc955d6SMark F. Adams 
1166ffc955d6SMark F. Adams 
1167ffc955d6SMark F. Adams    Options Database Key:
1168ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1169ffc955d6SMark F. Adams 
1170ffc955d6SMark F. Adams    Level: intermediate
1171ffc955d6SMark F. Adams 
1172ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1173ffc955d6SMark F. Adams 
1174ffc955d6SMark F. Adams .seealso: ()
1175ffc955d6SMark F. Adams @*/
1176ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1177ffc955d6SMark F. Adams {
1178ffc955d6SMark F. Adams   PetscErrorCode ierr;
1179ffc955d6SMark F. Adams 
1180ffc955d6SMark F. Adams   PetscFunctionBegin;
1181ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1182ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1183ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1184ffc955d6SMark F. Adams }
1185ffc955d6SMark F. Adams 
1186ffc955d6SMark F. Adams EXTERN_C_BEGIN
1187ffc955d6SMark F. Adams #undef __FUNCT__
1188ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
1189ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1190ffc955d6SMark F. Adams {
1191ffc955d6SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1192ffc955d6SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1193ffc955d6SMark F. Adams 
1194ffc955d6SMark F. Adams   PetscFunctionBegin;
1195ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1196ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1197ffc955d6SMark F. Adams }
1198ffc955d6SMark F. Adams EXTERN_C_END
1199ffc955d6SMark F. Adams 
1200ffc955d6SMark F. Adams #undef __FUNCT__
12014ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
12024ef23d27SMark F. Adams /*@
12034ef23d27SMark F. Adams    PCGAMGSetNlevels -
12044ef23d27SMark F. Adams 
12054ef23d27SMark F. Adams    Not collective on PC
12064ef23d27SMark F. Adams 
12074ef23d27SMark F. Adams    Input Parameters:
12084ef23d27SMark F. Adams .  pc - the preconditioner context
12094ef23d27SMark F. Adams 
12104ef23d27SMark F. Adams 
12114ef23d27SMark F. Adams    Options Database Key:
12124ef23d27SMark F. Adams .  -pc_mg_levels
12134ef23d27SMark F. Adams 
12144ef23d27SMark F. Adams    Level: intermediate
12154ef23d27SMark F. Adams 
12164ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12174ef23d27SMark F. Adams 
12184ef23d27SMark F. Adams .seealso: ()
12194ef23d27SMark F. Adams @*/
12204ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12214ef23d27SMark F. Adams {
12224ef23d27SMark F. Adams   PetscErrorCode ierr;
12234ef23d27SMark F. Adams 
12244ef23d27SMark F. Adams   PetscFunctionBegin;
12254ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12264ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
12274ef23d27SMark F. Adams   PetscFunctionReturn(0);
12284ef23d27SMark F. Adams }
12294ef23d27SMark F. Adams 
12304ef23d27SMark F. Adams EXTERN_C_BEGIN
12314ef23d27SMark F. Adams #undef __FUNCT__
12324ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
12334ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12344ef23d27SMark F. Adams {
12354ef23d27SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
12364ef23d27SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
12374ef23d27SMark F. Adams 
12384ef23d27SMark F. Adams   PetscFunctionBegin;
12399d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12404ef23d27SMark F. Adams   PetscFunctionReturn(0);
12414ef23d27SMark F. Adams }
12424ef23d27SMark F. Adams EXTERN_C_END
12434ef23d27SMark F. Adams 
12444ef23d27SMark F. Adams #undef __FUNCT__
12453542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
12463542efc5SMark F. Adams /*@
12473542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12483542efc5SMark F. Adams 
12493542efc5SMark F. Adams    Not collective on PC
12503542efc5SMark F. Adams 
12513542efc5SMark F. Adams    Input Parameters:
12523542efc5SMark F. Adams .  pc - the preconditioner context
12533542efc5SMark F. Adams 
12543542efc5SMark F. Adams 
12553542efc5SMark F. Adams    Options Database Key:
12563542efc5SMark F. Adams .  -pc_gamg_threshold
12573542efc5SMark F. Adams 
12583542efc5SMark F. Adams    Level: intermediate
12593542efc5SMark F. Adams 
12603542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
12613542efc5SMark F. Adams 
12623542efc5SMark F. Adams .seealso: ()
12633542efc5SMark F. Adams @*/
12643542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
12653542efc5SMark F. Adams {
12663542efc5SMark F. Adams   PetscErrorCode ierr;
12673542efc5SMark F. Adams 
12683542efc5SMark F. Adams   PetscFunctionBegin;
12693542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12703542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
12713542efc5SMark F. Adams   PetscFunctionReturn(0);
12723542efc5SMark F. Adams }
12733542efc5SMark F. Adams 
12743542efc5SMark F. Adams EXTERN_C_BEGIN
12753542efc5SMark F. Adams #undef __FUNCT__
12763542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
12773542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
12783542efc5SMark F. Adams {
1279c20e4228SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1280c20e4228SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
12813542efc5SMark F. Adams 
12823542efc5SMark F. Adams   PetscFunctionBegin;
12839d5b6da9SMark F. Adams   pc_gamg->threshold = n;
12843542efc5SMark F. Adams   PetscFunctionReturn(0);
12853542efc5SMark F. Adams }
12863542efc5SMark F. Adams EXTERN_C_END
12873542efc5SMark F. Adams 
12883542efc5SMark F. Adams #undef __FUNCT__
12899d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1290676e1743SMark F. Adams /*@
12919d5b6da9SMark F. Adams    PCGAMGSetType - Set solution method - calls sub create method
1292676e1743SMark F. Adams 
1293676e1743SMark F. Adams    Collective on PC
1294676e1743SMark F. Adams 
1295676e1743SMark F. Adams    Input Parameters:
1296676e1743SMark F. Adams .  pc - the preconditioner context
1297676e1743SMark F. Adams 
1298676e1743SMark F. Adams 
1299676e1743SMark F. Adams    Options Database Key:
13003542efc5SMark F. Adams .  -pc_gamg_type
1301676e1743SMark F. Adams 
1302676e1743SMark F. Adams    Level: intermediate
1303676e1743SMark F. Adams 
1304676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1305676e1743SMark F. Adams 
1306676e1743SMark F. Adams .seealso: ()
1307676e1743SMark F. Adams @*/
130819fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1309676e1743SMark F. Adams {
1310676e1743SMark F. Adams   PetscErrorCode ierr;
1311676e1743SMark F. Adams 
1312676e1743SMark F. Adams   PetscFunctionBegin;
1313676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1314806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1315676e1743SMark F. Adams   PetscFunctionReturn(0);
1316676e1743SMark F. Adams }
1317676e1743SMark F. Adams 
1318676e1743SMark F. Adams EXTERN_C_BEGIN
1319676e1743SMark F. Adams #undef __FUNCT__
13209d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
132119fd82e9SBarry Smith PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1322676e1743SMark F. Adams {
13239d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
1324676e1743SMark F. Adams 
1325676e1743SMark F. Adams   PetscFunctionBegin;
1326806fa848SBarry Smith   ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r);CHKERRQ(ierr);
13279d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
13289d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1329676e1743SMark F. Adams   PetscFunctionReturn(0);
1330676e1743SMark F. Adams }
1331676e1743SMark F. Adams EXTERN_C_END
1332676e1743SMark F. Adams 
13335b89ad90SMark F. Adams #undef __FUNCT__
13345b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
13355b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
13365b89ad90SMark F. Adams {
1337676e1743SMark F. Adams   PetscErrorCode  ierr;
1338676e1743SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1339676e1743SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1340676e1743SMark F. Adams   PetscBool        flag;
13415e7c91beSJed Brown   PetscInt         two = 2;
1342b43b03e9SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
13435b89ad90SMark F. Adams 
13445b89ad90SMark F. Adams   PetscFunctionBegin;
1345676e1743SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
1346676e1743SMark F. Adams   {
134775b74bdaSMark F. Adams     /* -pc_gamg_verbose */
13489d5b6da9SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
13499d5b6da9SMark F. Adams                            "none", pc_gamg->verbose,
1350806fa848SBarry Smith                            &pc_gamg->verbose, PETSC_NULL);CHKERRQ(ierr);
13518263b398SMark F. Adams     /* -pc_gamg_repartition */
13528263b398SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_repartition",
13538263b398SMark F. Adams                             "Repartion coarse grids (false)",
13548263b398SMark F. Adams                             "PCGAMGRepartitioning",
13559d5b6da9SMark F. Adams                             pc_gamg->repart,
13569d5b6da9SMark F. Adams                             &pc_gamg->repart,
1357806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1358ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
1359ffc955d6SMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1360ffc955d6SMark F. Adams                             "Use aggregation agragates for GASM smoother (false)",
1361ffc955d6SMark F. Adams                             "PCGAMGUseASMAggs",
1362ffc955d6SMark F. Adams                             pc_gamg->use_aggs_in_gasm,
1363ffc955d6SMark F. Adams                             &pc_gamg->use_aggs_in_gasm,
1364806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1365c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
1366676e1743SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1367676e1743SMark F. Adams                            "Limit (goal) on number of equations per process on coarse grids",
1368676e1743SMark F. Adams                            "PCGAMGSetProcEqLim",
13699d5b6da9SMark F. Adams                            pc_gamg->min_eq_proc,
13709d5b6da9SMark F. Adams                            &pc_gamg->min_eq_proc,
1371806fa848SBarry Smith                            &flag); CHKERRQ(ierr);
1372389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
1373389730f3SMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1374389730f3SMark F. Adams                            "Limit on number of equations for the coarse grid",
1375389730f3SMark F. Adams                            "PCGAMGSetCoarseEqLim",
13769d5b6da9SMark F. Adams                            pc_gamg->coarse_eq_limit,
13779d5b6da9SMark F. Adams                            &pc_gamg->coarse_eq_limit,
1378806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1379c20e4228SMark F. Adams     /* -pc_gamg_threshold */
13803542efc5SMark F. Adams     ierr = PetscOptionsReal("-pc_gamg_threshold",
13813542efc5SMark F. Adams                             "Relative threshold to use for dropping edges in aggregation graph",
13823542efc5SMark F. Adams                             "PCGAMGSetThreshold",
13839d5b6da9SMark F. Adams                             pc_gamg->threshold,
13849d5b6da9SMark F. Adams                             &pc_gamg->threshold,
1385806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1386806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
1387806fa848SBarry Smith       ierr = PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1388806fa848SBarry Smith     }
13894ef23d27SMark F. Adams 
13905e7c91beSJed Brown     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,PETSC_NULL);CHKERRQ(ierr);
13914ef23d27SMark F. Adams     ierr = PetscOptionsInt("-pc_mg_levels",
13924ef23d27SMark F. Adams                            "Set number of MG levels",
13934ef23d27SMark F. Adams                            "PCGAMGSetNlevels",
13949d5b6da9SMark F. Adams                            pc_gamg->Nlevels,
13959d5b6da9SMark F. Adams                            &pc_gamg->Nlevels,
1396806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
1397676e1743SMark F. Adams   }
1398676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1399676e1743SMark F. Adams 
14005b89ad90SMark F. Adams   PetscFunctionReturn(0);
14015b89ad90SMark F. Adams }
14025b89ad90SMark F. Adams 
14035b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
14045b89ad90SMark F. Adams /*MC
1405856830bfSJed Brown      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
14069d5b6da9SMark F. Adams        - This is the entry point to GAMG, registered in pcregis.c
14075b89ad90SMark F. Adams 
1408280d9858SJed Brown    Options Database Keys:
14095b89ad90SMark F. Adams    Multigrid options(inherited)
1410280d9858SJed Brown +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1411280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1412280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
14138c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
14145b89ad90SMark F. Adams 
14155b89ad90SMark F. Adams   Level: intermediate
1416280d9858SJed Brown 
14175b89ad90SMark F. Adams   Concepts: multigrid
14185b89ad90SMark F. Adams 
14195b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1420280d9858SJed Brown            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
14215b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
14225b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
14235b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
14245b89ad90SMark F. Adams M*/
14255b89ad90SMark F. Adams EXTERN_C_BEGIN
14265b89ad90SMark F. Adams #undef __FUNCT__
14275b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
14285b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
14295b89ad90SMark F. Adams {
14305b89ad90SMark F. Adams   PetscErrorCode  ierr;
14315b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
14325b89ad90SMark F. Adams   PC_MG           *mg;
14330cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14349d5b6da9SMark F. Adams   static long count = 0;
14355ee9c036SSatish Balay #endif
14365b89ad90SMark F. Adams 
14375b89ad90SMark F. Adams   PetscFunctionBegin;
143860a6d8f8SMark F. Adams 
14395b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14405b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
14415b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14425b89ad90SMark F. Adams 
14435b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
14445b89ad90SMark F. Adams   ierr = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr);
14455b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
1446ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
14475b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
14485b89ad90SMark F. Adams 
14499d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14509d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14519d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14529d5b6da9SMark F. Adams   pc_gamg->data = 0;
14535b89ad90SMark F. Adams 
14549d5b6da9SMark F. Adams   /* register AMG type */
14559d5b6da9SMark F. Adams   if (!GAMGList){
14569d5b6da9SMark F. Adams     ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
14579d5b6da9SMark F. Adams     ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
14589d5b6da9SMark F. Adams   }
14599d5b6da9SMark F. Adams 
14609d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
14615b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
14625b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
14635b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
14645b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
14655b89ad90SMark F. Adams 
14665b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1467676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_C",
1468676e1743SMark F. Adams 					    "PCGAMGSetProcEqLim_GAMG",
1469806fa848SBarry Smith 					    PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1470676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1471389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_C",
1472389730f3SMark F. Adams 					    "PCGAMGSetCoarseEqLim_GAMG",
1473806fa848SBarry Smith 					    PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1474389730f3SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
14758263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_C",
14768263b398SMark F. Adams 					    "PCGAMGSetRepartitioning_GAMG",
1477806fa848SBarry Smith 					    PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
1478676e1743SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1479ffc955d6SMark F. Adams 					    "PCGAMGSetUseASMAggs_C",
1480ffc955d6SMark F. Adams 					    "PCGAMGSetUseASMAggs_GAMG",
1481806fa848SBarry Smith 					    PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1482ffc955d6SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
14833542efc5SMark F. Adams 					    "PCGAMGSetThreshold_C",
14843542efc5SMark F. Adams 					    "PCGAMGSetThreshold_GAMG",
1485806fa848SBarry Smith 					    PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
14863542efc5SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
14879d5b6da9SMark F. Adams 					    "PCGAMGSetType_C",
14889d5b6da9SMark F. Adams 					    "PCGAMGSetType_GAMG",
1489806fa848SBarry Smith 					    PCGAMGSetType_GAMG);CHKERRQ(ierr);
14909d5b6da9SMark F. Adams   pc_gamg->repart = PETSC_FALSE;
1491ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
14929d5b6da9SMark F. Adams   pc_gamg->min_eq_proc = 100;
1493c8b0795cSMark F. Adams   pc_gamg->coarse_eq_limit = 800;
1494a2f3521dSMark F. Adams   pc_gamg->threshold = 0.001;
14959d5b6da9SMark F. Adams   pc_gamg->Nlevels = GAMG_MAXLEVELS;
14969d5b6da9SMark F. Adams   pc_gamg->verbose = 0;
14979d5b6da9SMark F. Adams   pc_gamg->emax_id = -1;
14985e7c91beSJed Brown   pc_gamg->eigtarget[0] = 0.05;
14995e7c91beSJed Brown   pc_gamg->eigtarget[1] = 1.05;
15009d5b6da9SMark F. Adams 
15010cbbd2e1SMark F. Adams   /* private events */
15020cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1503785cba28SMark F. Adams   if (count++ == 0) {
1504806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1505806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
15060cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
15070cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
15080cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1509806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1510806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1511806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1512806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1513806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1514806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1515806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1516806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1517806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1518806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1519806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1520806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1521f852f58cSMark F. Adams 
15220cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
15230cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
15240cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1525b4fbaa2aSMark F. Adams     /* create timer stages */
1526b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1527b4fbaa2aSMark F. Adams     {
1528b4fbaa2aSMark F. Adams       char     str[32];
1529b4fbaa2aSMark F. Adams       PetscInt lidx;
1530806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1531806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1532b4fbaa2aSMark F. Adams       for (lidx=1;lidx<9;lidx++){
1533b4fbaa2aSMark F. Adams 	sprintf(str,"MG Level %d",lidx);
1534806fa848SBarry Smith 	ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1535b4fbaa2aSMark F. Adams       }
1536b4fbaa2aSMark F. Adams     }
1537b4fbaa2aSMark F. Adams #endif
1538b4fbaa2aSMark F. Adams   }
1539b4fbaa2aSMark F. Adams #endif
15400cbbd2e1SMark F. Adams   /* general events */
15410cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
1542806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1543806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1544806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1545806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1546806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1547806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1548806fa848SBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1549806fa848SBarry Smith   ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr);
15500cbbd2e1SMark F. Adams #endif
15510cbbd2e1SMark F. Adams 
155260a6d8f8SMark F. Adams   /* instantiate derived type */
155360a6d8f8SMark F. Adams   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
155460a6d8f8SMark F. Adams   {
155560a6d8f8SMark F. Adams     char tname[256] = GAMGAGG;
1556806fa848SBarry Smith     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), PETSC_NULL);CHKERRQ(ierr);
155760a6d8f8SMark F. Adams     ierr = PCGAMGSetType(pc, tname);CHKERRQ(ierr);
155860a6d8f8SMark F. Adams   }
155960a6d8f8SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
15605b89ad90SMark F. Adams   PetscFunctionReturn(0);
15615b89ad90SMark F. Adams }
15625b89ad90SMark F. Adams EXTERN_C_END
1563