xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 3cb8563fbfdf202cd42eb0238066e5df626afb3a)
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,
75*3cb8563fSToby Isaac                                   Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,
76*3cb8563fSToby Isaac                                   IS * Pcolumnperm)
775b89ad90SMark F. Adams {
78a2f3521dSMark F. Adams   PetscErrorCode  ierr;
799d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
80486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
819d5b6da9SMark F. Adams   const PetscBool repart      = pc_gamg->repart;
829d5b6da9SMark F. Adams   const PetscInt  min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
83a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
843b4367a7SBarry Smith   MPI_Comm        comm;
85c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
863ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
875b89ad90SMark F. Adams 
885b89ad90SMark F. Adams   PetscFunctionBegin;
893b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
903b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
913b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
92c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
9311e60469SMark F. Adams   /* RAP */
949d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
95038e3b61SMark F. Adams 
963ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
970298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
983ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
993ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
1003ae0bb68SMark Adams   }
1013ae0bb68SMark Adams   else {
1023ae0bb68SMark Adams     PetscInt  bs;
1033ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
1043ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
1053ae0bb68SMark Adams   }
106a2f3521dSMark F. Adams 
107c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
108a2f3521dSMark F. Adams   {
109472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
1100298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
111c5df96a5SBarry Smith     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
112c5df96a5SBarry Smith     if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1;
113c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
114a2f3521dSMark F. Adams   }
115f852f58cSMark F. Adams 
116*3cb8563fSToby Isaac   if (Pcolumnperm) *Pcolumnperm = NULL;
117*3cb8563fSToby Isaac 
1182fa5cd67SKarl Rupp   if (!repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1192fa5cd67SKarl Rupp   else {
1203ae0bb68SMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
121885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
122e33ef3b1SMark F. Adams 
12371959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
12471959b99SBarry 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);
1250cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1260cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
127b4fbaa2aSMark F. Adams #endif
128a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
129785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
1307700e67bSMark Adams     if (repart) {
131a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1325a9b9e01SMark F. Adams       Mat adj;
1335a9b9e01SMark F. Adams 
134a2f3521dSMark F. Adams       if (pc_gamg->verbose>0) {
1353b4367a7SBarry 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);
136a2f3521dSMark F. Adams         else {
137a2f3521dSMark F. Adams           PetscInt n;
1383b4367a7SBarry Smith           ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
1393b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n);
140a2f3521dSMark F. Adams         }
141a2f3521dSMark F. Adams       }
1425a9b9e01SMark F. Adams 
143a2f3521dSMark F. Adams       /* get 'adj' */
144c5bfad50SMark F. Adams       if (cr_bs == 1) {
145038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
146806fa848SBarry Smith       } else {
147a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
148eb07cef2SMark F. Adams         Mat               tMat;
149a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
150b4fbaa2aSMark F. Adams         const PetscScalar *vals;
151b4fbaa2aSMark F. Adams         const PetscInt    *idx;
152a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1539057884aSMark F. Adams         static PetscInt   llev = 0;
154b4fbaa2aSMark F. Adams 
155578f55a3SPeter Brune         ierr = PetscMalloc1(ncrs, &d_nnz);CHKERRQ(ierr);
156578f55a3SPeter Brune         ierr = PetscMalloc1(ncrs, &o_nnz);CHKERRQ(ierr);
157a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
158a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
159c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
16058471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
161c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
162c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
16358471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1643ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1653ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
16658471d46SMark F. Adams         }
1676876a03eSMark F. Adams 
1683b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1693ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
170a2f3521dSMark F. Adams         ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr);
171a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
172a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
17358471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
1745f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
175eb07cef2SMark F. Adams 
176a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
177c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
17822063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
179eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
180c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
181eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
182eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
183eb07cef2SMark F. Adams           }
18422063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
185eb07cef2SMark F. Adams         }
186eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
187eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
188eb07cef2SMark F. Adams 
189b4fbaa2aSMark F. Adams         if (llev++ == -1) {
190b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
1918caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
1923b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
193b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
1943bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
195b4fbaa2aSMark F. Adams         }
196b4fbaa2aSMark F. Adams 
197eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
198eb07cef2SMark F. Adams 
199eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
200a2f3521dSMark F. Adams       } /* create 'adj' */
201f150b916SMark F. Adams 
202a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2035a9b9e01SMark F. Adams         char            prefix[256];
2045a9b9e01SMark F. Adams         const char      *pcpre;
205b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
206b4fbaa2aSMark F. Adams         MatPartitioning mpart;
207a4b7d37bSMark F. Adams         IS              proc_is;
208a2f3521dSMark F. Adams         PetscInt        targetPE;
2092f03bc48SMark F. Adams 
2103b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2115ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2129d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2138caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
21459a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
21511e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
216c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
217a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
21811e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2195a9b9e01SMark F. Adams 
2205ef31b24SMark F. Adams         /* collect IS info */
221785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
222a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
223a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
224c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
225a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
226c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
227a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
228eb07cef2SMark F. Adams           }
2295ef31b24SMark F. Adams         }
230a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
231a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2325ef31b24SMark F. Adams       }
2335ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2345a9b9e01SMark F. Adams 
2353b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2368263b398SMark F. Adams       if (newproc_idx != 0) {
2378263b398SMark F. Adams         ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
2385ef31b24SMark F. Adams       }
239806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
240a2f3521dSMark F. Adams 
241a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
2425a9b9e01SMark F. Adams       /* find factor */
243c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2445a9b9e01SMark F. Adams       else {
2455a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2465a9b9e01SMark F. Adams         jj = -1;
247c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
248c5df96a5SBarry Smith           if (size%kk==0) { /* a candidate */
249c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2505a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2515a9b9e01SMark F. Adams             if (fact > best_fact) {
2525a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
2535a9b9e01SMark F. Adams             }
2545a9b9e01SMark F. Adams           }
2555a9b9e01SMark F. Adams         }
2565a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
257a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
2585a9b9e01SMark F. Adams       }
259c5df96a5SBarry Smith       new_size = size/rfactor;
2605a9b9e01SMark F. Adams 
261c5df96a5SBarry Smith       if (new_size==nactive) {
262a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2635a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
264a2f3521dSMark F. Adams         if (pc_gamg->verbose>0) {
2653b4367a7SBarry Smith           PetscPrintf(comm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq);
266a2f3521dSMark F. Adams         }
2675a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2685a9b9e01SMark F. Adams       }
2695a9b9e01SMark F. Adams 
2703b4367a7SBarry Smith       if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq);
271c5df96a5SBarry Smith       targetPE = rank/rfactor;
2723b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
273a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
274e33ef3b1SMark F. Adams 
27511e60469SMark F. Adams     /*
276a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
27711e60469SMark F. Adams      */
278a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2797700e67bSMark Adams     is_eq_num_prim = is_eq_num;
28011e60469SMark F. Adams     /*
281a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
28211e60469SMark F. Adams      */
283c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
284c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
285a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
2863ae0bb68SMark Adams     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */
287a2f3521dSMark F. Adams 
288a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2890cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2900cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
291b4fbaa2aSMark F. Adams #endif
292885364a3SMark 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 */
293885364a3SMark Adams     {
294885364a3SMark Adams     Vec            src_crd, dest_crd;
295885364a3SMark 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;
296885364a3SMark Adams     VecScatter     vecscat;
297885364a3SMark Adams     PetscScalar    *array;
298885364a3SMark Adams     IS isscat;
299a2f3521dSMark F. Adams 
300a2f3521dSMark F. Adams     /* move data (for primal equations only) */
30122063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
3023b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
3033ae0bb68SMark Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
304c0dedaeaSBarry Smith     ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
30511e60469SMark F. Adams     /*
3069d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
307c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
30811e60469SMark F. Adams      */
309854ce69bSBarry Smith     ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
310a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3113ae0bb68SMark Adams     for (ii=0,jj=0; ii<ncrs; ii++) {
312c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
313a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
31411e60469SMark F. Adams     }
315a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3163ae0bb68SMark Adams     ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
31792a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
31811e60469SMark F. Adams     /*
31911e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
32011e60469SMark F. Adams      */
3213ae0bb68SMark Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
3229d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3233ae0bb68SMark Adams       const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3243ae0bb68SMark Adams       for (ii=0; ii<ncrs; ii++) {
3259d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
326a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
327c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
328676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
329d3d6bff4SMark F. Adams         }
330038e3b61SMark F. Adams       }
331eb07cef2SMark F. Adams     }
332eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
333eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
33411e60469SMark F. Adams     /*
33511e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
33611e60469SMark F. Adams       to the correct processor
33711e60469SMark F. Adams     */
3380298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
33911e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
34011e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
34111e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
34211e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
34311e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
34411e60469SMark F. Adams     /*
34511e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
34611e60469SMark F. Adams     */
347c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
348578f55a3SPeter Brune     ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3492fa5cd67SKarl Rupp 
3503ae0bb68SMark Adams     pc_gamg->data_sz = node_data_sz*ncrs_new;
3513ae0bb68SMark Adams     strideNew        = ncrs_new*ndata_rows;
3522fa5cd67SKarl Rupp 
35311e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3549d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3553ae0bb68SMark Adams       for (ii=0; ii<ncrs_new; ii++) {
3569d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
357a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
358c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
359d3d6bff4SMark F. Adams         }
360038e3b61SMark F. Adams       }
361038e3b61SMark F. Adams     }
36211e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
36311e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
364885364a3SMark Adams     }
365a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3660cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3670cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
368ed3f9983SMark F. Adams #endif
369a2f3521dSMark F. Adams 
37011e60469SMark F. Adams     /*
37111e60469SMark F. Adams       Invert for MatGetSubMatrix
37211e60469SMark F. Adams     */
373a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
374a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
375c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
376a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
377a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
378a2f3521dSMark F. Adams     }
379*3cb8563fSToby Isaac     if (Pcolumnperm) {
380*3cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
381*3cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
382*3cb8563fSToby Isaac     }
383a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3840cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3850cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3860cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
387ed3f9983SMark F. Adams #endif
388a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
389a2f3521dSMark F. Adams     {
390a2f3521dSMark F. Adams       Mat mat;
391806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
392a2f3521dSMark F. Adams       *a_Amat_crs = mat;
393c5bfad50SMark F. Adams 
394c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
395c5bfad50SMark F. Adams         PetscInt cbs, rbs;
396c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
397c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
398c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
399c5df96a5SBarry 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);
400c5bfad50SMark F. Adams       }
401a2f3521dSMark F. Adams     }
402038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
403a2f3521dSMark F. Adams 
4040cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4050cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
406ed3f9983SMark F. Adams #endif
40711e60469SMark F. Adams     /* prolongator */
40811e60469SMark F. Adams     {
40911e60469SMark F. Adams       IS       findices;
410a2f3521dSMark F. Adams       PetscInt Istart,Iend;
411a2f3521dSMark F. Adams       Mat      Pnew;
412a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4140cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
415ed3f9983SMark F. Adams #endif
4163b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
417c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
418806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
41911e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
420c5bfad50SMark F. Adams 
421c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
422c5bfad50SMark F. Adams         PetscInt cbs, rbs;
423c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
424c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
425c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
426c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
427c5bfad50SMark F. Adams       }
4280cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4290cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
430ed3f9983SMark F. Adams #endif
4313530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
432a2f3521dSMark F. Adams 
433a2f3521dSMark F. Adams       /* output - repartitioned */
434a2f3521dSMark F. Adams       *a_P_inout = Pnew;
435e33ef3b1SMark F. Adams     }
436a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4375b89ad90SMark F. Adams 
438c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
439a2f3521dSMark F. Adams   }
4405a9b9e01SMark F. Adams 
441a2f3521dSMark F. Adams   /* outout matrix data */
442c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
443c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
444c8b0795cSMark F. Adams     if (llev==0) {
445c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
4463b4367a7SBarry Smith       PetscViewerASCIIOpen(comm,fname,&viewer);
447c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
448c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
449c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
450c8b0795cSMark F. Adams     }
451c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
4523b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
453c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
454c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
455c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
456c8b0795cSMark F. Adams   }
4575b89ad90SMark F. Adams   PetscFunctionReturn(0);
4585b89ad90SMark F. Adams }
4595b89ad90SMark F. Adams 
4605b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4615b89ad90SMark F. Adams /*
4625b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4635b89ad90SMark F. Adams                     by setting data structures and options.
4645b89ad90SMark F. Adams 
4655b89ad90SMark F. Adams    Input Parameter:
4665b89ad90SMark F. Adams .  pc - the preconditioner context
4675b89ad90SMark F. Adams 
4685b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
4695b89ad90SMark F. Adams 
4705b89ad90SMark F. Adams    Notes:
4715b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
4725b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
4735b89ad90SMark F. Adams */
4745b89ad90SMark F. Adams #undef __FUNCT__
4755b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4769d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4775b89ad90SMark F. Adams {
4785b89ad90SMark F. Adams   PetscErrorCode ierr;
4799d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4805b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4812adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
482a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
4833b4367a7SBarry Smith   MPI_Comm       comm;
484c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
485587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
486c8b0795cSMark F. Adams   PetscReal      emaxs[GAMG_MAXLEVELS];
487e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
488a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
489737a81a9SMark F. Adams   MatInfo        info;
4907700e67bSMark Adams   PetscBool      redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);
4915ef31b24SMark F. Adams 
4925b89ad90SMark F. Adams   PetscFunctionBegin;
4933b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4943b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4953b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
496dfd5c07aSMark F. Adams 
4973b4367a7SBarry 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);
498dfd5c07aSMark F. Adams 
49984d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
500878e152fSMark F. Adams     if (redo_mesh_setup) {
501878e152fSMark F. Adams       /* reset everything */
502878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
503878e152fSMark F. Adams       pc->setupcalled = 0;
504806fa848SBarry Smith     } else {
50584d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
50603a628feSMark F. Adams       /* just do Galerkin grids */
50758471d46SMark F. Adams       Mat          B,dA,dB;
50858471d46SMark F. Adams 
50971959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
5109d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
51158471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
51223ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
51358471d46SMark F. Adams         /* (re)set to get dirty flag */
51423ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
51558471d46SMark F. Adams 
5162fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
51703a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
5180a97e771SToby Isaac           if (pc_gamg->setup_count==2) {
51903a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
520084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
5212fa5cd67SKarl Rupp 
52203a628feSMark F. Adams             mglevels[level]->A = B;
523806fa848SBarry Smith           } else {
52423ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
52558471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
52603a628feSMark F. Adams           }
52723ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
52858471d46SMark F. Adams           dB   = B;
52958471d46SMark F. Adams         }
5305f8cf99dSMark F. Adams       }
531d5280255SMark F. Adams 
532d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
533d5280255SMark F. Adams 
53458471d46SMark F. Adams       PetscFunctionReturn(0);
535eb07cef2SMark F. Adams     }
536878e152fSMark F. Adams   }
537f6536408SMark F. Adams 
538878e152fSMark F. Adams   if (!pc_gamg->data) {
539878e152fSMark F. Adams     if (pc_gamg->orig_data) {
540878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5410298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5422fa5cd67SKarl Rupp 
543878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
544878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
545878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5462fa5cd67SKarl Rupp 
547785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
548878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
549806fa848SBarry Smith     } else {
5501ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5517700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5529d5b6da9SMark F. Adams     }
553878e152fSMark F. Adams   }
554878e152fSMark F. Adams 
555878e152fSMark F. Adams   /* cache original data for reuse */
556878e152fSMark F. Adams   if (!pc_gamg->orig_data && redo_mesh_setup) {
557785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
558878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
559878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
560878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
561878e152fSMark F. Adams   }
562038e3b61SMark F. Adams 
563302f38e8SMark F. Adams   /* get basic dims */
564302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
565a2f3521dSMark F. Adams 
566a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
567c8b0795cSMark F. Adams   if (pc_gamg->verbose) {
56884f9421dSMark F. Adams     PetscInt NN = M;
56984f9421dSMark F. Adams     if (pc_gamg->verbose==1) {
57084f9421dSMark F. Adams       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
5713bf036e2SBarry Smith       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
572fe06f982SMark Adams       if (!NN) NN=1;
573806fa848SBarry Smith     } else {
574806fa848SBarry Smith       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
57584f9421dSMark F. Adams     }
576b2a4f308SMark F. Adams     nnz0   = info.nz_used;
577b2a4f308SMark F. Adams     nnztot = info.nz_used;
5783b4367a7SBarry 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",
579c5df96a5SBarry Smith                          rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
580c5df96a5SBarry Smith                          (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr);
581c8b0795cSMark F. Adams   }
58284d3f75bSMark F. Adams 
583a2f3521dSMark F. Adams   /* Get A_i and R_i */
584c5df96a5SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
585c5df96a5SBarry Smith        level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit);  /* && (size==1 || nactivepe>1); */
5860205a208SMark F. Adams        level++) {
58757d29afaSToby Isaac     pc_gamg->firstCoarsen = (level ? PETSC_FALSE : PETSC_TRUE);
5885b89ad90SMark F. Adams     level1 = level + 1;
5890cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5900cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
591a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
592a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
593b4fbaa2aSMark F. Adams #endif
594a2f3521dSMark F. Adams #endif
595c8b0795cSMark F. Adams     { /* construct prolongator */
596725b86d8SJed Brown       Mat              Gmat;
5970cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5987700e67bSMark Adams       Mat              Prol11;
599c8b0795cSMark F. Adams 
6007700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
6011ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
6027700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
603c8b0795cSMark F. Adams 
604a2f3521dSMark F. Adams       /* could have failed to create new level */
605a2f3521dSMark F. Adams       if (Prol11) {
6069d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6070298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
608a2f3521dSMark F. Adams 
6091ab5ffc9SJed Brown         if (pc_gamg->ops->optprol) {
610c8b0795cSMark F. Adams           /* smooth */
6117700e67bSMark Adams           ierr = pc_gamg->ops->optprol(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
612c8b0795cSMark F. Adams         }
613c8b0795cSMark F. Adams 
6147700e67bSMark Adams         Parr[level1] = Prol11;
6150298fd71SBarry Smith       } else Parr[level1] = NULL;
616ffc955d6SMark F. Adams 
617ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
6181b18a24aSMark Adams         PetscInt bs;
6191b18a24aSMark Adams         ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
620806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
621ffc955d6SMark F. Adams       }
622ffc955d6SMark F. Adams 
623a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
62441b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
625a2f3521dSMark F. Adams     } /* construct prolongator scope */
6260cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6270cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
628c8b0795cSMark F. Adams #endif
6299d5b6da9SMark F. Adams     /* cache eigen estimate */
6309d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
6319d5b6da9SMark F. Adams       PetscBool flag;
6327700e67bSMark Adams       ierr = PetscObjectComposedDataGetReal((PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
6339d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
634806fa848SBarry Smith     } else emaxs[level] = -1.;
6352adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
636c8b0795cSMark F. Adams     if (!Parr[level1]) {
637806fa848SBarry Smith       if (pc_gamg->verbose) {
6383b4367a7SBarry Smith         ierr =  PetscPrintf(comm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr);
639806fa848SBarry Smith       }
640c8b0795cSMark F. Adams       break;
641c8b0795cSMark F. Adams     }
6420cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6430cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
644b4fbaa2aSMark F. Adams #endif
645a2f3521dSMark F. Adams 
646c238b0ebSToby Isaac     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs,
647*3cb8563fSToby Isaac                        &Parr[level1], &Aarr[level1], &nactivepe, NULL);CHKERRQ(ierr);
648a2f3521dSMark F. Adams 
6490cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6500cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
651b4fbaa2aSMark F. Adams #endif
652a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
653a2f3521dSMark F. Adams 
654a2f3521dSMark F. Adams     if (pc_gamg->verbose > 0) {
6550cbbd2e1SMark F. Adams       PetscInt NN = M;
6560cbbd2e1SMark F. Adams       if (pc_gamg->verbose==1) {
657a2f3521dSMark F. Adams         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
6583bf036e2SBarry Smith         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
659fe06f982SMark Adams         if (!NN) NN=1;
660806fa848SBarry Smith       } else {
661806fa848SBarry Smith         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
6620cbbd2e1SMark F. Adams       }
663a2f3521dSMark F. Adams 
6640cbbd2e1SMark F. Adams       nnztot += info.nz_used;
6653b4367a7SBarry Smith       ierr    = PetscPrintf(comm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
666c5df96a5SBarry Smith                             rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
667806fa848SBarry Smith                             (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
668c8b0795cSMark F. Adams     }
669a2f3521dSMark F. Adams 
67025a145a7SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
67125a145a7SMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M < 2)
67225a145a7SMark Adams          || nactivepe==1 ) {
67381708dccSMark F. Adams       level++;
67481708dccSMark F. Adams       break;
67581708dccSMark F. Adams     }
6760cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
677b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
678b4fbaa2aSMark F. Adams #endif
679c8b0795cSMark F. Adams   } /* levels */
68057d29afaSToby Isaac   pc_gamg->firstCoarsen = PETSC_FALSE;
681c8b0795cSMark F. Adams 
682c8b0795cSMark F. Adams   if (pc_gamg->data) {
683c8b0795cSMark F. Adams     ierr          = PetscFree(pc_gamg->data);CHKERRQ(ierr);
6840298fd71SBarry Smith     pc_gamg->data = NULL;
6855b89ad90SMark F. Adams   }
686c8b0795cSMark F. Adams 
6873b4367a7SBarry Smith   if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
6889d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6895b89ad90SMark F. Adams   fine_level       = level;
6900298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
6915b89ad90SMark F. Adams 
69284d3f75bSMark F. Adams   /* simple setup */
69384d3f75bSMark F. Adams   if (!PETSC_TRUE) {
69484d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
69584d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
69684d3f75bSMark F. Adams          lidx<fine_level;
69784d3f75bSMark F. Adams          lidx++, level--) {
69884d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
69923ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level]);CHKERRQ(ierr);
70084d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
70184d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
70284d3f75bSMark F. Adams     }
70323ee1639SBarry Smith     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0]);CHKERRQ(ierr);
70484d3f75bSMark F. Adams 
70584d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
706806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
707d5280255SMark F. Adams     /* set default smoothers & set operators */
7089d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
709587fa25dSMark F. Adams          lidx <= fine_level;
710587fa25dSMark F. Adams          lidx++, level--) {
711ffc955d6SMark F. Adams       KSP smoother;
712ffc955d6SMark F. Adams       PC  subpc;
713a2f3521dSMark F. Adams 
7149d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
715f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
716ffc955d6SMark F. Adams 
717a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
718a2f3521dSMark F. Adams       /* set ops */
71923ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
720a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
721a2f3521dSMark F. Adams 
722a2f3521dSMark F. Adams       /* set defaults */
7236c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
724a2f3521dSMark F. Adams 
7251b18a24aSMark Adams       /* set blocks for GASM smoother that uses the 'aggregates' */
726ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
7272d3561bbSSatish Balay         PetscInt sz;
7282d3561bbSSatish Balay         IS       *is;
729a2f3521dSMark F. Adams 
7302d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
7312d3561bbSSatish Balay         is   = ASMLocalIDsArr[level];
732ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
7331b18a24aSMark Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
734ffc955d6SMark F. Adams         if (sz==0) {
735ffc955d6SMark F. Adams           IS       is;
736ffc955d6SMark F. Adams           PetscInt my0,kk;
737ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
738ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
7390298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr);
740a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
741806fa848SBarry Smith         } else {
742a94c3b12SMark F. Adams           PetscInt kk;
7430298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr);
744a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
745a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
746a94c3b12SMark F. Adams           }
747ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
748ffc955d6SMark F. Adams         }
7490298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
750ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
751ffc955d6SMark F. Adams         ierr                  = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
752806fa848SBarry Smith       } else {
753890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
754ffc955d6SMark F. Adams       }
755d5280255SMark F. Adams     }
756d5280255SMark F. Adams     {
757d5280255SMark F. Adams       /* coarse grid */
758d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
759d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
760d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
76123ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
762d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
763d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
764d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
765d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
76671959b99SBarry Smith       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
76771959b99SBarry Smith       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
768d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
769d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
7709dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
7712fb0b348SMark F. Adams       ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
7725b42dca8SJed Brown       /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
7735b42dca8SJed Brown        * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
7745b42dca8SJed Brown        * view every subdomain as though they were different. */
7755b42dca8SJed Brown       ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
776d5280255SMark F. Adams     }
777d5280255SMark F. Adams 
778d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
779d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
780e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
781d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
7823b4367a7SBarry Smith     if (mg->galerkin != 2) SETERRQ(comm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
783d5280255SMark F. Adams 
784d5280255SMark F. Adams     /* create cheby smoothers */
785d5280255SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
786d5280255SMark F. Adams          lidx <= fine_level;
787d5280255SMark F. Adams          lidx++, level--) {
788d5280255SMark F. Adams       KSP       smoother;
789890ffe84SMark Adams       PetscBool flag,flag2;
790d5280255SMark F. Adams       PC        subpc;
791d5280255SMark F. Adams 
792ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
793a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
794a2f3521dSMark F. Adams 
795ffc955d6SMark F. Adams       /* do my own cheby */
7966c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
797ffc955d6SMark F. Adams       if (flag) {
798ffc955d6SMark F. Adams         PetscReal emax, emin;
799251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
800890ffe84SMark Adams         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCSOR, &flag2);CHKERRQ(ierr);
801890ffe84SMark 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) */
802890ffe84SMark Adams         else { /* eigen estimate 'emax' -- this is done in cheby */
803e696ed0bSMark F. Adams           KSP eksp;
804e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
805b2a4f308SMark F. Adams           Vec bb, xx;
806038e3b61SMark F. Adams 
8072a7a6963SBarry Smith           ierr = MatCreateVecs(Lmat, &bb, 0);CHKERRQ(ierr);
8082a7a6963SBarry Smith           ierr = MatCreateVecs(Lmat, &xx, 0);CHKERRQ(ierr);
809fc4362bfSMark F. Adams           {
810fc4362bfSMark F. Adams             PetscRandom rctx;
8113b4367a7SBarry Smith             ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
812fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
813fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
814fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
8155b89ad90SMark F. Adams           }
816a2f3521dSMark F. Adams 
817e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
818e696ed0bSMark F. Adams           {
819e696ed0bSMark F. Adams             PetscInt    Istart,Iend,ncols,jj,Ii;
820e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
821e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
822e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
823e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
824e696ed0bSMark F. Adams               if (ncols <= 1) {
825e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
826a94c3b12SMark F. Adams               }
827e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
828a94c3b12SMark F. Adams             }
829a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
830a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
831a94c3b12SMark F. Adams           }
832e696ed0bSMark F. Adams 
8333b4367a7SBarry Smith           ierr = KSPCreate(comm, &eksp);CHKERRQ(ierr);
834806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
835fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
8361a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
8371a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
838f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
839f6536408SMark F. Adams 
840f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
84123ee1639SBarry Smith           ierr = KSPSetOperators(eksp, Lmat, Lmat);CHKERRQ(ierr);
842fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
8435a9b9e01SMark F. Adams 
844d3d0db20SJed Brown           /* set PC type to be same as smoother */
845ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
846b2a4f308SMark F. Adams 
8475a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
8485a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
8495a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
850fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
8515a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
8525a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
8535a9b9e01SMark F. Adams 
854fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
8555a9b9e01SMark F. Adams 
856fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
857fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
858fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
859f6536408SMark F. Adams 
860ffc955d6SMark F. Adams           if (pc_gamg->verbose > 0) {
861a94c3b12SMark F. Adams             PetscInt N1, tt;
862a94c3b12SMark F. Adams             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
8633b4367a7SBarry 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);
864f6536408SMark F. Adams           }
865fc4362bfSMark F. Adams         }
866038e3b61SMark F. Adams         {
867c5bfad50SMark F. Adams           PetscInt N1, N0;
8680298fd71SBarry Smith           ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr);
8690298fd71SBarry Smith           ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr);
870f6536408SMark F. Adams           /* heuristic - is this crap? */
871b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
8725e7c91beSJed Brown           emin  = emax * pc_gamg->eigtarget[0];
8735e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
874038e3b61SMark F. Adams         }
8756c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
876ffc955d6SMark F. Adams       } /* setup checby flag */
877ffc955d6SMark F. Adams     } /* non-coarse levels */
878737a81a9SMark F. Adams 
879d5280255SMark F. Adams     /* clean up */
880d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
881587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
882587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
8835b89ad90SMark F. Adams     }
8840cbbd2e1SMark F. Adams 
8850cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
886806fa848SBarry Smith   } else {
8875f8cf99dSMark F. Adams     KSP smoother;
8883b4367a7SBarry 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__);
8899d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
89023ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
8915f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
8929d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8935f8cf99dSMark F. Adams   }
8945b89ad90SMark F. Adams   PetscFunctionReturn(0);
8955b89ad90SMark F. Adams }
8965b89ad90SMark F. Adams 
897eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8985b89ad90SMark F. Adams /*
8995b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
9005b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
9015b89ad90SMark F. Adams 
9025b89ad90SMark F. Adams    Input Parameter:
9035b89ad90SMark F. Adams .  pc - the preconditioner context
9045b89ad90SMark F. Adams 
9055b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
9065b89ad90SMark F. Adams */
9075b89ad90SMark F. Adams #undef __FUNCT__
9085b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
9095b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
9105b89ad90SMark F. Adams {
9115b89ad90SMark F. Adams   PetscErrorCode ierr;
9125b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
9135b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
9145b89ad90SMark F. Adams 
9155b89ad90SMark F. Adams   PetscFunctionBegin;
9165b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
9179b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
9189b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
9199b8ffb57SJed Brown   }
9201ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
9211ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
9225b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
9235b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
9245b89ad90SMark F. Adams   PetscFunctionReturn(0);
9255b89ad90SMark F. Adams }
9265b89ad90SMark F. Adams 
927676e1743SMark F. Adams 
928676e1743SMark F. Adams #undef __FUNCT__
929676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
930676e1743SMark F. Adams /*@
9311cc46a46SBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via processor reduction.
932676e1743SMark F. Adams 
9331cc46a46SBarry Smith    Logically Collective on PC
934676e1743SMark F. Adams 
935676e1743SMark F. Adams    Input Parameters:
9361cc46a46SBarry Smith +  pc - the preconditioner context
9371cc46a46SBarry Smith -  n - the number of equations
938676e1743SMark F. Adams 
939676e1743SMark F. Adams 
940676e1743SMark F. Adams    Options Database Key:
9411cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
942676e1743SMark F. Adams 
943676e1743SMark F. Adams    Level: intermediate
944676e1743SMark F. Adams 
945676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
946676e1743SMark F. Adams 
947676e1743SMark F. Adams .seealso: ()
948676e1743SMark F. Adams @*/
949676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
950676e1743SMark F. Adams {
951676e1743SMark F. Adams   PetscErrorCode ierr;
952676e1743SMark F. Adams 
953676e1743SMark F. Adams   PetscFunctionBegin;
954676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
955676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
956676e1743SMark F. Adams   PetscFunctionReturn(0);
957676e1743SMark F. Adams }
958676e1743SMark F. Adams 
959676e1743SMark F. Adams #undef __FUNCT__
960676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
9611e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
962676e1743SMark F. Adams {
963c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
964c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
965676e1743SMark F. Adams 
966676e1743SMark F. Adams   PetscFunctionBegin;
9679d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
968676e1743SMark F. Adams   PetscFunctionReturn(0);
969676e1743SMark F. Adams }
970676e1743SMark F. Adams 
971676e1743SMark F. Adams #undef __FUNCT__
972389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
973389730f3SMark F. Adams /*@
974389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
975389730f3SMark F. Adams 
976389730f3SMark F. Adams  Collective on PC
977389730f3SMark F. Adams 
978389730f3SMark F. Adams    Input Parameters:
9791cc46a46SBarry Smith +  pc - the preconditioner context
9801cc46a46SBarry Smith -  n - maximum number of equations to aim for
981389730f3SMark F. Adams 
982389730f3SMark F. Adams    Options Database Key:
9831cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
984389730f3SMark F. Adams 
985389730f3SMark F. Adams    Level: intermediate
986389730f3SMark F. Adams 
987389730f3SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
988389730f3SMark F. Adams 
989389730f3SMark F. Adams @*/
990389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
991389730f3SMark F. Adams {
992389730f3SMark F. Adams   PetscErrorCode ierr;
993389730f3SMark F. Adams 
994389730f3SMark F. Adams   PetscFunctionBegin;
995389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
996389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
997389730f3SMark F. Adams   PetscFunctionReturn(0);
998389730f3SMark F. Adams }
999389730f3SMark F. Adams 
1000389730f3SMark F. Adams #undef __FUNCT__
1001389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
10021e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1003389730f3SMark F. Adams {
1004389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1005389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1006389730f3SMark F. Adams 
1007389730f3SMark F. Adams   PetscFunctionBegin;
10089d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
1009389730f3SMark F. Adams   PetscFunctionReturn(0);
1010389730f3SMark F. Adams }
1011389730f3SMark F. Adams 
1012389730f3SMark F. Adams #undef __FUNCT__
10138263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
1014676e1743SMark F. Adams /*@
10158263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
1016676e1743SMark F. Adams 
1017676e1743SMark F. Adams    Collective on PC
1018676e1743SMark F. Adams 
1019676e1743SMark F. Adams    Input Parameters:
10201cc46a46SBarry Smith +  pc - the preconditioner context
10211cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1022676e1743SMark F. Adams 
1023676e1743SMark F. Adams    Options Database Key:
10241cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
1025676e1743SMark F. Adams 
1026676e1743SMark F. Adams    Level: intermediate
1027676e1743SMark F. Adams 
1028676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1029676e1743SMark F. Adams 
1030676e1743SMark F. Adams .seealso: ()
1031676e1743SMark F. Adams @*/
10328263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1033676e1743SMark F. Adams {
1034676e1743SMark F. Adams   PetscErrorCode ierr;
1035676e1743SMark F. Adams 
1036676e1743SMark F. Adams   PetscFunctionBegin;
1037676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10388263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1039676e1743SMark F. Adams   PetscFunctionReturn(0);
1040676e1743SMark F. Adams }
1041676e1743SMark F. Adams 
1042676e1743SMark F. Adams #undef __FUNCT__
10438263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
10441e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1045676e1743SMark F. Adams {
1046c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1047c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1048676e1743SMark F. Adams 
1049676e1743SMark F. Adams   PetscFunctionBegin;
10509d5b6da9SMark F. Adams   pc_gamg->repart = n;
1051676e1743SMark F. Adams   PetscFunctionReturn(0);
1052676e1743SMark F. Adams }
1053676e1743SMark F. Adams 
1054676e1743SMark F. Adams #undef __FUNCT__
10551cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation"
1056dfd5c07aSMark F. Adams /*@
10571cc46a46SBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner
1058dfd5c07aSMark F. Adams 
1059dfd5c07aSMark F. Adams    Collective on PC
1060dfd5c07aSMark F. Adams 
1061dfd5c07aSMark F. Adams    Input Parameters:
10621cc46a46SBarry Smith +  pc - the preconditioner context
10631cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1064dfd5c07aSMark F. Adams 
1065dfd5c07aSMark F. Adams    Options Database Key:
10661cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
1067dfd5c07aSMark F. Adams 
1068dfd5c07aSMark F. Adams    Level: intermediate
1069dfd5c07aSMark F. Adams 
1070dfd5c07aSMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1071dfd5c07aSMark F. Adams 
1072dfd5c07aSMark F. Adams .seealso: ()
1073dfd5c07aSMark F. Adams @*/
10741cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1075dfd5c07aSMark F. Adams {
1076dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1077dfd5c07aSMark F. Adams 
1078dfd5c07aSMark F. Adams   PetscFunctionBegin;
1079dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10801cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1081dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1082dfd5c07aSMark F. Adams }
1083dfd5c07aSMark F. Adams 
1084dfd5c07aSMark F. Adams #undef __FUNCT__
10851cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG"
10861cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1087dfd5c07aSMark F. Adams {
1088dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1089dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1090dfd5c07aSMark F. Adams 
1091dfd5c07aSMark F. Adams   PetscFunctionBegin;
1092dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1093dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1094dfd5c07aSMark F. Adams }
1095dfd5c07aSMark F. Adams 
1096dfd5c07aSMark F. Adams #undef __FUNCT__
1097ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1098ffc955d6SMark F. Adams /*@
1099ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1100ffc955d6SMark F. Adams 
1101ffc955d6SMark F. Adams    Collective on PC
1102ffc955d6SMark F. Adams 
1103ffc955d6SMark F. Adams    Input Parameters:
1104ffc955d6SMark F. Adams .  pc - the preconditioner context
1105ffc955d6SMark F. Adams 
1106ffc955d6SMark F. Adams 
1107ffc955d6SMark F. Adams    Options Database Key:
1108ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1109ffc955d6SMark F. Adams 
1110ffc955d6SMark F. Adams    Level: intermediate
1111ffc955d6SMark F. Adams 
1112ffc955d6SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1113ffc955d6SMark F. Adams 
1114ffc955d6SMark F. Adams .seealso: ()
1115ffc955d6SMark F. Adams @*/
1116ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1117ffc955d6SMark F. Adams {
1118ffc955d6SMark F. Adams   PetscErrorCode ierr;
1119ffc955d6SMark F. Adams 
1120ffc955d6SMark F. Adams   PetscFunctionBegin;
1121ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1122ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1123ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1124ffc955d6SMark F. Adams }
1125ffc955d6SMark F. Adams 
1126ffc955d6SMark F. Adams #undef __FUNCT__
1127ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
11281e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1129ffc955d6SMark F. Adams {
1130ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1131ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1132ffc955d6SMark F. Adams 
1133ffc955d6SMark F. Adams   PetscFunctionBegin;
1134ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1135ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1136ffc955d6SMark F. Adams }
1137ffc955d6SMark F. Adams 
1138ffc955d6SMark F. Adams #undef __FUNCT__
11394ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
11404ef23d27SMark F. Adams /*@
11411cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
11424ef23d27SMark F. Adams 
11434ef23d27SMark F. Adams    Not collective on PC
11444ef23d27SMark F. Adams 
11454ef23d27SMark F. Adams    Input Parameters:
11461cc46a46SBarry Smith +  pc - the preconditioner
11471cc46a46SBarry Smith -  n - the maximum number of levels to use
11484ef23d27SMark F. Adams 
11494ef23d27SMark F. Adams    Options Database Key:
11504ef23d27SMark F. Adams .  -pc_mg_levels
11514ef23d27SMark F. Adams 
11524ef23d27SMark F. Adams    Level: intermediate
11534ef23d27SMark F. Adams 
11544ef23d27SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
11554ef23d27SMark F. Adams 
11564ef23d27SMark F. Adams .seealso: ()
11574ef23d27SMark F. Adams @*/
11584ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
11594ef23d27SMark F. Adams {
11604ef23d27SMark F. Adams   PetscErrorCode ierr;
11614ef23d27SMark F. Adams 
11624ef23d27SMark F. Adams   PetscFunctionBegin;
11634ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11644ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
11654ef23d27SMark F. Adams   PetscFunctionReturn(0);
11664ef23d27SMark F. Adams }
11674ef23d27SMark F. Adams 
11684ef23d27SMark F. Adams #undef __FUNCT__
11694ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
11701e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
11714ef23d27SMark F. Adams {
11724ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
11734ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
11744ef23d27SMark F. Adams 
11754ef23d27SMark F. Adams   PetscFunctionBegin;
11769d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
11774ef23d27SMark F. Adams   PetscFunctionReturn(0);
11784ef23d27SMark F. Adams }
11794ef23d27SMark F. Adams 
11804ef23d27SMark F. Adams #undef __FUNCT__
11813542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
11823542efc5SMark F. Adams /*@
11833542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
11843542efc5SMark F. Adams 
11853542efc5SMark F. Adams    Not collective on PC
11863542efc5SMark F. Adams 
11873542efc5SMark F. Adams    Input Parameters:
11881cc46a46SBarry Smith +  pc - the preconditioner context
11891cc46a46SBarry Smith -  threshold - the threshold value, 0.0 is the lowest possible
11903542efc5SMark F. Adams 
11913542efc5SMark F. Adams    Options Database Key:
11921cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
11933542efc5SMark F. Adams 
11943542efc5SMark F. Adams    Level: intermediate
11953542efc5SMark F. Adams 
11963542efc5SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
11973542efc5SMark F. Adams 
11983542efc5SMark F. Adams .seealso: ()
11993542efc5SMark F. Adams @*/
12003542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
12013542efc5SMark F. Adams {
12023542efc5SMark F. Adams   PetscErrorCode ierr;
12033542efc5SMark F. Adams 
12043542efc5SMark F. Adams   PetscFunctionBegin;
12053542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12063542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
12073542efc5SMark F. Adams   PetscFunctionReturn(0);
12083542efc5SMark F. Adams }
12093542efc5SMark F. Adams 
12103542efc5SMark F. Adams #undef __FUNCT__
12113542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
12121e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
12133542efc5SMark F. Adams {
1214c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1215c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12163542efc5SMark F. Adams 
12173542efc5SMark F. Adams   PetscFunctionBegin;
12189d5b6da9SMark F. Adams   pc_gamg->threshold = n;
12193542efc5SMark F. Adams   PetscFunctionReturn(0);
12203542efc5SMark F. Adams }
12213542efc5SMark F. Adams 
12223542efc5SMark F. Adams #undef __FUNCT__
12239d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1224676e1743SMark F. Adams /*@
1225c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1226676e1743SMark F. Adams 
1227676e1743SMark F. Adams    Collective on PC
1228676e1743SMark F. Adams 
1229676e1743SMark F. Adams    Input Parameters:
1230c60c7ad4SBarry Smith +  pc - the preconditioner context
1231c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1232676e1743SMark F. Adams 
1233676e1743SMark F. Adams    Options Database Key:
1234c60c7ad4SBarry Smith .  -pc_gamg_type <agg,geo,classical>
1235676e1743SMark F. Adams 
1236676e1743SMark F. Adams    Level: intermediate
1237676e1743SMark F. Adams 
1238676e1743SMark F. Adams    Concepts: Unstructured multrigrid preconditioner
1239676e1743SMark F. Adams 
1240c60c7ad4SBarry Smith .seealso: PCGAMGGetType(), PCGAMG
1241676e1743SMark F. Adams @*/
124219fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1243676e1743SMark F. Adams {
1244676e1743SMark F. Adams   PetscErrorCode ierr;
1245676e1743SMark F. Adams 
1246676e1743SMark F. Adams   PetscFunctionBegin;
1247676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1248806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1249676e1743SMark F. Adams   PetscFunctionReturn(0);
1250676e1743SMark F. Adams }
1251676e1743SMark F. Adams 
1252676e1743SMark F. Adams #undef __FUNCT__
1253c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType"
1254c60c7ad4SBarry Smith /*@
1255c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1256c60c7ad4SBarry Smith 
1257c60c7ad4SBarry Smith    Collective on PC
1258c60c7ad4SBarry Smith 
1259c60c7ad4SBarry Smith    Input Parameter:
1260c60c7ad4SBarry Smith .  pc - the preconditioner context
1261c60c7ad4SBarry Smith 
1262c60c7ad4SBarry Smith    Output Parameter:
1263c60c7ad4SBarry Smith .  type - the type of algorithm used
1264c60c7ad4SBarry Smith 
1265c60c7ad4SBarry Smith    Level: intermediate
1266c60c7ad4SBarry Smith 
1267c60c7ad4SBarry Smith    Concepts: Unstructured multrigrid preconditioner
1268c60c7ad4SBarry Smith 
1269c60c7ad4SBarry Smith .seealso: PCGAMGSetType()
1270c60c7ad4SBarry Smith @*/
1271c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1272c60c7ad4SBarry Smith {
1273c60c7ad4SBarry Smith   PetscErrorCode ierr;
1274c60c7ad4SBarry Smith 
1275c60c7ad4SBarry Smith   PetscFunctionBegin;
1276c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1277c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1278c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1279c60c7ad4SBarry Smith }
1280c60c7ad4SBarry Smith 
1281c60c7ad4SBarry Smith #undef __FUNCT__
1282c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG"
1283c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1284c60c7ad4SBarry Smith {
1285c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1286c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1287c60c7ad4SBarry Smith 
1288c60c7ad4SBarry Smith   PetscFunctionBegin;
1289c60c7ad4SBarry Smith   *type = pc_gamg->type;
1290c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1291c60c7ad4SBarry Smith }
1292c60c7ad4SBarry Smith 
1293c60c7ad4SBarry Smith #undef __FUNCT__
12949d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
12951e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1296676e1743SMark F. Adams {
12979d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
12981ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
12991ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1300676e1743SMark F. Adams 
1301676e1743SMark F. Adams   PetscFunctionBegin;
1302c60c7ad4SBarry Smith   pc_gamg->type = type;
13031c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
13049d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
13051ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
13063ae0bb68SMark Adams     /* there was something here - kill it */
13071ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
13081ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1309e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
13103ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
13113ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
13123ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
13133ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
13143ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
13153ae0bb68SMark Adams     if (pc_gamg->data_sz) {
13163ae0bb68SMark Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
13173ae0bb68SMark Adams       pc_gamg->data_sz = 0;
1318c60c7ad4SBarry Smith     } else if (pc_gamg->data) {
13193ae0bb68SMark Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); /* can this happen ? */
13203ae0bb68SMark Adams     }
13211ab5ffc9SJed Brown   }
13221ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
13231ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
13249d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1325676e1743SMark F. Adams   PetscFunctionReturn(0);
1326676e1743SMark F. Adams }
1327676e1743SMark F. Adams 
13285b89ad90SMark F. Adams #undef __FUNCT__
13295b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
13308c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptions *PetscOptionsObject,PC pc)
13315b89ad90SMark F. Adams {
1332676e1743SMark F. Adams   PetscErrorCode ierr;
1333676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1334676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1335676e1743SMark F. Adams   PetscBool      flag;
13365e7c91beSJed Brown   PetscInt       two   = 2;
13373b4367a7SBarry Smith   MPI_Comm       comm;
13385b89ad90SMark F. Adams 
13395b89ad90SMark F. Adams   PetscFunctionBegin;
13403b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1341e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1342676e1743SMark F. Adams   {
1343b7cbab4eSMark Adams     /* -pc_gamg_type */
1344b7cbab4eSMark Adams     {
1345bd94a7aaSJed Brown       char tname[256];
13461a1c1e04SBarry Smith       ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1347bd94a7aaSJed Brown       if (flag) {
1348bd94a7aaSJed Brown         ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
13491ab5ffc9SJed Brown       }
1350b7cbab4eSMark Adams     }
135175b74bdaSMark F. Adams     /* -pc_gamg_verbose */
135294ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG","none", pc_gamg->verbose,&pc_gamg->verbose, NULL);CHKERRQ(ierr);
13538263b398SMark F. Adams     /* -pc_gamg_repartition */
135494ae4db5SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
1355dfd5c07aSMark F. Adams     /* -pc_gamg_reuse_interpolation */
13561cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1357ffc955d6SMark F. Adams     /* -pc_gamg_use_agg_gasm */
135894ae4db5SBarry 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);
1359c20e4228SMark F. Adams     /* -pc_gamg_process_eq_limit */
136094ae4db5SBarry 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);
1361389730f3SMark F. Adams     /* -pc_gamg_coarse_eq_limit */
136294ae4db5SBarry 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);
1363c20e4228SMark F. Adams     /* -pc_gamg_threshold */
136494ae4db5SBarry 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);
1365806fa848SBarry Smith     if (flag && pc_gamg->verbose) {
13663b4367a7SBarry Smith       ierr = PetscPrintf(comm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1367806fa848SBarry Smith     }
1368b7cbab4eSMark Adams     /* -pc_gamg_eigtarget */
13690298fd71SBarry Smith     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr);
137094ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1371b7cbab4eSMark Adams 
1372b7cbab4eSMark Adams     /* set options for subtype */
1373e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1374676e1743SMark F. Adams   }
1375676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
13765b89ad90SMark F. Adams   PetscFunctionReturn(0);
13775b89ad90SMark F. Adams }
13785b89ad90SMark F. Adams 
13795b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
13805b89ad90SMark F. Adams /*MC
13811cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
13825b89ad90SMark F. Adams 
1383280d9858SJed Brown    Options Database Keys:
13845b89ad90SMark F. Adams    Multigrid options(inherited)
13851cc46a46SBarry Smith +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1386280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1387280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
13888c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
13895b89ad90SMark F. Adams 
13901cc46a46SBarry Smith 
13911cc46a46SBarry Smith   Notes: In order to obtain good performance for PCGAMG for vector valued problems you must
13921cc46a46SBarry Smith $       Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
13931cc46a46SBarry Smith $       Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
13941cc46a46SBarry Smith $       See the Users Manual Chapter 4 for more details
13951cc46a46SBarry Smith 
13965b89ad90SMark F. Adams   Level: intermediate
1397280d9858SJed Brown 
13981cc46a46SBarry Smith   Concepts: algebraic multigrid
13995b89ad90SMark F. Adams 
14001cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
14011cc46a46SBarry Smith            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartitioning(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseASMAggs(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType()
14025b89ad90SMark F. Adams M*/
1403b2573a8aSBarry Smith 
14045b89ad90SMark F. Adams #undef __FUNCT__
14055b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
14068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
14075b89ad90SMark F. Adams {
14085b89ad90SMark F. Adams   PetscErrorCode ierr;
14095b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
14105b89ad90SMark F. Adams   PC_MG          *mg;
14110cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14129d5b6da9SMark F. Adams   static long count = 0;
14135ee9c036SSatish Balay #endif
14145b89ad90SMark F. Adams 
14155b89ad90SMark F. Adams   PetscFunctionBegin;
14165b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14175b89ad90SMark F. Adams   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
14185b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14195b89ad90SMark F. Adams 
14205b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1421b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
14225b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
1423ce4cda84SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
14245b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
14255b89ad90SMark F. Adams 
1426b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
14271ab5ffc9SJed Brown 
14289d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14299d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14309d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14319d5b6da9SMark F. Adams   pc_gamg->data    = 0;
14325b89ad90SMark F. Adams 
14339d5b6da9SMark F. Adams   /* register AMG type */
14343e3471ccSMark Adams   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
14359d5b6da9SMark F. Adams 
14369d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
14375b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
14385b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
14395b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
14405b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
14415b89ad90SMark F. Adams 
1442bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1443bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1444bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
14451cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1446bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1447bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1448bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1449c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1450bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
14519d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1452d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
1453ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1454038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
145525a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1456d3042614SMark Adams   pc_gamg->threshold        = 0.;
14579d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
14589d5b6da9SMark F. Adams   pc_gamg->verbose          = 0;
14599d5b6da9SMark F. Adams   pc_gamg->emax_id          = -1;
146057d29afaSToby Isaac   pc_gamg->firstCoarsen     = PETSC_FALSE;
14615e7c91beSJed Brown   pc_gamg->eigtarget[0]     = 0.05;
14625e7c91beSJed Brown   pc_gamg->eigtarget[1]     = 1.05;
1463c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
14649d5b6da9SMark F. Adams 
14650cbbd2e1SMark F. Adams   /* private events */
14660cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1467785cba28SMark F. Adams   if (count++ == 0) {
1468806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1469806fa848SBarry Smith     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
14700cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
14710cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
14720cbbd2e1SMark F. Adams     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1473806fa848SBarry Smith     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1474806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1475806fa848SBarry Smith     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1476806fa848SBarry Smith     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1477806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1478806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1479806fa848SBarry Smith     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1480806fa848SBarry Smith     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1481806fa848SBarry Smith     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1482806fa848SBarry Smith     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1483806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1484806fa848SBarry Smith     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1485f852f58cSMark F. Adams 
14860cbbd2e1SMark F. Adams     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
14870cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
14880cbbd2e1SMark F. Adams     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1489b4fbaa2aSMark F. Adams     /* create timer stages */
1490b4fbaa2aSMark F. Adams #if defined GAMG_STAGES
1491b4fbaa2aSMark F. Adams     {
1492b4fbaa2aSMark F. Adams       char     str[32];
1493b4fbaa2aSMark F. Adams       PetscInt lidx;
1494806fa848SBarry Smith       sprintf(str,"MG Level %d (finest)",0);
1495806fa848SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1496b4fbaa2aSMark F. Adams       for (lidx=1; lidx<9; lidx++) {
1497b4fbaa2aSMark F. Adams         sprintf(str,"MG Level %d",lidx);
1498806fa848SBarry Smith         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1499b4fbaa2aSMark F. Adams       }
1500b4fbaa2aSMark F. Adams     }
1501b4fbaa2aSMark F. Adams #endif
1502b4fbaa2aSMark F. Adams   }
1503b4fbaa2aSMark F. Adams #endif
1504bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1505bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
15065b89ad90SMark F. Adams   PetscFunctionReturn(0);
15075b89ad90SMark F. Adams }
15083e3471ccSMark Adams 
15093e3471ccSMark Adams #undef __FUNCT__
15103e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage"
15113e3471ccSMark Adams /*@C
15123e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
15133e3471ccSMark Adams  from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
15143e3471ccSMark Adams  when using static libraries.
15153e3471ccSMark Adams 
15163e3471ccSMark Adams  Level: developer
15173e3471ccSMark Adams 
15183e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
15193e3471ccSMark Adams  .seealso: PetscInitialize()
15203e3471ccSMark Adams @*/
15213e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
15223e3471ccSMark Adams {
15233e3471ccSMark Adams   PetscErrorCode ierr;
15243e3471ccSMark Adams 
15253e3471ccSMark Adams   PetscFunctionBegin;
15263e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
15273e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
15283e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
15293e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
15308e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
15313e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1532c1c463dbSMark Adams 
1533c1c463dbSMark Adams   /* general events */
1534c1c463dbSMark Adams #if defined PETSC_USE_LOG
1535c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1536c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1537c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1538c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1539c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1540c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1541c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1542c1c463dbSMark Adams #endif
1543c1c463dbSMark Adams 
15443e3471ccSMark Adams   PetscFunctionReturn(0);
15453e3471ccSMark Adams }
15463e3471ccSMark Adams 
15473e3471ccSMark Adams #undef __FUNCT__
15483e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage"
15493e3471ccSMark Adams /*@C
15503e3471ccSMark Adams  PCGAMGFinalizePackage - This function destroys everything in the PCGAMG package. It is
15513e3471ccSMark Adams  called from PetscFinalize().
15523e3471ccSMark Adams 
15533e3471ccSMark Adams  Level: developer
15543e3471ccSMark Adams 
15553e3471ccSMark Adams  .keywords: Petsc, destroy, package
15563e3471ccSMark Adams  .seealso: PetscFinalize()
15573e3471ccSMark Adams @*/
15583e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
15593e3471ccSMark Adams {
15603e3471ccSMark Adams   PetscErrorCode ierr;
15613e3471ccSMark Adams 
15623e3471ccSMark Adams   PetscFunctionBegin;
15633e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
15643e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
15653e3471ccSMark Adams   PetscFunctionReturn(0);
15663e3471ccSMark Adams }
1567a36cf38bSToby Isaac 
1568a36cf38bSToby Isaac #undef __FUNCT__
1569a36cf38bSToby Isaac #define __FUNCT__ "PCGAMGRegister"
1570a36cf38bSToby Isaac /*@C
1571a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1572a36cf38bSToby Isaac 
1573a36cf38bSToby Isaac  Input Parameters:
1574a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1575a36cf38bSToby Isaac  - create - function for creating the gamg context.
1576a36cf38bSToby Isaac 
1577a36cf38bSToby Isaac   Level: advanced
1578a36cf38bSToby Isaac 
1579a36cf38bSToby Isaac  .seealso: PCGAMGGetContext(), PCGAMGSetContext()
1580a36cf38bSToby Isaac @*/
1581a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1582a36cf38bSToby Isaac {
1583a36cf38bSToby Isaac   PetscErrorCode ierr;
1584a36cf38bSToby Isaac 
1585a36cf38bSToby Isaac   PetscFunctionBegin;
1586a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1587a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1588a36cf38bSToby Isaac   PetscFunctionReturn(0);
1589a36cf38bSToby Isaac }
1590a36cf38bSToby Isaac 
1591