xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 7a28f3e5b83784e012fb5372ca74592ce67c00ce)
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 
22b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
23b4fbaa2aSMark F. Adams 
24b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
250cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
26b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
27b4fbaa2aSMark F. Adams #endif
28f96513f1SMatthew G Knepley 
29140e18c1SBarry Smith static PetscFunctionList GAMGList = 0;
303e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
319d5b6da9SMark F. Adams 
32d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
33d3d6bff4SMark F. Adams #undef __FUNCT__
34d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
35d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
36d3d6bff4SMark F. Adams {
37d3d6bff4SMark F. Adams   PetscErrorCode ierr;
38d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
39d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
40d3d6bff4SMark F. Adams 
41d3d6bff4SMark F. Adams   PetscFunctionBegin;
4262294041SBarry Smith   if (pc_gamg->data) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen, cleaned up in SetUp\n");
431c1aac46SBarry Smith   pc_gamg->data_sz = 0;
44878e152fSMark F. Adams   ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
45a2f3521dSMark F. Adams   PetscFunctionReturn(0);
46a2f3521dSMark F. Adams }
47a2f3521dSMark F. Adams 
485b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
495b89ad90SMark F. Adams /*
50c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
51a147abb0SMark F. Adams      of active processors.
525b89ad90SMark F. Adams 
535b89ad90SMark F. Adams    Input Parameter:
54a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
55a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
569d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
57c5bfad50SMark F. Adams    . cr_bs - coarse block size
583530afc2SMark F. Adams    In/Output Parameter:
59a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
60afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
6111e60469SMark F. Adams    Output Parameter:
623530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
635b89ad90SMark F. Adams */
645cb416c2SMark F. Adams 
655b89ad90SMark F. Adams #undef __FUNCT__
66c238b0ebSToby Isaac #define __FUNCT__ "PCGAMGCreateLevel_GAMG"
6762294041SBarry Smith 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)
685b89ad90SMark F. Adams {
69a2f3521dSMark F. Adams   PetscErrorCode  ierr;
709d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
71486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
72a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
733b4367a7SBarry Smith   MPI_Comm        comm;
74c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
753ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
765b89ad90SMark F. Adams 
775b89ad90SMark F. Adams   PetscFunctionBegin;
783b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
793b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
803b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
81c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
829d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
83038e3b61SMark F. Adams 
843ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
850298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
863ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
873ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
8873911c69SBarry Smith   } else {
893ae0bb68SMark Adams     PetscInt  bs;
903ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
913ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
923ae0bb68SMark Adams   }
93a2f3521dSMark F. Adams 
94c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
95a2f3521dSMark F. Adams   {
96472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
970298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
98a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
997f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
100c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
101a2f3521dSMark F. Adams   }
102f852f58cSMark F. Adams 
1033cb8563fSToby Isaac   if (Pcolumnperm) *Pcolumnperm = NULL;
1043cb8563fSToby Isaac 
105a90e85d9SMark Adams   if (!pc_gamg->repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1062fa5cd67SKarl Rupp   else {
1073ae0bb68SMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
108885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
109e33ef3b1SMark F. Adams 
11071959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
11171959b99SBarry 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);
1120cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1130cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
114b4fbaa2aSMark F. Adams #endif
115a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
116785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
117a90e85d9SMark Adams     if (pc_gamg->repart) {
118a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1195a9b9e01SMark F. Adams       Mat adj;
1205a9b9e01SMark F. Adams 
121302440fdSBarry Smith      ierr = PetscInfo3(pc,"Repartition: size (active): %D --> %D, neq = %D\n",*a_nactive_proc,new_size,ncrs_eq);CHKERRQ(ierr);
1225a9b9e01SMark F. Adams 
123a2f3521dSMark F. Adams       /* get 'adj' */
124c5bfad50SMark F. Adams       if (cr_bs == 1) {
125038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
126806fa848SBarry Smith       } else {
127a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
128eb07cef2SMark F. Adams         Mat               tMat;
129a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
130b4fbaa2aSMark F. Adams         const PetscScalar *vals;
131b4fbaa2aSMark F. Adams         const PetscInt    *idx;
132a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1339057884aSMark F. Adams         static PetscInt   llev = 0;
134d9558ea9SBarry Smith         MatType           mtype;
135b4fbaa2aSMark F. Adams 
136e632b94dSBarry Smith         ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr);
137a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
138a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
139c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
14058471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
141c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
142c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
14358471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1443ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1453ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
14658471d46SMark F. Adams         }
1476876a03eSMark F. Adams 
148d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
1493b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1503ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
151d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
152a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
153a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
154e632b94dSBarry Smith         ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
155eb07cef2SMark F. Adams 
156a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
157c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
15822063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
159eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
160c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
161eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
162eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
163eb07cef2SMark F. Adams           }
16422063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
165eb07cef2SMark F. Adams         }
166eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168eb07cef2SMark F. Adams 
169b4fbaa2aSMark F. Adams         if (llev++ == -1) {
170b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
1718caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
1723b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
173b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
1743bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
175b4fbaa2aSMark F. Adams         }
176eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
177eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
178a2f3521dSMark F. Adams       } /* create 'adj' */
179f150b916SMark F. Adams 
180a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
1815a9b9e01SMark F. Adams         char            prefix[256];
1825a9b9e01SMark F. Adams         const char      *pcpre;
183b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
184b4fbaa2aSMark F. Adams         MatPartitioning mpart;
185a4b7d37bSMark F. Adams         IS              proc_is;
186a2f3521dSMark F. Adams         PetscInt        targetPE;
1872f03bc48SMark F. Adams 
1883b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
1895ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
1909d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
1918caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
19259a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
19311e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
194c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
195a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
19611e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1975a9b9e01SMark F. Adams 
1985ef31b24SMark F. Adams         /* collect IS info */
199785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
200a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
201a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
202c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
203a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
204c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
205a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
206eb07cef2SMark F. Adams           }
2075ef31b24SMark F. Adams         }
208a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
209a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2105ef31b24SMark F. Adams       }
2115ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2125a9b9e01SMark F. Adams 
2133b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2148263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
215806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
216a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
21762294041SBarry Smith 
2185a9b9e01SMark F. Adams       /* find factor */
219c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2205a9b9e01SMark F. Adams       else {
2215a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2225a9b9e01SMark F. Adams         jj = -1;
223c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
2242a82295dSMark Adams           if (!(size%kk)) { /* a candidate */
225c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2265a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2275a9b9e01SMark F. Adams             if (fact > best_fact) {
2285a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
2295a9b9e01SMark F. Adams             }
2305a9b9e01SMark F. Adams           }
2315a9b9e01SMark F. Adams         }
2325a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
233a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
2345a9b9e01SMark F. Adams       }
235c5df96a5SBarry Smith       new_size = size/rfactor;
2365a9b9e01SMark F. Adams 
237c5df96a5SBarry Smith       if (new_size==nactive) {
238a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2395a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
240302440fdSBarry Smith         ierr = PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr);
241fe682e87SBarry Smith #if defined PETSC_GAMG_USE_LOG
242fe682e87SBarry Smith         ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
243fe682e87SBarry Smith #endif
2445a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2455a9b9e01SMark F. Adams       }
2465a9b9e01SMark F. Adams 
247302440fdSBarry Smith       ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr);
248c5df96a5SBarry Smith       targetPE = rank/rfactor;
2493b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
250a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
251e33ef3b1SMark F. Adams 
25211e60469SMark F. Adams     /*
253a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
25411e60469SMark F. Adams      */
255a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2567700e67bSMark Adams     is_eq_num_prim = is_eq_num;
25711e60469SMark F. Adams     /*
258a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
25911e60469SMark F. Adams      */
260c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
261c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
262a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
2633ae0bb68SMark Adams     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */
264a2f3521dSMark F. Adams 
265a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2660cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2670cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
268b4fbaa2aSMark F. Adams #endif
269885364a3SMark 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 */
270885364a3SMark Adams     {
271885364a3SMark Adams     Vec            src_crd, dest_crd;
272885364a3SMark 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;
273885364a3SMark Adams     VecScatter     vecscat;
274885364a3SMark Adams     PetscScalar    *array;
275885364a3SMark Adams     IS isscat;
276a2f3521dSMark F. Adams 
277a2f3521dSMark F. Adams     /* move data (for primal equations only) */
27822063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
2793b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
2803ae0bb68SMark Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
281c0dedaeaSBarry Smith     ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
28211e60469SMark F. Adams     /*
2839d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
284c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
28511e60469SMark F. Adams      */
286854ce69bSBarry Smith     ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
287a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2883ae0bb68SMark Adams     for (ii=0,jj=0; ii<ncrs; ii++) {
289c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
290a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
29111e60469SMark F. Adams     }
292a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2933ae0bb68SMark Adams     ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
29492a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
29511e60469SMark F. Adams     /*
29611e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
29711e60469SMark F. Adams      */
2983ae0bb68SMark Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
2999d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3003ae0bb68SMark Adams       const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3013ae0bb68SMark Adams       for (ii=0; ii<ncrs; ii++) {
3029d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
303a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
304c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
305676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
306d3d6bff4SMark F. Adams         }
307038e3b61SMark F. Adams       }
308eb07cef2SMark F. Adams     }
309eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
310eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
31111e60469SMark F. Adams     /*
31211e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
31311e60469SMark F. Adams       to the correct processor
31411e60469SMark F. Adams     */
3150298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
31611e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
31711e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31811e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31911e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
32011e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
32111e60469SMark F. Adams     /*
32211e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
32311e60469SMark F. Adams     */
324c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
325578f55a3SPeter Brune     ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3262fa5cd67SKarl Rupp 
3273ae0bb68SMark Adams     pc_gamg->data_sz = node_data_sz*ncrs_new;
3283ae0bb68SMark Adams     strideNew        = ncrs_new*ndata_rows;
3292fa5cd67SKarl Rupp 
33011e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3319d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3323ae0bb68SMark Adams       for (ii=0; ii<ncrs_new; ii++) {
3339d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
334a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
335c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
336d3d6bff4SMark F. Adams         }
337038e3b61SMark F. Adams       }
338038e3b61SMark F. Adams     }
33911e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
34011e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
341885364a3SMark Adams     }
342a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3430cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3440cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
345ed3f9983SMark F. Adams #endif
346a2f3521dSMark F. Adams 
34711e60469SMark F. Adams     /*
34811e60469SMark F. Adams       Invert for MatGetSubMatrix
34911e60469SMark F. Adams     */
350a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
351a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
352c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
353a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
354a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
355a2f3521dSMark F. Adams     }
3563cb8563fSToby Isaac     if (Pcolumnperm) {
3573cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3583cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3593cb8563fSToby Isaac     }
360a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3610cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3620cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3630cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
364ed3f9983SMark F. Adams #endif
365a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
366a2f3521dSMark F. Adams     {
367a2f3521dSMark F. Adams       Mat mat;
368806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
369a2f3521dSMark F. Adams       *a_Amat_crs = mat;
370a2f3521dSMark F. Adams     }
371038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
372a2f3521dSMark F. Adams 
3730cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3740cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
375ed3f9983SMark F. Adams #endif
37611e60469SMark F. Adams     /* prolongator */
37711e60469SMark F. Adams     {
37811e60469SMark F. Adams       IS       findices;
379a2f3521dSMark F. Adams       PetscInt Istart,Iend;
380a2f3521dSMark F. Adams       Mat      Pnew;
38162294041SBarry Smith 
382a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
3830cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3840cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
385ed3f9983SMark F. Adams #endif
3863b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
387c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
388806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
38911e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
390c5bfad50SMark F. Adams 
3910cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3920cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
393ed3f9983SMark F. Adams #endif
3943530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
395a2f3521dSMark F. Adams 
396a2f3521dSMark F. Adams       /* output - repartitioned */
397a2f3521dSMark F. Adams       *a_P_inout = Pnew;
398e33ef3b1SMark F. Adams     }
399a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4005b89ad90SMark F. Adams 
401c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
402a2f3521dSMark F. Adams   }
4035b89ad90SMark F. Adams   PetscFunctionReturn(0);
4045b89ad90SMark F. Adams }
4055b89ad90SMark F. Adams 
4065b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4075b89ad90SMark F. Adams /*
4085b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4095b89ad90SMark F. Adams                     by setting data structures and options.
4105b89ad90SMark F. Adams 
4115b89ad90SMark F. Adams    Input Parameter:
4125b89ad90SMark F. Adams .  pc - the preconditioner context
4135b89ad90SMark F. Adams 
4145b89ad90SMark F. Adams */
4155b89ad90SMark F. Adams #undef __FUNCT__
4165b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4179d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4185b89ad90SMark F. Adams {
4195b89ad90SMark F. Adams   PetscErrorCode ierr;
4209d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4215b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4222adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
423a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
4243b4367a7SBarry Smith   MPI_Comm       comm;
425c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
426587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
427e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
428a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
429569f4572SMark Adams   MatInfo        info;
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);
501a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);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;
5380298fd71SBarry Smith       } else Parr[level1] = NULL;
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 */
553c8b0795cSMark F. Adams     if (!Parr[level1]) {
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
563a2f3521dSMark F. Adams 
5641c1aac46SBarry Smith     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs,&Parr[level1], &Aarr[level1], &nactivepe, NULL);CHKERRQ(ierr);
565a2f3521dSMark F. Adams 
5660cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5670cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
568b4fbaa2aSMark F. Adams #endif
569a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
570569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
571569f4572SMark Adams     nnztot += info.nz_used;
5721d5b2942SMark 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);
573569f4572SMark Adams 
5740cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
575b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
576b4fbaa2aSMark F. Adams #endif
577a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
5789ab59c8bSMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) {
5799ab59c8bSMark Adams       ierr =  PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr);
580a90e85d9SMark Adams       level++;
581a90e85d9SMark Adams       break;
582a90e85d9SMark Adams     }
583c8b0795cSMark F. Adams   } /* levels */
584c8b0795cSMark F. Adams   ierr                  = PetscFree(pc_gamg->data);CHKERRQ(ierr);
585c8b0795cSMark F. Adams 
586569f4572SMark Adams   ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr);
5879d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
5885b89ad90SMark F. Adams   fine_level       = level;
5890298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
5905b89ad90SMark F. Adams 
59162294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
592d5280255SMark F. Adams     /* set default smoothers & set operators */
59362294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
594ffc955d6SMark F. Adams       KSP smoother;
595ffc955d6SMark F. Adams       PC  subpc;
596a2f3521dSMark F. Adams 
5979d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
598f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
599ffc955d6SMark F. Adams 
600a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
601a2f3521dSMark F. Adams       /* set ops */
60223ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
603a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
604a2f3521dSMark F. Adams 
605a2f3521dSMark F. Adams       /* set defaults */
6066c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
607a2f3521dSMark F. Adams 
6080c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
6090c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
6102d3561bbSSatish Balay         PetscInt sz;
611*7a28f3e5SMark Adams         IS       *iss;
612a2f3521dSMark F. Adams 
6132d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
614*7a28f3e5SMark Adams         iss   = ASMLocalIDsArr[level];
6150c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
6160a3c815dSMark Adams         ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr);
6170c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
6187f66b68fSMark Adams         if (!sz) {
619ffc955d6SMark F. Adams           IS       is;
6200a3c815dSMark Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
621*7a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
622a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
623806fa848SBarry Smith         } else {
624a94c3b12SMark F. Adams           PetscInt kk;
625*7a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr);
626a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
627*7a28f3e5SMark Adams             ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr);
628a94c3b12SMark F. Adams           }
629*7a28f3e5SMark Adams           ierr = PetscFree(iss);CHKERRQ(ierr);
630ffc955d6SMark F. Adams         }
6310298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
632ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
633806fa848SBarry Smith       } else {
634890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
635ffc955d6SMark F. Adams       }
636d5280255SMark F. Adams     }
637d5280255SMark F. Adams     {
638d5280255SMark F. Adams       /* coarse grid */
639d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
640d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
641d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
64223ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
643d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
644d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
645d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
646d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
64771959b99SBarry Smith       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
64871959b99SBarry Smith       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
649d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
650d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
6519dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
6522fb0b348SMark F. Adams       ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
65308e36f19SMark Adams       ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
6545b42dca8SJed Brown       /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
6555b42dca8SJed Brown        * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
6565b42dca8SJed Brown        * view every subdomain as though they were different. */
6575b42dca8SJed Brown       ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
658d5280255SMark F. Adams     }
659d5280255SMark F. Adams 
660d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
661d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
662e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
663d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
6641c1aac46SBarry Smith     if (!mg->galerkin) SETERRQ(comm,PETSC_ERR_USER,"PCGAMG must use Galerkin for coarse operators.");
6651c1aac46SBarry Smith     if (mg->galerkin == 1) mg->galerkin = 2;
666d5280255SMark F. Adams 
667d5280255SMark F. Adams     /* clean up */
668d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
669587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
670587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
6715b89ad90SMark F. Adams     }
6720cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
673806fa848SBarry Smith   } else {
6745f8cf99dSMark F. Adams     KSP smoother;
675302440fdSBarry Smith     ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
6769d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
67723ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
6785f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
6799d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
6805f8cf99dSMark F. Adams   }
6815b89ad90SMark F. Adams   PetscFunctionReturn(0);
6825b89ad90SMark F. Adams }
6835b89ad90SMark F. Adams 
684eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
6855b89ad90SMark F. Adams /*
6865b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
6875b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
6885b89ad90SMark F. Adams 
6895b89ad90SMark F. Adams    Input Parameter:
6905b89ad90SMark F. Adams .  pc - the preconditioner context
6915b89ad90SMark F. Adams 
6925b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
6935b89ad90SMark F. Adams */
6945b89ad90SMark F. Adams #undef __FUNCT__
6955b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
6965b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
6975b89ad90SMark F. Adams {
6985b89ad90SMark F. Adams   PetscErrorCode ierr;
6995b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
7005b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
7015b89ad90SMark F. Adams 
7025b89ad90SMark F. Adams   PetscFunctionBegin;
7035b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7049b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
7059b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
7069b8ffb57SJed Brown   }
70714a9496bSBarry Smith   ierr = PetscRandomDestroy(&pc_gamg->random);CHKERRQ(ierr);
7081ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
7091ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
7105b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7115b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7125b89ad90SMark F. Adams   PetscFunctionReturn(0);
7135b89ad90SMark F. Adams }
7145b89ad90SMark F. Adams 
715676e1743SMark F. Adams #undef __FUNCT__
716676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
717676e1743SMark F. Adams /*@
718cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
719676e1743SMark F. Adams 
7201cc46a46SBarry Smith    Logically Collective on PC
721676e1743SMark F. Adams 
722676e1743SMark F. Adams    Input Parameters:
7231cc46a46SBarry Smith +  pc - the preconditioner context
7241cc46a46SBarry Smith -  n - the number of equations
725676e1743SMark F. Adams 
726676e1743SMark F. Adams 
727676e1743SMark F. Adams    Options Database Key:
7281cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
729676e1743SMark F. Adams 
730cab9ed1eSBarry 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
731cab9ed1eSBarry Smith           that has degrees of freedom
732cab9ed1eSBarry Smith 
733676e1743SMark F. Adams    Level: intermediate
734676e1743SMark F. Adams 
7351c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
736676e1743SMark F. Adams 
7371c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim()
738676e1743SMark F. Adams @*/
739676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
740676e1743SMark F. Adams {
741676e1743SMark F. Adams   PetscErrorCode ierr;
742676e1743SMark F. Adams 
743676e1743SMark F. Adams   PetscFunctionBegin;
744676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
745676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
746676e1743SMark F. Adams   PetscFunctionReturn(0);
747676e1743SMark F. Adams }
748676e1743SMark F. Adams 
749676e1743SMark F. Adams #undef __FUNCT__
750676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
7511e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
752676e1743SMark F. Adams {
753c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
754c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
755676e1743SMark F. Adams 
756676e1743SMark F. Adams   PetscFunctionBegin;
7579d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
758676e1743SMark F. Adams   PetscFunctionReturn(0);
759676e1743SMark F. Adams }
760676e1743SMark F. Adams 
761676e1743SMark F. Adams #undef __FUNCT__
762389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
763389730f3SMark F. Adams /*@
764cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
765389730f3SMark F. Adams 
766389730f3SMark F. Adams  Collective on PC
767389730f3SMark F. Adams 
768389730f3SMark F. Adams    Input Parameters:
7691cc46a46SBarry Smith +  pc - the preconditioner context
7701cc46a46SBarry Smith -  n - maximum number of equations to aim for
771389730f3SMark F. Adams 
772389730f3SMark F. Adams    Options Database Key:
7731cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
774389730f3SMark F. Adams 
775389730f3SMark F. Adams    Level: intermediate
776389730f3SMark F. Adams 
7771c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
778389730f3SMark F. Adams 
7791c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim()
780389730f3SMark F. Adams @*/
781389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
782389730f3SMark F. Adams {
783389730f3SMark F. Adams   PetscErrorCode ierr;
784389730f3SMark F. Adams 
785389730f3SMark F. Adams   PetscFunctionBegin;
786389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
787389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
788389730f3SMark F. Adams   PetscFunctionReturn(0);
789389730f3SMark F. Adams }
790389730f3SMark F. Adams 
791389730f3SMark F. Adams #undef __FUNCT__
792389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
7931e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
794389730f3SMark F. Adams {
795389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
796389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
797389730f3SMark F. Adams 
798389730f3SMark F. Adams   PetscFunctionBegin;
7999d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
800389730f3SMark F. Adams   PetscFunctionReturn(0);
801389730f3SMark F. Adams }
802389730f3SMark F. Adams 
803389730f3SMark F. Adams #undef __FUNCT__
804cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGSetRepartition"
805676e1743SMark F. Adams /*@
806cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
807676e1743SMark F. Adams 
808676e1743SMark F. Adams    Collective on PC
809676e1743SMark F. Adams 
810676e1743SMark F. Adams    Input Parameters:
8111cc46a46SBarry Smith +  pc - the preconditioner context
8121cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
813676e1743SMark F. Adams 
814676e1743SMark F. Adams    Options Database Key:
8151cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
816676e1743SMark F. Adams 
817cab9ed1eSBarry Smith    Notes: this will generally improve the loading balancing of the work on each level
818cab9ed1eSBarry Smith 
819676e1743SMark F. Adams    Level: intermediate
820676e1743SMark F. Adams 
8211c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
822676e1743SMark F. Adams 
823676e1743SMark F. Adams .seealso: ()
824676e1743SMark F. Adams @*/
825cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
826676e1743SMark F. Adams {
827676e1743SMark F. Adams   PetscErrorCode ierr;
828676e1743SMark F. Adams 
829676e1743SMark F. Adams   PetscFunctionBegin;
830676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
831cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
832676e1743SMark F. Adams   PetscFunctionReturn(0);
833676e1743SMark F. Adams }
834676e1743SMark F. Adams 
835676e1743SMark F. Adams #undef __FUNCT__
836cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGSetRepartition_GAMG"
837cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
838676e1743SMark F. Adams {
839c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
840c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
841676e1743SMark F. Adams 
842676e1743SMark F. Adams   PetscFunctionBegin;
8439d5b6da9SMark F. Adams   pc_gamg->repart = n;
844676e1743SMark F. Adams   PetscFunctionReturn(0);
845676e1743SMark F. Adams }
846676e1743SMark F. Adams 
847676e1743SMark F. Adams #undef __FUNCT__
8481cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation"
849dfd5c07aSMark F. Adams /*@
850cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
851dfd5c07aSMark F. Adams 
852dfd5c07aSMark F. Adams    Collective on PC
853dfd5c07aSMark F. Adams 
854dfd5c07aSMark F. Adams    Input Parameters:
8551cc46a46SBarry Smith +  pc - the preconditioner context
8561cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
857dfd5c07aSMark F. Adams 
858dfd5c07aSMark F. Adams    Options Database Key:
8591cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
860dfd5c07aSMark F. Adams 
861dfd5c07aSMark F. Adams    Level: intermediate
862dfd5c07aSMark F. Adams 
863cab9ed1eSBarry 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
864cab9ed1eSBarry Smith           rebuilding the preconditioner quicker.
865cab9ed1eSBarry Smith 
8661c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
867dfd5c07aSMark F. Adams 
868dfd5c07aSMark F. Adams .seealso: ()
869dfd5c07aSMark F. Adams @*/
8701cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
871dfd5c07aSMark F. Adams {
872dfd5c07aSMark F. Adams   PetscErrorCode ierr;
873dfd5c07aSMark F. Adams 
874dfd5c07aSMark F. Adams   PetscFunctionBegin;
875dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8761cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
877dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
878dfd5c07aSMark F. Adams }
879dfd5c07aSMark F. Adams 
880dfd5c07aSMark F. Adams #undef __FUNCT__
8811cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG"
8821cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
883dfd5c07aSMark F. Adams {
884dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
885dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
886dfd5c07aSMark F. Adams 
887dfd5c07aSMark F. Adams   PetscFunctionBegin;
888dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
889dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
890dfd5c07aSMark F. Adams }
891dfd5c07aSMark F. Adams 
892dfd5c07aSMark F. Adams #undef __FUNCT__
893cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGASMSetUseAggs"
894ffc955d6SMark F. Adams /*@
895cab9ed1eSBarry 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.
896ffc955d6SMark F. Adams 
897ffc955d6SMark F. Adams    Collective on PC
898ffc955d6SMark F. Adams 
899ffc955d6SMark F. Adams    Input Parameters:
900cab9ed1eSBarry Smith +  pc - the preconditioner context
901cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
902ffc955d6SMark F. Adams 
903ffc955d6SMark F. Adams    Options Database Key:
904cab9ed1eSBarry Smith .  -pc_gamg_asm_use_agg
905ffc955d6SMark F. Adams 
906ffc955d6SMark F. Adams    Level: intermediate
907ffc955d6SMark F. Adams 
9081c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
909ffc955d6SMark F. Adams 
910ffc955d6SMark F. Adams .seealso: ()
911ffc955d6SMark F. Adams @*/
912cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
913ffc955d6SMark F. Adams {
914ffc955d6SMark F. Adams   PetscErrorCode ierr;
915ffc955d6SMark F. Adams 
916ffc955d6SMark F. Adams   PetscFunctionBegin;
917ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
918cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
919ffc955d6SMark F. Adams   PetscFunctionReturn(0);
920ffc955d6SMark F. Adams }
921ffc955d6SMark F. Adams 
922ffc955d6SMark F. Adams #undef __FUNCT__
923cab9ed1eSBarry Smith #define __FUNCT__ "PCGAMGASMSetUseAggs_GAMG"
924cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
925ffc955d6SMark F. Adams {
926ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
927ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
928ffc955d6SMark F. Adams 
929ffc955d6SMark F. Adams   PetscFunctionBegin;
930cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
931ffc955d6SMark F. Adams   PetscFunctionReturn(0);
932ffc955d6SMark F. Adams }
933ffc955d6SMark F. Adams 
934ffc955d6SMark F. Adams #undef __FUNCT__
9354ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
9364ef23d27SMark F. Adams /*@
9371cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
9384ef23d27SMark F. Adams 
9394ef23d27SMark F. Adams    Not collective on PC
9404ef23d27SMark F. Adams 
9414ef23d27SMark F. Adams    Input Parameters:
9421cc46a46SBarry Smith +  pc - the preconditioner
9431cc46a46SBarry Smith -  n - the maximum number of levels to use
9444ef23d27SMark F. Adams 
9454ef23d27SMark F. Adams    Options Database Key:
9464ef23d27SMark F. Adams .  -pc_mg_levels
9474ef23d27SMark F. Adams 
9484ef23d27SMark F. Adams    Level: intermediate
9494ef23d27SMark F. Adams 
9501c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
9514ef23d27SMark F. Adams 
9524ef23d27SMark F. Adams .seealso: ()
9534ef23d27SMark F. Adams @*/
9544ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
9554ef23d27SMark F. Adams {
9564ef23d27SMark F. Adams   PetscErrorCode ierr;
9574ef23d27SMark F. Adams 
9584ef23d27SMark F. Adams   PetscFunctionBegin;
9594ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9604ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
9614ef23d27SMark F. Adams   PetscFunctionReturn(0);
9624ef23d27SMark F. Adams }
9634ef23d27SMark F. Adams 
9644ef23d27SMark F. Adams #undef __FUNCT__
9654ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
9661e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
9674ef23d27SMark F. Adams {
9684ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
9694ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
9704ef23d27SMark F. Adams 
9714ef23d27SMark F. Adams   PetscFunctionBegin;
9729d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
9734ef23d27SMark F. Adams   PetscFunctionReturn(0);
9744ef23d27SMark F. Adams }
9754ef23d27SMark F. Adams 
9764ef23d27SMark F. Adams #undef __FUNCT__
9773542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
9783542efc5SMark F. Adams /*@
9793542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
9803542efc5SMark F. Adams 
9813542efc5SMark F. Adams    Not collective on PC
9823542efc5SMark F. Adams 
9833542efc5SMark F. Adams    Input Parameters:
9841cc46a46SBarry Smith +  pc - the preconditioner context
985b001cb0fSBarry Smith -  threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
9863542efc5SMark F. Adams 
9873542efc5SMark F. Adams    Options Database Key:
9881cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
9893542efc5SMark F. Adams 
990cab9ed1eSBarry Smith    Notes: Before aggregating the graph GAMG will remove small values from the graph thus reducing the coupling in the graph and a different
991cab9ed1eSBarry Smith     (perhaps better) coarser set of points.
992cab9ed1eSBarry Smith 
9933542efc5SMark F. Adams    Level: intermediate
9943542efc5SMark F. Adams 
9951c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
9963542efc5SMark F. Adams 
9973542efc5SMark F. Adams .seealso: ()
9983542efc5SMark F. Adams @*/
9993542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
10003542efc5SMark F. Adams {
10013542efc5SMark F. Adams   PetscErrorCode ierr;
10023542efc5SMark F. Adams 
10033542efc5SMark F. Adams   PetscFunctionBegin;
10043542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10053542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
10063542efc5SMark F. Adams   PetscFunctionReturn(0);
10073542efc5SMark F. Adams }
10083542efc5SMark F. Adams 
10093542efc5SMark F. Adams #undef __FUNCT__
10103542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
10111e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
10123542efc5SMark F. Adams {
1013c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1014c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
10153542efc5SMark F. Adams 
10163542efc5SMark F. Adams   PetscFunctionBegin;
10179d5b6da9SMark F. Adams   pc_gamg->threshold = n;
10183542efc5SMark F. Adams   PetscFunctionReturn(0);
10193542efc5SMark F. Adams }
10203542efc5SMark F. Adams 
10213542efc5SMark F. Adams #undef __FUNCT__
10229d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1023676e1743SMark F. Adams /*@
1024c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1025676e1743SMark F. Adams 
1026676e1743SMark F. Adams    Collective on PC
1027676e1743SMark F. Adams 
1028676e1743SMark F. Adams    Input Parameters:
1029c60c7ad4SBarry Smith +  pc - the preconditioner context
1030c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1031676e1743SMark F. Adams 
1032676e1743SMark F. Adams    Options Database Key:
1033cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1034676e1743SMark F. Adams 
1035676e1743SMark F. Adams    Level: intermediate
1036676e1743SMark F. Adams 
10371c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1038676e1743SMark F. Adams 
1039cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1040676e1743SMark F. Adams @*/
104119fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1042676e1743SMark F. Adams {
1043676e1743SMark F. Adams   PetscErrorCode ierr;
1044676e1743SMark F. Adams 
1045676e1743SMark F. Adams   PetscFunctionBegin;
1046676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1047806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1048676e1743SMark F. Adams   PetscFunctionReturn(0);
1049676e1743SMark F. Adams }
1050676e1743SMark F. Adams 
1051676e1743SMark F. Adams #undef __FUNCT__
1052c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType"
1053c60c7ad4SBarry Smith /*@
1054c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1055c60c7ad4SBarry Smith 
1056c60c7ad4SBarry Smith    Collective on PC
1057c60c7ad4SBarry Smith 
1058c60c7ad4SBarry Smith    Input Parameter:
1059c60c7ad4SBarry Smith .  pc - the preconditioner context
1060c60c7ad4SBarry Smith 
1061c60c7ad4SBarry Smith    Output Parameter:
1062c60c7ad4SBarry Smith .  type - the type of algorithm used
1063c60c7ad4SBarry Smith 
1064c60c7ad4SBarry Smith    Level: intermediate
1065c60c7ad4SBarry Smith 
10661c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1067c60c7ad4SBarry Smith 
10681c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1069c60c7ad4SBarry Smith @*/
1070c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1071c60c7ad4SBarry Smith {
1072c60c7ad4SBarry Smith   PetscErrorCode ierr;
1073c60c7ad4SBarry Smith 
1074c60c7ad4SBarry Smith   PetscFunctionBegin;
1075c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1076c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1077c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1078c60c7ad4SBarry Smith }
1079c60c7ad4SBarry Smith 
1080c60c7ad4SBarry Smith #undef __FUNCT__
1081c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG"
1082c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1083c60c7ad4SBarry Smith {
1084c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1085c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1086c60c7ad4SBarry Smith 
1087c60c7ad4SBarry Smith   PetscFunctionBegin;
1088c60c7ad4SBarry Smith   *type = pc_gamg->type;
1089c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1090c60c7ad4SBarry Smith }
1091c60c7ad4SBarry Smith 
1092c60c7ad4SBarry Smith #undef __FUNCT__
10939d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
10941e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1095676e1743SMark F. Adams {
10969d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
10971ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
10981ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1099676e1743SMark F. Adams 
1100676e1743SMark F. Adams   PetscFunctionBegin;
1101c60c7ad4SBarry Smith   pc_gamg->type = type;
11021c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
11039d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
11041ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
11051ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
11061ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1107e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
11083ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
11093ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
11103ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
11113ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
11123ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
11133ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
11143ae0bb68SMark Adams     pc_gamg->data_sz = 0;
11151ab5ffc9SJed Brown   }
11161ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
11171ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
11189d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1119676e1743SMark F. Adams   PetscFunctionReturn(0);
1120676e1743SMark F. Adams }
1121676e1743SMark F. Adams 
11225b89ad90SMark F. Adams #undef __FUNCT__
11235adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG"
11245adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
11255adeb434SBarry Smith {
11265adeb434SBarry Smith   PetscErrorCode ierr;
11275adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
11285adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
11295adeb434SBarry Smith 
11305adeb434SBarry Smith   PetscFunctionBegin;
11315adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1132b001cb0fSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values from graph %g\n",(double)pc_gamg->threshold);CHKERRQ(ierr);
1133cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
1134cab9ed1eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr);
1135cab9ed1eSBarry Smith   }
11365adeb434SBarry Smith   if (pc_gamg->ops->view) {
11375adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
11385adeb434SBarry Smith   }
11395adeb434SBarry Smith   PetscFunctionReturn(0);
11405adeb434SBarry Smith }
11415adeb434SBarry Smith 
11425adeb434SBarry Smith #undef __FUNCT__
11435b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
11444416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
11455b89ad90SMark F. Adams {
1146676e1743SMark F. Adams   PetscErrorCode ierr;
1147676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1148676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1149676e1743SMark F. Adams   PetscBool      flag;
11503b4367a7SBarry Smith   MPI_Comm       comm;
115114a9496bSBarry Smith   char           prefix[256];
115214a9496bSBarry Smith   const char     *pcpre;
11535b89ad90SMark F. Adams 
11545b89ad90SMark F. Adams   PetscFunctionBegin;
11553b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1156e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1157676e1743SMark F. Adams   {
1158bd94a7aaSJed Brown     char tname[256];
11591a1c1e04SBarry Smith     ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1160bd94a7aaSJed Brown     if (flag) {
1161bd94a7aaSJed Brown       ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
11621ab5ffc9SJed Brown     }
1163cab9ed1eSBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
11641cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1165cab9ed1eSBarry 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);
116694ae4db5SBarry 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);
116794ae4db5SBarry 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);
116894ae4db5SBarry Smith     ierr = PetscOptionsReal("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&pc_gamg->threshold,&flag);CHKERRQ(ierr);
116994ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1170b7cbab4eSMark Adams 
1171b7cbab4eSMark Adams     /* set options for subtype */
1172e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1173676e1743SMark F. Adams   }
117414a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
117514a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
117614a9496bSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)pc_gamg->random,prefix);CHKERRQ(ierr);
117714a9496bSBarry Smith   ierr = PetscRandomSetFromOptions(pc_gamg->random);CHKERRQ(ierr);
1178676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
11795b89ad90SMark F. Adams   PetscFunctionReturn(0);
11805b89ad90SMark F. Adams }
11815b89ad90SMark F. Adams 
11825b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
11835b89ad90SMark F. Adams /*MC
11841cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
11855b89ad90SMark F. Adams 
1186280d9858SJed Brown    Options Database Keys:
1187cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1188cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1189cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1190cab9ed1eSBarry 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
1191cab9ed1eSBarry 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>
1192cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1193cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
1194cab9ed1eSBarry Smith -   -pc_gamg_threshold <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph thus reducing the coupling in the graph and a different
1195cab9ed1eSBarry Smith 
1196cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1197cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1198cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1199cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1200cab9ed1eSBarry Smith 
1201cab9ed1eSBarry Smith    Multigrid options(inherited):
12021cc46a46SBarry Smith +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1203280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1204280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1205cab9ed1eSBarry Smith .  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
1206cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
12075b89ad90SMark F. Adams 
12081cc46a46SBarry Smith 
12091cc46a46SBarry Smith   Notes: In order to obtain good performance for PCGAMG for vector valued problems you must
12101cc46a46SBarry Smith $       Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
12111cc46a46SBarry Smith $       Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
12121cc46a46SBarry Smith $       See the Users Manual Chapter 4 for more details
12131cc46a46SBarry Smith 
12145b89ad90SMark F. Adams   Level: intermediate
1215280d9858SJed Brown 
12161cc46a46SBarry Smith   Concepts: algebraic multigrid
12175b89ad90SMark F. Adams 
12181cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1219cab9ed1eSBarry Smith            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation()
12205b89ad90SMark F. Adams M*/
1221b2573a8aSBarry Smith 
12225b89ad90SMark F. Adams #undef __FUNCT__
12235b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
12248cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
12255b89ad90SMark F. Adams {
12265b89ad90SMark F. Adams   PetscErrorCode ierr;
12275b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
12285b89ad90SMark F. Adams   PC_MG          *mg;
12295b89ad90SMark F. Adams 
12305b89ad90SMark F. Adams   PetscFunctionBegin;
12311c1aac46SBarry Smith    /* register AMG type */
12321c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
12331c1aac46SBarry Smith 
12345b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
12351c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
12365b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
12375b89ad90SMark F. Adams 
12385b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1239b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
12405b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
12411c1aac46SBarry Smith   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally from PCMG by GAMG code */
12425b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
12435b89ad90SMark F. Adams 
1244b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
12451ab5ffc9SJed Brown 
12469d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
12479d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
12489d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
12499d5b6da9SMark F. Adams   pc_gamg->data    = 0;
12505b89ad90SMark F. Adams 
12519d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
12525b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
12535b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
12545b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
12555b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
12565adeb434SBarry Smith   mg->view                = PCView_GAMG;
12575b89ad90SMark F. Adams 
1258bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1259bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1260cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr);
12611cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1262cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr);
1263bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1264bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1265c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1266bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
12679d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1268d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
12690c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1270038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
127125a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1272d3042614SMark Adams   pc_gamg->threshold        = 0.;
12739d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
12749ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
1275c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
12769d5b6da9SMark F. Adams 
127714a9496bSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)pc),&pc_gamg->random);CHKERRQ(ierr);
127814a9496bSBarry Smith 
1279bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1280bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
12815b89ad90SMark F. Adams   PetscFunctionReturn(0);
12825b89ad90SMark F. Adams }
12833e3471ccSMark Adams 
12843e3471ccSMark Adams #undef __FUNCT__
12853e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage"
12863e3471ccSMark Adams /*@C
12873e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
12883e3471ccSMark Adams     from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
12893e3471ccSMark Adams     when using static libraries.
12903e3471ccSMark Adams 
12913e3471ccSMark Adams  Level: developer
12923e3471ccSMark Adams 
12933e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
12943e3471ccSMark Adams  .seealso: PetscInitialize()
12953e3471ccSMark Adams @*/
12963e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
12973e3471ccSMark Adams {
12983e3471ccSMark Adams   PetscErrorCode ierr;
12993e3471ccSMark Adams 
13003e3471ccSMark Adams   PetscFunctionBegin;
13013e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
13023e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
13033e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
13043e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
13058e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
13063e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1307c1c463dbSMark Adams 
1308c1c463dbSMark Adams   /* general events */
1309fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1310fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1311fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1312fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1313c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1314c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1315fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1316c1c463dbSMark Adams 
13175b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13185b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
13195b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
13205b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
13215b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
13225b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
13235b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
13245b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
13255b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1326bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
13275b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
13285b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
13295b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
13305b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
13315b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
13325b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
13335b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
13345b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
13355b89ad90SMark F. Adams 
13365b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
13375b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
13385b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
13395b89ad90SMark F. Adams   /* create timer stages */
13405b89ad90SMark F. Adams #if defined GAMG_STAGES
13415b89ad90SMark F. Adams   {
13425b89ad90SMark F. Adams     char     str[32];
13435b89ad90SMark F. Adams     PetscInt lidx;
13445b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
13455b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
13465b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
13475b89ad90SMark F. Adams       sprintf(str,"MG Level %d",lidx);
13485b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
13495b89ad90SMark F. Adams     }
13505b89ad90SMark F. Adams   }
13515b89ad90SMark F. Adams #endif
13525b89ad90SMark F. Adams #endif
13533e3471ccSMark Adams   PetscFunctionReturn(0);
13543e3471ccSMark Adams }
13553e3471ccSMark Adams 
13563e3471ccSMark Adams #undef __FUNCT__
13573e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage"
13583e3471ccSMark Adams /*@C
13591c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
13601c1aac46SBarry Smith     called from PetscFinalize() automatically.
13613e3471ccSMark Adams 
13623e3471ccSMark Adams  Level: developer
13633e3471ccSMark Adams 
13643e3471ccSMark Adams  .keywords: Petsc, destroy, package
13653e3471ccSMark Adams  .seealso: PetscFinalize()
13663e3471ccSMark Adams @*/
13673e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
13683e3471ccSMark Adams {
13693e3471ccSMark Adams   PetscErrorCode ierr;
13703e3471ccSMark Adams 
13713e3471ccSMark Adams   PetscFunctionBegin;
13723e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
13733e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
13743e3471ccSMark Adams   PetscFunctionReturn(0);
13753e3471ccSMark Adams }
1376a36cf38bSToby Isaac 
1377a36cf38bSToby Isaac #undef __FUNCT__
1378a36cf38bSToby Isaac #define __FUNCT__ "PCGAMGRegister"
1379a36cf38bSToby Isaac /*@C
1380a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1381a36cf38bSToby Isaac 
1382a36cf38bSToby Isaac  Input Parameters:
1383a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1384a36cf38bSToby Isaac  - create - function for creating the gamg context.
1385a36cf38bSToby Isaac 
1386a36cf38bSToby Isaac   Level: advanced
1387a36cf38bSToby Isaac 
13881c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1389a36cf38bSToby Isaac @*/
1390a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1391a36cf38bSToby Isaac {
1392a36cf38bSToby Isaac   PetscErrorCode ierr;
1393a36cf38bSToby Isaac 
1394a36cf38bSToby Isaac   PetscFunctionBegin;
1395a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1396a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1397a36cf38bSToby Isaac   PetscFunctionReturn(0);
1398a36cf38bSToby Isaac }
1399a36cf38bSToby Isaac 
1400