xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision d9558ea91a0ca1fa4633d87e21172ffbe9a39015)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4aaa7dc30SBarry Smith #include <petsc-private/matimpl.h>
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
6b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
75b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */
8f96513f1SMatthew G Knepley 
90cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
100cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
11b4fbaa2aSMark F. Adams #endif
120cbbd2e1SMark F. Adams 
130cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
14fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG;
15fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO;
160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
20fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG;
210cbbd2e1SMark F. Adams #endif
220cbbd2e1SMark F. Adams 
23b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30
24b4fbaa2aSMark F. Adams 
25b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
260cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
27b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
28b4fbaa2aSMark F. Adams #endif
29f96513f1SMatthew G Knepley 
30140e18c1SBarry Smith static PetscFunctionList GAMGList = 0;
313e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
329d5b6da9SMark F. Adams 
33d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
34d3d6bff4SMark F. Adams #undef __FUNCT__
35d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG"
36d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
37d3d6bff4SMark F. Adams {
38d3d6bff4SMark F. Adams   PetscErrorCode ierr;
39d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
40d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
41d3d6bff4SMark F. Adams 
42d3d6bff4SMark F. Adams   PetscFunctionBegin;
43a2f3521dSMark F. Adams   if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */
44ce94432eSBarry Smith     PetscPrintf(PetscObjectComm((PetscObject)pc),"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__);
459d5b6da9SMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
4658471d46SMark F. Adams   }
471c1aac46SBarry Smith   pc_gamg->data_sz = 0;
48878e152fSMark F. Adams   ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
49a2f3521dSMark F. Adams   PetscFunctionReturn(0);
50a2f3521dSMark F. Adams }
51a2f3521dSMark F. Adams 
525b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
535b89ad90SMark F. Adams /*
54c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
55a147abb0SMark F. Adams      of active processors.
565b89ad90SMark F. Adams 
575b89ad90SMark F. Adams    Input Parameter:
58a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
59a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
609d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
61c5bfad50SMark F. Adams    . cr_bs - coarse block size
623530afc2SMark F. Adams    In/Output Parameter:
63a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
64afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
6511e60469SMark F. Adams    Output Parameter:
663530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
675b89ad90SMark F. Adams */
685cb416c2SMark F. Adams 
695b89ad90SMark F. Adams #undef __FUNCT__
70c238b0ebSToby Isaac #define __FUNCT__ "PCGAMGCreateLevel_GAMG"
71b34066adSToby Isaac static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,
723cb8563fSToby Isaac                                   Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,
733cb8563fSToby Isaac                                   IS * Pcolumnperm)
745b89ad90SMark F. Adams {
75a2f3521dSMark F. Adams   PetscErrorCode  ierr;
769d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
77486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
78a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
793b4367a7SBarry Smith   MPI_Comm        comm;
80c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
813ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
825b89ad90SMark F. Adams 
835b89ad90SMark F. Adams   PetscFunctionBegin;
843b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
853b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
863b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
87c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
889d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
89038e3b61SMark F. Adams 
903ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
910298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
923ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
933ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
9473911c69SBarry Smith   } else {
953ae0bb68SMark Adams     PetscInt  bs;
963ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
973ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
983ae0bb68SMark Adams   }
99a2f3521dSMark F. Adams 
100c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
101a2f3521dSMark F. Adams   {
102472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
1030298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
104a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
105a90e85d9SMark Adams     if (new_size == 0) new_size = 1; /* not likely, posible? */
106c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
107a2f3521dSMark F. Adams   }
108f852f58cSMark F. Adams 
1093cb8563fSToby Isaac   if (Pcolumnperm) *Pcolumnperm = NULL;
1103cb8563fSToby Isaac 
111a90e85d9SMark Adams   if (!pc_gamg->repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
1122fa5cd67SKarl Rupp   else {
1133ae0bb68SMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
114885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
115e33ef3b1SMark F. Adams 
11671959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
11771959b99SBarry 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);
1180cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
1190cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
120b4fbaa2aSMark F. Adams #endif
121a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
122785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
123a90e85d9SMark Adams     if (pc_gamg->repart) {
124a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1255a9b9e01SMark F. Adams       Mat adj;
1265a9b9e01SMark F. Adams 
127bf4339c2SBarry Smith      ierr = PetscInfo3(pc,"\t repartition: size (active): %D --> %D, neq = %D\n",*a_nactive_proc,new_size,ncrs_eq);
1285a9b9e01SMark F. Adams 
129a2f3521dSMark F. Adams       /* get 'adj' */
130c5bfad50SMark F. Adams       if (cr_bs == 1) {
131038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
132806fa848SBarry Smith       } else {
133a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
134eb07cef2SMark F. Adams         Mat               tMat;
135a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
136b4fbaa2aSMark F. Adams         const PetscScalar *vals;
137b4fbaa2aSMark F. Adams         const PetscInt    *idx;
138a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
1399057884aSMark F. Adams         static PetscInt   llev = 0;
140*d9558ea9SBarry Smith         MatType           mtype;
141b4fbaa2aSMark F. Adams 
142578f55a3SPeter Brune         ierr = PetscMalloc1(ncrs, &d_nnz);CHKERRQ(ierr);
143578f55a3SPeter Brune         ierr = PetscMalloc1(ncrs, &o_nnz);CHKERRQ(ierr);
144a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
145a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
146c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
14758471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
148c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
149c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
15058471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1513ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1523ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
15358471d46SMark F. Adams         }
1546876a03eSMark F. Adams 
155*d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
1563b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1573ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
158*d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
159a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
160a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
16158471d46SMark F. Adams         ierr = PetscFree(d_nnz);CHKERRQ(ierr);
1625f8cf99dSMark F. Adams         ierr = PetscFree(o_nnz);CHKERRQ(ierr);
163eb07cef2SMark F. Adams 
164a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
165c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
16622063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
167eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
168c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
169eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
170eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
171eb07cef2SMark F. Adams           }
17222063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
173eb07cef2SMark F. Adams         }
174eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
175eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
176eb07cef2SMark F. Adams 
177b4fbaa2aSMark F. Adams         if (llev++ == -1) {
178b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
1798caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
1803b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
181b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
1823bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
183b4fbaa2aSMark F. Adams         }
184b4fbaa2aSMark F. Adams 
185eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
186eb07cef2SMark F. Adams 
187eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
188a2f3521dSMark F. Adams       } /* create 'adj' */
189f150b916SMark F. Adams 
190a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
1915a9b9e01SMark F. Adams         char            prefix[256];
1925a9b9e01SMark F. Adams         const char      *pcpre;
193b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
194b4fbaa2aSMark F. Adams         MatPartitioning mpart;
195a4b7d37bSMark F. Adams         IS              proc_is;
196a2f3521dSMark F. Adams         PetscInt        targetPE;
1972f03bc48SMark F. Adams 
1983b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
1995ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2009d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2018caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
20259a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
20311e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
204c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
205a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
20611e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2075a9b9e01SMark F. Adams 
2085ef31b24SMark F. Adams         /* collect IS info */
209785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
210a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
211a2f3521dSMark F. Adams         targetPE = 1; /* bring to "front" of machine */
212c5df96a5SBarry Smith         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
213a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
214c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
215a2f3521dSMark F. Adams             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
216eb07cef2SMark F. Adams           }
2175ef31b24SMark F. Adams         }
218a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
219a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2205ef31b24SMark F. Adams       }
2215ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2225a9b9e01SMark F. Adams 
2233b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2248263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
225806fa848SBarry Smith     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
226a2f3521dSMark F. Adams 
227a2f3521dSMark F. Adams       PetscInt rfactor,targetPE;
2285a9b9e01SMark F. Adams       /* find factor */
229c5df96a5SBarry Smith       if (new_size == 1) rfactor = size; /* easy */
2305a9b9e01SMark F. Adams       else {
2315a9b9e01SMark F. Adams         PetscReal best_fact = 0.;
2325a9b9e01SMark F. Adams         jj = -1;
233c5df96a5SBarry Smith         for (kk = 1 ; kk <= size ; kk++) {
234c5df96a5SBarry Smith           if (size%kk==0) { /* a candidate */
235c5df96a5SBarry Smith             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
2365a9b9e01SMark F. Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
2375a9b9e01SMark F. Adams             if (fact > best_fact) {
2385a9b9e01SMark F. Adams               best_fact = fact; jj = kk;
2395a9b9e01SMark F. Adams             }
2405a9b9e01SMark F. Adams           }
2415a9b9e01SMark F. Adams         }
2425a9b9e01SMark F. Adams         if (jj != -1) rfactor = jj;
243a2f3521dSMark F. Adams         else rfactor = 1; /* does this happen .. a prime */
2445a9b9e01SMark F. Adams       }
245c5df96a5SBarry Smith       new_size = size/rfactor;
2465a9b9e01SMark F. Adams 
247c5df96a5SBarry Smith       if (new_size==nactive) {
248a2f3521dSMark F. Adams         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
2495a9b9e01SMark F. Adams         ierr        = PetscFree(counts);CHKERRQ(ierr);
250bf4339c2SBarry Smith         ierr = PetscInfo2(pc,"\t aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);
2515a9b9e01SMark F. Adams         PetscFunctionReturn(0);
2525a9b9e01SMark F. Adams       }
2535a9b9e01SMark F. Adams 
254bf4339c2SBarry Smith       ierr = PetscInfo1(pc,"\t number of equations (loc) %D with simple aggregation\n",ncrs_eq);
255c5df96a5SBarry Smith       targetPE = rank/rfactor;
2563b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
257a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
258e33ef3b1SMark F. Adams 
25911e60469SMark F. Adams     /*
260a2f3521dSMark F. Adams      Create an index set from the is_eq_newproc index set to indicate the mapping TO
26111e60469SMark F. Adams      */
262a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2637700e67bSMark Adams     is_eq_num_prim = is_eq_num;
26411e60469SMark F. Adams     /*
265a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
26611e60469SMark F. Adams      */
267c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
268c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
269a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
2703ae0bb68SMark Adams     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */
271a2f3521dSMark F. Adams 
272a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
2730cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2740cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
275b4fbaa2aSMark F. Adams #endif
276885364a3SMark 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 */
277885364a3SMark Adams     {
278885364a3SMark Adams     Vec            src_crd, dest_crd;
279885364a3SMark 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;
280885364a3SMark Adams     VecScatter     vecscat;
281885364a3SMark Adams     PetscScalar    *array;
282885364a3SMark Adams     IS isscat;
283a2f3521dSMark F. Adams 
284a2f3521dSMark F. Adams     /* move data (for primal equations only) */
28522063be5SMark F. Adams     /* Create a vector to contain the newly ordered element information */
2863b4367a7SBarry Smith     ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
2873ae0bb68SMark Adams     ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
288c0dedaeaSBarry Smith     ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
28911e60469SMark F. Adams     /*
2909d5b6da9SMark F. Adams      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
291c5bfad50SMark F. Adams      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
29211e60469SMark F. Adams      */
293854ce69bSBarry Smith     ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
294a2f3521dSMark F. Adams     ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2953ae0bb68SMark Adams     for (ii=0,jj=0; ii<ncrs; ii++) {
296c5bfad50SMark F. Adams       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
297a2f3521dSMark F. Adams       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
29811e60469SMark F. Adams     }
299a2f3521dSMark F. Adams     ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3003ae0bb68SMark Adams     ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
30192a756f0SMark F. Adams     ierr = PetscFree(tidx);CHKERRQ(ierr);
30211e60469SMark F. Adams     /*
30311e60469SMark F. Adams      Create a vector to contain the original vertex information for each element
30411e60469SMark F. Adams      */
3053ae0bb68SMark Adams     ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
3069d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3073ae0bb68SMark Adams       const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3083ae0bb68SMark Adams       for (ii=0; ii<ncrs; ii++) {
3099d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
310a2f3521dSMark F. Adams           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
311c8b0795cSMark F. Adams           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
312676e1743SMark F. Adams           ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
313d3d6bff4SMark F. Adams         }
314038e3b61SMark F. Adams       }
315eb07cef2SMark F. Adams     }
316eb07cef2SMark F. Adams     ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
317eb07cef2SMark F. Adams     ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
31811e60469SMark F. Adams     /*
31911e60469SMark F. Adams       Scatter the element vertex information (still in the original vertex ordering)
32011e60469SMark F. Adams       to the correct processor
32111e60469SMark F. Adams     */
3220298fd71SBarry Smith     ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
32311e60469SMark F. Adams     ierr = ISDestroy(&isscat);CHKERRQ(ierr);
32411e60469SMark F. Adams     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32511e60469SMark F. Adams     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32611e60469SMark F. Adams     ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
32711e60469SMark F. Adams     ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
32811e60469SMark F. Adams     /*
32911e60469SMark F. Adams       Put the element vertex data into a new allocation of the gdata->ele
33011e60469SMark F. Adams     */
331c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
332578f55a3SPeter Brune     ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3332fa5cd67SKarl Rupp 
3343ae0bb68SMark Adams     pc_gamg->data_sz = node_data_sz*ncrs_new;
3353ae0bb68SMark Adams     strideNew        = ncrs_new*ndata_rows;
3362fa5cd67SKarl Rupp 
33711e60469SMark F. Adams     ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3389d5b6da9SMark F. Adams     for (jj=0; jj<ndata_cols; jj++) {
3393ae0bb68SMark Adams       for (ii=0; ii<ncrs_new; ii++) {
3409d5b6da9SMark F. Adams         for (kk=0; kk<ndata_rows; kk++) {
341a2f3521dSMark F. Adams           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
342c8b0795cSMark F. Adams           pc_gamg->data[ix] = PetscRealPart(array[jx]);
343d3d6bff4SMark F. Adams         }
344038e3b61SMark F. Adams       }
345038e3b61SMark F. Adams     }
34611e60469SMark F. Adams     ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
34711e60469SMark F. Adams     ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
348885364a3SMark Adams     }
349a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3500cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3510cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
352ed3f9983SMark F. Adams #endif
353a2f3521dSMark F. Adams 
35411e60469SMark F. Adams     /*
35511e60469SMark F. Adams       Invert for MatGetSubMatrix
35611e60469SMark F. Adams     */
357a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
358a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
359c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
360a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
361a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
362a2f3521dSMark F. Adams     }
3633cb8563fSToby Isaac     if (Pcolumnperm) {
3643cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3653cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3663cb8563fSToby Isaac     }
367a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3680cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3690cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3700cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
371ed3f9983SMark F. Adams #endif
372a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
373a2f3521dSMark F. Adams     {
374a2f3521dSMark F. Adams       Mat mat;
375806fa848SBarry Smith       ierr        = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
376a2f3521dSMark F. Adams       *a_Amat_crs = mat;
377c5bfad50SMark F. Adams 
378c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
379c5bfad50SMark F. Adams         PetscInt cbs, rbs;
380c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr);
381c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
382c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr);
383c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",rank,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr);
384c5bfad50SMark F. Adams       }
385a2f3521dSMark F. Adams     }
386038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
387a2f3521dSMark F. Adams 
3880cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3890cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
390ed3f9983SMark F. Adams #endif
39111e60469SMark F. Adams     /* prolongator */
39211e60469SMark F. Adams     {
39311e60469SMark F. Adams       IS       findices;
394a2f3521dSMark F. Adams       PetscInt Istart,Iend;
395a2f3521dSMark F. Adams       Mat      Pnew;
396a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
3970cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3980cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
399ed3f9983SMark F. Adams #endif
4003b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
401c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
402806fa848SBarry Smith       ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
40311e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
404c5bfad50SMark F. Adams 
405c5bfad50SMark F. Adams       if (!PETSC_TRUE) {
406c5bfad50SMark F. Adams         PetscInt cbs, rbs;
407c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr);
408c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
409c5bfad50SMark F. Adams         ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr);
410c5df96a5SBarry Smith         ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr);
411c5bfad50SMark F. Adams       }
4120cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4130cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
414ed3f9983SMark F. Adams #endif
4153530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
416a2f3521dSMark F. Adams 
417a2f3521dSMark F. Adams       /* output - repartitioned */
418a2f3521dSMark F. Adams       *a_P_inout = Pnew;
419e33ef3b1SMark F. Adams     }
420a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4215b89ad90SMark F. Adams 
422c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
423a2f3521dSMark F. Adams   }
4245a9b9e01SMark F. Adams 
425a2f3521dSMark F. Adams   /* outout matrix data */
426c8b0795cSMark F. Adams   if (!PETSC_TRUE) {
427c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
428c8b0795cSMark F. Adams     if (llev==0) {
429c8b0795cSMark F. Adams       sprintf(fname,"Cmat_%d.m",llev++);
4303b4367a7SBarry Smith       PetscViewerASCIIOpen(comm,fname,&viewer);
431c8b0795cSMark F. Adams       ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
432c8b0795cSMark F. Adams       ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr);
433c8b0795cSMark F. Adams       ierr = PetscViewerDestroy(&viewer);
434c8b0795cSMark F. Adams     }
435c8b0795cSMark F. Adams     sprintf(fname,"Cmat_%d.m",llev++);
4363b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
437c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
438c8b0795cSMark F. Adams     ierr = MatView(Cmat, viewer);CHKERRQ(ierr);
439c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
440c8b0795cSMark F. Adams   }
4415b89ad90SMark F. Adams   PetscFunctionReturn(0);
4425b89ad90SMark F. Adams }
4435b89ad90SMark F. Adams 
4445b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4455b89ad90SMark F. Adams /*
4465b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4475b89ad90SMark F. Adams                     by setting data structures and options.
4485b89ad90SMark F. Adams 
4495b89ad90SMark F. Adams    Input Parameter:
4505b89ad90SMark F. Adams .  pc - the preconditioner context
4515b89ad90SMark F. Adams 
4525b89ad90SMark F. Adams */
4535b89ad90SMark F. Adams #undef __FUNCT__
4545b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
4559d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4565b89ad90SMark F. Adams {
4575b89ad90SMark F. Adams   PetscErrorCode ierr;
4589d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4595b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4602adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
461a2f3521dSMark F. Adams   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
4623b4367a7SBarry Smith   MPI_Comm       comm;
463c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
464587fa25dSMark F. Adams   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
465c8b0795cSMark F. Adams   PetscReal      emaxs[GAMG_MAXLEVELS];
466e696ed0bSMark F. Adams   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
467a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
4685ef31b24SMark F. Adams 
4695b89ad90SMark F. Adams   PetscFunctionBegin;
4703b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4713b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4723b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
473dfd5c07aSMark F. Adams 
47484d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
4751c1aac46SBarry Smith     if ((PetscBool)(!pc_gamg->reuse_prol)) {
476878e152fSMark F. Adams       /* reset everything */
477878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
478878e152fSMark F. Adams       pc->setupcalled = 0;
479806fa848SBarry Smith     } else {
48084d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
48103a628feSMark F. Adams       /* just do Galerkin grids */
48258471d46SMark F. Adams       Mat          B,dA,dB;
48358471d46SMark F. Adams 
48471959b99SBarry Smith      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
4859d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
48658471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
48723ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
48858471d46SMark F. Adams         /* (re)set to get dirty flag */
48923ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
49058471d46SMark F. Adams 
4912fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
49203a628feSMark F. Adams           /* the first time through the matrix structure has changed from repartitioning */
4930a97e771SToby Isaac           if (pc_gamg->setup_count==2) {
49403a628feSMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
495084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
4962fa5cd67SKarl Rupp 
49703a628feSMark F. Adams             mglevels[level]->A = B;
498806fa848SBarry Smith           } else {
49923ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
50058471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
50103a628feSMark F. Adams           }
50223ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
50358471d46SMark F. Adams           dB   = B;
50458471d46SMark F. Adams         }
5055f8cf99dSMark F. Adams       }
506d5280255SMark F. Adams 
507d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
508d5280255SMark F. Adams 
50958471d46SMark F. Adams       PetscFunctionReturn(0);
510eb07cef2SMark F. Adams     }
511878e152fSMark F. Adams   }
512f6536408SMark F. Adams 
513878e152fSMark F. Adams   if (!pc_gamg->data) {
514878e152fSMark F. Adams     if (pc_gamg->orig_data) {
515878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5160298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5172fa5cd67SKarl Rupp 
518878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
519878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
520878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5212fa5cd67SKarl Rupp 
522785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
523878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
524806fa848SBarry Smith     } else {
5251ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5267700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5279d5b6da9SMark F. Adams     }
528878e152fSMark F. Adams   }
529878e152fSMark F. Adams 
530878e152fSMark F. Adams   /* cache original data for reuse */
5311c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
532785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
533878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
534878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
535878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
536878e152fSMark F. Adams   }
537038e3b61SMark F. Adams 
538302f38e8SMark F. Adams   /* get basic dims */
539302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
540a2f3521dSMark F. Adams   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
54184d3f75bSMark F. Adams 
542a2f3521dSMark F. Adams   /* Get A_i and R_i */
543c5df96a5SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
544a90e85d9SMark Adams        level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit);
5450205a208SMark F. Adams        level++) {
54657d29afaSToby Isaac     pc_gamg->firstCoarsen = (level ? PETSC_FALSE : PETSC_TRUE);
5475b89ad90SMark F. Adams     level1 = level + 1;
5480cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5490cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
550a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
551a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
552b4fbaa2aSMark F. Adams #endif
553a2f3521dSMark F. Adams #endif
554c8b0795cSMark F. Adams     { /* construct prolongator */
555725b86d8SJed Brown       Mat              Gmat;
5560cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5577700e67bSMark Adams       Mat              Prol11;
558c8b0795cSMark F. Adams 
5597700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5601ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5617700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
562c8b0795cSMark F. Adams 
563a2f3521dSMark F. Adams       /* could have failed to create new level */
564a2f3521dSMark F. Adams       if (Prol11) {
5659d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5660298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
567a2f3521dSMark F. Adams 
568fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
569c8b0795cSMark F. Adams           /* smooth */
570fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
571c8b0795cSMark F. Adams         }
572c8b0795cSMark F. Adams 
5737700e67bSMark Adams         Parr[level1] = Prol11;
5740298fd71SBarry Smith       } else Parr[level1] = NULL;
575ffc955d6SMark F. Adams 
576ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
5771b18a24aSMark Adams         PetscInt bs;
5781b18a24aSMark Adams         ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
579806fa848SBarry Smith         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
580ffc955d6SMark F. Adams       }
581ffc955d6SMark F. Adams 
582a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
58341b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
584a2f3521dSMark F. Adams     } /* construct prolongator scope */
5850cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5860cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
587c8b0795cSMark F. Adams #endif
5889d5b6da9SMark F. Adams     /* cache eigen estimate */
5899d5b6da9SMark F. Adams     if (pc_gamg->emax_id != -1) {
5909d5b6da9SMark F. Adams       PetscBool flag;
5917700e67bSMark Adams       ierr = PetscObjectComposedDataGetReal((PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
5929d5b6da9SMark F. Adams       if (!flag) emaxs[level] = -1.;
593806fa848SBarry Smith     } else emaxs[level] = -1.;
5942adcac29SMark F. Adams     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
595c8b0795cSMark F. Adams     if (!Parr[level1]) {
596bf4339c2SBarry Smith       ierr =  PetscInfo1(pc,"\t stop gridding, level %D\n",level);CHKERRQ(ierr);
597a90e85d9SMark Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
598a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
599a90e85d9SMark Adams #endif
600c8b0795cSMark F. Adams       break;
601c8b0795cSMark F. Adams     }
6020cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6030cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
604b4fbaa2aSMark F. Adams #endif
605a2f3521dSMark F. Adams 
6061c1aac46SBarry Smith     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs,&Parr[level1], &Aarr[level1], &nactivepe, NULL);CHKERRQ(ierr);
607a2f3521dSMark F. Adams 
6080cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6090cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
610b4fbaa2aSMark F. Adams #endif
611a2f3521dSMark F. Adams     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
6120cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
613b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
614b4fbaa2aSMark F. Adams #endif
615a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
616a90e85d9SMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M < 2) ) {
617bf4339c2SBarry Smith       ierr =  PetscInfo1(pc,"\t HARD stop of coarsening ?????????, level %D\n",level);CHKERRQ(ierr);
618a90e85d9SMark Adams       level++;
619a90e85d9SMark Adams       break;
620a90e85d9SMark Adams     }
621c8b0795cSMark F. Adams   } /* levels */
62257d29afaSToby Isaac   pc_gamg->firstCoarsen = PETSC_FALSE;
623c8b0795cSMark F. Adams   ierr                  = PetscFree(pc_gamg->data);CHKERRQ(ierr);
624c8b0795cSMark F. Adams 
625bf4339c2SBarry Smith   ierr = PetscInfo2(pc,"\t %D levels, grid complexity = %g\n",level+1,((double)nnztot)/nnz0);CHKERRQ(ierr);
6269d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6275b89ad90SMark F. Adams   fine_level       = level;
6280298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
6295b89ad90SMark F. Adams 
63084d3f75bSMark F. Adams   /* simple setup */
63184d3f75bSMark F. Adams   if (!PETSC_TRUE) {
63284d3f75bSMark F. Adams     PC_MG_Levels **mglevels = mg->levels;
63384d3f75bSMark F. Adams     for (lidx=0,level=pc_gamg->Nlevels-1;
63484d3f75bSMark F. Adams          lidx<fine_level;
63584d3f75bSMark F. Adams          lidx++, level--) {
63684d3f75bSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
63723ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level]);CHKERRQ(ierr);
63884d3f75bSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
63984d3f75bSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
64084d3f75bSMark F. Adams     }
64123ee1639SBarry Smith     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0]);CHKERRQ(ierr);
64284d3f75bSMark F. Adams 
64384d3f75bSMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
644806fa848SBarry Smith   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
645d5280255SMark F. Adams     /* set default smoothers & set operators */
6469d5b6da9SMark F. Adams     for (lidx = 1, level = pc_gamg->Nlevels-2;
647587fa25dSMark F. Adams          lidx <= fine_level;
648587fa25dSMark F. Adams          lidx++, level--) {
649ffc955d6SMark F. Adams       KSP smoother;
650ffc955d6SMark F. Adams       PC  subpc;
651a2f3521dSMark F. Adams 
6529d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
653f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
654ffc955d6SMark F. Adams 
655a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
656a2f3521dSMark F. Adams       /* set ops */
65723ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
658a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
659a2f3521dSMark F. Adams 
660a2f3521dSMark F. Adams       /* set defaults */
6616c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
662a2f3521dSMark F. Adams 
6631b18a24aSMark Adams       /* set blocks for GASM smoother that uses the 'aggregates' */
664ffc955d6SMark F. Adams       if (pc_gamg->use_aggs_in_gasm) {
6652d3561bbSSatish Balay         PetscInt sz;
6662d3561bbSSatish Balay         IS       *is;
667a2f3521dSMark F. Adams 
6682d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
6692d3561bbSSatish Balay         is   = ASMLocalIDsArr[level];
670ffc955d6SMark F. Adams         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
6711b18a24aSMark Adams         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
672ffc955d6SMark F. Adams         if (sz==0) {
673ffc955d6SMark F. Adams           IS       is;
674ffc955d6SMark F. Adams           PetscInt my0,kk;
675ffc955d6SMark F. Adams           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
676ffc955d6SMark F. Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
6770298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr);
678a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
679806fa848SBarry Smith         } else {
680a94c3b12SMark F. Adams           PetscInt kk;
6810298fd71SBarry Smith           ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr);
682a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
683a94c3b12SMark F. Adams             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
684a94c3b12SMark F. Adams           }
685ffc955d6SMark F. Adams           ierr = PetscFree(is);CHKERRQ(ierr);
686ffc955d6SMark F. Adams         }
6870298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
688ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
689ffc955d6SMark F. Adams         ierr                  = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
690806fa848SBarry Smith       } else {
691890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
692ffc955d6SMark F. Adams       }
693d5280255SMark F. Adams     }
694d5280255SMark F. Adams     {
695d5280255SMark F. Adams       /* coarse grid */
696d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
697d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
698d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
69923ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
700d5280255SMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
701d5280255SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
702d5280255SMark F. Adams       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
703d5280255SMark F. Adams       ierr = PCSetUp(subpc);CHKERRQ(ierr);
70471959b99SBarry Smith       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
70571959b99SBarry Smith       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
706d5280255SMark F. Adams       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
707d5280255SMark F. Adams       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
7089dbfc187SHong Zhang       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
7092fb0b348SMark F. Adams       ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
7105b42dca8SJed Brown       /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
7115b42dca8SJed Brown        * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
7125b42dca8SJed Brown        * view every subdomain as though they were different. */
7135b42dca8SJed Brown       ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
714d5280255SMark F. Adams     }
715d5280255SMark F. Adams 
716d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
717d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
718e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
719d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
7201c1aac46SBarry Smith     if (!mg->galerkin) SETERRQ(comm,PETSC_ERR_USER,"PCGAMG must use Galerkin for coarse operators.");
7211c1aac46SBarry Smith     if (mg->galerkin == 1) mg->galerkin = 2;
722d5280255SMark F. Adams 
723d5280255SMark F. Adams     /* create cheby smoothers */
7241c1aac46SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
725d5280255SMark F. Adams       KSP       smoother;
726890ffe84SMark Adams       PetscBool flag,flag2;
727d5280255SMark F. Adams       PC        subpc;
728d5280255SMark F. Adams 
729ffc955d6SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
730a2f3521dSMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
731a2f3521dSMark F. Adams 
732ffc955d6SMark F. Adams       /* do my own cheby */
7336c9de887SHong Zhang       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
7341c1aac46SBarry Smith       if (0 && flag) {
735ffc955d6SMark F. Adams         PetscReal emax, emin;
736251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
737890ffe84SMark Adams         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCSOR, &flag2);CHKERRQ(ierr);
738890ffe84SMark Adams         if ((flag||flag2) && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC but lets acccept SOR because it is close and safe (always lower) */
739890ffe84SMark Adams         else { /* eigen estimate 'emax' -- this is done in cheby */
740e696ed0bSMark F. Adams           KSP eksp;
741e696ed0bSMark F. Adams           Mat Lmat = Aarr[level];
742b2a4f308SMark F. Adams           Vec bb, xx;
743038e3b61SMark F. Adams 
7442a7a6963SBarry Smith           ierr = MatCreateVecs(Lmat, &bb, 0);CHKERRQ(ierr);
7452a7a6963SBarry Smith           ierr = MatCreateVecs(Lmat, &xx, 0);CHKERRQ(ierr);
746fc4362bfSMark F. Adams           {
747fc4362bfSMark F. Adams             PetscRandom rctx;
7483b4367a7SBarry Smith             ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
749fc4362bfSMark F. Adams             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
750fc4362bfSMark F. Adams             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
751fc4362bfSMark F. Adams             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
7525b89ad90SMark F. Adams           }
753a2f3521dSMark F. Adams 
754e696ed0bSMark F. Adams           /* zeroing out BC rows -- needed for crazy matrices */
755e696ed0bSMark F. Adams           {
756e696ed0bSMark F. Adams             PetscInt    Istart,Iend,ncols,jj,Ii;
757e696ed0bSMark F. Adams             PetscScalar zero = 0.0;
758e696ed0bSMark F. Adams             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
759e696ed0bSMark F. Adams             for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
760e696ed0bSMark F. Adams               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
761e696ed0bSMark F. Adams               if (ncols <= 1) {
762e696ed0bSMark F. Adams                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
763a94c3b12SMark F. Adams               }
764e696ed0bSMark F. Adams               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
765a94c3b12SMark F. Adams             }
766a94c3b12SMark F. Adams             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
767a94c3b12SMark F. Adams             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
768a94c3b12SMark F. Adams           }
769e696ed0bSMark F. Adams 
7703b4367a7SBarry Smith           ierr = KSPCreate(comm, &eksp);CHKERRQ(ierr);
771806fa848SBarry Smith           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
772fc4362bfSMark F. Adams           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
7731a166f3bSJed Brown           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
7741a166f3bSJed Brown           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
775f6536408SMark F. Adams           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
776f6536408SMark F. Adams 
777f6536408SMark F. Adams           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
77823ee1639SBarry Smith           ierr = KSPSetOperators(eksp, Lmat, Lmat);CHKERRQ(ierr);
779fc4362bfSMark F. Adams           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
7805a9b9e01SMark F. Adams 
781d3d0db20SJed Brown           /* set PC type to be same as smoother */
782ffc955d6SMark F. Adams           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
783b2a4f308SMark F. Adams 
7845a9b9e01SMark F. Adams           /* solve - keep stuff out of logging */
7855a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
7865a9b9e01SMark F. Adams           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
787fc4362bfSMark F. Adams           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
7885a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
7895a9b9e01SMark F. Adams           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
7905a9b9e01SMark F. Adams 
791fc4362bfSMark F. Adams           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
7925a9b9e01SMark F. Adams 
793fc4362bfSMark F. Adams           ierr = VecDestroy(&xx);CHKERRQ(ierr);
794fc4362bfSMark F. Adams           ierr = VecDestroy(&bb);CHKERRQ(ierr);
795fc4362bfSMark F. Adams           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
796f6536408SMark F. Adams 
797bf4339c2SBarry Smith #if defined(PETSC_USE_INFO)
798bf4339c2SBarry Smith           {
799bf4339c2SBarry Smith             PetscInt N1;
800bf4339c2SBarry Smith             ierr = MatGetSize(Aarr[level], &N1,NULL);CHKERRQ(ierr);
801bf4339c2SBarry Smith             ierr = PetscInfo4(pc,"\t\t\t PC setup max eigen=%e min=%e on level %D (N=%D)\n",(double)emax,(double)emin,lidx,N1);CHKERRQ(ierr);
802f6536408SMark F. Adams           }
803bf4339c2SBarry Smith #endif
804fc4362bfSMark F. Adams         }
805038e3b61SMark F. Adams         {
806c5bfad50SMark F. Adams           PetscInt N1, N0;
8070298fd71SBarry Smith           ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr);
8080298fd71SBarry Smith           ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr);
809f6536408SMark F. Adams           /* heuristic - is this crap? */
810b4ec6429SMark F. Adams           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
8115e7c91beSJed Brown           emin  = emax * pc_gamg->eigtarget[0];
8125e7c91beSJed Brown           emax *= pc_gamg->eigtarget[1];
813038e3b61SMark F. Adams         }
8146c9de887SHong Zhang         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
815ffc955d6SMark F. Adams       } /* setup checby flag */
816ffc955d6SMark F. Adams     } /* non-coarse levels */
817737a81a9SMark F. Adams 
818d5280255SMark F. Adams     /* clean up */
819d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
820587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
821587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
8225b89ad90SMark F. Adams     }
8230cbbd2e1SMark F. Adams 
8240cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
825806fa848SBarry Smith   } else {
8265f8cf99dSMark F. Adams     KSP smoother;
827bf4339c2SBarry Smith     ierr = PetscInfo(pc,"\tone level solver used (system is seen as DD). Using default solver.\n");
8289d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
82923ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
8305f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
8319d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8325f8cf99dSMark F. Adams   }
8335b89ad90SMark F. Adams   PetscFunctionReturn(0);
8345b89ad90SMark F. Adams }
8355b89ad90SMark F. Adams 
836eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8375b89ad90SMark F. Adams /*
8385b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8395b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8405b89ad90SMark F. Adams 
8415b89ad90SMark F. Adams    Input Parameter:
8425b89ad90SMark F. Adams .  pc - the preconditioner context
8435b89ad90SMark F. Adams 
8445b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8455b89ad90SMark F. Adams */
8465b89ad90SMark F. Adams #undef __FUNCT__
8475b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
8485b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8495b89ad90SMark F. Adams {
8505b89ad90SMark F. Adams   PetscErrorCode ierr;
8515b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
8525b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
8535b89ad90SMark F. Adams 
8545b89ad90SMark F. Adams   PetscFunctionBegin;
8555b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8569b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
8579b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
8589b8ffb57SJed Brown   }
8591ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
8601ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
8615b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8625b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8635b89ad90SMark F. Adams   PetscFunctionReturn(0);
8645b89ad90SMark F. Adams }
8655b89ad90SMark F. Adams 
866676e1743SMark F. Adams 
867676e1743SMark F. Adams #undef __FUNCT__
868676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim"
869676e1743SMark F. Adams /*@
8701cc46a46SBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via processor reduction.
871676e1743SMark F. Adams 
8721cc46a46SBarry Smith    Logically Collective on PC
873676e1743SMark F. Adams 
874676e1743SMark F. Adams    Input Parameters:
8751cc46a46SBarry Smith +  pc - the preconditioner context
8761cc46a46SBarry Smith -  n - the number of equations
877676e1743SMark F. Adams 
878676e1743SMark F. Adams 
879676e1743SMark F. Adams    Options Database Key:
8801cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
881676e1743SMark F. Adams 
882676e1743SMark F. Adams    Level: intermediate
883676e1743SMark F. Adams 
8841c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
885676e1743SMark F. Adams 
8861c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim()
887676e1743SMark F. Adams @*/
888676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
889676e1743SMark F. Adams {
890676e1743SMark F. Adams   PetscErrorCode ierr;
891676e1743SMark F. Adams 
892676e1743SMark F. Adams   PetscFunctionBegin;
893676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
894676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
895676e1743SMark F. Adams   PetscFunctionReturn(0);
896676e1743SMark F. Adams }
897676e1743SMark F. Adams 
898676e1743SMark F. Adams #undef __FUNCT__
899676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
9001e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
901676e1743SMark F. Adams {
902c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
903c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
904676e1743SMark F. Adams 
905676e1743SMark F. Adams   PetscFunctionBegin;
9069d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
907676e1743SMark F. Adams   PetscFunctionReturn(0);
908676e1743SMark F. Adams }
909676e1743SMark F. Adams 
910676e1743SMark F. Adams #undef __FUNCT__
911389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim"
912389730f3SMark F. Adams /*@
913389730f3SMark F. Adams    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
914389730f3SMark F. Adams 
915389730f3SMark F. Adams  Collective on PC
916389730f3SMark F. Adams 
917389730f3SMark F. Adams    Input Parameters:
9181cc46a46SBarry Smith +  pc - the preconditioner context
9191cc46a46SBarry Smith -  n - maximum number of equations to aim for
920389730f3SMark F. Adams 
921389730f3SMark F. Adams    Options Database Key:
9221cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
923389730f3SMark F. Adams 
924389730f3SMark F. Adams    Level: intermediate
925389730f3SMark F. Adams 
9261c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
927389730f3SMark F. Adams 
9281c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim()
929389730f3SMark F. Adams @*/
930389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
931389730f3SMark F. Adams {
932389730f3SMark F. Adams   PetscErrorCode ierr;
933389730f3SMark F. Adams 
934389730f3SMark F. Adams   PetscFunctionBegin;
935389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
936389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
937389730f3SMark F. Adams   PetscFunctionReturn(0);
938389730f3SMark F. Adams }
939389730f3SMark F. Adams 
940389730f3SMark F. Adams #undef __FUNCT__
941389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
9421e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
943389730f3SMark F. Adams {
944389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
945389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
946389730f3SMark F. Adams 
947389730f3SMark F. Adams   PetscFunctionBegin;
9489d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
949389730f3SMark F. Adams   PetscFunctionReturn(0);
950389730f3SMark F. Adams }
951389730f3SMark F. Adams 
952389730f3SMark F. Adams #undef __FUNCT__
9538263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning"
954676e1743SMark F. Adams /*@
9558263b398SMark F. Adams    PCGAMGSetRepartitioning - Repartition the coarse grids
956676e1743SMark F. Adams 
957676e1743SMark F. Adams    Collective on PC
958676e1743SMark F. Adams 
959676e1743SMark F. Adams    Input Parameters:
9601cc46a46SBarry Smith +  pc - the preconditioner context
9611cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
962676e1743SMark F. Adams 
963676e1743SMark F. Adams    Options Database Key:
9641cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
965676e1743SMark F. Adams 
966676e1743SMark F. Adams    Level: intermediate
967676e1743SMark F. Adams 
9681c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
969676e1743SMark F. Adams 
970676e1743SMark F. Adams .seealso: ()
971676e1743SMark F. Adams @*/
9728263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
973676e1743SMark F. Adams {
974676e1743SMark F. Adams   PetscErrorCode ierr;
975676e1743SMark F. Adams 
976676e1743SMark F. Adams   PetscFunctionBegin;
977676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9788263b398SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
979676e1743SMark F. Adams   PetscFunctionReturn(0);
980676e1743SMark F. Adams }
981676e1743SMark F. Adams 
982676e1743SMark F. Adams #undef __FUNCT__
9838263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
9841e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
985676e1743SMark F. Adams {
986c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
987c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
988676e1743SMark F. Adams 
989676e1743SMark F. Adams   PetscFunctionBegin;
9909d5b6da9SMark F. Adams   pc_gamg->repart = n;
991676e1743SMark F. Adams   PetscFunctionReturn(0);
992676e1743SMark F. Adams }
993676e1743SMark F. Adams 
994676e1743SMark F. Adams #undef __FUNCT__
9951cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation"
996dfd5c07aSMark F. Adams /*@
9971cc46a46SBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner
998dfd5c07aSMark F. Adams 
999dfd5c07aSMark F. Adams    Collective on PC
1000dfd5c07aSMark F. Adams 
1001dfd5c07aSMark F. Adams    Input Parameters:
10021cc46a46SBarry Smith +  pc - the preconditioner context
10031cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1004dfd5c07aSMark F. Adams 
1005dfd5c07aSMark F. Adams    Options Database Key:
10061cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
1007dfd5c07aSMark F. Adams 
1008dfd5c07aSMark F. Adams    Level: intermediate
1009dfd5c07aSMark F. Adams 
10101c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1011dfd5c07aSMark F. Adams 
1012dfd5c07aSMark F. Adams .seealso: ()
1013dfd5c07aSMark F. Adams @*/
10141cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1015dfd5c07aSMark F. Adams {
1016dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1017dfd5c07aSMark F. Adams 
1018dfd5c07aSMark F. Adams   PetscFunctionBegin;
1019dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10201cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1021dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1022dfd5c07aSMark F. Adams }
1023dfd5c07aSMark F. Adams 
1024dfd5c07aSMark F. Adams #undef __FUNCT__
10251cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG"
10261cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1027dfd5c07aSMark F. Adams {
1028dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1029dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1030dfd5c07aSMark F. Adams 
1031dfd5c07aSMark F. Adams   PetscFunctionBegin;
1032dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1033dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1034dfd5c07aSMark F. Adams }
1035dfd5c07aSMark F. Adams 
1036dfd5c07aSMark F. Adams #undef __FUNCT__
1037ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs"
1038ffc955d6SMark F. Adams /*@
1039ffc955d6SMark F. Adams    PCGAMGSetUseASMAggs -
1040ffc955d6SMark F. Adams 
1041ffc955d6SMark F. Adams    Collective on PC
1042ffc955d6SMark F. Adams 
1043ffc955d6SMark F. Adams    Input Parameters:
1044ffc955d6SMark F. Adams .  pc - the preconditioner context
1045ffc955d6SMark F. Adams 
1046ffc955d6SMark F. Adams 
1047ffc955d6SMark F. Adams    Options Database Key:
1048ffc955d6SMark F. Adams .  -pc_gamg_use_agg_gasm
1049ffc955d6SMark F. Adams 
1050ffc955d6SMark F. Adams    Level: intermediate
1051ffc955d6SMark F. Adams 
10521c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1053ffc955d6SMark F. Adams 
1054ffc955d6SMark F. Adams .seealso: ()
1055ffc955d6SMark F. Adams @*/
1056ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1057ffc955d6SMark F. Adams {
1058ffc955d6SMark F. Adams   PetscErrorCode ierr;
1059ffc955d6SMark F. Adams 
1060ffc955d6SMark F. Adams   PetscFunctionBegin;
1061ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1062ffc955d6SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1063ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1064ffc955d6SMark F. Adams }
1065ffc955d6SMark F. Adams 
1066ffc955d6SMark F. Adams #undef __FUNCT__
1067ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
10681e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1069ffc955d6SMark F. Adams {
1070ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1071ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1072ffc955d6SMark F. Adams 
1073ffc955d6SMark F. Adams   PetscFunctionBegin;
1074ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = n;
1075ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1076ffc955d6SMark F. Adams }
1077ffc955d6SMark F. Adams 
1078ffc955d6SMark F. Adams #undef __FUNCT__
10794ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels"
10804ef23d27SMark F. Adams /*@
10811cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
10824ef23d27SMark F. Adams 
10834ef23d27SMark F. Adams    Not collective on PC
10844ef23d27SMark F. Adams 
10854ef23d27SMark F. Adams    Input Parameters:
10861cc46a46SBarry Smith +  pc - the preconditioner
10871cc46a46SBarry Smith -  n - the maximum number of levels to use
10884ef23d27SMark F. Adams 
10894ef23d27SMark F. Adams    Options Database Key:
10904ef23d27SMark F. Adams .  -pc_mg_levels
10914ef23d27SMark F. Adams 
10924ef23d27SMark F. Adams    Level: intermediate
10934ef23d27SMark F. Adams 
10941c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
10954ef23d27SMark F. Adams 
10964ef23d27SMark F. Adams .seealso: ()
10974ef23d27SMark F. Adams @*/
10984ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
10994ef23d27SMark F. Adams {
11004ef23d27SMark F. Adams   PetscErrorCode ierr;
11014ef23d27SMark F. Adams 
11024ef23d27SMark F. Adams   PetscFunctionBegin;
11034ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11044ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
11054ef23d27SMark F. Adams   PetscFunctionReturn(0);
11064ef23d27SMark F. Adams }
11074ef23d27SMark F. Adams 
11084ef23d27SMark F. Adams #undef __FUNCT__
11094ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
11101e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
11114ef23d27SMark F. Adams {
11124ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
11134ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
11144ef23d27SMark F. Adams 
11154ef23d27SMark F. Adams   PetscFunctionBegin;
11169d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
11174ef23d27SMark F. Adams   PetscFunctionReturn(0);
11184ef23d27SMark F. Adams }
11194ef23d27SMark F. Adams 
11204ef23d27SMark F. Adams #undef __FUNCT__
11213542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold"
11223542efc5SMark F. Adams /*@
11233542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
11243542efc5SMark F. Adams 
11253542efc5SMark F. Adams    Not collective on PC
11263542efc5SMark F. Adams 
11273542efc5SMark F. Adams    Input Parameters:
11281cc46a46SBarry Smith +  pc - the preconditioner context
1129b001cb0fSBarry Smith -  threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
11303542efc5SMark F. Adams 
11313542efc5SMark F. Adams    Options Database Key:
11321cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
11333542efc5SMark F. Adams 
11343542efc5SMark F. Adams    Level: intermediate
11353542efc5SMark F. Adams 
11361c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
11373542efc5SMark F. Adams 
11383542efc5SMark F. Adams .seealso: ()
11393542efc5SMark F. Adams @*/
11403542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
11413542efc5SMark F. Adams {
11423542efc5SMark F. Adams   PetscErrorCode ierr;
11433542efc5SMark F. Adams 
11443542efc5SMark F. Adams   PetscFunctionBegin;
11453542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11463542efc5SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
11473542efc5SMark F. Adams   PetscFunctionReturn(0);
11483542efc5SMark F. Adams }
11493542efc5SMark F. Adams 
11503542efc5SMark F. Adams #undef __FUNCT__
11513542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
11521e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
11533542efc5SMark F. Adams {
1154c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1155c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
11563542efc5SMark F. Adams 
11573542efc5SMark F. Adams   PetscFunctionBegin;
11589d5b6da9SMark F. Adams   pc_gamg->threshold = n;
11593542efc5SMark F. Adams   PetscFunctionReturn(0);
11603542efc5SMark F. Adams }
11613542efc5SMark F. Adams 
11623542efc5SMark F. Adams #undef __FUNCT__
11639d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType"
1164676e1743SMark F. Adams /*@
1165c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1166676e1743SMark F. Adams 
1167676e1743SMark F. Adams    Collective on PC
1168676e1743SMark F. Adams 
1169676e1743SMark F. Adams    Input Parameters:
1170c60c7ad4SBarry Smith +  pc - the preconditioner context
1171c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1172676e1743SMark F. Adams 
1173676e1743SMark F. Adams    Options Database Key:
1174c60c7ad4SBarry Smith .  -pc_gamg_type <agg,geo,classical>
1175676e1743SMark F. Adams 
1176676e1743SMark F. Adams    Level: intermediate
1177676e1743SMark F. Adams 
11781c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1179676e1743SMark F. Adams 
1180c60c7ad4SBarry Smith .seealso: PCGAMGGetType(), PCGAMG
1181676e1743SMark F. Adams @*/
118219fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1183676e1743SMark F. Adams {
1184676e1743SMark F. Adams   PetscErrorCode ierr;
1185676e1743SMark F. Adams 
1186676e1743SMark F. Adams   PetscFunctionBegin;
1187676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1188806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1189676e1743SMark F. Adams   PetscFunctionReturn(0);
1190676e1743SMark F. Adams }
1191676e1743SMark F. Adams 
1192676e1743SMark F. Adams #undef __FUNCT__
1193c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType"
1194c60c7ad4SBarry Smith /*@
1195c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1196c60c7ad4SBarry Smith 
1197c60c7ad4SBarry Smith    Collective on PC
1198c60c7ad4SBarry Smith 
1199c60c7ad4SBarry Smith    Input Parameter:
1200c60c7ad4SBarry Smith .  pc - the preconditioner context
1201c60c7ad4SBarry Smith 
1202c60c7ad4SBarry Smith    Output Parameter:
1203c60c7ad4SBarry Smith .  type - the type of algorithm used
1204c60c7ad4SBarry Smith 
1205c60c7ad4SBarry Smith    Level: intermediate
1206c60c7ad4SBarry Smith 
12071c1aac46SBarry Smith    Concepts: Unstructured multigrid preconditioner
1208c60c7ad4SBarry Smith 
12091c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1210c60c7ad4SBarry Smith @*/
1211c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1212c60c7ad4SBarry Smith {
1213c60c7ad4SBarry Smith   PetscErrorCode ierr;
1214c60c7ad4SBarry Smith 
1215c60c7ad4SBarry Smith   PetscFunctionBegin;
1216c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1217c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1218c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1219c60c7ad4SBarry Smith }
1220c60c7ad4SBarry Smith 
1221c60c7ad4SBarry Smith #undef __FUNCT__
1222c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG"
1223c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1224c60c7ad4SBarry Smith {
1225c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1226c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1227c60c7ad4SBarry Smith 
1228c60c7ad4SBarry Smith   PetscFunctionBegin;
1229c60c7ad4SBarry Smith   *type = pc_gamg->type;
1230c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1231c60c7ad4SBarry Smith }
1232c60c7ad4SBarry Smith 
1233c60c7ad4SBarry Smith #undef __FUNCT__
12349d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG"
12351e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1236676e1743SMark F. Adams {
12379d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
12381ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
12391ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1240676e1743SMark F. Adams 
1241676e1743SMark F. Adams   PetscFunctionBegin;
1242c60c7ad4SBarry Smith   pc_gamg->type = type;
12431c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
12449d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
12451ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
12461ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
12471ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1248e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
12493ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
12503ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
12513ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
12523ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
12533ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
12543ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
12553ae0bb68SMark Adams     pc_gamg->data_sz = 0;
12561ab5ffc9SJed Brown   }
12571ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
12581ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
12599d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1260676e1743SMark F. Adams   PetscFunctionReturn(0);
1261676e1743SMark F. Adams }
1262676e1743SMark F. Adams 
12635b89ad90SMark F. Adams #undef __FUNCT__
12645adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG"
12655adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
12665adeb434SBarry Smith {
12675adeb434SBarry Smith   PetscErrorCode ierr;
12685adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
12695adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
12705adeb434SBarry Smith 
12715adeb434SBarry Smith   PetscFunctionBegin;
12725adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1273b001cb0fSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values from graph %g\n",(double)pc_gamg->threshold);CHKERRQ(ierr);
12745adeb434SBarry Smith   if (pc_gamg->ops->view) {
12755adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
12765adeb434SBarry Smith   }
12775adeb434SBarry Smith   PetscFunctionReturn(0);
12785adeb434SBarry Smith }
12795adeb434SBarry Smith 
12805adeb434SBarry Smith #undef __FUNCT__
12815b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
12828c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptions *PetscOptionsObject,PC pc)
12835b89ad90SMark F. Adams {
1284676e1743SMark F. Adams   PetscErrorCode ierr;
1285676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1286676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1287676e1743SMark F. Adams   PetscBool      flag;
12885e7c91beSJed Brown   PetscInt       two   = 2;
12893b4367a7SBarry Smith   MPI_Comm       comm;
12905b89ad90SMark F. Adams 
12915b89ad90SMark F. Adams   PetscFunctionBegin;
12923b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1293e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1294676e1743SMark F. Adams   {
1295bd94a7aaSJed Brown     char tname[256];
12961a1c1e04SBarry Smith     ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1297bd94a7aaSJed Brown     if (flag) {
1298bd94a7aaSJed Brown       ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
12991ab5ffc9SJed Brown     }
130094ae4db5SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
13011cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
130294ae4db5SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm","Use aggregation agragates for GASM smoother","PCGAMGUseASMAggs",pc_gamg->use_aggs_in_gasm,&pc_gamg->use_aggs_in_gasm,NULL);CHKERRQ(ierr);
130394ae4db5SBarry 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);
130494ae4db5SBarry 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);
130594ae4db5SBarry 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);
13060298fd71SBarry Smith     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr);
130794ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1308b7cbab4eSMark Adams 
1309b7cbab4eSMark Adams     /* set options for subtype */
1310e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1311676e1743SMark F. Adams   }
1312676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
13135b89ad90SMark F. Adams   PetscFunctionReturn(0);
13145b89ad90SMark F. Adams }
13155b89ad90SMark F. Adams 
13165b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
13175b89ad90SMark F. Adams /*MC
13181cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
13195b89ad90SMark F. Adams 
1320280d9858SJed Brown    Options Database Keys:
13215b89ad90SMark F. Adams    Multigrid options(inherited)
13221cc46a46SBarry Smith +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1323280d9858SJed Brown .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1324280d9858SJed Brown .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
13258c1c2452SJed Brown -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
13265b89ad90SMark F. Adams 
13271cc46a46SBarry Smith 
13281cc46a46SBarry Smith   Notes: In order to obtain good performance for PCGAMG for vector valued problems you must
13291cc46a46SBarry Smith $       Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
13301cc46a46SBarry Smith $       Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
13311cc46a46SBarry Smith $       See the Users Manual Chapter 4 for more details
13321cc46a46SBarry Smith 
13335b89ad90SMark F. Adams   Level: intermediate
1334280d9858SJed Brown 
13351cc46a46SBarry Smith   Concepts: algebraic multigrid
13365b89ad90SMark F. Adams 
13371cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
13381cc46a46SBarry Smith            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartitioning(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseASMAggs(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType()
13395b89ad90SMark F. Adams M*/
1340b2573a8aSBarry Smith 
13415b89ad90SMark F. Adams #undef __FUNCT__
13425b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
13438cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
13445b89ad90SMark F. Adams {
13455b89ad90SMark F. Adams   PetscErrorCode ierr;
13465b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
13475b89ad90SMark F. Adams   PC_MG          *mg;
13485b89ad90SMark F. Adams 
13495b89ad90SMark F. Adams   PetscFunctionBegin;
13501c1aac46SBarry Smith   /* register AMG type */
13511c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
13521c1aac46SBarry Smith 
13535b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
13541c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
13555b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
13565b89ad90SMark F. Adams 
13575b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1358b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
13595b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
13601c1aac46SBarry Smith   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally from PCMG by GAMG code */
13615b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
13625b89ad90SMark F. Adams 
1363b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
13641ab5ffc9SJed Brown 
13659d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
13669d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
13679d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
13689d5b6da9SMark F. Adams   pc_gamg->data    = 0;
13695b89ad90SMark F. Adams 
13709d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
13715b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
13725b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
13735b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
13745b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
13755adeb434SBarry Smith   mg->view                = PCView_GAMG;
13765b89ad90SMark F. Adams 
1377bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1378bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1379bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
13801cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1381bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1382bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1383bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1384c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1385bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
13869d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1387d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
1388ffc955d6SMark F. Adams   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1389038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
139025a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1391d3042614SMark Adams   pc_gamg->threshold        = 0.;
13929d5b6da9SMark F. Adams   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
13939d5b6da9SMark F. Adams   pc_gamg->emax_id          = -1;
139457d29afaSToby Isaac   pc_gamg->firstCoarsen     = PETSC_FALSE;
13955e7c91beSJed Brown   pc_gamg->eigtarget[0]     = 0.05;
13965e7c91beSJed Brown   pc_gamg->eigtarget[1]     = 1.05;
1397c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
13989d5b6da9SMark F. Adams 
1399bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1400bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
14015b89ad90SMark F. Adams   PetscFunctionReturn(0);
14025b89ad90SMark F. Adams }
14033e3471ccSMark Adams 
14043e3471ccSMark Adams #undef __FUNCT__
14053e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage"
14063e3471ccSMark Adams /*@C
14073e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
14083e3471ccSMark Adams     from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
14093e3471ccSMark Adams     when using static libraries.
14103e3471ccSMark Adams 
14113e3471ccSMark Adams  Level: developer
14123e3471ccSMark Adams 
14133e3471ccSMark Adams  .keywords: PC, PCGAMG, initialize, package
14143e3471ccSMark Adams  .seealso: PetscInitialize()
14153e3471ccSMark Adams @*/
14163e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
14173e3471ccSMark Adams {
14183e3471ccSMark Adams   PetscErrorCode ierr;
14193e3471ccSMark Adams 
14203e3471ccSMark Adams   PetscFunctionBegin;
14213e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
14223e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
14233e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
14243e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
14258e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
14263e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1427c1c463dbSMark Adams 
1428c1c463dbSMark Adams   /* general events */
1429fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1430fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1431fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1432fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1433c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1434c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1435fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1436c1c463dbSMark Adams 
14371c1aac46SBarry Smith #if defined PETSC_GAMG_USE_LOG
14381c1aac46SBarry Smith   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
14391c1aac46SBarry Smith   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
14401c1aac46SBarry Smith   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
14411c1aac46SBarry Smith   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
14421c1aac46SBarry Smith   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
14431c1aac46SBarry Smith   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
14441c1aac46SBarry Smith   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
14451c1aac46SBarry Smith   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
14461c1aac46SBarry Smith   ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
14471c1aac46SBarry Smith   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
14481c1aac46SBarry Smith   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
14491c1aac46SBarry Smith   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
14501c1aac46SBarry Smith   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
14511c1aac46SBarry Smith   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
14521c1aac46SBarry Smith   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
14531c1aac46SBarry Smith   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
14541c1aac46SBarry Smith   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
14551c1aac46SBarry Smith 
14561c1aac46SBarry Smith   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
14571c1aac46SBarry Smith   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
14581c1aac46SBarry Smith   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
14591c1aac46SBarry Smith   /* create timer stages */
14601c1aac46SBarry Smith #if defined GAMG_STAGES
14611c1aac46SBarry Smith   {
14621c1aac46SBarry Smith     char     str[32];
14631c1aac46SBarry Smith     PetscInt lidx;
14641c1aac46SBarry Smith     sprintf(str,"MG Level %d (finest)",0);
14651c1aac46SBarry Smith     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
14661c1aac46SBarry Smith     for (lidx=1; lidx<9; lidx++) {
14671c1aac46SBarry Smith       sprintf(str,"MG Level %d",lidx);
14681c1aac46SBarry Smith       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
14691c1aac46SBarry Smith     }
14701c1aac46SBarry Smith   }
14711c1aac46SBarry Smith #endif
14721c1aac46SBarry Smith #endif
14733e3471ccSMark Adams   PetscFunctionReturn(0);
14743e3471ccSMark Adams }
14753e3471ccSMark Adams 
14763e3471ccSMark Adams #undef __FUNCT__
14773e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage"
14783e3471ccSMark Adams /*@C
14791c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
14801c1aac46SBarry Smith     called from PetscFinalize() automatically.
14813e3471ccSMark Adams 
14823e3471ccSMark Adams  Level: developer
14833e3471ccSMark Adams 
14843e3471ccSMark Adams  .keywords: Petsc, destroy, package
14853e3471ccSMark Adams  .seealso: PetscFinalize()
14863e3471ccSMark Adams @*/
14873e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
14883e3471ccSMark Adams {
14893e3471ccSMark Adams   PetscErrorCode ierr;
14903e3471ccSMark Adams 
14913e3471ccSMark Adams   PetscFunctionBegin;
14923e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
14933e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
14943e3471ccSMark Adams   PetscFunctionReturn(0);
14953e3471ccSMark Adams }
1496a36cf38bSToby Isaac 
1497a36cf38bSToby Isaac #undef __FUNCT__
1498a36cf38bSToby Isaac #define __FUNCT__ "PCGAMGRegister"
1499a36cf38bSToby Isaac /*@C
1500a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1501a36cf38bSToby Isaac 
1502a36cf38bSToby Isaac  Input Parameters:
1503a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1504a36cf38bSToby Isaac  - create - function for creating the gamg context.
1505a36cf38bSToby Isaac 
1506a36cf38bSToby Isaac   Level: advanced
1507a36cf38bSToby Isaac 
15081c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1509a36cf38bSToby Isaac @*/
1510a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1511a36cf38bSToby Isaac {
1512a36cf38bSToby Isaac   PetscErrorCode ierr;
1513a36cf38bSToby Isaac 
1514a36cf38bSToby Isaac   PetscFunctionBegin;
1515a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1516a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1517a36cf38bSToby Isaac   PetscFunctionReturn(0);
1518a36cf38bSToby Isaac }
1519a36cf38bSToby Isaac 
1520