xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision c1eae691a6266ae3a8f3c32c4f419f5f45c3ffa8)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4af0996ceSBarry Smith #include <petsc/private/matimpl.h>
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
65b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */
7f96513f1SMatthew G Knepley 
80cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
90cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
10b4fbaa2aSMark F. Adams #endif
110cbbd2e1SMark F. Adams 
120cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
13fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG;
14fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO;
150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
19fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG;
200cbbd2e1SMark F. Adams #endif
210cbbd2e1SMark F. Adams 
22b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
230cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
24*c1eae691SMark Adams static PetscLogStage gamg_stages[PETSC_GAMG_MAXLEVELS];
25b4fbaa2aSMark F. Adams #endif
26f96513f1SMatthew G Knepley 
27140e18c1SBarry Smith static PetscFunctionList GAMGList = 0;
283e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
299d5b6da9SMark F. Adams 
30d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
31d3d6bff4SMark F. Adams #undef __FUNCT__
32d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
33d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
34d3d6bff4SMark F. Adams {
35d3d6bff4SMark F. Adams   PetscErrorCode ierr;
36d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
37d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
38d3d6bff4SMark F. Adams 
39d3d6bff4SMark F. Adams   PetscFunctionBegin;
4062294041SBarry Smith   if (pc_gamg->data) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen, cleaned up in SetUp\n");
411c1aac46SBarry Smith   pc_gamg->data_sz = 0;
42878e152fSMark F. Adams   ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
43a2f3521dSMark F. Adams   PetscFunctionReturn(0);
44a2f3521dSMark F. Adams }
45a2f3521dSMark F. Adams 
465b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
475b89ad90SMark F. Adams /*
48c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
49a147abb0SMark F. Adams      of active processors.
505b89ad90SMark F. Adams 
515b89ad90SMark F. Adams    Input Parameter:
52a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
53a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
549d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
55c5bfad50SMark F. Adams    . cr_bs - coarse block size
563530afc2SMark F. Adams    In/Output Parameter:
57a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
58afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
5911e60469SMark F. Adams    Output Parameter:
603530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
615b89ad90SMark F. Adams */
625cb416c2SMark F. Adams 
635b89ad90SMark F. Adams #undef __FUNCT__
64c238b0ebSToby Isaac #define __FUNCT__ "PCGAMGCreateLevel_GAMG"
65171cca9aSMark Adams static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,IS * Pcolumnperm, PetscBool is_last)
665b89ad90SMark F. Adams {
67a2f3521dSMark F. Adams   PetscErrorCode  ierr;
689d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
69486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
70a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
713b4367a7SBarry Smith   MPI_Comm        comm;
72c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
733ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
745b89ad90SMark F. Adams 
755b89ad90SMark F. Adams   PetscFunctionBegin;
763b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
773b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
783b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
79c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
809d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
81038e3b61SMark F. Adams 
823ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
830298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
843ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
853ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
8673911c69SBarry Smith   } else {
873ae0bb68SMark Adams     PetscInt  bs;
883ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
893ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
903ae0bb68SMark Adams   }
91a2f3521dSMark F. Adams 
92c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
93171cca9aSMark Adams   if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) new_size = 1;
94171cca9aSMark Adams   else {
95472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
960298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
97a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
987f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
99c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
100a2f3521dSMark F. Adams   }
101f852f58cSMark F. Adams 
1023cb8563fSToby Isaac   if (Pcolumnperm) *Pcolumnperm = NULL;
1033cb8563fSToby Isaac 
104a90e85d9SMark Adams   if (!pc_gamg->repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1052fa5cd67SKarl Rupp   else {
1063ae0bb68SMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
107885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
108e33ef3b1SMark F. Adams 
10971959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
11071959b99SBarry 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);
1110cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1120cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
113b4fbaa2aSMark F. Adams #endif
114a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
115785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
116a90e85d9SMark Adams     if (pc_gamg->repart) {
117a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1185a9b9e01SMark F. Adams       Mat adj;
1195a9b9e01SMark F. Adams 
120cf8ae1d3SMark Adams      ierr = PetscInfo3(pc,"Repartition: size (active): %D --> %D, %D local equations\n",*a_nactive_proc,new_size,ncrs_eq);CHKERRQ(ierr);
1215a9b9e01SMark F. Adams 
122a2f3521dSMark F. Adams       /* get 'adj' */
123c5bfad50SMark F. Adams       if (cr_bs == 1) {
124038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
125806fa848SBarry Smith       } else {
126a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
127eb07cef2SMark F. Adams         Mat               tMat;
128a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
129b4fbaa2aSMark F. Adams         const PetscScalar *vals;
130b4fbaa2aSMark F. Adams         const PetscInt    *idx;
131a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1329057884aSMark F. Adams         static PetscInt   llev = 0;
133d9558ea9SBarry Smith         MatType           mtype;
134b4fbaa2aSMark F. Adams 
135e632b94dSBarry Smith         ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr);
136a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
137a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
138c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
13958471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
140c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
141c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
14258471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1433ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1443ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
14558471d46SMark F. Adams         }
1466876a03eSMark F. Adams 
147d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
1483b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1493ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
150d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
151a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
152a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
153e632b94dSBarry Smith         ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
154eb07cef2SMark F. Adams 
155a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
156c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
15722063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
158eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
159c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
160eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
161eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
162eb07cef2SMark F. Adams           }
16322063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
164eb07cef2SMark F. Adams         }
165eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167eb07cef2SMark F. Adams 
168b4fbaa2aSMark F. Adams         if (llev++ == -1) {
169b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
1708caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
1713b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
172b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
1733bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
174b4fbaa2aSMark F. Adams         }
175eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
176eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
177a2f3521dSMark F. Adams       } /* create 'adj' */
178f150b916SMark F. Adams 
179a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
1805a9b9e01SMark F. Adams         char            prefix[256];
1815a9b9e01SMark F. Adams         const char      *pcpre;
182b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
183b4fbaa2aSMark F. Adams         MatPartitioning mpart;
184a4b7d37bSMark F. Adams         IS              proc_is;
185a2f3521dSMark F. Adams         PetscInt        targetPE;
1862f03bc48SMark F. Adams 
1873b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
1885ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
1899d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
1908caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
19159a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
19211e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
193c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
194a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
19511e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1965a9b9e01SMark F. Adams 
1975ef31b24SMark F. Adams         /* collect IS info */
198785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
199a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
200a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
201c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
202a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
203c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
204a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
205eb07cef2SMark F. Adams           }
2065ef31b24SMark F. Adams         }
207a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
208a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2095ef31b24SMark F. Adams       }
2105ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2115a9b9e01SMark F. Adams 
2123b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2138263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
214806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
215a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
21662294041SBarry Smith 
2175a9b9e01SMark F. Adams       /* find factor */
218c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2195a9b9e01SMark F. Adams       else {
2205a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2215a9b9e01SMark F. Adams         jj = -1;
222c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
2232a82295dSMark Adams           if (!(size%kk)) { /* a candidate */
224c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2255a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2265a9b9e01SMark F. Adams             if (fact > best_fact) {
2275a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
2285a9b9e01SMark F. Adams             }
2295a9b9e01SMark F. Adams           }
2305a9b9e01SMark F. Adams         }
2315a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
232a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
2335a9b9e01SMark F. Adams       }
234c5df96a5SBarry Smith       new_size = size/rfactor;
2355a9b9e01SMark F. Adams 
236c5df96a5SBarry Smith       if (new_size==nactive) {
237a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2385a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
239302440fdSBarry Smith         ierr = PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr);
240fe682e87SBarry Smith #if defined PETSC_GAMG_USE_LOG
241fe682e87SBarry Smith         ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
242fe682e87SBarry Smith #endif
2435a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2445a9b9e01SMark F. Adams       }
2455a9b9e01SMark F. Adams 
246302440fdSBarry Smith       ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr);
247c5df96a5SBarry Smith       targetPE = rank/rfactor;
2483b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
249a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
250e33ef3b1SMark F. Adams 
25111e60469SMark F. Adams     /*
252a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
25311e60469SMark F. Adams      */
254a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2557700e67bSMark Adams     is_eq_num_prim = is_eq_num;
25611e60469SMark F. Adams     /*
257a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
25811e60469SMark F. Adams      */
259c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
260c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
261a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
2623ae0bb68SMark Adams     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */
263a2f3521dSMark F. Adams 
264a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2650cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2660cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
267b4fbaa2aSMark F. Adams #endif
268885364a3SMark 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 */
269885364a3SMark Adams     {
270885364a3SMark Adams     Vec            src_crd, dest_crd;
271885364a3SMark 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;
272885364a3SMark Adams     VecScatter     vecscat;
273885364a3SMark Adams     PetscScalar    *array;
274885364a3SMark Adams     IS isscat;
275a2f3521dSMark F. Adams 
276a2f3521dSMark F. Adams     /* move data (for primal equations only) */
27722063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
2783b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
2793ae0bb68SMark Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
280c0dedaeaSBarry Smith     ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
28111e60469SMark F. Adams     /*
2829d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
283c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
28411e60469SMark F. Adams      */
285854ce69bSBarry Smith     ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
286a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2873ae0bb68SMark Adams     for (ii=0,jj=0; ii<ncrs; ii++) {
288c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
289a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
29011e60469SMark F. Adams     }
291a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2923ae0bb68SMark Adams     ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
29392a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
29411e60469SMark F. Adams     /*
29511e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
29611e60469SMark F. Adams      */
2973ae0bb68SMark Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
2989d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
2993ae0bb68SMark Adams       const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3003ae0bb68SMark Adams       for (ii=0; ii<ncrs; ii++) {
3019d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
302a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
303c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
304676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
305d3d6bff4SMark F. Adams         }
306038e3b61SMark F. Adams       }
307eb07cef2SMark F. Adams     }
308eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
309eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
31011e60469SMark F. Adams     /*
31111e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
31211e60469SMark F. Adams       to the correct processor
31311e60469SMark F. Adams     */
3140298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
31511e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
31611e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31711e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31811e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
31911e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
32011e60469SMark F. Adams     /*
32111e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
32211e60469SMark F. Adams     */
323c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
324578f55a3SPeter Brune     ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3252fa5cd67SKarl Rupp 
3263ae0bb68SMark Adams     pc_gamg->data_sz = node_data_sz*ncrs_new;
3273ae0bb68SMark Adams     strideNew        = ncrs_new*ndata_rows;
3282fa5cd67SKarl Rupp 
32911e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3309d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3313ae0bb68SMark Adams       for (ii=0; ii<ncrs_new; ii++) {
3329d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
333a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
334c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
335d3d6bff4SMark F. Adams         }
336038e3b61SMark F. Adams       }
337038e3b61SMark F. Adams     }
33811e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
33911e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
340885364a3SMark Adams     }
341a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3420cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3430cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
344ed3f9983SMark F. Adams #endif
345a2f3521dSMark F. Adams 
34611e60469SMark F. Adams     /*
34711e60469SMark F. Adams       Invert for MatGetSubMatrix
34811e60469SMark F. Adams     */
349a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
350a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
351c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
352a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
353a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
354a2f3521dSMark F. Adams     }
3553cb8563fSToby Isaac     if (Pcolumnperm) {
3563cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3573cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3583cb8563fSToby Isaac     }
359a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3610cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3620cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
363ed3f9983SMark F. Adams #endif
364a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
365a2f3521dSMark F. Adams     {
366a2f3521dSMark F. Adams       Mat mat;
367806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
368a2f3521dSMark F. Adams       *a_Amat_crs = mat;
369a2f3521dSMark F. Adams     }
370038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
371a2f3521dSMark F. Adams 
3720cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3730cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
374ed3f9983SMark F. Adams #endif
37511e60469SMark F. Adams     /* prolongator */
37611e60469SMark F. Adams     {
37711e60469SMark F. Adams       IS       findices;
378a2f3521dSMark F. Adams       PetscInt Istart,Iend;
379a2f3521dSMark F. Adams       Mat      Pnew;
38062294041SBarry Smith 
381a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
3820cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3830cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
384ed3f9983SMark F. Adams #endif
3853b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
386c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
387806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
38811e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
389c5bfad50SMark F. Adams 
3900cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3910cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
392ed3f9983SMark F. Adams #endif
3933530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
394a2f3521dSMark F. Adams 
395a2f3521dSMark F. Adams       /* output - repartitioned */
396a2f3521dSMark F. Adams       *a_P_inout = Pnew;
397e33ef3b1SMark F. Adams     }
398a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
3995b89ad90SMark F. Adams 
400c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
401a2f3521dSMark F. Adams   }
4025b89ad90SMark F. Adams   PetscFunctionReturn(0);
4035b89ad90SMark F. Adams }
4045b89ad90SMark F. Adams 
4055b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4065b89ad90SMark F. Adams /*
4075b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4085b89ad90SMark F. Adams                     by setting data structures and options.
4095b89ad90SMark F. Adams 
4105b89ad90SMark F. Adams    Input Parameter:
4115b89ad90SMark F. Adams .  pc - the preconditioner context
4125b89ad90SMark F. Adams 
4135b89ad90SMark F. Adams */
4145b89ad90SMark F. Adams #undef __FUNCT__
4155b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4169d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4175b89ad90SMark F. Adams {
4185b89ad90SMark F. Adams   PetscErrorCode ierr;
4199d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4205b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4212adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
422*c1eae691SMark Adams   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_GAMG_MAXLEVELS];
4233b4367a7SBarry Smith   MPI_Comm       comm;
424c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
425*c1eae691SMark Adams   Mat            Aarr[PETSC_GAMG_MAXLEVELS],Parr[PETSC_GAMG_MAXLEVELS];
426*c1eae691SMark Adams   IS             *ASMLocalIDsArr[PETSC_GAMG_MAXLEVELS];
427a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
428569f4572SMark Adams   MatInfo        info;
429171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
4305ef31b24SMark F. Adams 
4315b89ad90SMark F. Adams   PetscFunctionBegin;
4323b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4333b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4343b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
435dfd5c07aSMark F. Adams 
43684d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
4371c1aac46SBarry Smith     if ((PetscBool)(!pc_gamg->reuse_prol)) {
438878e152fSMark F. Adams       /* reset everything */
439878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
440878e152fSMark F. Adams       pc->setupcalled = 0;
441806fa848SBarry Smith     } else {
44284d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
44303a628feSMark F. Adams       /* just do Galerkin grids */
44458471d46SMark F. Adams       Mat          B,dA,dB;
44558471d46SMark F. Adams 
44671959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
4479d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
44858471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
44923ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
45058471d46SMark F. Adams         /* (re)set to get dirty flag */
45123ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
45258471d46SMark F. Adams 
4532fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
45403a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
4550a97e771SToby Isaac           if (pc_gamg->setup_count==2) {
45603a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
457084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
4582fa5cd67SKarl Rupp 
45903a628feSMark F. Adams             mglevels[level]->A = B;
460806fa848SBarry Smith           } else {
46123ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
46258471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
46303a628feSMark F. Adams           }
46423ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
46558471d46SMark F. Adams           dB   = B;
46658471d46SMark F. Adams         }
4675f8cf99dSMark F. Adams       }
468d5280255SMark F. Adams 
469d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
47058471d46SMark F. Adams       PetscFunctionReturn(0);
471eb07cef2SMark F. Adams     }
472878e152fSMark F. Adams   }
473f6536408SMark F. Adams 
474878e152fSMark F. Adams   if (!pc_gamg->data) {
475878e152fSMark F. Adams     if (pc_gamg->orig_data) {
476878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
4770298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
4782fa5cd67SKarl Rupp 
479878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
480878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
481878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
4822fa5cd67SKarl Rupp 
483785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
484878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
485806fa848SBarry Smith     } else {
4861ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
4877700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
4889d5b6da9SMark F. Adams     }
489878e152fSMark F. Adams   }
490878e152fSMark F. Adams 
491878e152fSMark F. Adams   /* cache original data for reuse */
4921c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
493785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
494878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
495878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
496878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
497878e152fSMark F. Adams   }
498038e3b61SMark F. Adams 
499302f38e8SMark F. Adams   /* get basic dims */
500302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
501171cca9aSMark Adams   ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr);
50284d3f75bSMark F. Adams 
503569f4572SMark Adams   ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
504569f4572SMark Adams   nnz0   = info.nz_used;
505569f4572SMark Adams   nnztot = info.nz_used;
50662294041SBarry Smith   ierr = PetscInfo6(pc,"level %d) N=%D, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(int)(nnz0/(PetscReal)M+0.5),size);CHKERRQ(ierr);
507569f4572SMark Adams 
508a2f3521dSMark F. Adams   /* Get A_i and R_i */
50962294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
5109ab59c8bSMark Adams     pc_gamg->current_level = level;
5115b89ad90SMark F. Adams     level1 = level + 1;
5120cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5130cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
514a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
515a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
516b4fbaa2aSMark F. Adams #endif
517a2f3521dSMark F. Adams #endif
518c8b0795cSMark F. Adams     { /* construct prolongator */
519725b86d8SJed Brown       Mat              Gmat;
5200cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5217700e67bSMark Adams       Mat              Prol11;
522c8b0795cSMark F. Adams 
5237700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5241ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5257700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
526c8b0795cSMark F. Adams 
527a2f3521dSMark F. Adams       /* could have failed to create new level */
528a2f3521dSMark F. Adams       if (Prol11) {
5299d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5300298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
531a2f3521dSMark F. Adams 
532fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
533c8b0795cSMark F. Adams           /* smooth */
534fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
535c8b0795cSMark F. Adams         }
536c8b0795cSMark F. Adams 
5377700e67bSMark Adams         Parr[level1] = Prol11;
538171cca9aSMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
539ffc955d6SMark F. Adams 
5400c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
5411b18a24aSMark Adams         PetscInt bs;
5421b18a24aSMark Adams         ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
5430a3c815dSMark Adams         ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
544ffc955d6SMark F. Adams       }
545ffc955d6SMark F. Adams 
546a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
54741b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
548a2f3521dSMark F. Adams     } /* construct prolongator scope */
5490cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5500cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
551c8b0795cSMark F. Adams #endif
5527f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
553171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
554569f4572SMark Adams       ierr =  PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr);
55562294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
556a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
557a90e85d9SMark Adams #endif
558c8b0795cSMark F. Adams       break;
559c8b0795cSMark F. Adams     }
5600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5610cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
562b4fbaa2aSMark F. Adams #endif
563171cca9aSMark Adams     ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */
564171cca9aSMark Adams     if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ????????");
565171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
566171cca9aSMark Adams     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr);
567a2f3521dSMark F. Adams 
5680cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5690cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
570b4fbaa2aSMark F. Adams #endif
571171cca9aSMark Adams     ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */
572569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
573569f4572SMark Adams     nnztot += info.nz_used;
5741d5b2942SMark Adams     ierr = PetscInfo5(pc,"%d) N=%D, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",level1,M,pc_gamg->data_cell_cols,(int)(info.nz_used/(PetscReal)M),nactivepe);CHKERRQ(ierr);
575569f4572SMark Adams 
5760cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
577b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
578b4fbaa2aSMark F. Adams #endif
579a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
5809ab59c8bSMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) {
5819ab59c8bSMark Adams       ierr =  PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr);
582a90e85d9SMark Adams       level++;
583a90e85d9SMark Adams       break;
584a90e85d9SMark Adams     }
585c8b0795cSMark F. Adams   } /* levels */
586c8b0795cSMark F. Adams   ierr                  = PetscFree(pc_gamg->data);CHKERRQ(ierr);
587c8b0795cSMark F. Adams 
588569f4572SMark Adams   ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr);
5899d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
5905b89ad90SMark F. Adams   fine_level       = level;
5910298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
5925b89ad90SMark F. Adams 
59362294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
594d5280255SMark F. Adams     /* set default smoothers & set operators */
59562294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
596ffc955d6SMark F. Adams       KSP smoother;
597ffc955d6SMark F. Adams       PC  subpc;
598a2f3521dSMark F. Adams 
5999d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
600f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
601ffc955d6SMark F. Adams 
602a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
603a2f3521dSMark F. Adams       /* set ops */
60423ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
605a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
606a2f3521dSMark F. Adams 
607a2f3521dSMark F. Adams       /* set defaults */
6086c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
609a2f3521dSMark F. Adams 
6100c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
6110c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
6122d3561bbSSatish Balay         PetscInt sz;
6137a28f3e5SMark Adams         IS       *iss;
614a2f3521dSMark F. Adams 
6152d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
6167a28f3e5SMark Adams         iss   = ASMLocalIDsArr[level];
6170c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
6180a3c815dSMark Adams         ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr);
6190c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
6207f66b68fSMark Adams         if (!sz) {
621ffc955d6SMark F. Adams           IS       is;
6220a3c815dSMark Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
6237a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
624a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
625806fa848SBarry Smith         } else {
626a94c3b12SMark F. Adams           PetscInt kk;
6277a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr);
628a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
6297a28f3e5SMark Adams             ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr);
630a94c3b12SMark F. Adams           }
6317a28f3e5SMark Adams           ierr = PetscFree(iss);CHKERRQ(ierr);
632ffc955d6SMark F. Adams         }
6330298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
634ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
635806fa848SBarry Smith       } else {
636890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
637ffc955d6SMark F. Adams       }
638d5280255SMark F. Adams     }
639d5280255SMark F. Adams     {
640d5280255SMark F. Adams       /* coarse grid */
641d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
642d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
643d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
64423ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
645cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
646d5280255SMark F. Adams         ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
647d5280255SMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
648d5280255SMark F. Adams         ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
649d5280255SMark F. Adams         ierr = PCSetUp(subpc);CHKERRQ(ierr);
65071959b99SBarry Smith         ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
65171959b99SBarry Smith         if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
652d5280255SMark F. Adams         ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
653d5280255SMark F. Adams         ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
6549dbfc187SHong Zhang         ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
6552fb0b348SMark F. Adams         ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
65608e36f19SMark Adams         ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
6575b42dca8SJed Brown         /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
6585b42dca8SJed Brown          * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
6595b42dca8SJed Brown          * view every subdomain as though they were different. */
6605b42dca8SJed Brown         ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
661d5280255SMark F. Adams       }
662cf8ae1d3SMark Adams     }
663d5280255SMark F. Adams 
664d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
665d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
666e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
667d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
66869aca0b8SBarry Smith     ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
669d5280255SMark F. Adams 
670d5280255SMark F. Adams     /* clean up */
671d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
672587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
673587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
6745b89ad90SMark F. Adams     }
6750cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
676806fa848SBarry Smith   } else {
6775f8cf99dSMark F. Adams     KSP smoother;
678302440fdSBarry Smith     ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
6799d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
68023ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
6815f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
6829d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
6835f8cf99dSMark F. Adams   }
6845b89ad90SMark F. Adams   PetscFunctionReturn(0);
6855b89ad90SMark F. Adams }
6865b89ad90SMark F. Adams 
687eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
6885b89ad90SMark F. Adams /*
6895b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
6905b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
6915b89ad90SMark F. Adams 
6925b89ad90SMark F. Adams    Input Parameter:
6935b89ad90SMark F. Adams .  pc - the preconditioner context
6945b89ad90SMark F. Adams 
6955b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
6965b89ad90SMark F. Adams */
6975b89ad90SMark F. Adams #undef __FUNCT__
6985b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
6995b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
7005b89ad90SMark F. Adams {
7015b89ad90SMark F. Adams   PetscErrorCode ierr;
7025b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
7035b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
7045b89ad90SMark F. Adams 
7055b89ad90SMark F. Adams   PetscFunctionBegin;
7065b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7079b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
7089b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
7099b8ffb57SJed Brown   }
7101ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
7111ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
7125b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7135b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7145b89ad90SMark F. Adams   PetscFunctionReturn(0);
7155b89ad90SMark F. Adams }
7165b89ad90SMark F. Adams 
717676e1743SMark F. Adams #undef __FUNCT__
718676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
719676e1743SMark F. Adams /*@
720cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
721676e1743SMark F. Adams 
7221cc46a46SBarry Smith    Logically Collective on PC
723676e1743SMark F. Adams 
724676e1743SMark F. Adams    Input Parameters:
7251cc46a46SBarry Smith +  pc - the preconditioner context
7261cc46a46SBarry Smith -  n - the number of equations
727676e1743SMark F. Adams 
728676e1743SMark F. Adams 
729676e1743SMark F. Adams    Options Database Key:
7301cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
731676e1743SMark F. Adams 
732cab9ed1eSBarry Smith    Notes: GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
733cab9ed1eSBarry Smith           that has degrees of freedom
734cab9ed1eSBarry Smith 
735676e1743SMark F. Adams    Level: intermediate
736676e1743SMark F. Adams 
7371c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
738676e1743SMark F. Adams 
7391c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim()
740676e1743SMark F. Adams @*/
741676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
742676e1743SMark F. Adams {
743676e1743SMark F. Adams   PetscErrorCode ierr;
744676e1743SMark F. Adams 
745676e1743SMark F. Adams   PetscFunctionBegin;
746676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
747676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
748676e1743SMark F. Adams   PetscFunctionReturn(0);
749676e1743SMark F. Adams }
750676e1743SMark F. Adams 
751676e1743SMark F. Adams #undef __FUNCT__
752676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
7531e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
754676e1743SMark F. Adams {
755c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
756c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
757676e1743SMark F. Adams 
758676e1743SMark F. Adams   PetscFunctionBegin;
7599d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
760676e1743SMark F. Adams   PetscFunctionReturn(0);
761676e1743SMark F. Adams }
762676e1743SMark F. Adams 
763676e1743SMark F. Adams #undef __FUNCT__
764389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
765389730f3SMark F. Adams /*@
766cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
767389730f3SMark F. Adams 
768389730f3SMark F. Adams  Collective on PC
769389730f3SMark F. Adams 
770389730f3SMark F. Adams    Input Parameters:
7711cc46a46SBarry Smith +  pc - the preconditioner context
7721cc46a46SBarry Smith -  n - maximum number of equations to aim for
773389730f3SMark F. Adams 
774389730f3SMark F. Adams    Options Database Key:
7751cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
776389730f3SMark F. Adams 
777389730f3SMark F. Adams    Level: intermediate
778389730f3SMark F. Adams 
7791c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
780389730f3SMark F. Adams 
7811c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim()
782389730f3SMark F. Adams @*/
783389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
784389730f3SMark F. Adams {
785389730f3SMark F. Adams   PetscErrorCode ierr;
786389730f3SMark F. Adams 
787389730f3SMark F. Adams   PetscFunctionBegin;
788389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
789389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
790389730f3SMark F. Adams   PetscFunctionReturn(0);
791389730f3SMark F. Adams }
792389730f3SMark F. Adams 
793389730f3SMark F. Adams #undef __FUNCT__
794389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
7951e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
796389730f3SMark F. Adams {
797389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
798389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
799389730f3SMark F. Adams 
800389730f3SMark F. Adams   PetscFunctionBegin;
8019d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
802389730f3SMark F. Adams   PetscFunctionReturn(0);
803389730f3SMark F. Adams }
804389730f3SMark F. Adams 
805389730f3SMark F. Adams #undef __FUNCT__
806cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGSetRepartition"
807676e1743SMark F. Adams /*@
808cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
809676e1743SMark F. Adams 
810676e1743SMark F. Adams    Collective on PC
811676e1743SMark F. Adams 
812676e1743SMark F. Adams    Input Parameters:
8131cc46a46SBarry Smith +  pc - the preconditioner context
8141cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
815676e1743SMark F. Adams 
816676e1743SMark F. Adams    Options Database Key:
8171cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
818676e1743SMark F. Adams 
819cab9ed1eSBarry Smith    Notes: this will generally improve the loading balancing of the work on each level
820cab9ed1eSBarry Smith 
821676e1743SMark F. Adams    Level: intermediate
822676e1743SMark F. Adams 
8231c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
824676e1743SMark F. Adams 
825676e1743SMark F. Adams .seealso: ()
826676e1743SMark F. Adams @*/
827cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
828676e1743SMark F. Adams {
829676e1743SMark F. Adams   PetscErrorCode ierr;
830676e1743SMark F. Adams 
831676e1743SMark F. Adams   PetscFunctionBegin;
832676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
833cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
834676e1743SMark F. Adams   PetscFunctionReturn(0);
835676e1743SMark F. Adams }
836676e1743SMark F. Adams 
837676e1743SMark F. Adams #undef __FUNCT__
838cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGSetRepartition_GAMG"
839cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
840676e1743SMark F. Adams {
841c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
842c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
843676e1743SMark F. Adams 
844676e1743SMark F. Adams   PetscFunctionBegin;
8459d5b6da9SMark F. Adams   pc_gamg->repart = n;
846676e1743SMark F. Adams   PetscFunctionReturn(0);
847676e1743SMark F. Adams }
848676e1743SMark F. Adams 
849676e1743SMark F. Adams #undef __FUNCT__
8501cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation"
851dfd5c07aSMark F. Adams /*@
852cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
853dfd5c07aSMark F. Adams 
854dfd5c07aSMark F. Adams    Collective on PC
855dfd5c07aSMark F. Adams 
856dfd5c07aSMark F. Adams    Input Parameters:
8571cc46a46SBarry Smith +  pc - the preconditioner context
8581cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
859dfd5c07aSMark F. Adams 
860dfd5c07aSMark F. Adams    Options Database Key:
8611cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
862dfd5c07aSMark F. Adams 
863dfd5c07aSMark F. Adams    Level: intermediate
864dfd5c07aSMark F. Adams 
865cab9ed1eSBarry Smith    Notes: this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
866cab9ed1eSBarry Smith           rebuilding the preconditioner quicker.
867cab9ed1eSBarry Smith 
8681c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
869dfd5c07aSMark F. Adams 
870dfd5c07aSMark F. Adams .seealso: ()
871dfd5c07aSMark F. Adams @*/
8721cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
873dfd5c07aSMark F. Adams {
874dfd5c07aSMark F. Adams   PetscErrorCode ierr;
875dfd5c07aSMark F. Adams 
876dfd5c07aSMark F. Adams   PetscFunctionBegin;
877dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8781cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
879dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
880dfd5c07aSMark F. Adams }
881dfd5c07aSMark F. Adams 
882dfd5c07aSMark F. Adams #undef __FUNCT__
8831cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG"
8841cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
885dfd5c07aSMark F. Adams {
886dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
887dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
888dfd5c07aSMark F. Adams 
889dfd5c07aSMark F. Adams   PetscFunctionBegin;
890dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
891dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
892dfd5c07aSMark F. Adams }
893dfd5c07aSMark F. Adams 
894dfd5c07aSMark F. Adams #undef __FUNCT__
895cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGASMSetUseAggs"
896ffc955d6SMark F. Adams /*@
897cab9ed1eSBarry Smith    PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner.
898ffc955d6SMark F. Adams 
899ffc955d6SMark F. Adams    Collective on PC
900ffc955d6SMark F. Adams 
901ffc955d6SMark F. Adams    Input Parameters:
902cab9ed1eSBarry Smith +  pc - the preconditioner context
903cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
904ffc955d6SMark F. Adams 
905ffc955d6SMark F. Adams    Options Database Key:
906cab9ed1eSBarry Smith .  -pc_gamg_asm_use_agg
907ffc955d6SMark F. Adams 
908ffc955d6SMark F. Adams    Level: intermediate
909ffc955d6SMark F. Adams 
9101c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
911ffc955d6SMark F. Adams 
912ffc955d6SMark F. Adams .seealso: ()
913ffc955d6SMark F. Adams @*/
914cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
915ffc955d6SMark F. Adams {
916ffc955d6SMark F. Adams   PetscErrorCode ierr;
917ffc955d6SMark F. Adams 
918ffc955d6SMark F. Adams   PetscFunctionBegin;
919ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
920cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
921ffc955d6SMark F. Adams   PetscFunctionReturn(0);
922ffc955d6SMark F. Adams }
923ffc955d6SMark F. Adams 
924ffc955d6SMark F. Adams #undef __FUNCT__
925cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGASMSetUseAggs_GAMG"
926cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
927ffc955d6SMark F. Adams {
928ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
929ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
930ffc955d6SMark F. Adams 
931ffc955d6SMark F. Adams   PetscFunctionBegin;
932cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
933ffc955d6SMark F. Adams   PetscFunctionReturn(0);
934ffc955d6SMark F. Adams }
935ffc955d6SMark F. Adams 
936ffc955d6SMark F. Adams #undef __FUNCT__
937171cca9aSMark Adams #define __FUNCT__ "PCGAMGSetUseParallelCoarseGridSolve"
938171cca9aSMark Adams /*@
939cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
940171cca9aSMark Adams 
941171cca9aSMark Adams    Collective on PC
942171cca9aSMark Adams 
943171cca9aSMark Adams    Input Parameters:
944171cca9aSMark Adams +  pc - the preconditioner context
945cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
946171cca9aSMark Adams 
947171cca9aSMark Adams    Options Database Key:
948cf8ae1d3SMark Adams .  -pc_gamg_use_parallel_coarse_grid_solver
949171cca9aSMark Adams 
950171cca9aSMark Adams    Level: intermediate
951171cca9aSMark Adams 
952171cca9aSMark Adams    Concepts: Unstructured multigrid preconditioner
953171cca9aSMark Adams 
954171cca9aSMark Adams .seealso: ()
955171cca9aSMark Adams @*/
956171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
957171cca9aSMark Adams {
958171cca9aSMark Adams   PetscErrorCode ierr;
959171cca9aSMark Adams 
960171cca9aSMark Adams   PetscFunctionBegin;
961171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
962171cca9aSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
963171cca9aSMark Adams   PetscFunctionReturn(0);
964171cca9aSMark Adams }
965171cca9aSMark Adams 
966171cca9aSMark Adams #undef __FUNCT__
967171cca9aSMark Adams #define __FUNCT__ "PCGAMGSetUseParallelCoarseGridSolve_GAMG"
968171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
969171cca9aSMark Adams {
970171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
971171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
972171cca9aSMark Adams 
973171cca9aSMark Adams   PetscFunctionBegin;
974171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
975676e1743SMark F. Adams   PetscFunctionReturn(0);
976676e1743SMark F. Adams }
977676e1743SMark F. Adams 
978676e1743SMark F. Adams #undef __FUNCT__
9794ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
9804ef23d27SMark F. Adams /*@
9811cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
9824ef23d27SMark F. Adams 
9834ef23d27SMark F. Adams    Not collective on PC
9844ef23d27SMark F. Adams 
9854ef23d27SMark F. Adams    Input Parameters:
9861cc46a46SBarry Smith +  pc - the preconditioner
9871cc46a46SBarry Smith -  n - the maximum number of levels to use
9884ef23d27SMark F. Adams 
9894ef23d27SMark F. Adams    Options Database Key:
9904ef23d27SMark F. Adams .  -pc_mg_levels
9914ef23d27SMark F. Adams 
9924ef23d27SMark F. Adams    Level: intermediate
9934ef23d27SMark F. Adams 
9941c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
9954ef23d27SMark F. Adams 
9964ef23d27SMark F. Adams .seealso: ()
9974ef23d27SMark F. Adams @*/
9984ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
9994ef23d27SMark F. Adams {
10004ef23d27SMark F. Adams   PetscErrorCode ierr;
10014ef23d27SMark F. Adams 
10024ef23d27SMark F. Adams   PetscFunctionBegin;
10034ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10044ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
10054ef23d27SMark F. Adams   PetscFunctionReturn(0);
10064ef23d27SMark F. Adams }
10074ef23d27SMark F. Adams 
10084ef23d27SMark F. Adams #undef __FUNCT__
10094ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
10101e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
10114ef23d27SMark F. Adams {
10124ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
10134ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
10144ef23d27SMark F. Adams 
10154ef23d27SMark F. Adams   PetscFunctionBegin;
10169d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
10174ef23d27SMark F. Adams   PetscFunctionReturn(0);
10184ef23d27SMark F. Adams }
10194ef23d27SMark F. Adams 
10204ef23d27SMark F. Adams #undef __FUNCT__
10213542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
10223542efc5SMark F. Adams /*@
10233542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
10243542efc5SMark F. Adams 
10253542efc5SMark F. Adams    Not collective on PC
10263542efc5SMark F. Adams 
10273542efc5SMark F. Adams    Input Parameters:
10281cc46a46SBarry Smith +  pc - the preconditioner context
1029b001cb0fSBarry Smith -  threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
10303542efc5SMark F. Adams 
10313542efc5SMark F. Adams    Options Database Key:
10321cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
10333542efc5SMark F. Adams 
1034cab9ed1eSBarry Smith    Notes: Before aggregating the graph GAMG will remove small values from the graph thus reducing the coupling in the graph and a different
1035cab9ed1eSBarry Smith     (perhaps better) coarser set of points.
1036cab9ed1eSBarry Smith 
10373542efc5SMark F. Adams    Level: intermediate
10383542efc5SMark F. Adams 
10391c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
10403542efc5SMark F. Adams 
10413542efc5SMark F. Adams .seealso: ()
10423542efc5SMark F. Adams @*/
1043*c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
10443542efc5SMark F. Adams {
10453542efc5SMark F. Adams   PetscErrorCode ierr;
10463542efc5SMark F. Adams 
10473542efc5SMark F. Adams   PetscFunctionBegin;
10483542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1049*c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr);
10503542efc5SMark F. Adams   PetscFunctionReturn(0);
10513542efc5SMark F. Adams }
10523542efc5SMark F. Adams 
10533542efc5SMark F. Adams #undef __FUNCT__
10543542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
1055*c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
10563542efc5SMark F. Adams {
1057c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1058c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1059*c1eae691SMark Adams   PetscInt i;
1060*c1eae691SMark Adams   PetscFunctionBegin;
1061*c1eae691SMark Adams   for (i=0;i<n;i++) pc_gamg->threshold[i] = v[i];
1062*c1eae691SMark Adams   do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1063*c1eae691SMark Adams   PetscFunctionReturn(0);
1064*c1eae691SMark Adams }
1065*c1eae691SMark Adams 
1066*c1eae691SMark Adams #undef __FUNCT__
1067*c1eae691SMark Adams #define __FUNCT__ "PCGAMGSetThresholdScale"
1068*c1eae691SMark Adams /*@
1069*c1eae691SMark Adams    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1070*c1eae691SMark Adams 
1071*c1eae691SMark Adams    Not collective on PC
1072*c1eae691SMark Adams 
1073*c1eae691SMark Adams    Input Parameters:
1074*c1eae691SMark Adams +  pc - the preconditioner context
1075*c1eae691SMark Adams -  scale - the threshold value reduction, ussually < 1.0
1076*c1eae691SMark Adams 
1077*c1eae691SMark Adams    Options Database Key:
1078*c1eae691SMark Adams .  -pc_gamg_threshold_scale <v>
1079*c1eae691SMark Adams 
1080*c1eae691SMark Adams    Level: advanced
1081*c1eae691SMark Adams 
1082*c1eae691SMark Adams    Concepts: Unstructured multigrid preconditioner
1083*c1eae691SMark Adams 
1084*c1eae691SMark Adams .seealso: ()
1085*c1eae691SMark Adams @*/
1086*c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1087*c1eae691SMark Adams {
1088*c1eae691SMark Adams   PetscErrorCode ierr;
10893542efc5SMark F. Adams 
10903542efc5SMark F. Adams   PetscFunctionBegin;
1091*c1eae691SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1092*c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr);
1093*c1eae691SMark Adams   PetscFunctionReturn(0);
1094*c1eae691SMark Adams }
1095*c1eae691SMark Adams 
1096*c1eae691SMark Adams #undef __FUNCT__
1097*c1eae691SMark Adams #define __FUNCT__ "PCGAMGSetThresholdScale_GAMG"
1098*c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1099*c1eae691SMark Adams {
1100*c1eae691SMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1101*c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1102*c1eae691SMark Adams   PetscFunctionBegin;
1103*c1eae691SMark Adams   pc_gamg->threshold_scale = v;
11043542efc5SMark F. Adams   PetscFunctionReturn(0);
11053542efc5SMark F. Adams }
11063542efc5SMark F. Adams 
11073542efc5SMark F. Adams #undef __FUNCT__
11089d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1109e20c40e8SBarry Smith /*@C
1110c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1111676e1743SMark F. Adams 
1112676e1743SMark F. Adams    Collective on PC
1113676e1743SMark F. Adams 
1114676e1743SMark F. Adams    Input Parameters:
1115c60c7ad4SBarry Smith +  pc - the preconditioner context
1116c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1117676e1743SMark F. Adams 
1118676e1743SMark F. Adams    Options Database Key:
1119cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1120676e1743SMark F. Adams 
1121676e1743SMark F. Adams    Level: intermediate
1122676e1743SMark F. Adams 
11231c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1124676e1743SMark F. Adams 
1125cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1126676e1743SMark F. Adams @*/
112719fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1128676e1743SMark F. Adams {
1129676e1743SMark F. Adams   PetscErrorCode ierr;
1130676e1743SMark F. Adams 
1131676e1743SMark F. Adams   PetscFunctionBegin;
1132676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1133806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1134676e1743SMark F. Adams   PetscFunctionReturn(0);
1135676e1743SMark F. Adams }
1136676e1743SMark F. Adams 
1137676e1743SMark F. Adams #undef __FUNCT__
1138c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType"
1139e20c40e8SBarry Smith /*@C
1140c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1141c60c7ad4SBarry Smith 
1142c60c7ad4SBarry Smith    Collective on PC
1143c60c7ad4SBarry Smith 
1144c60c7ad4SBarry Smith    Input Parameter:
1145c60c7ad4SBarry Smith .  pc - the preconditioner context
1146c60c7ad4SBarry Smith 
1147c60c7ad4SBarry Smith    Output Parameter:
1148c60c7ad4SBarry Smith .  type - the type of algorithm used
1149c60c7ad4SBarry Smith 
1150c60c7ad4SBarry Smith    Level: intermediate
1151c60c7ad4SBarry Smith 
11521c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1153c60c7ad4SBarry Smith 
11541c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1155c60c7ad4SBarry Smith @*/
1156c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1157c60c7ad4SBarry Smith {
1158c60c7ad4SBarry Smith   PetscErrorCode ierr;
1159c60c7ad4SBarry Smith 
1160c60c7ad4SBarry Smith   PetscFunctionBegin;
1161c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1162c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1163c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1164c60c7ad4SBarry Smith }
1165c60c7ad4SBarry Smith 
1166c60c7ad4SBarry Smith #undef __FUNCT__
1167c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG"
1168c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1169c60c7ad4SBarry Smith {
1170c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1171c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1172c60c7ad4SBarry Smith 
1173c60c7ad4SBarry Smith   PetscFunctionBegin;
1174c60c7ad4SBarry Smith   *type = pc_gamg->type;
1175c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1176c60c7ad4SBarry Smith }
1177c60c7ad4SBarry Smith 
1178c60c7ad4SBarry Smith #undef __FUNCT__
11799d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
11801e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1181676e1743SMark F. Adams {
11829d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
11831ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
11841ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1185676e1743SMark F. Adams 
1186676e1743SMark F. Adams   PetscFunctionBegin;
1187c60c7ad4SBarry Smith   pc_gamg->type = type;
11881c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
11899d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
11901ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
11911ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
11921ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1193e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
11943ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
11953ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
11963ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
11973ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
11983ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
11993ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
12003ae0bb68SMark Adams     pc_gamg->data_sz = 0;
12011ab5ffc9SJed Brown   }
12021ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
12031ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
12049d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1205676e1743SMark F. Adams   PetscFunctionReturn(0);
1206676e1743SMark F. Adams }
1207676e1743SMark F. Adams 
12085b89ad90SMark F. Adams #undef __FUNCT__
12095adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG"
12105adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
12115adeb434SBarry Smith {
1212*c1eae691SMark Adams   PetscErrorCode ierr,i;
12135adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
12145adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
12155adeb434SBarry Smith 
12165adeb434SBarry Smith   PetscFunctionBegin;
12175adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1218*c1eae691SMark Adams   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");
1219*c1eae691SMark Adams   for (i=0;i<pc_gamg->current_level;i++) {
1220*c1eae691SMark Adams     ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr);
1221*c1eae691SMark Adams   }
1222*c1eae691SMark Adams   ierr = PetscViewerASCIIPrintf(viewer,"\n");
1223*c1eae691SMark Adams   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g",(double)pc_gamg->threshold_scale);
1224cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
1225cab9ed1eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr);
1226cab9ed1eSBarry Smith   }
1227171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
1228171cca9aSMark Adams     ierr = PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr);
1229171cca9aSMark Adams   }
12305adeb434SBarry Smith   if (pc_gamg->ops->view) {
12315adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
12325adeb434SBarry Smith   }
12335adeb434SBarry Smith   PetscFunctionReturn(0);
12345adeb434SBarry Smith }
12355adeb434SBarry Smith 
12365adeb434SBarry Smith #undef __FUNCT__
12375b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
12384416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
12395b89ad90SMark F. Adams {
1240676e1743SMark F. Adams   PetscErrorCode ierr;
1241676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1242676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1243676e1743SMark F. Adams   PetscBool      flag;
12443b4367a7SBarry Smith   MPI_Comm       comm;
124514a9496bSBarry Smith   char           prefix[256];
1246*c1eae691SMark Adams   PetscInt       i,n;
124714a9496bSBarry Smith   const char     *pcpre;
12485b89ad90SMark F. Adams 
12495b89ad90SMark F. Adams   PetscFunctionBegin;
12503b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1251e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1252676e1743SMark F. Adams   {
1253bd94a7aaSJed Brown     char tname[256];
12541a1c1e04SBarry Smith     ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1255bd94a7aaSJed Brown     if (flag) {
1256bd94a7aaSJed Brown       ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
12571ab5ffc9SJed Brown     }
1258cab9ed1eSBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
12591cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1260cab9ed1eSBarry Smith     ierr = PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation agragates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL);CHKERRQ(ierr);
1261cf8ae1d3SMark Adams     ierr = PetscOptionsBool("-pc_gamg_use_parallel_coarse_grid_solver","Use parallel coarse grid solver (otherwise put last grid on one process)","PCGAMGSetUseParallelCoarseGridSolve",pc_gamg->use_parallel_coarse_grid_solver,&pc_gamg->use_parallel_coarse_grid_solver,NULL);CHKERRQ(ierr);
126294ae4db5SBarry 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);
126394ae4db5SBarry 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);
1264*c1eae691SMark Adams     n = PETSC_GAMG_MAXLEVELS;
1265*c1eae691SMark Adams     ierr = PetscOptionsReal("-pc_gamg_threshold_scale","scaleing of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,&flag);CHKERRQ(ierr);
1266*c1eae691SMark Adams     ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr);
1267*c1eae691SMark Adams     if (flag && n < PETSC_GAMG_MAXLEVELS) {
1268*c1eae691SMark Adams       i = n;
1269*c1eae691SMark Adams       do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1270*c1eae691SMark Adams     }
127194ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1272b7cbab4eSMark Adams 
1273b7cbab4eSMark Adams     /* set options for subtype */
1274e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1275676e1743SMark F. Adams   }
127614a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
127714a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
1278676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
12795b89ad90SMark F. Adams   PetscFunctionReturn(0);
12805b89ad90SMark F. Adams }
12815b89ad90SMark F. Adams 
12825b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
12835b89ad90SMark F. Adams /*MC
12841cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
12855b89ad90SMark F. Adams 
1286280d9858SJed Brown    Options Database Keys:
1287cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1288cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1289cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1290cab9ed1eSBarry Smith .   -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother
1291cab9ed1eSBarry Smith .   -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit>
1292cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1293cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
1294*c1eae691SMark Adams -   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1295*c1eae691SMark Adams -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1296cab9ed1eSBarry Smith 
1297cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1298cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1299cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1300cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1301cab9ed1eSBarry Smith 
1302cab9ed1eSBarry Smith    Multigrid options(inherited):
13031cc46a46SBarry Smith +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1304280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1305280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1306cab9ed1eSBarry Smith .  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
1307cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
13085b89ad90SMark F. Adams 
13091cc46a46SBarry Smith 
13101cc46a46SBarry Smith   Notes: In order to obtain good performance for PCGAMG for vector valued problems you must
13111cc46a46SBarry Smith $       Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
13121cc46a46SBarry Smith $       Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
13131cc46a46SBarry Smith $       See the Users Manual Chapter 4 for more details
13141cc46a46SBarry Smith 
13155b89ad90SMark F. Adams   Level: intermediate
1316280d9858SJed Brown 
13171cc46a46SBarry Smith   Concepts: algebraic multigrid
13185b89ad90SMark F. Adams 
13191cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1320171cca9aSMark Adams            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation()
13215b89ad90SMark F. Adams M*/
1322b2573a8aSBarry Smith 
13235b89ad90SMark F. Adams #undef __FUNCT__
13245b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
13258cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
13265b89ad90SMark F. Adams {
1327*c1eae691SMark Adams   PetscErrorCode ierr,i;
13285b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
13295b89ad90SMark F. Adams   PC_MG          *mg;
13305b89ad90SMark F. Adams 
13315b89ad90SMark F. Adams   PetscFunctionBegin;
13321c1aac46SBarry Smith    /* register AMG type */
13331c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
13341c1aac46SBarry Smith 
13355b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
13361c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
13375b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
13385b89ad90SMark F. Adams 
13395b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1340b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
134169aca0b8SBarry Smith   ierr         = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
13425b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
13435b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
13445b89ad90SMark F. Adams 
1345b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
13461ab5ffc9SJed Brown 
13479d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
13489d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
13499d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
13509d5b6da9SMark F. Adams   pc_gamg->data    = 0;
13515b89ad90SMark F. Adams 
13529d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
13535b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
13545b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
13555b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
13565b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
13575adeb434SBarry Smith   mg->view                = PCView_GAMG;
13585b89ad90SMark F. Adams 
1359bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1360bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1361cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr);
13621cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1363cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr);
1364171cca9aSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr);
1365bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1366*c1eae691SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr);
1367bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1368c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1369bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
13709d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1371d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
13720c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1373171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1374038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
137525a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1376*c1eae691SMark Adams   for (i=0;i<PETSC_GAMG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1377*c1eae691SMark Adams   pc_gamg->threshold_scale = 1.;
1378*c1eae691SMark Adams   pc_gamg->Nlevels          = PETSC_GAMG_MAXLEVELS;
13799ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
1380c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
13819d5b6da9SMark F. Adams 
1382bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1383bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
13845b89ad90SMark F. Adams   PetscFunctionReturn(0);
13855b89ad90SMark F. Adams }
13863e3471ccSMark Adams 
13873e3471ccSMark Adams #undef __FUNCT__
13883e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage"
13893e3471ccSMark Adams /*@C
13903e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
13913e3471ccSMark Adams     from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
13923e3471ccSMark Adams     when using static libraries.
13933e3471ccSMark Adams 
13943e3471ccSMark Adams  Level: developer
13953e3471ccSMark Adams 
13963e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
13973e3471ccSMark Adams  .seealso: PetscInitialize()
13983e3471ccSMark Adams @*/
13993e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
14003e3471ccSMark Adams {
14013e3471ccSMark Adams   PetscErrorCode ierr;
14023e3471ccSMark Adams 
14033e3471ccSMark Adams   PetscFunctionBegin;
14043e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
14053e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
14063e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
14073e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
14088e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
14093e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1410c1c463dbSMark Adams 
1411c1c463dbSMark Adams   /* general events */
1412fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1413fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1414fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1415fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1416c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1417c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1418fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1419c1c463dbSMark Adams 
14205b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14215b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
14225b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
14235b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
14245b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
14255b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
14265b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
14275b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
14285b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1429bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
14305b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
14315b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
14325b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
14335b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
14345b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
14355b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
14365b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
14375b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
14385b89ad90SMark F. Adams 
14395b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
14405b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
14415b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
14425b89ad90SMark F. Adams   /* create timer stages */
14435b89ad90SMark F. Adams #if defined GAMG_STAGES
14445b89ad90SMark F. Adams   {
14455b89ad90SMark F. Adams     char     str[32];
14465b89ad90SMark F. Adams     PetscInt lidx;
14475b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
14485b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
14495b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
14505b89ad90SMark F. Adams       sprintf(str,"MG Level %d",lidx);
14515b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
14525b89ad90SMark F. Adams     }
14535b89ad90SMark F. Adams   }
14545b89ad90SMark F. Adams #endif
14555b89ad90SMark F. Adams #endif
14563e3471ccSMark Adams   PetscFunctionReturn(0);
14573e3471ccSMark Adams }
14583e3471ccSMark Adams 
14593e3471ccSMark Adams #undef __FUNCT__
14603e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage"
14613e3471ccSMark Adams /*@C
14621c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
14631c1aac46SBarry Smith     called from PetscFinalize() automatically.
14643e3471ccSMark Adams 
14653e3471ccSMark Adams  Level: developer
14663e3471ccSMark Adams 
14673e3471ccSMark Adams  .keywords: Petsc, destroy, package
14683e3471ccSMark Adams  .seealso: PetscFinalize()
14693e3471ccSMark Adams @*/
14703e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
14713e3471ccSMark Adams {
14723e3471ccSMark Adams   PetscErrorCode ierr;
14733e3471ccSMark Adams 
14743e3471ccSMark Adams   PetscFunctionBegin;
14753e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
14763e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
14773e3471ccSMark Adams   PetscFunctionReturn(0);
14783e3471ccSMark Adams }
1479a36cf38bSToby Isaac 
1480a36cf38bSToby Isaac #undef __FUNCT__
1481a36cf38bSToby Isaac #define __FUNCT__ "PCGAMGRegister"
1482a36cf38bSToby Isaac /*@C
1483a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1484a36cf38bSToby Isaac 
1485a36cf38bSToby Isaac  Input Parameters:
1486a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1487a36cf38bSToby Isaac  - create - function for creating the gamg context.
1488a36cf38bSToby Isaac 
1489a36cf38bSToby Isaac   Level: advanced
1490a36cf38bSToby Isaac 
14911c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1492a36cf38bSToby Isaac @*/
1493a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1494a36cf38bSToby Isaac {
1495a36cf38bSToby Isaac   PetscErrorCode ierr;
1496a36cf38bSToby Isaac 
1497a36cf38bSToby Isaac   PetscFunctionBegin;
1498a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1499a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1500a36cf38bSToby Isaac   PetscFunctionReturn(0);
1501a36cf38bSToby Isaac }
1502a36cf38bSToby Isaac 
1503