xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 1cc46a469a872222aebe168517d4f255f409e0aa)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4aaa7dc30SBarry Smith #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>
75b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */
8f96513f1SMatthew G Knepley 
90cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
100cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
11b4fbaa2aSMark F. Adams #endif
120cbbd2e1SMark F. Adams 
130cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
140cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_AGG;
150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGGgraph_GEO;
160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
200cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGOptprol_AGG;
210cbbd2e1SMark F. Adams #endif
220cbbd2e1SMark F. Adams 
23b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
24b4fbaa2aSMark F. Adams 
25b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
260cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
27b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
28b4fbaa2aSMark F. Adams #endif
29f96513f1SMatthew G Knepley 
30140e18c1SBarry Smith static PetscFunctionList GAMGList = 0;
313e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
329d5b6da9SMark F. Adams 
33d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
34d3d6bff4SMark F. Adams #undef __FUNCT__
35d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
36d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
37d3d6bff4SMark F. Adams {
38d3d6bff4SMark F. Adams   PetscErrorCode ierr;
39d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
40d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
41d3d6bff4SMark F. Adams 
42d3d6bff4SMark F. Adams   PetscFunctionBegin;
43a2f3521dSMark F. Adams   if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */
44ce94432eSBarry Smith     PetscPrintf(PetscObjectComm((PetscObject)pc),"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__);
459d5b6da9SMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
4658471d46SMark F. Adams   }
470298fd71SBarry Smith   pc_gamg->data = NULL; pc_gamg->data_sz = 0;
48878e152fSMark F. Adams 
49878e152fSMark F. Adams   if (pc_gamg->orig_data) {
50878e152fSMark F. Adams     ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
51878e152fSMark F. Adams   }
52a2f3521dSMark F. Adams   PetscFunctionReturn(0);
53a2f3521dSMark F. Adams }
54a2f3521dSMark F. Adams 
555b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
565b89ad90SMark F. Adams /*
57c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
58a147abb0SMark F. Adams      of active processors.
595b89ad90SMark F. Adams 
605b89ad90SMark F. Adams    Input Parameter:
61a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
62a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
639d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
64c5bfad50SMark F. Adams    . cr_bs - coarse block size
653530afc2SMark F. Adams    In/Output Parameter:
66a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
67afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
6811e60469SMark F. Adams    Output Parameter:
693530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
705b89ad90SMark F. Adams */
715cb416c2SMark F. Adams 
725b89ad90SMark F. Adams #undef __FUNCT__
73c238b0ebSToby Isaac #define __FUNCT__ "PCGAMGCreateLevel_GAMG"
74b34066adSToby Isaac static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,
757700e67bSMark Adams                                   Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc)
765b89ad90SMark F. Adams {
77a2f3521dSMark F. Adams   PetscErrorCode  ierr;
789d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
79486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
809d5b6da9SMark F. Adams   const PetscBool repart      = pc_gamg->repart;
819d5b6da9SMark F. Adams   const PetscInt  min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
82a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
833b4367a7SBarry Smith   MPI_Comm        comm;
84c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
853ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
865b89ad90SMark F. Adams 
875b89ad90SMark F. Adams   PetscFunctionBegin;
883b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
893b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
903b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
91c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
9211e60469SMark F. Adams   /* RAP */
939d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
94038e3b61SMark F. Adams 
953ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
960298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
973ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
983ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
993ae0bb68SMark Adams   }
1003ae0bb68SMark Adams   else {
1013ae0bb68SMark Adams     PetscInt  bs;
1023ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
1033ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
1043ae0bb68SMark Adams   }
105a2f3521dSMark F. Adams 
106c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
107a2f3521dSMark F. Adams   {
108472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
1090298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
110c5df96a5SBarry Smith     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
111c5df96a5SBarry Smith     if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1;
112c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
113a2f3521dSMark F. Adams   }
114f852f58cSMark F. Adams 
1152fa5cd67SKarl Rupp   if (!repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1162fa5cd67SKarl Rupp   else {
1173ae0bb68SMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
118885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
119e33ef3b1SMark F. Adams 
12071959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
12171959b99SBarry Smith     if (ncrs_eq % cr_bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs);
1220cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1230cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
124b4fbaa2aSMark F. Adams #endif
125a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
126785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
1277700e67bSMark Adams     if (repart) {
128a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1295a9b9e01SMark F. Adams       Mat adj;
1305a9b9e01SMark F. Adams 
131a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
1323b4367a7SBarry Smith         if (pc_gamg->verbose==1) PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,ncrs_eq);
133a2f3521dSMark F. Adams         else {
134a2f3521dSMark F. Adams           PetscInt n;
1353b4367a7SBarry Smith           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
1363b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n);
137a2f3521dSMark F. Adams         }
138a2f3521dSMark F. Adams       }
1395a9b9e01SMark F. Adams 
140a2f3521dSMark F. Adams       /* get 'adj' */
141c5bfad50SMark F. Adams       if (cr_bs == 1) {
142038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
143806fa848SBarry Smith       } else {
144a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
145eb07cef2SMark F. Adams         Mat               tMat;
146a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
147b4fbaa2aSMark F. Adams         const PetscScalar *vals;
148b4fbaa2aSMark F. Adams         const PetscInt    *idx;
149a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1509057884aSMark F. Adams         static PetscInt   llev = 0;
151b4fbaa2aSMark F. Adams 
152578f55a3SPeter Brune         ierr = PetscMalloc1(ncrs, &d_nnz);CHKERRQ(ierr);
153578f55a3SPeter Brune         ierr = PetscMalloc1(ncrs, &o_nnz);CHKERRQ(ierr);
154a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
155a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
156c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
15758471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
158c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
159c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
16058471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1613ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1623ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
16358471d46SMark F. Adams         }
1646876a03eSMark F. Adams 
1653b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1663ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
167a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
168a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
169a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
17058471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
1715f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
172eb07cef2SMark F. Adams 
173a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
174c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
17522063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
176eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
177c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
178eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
179eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
180eb07cef2SMark F. Adams           }
18122063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
182eb07cef2SMark F. Adams         }
183eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
184eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
185eb07cef2SMark F. Adams 
186b4fbaa2aSMark F. Adams         if (llev++ == -1) {
187b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
1888caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
1893b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
190b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
1913bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
192b4fbaa2aSMark F. Adams         }
193b4fbaa2aSMark F. Adams 
194eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
195eb07cef2SMark F. Adams 
196eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
197a2f3521dSMark F. Adams       } /* create 'adj' */
198f150b916SMark F. Adams 
199a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2005a9b9e01SMark F. Adams         char            prefix[256];
2015a9b9e01SMark F. Adams         const char      *pcpre;
202b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
203b4fbaa2aSMark F. Adams         MatPartitioning mpart;
204a4b7d37bSMark F. Adams         IS              proc_is;
205a2f3521dSMark F. Adams         PetscInt        targetPE;
2062f03bc48SMark F. Adams 
2073b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2085ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2099d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2108caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
21159a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
21211e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
213c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
214a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
21511e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2165a9b9e01SMark F. Adams 
2175ef31b24SMark F. Adams         /* collect IS info */
218785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
219a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
220a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
221c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
222a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
223c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
224a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
225eb07cef2SMark F. Adams           }
2265ef31b24SMark F. Adams         }
227a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
228a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2295ef31b24SMark F. Adams       }
2305ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2315a9b9e01SMark F. Adams 
2323b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2338263b398SMark F. Adams       if (newproc_idx != 0) {
2348263b398SMark F. Adams         ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2355ef31b24SMark F. Adams       }
236806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
237a2f3521dSMark F. Adams 
238a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
2395a9b9e01SMark F. Adams       /* find factor */
240c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2415a9b9e01SMark F. Adams       else {
2425a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2435a9b9e01SMark F. Adams         jj = -1;
244c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
245c5df96a5SBarry Smith           if (size%kk==0) { /* a candidate */
246c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2475a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2485a9b9e01SMark F. Adams             if (fact > best_fact) {
2495a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
2505a9b9e01SMark F. Adams             }
2515a9b9e01SMark F. Adams           }
2525a9b9e01SMark F. Adams         }
2535a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
254a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
2555a9b9e01SMark F. Adams       }
256c5df96a5SBarry Smith       new_size = size/rfactor;
2575a9b9e01SMark F. Adams 
258c5df96a5SBarry Smith       if (new_size==nactive) {
259a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2605a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
261a2f3521dSMark F. Adams         if (pc_gamg->verbose>0) {
2623b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq);
263a2f3521dSMark F. Adams         }
2645a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2655a9b9e01SMark F. Adams       }
2665a9b9e01SMark F. Adams 
2673b4367a7SBarry Smith       if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq);
268c5df96a5SBarry Smith       targetPE = rank/rfactor;
2693b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
270a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
271e33ef3b1SMark F. Adams 
27211e60469SMark F. Adams     /*
273a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
27411e60469SMark F. Adams      */
275a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2767700e67bSMark Adams     is_eq_num_prim = is_eq_num;
27711e60469SMark F. Adams     /*
278a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
27911e60469SMark F. Adams      */
280c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
281c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
282a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
2833ae0bb68SMark Adams     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */
284a2f3521dSMark F. Adams 
285a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2860cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2870cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
288b4fbaa2aSMark F. Adams #endif
289885364a3SMark Adams     /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxilary data into some complex abstracted thing */
290885364a3SMark Adams     {
291885364a3SMark Adams     Vec            src_crd, dest_crd;
292885364a3SMark 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;
293885364a3SMark Adams     VecScatter     vecscat;
294885364a3SMark Adams     PetscScalar    *array;
295885364a3SMark Adams     IS isscat;
296a2f3521dSMark F. Adams 
297a2f3521dSMark F. Adams     /* move data (for primal equations only) */
29822063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
2993b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
3003ae0bb68SMark Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
301c0dedaeaSBarry Smith     ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
30211e60469SMark F. Adams     /*
3039d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
304c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
30511e60469SMark F. Adams      */
306854ce69bSBarry Smith     ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
307a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3083ae0bb68SMark Adams     for (ii=0,jj=0; ii<ncrs; ii++) {
309c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
310a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
31111e60469SMark F. Adams     }
312a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3133ae0bb68SMark Adams     ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
31492a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
31511e60469SMark F. Adams     /*
31611e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
31711e60469SMark F. Adams      */
3183ae0bb68SMark Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
3199d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3203ae0bb68SMark Adams       const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3213ae0bb68SMark Adams       for (ii=0; ii<ncrs; ii++) {
3229d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
323a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
324c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
325676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
326d3d6bff4SMark F. Adams         }
327038e3b61SMark F. Adams       }
328eb07cef2SMark F. Adams     }
329eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
330eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
33111e60469SMark F. Adams     /*
33211e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
33311e60469SMark F. Adams       to the correct processor
33411e60469SMark F. Adams     */
3350298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
33611e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
33711e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
33811e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
33911e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
34011e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
34111e60469SMark F. Adams     /*
34211e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
34311e60469SMark F. Adams     */
344c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
345578f55a3SPeter Brune     ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3462fa5cd67SKarl Rupp 
3473ae0bb68SMark Adams     pc_gamg->data_sz = node_data_sz*ncrs_new;
3483ae0bb68SMark Adams     strideNew        = ncrs_new*ndata_rows;
3492fa5cd67SKarl Rupp 
35011e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3519d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3523ae0bb68SMark Adams       for (ii=0; ii<ncrs_new; ii++) {
3539d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
354a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
355c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
356d3d6bff4SMark F. Adams         }
357038e3b61SMark F. Adams       }
358038e3b61SMark F. Adams     }
35911e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
36011e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
361885364a3SMark Adams     }
362a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3630cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3640cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
365ed3f9983SMark F. Adams #endif
366a2f3521dSMark F. Adams 
36711e60469SMark F. Adams     /*
36811e60469SMark F. Adams       Invert for MatGetSubMatrix
36911e60469SMark F. Adams     */
370a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
371a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
372c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
373a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
374a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
375a2f3521dSMark F. Adams     }
376a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3770cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3780cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3790cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
380ed3f9983SMark F. Adams #endif
381a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
382a2f3521dSMark F. Adams     {
383a2f3521dSMark F. Adams       Mat mat;
384806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
385a2f3521dSMark F. Adams       *a_Amat_crs = mat;
386c5bfad50SMark F. Adams 
387c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
388c5bfad50SMark F. Adams         PetscInt cbs, rbs;
389c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
390c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
391c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
392c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",rank,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr);
393c5bfad50SMark F. Adams       }
394a2f3521dSMark F. Adams     }
395038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
396a2f3521dSMark F. Adams 
3970cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3980cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
399ed3f9983SMark F. Adams #endif
40011e60469SMark F. Adams     /* prolongator */
40111e60469SMark F. Adams     {
40211e60469SMark F. Adams       IS       findices;
403a2f3521dSMark F. Adams       PetscInt Istart,Iend;
404a2f3521dSMark F. Adams       Mat      Pnew;
405a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4060cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4070cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
408ed3f9983SMark F. Adams #endif
4093b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
410c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
411806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
41211e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
413c5bfad50SMark F. Adams 
414c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
415c5bfad50SMark F. Adams         PetscInt cbs, rbs;
416c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
417c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
418c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
419c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
420c5bfad50SMark F. Adams       }
4210cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4220cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
423ed3f9983SMark F. Adams #endif
4243530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
425a2f3521dSMark F. Adams 
426a2f3521dSMark F. Adams       /* output - repartitioned */
427a2f3521dSMark F. Adams       *a_P_inout = Pnew;
428e33ef3b1SMark F. Adams     }
429a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4305b89ad90SMark F. Adams 
431c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
432a2f3521dSMark F. Adams   }
4335a9b9e01SMark F. Adams 
434a2f3521dSMark F. Adams   /* outout matrix data */
435c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
436c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
437c8b0795cSMark F. Adams     if (llev==0) {
438c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
4393b4367a7SBarry Smith       PetscViewerASCIIOpen(comm,fname,&viewer);
440c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
441c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
442c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
443c8b0795cSMark F. Adams     }
444c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
4453b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
446c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
447c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
448c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
449c8b0795cSMark F. Adams   }
4505b89ad90SMark F. Adams   PetscFunctionReturn(0);
4515b89ad90SMark F. Adams }
4525b89ad90SMark F. Adams 
4535b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4545b89ad90SMark F. Adams /*
4555b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4565b89ad90SMark F. Adams                     by setting data structures and options.
4575b89ad90SMark F. Adams 
4585b89ad90SMark F. Adams    Input Parameter:
4595b89ad90SMark F. Adams .  pc - the preconditioner context
4605b89ad90SMark F. Adams 
4615b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4625b89ad90SMark F. Adams 
4635b89ad90SMark F. Adams    Notes:
4645b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4655b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4665b89ad90SMark F. Adams */
4675b89ad90SMark F. Adams #undef __FUNCT__
4685b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4699d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4705b89ad90SMark F. Adams {
4715b89ad90SMark F. Adams   PetscErrorCode ierr;
4729d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4735b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4742adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
475a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
4763b4367a7SBarry Smith   MPI_Comm       comm;
477c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
478587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
479c8b0795cSMark F. Adams   PetscReal      emaxs[GAMG_MAXLEVELS];
480e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
481a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
482737a81a9SMark F. Adams   MatInfo        info;
4837700e67bSMark Adams   PetscBool      redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);
4845ef31b24SMark F. Adams 
4855b89ad90SMark F. Adams   PetscFunctionBegin;
4863b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4873b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4883b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
489dfd5c07aSMark F. Adams 
4903b4367a7SBarry Smith   if (pc_gamg->verbose>2) PetscPrintf(comm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",rank,__FUNCT__,pc_gamg->setup_count,pc->setupcalled);
491dfd5c07aSMark F. Adams 
49284d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
493878e152fSMark F. Adams     if (redo_mesh_setup) {
494878e152fSMark F. Adams       /* reset everything */
495878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
496878e152fSMark F. Adams       pc->setupcalled = 0;
497806fa848SBarry Smith     } else {
49884d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
49903a628feSMark F. Adams       /* just do Galerkin grids */
50058471d46SMark F. Adams       Mat          B,dA,dB;
50158471d46SMark F. Adams 
50271959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
5039d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
50458471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
50523ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
50658471d46SMark F. Adams         /* (re)set to get dirty flag */
50723ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
50858471d46SMark F. Adams 
5092fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
51003a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
5110a97e771SToby Isaac           if (pc_gamg->setup_count==2) {
51203a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
513084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
5142fa5cd67SKarl Rupp 
51503a628feSMark F. Adams             mglevels[level]->A = B;
516806fa848SBarry Smith           } else {
51723ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
51858471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
51903a628feSMark F. Adams           }
52023ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
52158471d46SMark F. Adams           dB   = B;
52258471d46SMark F. Adams         }
5235f8cf99dSMark F. Adams       }
524d5280255SMark F. Adams 
525d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
526d5280255SMark F. Adams 
52758471d46SMark F. Adams       PetscFunctionReturn(0);
528eb07cef2SMark F. Adams     }
529878e152fSMark F. Adams   }
530f6536408SMark F. Adams 
531878e152fSMark F. Adams   if (!pc_gamg->data) {
532878e152fSMark F. Adams     if (pc_gamg->orig_data) {
533878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5340298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5352fa5cd67SKarl Rupp 
536878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
537878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
538878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5392fa5cd67SKarl Rupp 
540785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
541878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
542806fa848SBarry Smith     } else {
5431ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5447700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5459d5b6da9SMark F. Adams     }
546878e152fSMark F. Adams   }
547878e152fSMark F. Adams 
548878e152fSMark F. Adams   /* cache original data for reuse */
549878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
550785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
551878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
552878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
553878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
554878e152fSMark F. Adams   }
555038e3b61SMark F. Adams 
556302f38e8SMark F. Adams   /* get basic dims */
557302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
558a2f3521dSMark F. Adams 
559a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
560c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
56184f9421dSMark F. Adams     PetscInt NN = M;
56284f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
56384f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
5643bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
565fe06f982SMark Adams       if (!NN) NN=1;
566806fa848SBarry Smith     } else {
567806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
56884f9421dSMark F. Adams     }
569b2a4f308SMark F. Adams     nnz0   = info.nz_used;
570b2a4f308SMark F. Adams     nnztot = info.nz_used;
5713b4367a7SBarry Smith     ierr   = PetscPrintf(comm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
572c5df96a5SBarry Smith                          rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
573c5df96a5SBarry Smith                          (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr);
574c8b0795cSMark F. Adams   }
57584d3f75bSMark F. Adams 
576a2f3521dSMark F. Adams   /* Get A_i and R_i */
577c5df96a5SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
578c5df96a5SBarry Smith        level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit);  /* && (size==1 || nactivepe>1); */
5790205a208SMark F. Adams        level++) {
5805b89ad90SMark F. Adams     level1 = level + 1;
5810cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5820cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
583a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
584a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
585b4fbaa2aSMark F. Adams #endif
586a2f3521dSMark F. Adams #endif
587c8b0795cSMark F. Adams     { /* construct prolongator */
588725b86d8SJed Brown       Mat              Gmat;
5890cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5907700e67bSMark Adams       Mat              Prol11;
591c8b0795cSMark F. Adams 
5927700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5931ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5947700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
595c8b0795cSMark F. Adams 
596a2f3521dSMark F. Adams       /* could have failed to create new level */
597a2f3521dSMark F. Adams       if (Prol11) {
5989d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5990298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
600a2f3521dSMark F. Adams 
6011ab5ffc9SJed Brown         if (pc_gamg->ops->optprol) {
602c8b0795cSMark F. Adams           /* smooth */
6037700e67bSMark Adams           ierr = pc_gamg->ops->optprol(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
604c8b0795cSMark F. Adams         }
605c8b0795cSMark F. Adams 
6067700e67bSMark Adams         Parr[level1] = Prol11;
6070298fd71SBarry Smith       } else Parr[level1] = NULL;
608ffc955d6SMark F. Adams 
609ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
6101b18a24aSMark Adams         PetscInt bs;
6111b18a24aSMark Adams         ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
612806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
613ffc955d6SMark F. Adams       }
614ffc955d6SMark F. Adams 
615a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
61641b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
617a2f3521dSMark F. Adams     } /* construct prolongator scope */
6180cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6190cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
620c8b0795cSMark F. Adams #endif
6219d5b6da9SMark F. Adams     /* cache eigen estimate */
6229d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
6239d5b6da9SMark F. Adams       PetscBool flag;
6247700e67bSMark Adams       ierr = PetscObjectComposedDataGetReal((PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
6259d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
626806fa848SBarry Smith     } else emaxs[level] = -1.;
6272adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
628c8b0795cSMark F. Adams     if (!Parr[level1]) {
629806fa848SBarry Smith       if (pc_gamg->verbose) {
6303b4367a7SBarry Smith         ierr =  PetscPrintf(comm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr);
631806fa848SBarry Smith       }
632c8b0795cSMark F. Adams       break;
633c8b0795cSMark F. Adams     }
6340cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6350cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
636b4fbaa2aSMark F. Adams #endif
637a2f3521dSMark F. Adams 
638c238b0ebSToby Isaac     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs,
6397700e67bSMark Adams                        &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr);
640a2f3521dSMark F. Adams 
6410cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6420cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
643b4fbaa2aSMark F. Adams #endif
644a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
645a2f3521dSMark F. Adams 
646a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0) {
6470cbbd2e1SMark F. Adams       PetscInt NN = M;
6480cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
649a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
6503bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
651fe06f982SMark Adams         if (!NN) NN=1;
652806fa848SBarry Smith       } else {
653806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
6540cbbd2e1SMark F. Adams       }
655a2f3521dSMark F. Adams 
6560cbbd2e1SMark F. Adams       nnztot += info.nz_used;
6573b4367a7SBarry Smith       ierr    = PetscPrintf(comm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
658c5df96a5SBarry Smith                             rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
659806fa848SBarry Smith                             (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
660c8b0795cSMark F. Adams     }
661a2f3521dSMark F. Adams 
66225a145a7SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
66325a145a7SMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M < 2)
66425a145a7SMark Adams          || nactivepe==1 ) {
66581708dccSMark F. Adams       level++;
66681708dccSMark F. Adams       break;
66781708dccSMark F. Adams     }
6680cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
669b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
670b4fbaa2aSMark F. Adams #endif
671c8b0795cSMark F. Adams   } /* levels */
672c8b0795cSMark F. Adams 
673c8b0795cSMark F. Adams   if (pc_gamg->data) {
674c8b0795cSMark F. Adams     ierr          = PetscFree(pc_gamg->data);CHKERRQ(ierr);
6750298fd71SBarry Smith     pc_gamg->data = NULL;
6765b89ad90SMark F. Adams   }
677c8b0795cSMark F. Adams 
6783b4367a7SBarry Smith   if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
6799d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6805b89ad90SMark F. Adams   fine_level       = level;
6810298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
6825b89ad90SMark F. Adams 
68384d3f75bSMark F. Adams   /* simple setup */
68484d3f75bSMark F. Adams   if (!PETSC_TRUE) {
68584d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
68684d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
68784d3f75bSMark F. Adams          lidx<fine_level;
68884d3f75bSMark F. Adams          lidx++, level--) {
68984d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
69023ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level]);CHKERRQ(ierr);
69184d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
69284d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
69384d3f75bSMark F. Adams     }
69423ee1639SBarry Smith     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0]);CHKERRQ(ierr);
69584d3f75bSMark F. Adams 
69684d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
697806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
698d5280255SMark F. Adams     /* set default smoothers & set operators */
6999d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
700587fa25dSMark F. Adams          lidx <= fine_level;
701587fa25dSMark F. Adams          lidx++, level--) {
702ffc955d6SMark F. Adams       KSP smoother;
703ffc955d6SMark F. Adams       PC  subpc;
704a2f3521dSMark F. Adams 
7059d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
706f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
707ffc955d6SMark F. Adams 
708a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
709a2f3521dSMark F. Adams       /* set ops */
71023ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
711a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
712a2f3521dSMark F. Adams 
713a2f3521dSMark F. Adams       /* set defaults */
7146c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
715a2f3521dSMark F. Adams 
7161b18a24aSMark Adams       /* set blocks for GASM smoother that uses the 'aggregates' */
717ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
7182d3561bbSSatish Balay         PetscInt sz;
7192d3561bbSSatish Balay         IS       *is;
720a2f3521dSMark F. Adams 
7212d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
7222d3561bbSSatish Balay         is   = ASMLocalIDsArr[level];
723ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
7241b18a24aSMark Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
725ffc955d6SMark F. Adams         if (sz==0) {
726ffc955d6SMark F. Adams           IS       is;
727ffc955d6SMark F. Adams           PetscInt my0,kk;
728ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
729ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
7300298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr);
731a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
732806fa848SBarry Smith         } else {
733a94c3b12SMark F. Adams           PetscInt kk;
7340298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr);
735a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
736a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
737a94c3b12SMark F. Adams           }
738ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
739ffc955d6SMark F. Adams         }
7400298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
741ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
742ffc955d6SMark F. Adams         ierr                  = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
743806fa848SBarry Smith       } else {
744890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
745ffc955d6SMark F. Adams       }
746d5280255SMark F. Adams     }
747d5280255SMark F. Adams     {
748d5280255SMark F. Adams       /* coarse grid */
749d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
750d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
751d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
75223ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
753d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
754d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
755d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
756d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
75771959b99SBarry Smith       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
75871959b99SBarry Smith       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
759d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
760d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
7619dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
7622fb0b348SMark F. Adams       ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
7635b42dca8SJed Brown       /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
7645b42dca8SJed Brown        * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
7655b42dca8SJed Brown        * view every subdomain as though they were different. */
7665b42dca8SJed Brown       ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
767d5280255SMark F. Adams     }
768d5280255SMark F. Adams 
769d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
770d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
771e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
772d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
7733b4367a7SBarry Smith     if (mg->galerkin != 2) SETERRQ(comm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
774d5280255SMark F. Adams 
775d5280255SMark F. Adams     /* create cheby smoothers */
776d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
777d5280255SMark F. Adams          lidx <= fine_level;
778d5280255SMark F. Adams          lidx++, level--) {
779d5280255SMark F. Adams       KSP       smoother;
780890ffe84SMark Adams       PetscBool flag,flag2;
781d5280255SMark F. Adams       PC        subpc;
782d5280255SMark F. Adams 
783ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
784a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
785a2f3521dSMark F. Adams 
786ffc955d6SMark F. Adams       /* do my own cheby */
7876c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
788ffc955d6SMark F. Adams       if (flag) {
789ffc955d6SMark F. Adams         PetscReal emax, emin;
790251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
791890ffe84SMark Adams         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCSOR, &flag2);CHKERRQ(ierr);
792890ffe84SMark Adams         if ((flag||flag2) && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC but lets acccept SOR because it is close and safe (always lower) */
793890ffe84SMark Adams         else { /* eigen estimate 'emax' -- this is done in cheby */
794e696ed0bSMark F. Adams           KSP eksp;
795e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
796b2a4f308SMark F. Adams           Vec bb, xx;
797038e3b61SMark F. Adams 
7982a7a6963SBarry Smith           ierr = MatCreateVecs(Lmat, &bb, 0);CHKERRQ(ierr);
7992a7a6963SBarry Smith           ierr = MatCreateVecs(Lmat, &xx, 0);CHKERRQ(ierr);
800fc4362bfSMark F. Adams           {
801fc4362bfSMark F. Adams             PetscRandom rctx;
8023b4367a7SBarry Smith             ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
803fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
804fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
805fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
8065b89ad90SMark F. Adams           }
807a2f3521dSMark F. Adams 
808e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
809e696ed0bSMark F. Adams           {
810e696ed0bSMark F. Adams             PetscInt    Istart,Iend,ncols,jj,Ii;
811e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
812e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
813e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
814e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
815e696ed0bSMark F. Adams               if (ncols <= 1) {
816e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
817a94c3b12SMark F. Adams               }
818e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
819a94c3b12SMark F. Adams             }
820a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
821a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
822a94c3b12SMark F. Adams           }
823e696ed0bSMark F. Adams 
8243b4367a7SBarry Smith           ierr = KSPCreate(comm, &eksp);CHKERRQ(ierr);
825806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
826fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
8271a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
8281a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
829f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
830f6536408SMark F. Adams 
831f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
83223ee1639SBarry Smith           ierr = KSPSetOperators(eksp, Lmat, Lmat);CHKERRQ(ierr);
833fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
8345a9b9e01SMark F. Adams 
835d3d0db20SJed Brown           /* set PC type to be same as smoother */
836ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
837b2a4f308SMark F. Adams 
8385a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
8395a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
8405a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
841fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
8425a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
8435a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
8445a9b9e01SMark F. Adams 
845fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
8465a9b9e01SMark F. Adams 
847fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
848fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
849fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
850f6536408SMark F. Adams 
851ffc955d6SMark F. Adams           if (pc_gamg->verbose > 0) {
852a94c3b12SMark F. Adams             PetscInt N1, tt;
853a94c3b12SMark F. Adams             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
8543b4367a7SBarry Smith             PetscPrintf(comm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1);
855f6536408SMark F. Adams           }
856fc4362bfSMark F. Adams         }
857038e3b61SMark F. Adams         {
858c5bfad50SMark F. Adams           PetscInt N1, N0;
8590298fd71SBarry Smith           ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr);
8600298fd71SBarry Smith           ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr);
861f6536408SMark F. Adams           /* heuristic - is this crap? */
862b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
8635e7c91beSJed Brown           emin  = emax * pc_gamg->eigtarget[0];
8645e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
865038e3b61SMark F. Adams         }
8666c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
867ffc955d6SMark F. Adams       } /* setup checby flag */
868ffc955d6SMark F. Adams     } /* non-coarse levels */
869737a81a9SMark F. Adams 
870d5280255SMark F. Adams     /* clean up */
871d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
872587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
873587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
8745b89ad90SMark F. Adams     }
8750cbbd2e1SMark F. Adams 
8760cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
877806fa848SBarry Smith   } else {
8785f8cf99dSMark F. Adams     KSP smoother;
8793b4367a7SBarry Smith     if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",rank,__FUNCT__);
8809d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
88123ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
8825f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
8839d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8845f8cf99dSMark F. Adams   }
8855b89ad90SMark F. Adams   PetscFunctionReturn(0);
8865b89ad90SMark F. Adams }
8875b89ad90SMark F. Adams 
888eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8895b89ad90SMark F. Adams /*
8905b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8915b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8925b89ad90SMark F. Adams 
8935b89ad90SMark F. Adams    Input Parameter:
8945b89ad90SMark F. Adams .  pc - the preconditioner context
8955b89ad90SMark F. Adams 
8965b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8975b89ad90SMark F. Adams */
8985b89ad90SMark F. Adams #undef __FUNCT__
8995b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
9005b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
9015b89ad90SMark F. Adams {
9025b89ad90SMark F. Adams   PetscErrorCode ierr;
9035b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
9045b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
9055b89ad90SMark F. Adams 
9065b89ad90SMark F. Adams   PetscFunctionBegin;
9075b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
9089b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
9099b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
9109b8ffb57SJed Brown   }
9111ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
9121ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
9135b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
9145b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
9155b89ad90SMark F. Adams   PetscFunctionReturn(0);
9165b89ad90SMark F. Adams }
9175b89ad90SMark F. Adams 
918676e1743SMark F. Adams 
919676e1743SMark F. Adams #undef __FUNCT__
920676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
921676e1743SMark F. Adams /*@
922*1cc46a46SBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via processor reduction.
923676e1743SMark F. Adams 
924*1cc46a46SBarry Smith    Logically Collective on PC
925676e1743SMark F. Adams 
926676e1743SMark F. Adams    Input Parameters:
927*1cc46a46SBarry Smith +  pc - the preconditioner context
928*1cc46a46SBarry Smith -  n - the number of equations
929676e1743SMark F. Adams 
930676e1743SMark F. Adams 
931676e1743SMark F. Adams    Options Database Key:
932*1cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
933676e1743SMark F. Adams 
934676e1743SMark F. Adams    Level: intermediate
935676e1743SMark F. Adams 
936676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
937676e1743SMark F. Adams 
938676e1743SMark F. Adams .seealso: ()
939676e1743SMark F. Adams @*/
940676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
941676e1743SMark F. Adams {
942676e1743SMark F. Adams   PetscErrorCode ierr;
943676e1743SMark F. Adams 
944676e1743SMark F. Adams   PetscFunctionBegin;
945676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
946676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
947676e1743SMark F. Adams   PetscFunctionReturn(0);
948676e1743SMark F. Adams }
949676e1743SMark F. Adams 
950676e1743SMark F. Adams #undef __FUNCT__
951676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
9521e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
953676e1743SMark F. Adams {
954c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
955c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
956676e1743SMark F. Adams 
957676e1743SMark F. Adams   PetscFunctionBegin;
9589d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
959676e1743SMark F. Adams   PetscFunctionReturn(0);
960676e1743SMark F. Adams }
961676e1743SMark F. Adams 
962676e1743SMark F. Adams #undef __FUNCT__
963389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
964389730f3SMark F. Adams /*@
965389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
966389730f3SMark F. Adams 
967389730f3SMark F. Adams  Collective on PC
968389730f3SMark F. Adams 
969389730f3SMark F. Adams    Input Parameters:
970*1cc46a46SBarry Smith +  pc - the preconditioner context
971*1cc46a46SBarry Smith -  n - maximum number of equations to aim for
972389730f3SMark F. Adams 
973389730f3SMark F. Adams    Options Database Key:
974*1cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
975389730f3SMark F. Adams 
976389730f3SMark F. Adams    Level: intermediate
977389730f3SMark F. Adams 
978389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
979389730f3SMark F. Adams 
980389730f3SMark F. Adams @*/
981389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
982389730f3SMark F. Adams {
983389730f3SMark F. Adams   PetscErrorCode ierr;
984389730f3SMark F. Adams 
985389730f3SMark F. Adams   PetscFunctionBegin;
986389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
987389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
988389730f3SMark F. Adams   PetscFunctionReturn(0);
989389730f3SMark F. Adams }
990389730f3SMark F. Adams 
991389730f3SMark F. Adams #undef __FUNCT__
992389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
9931e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
994389730f3SMark F. Adams {
995389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
996389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
997389730f3SMark F. Adams 
998389730f3SMark F. Adams   PetscFunctionBegin;
9999d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
1000389730f3SMark F. Adams   PetscFunctionReturn(0);
1001389730f3SMark F. Adams }
1002389730f3SMark F. Adams 
1003389730f3SMark F. Adams #undef __FUNCT__
10048263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1005676e1743SMark F. Adams /*@
10068263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1007676e1743SMark F. Adams 
1008676e1743SMark F. Adams    Collective on PC
1009676e1743SMark F. Adams 
1010676e1743SMark F. Adams    Input Parameters:
1011*1cc46a46SBarry Smith +  pc - the preconditioner context
1012*1cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1013676e1743SMark F. Adams 
1014676e1743SMark F. Adams    Options Database Key:
1015*1cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
1016676e1743SMark F. Adams 
1017676e1743SMark F. Adams    Level: intermediate
1018676e1743SMark F. Adams 
1019676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1020676e1743SMark F. Adams 
1021676e1743SMark F. Adams .seealso: ()
1022676e1743SMark F. Adams @*/
10238263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1024676e1743SMark F. Adams {
1025676e1743SMark F. Adams   PetscErrorCode ierr;
1026676e1743SMark F. Adams 
1027676e1743SMark F. Adams   PetscFunctionBegin;
1028676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10298263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1030676e1743SMark F. Adams   PetscFunctionReturn(0);
1031676e1743SMark F. Adams }
1032676e1743SMark F. Adams 
1033676e1743SMark F. Adams #undef __FUNCT__
10348263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
10351e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1036676e1743SMark F. Adams {
1037c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1038c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1039676e1743SMark F. Adams 
1040676e1743SMark F. Adams   PetscFunctionBegin;
10419d5b6da9SMark F. Adams   pc_gamg->repart = n;
1042676e1743SMark F. Adams   PetscFunctionReturn(0);
1043676e1743SMark F. Adams }
1044676e1743SMark F. Adams 
1045676e1743SMark F. Adams #undef __FUNCT__
1046*1cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation"
1047dfd5c07aSMark F. Adams /*@
1048*1cc46a46SBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner
1049dfd5c07aSMark F. Adams 
1050dfd5c07aSMark F. Adams    Collective on PC
1051dfd5c07aSMark F. Adams 
1052dfd5c07aSMark F. Adams    Input Parameters:
1053*1cc46a46SBarry Smith +  pc - the preconditioner context
1054*1cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1055dfd5c07aSMark F. Adams 
1056dfd5c07aSMark F. Adams    Options Database Key:
1057*1cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
1058dfd5c07aSMark F. Adams 
1059dfd5c07aSMark F. Adams    Level: intermediate
1060dfd5c07aSMark F. Adams 
1061dfd5c07aSMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1062dfd5c07aSMark F. Adams 
1063dfd5c07aSMark F. Adams .seealso: ()
1064dfd5c07aSMark F. Adams @*/
1065*1cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1066dfd5c07aSMark F. Adams {
1067dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1068dfd5c07aSMark F. Adams 
1069dfd5c07aSMark F. Adams   PetscFunctionBegin;
1070dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1071*1cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1072dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1073dfd5c07aSMark F. Adams }
1074dfd5c07aSMark F. Adams 
1075dfd5c07aSMark F. Adams #undef __FUNCT__
1076*1cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG"
1077*1cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1078dfd5c07aSMark F. Adams {
1079dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1080dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1081dfd5c07aSMark F. Adams 
1082dfd5c07aSMark F. Adams   PetscFunctionBegin;
1083dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1084dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1085dfd5c07aSMark F. Adams }
1086dfd5c07aSMark F. Adams 
1087dfd5c07aSMark F. Adams #undef __FUNCT__
1088ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1089ffc955d6SMark F. Adams /*@
1090ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1091ffc955d6SMark F. Adams 
1092ffc955d6SMark F. Adams    Collective on PC
1093ffc955d6SMark F. Adams 
1094ffc955d6SMark F. Adams    Input Parameters:
1095ffc955d6SMark F. Adams .  pc - the preconditioner context
1096ffc955d6SMark F. Adams 
1097ffc955d6SMark F. Adams 
1098ffc955d6SMark F. Adams    Options Database Key:
1099ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1100ffc955d6SMark F. Adams 
1101ffc955d6SMark F. Adams    Level: intermediate
1102ffc955d6SMark F. Adams 
1103ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1104ffc955d6SMark F. Adams 
1105ffc955d6SMark F. Adams .seealso: ()
1106ffc955d6SMark F. Adams @*/
1107ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1108ffc955d6SMark F. Adams {
1109ffc955d6SMark F. Adams   PetscErrorCode ierr;
1110ffc955d6SMark F. Adams 
1111ffc955d6SMark F. Adams   PetscFunctionBegin;
1112ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1113ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1114ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1115ffc955d6SMark F. Adams }
1116ffc955d6SMark F. Adams 
1117ffc955d6SMark F. Adams #undef __FUNCT__
1118ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
11191e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1120ffc955d6SMark F. Adams {
1121ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1122ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1123ffc955d6SMark F. Adams 
1124ffc955d6SMark F. Adams   PetscFunctionBegin;
1125ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1126ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1127ffc955d6SMark F. Adams }
1128ffc955d6SMark F. Adams 
1129ffc955d6SMark F. Adams #undef __FUNCT__
11304ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
11314ef23d27SMark F. Adams /*@
1132*1cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
11334ef23d27SMark F. Adams 
11344ef23d27SMark F. Adams    Not collective on PC
11354ef23d27SMark F. Adams 
11364ef23d27SMark F. Adams    Input Parameters:
1137*1cc46a46SBarry Smith +  pc - the preconditioner
1138*1cc46a46SBarry Smith -  n - the maximum number of levels to use
11394ef23d27SMark F. Adams 
11404ef23d27SMark F. Adams    Options Database Key:
11414ef23d27SMark F. Adams .  -pc_mg_levels
11424ef23d27SMark F. Adams 
11434ef23d27SMark F. Adams    Level: intermediate
11444ef23d27SMark F. Adams 
11454ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
11464ef23d27SMark F. Adams 
11474ef23d27SMark F. Adams .seealso: ()
11484ef23d27SMark F. Adams @*/
11494ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
11504ef23d27SMark F. Adams {
11514ef23d27SMark F. Adams   PetscErrorCode ierr;
11524ef23d27SMark F. Adams 
11534ef23d27SMark F. Adams   PetscFunctionBegin;
11544ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11554ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
11564ef23d27SMark F. Adams   PetscFunctionReturn(0);
11574ef23d27SMark F. Adams }
11584ef23d27SMark F. Adams 
11594ef23d27SMark F. Adams #undef __FUNCT__
11604ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
11611e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
11624ef23d27SMark F. Adams {
11634ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
11644ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
11654ef23d27SMark F. Adams 
11664ef23d27SMark F. Adams   PetscFunctionBegin;
11679d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
11684ef23d27SMark F. Adams   PetscFunctionReturn(0);
11694ef23d27SMark F. Adams }
11704ef23d27SMark F. Adams 
11714ef23d27SMark F. Adams #undef __FUNCT__
11723542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
11733542efc5SMark F. Adams /*@
11743542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
11753542efc5SMark F. Adams 
11763542efc5SMark F. Adams    Not collective on PC
11773542efc5SMark F. Adams 
11783542efc5SMark F. Adams    Input Parameters:
1179*1cc46a46SBarry Smith +  pc - the preconditioner context
1180*1cc46a46SBarry Smith -  threshold - the threshold value, 0.0 is the lowest possible
11813542efc5SMark F. Adams 
11823542efc5SMark F. Adams    Options Database Key:
1183*1cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
11843542efc5SMark F. Adams 
11853542efc5SMark F. Adams    Level: intermediate
11863542efc5SMark F. Adams 
11873542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
11883542efc5SMark F. Adams 
11893542efc5SMark F. Adams .seealso: ()
11903542efc5SMark F. Adams @*/
11913542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
11923542efc5SMark F. Adams {
11933542efc5SMark F. Adams   PetscErrorCode ierr;
11943542efc5SMark F. Adams 
11953542efc5SMark F. Adams   PetscFunctionBegin;
11963542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11973542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
11983542efc5SMark F. Adams   PetscFunctionReturn(0);
11993542efc5SMark F. Adams }
12003542efc5SMark F. Adams 
12013542efc5SMark F. Adams #undef __FUNCT__
12023542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
12031e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
12043542efc5SMark F. Adams {
1205c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1206c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12073542efc5SMark F. Adams 
12083542efc5SMark F. Adams   PetscFunctionBegin;
12099d5b6da9SMark F. Adams   pc_gamg->threshold = n;
12103542efc5SMark F. Adams   PetscFunctionReturn(0);
12113542efc5SMark F. Adams }
12123542efc5SMark F. Adams 
12133542efc5SMark F. Adams #undef __FUNCT__
12149d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1215676e1743SMark F. Adams /*@
1216c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1217676e1743SMark F. Adams 
1218676e1743SMark F. Adams    Collective on PC
1219676e1743SMark F. Adams 
1220676e1743SMark F. Adams    Input Parameters:
1221c60c7ad4SBarry Smith +  pc - the preconditioner context
1222c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1223676e1743SMark F. Adams 
1224676e1743SMark F. Adams    Options Database Key:
1225c60c7ad4SBarry Smith .  -pc_gamg_type <agg,geo,classical>
1226676e1743SMark F. Adams 
1227676e1743SMark F. Adams    Level: intermediate
1228676e1743SMark F. Adams 
1229676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1230676e1743SMark F. Adams 
1231c60c7ad4SBarry Smith .seealso: PCGAMGGetType(), PCGAMG
1232676e1743SMark F. Adams @*/
123319fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1234676e1743SMark F. Adams {
1235676e1743SMark F. Adams   PetscErrorCode ierr;
1236676e1743SMark F. Adams 
1237676e1743SMark F. Adams   PetscFunctionBegin;
1238676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1239806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1240676e1743SMark F. Adams   PetscFunctionReturn(0);
1241676e1743SMark F. Adams }
1242676e1743SMark F. Adams 
1243676e1743SMark F. Adams #undef __FUNCT__
1244c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType"
1245c60c7ad4SBarry Smith /*@
1246c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1247c60c7ad4SBarry Smith 
1248c60c7ad4SBarry Smith    Collective on PC
1249c60c7ad4SBarry Smith 
1250c60c7ad4SBarry Smith    Input Parameter:
1251c60c7ad4SBarry Smith .  pc - the preconditioner context
1252c60c7ad4SBarry Smith 
1253c60c7ad4SBarry Smith    Output Parameter:
1254c60c7ad4SBarry Smith .  type - the type of algorithm used
1255c60c7ad4SBarry Smith 
1256c60c7ad4SBarry Smith    Level: intermediate
1257c60c7ad4SBarry Smith 
1258c60c7ad4SBarry Smith    Concepts: Unstructured multrigrid preconditioner
1259c60c7ad4SBarry Smith 
1260c60c7ad4SBarry Smith .seealso: PCGAMGSetType()
1261c60c7ad4SBarry Smith @*/
1262c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1263c60c7ad4SBarry Smith {
1264c60c7ad4SBarry Smith   PetscErrorCode ierr;
1265c60c7ad4SBarry Smith 
1266c60c7ad4SBarry Smith   PetscFunctionBegin;
1267c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1268c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1269c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1270c60c7ad4SBarry Smith }
1271c60c7ad4SBarry Smith 
1272c60c7ad4SBarry Smith #undef __FUNCT__
1273c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG"
1274c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1275c60c7ad4SBarry Smith {
1276c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1277c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1278c60c7ad4SBarry Smith 
1279c60c7ad4SBarry Smith   PetscFunctionBegin;
1280c60c7ad4SBarry Smith   *type = pc_gamg->type;
1281c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1282c60c7ad4SBarry Smith }
1283c60c7ad4SBarry Smith 
1284c60c7ad4SBarry Smith #undef __FUNCT__
12859d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
12861e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1287676e1743SMark F. Adams {
12889d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
12891ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
12901ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1291676e1743SMark F. Adams 
1292676e1743SMark F. Adams   PetscFunctionBegin;
1293c60c7ad4SBarry Smith   pc_gamg->type = type;
12941c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
12959d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
12961ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
12973ae0bb68SMark Adams     /* there was something here - kill it */
12981ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
12991ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1300e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
13013ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
13023ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
13033ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
13043ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
13053ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
13063ae0bb68SMark Adams     if (pc_gamg->data_sz) {
13073ae0bb68SMark Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
13083ae0bb68SMark Adams       pc_gamg->data_sz = 0;
1309c60c7ad4SBarry Smith     } else if (pc_gamg->data) {
13103ae0bb68SMark Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); /* can this happen ? */
13113ae0bb68SMark Adams     }
13121ab5ffc9SJed Brown   }
13131ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
13141ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
13159d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1316676e1743SMark F. Adams   PetscFunctionReturn(0);
1317676e1743SMark F. Adams }
1318676e1743SMark F. Adams 
13195b89ad90SMark F. Adams #undef __FUNCT__
13205b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
13218c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptions *PetscOptionsObject,PC pc)
13225b89ad90SMark F. Adams {
1323676e1743SMark F. Adams   PetscErrorCode ierr;
1324676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1325676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1326676e1743SMark F. Adams   PetscBool      flag;
13275e7c91beSJed Brown   PetscInt       two   = 2;
13283b4367a7SBarry Smith   MPI_Comm       comm;
13295b89ad90SMark F. Adams 
13305b89ad90SMark F. Adams   PetscFunctionBegin;
13313b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1332e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1333676e1743SMark F. Adams   {
1334b7cbab4eSMark Adams     /* -pc_gamg_type */
1335b7cbab4eSMark Adams     {
1336bd94a7aaSJed Brown       char tname[256];
13371a1c1e04SBarry Smith       ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1338bd94a7aaSJed Brown       if (flag) {
1339bd94a7aaSJed Brown         ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
13401ab5ffc9SJed Brown       }
1341b7cbab4eSMark Adams     }
134275b74bdaSMark F. Adams     /* -pc_gamg_verbose */
134394ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG","none", pc_gamg->verbose,&pc_gamg->verbose, NULL);CHKERRQ(ierr);
13448263b398SMark F. Adams     /* -pc_gamg_repartition */
134594ae4db5SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
1346dfd5c07aSMark F. Adams     /* -pc_gamg_reuse_interpolation */
1347*1cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1348ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
134994ae4db5SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm","Use aggregation agragates for GASM smoother","PCGAMGUseASMAggs",pc_gamg->use_aggs_in_gasm,&pc_gamg->use_aggs_in_gasm,NULL);CHKERRQ(ierr);
1350c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
135194ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit","Limit (goal) on number of equations per process on coarse grids","PCGAMGSetProcEqLim",pc_gamg->min_eq_proc,&pc_gamg->min_eq_proc,NULL);CHKERRQ(ierr);
1352389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
135394ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit","Limit on number of equations for the coarse grid","PCGAMGSetCoarseEqLim",pc_gamg->coarse_eq_limit,&pc_gamg->coarse_eq_limit,NULL);CHKERRQ(ierr);
1354c20e4228SMark F. Adams     /* -pc_gamg_threshold */
135594ae4db5SBarry Smith     ierr = PetscOptionsReal("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&pc_gamg->threshold,&flag);CHKERRQ(ierr);
1356806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
13573b4367a7SBarry Smith       ierr = PetscPrintf(comm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1358806fa848SBarry Smith     }
1359b7cbab4eSMark Adams     /* -pc_gamg_eigtarget */
13600298fd71SBarry Smith     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr);
136194ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1362b7cbab4eSMark Adams 
1363b7cbab4eSMark Adams     /* set options for subtype */
1364e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1365676e1743SMark F. Adams   }
1366676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
13675b89ad90SMark F. Adams   PetscFunctionReturn(0);
13685b89ad90SMark F. Adams }
13695b89ad90SMark F. Adams 
13705b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
13715b89ad90SMark F. Adams /*MC
1372*1cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
13735b89ad90SMark F. Adams 
1374280d9858SJed Brown    Options Database Keys:
13755b89ad90SMark F. Adams    Multigrid options(inherited)
1376*1cc46a46SBarry Smith +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1377280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1378280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
13798c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
13805b89ad90SMark F. Adams 
1381*1cc46a46SBarry Smith 
1382*1cc46a46SBarry Smith   Notes: In order to obtain good performance for PCGAMG for vector valued problems you must
1383*1cc46a46SBarry Smith $       Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1384*1cc46a46SBarry Smith $       Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1385*1cc46a46SBarry Smith $       See the Users Manual Chapter 4 for more details
1386*1cc46a46SBarry Smith 
13875b89ad90SMark F. Adams   Level: intermediate
1388280d9858SJed Brown 
1389*1cc46a46SBarry Smith   Concepts: algebraic multigrid
13905b89ad90SMark F. Adams 
1391*1cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1392*1cc46a46SBarry Smith            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartitioning(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseASMAggs(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType()
13935b89ad90SMark F. Adams M*/
1394b2573a8aSBarry Smith 
13955b89ad90SMark F. Adams #undef __FUNCT__
13965b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
13978cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
13985b89ad90SMark F. Adams {
13995b89ad90SMark F. Adams   PetscErrorCode ierr;
14005b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
14015b89ad90SMark F. Adams   PC_MG          *mg;
14020cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14039d5b6da9SMark F. Adams   static long count = 0;
14045ee9c036SSatish Balay #endif
14055b89ad90SMark F. Adams 
14065b89ad90SMark F. Adams   PetscFunctionBegin;
14075b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14085b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
14095b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14105b89ad90SMark F. Adams 
14115b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1412b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
14135b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
1414ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
14155b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
14165b89ad90SMark F. Adams 
1417b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
14181ab5ffc9SJed Brown 
14199d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14209d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14219d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14229d5b6da9SMark F. Adams   pc_gamg->data    = 0;
14235b89ad90SMark F. Adams 
14249d5b6da9SMark F. Adams   /* register AMG type */
14253e3471ccSMark Adams   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
14269d5b6da9SMark F. Adams 
14279d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
14285b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
14295b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
14305b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
14315b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
14325b89ad90SMark F. Adams 
1433bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1434bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1435bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
1436*1cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1437bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1438bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1439bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1440c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1441bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
14429d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1443d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
1444ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1445038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
144625a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1447d3042614SMark Adams   pc_gamg->threshold        = 0.;
14489d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
14499d5b6da9SMark F. Adams   pc_gamg->verbose          = 0;
14509d5b6da9SMark F. Adams   pc_gamg->emax_id          = -1;
14515e7c91beSJed Brown   pc_gamg->eigtarget[0]     = 0.05;
14525e7c91beSJed Brown   pc_gamg->eigtarget[1]     = 1.05;
1453c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
14549d5b6da9SMark F. Adams 
14550cbbd2e1SMark F. Adams   /* private events */
14560cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1457785cba28SMark F. Adams   if (count++ == 0) {
1458806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1459806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
14600cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
14610cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
14620cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1463806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1464806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1465806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1466806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1467806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1468806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1469806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1470806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1471806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1472806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1473806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1474806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1475f852f58cSMark F. Adams 
14760cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
14770cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
14780cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1479b4fbaa2aSMark F. Adams     /* create timer stages */
1480b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1481b4fbaa2aSMark F. Adams     {
1482b4fbaa2aSMark F. Adams       char     str[32];
1483b4fbaa2aSMark F. Adams       PetscInt lidx;
1484806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1485806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1486b4fbaa2aSMark F. Adams       for (lidx=1; lidx<9; lidx++) {
1487b4fbaa2aSMark F. Adams         sprintf(str,"MG Level %d",lidx);
1488806fa848SBarry Smith         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1489b4fbaa2aSMark F. Adams       }
1490b4fbaa2aSMark F. Adams     }
1491b4fbaa2aSMark F. Adams #endif
1492b4fbaa2aSMark F. Adams   }
1493b4fbaa2aSMark F. Adams #endif
1494bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1495bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
14965b89ad90SMark F. Adams   PetscFunctionReturn(0);
14975b89ad90SMark F. Adams }
14983e3471ccSMark Adams 
14993e3471ccSMark Adams #undef __FUNCT__
15003e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage"
15013e3471ccSMark Adams /*@C
15023e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
15033e3471ccSMark Adams  from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
15043e3471ccSMark Adams  when using static libraries.
15053e3471ccSMark Adams 
15063e3471ccSMark Adams  Level: developer
15073e3471ccSMark Adams 
15083e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
15093e3471ccSMark Adams  .seealso: PetscInitialize()
15103e3471ccSMark Adams @*/
15113e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
15123e3471ccSMark Adams {
15133e3471ccSMark Adams   PetscErrorCode ierr;
15143e3471ccSMark Adams 
15153e3471ccSMark Adams   PetscFunctionBegin;
15163e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
15173e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
15183e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
15193e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
15208e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
15213e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1522c1c463dbSMark Adams 
1523c1c463dbSMark Adams   /* general events */
1524c1c463dbSMark Adams #if defined PETSC_USE_LOG
1525c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1526c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1527c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1528c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1529c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1530c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1531c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1532c1c463dbSMark Adams #endif
1533c1c463dbSMark Adams 
15343e3471ccSMark Adams   PetscFunctionReturn(0);
15353e3471ccSMark Adams }
15363e3471ccSMark Adams 
15373e3471ccSMark Adams #undef __FUNCT__
15383e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage"
15393e3471ccSMark Adams /*@C
15403e3471ccSMark Adams  PCGAMGFinalizePackage - This function destroys everything in the PCGAMG package. It is
15413e3471ccSMark Adams  called from PetscFinalize().
15423e3471ccSMark Adams 
15433e3471ccSMark Adams  Level: developer
15443e3471ccSMark Adams 
15453e3471ccSMark Adams  .keywords: Petsc, destroy, package
15463e3471ccSMark Adams  .seealso: PetscFinalize()
15473e3471ccSMark Adams @*/
15483e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
15493e3471ccSMark Adams {
15503e3471ccSMark Adams   PetscErrorCode ierr;
15513e3471ccSMark Adams 
15523e3471ccSMark Adams   PetscFunctionBegin;
15533e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
15543e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
15553e3471ccSMark Adams   PetscFunctionReturn(0);
15563e3471ccSMark Adams }
1557a36cf38bSToby Isaac 
1558a36cf38bSToby Isaac #undef __FUNCT__
1559a36cf38bSToby Isaac #define __FUNCT__ "PCGAMGRegister"
1560a36cf38bSToby Isaac /*@C
1561a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1562a36cf38bSToby Isaac 
1563a36cf38bSToby Isaac  Input Parameters:
1564a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1565a36cf38bSToby Isaac  - create - function for creating the gamg context.
1566a36cf38bSToby Isaac 
1567a36cf38bSToby Isaac   Level: advanced
1568a36cf38bSToby Isaac 
1569a36cf38bSToby Isaac  .seealso: PCGAMGGetContext(), PCGAMGSetContext()
1570a36cf38bSToby Isaac @*/
1571a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1572a36cf38bSToby Isaac {
1573a36cf38bSToby Isaac   PetscErrorCode ierr;
1574a36cf38bSToby Isaac 
1575a36cf38bSToby Isaac   PetscFunctionBegin;
1576a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1577a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1578a36cf38bSToby Isaac   PetscFunctionReturn(0);
1579a36cf38bSToby Isaac }
1580a36cf38bSToby Isaac 
1581