xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 5d5a09d9ef0b61c61936e6f45bee36e62c8236a6)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4af0996ceSBarry Smith #include <petsc/private/matimpl.h>
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
618c3aa7eSMark #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h>    /*I "petscksp.h" I*/
7f96513f1SMatthew G Knepley 
8c9567895SMark #if defined(PETSC_HAVE_CUDA)
9c9567895SMark   #include <cuda_runtime.h>
10c9567895SMark #endif
11c9567895SMark 
12c9567895SMark #if defined(PETSC_HAVE_HIP)
13c9567895SMark   #include <hip/hip_runtime.h>
14c9567895SMark #endif
15c9567895SMark 
16849bee69SMark Adams PetscLogEvent petsc_gamg_setup_events[GAMG_NUM_SET];
174555aa8cSStefano Zampini PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3];
180cbbd2e1SMark F. Adams 
19849bee69SMark Adams // #define GAMG_STAGES
204555aa8cSStefano Zampini #if defined(GAMG_STAGES)
2118c3aa7eSMark static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
22b4fbaa2aSMark F. Adams #endif
23f96513f1SMatthew G Knepley 
240a545947SLisandro Dalcin static PetscFunctionList GAMGList = NULL;
253e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
269d5b6da9SMark F. Adams 
27d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
28d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
29d3d6bff4SMark F. Adams {
30d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
31d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
32d3d6bff4SMark F. Adams 
33d3d6bff4SMark F. Adams   PetscFunctionBegin;
349566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->data));
351c1aac46SBarry Smith   pc_gamg->data_sz = 0;
369566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->orig_data));
375f80ce2aSJacob Faibussowitsch   for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS ; level++) {
3818c3aa7eSMark     mg->min_eigen_DinvA[level] = 0;
3918c3aa7eSMark     mg->max_eigen_DinvA[level] = 0;
4018c3aa7eSMark   }
4118c3aa7eSMark   pc_gamg->emin = 0;
4218c3aa7eSMark   pc_gamg->emax = 0;
43a2f3521dSMark F. Adams   PetscFunctionReturn(0);
44a2f3521dSMark F. Adams }
45a2f3521dSMark F. Adams 
465b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
475b89ad90SMark F. Adams /*
48c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
49a147abb0SMark F. Adams      of active processors.
505b89ad90SMark F. Adams 
515b89ad90SMark F. Adams    Input Parameter:
52a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
53a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
549d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
55c5bfad50SMark F. Adams    . cr_bs - coarse block size
563530afc2SMark F. Adams    In/Output Parameter:
57a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
58afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
5911e60469SMark F. Adams    Output Parameter:
603530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
615b89ad90SMark F. Adams */
62171cca9aSMark Adams static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,IS * Pcolumnperm, PetscBool is_last)
635b89ad90SMark F. Adams {
649d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
65486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
66a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
673b4367a7SBarry Smith   MPI_Comm        comm;
68c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
693ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
705b89ad90SMark F. Adams 
715b89ad90SMark F. Adams   PetscFunctionBegin;
729566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat_fine,&comm));
739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
759566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat_fine, &f_bs));
76849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP],0,0,0,0));
779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0));
789566063dSJacob Faibussowitsch   PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat));
799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0));
80849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP],0,0,0,0));
81038e3b61SMark F. Adams 
82ce7c7f2fSMark Adams   if (Pcolumnperm) *Pcolumnperm = NULL;
83ce7c7f2fSMark Adams 
843ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
859566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL));
863ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
873ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
8873911c69SBarry Smith   } else {
893ae0bb68SMark Adams     PetscInt  bs;
909566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(Cmat, &bs));
913ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
923ae0bb68SMark Adams   }
93c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
94c9567895SMark   if (pc_gamg->level_reduction_factors[pc_gamg->current_level] == 0 && PetscDefined(HAVE_CUDA) && pc_gamg->current_level==0) { /* 0 turns reducing to 1 process/device on; do for HIP, etc. */
95c9567895SMark #if defined(PETSC_HAVE_CUDA)
96c9567895SMark     PetscShmComm pshmcomm;
97c9567895SMark     PetscMPIInt  locrank;
98c9567895SMark     MPI_Comm     loccomm;
99c9567895SMark     PetscInt     s_nnodes,r_nnodes, new_new_size;
100c9567895SMark     cudaError_t  cerr;
101c9567895SMark     int          devCount;
1029566063dSJacob Faibussowitsch     PetscCall(PetscShmCommGet(comm,&pshmcomm));
1039566063dSJacob Faibussowitsch     PetscCall(PetscShmCommGetMpiShmComm(pshmcomm,&loccomm));
1049566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(loccomm, &locrank));
105c9567895SMark     s_nnodes = !locrank;
1069566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm));
10763a3b9bcSJacob Faibussowitsch     PetscCheck((size%r_nnodes) == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of nodes np=%d nnodes%" PetscInt_FMT,size,r_nnodes);
108c9567895SMark     devCount = 0;
109c9567895SMark     cerr = cudaGetDeviceCount(&devCount);
110c9567895SMark     cudaGetLastError(); /* Reset the last error */
111c9567895SMark     if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */
112c9567895SMark       new_new_size = r_nnodes * devCount;
113c9567895SMark       new_size = new_new_size;
11463a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: Fine grid with Cuda. %" PetscInt_FMT " nodes. Change new active set size %d --> %d (devCount=%d #nodes=%" PetscInt_FMT ")\n",((PetscObject)pc)->prefix,r_nnodes,nactive,new_size,devCount,r_nnodes));
115c9567895SMark     } else {
11663a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: With Cuda but no device. Use heuristics.",((PetscObject)pc)->prefix));
117c9567895SMark       goto HEURISTIC;
118c9567895SMark     }
119c9567895SMark #else
120c9567895SMark     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here");
121c9567895SMark #endif
122c9567895SMark   } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) {
12363a3b9bcSJacob Faibussowitsch     PetscCheck(nactive%pc_gamg->level_reduction_factors[pc_gamg->current_level] == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of active process %d wrt reduction factor %" PetscInt_FMT,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]);
124c9567895SMark     new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level];
12563a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"%s: Manually setting reduction to %d active processes (%d/%" PetscInt_FMT ")\n",((PetscObject)pc)->prefix,new_size,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]));
126c9567895SMark   } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) {
127c9567895SMark     new_size = 1;
1289566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"%s: Force coarsest grid reduction to %d active processes\n",((PetscObject)pc)->prefix,new_size));
129c9567895SMark   } else {
130472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
131c9567895SMark #if defined(PETSC_HAVE_CUDA)
132c9567895SMark     HEURISTIC:
133c9567895SMark #endif
1349566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL));
135a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
1367f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
137c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
1389566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"%s: Coarse grid reduction from %d to %d active processes\n",((PetscObject)pc)->prefix,nactive,new_size));
139a2f3521dSMark F. Adams   }
1402e3501ffSMark Adams   if (new_size==nactive) {
141ef3f0257SMark Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
142ce7c7f2fSMark Adams     if (new_size < size) {
143ce7c7f2fSMark Adams       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
1449566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",((PetscObject)pc)->prefix,nactive));
145ce7c7f2fSMark Adams       if (pc_gamg->cpu_pin_coarse_grids) {
1469566063dSJacob Faibussowitsch         PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE));
1479566063dSJacob Faibussowitsch         PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE));
148ce7c7f2fSMark Adams       }
149ce7c7f2fSMark Adams     }
150ef3f0257SMark Adams     /* we know that the grid structure can be reused in MatPtAP */
1512e3501ffSMark Adams   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
152192c0e8bSMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1;
153885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
154849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REDUCE],0,0,0,0));
15571959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
15663a3b9bcSJacob Faibussowitsch     PetscCheck(ncrs_eq % cr_bs == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %" PetscInt_FMT " not divisible by cr_bs %" PetscInt_FMT,ncrs_eq,cr_bs);
157ce7c7f2fSMark Adams     /* get new_size and rfactor */
158ce7c7f2fSMark Adams     if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
159ce7c7f2fSMark Adams       /* find factor */
160ce7c7f2fSMark Adams       if (new_size == 1) rfactor = size; /* don't modify */
161ce7c7f2fSMark Adams       else {
162ce7c7f2fSMark Adams         PetscReal best_fact = 0.;
163ce7c7f2fSMark Adams         jj = -1;
164ce7c7f2fSMark Adams         for (kk = 1 ; kk <= size ; kk++) {
165ce7c7f2fSMark Adams           if (!(size%kk)) { /* a candidate */
166ce7c7f2fSMark Adams             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
167ce7c7f2fSMark Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
168ce7c7f2fSMark Adams             if (fact > best_fact) {
169ce7c7f2fSMark Adams               best_fact = fact; jj = kk;
170ce7c7f2fSMark Adams             }
171ce7c7f2fSMark Adams           }
172ce7c7f2fSMark Adams         }
173ce7c7f2fSMark Adams         if (jj != -1) rfactor = jj;
174ce7c7f2fSMark Adams         else rfactor = 1; /* a prime */
175ce7c7f2fSMark Adams         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
176ce7c7f2fSMark Adams         else expand_factor = rfactor;
177ce7c7f2fSMark Adams       }
178ce7c7f2fSMark Adams       new_size = size/rfactor; /* make new size one that is factor */
1794cdfd227SMark       if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
1804cdfd227SMark         *a_Amat_crs = Cmat;
18163a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(pc,"%s: Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%" PetscInt_FMT "\n",((PetscObject)pc)->prefix,new_size,ncrs_eq));
182849bee69SMark Adams         PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE],0,0,0,0));
183ce7c7f2fSMark Adams         PetscFunctionReturn(0);
184ce7c7f2fSMark Adams       }
185ce7c7f2fSMark Adams     }
186a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
1879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &counts));
188849bee69SMark Adams     if (pc_gamg->repart) { /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */
1895a9b9e01SMark F. Adams       Mat      adj;
190849bee69SMark Adams       PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REPART],0,0,0,0));
19163a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: Repartition: size (active): %d --> %d, %" PetscInt_FMT " local equations, using %s process layout\n",((PetscObject)pc)->prefix,*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread"));
192a2f3521dSMark F. Adams       /* get 'adj' */
193c5bfad50SMark F. Adams       if (cr_bs == 1) {
1949566063dSJacob Faibussowitsch         PetscCall(MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
195806fa848SBarry Smith       } else {
196a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
197eb07cef2SMark F. Adams         Mat               tMat;
198a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
199b4fbaa2aSMark F. Adams         const PetscScalar *vals;
200b4fbaa2aSMark F. Adams         const PetscInt    *idx;
201a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
20239d09545SMark Adams         static PetscInt   llev = 0; /* ugly but just used for debugging */
203d9558ea9SBarry Smith         MatType           mtype;
204b4fbaa2aSMark F. Adams 
2059566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz));
2069566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs));
2079566063dSJacob Faibussowitsch         PetscCall(MatGetSize(Cmat, &M, &N));
208c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
2099566063dSJacob Faibussowitsch           PetscCall(MatGetRow(Cmat,Ii,&ncols,NULL,NULL));
210c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
211c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
2129566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL));
2133ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
2143ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
21558471d46SMark F. Adams         }
2166876a03eSMark F. Adams 
2179566063dSJacob Faibussowitsch         PetscCall(MatGetType(Amat_fine,&mtype));
2189566063dSJacob Faibussowitsch         PetscCall(MatCreate(comm, &tMat));
2199566063dSJacob Faibussowitsch         PetscCall(MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE));
2209566063dSJacob Faibussowitsch         PetscCall(MatSetType(tMat,mtype));
2219566063dSJacob Faibussowitsch         PetscCall(MatSeqAIJSetPreallocation(tMat,0,d_nnz));
2229566063dSJacob Faibussowitsch         PetscCall(MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz));
2239566063dSJacob Faibussowitsch         PetscCall(PetscFree2(d_nnz,o_nnz));
224eb07cef2SMark F. Adams 
225a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
226c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
2279566063dSJacob Faibussowitsch           PetscCall(MatGetRow(Cmat,ii,&ncols,&idx,&vals));
228eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
229c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
230eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
2319566063dSJacob Faibussowitsch             PetscCall(MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES));
232eb07cef2SMark F. Adams           }
2339566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(Cmat,ii,&ncols,&idx,&vals));
234eb07cef2SMark F. Adams         }
2359566063dSJacob Faibussowitsch         PetscCall(MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY));
2369566063dSJacob Faibussowitsch         PetscCall(MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY));
237eb07cef2SMark F. Adams 
238b4fbaa2aSMark F. Adams         if (llev++ == -1) {
239b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
24063a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(fname,sizeof(fname),"part_mat_%" PetscInt_FMT ".mat",llev));
2413b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
2429566063dSJacob Faibussowitsch           PetscCall(MatView(tMat, viewer));
2439566063dSJacob Faibussowitsch           PetscCall(PetscViewerDestroy(&viewer));
244b4fbaa2aSMark F. Adams         }
2459566063dSJacob Faibussowitsch         PetscCall(MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
2469566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tMat));
247a2f3521dSMark F. Adams       } /* create 'adj' */
248f150b916SMark F. Adams 
249a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2505a9b9e01SMark F. Adams         char            prefix[256];
2515a9b9e01SMark F. Adams         const char      *pcpre;
252b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
253b4fbaa2aSMark F. Adams         MatPartitioning mpart;
254a4b7d37bSMark F. Adams         IS              proc_is;
2552f03bc48SMark F. Adams 
2569566063dSJacob Faibussowitsch         PetscCall(MatPartitioningCreate(comm, &mpart));
2579566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetAdjacency(mpart, adj));
2589566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &pcpre));
2599566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : ""));
2609566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix));
2619566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetFromOptions(mpart));
2629566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetNParts(mpart, new_size));
2639566063dSJacob Faibussowitsch         PetscCall(MatPartitioningApply(mpart, &proc_is));
2649566063dSJacob Faibussowitsch         PetscCall(MatPartitioningDestroy(&mpart));
2655a9b9e01SMark F. Adams 
2665ef31b24SMark F. Adams         /* collect IS info */
2679566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(ncrs_eq, &newproc_idx));
2689566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(proc_is, &is_idx));
269a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
270c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
271ce7c7f2fSMark Adams             newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
272eb07cef2SMark F. Adams           }
2735ef31b24SMark F. Adams         }
2749566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(proc_is, &is_idx));
2759566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&proc_is));
2765ef31b24SMark F. Adams       }
2779566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&adj));
2785a9b9e01SMark F. Adams 
2799566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc));
2809566063dSJacob Faibussowitsch       PetscCall(PetscFree(newproc_idx));
281849bee69SMark Adams       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART],0,0,0,0));
28231cb4603SMark Adams     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
283ce7c7f2fSMark Adams       PetscInt targetPE;
28408401ef6SPierre Jolivet       PetscCheck(new_size!=nactive,PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen");
28563a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: Number of equations (loc) %" PetscInt_FMT " with simple aggregation\n",((PetscObject)pc)->prefix,ncrs_eq));
286ce7c7f2fSMark Adams       targetPE = (rank/rfactor)*expand_factor;
2879566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc));
288a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
289e33ef3b1SMark F. Adams 
29011e60469SMark F. Adams     /*
291a2f3521dSMark F. Adams       Create an index set from the is_eq_newproc index set to indicate the mapping TO
29211e60469SMark F. Adams     */
2939566063dSJacob Faibussowitsch     PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num));
2947700e67bSMark Adams     is_eq_num_prim = is_eq_num;
29511e60469SMark F. Adams     /*
296a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
29711e60469SMark F. Adams     */
2989566063dSJacob Faibussowitsch     PetscCall(ISPartitioningCount(is_eq_newproc, size, counts));
299c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
3009566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_eq_newproc));
301ce7c7f2fSMark Adams     ncrs_new = ncrs_eq_new/cr_bs;
302a2f3521dSMark F. Adams 
3039566063dSJacob Faibussowitsch     PetscCall(PetscFree(counts));
3046aad120cSJose E. Roman     /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxiliary data into some complex abstracted thing */
305885364a3SMark Adams     {
306885364a3SMark Adams       Vec            src_crd, dest_crd;
307885364a3SMark 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;
308885364a3SMark Adams       VecScatter     vecscat;
309885364a3SMark Adams       PetscScalar    *array;
310885364a3SMark Adams       IS isscat;
311a2f3521dSMark F. Adams       /* move data (for primal equations only) */
31222063be5SMark F. Adams       /* Create a vector to contain the newly ordered element information */
3139566063dSJacob Faibussowitsch       PetscCall(VecCreate(comm, &dest_crd));
3149566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE));
3159566063dSJacob Faibussowitsch       PetscCall(VecSetType(dest_crd,VECSTANDARD)); /* this is needed! */
31611e60469SMark F. Adams       /*
3179d5b6da9SMark F. Adams         There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
318c5bfad50SMark F. Adams         a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
31911e60469SMark F. Adams       */
3209566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(ncrs*node_data_sz, &tidx));
3219566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(is_eq_num_prim, &idx));
3223ae0bb68SMark Adams       for (ii=0,jj=0; ii<ncrs; ii++) {
323c5bfad50SMark F. Adams         PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
324a2f3521dSMark F. Adams         for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
32511e60469SMark F. Adams       }
3269566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(is_eq_num_prim, &idx));
3279566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat));
3289566063dSJacob Faibussowitsch       PetscCall(PetscFree(tidx));
32911e60469SMark F. Adams       /*
33011e60469SMark F. Adams         Create a vector to contain the original vertex information for each element
33111e60469SMark F. Adams       */
3329566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd));
3339d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3343ae0bb68SMark Adams         const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3353ae0bb68SMark Adams         for (ii=0; ii<ncrs; ii++) {
3369d5b6da9SMark F. Adams           for (kk=0; kk<ndata_rows; kk++) {
337a2f3521dSMark F. Adams             PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
338c8b0795cSMark F. Adams             PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
3399566063dSJacob Faibussowitsch             PetscCall(VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES));
340d3d6bff4SMark F. Adams           }
341038e3b61SMark F. Adams         }
342eb07cef2SMark F. Adams       }
3439566063dSJacob Faibussowitsch       PetscCall(VecAssemblyBegin(src_crd));
3449566063dSJacob Faibussowitsch       PetscCall(VecAssemblyEnd(src_crd));
34511e60469SMark F. Adams       /*
34611e60469SMark F. Adams         Scatter the element vertex information (still in the original vertex ordering)
34711e60469SMark F. Adams         to the correct processor
34811e60469SMark F. Adams       */
3499566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat));
3509566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&isscat));
3519566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD));
3529566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD));
3539566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&vecscat));
3549566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&src_crd));
35511e60469SMark F. Adams       /*
35611e60469SMark F. Adams         Put the element vertex data into a new allocation of the gdata->ele
35711e60469SMark F. Adams       */
3589566063dSJacob Faibussowitsch       PetscCall(PetscFree(pc_gamg->data));
3599566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data));
3602fa5cd67SKarl Rupp 
3613ae0bb68SMark Adams       pc_gamg->data_sz = node_data_sz*ncrs_new;
3623ae0bb68SMark Adams       strideNew        = ncrs_new*ndata_rows;
3632fa5cd67SKarl Rupp 
3649566063dSJacob Faibussowitsch       PetscCall(VecGetArray(dest_crd, &array));
3659d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3663ae0bb68SMark Adams         for (ii=0; ii<ncrs_new; ii++) {
3679d5b6da9SMark F. Adams           for (kk=0; kk<ndata_rows; kk++) {
368a2f3521dSMark F. Adams             PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
369c8b0795cSMark F. Adams             pc_gamg->data[ix] = PetscRealPart(array[jx]);
370d3d6bff4SMark F. Adams           }
371038e3b61SMark F. Adams         }
372038e3b61SMark F. Adams       }
3739566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(dest_crd, &array));
3749566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&dest_crd));
375885364a3SMark Adams     }
376a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
377849bee69SMark Adams     /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0)); */
37811e60469SMark F. Adams     /*
3797dae84e0SHong Zhang       Invert for MatCreateSubMatrix
38011e60469SMark F. Adams     */
3819566063dSJacob Faibussowitsch     PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices));
3829566063dSJacob Faibussowitsch     PetscCall(ISSort(new_eq_indices)); /* is this needed? */
3839566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize(new_eq_indices, cr_bs));
384a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
3859566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */
386a2f3521dSMark F. Adams     }
3873cb8563fSToby Isaac     if (Pcolumnperm) {
3889566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)new_eq_indices));
3893cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3903cb8563fSToby Isaac     }
3919566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_eq_num));
392849bee69SMark Adams     /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0)); */
393849bee69SMark Adams     /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0)); */
394849bee69SMark Adams 
395a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
396a2f3521dSMark F. Adams     {
397a2f3521dSMark F. Adams       Mat       mat;
39890db8557SMark Adams       PetscBool flg;
3999566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat));
400*5d5a09d9SMark Adams       PetscCall(MatGetOption(Cmat, MAT_SPD, &flg)); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ?
40190db8557SMark Adams       if (flg) {
4029566063dSJacob Faibussowitsch         PetscCall(MatSetOption(mat, MAT_SPD,PETSC_TRUE));
40390db8557SMark Adams       } else {
4049566063dSJacob Faibussowitsch         PetscCall(MatGetOption(Cmat, MAT_HERMITIAN, &flg));
40590db8557SMark Adams         if (flg) {
4069566063dSJacob Faibussowitsch           PetscCall(MatSetOption(mat, MAT_HERMITIAN,PETSC_TRUE));
40790db8557SMark Adams         } else {
40890db8557SMark Adams #if !defined(PETSC_USE_COMPLEX)
4099566063dSJacob Faibussowitsch           PetscCall(MatGetOption(Cmat, MAT_SYMMETRIC, &flg));
41090db8557SMark Adams           if (flg) {
4119566063dSJacob Faibussowitsch             PetscCall(MatSetOption(mat, MAT_SYMMETRIC,PETSC_TRUE));
41290db8557SMark Adams           }
41390db8557SMark Adams #endif
41490db8557SMark Adams         }
41590db8557SMark Adams       }
416a2f3521dSMark F. Adams       *a_Amat_crs = mat;
417a2f3521dSMark F. Adams     }
4189566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Cmat));
419a2f3521dSMark F. Adams 
420849bee69SMark Adams     /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0)); */
42111e60469SMark F. Adams     /* prolongator */
42211e60469SMark F. Adams     {
42311e60469SMark F. Adams       IS       findices;
424a2f3521dSMark F. Adams       PetscInt Istart,Iend;
425a2f3521dSMark F. Adams       Mat      Pnew;
42662294041SBarry Smith 
4279566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend));
428849bee69SMark Adams       /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0)); */
4299566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(comm,Iend-Istart,Istart,1,&findices));
4309566063dSJacob Faibussowitsch       PetscCall(ISSetBlockSize(findices,f_bs));
4319566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew));
4329566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&findices));
4339566063dSJacob Faibussowitsch       PetscCall(MatSetOption(Pnew,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE));
434c5bfad50SMark F. Adams 
435849bee69SMark Adams       /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0)); */
4369566063dSJacob Faibussowitsch       PetscCall(MatDestroy(a_P_inout));
437a2f3521dSMark F. Adams 
438a2f3521dSMark F. Adams       /* output - repartitioned */
439a2f3521dSMark F. Adams       *a_P_inout = Pnew;
440e33ef3b1SMark F. Adams     }
4419566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&new_eq_indices));
4425b89ad90SMark F. Adams 
443c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
444ce7c7f2fSMark Adams 
445ce7c7f2fSMark Adams     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
446ce7c7f2fSMark Adams     if (pc_gamg->cpu_pin_coarse_grids) {
447ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
4488bca76a6SMark Adams       static PetscInt llev = 2;
44963a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: Pinning level %" PetscInt_FMT " to the CPU\n",((PetscObject)pc)->prefix,llev++));
450ce7c7f2fSMark Adams #endif
4519566063dSJacob Faibussowitsch       PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE));
4529566063dSJacob Faibussowitsch       PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE));
453adf5291fSStefano Zampini       if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
454ce7c7f2fSMark Adams         Mat         A = *a_Amat_crs, P = *a_P_inout;
455ce7c7f2fSMark Adams         PetscMPIInt size;
4569566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
457ce7c7f2fSMark Adams         if (size > 1) {
458ce7c7f2fSMark Adams           Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data;
4599566063dSJacob Faibussowitsch           PetscCall(VecBindToCPU(a->lvec,PETSC_TRUE));
4609566063dSJacob Faibussowitsch           PetscCall(VecBindToCPU(p->lvec,PETSC_TRUE));
461ce7c7f2fSMark Adams         }
462ce7c7f2fSMark Adams       }
463ce7c7f2fSMark Adams     }
464849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE],0,0,0,0));
465849bee69SMark Adams   } // processor reduce
4665b89ad90SMark F. Adams   PetscFunctionReturn(0);
4675b89ad90SMark F. Adams }
4685b89ad90SMark F. Adams 
4694b1575e2SStefano Zampini PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2)
4704b1575e2SStefano Zampini {
4714b1575e2SStefano Zampini   const char     *prefix;
4724b1575e2SStefano Zampini   char           addp[32];
4734b1575e2SStefano Zampini   PC_MG          *mg      = (PC_MG*)a_pc->data;
4744b1575e2SStefano Zampini   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4754b1575e2SStefano Zampini 
4764b1575e2SStefano Zampini   PetscFunctionBegin;
4779566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(a_pc,&prefix));
47863a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(a_pc,"%s: Square Graph on level %" PetscInt_FMT "\n",((PetscObject)a_pc)->prefix,pc_gamg->current_level+1));
4799566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(Gmat1,Gmat1,NULL,Gmat2));
4809566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(*Gmat2,prefix));
48163a3b9bcSJacob Faibussowitsch   PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%" PetscInt_FMT "_",pc_gamg->current_level));
4829566063dSJacob Faibussowitsch   PetscCall(MatAppendOptionsPrefix(*Gmat2,addp));
483b4da3a1bSStefano Zampini   if ((*Gmat2)->structurally_symmetric) {
4849566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AB));
485b4da3a1bSStefano Zampini   } else {
4869566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Gmat1,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE));
4879566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AtB));
488b4da3a1bSStefano Zampini   }
4899566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(*Gmat2));
4909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0));
4919566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(*Gmat2));
4929566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0));
4939566063dSJacob Faibussowitsch   PetscCall(MatProductClear(*Gmat2));
4944b1575e2SStefano Zampini   /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
4954b1575e2SStefano Zampini   (*Gmat2)->assembled = PETSC_TRUE;
4964b1575e2SStefano Zampini   PetscFunctionReturn(0);
4974b1575e2SStefano Zampini }
4984b1575e2SStefano Zampini 
4995b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5005b89ad90SMark F. Adams /*
5015b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5025b89ad90SMark F. Adams                     by setting data structures and options.
5035b89ad90SMark F. Adams 
5045b89ad90SMark F. Adams    Input Parameter:
5055b89ad90SMark F. Adams .  pc - the preconditioner context
5065b89ad90SMark F. Adams 
5075b89ad90SMark F. Adams */
5089d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5095b89ad90SMark F. Adams {
5109d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
5115b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
5122adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
51318c3aa7eSMark   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS];
5143b4367a7SBarry Smith   MPI_Comm       comm;
515c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
51618c3aa7eSMark   Mat            Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS];
51718c3aa7eSMark   IS             *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
518a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
519569f4572SMark Adams   MatInfo        info;
520171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
5215ef31b24SMark F. Adams 
5225b89ad90SMark F. Adams   PetscFunctionBegin;
5239566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
5249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
5259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
526849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0));
5278abdc6daSStefano Zampini   if (pc->setupcalled) {
5288abdc6daSStefano Zampini     if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) {
529878e152fSMark F. Adams       /* reset everything */
5309566063dSJacob Faibussowitsch       PetscCall(PCReset_MG(pc));
531878e152fSMark F. Adams       pc->setupcalled = 0;
532806fa848SBarry Smith     } else {
53384d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
53403a628feSMark F. Adams       /* just do Galerkin grids */
53558471d46SMark F. Adams       Mat          B,dA,dB;
5369d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
5374555aa8cSStefano Zampini         PetscInt gl;
53858471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5399566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB));
54058471d46SMark F. Adams         /* (re)set to get dirty flag */
5419566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB));
54258471d46SMark F. Adams 
5434555aa8cSStefano Zampini         for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) {
5448abdc6daSStefano Zampini           MatReuse reuse = MAT_INITIAL_MATRIX ;
545849bee69SMark Adams #if defined(GAMG_STAGES)
546849bee69SMark Adams           PetscCall(PetscLogStagePush(gamg_stages[gl]));
547849bee69SMark Adams #endif
5488abdc6daSStefano Zampini           /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
5499566063dSJacob Faibussowitsch           PetscCall(KSPGetOperators(mglevels[level]->smoothd,NULL,&B));
5508abdc6daSStefano Zampini           if (B->product) {
5518abdc6daSStefano Zampini             if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) {
5528abdc6daSStefano Zampini               reuse = MAT_REUSE_MATRIX;
55303a628feSMark F. Adams             }
5548abdc6daSStefano Zampini           }
5559566063dSJacob Faibussowitsch           if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A));
5568abdc6daSStefano Zampini           if (reuse == MAT_REUSE_MATRIX) {
55763a3b9bcSJacob Faibussowitsch             PetscCall(PetscInfo(pc,"%s: RAP after first solve, reuse matrix level %" PetscInt_FMT "\n",((PetscObject)pc)->prefix,level));
5588abdc6daSStefano Zampini           } else {
55963a3b9bcSJacob Faibussowitsch             PetscCall(PetscInfo(pc,"%s: RAP after first solve, new matrix level %" PetscInt_FMT "\n",((PetscObject)pc)->prefix,level));
5608abdc6daSStefano Zampini           }
5619566063dSJacob Faibussowitsch           PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0));
5629566063dSJacob Faibussowitsch           PetscCall(MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B));
5639566063dSJacob Faibussowitsch           PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0));
56463b77682SMark Adams           if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B;
5659566063dSJacob Faibussowitsch           PetscCall(KSPSetOperators(mglevels[level]->smoothd,B,B));
56658471d46SMark F. Adams           dB   = B;
567849bee69SMark Adams #if defined(GAMG_STAGES)
568849bee69SMark Adams           PetscCall(PetscLogStagePop());
569849bee69SMark Adams #endif
57058471d46SMark F. Adams         }
5715f8cf99dSMark F. Adams       }
572d5280255SMark F. Adams 
5739566063dSJacob Faibussowitsch       PetscCall(PCSetUp_MG(pc));
574849bee69SMark Adams       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0));
57558471d46SMark F. Adams       PetscFunctionReturn(0);
576eb07cef2SMark F. Adams     }
577878e152fSMark F. Adams   }
578f6536408SMark F. Adams 
579878e152fSMark F. Adams   if (!pc_gamg->data) {
580878e152fSMark F. Adams     if (pc_gamg->orig_data) {
5819566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSize(Pmat, &bs));
5829566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(Pmat, &qq, NULL));
5832fa5cd67SKarl Rupp 
584878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
585878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
586878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5872fa5cd67SKarl Rupp 
5889566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data));
589878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
590806fa848SBarry Smith     } else {
5915f80ce2aSJacob Faibussowitsch       PetscCheck(pc_gamg->ops->createdefaultdata,comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5929566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->createdefaultdata(pc,Pmat));
5939d5b6da9SMark F. Adams     }
594878e152fSMark F. Adams   }
595878e152fSMark F. Adams 
596878e152fSMark F. Adams   /* cache original data for reuse */
5971c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
5989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data));
599878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
600878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
601878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
602878e152fSMark F. Adams   }
603038e3b61SMark F. Adams 
604302f38e8SMark F. Adams   /* get basic dims */
6059566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Pmat, &bs));
6069566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Pmat, &M, &N));
60784d3f75bSMark F. Adams 
6089566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info)); /* global reduction */
609569f4572SMark Adams   nnz0   = info.nz_used;
610569f4572SMark Adams   nnztot = info.nz_used;
61163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc,"%s: level %" PetscInt_FMT ") N=0, n data rows=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", np=%d\n",((PetscObject)pc)->prefix,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(PetscInt)(nnz0/(PetscReal)M+0.5),size));
612569f4572SMark Adams 
613a2f3521dSMark F. Adams   /* Get A_i and R_i */
61462294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
6159ab59c8bSMark Adams     pc_gamg->current_level = level;
61663a3b9bcSJacob Faibussowitsch     PetscCheck(level < PETSC_MG_MAXLEVELS,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %" PetscInt_FMT,level);
6175b89ad90SMark F. Adams     level1 = level + 1;
6184555aa8cSStefano Zampini #if defined(GAMG_STAGES)
619849bee69SMark Adams     if (!gamg_stages[level]) {
620849bee69SMark Adams       char     str[32];
621849bee69SMark Adams       sprintf(str,"GAMG Level %d",(int)level);
622849bee69SMark Adams       PetscCall(PetscLogStageRegister(str, &gamg_stages[level]));
623849bee69SMark Adams     }
6249566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(gamg_stages[level]));
625b4fbaa2aSMark F. Adams #endif
626c8b0795cSMark F. Adams     { /* construct prolongator */
627725b86d8SJed Brown       Mat              Gmat;
6280cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
6297700e67bSMark Adams       Mat              Prol11;
630c8b0795cSMark F. Adams 
6319566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->graph(pc,Aarr[level], &Gmat));
6329566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists));
6339566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11));
634c8b0795cSMark F. Adams 
635a2f3521dSMark F. Adams       /* could have failed to create new level */
636a2f3521dSMark F. Adams       if (Prol11) {
637f7df55f0SStefano Zampini         const char *prefix;
638f7df55f0SStefano Zampini         char       addp[32];
639f7df55f0SStefano Zampini 
6409d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6419566063dSJacob Faibussowitsch         PetscCall(MatGetBlockSizes(Prol11, NULL, &bs));
642a2f3521dSMark F. Adams 
643fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
644c8b0795cSMark F. Adams           /* smooth */
6459566063dSJacob Faibussowitsch           PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11));
646c8b0795cSMark F. Adams         }
647c8b0795cSMark F. Adams 
6480c3bc534SBarry Smith         if (pc_gamg->use_aggs_in_asm) {
6491b18a24aSMark Adams           PetscInt bs;
650849bee69SMark Adams           PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // not timed directly, ugly, could remove, but good ASM method
6519566063dSJacob Faibussowitsch           PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]));
652ffc955d6SMark F. Adams         }
653ffc955d6SMark F. Adams 
6549566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc,&prefix));
6559566063dSJacob Faibussowitsch         PetscCall(MatSetOptionsPrefix(Prol11,prefix));
6569566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level));
6579566063dSJacob Faibussowitsch         PetscCall(MatAppendOptionsPrefix(Prol11,addp));
65891f31d3dSStefano Zampini         /* Always generate the transpose with CUDA
659f7df55f0SStefano Zampini            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
6609566063dSJacob Faibussowitsch         PetscCall(MatSetOption(Prol11,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE));
6619566063dSJacob Faibussowitsch         PetscCall(MatSetFromOptions(Prol11));
6624bde40a0SMark Adams         Parr[level1] = Prol11;
6634bde40a0SMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
6644bde40a0SMark Adams 
6659566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Gmat));
6669566063dSJacob Faibussowitsch       PetscCall(PetscCDDestroy(agg_lists));
667a2f3521dSMark F. Adams     } /* construct prolongator scope */
6687f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
669171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
67063a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: Stop gridding, level %" PetscInt_FMT "\n",((PetscObject)pc)->prefix,level));
6714555aa8cSStefano Zampini #if defined(GAMG_STAGES)
6729566063dSJacob Faibussowitsch       PetscCall(PetscLogStagePop());
673a90e85d9SMark Adams #endif
674c8b0795cSMark F. Adams       break;
675c8b0795cSMark F. Adams     }
6769566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */
6772472a847SBarry Smith     PetscCheck(!is_last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?");
678171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
6790e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
680849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL],0,0,0,0));
6819566063dSJacob Faibussowitsch     PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last));
682849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL],0,0,0,0));
683a2f3521dSMark F. Adams 
6849566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */
6859566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info));
686569f4572SMark Adams     nnztot += info.nz_used;
68763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"%s: %" PetscInt_FMT ") N=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", %d active pes\n",((PetscObject)pc)->prefix,level1,M,pc_gamg->data_cell_cols,(PetscInt)(info.nz_used/(PetscReal)M),nactivepe));
688569f4572SMark Adams 
6894555aa8cSStefano Zampini #if defined(GAMG_STAGES)
6909566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
691b4fbaa2aSMark F. Adams #endif
692a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
6939ab59c8bSMark Adams     if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) {
69463a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: HARD stop of coarsening on level %" PetscInt_FMT ".  Grid too small: %" PetscInt_FMT " block nodes\n",((PetscObject)pc)->prefix,level,M/bs));
695a90e85d9SMark Adams       level++;
696a90e85d9SMark Adams       break;
697a90e85d9SMark Adams     }
698c8b0795cSMark F. Adams   } /* levels */
6999566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->data));
700c8b0795cSMark F. Adams 
70163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc,"%s: %" PetscInt_FMT " levels, grid complexity = %g\n",((PetscObject)pc)->prefix,level+1,nnztot/nnz0));
7029d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7035b89ad90SMark F. Adams   fine_level       = level;
7049566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc,pc_gamg->Nlevels,NULL));
7055b89ad90SMark F. Adams 
70662294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
7070ed2132dSStefano Zampini 
708d5280255SMark F. Adams     /* set default smoothers & set operators */
70962294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
710ffc955d6SMark F. Adams       KSP smoother;
711ffc955d6SMark F. Adams       PC  subpc;
712a2f3521dSMark F. Adams 
7139566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
7149566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(smoother, &subpc));
715ffc955d6SMark F. Adams 
7169566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
717a2f3521dSMark F. Adams       /* set ops */
7189566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level]));
7199566063dSJacob Faibussowitsch       PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level+1]));
720a2f3521dSMark F. Adams 
721a2f3521dSMark F. Adams       /* set defaults */
7229566063dSJacob Faibussowitsch       PetscCall(KSPSetType(smoother, KSPCHEBYSHEV));
723a2f3521dSMark F. Adams 
7240c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
7250c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
7262d3561bbSSatish Balay         PetscInt sz;
7277a28f3e5SMark Adams         IS       *iss;
728a2f3521dSMark F. Adams 
7292d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
7307a28f3e5SMark Adams         iss  = ASMLocalIDsArr[level];
7319566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCASM));
7329566063dSJacob Faibussowitsch         PetscCall(PCASMSetOverlap(subpc, 0));
7339566063dSJacob Faibussowitsch         PetscCall(PCASMSetType(subpc,PC_ASM_BASIC));
7347f66b68fSMark Adams         if (!sz) {
735ffc955d6SMark F. Adams           IS       is;
7369566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is));
7379566063dSJacob Faibussowitsch           PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is));
7389566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&is));
739806fa848SBarry Smith         } else {
740a94c3b12SMark F. Adams           PetscInt kk;
7419566063dSJacob Faibussowitsch           PetscCall(PCASMSetLocalSubdomains(subpc, sz, NULL, iss));
742a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
7439566063dSJacob Faibussowitsch             PetscCall(ISDestroy(&iss[kk]));
744a94c3b12SMark F. Adams           }
7459566063dSJacob Faibussowitsch           PetscCall(PetscFree(iss));
746ffc955d6SMark F. Adams         }
7470298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
748ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
749806fa848SBarry Smith       } else {
7509566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCJACOBI));
751ffc955d6SMark F. Adams       }
752d5280255SMark F. Adams     }
753d5280255SMark F. Adams     {
754d5280255SMark F. Adams       /* coarse grid */
755d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
756d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
7570ed2132dSStefano Zampini 
7589566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
7599566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(smoother, Lmat, Lmat));
760cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
7619566063dSJacob Faibussowitsch         PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
7629566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(smoother, &subpc));
7639566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCBJACOBI));
7649566063dSJacob Faibussowitsch         PetscCall(PCSetUp(subpc));
7659566063dSJacob Faibussowitsch         PetscCall(PCBJacobiGetSubKSP(subpc,&ii,&first,&k2));
76663a3b9bcSJacob Faibussowitsch         PetscCheck(ii == 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %" PetscInt_FMT " is not one",ii);
7679566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(k2[0],&pc2));
7689566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc2, PCLU));
7699566063dSJacob Faibussowitsch         PetscCall(PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS));
7709566063dSJacob Faibussowitsch         PetscCall(KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1));
7719566063dSJacob Faibussowitsch         PetscCall(KSPSetType(k2[0], KSPPREONLY));
772d5280255SMark F. Adams       }
773cf8ae1d3SMark Adams     }
774d5280255SMark F. Adams 
775d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
776d0609cedSBarry Smith     PetscObjectOptionsBegin((PetscObject)pc);
7779566063dSJacob Faibussowitsch     PetscCall(PCSetFromOptions_MG(PetscOptionsObject,pc));
778d0609cedSBarry Smith     PetscOptionsEnd();
7799566063dSJacob Faibussowitsch     PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL));
780d5280255SMark F. Adams 
78118c3aa7eSMark     /* setup cheby eigen estimates from SA */
7827e6512fdSJed Brown     if (pc_gamg->use_sa_esteig) {
78318c3aa7eSMark       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
78418c3aa7eSMark         KSP       smoother;
78518c3aa7eSMark         PetscBool ischeb;
7860ed2132dSStefano Zampini 
7879566063dSJacob Faibussowitsch         PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
7889566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb));
78918c3aa7eSMark         if (ischeb) {
79018c3aa7eSMark           KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data;
7910ed2132dSStefano Zampini 
7922de708cbSJed Brown           // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG
793efe053fcSJed Brown           if (mg->max_eigen_DinvA[level] > 0) {
7947e6512fdSJed Brown             // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it.
7957e6512fdSJed Brown             // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.)
79618c3aa7eSMark             PetscReal emax,emin;
7970ed2132dSStefano Zampini 
79818c3aa7eSMark             emin = mg->min_eigen_DinvA[level];
79918c3aa7eSMark             emax = mg->max_eigen_DinvA[level];
80063a3b9bcSJacob Faibussowitsch             PetscCall(PetscInfo(pc,"%s: PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %" PetscInt_FMT " (N=%" PetscInt_FMT ") with emax = %g emin = %g\n",((PetscObject)pc)->prefix,level,Aarr[level]->rmap->N,(double)emax,(double)emin));
8010a94a983SJed Brown             cheb->emin_provided = emin;
8020a94a983SJed Brown             cheb->emax_provided = emax;
80318c3aa7eSMark           }
80418c3aa7eSMark         }
80518c3aa7eSMark       }
80618c3aa7eSMark     }
8070ed2132dSStefano Zampini 
8089566063dSJacob Faibussowitsch     PetscCall(PCSetUp_MG(pc));
8090ed2132dSStefano Zampini 
810d5280255SMark F. Adams     /* clean up */
811d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
8129566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Parr[level]));
8139566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Aarr[level]));
8145b89ad90SMark F. Adams     }
815806fa848SBarry Smith   } else {
8165f8cf99dSMark F. Adams     KSP smoother;
8170ed2132dSStefano Zampini 
81863a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"%s: One level solver used (system is seen as DD). Using default solver.\n",((PetscObject)pc)->prefix));
8199566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmoother(pc, 0, &smoother));
8209566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0]));
8219566063dSJacob Faibussowitsch     PetscCall(KSPSetType(smoother, KSPPREONLY));
8229566063dSJacob Faibussowitsch     PetscCall(PCSetUp_MG(pc));
8235f8cf99dSMark F. Adams   }
824849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0));
8255b89ad90SMark F. Adams   PetscFunctionReturn(0);
8265b89ad90SMark F. Adams }
8275b89ad90SMark F. Adams 
828eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8295b89ad90SMark F. Adams /*
8305b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8315b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8325b89ad90SMark F. Adams 
8335b89ad90SMark F. Adams    Input Parameter:
8345b89ad90SMark F. Adams .  pc - the preconditioner context
8355b89ad90SMark F. Adams 
8365b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8375b89ad90SMark F. Adams */
8385b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8395b89ad90SMark F. Adams {
8405b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
8415b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
8425b89ad90SMark F. Adams 
8435b89ad90SMark F. Adams   PetscFunctionBegin;
8449566063dSJacob Faibussowitsch   PetscCall(PCReset_GAMG(pc));
8459b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
8469566063dSJacob Faibussowitsch     PetscCall((*pc_gamg->ops->destroy)(pc));
8479b8ffb57SJed Brown   }
8489566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->ops));
8499566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->gamg_type_name));
8509566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg));
8519566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
8525b89ad90SMark F. Adams   PetscFunctionReturn(0);
8535b89ad90SMark F. Adams }
8545b89ad90SMark F. Adams 
855676e1743SMark F. Adams /*@
856cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
857676e1743SMark F. Adams 
8581cc46a46SBarry Smith    Logically Collective on PC
859676e1743SMark F. Adams 
860676e1743SMark F. Adams    Input Parameters:
8611cc46a46SBarry Smith +  pc - the preconditioner context
8621cc46a46SBarry Smith -  n - the number of equations
863676e1743SMark F. Adams 
864676e1743SMark F. Adams    Options Database Key:
865147403d9SBarry Smith .  -pc_gamg_process_eq_limit <limit> - set the limit
866676e1743SMark F. Adams 
86795452b02SPatrick Sanan    Notes:
86895452b02SPatrick Sanan     GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
869cab9ed1eSBarry Smith     that has degrees of freedom
870cab9ed1eSBarry Smith 
871676e1743SMark F. Adams    Level: intermediate
872676e1743SMark F. Adams 
873db781477SPatrick Sanan .seealso: `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()`
874676e1743SMark F. Adams @*/
875676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
876676e1743SMark F. Adams {
877676e1743SMark F. Adams   PetscFunctionBegin;
878676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
879cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));
880676e1743SMark F. Adams   PetscFunctionReturn(0);
881676e1743SMark F. Adams }
882676e1743SMark F. Adams 
8831e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
884676e1743SMark F. Adams {
885c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
886c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
887676e1743SMark F. Adams 
888676e1743SMark F. Adams   PetscFunctionBegin;
8899d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
890676e1743SMark F. Adams   PetscFunctionReturn(0);
891676e1743SMark F. Adams }
892676e1743SMark F. Adams 
893389730f3SMark F. Adams /*@
894cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
895389730f3SMark F. Adams 
896389730f3SMark F. Adams  Collective on PC
897389730f3SMark F. Adams 
898389730f3SMark F. Adams    Input Parameters:
8991cc46a46SBarry Smith +  pc - the preconditioner context
9001cc46a46SBarry Smith -  n - maximum number of equations to aim for
901389730f3SMark F. Adams 
902389730f3SMark F. Adams    Options Database Key:
903147403d9SBarry Smith .  -pc_gamg_coarse_eq_limit <limit> - set the limit
904389730f3SMark F. Adams 
90574329af1SBarry Smith    Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
90674329af1SBarry Smith      has less than 1000 unknowns.
90774329af1SBarry Smith 
908389730f3SMark F. Adams    Level: intermediate
909389730f3SMark F. Adams 
910db781477SPatrick Sanan .seealso: `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`
911389730f3SMark F. Adams @*/
912389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
913389730f3SMark F. Adams {
914389730f3SMark F. Adams   PetscFunctionBegin;
915389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
916cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
917389730f3SMark F. Adams   PetscFunctionReturn(0);
918389730f3SMark F. Adams }
919389730f3SMark F. Adams 
9201e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
921389730f3SMark F. Adams {
922389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
923389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
924389730f3SMark F. Adams 
925389730f3SMark F. Adams   PetscFunctionBegin;
9269d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
927389730f3SMark F. Adams   PetscFunctionReturn(0);
928389730f3SMark F. Adams }
929389730f3SMark F. Adams 
930676e1743SMark F. Adams /*@
931cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
932676e1743SMark F. Adams 
933676e1743SMark F. Adams    Collective on PC
934676e1743SMark F. Adams 
935676e1743SMark F. Adams    Input Parameters:
9361cc46a46SBarry Smith +  pc - the preconditioner context
9371cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
938676e1743SMark F. Adams 
939676e1743SMark F. Adams    Options Database Key:
940147403d9SBarry Smith .  -pc_gamg_repartition <true,false> - turn on the repartitioning
941676e1743SMark F. Adams 
94295452b02SPatrick Sanan    Notes:
94395452b02SPatrick Sanan     this will generally improve the loading balancing of the work on each level
944cab9ed1eSBarry Smith 
945676e1743SMark F. Adams    Level: intermediate
946676e1743SMark F. Adams 
947676e1743SMark F. Adams @*/
948cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
949676e1743SMark F. Adams {
950676e1743SMark F. Adams   PetscFunctionBegin;
951676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
952cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));
953676e1743SMark F. Adams   PetscFunctionReturn(0);
954676e1743SMark F. Adams }
955676e1743SMark F. Adams 
956cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
957676e1743SMark F. Adams {
958c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
959c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
960676e1743SMark F. Adams 
961676e1743SMark F. Adams   PetscFunctionBegin;
9629d5b6da9SMark F. Adams   pc_gamg->repart = n;
963676e1743SMark F. Adams   PetscFunctionReturn(0);
964676e1743SMark F. Adams }
965676e1743SMark F. Adams 
966dfd5c07aSMark F. Adams /*@
9677e6512fdSJed Brown    PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Chebyshev smoother
96818c3aa7eSMark 
96918c3aa7eSMark    Collective on PC
97018c3aa7eSMark 
97118c3aa7eSMark    Input Parameters:
97218c3aa7eSMark +  pc - the preconditioner context
97318c3aa7eSMark -  n - number of its
97418c3aa7eSMark 
97518c3aa7eSMark    Options Database Key:
976147403d9SBarry Smith .  -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate
97718c3aa7eSMark 
97818c3aa7eSMark    Notes:
9797e6512fdSJed Brown    Smoothed aggregation constructs the smoothed prolongator $P = (I - \omega D^{-1} A) T$ where $T$ is the tentative prolongator and $D$ is the diagonal of $A$.
9807e6512fdSJed Brown    Eigenvalue estimates (based on a few CG or GMRES iterations) are computed to choose $\omega$ so that this is a stable smoothing operation.
9817e6512fdSJed Brown    If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused.
982efe053fcSJed Brown    This option is only used when the smoother uses Jacobi, and should be turned off if a different PCJacobiType is used.
983efe053fcSJed Brown    It became default in PETSc 3.17.
98418c3aa7eSMark 
9857e6512fdSJed Brown    Level: advanced
98618c3aa7eSMark 
987db781477SPatrick Sanan .seealso: `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`
98818c3aa7eSMark @*/
98918c3aa7eSMark PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
99018c3aa7eSMark {
99118c3aa7eSMark   PetscFunctionBegin;
99218c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
993cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));
99418c3aa7eSMark   PetscFunctionReturn(0);
99518c3aa7eSMark }
99618c3aa7eSMark 
9970ed2132dSStefano Zampini static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
99818c3aa7eSMark {
99918c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
100018c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
100118c3aa7eSMark 
100218c3aa7eSMark   PetscFunctionBegin;
10037e6512fdSJed Brown   pc_gamg->use_sa_esteig = n;
100418c3aa7eSMark   PetscFunctionReturn(0);
100518c3aa7eSMark }
100618c3aa7eSMark 
100718c3aa7eSMark /*@
100818c3aa7eSMark    PCGAMGSetEigenvalues - Set eigenvalues
100918c3aa7eSMark 
101018c3aa7eSMark    Collective on PC
101118c3aa7eSMark 
101218c3aa7eSMark    Input Parameters:
101318c3aa7eSMark +  pc - the preconditioner context
101418c3aa7eSMark -  emax - max eigenvalue
101518c3aa7eSMark .  emin - min eigenvalue
101618c3aa7eSMark 
101718c3aa7eSMark    Options Database Key:
1018147403d9SBarry Smith .  -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues
101918c3aa7eSMark 
102018c3aa7eSMark    Level: intermediate
102118c3aa7eSMark 
1022db781477SPatrick Sanan .seealso: `PCGAMGSetUseSAEstEig()`
102318c3aa7eSMark @*/
102418c3aa7eSMark PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin)
102518c3aa7eSMark {
102618c3aa7eSMark   PetscFunctionBegin;
102718c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1028cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));
102918c3aa7eSMark   PetscFunctionReturn(0);
103018c3aa7eSMark }
103141ffd417SStefano Zampini 
103218c3aa7eSMark static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin)
103318c3aa7eSMark {
103418c3aa7eSMark   PC_MG          *mg      = (PC_MG*)pc->data;
103518c3aa7eSMark   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
103618c3aa7eSMark 
103718c3aa7eSMark   PetscFunctionBegin;
10385f80ce2aSJacob Faibussowitsch   PetscCheck(emax > emin,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin);
10395f80ce2aSJacob Faibussowitsch   PetscCheck(emax*emin > 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin);
104018c3aa7eSMark   pc_gamg->emax = emax;
104118c3aa7eSMark   pc_gamg->emin = emin;
104218c3aa7eSMark   PetscFunctionReturn(0);
104318c3aa7eSMark }
104418c3aa7eSMark 
104518c3aa7eSMark /*@
1046cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
1047dfd5c07aSMark F. Adams 
1048dfd5c07aSMark F. Adams    Collective on PC
1049dfd5c07aSMark F. Adams 
1050dfd5c07aSMark F. Adams    Input Parameters:
10511cc46a46SBarry Smith +  pc - the preconditioner context
10521cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1053dfd5c07aSMark F. Adams 
1054dfd5c07aSMark F. Adams    Options Database Key:
1055147403d9SBarry Smith .  -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation
1056dfd5c07aSMark F. Adams 
1057dfd5c07aSMark F. Adams    Level: intermediate
1058dfd5c07aSMark F. Adams 
105995452b02SPatrick Sanan    Notes:
1060147403d9SBarry Smith     May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1061cab9ed1eSBarry Smith     rebuilding the preconditioner quicker.
1062cab9ed1eSBarry Smith 
1063dfd5c07aSMark F. Adams @*/
10641cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1065dfd5c07aSMark F. Adams {
1066dfd5c07aSMark F. Adams   PetscFunctionBegin;
1067dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1068cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
1069dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1070dfd5c07aSMark F. Adams }
1071dfd5c07aSMark F. Adams 
10721cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1073dfd5c07aSMark F. Adams {
1074dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1075dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1076dfd5c07aSMark F. Adams 
1077dfd5c07aSMark F. Adams   PetscFunctionBegin;
1078dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1079dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1080dfd5c07aSMark F. Adams }
1081dfd5c07aSMark F. Adams 
1082ffc955d6SMark F. Adams /*@
1083cab9ed1eSBarry Smith    PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner.
1084ffc955d6SMark F. Adams 
1085ffc955d6SMark F. Adams    Collective on PC
1086ffc955d6SMark F. Adams 
1087ffc955d6SMark F. Adams    Input Parameters:
1088cab9ed1eSBarry Smith +  pc - the preconditioner context
1089cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
1090ffc955d6SMark F. Adams 
1091ffc955d6SMark F. Adams    Options Database Key:
1092147403d9SBarry Smith .  -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains
1093ffc955d6SMark F. Adams 
1094ffc955d6SMark F. Adams    Level: intermediate
1095ffc955d6SMark F. Adams 
1096ffc955d6SMark F. Adams @*/
1097cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1098ffc955d6SMark F. Adams {
1099ffc955d6SMark F. Adams   PetscFunctionBegin;
1100ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1101cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));
1102ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1103ffc955d6SMark F. Adams }
1104ffc955d6SMark F. Adams 
1105cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1106ffc955d6SMark F. Adams {
1107ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1108ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1109ffc955d6SMark F. Adams 
1110ffc955d6SMark F. Adams   PetscFunctionBegin;
1111cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
1112ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1113ffc955d6SMark F. Adams }
1114ffc955d6SMark F. Adams 
1115171cca9aSMark Adams /*@
1116cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
1117171cca9aSMark Adams 
1118171cca9aSMark Adams    Collective on PC
1119171cca9aSMark Adams 
1120171cca9aSMark Adams    Input Parameters:
1121171cca9aSMark Adams +  pc - the preconditioner context
1122cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
1123171cca9aSMark Adams 
1124171cca9aSMark Adams    Options Database Key:
1125147403d9SBarry Smith .  -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver
1126171cca9aSMark Adams 
1127171cca9aSMark Adams    Level: intermediate
1128171cca9aSMark Adams 
1129db781477SPatrick Sanan .seealso: `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()`
1130171cca9aSMark Adams @*/
1131171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1132171cca9aSMark Adams {
1133171cca9aSMark Adams   PetscFunctionBegin;
1134171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1135cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));
1136171cca9aSMark Adams   PetscFunctionReturn(0);
1137171cca9aSMark Adams }
1138171cca9aSMark Adams 
1139171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1140171cca9aSMark Adams {
1141171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1142171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1143171cca9aSMark Adams 
1144171cca9aSMark Adams   PetscFunctionBegin;
1145171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
1146ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1147ffc955d6SMark F. Adams }
1148ffc955d6SMark F. Adams 
11494ef23d27SMark F. Adams /*@
1150ce7c7f2fSMark Adams    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU
1151ce7c7f2fSMark Adams 
1152ce7c7f2fSMark Adams    Collective on PC
1153ce7c7f2fSMark Adams 
1154ce7c7f2fSMark Adams    Input Parameters:
1155ce7c7f2fSMark Adams +  pc - the preconditioner context
1156ce7c7f2fSMark Adams -  flg - PETSC_TRUE to pin coarse grids to CPU
1157ce7c7f2fSMark Adams 
1158ce7c7f2fSMark Adams    Options Database Key:
1159147403d9SBarry Smith .  -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU
1160ce7c7f2fSMark Adams 
1161ce7c7f2fSMark Adams    Level: intermediate
1162ce7c7f2fSMark Adams 
1163db781477SPatrick Sanan .seealso: `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetUseParallelCoarseGridSolve()`
1164ce7c7f2fSMark Adams @*/
1165ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1166ce7c7f2fSMark Adams {
1167ce7c7f2fSMark Adams   PetscFunctionBegin;
1168ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1169cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));
1170ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1171ce7c7f2fSMark Adams }
1172ce7c7f2fSMark Adams 
1173ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1174ce7c7f2fSMark Adams {
1175ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1176ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1177ce7c7f2fSMark Adams 
1178ce7c7f2fSMark Adams   PetscFunctionBegin;
1179ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = flg;
1180ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1181ce7c7f2fSMark Adams }
1182ce7c7f2fSMark Adams 
1183ce7c7f2fSMark Adams /*@
1184147403d9SBarry Smith    PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type)
1185ce7c7f2fSMark Adams 
1186ce7c7f2fSMark Adams    Collective on PC
1187ce7c7f2fSMark Adams 
1188ce7c7f2fSMark Adams    Input Parameters:
1189ce7c7f2fSMark Adams +  pc - the preconditioner context
1190ce7c7f2fSMark Adams -  flg - Layout type
1191ce7c7f2fSMark Adams 
1192ce7c7f2fSMark Adams    Options Database Key:
1193147403d9SBarry Smith .  -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering
1194ce7c7f2fSMark Adams 
1195ce7c7f2fSMark Adams    Level: intermediate
1196ce7c7f2fSMark Adams 
1197db781477SPatrick Sanan .seealso: `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()`
1198ce7c7f2fSMark Adams @*/
1199ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1200ce7c7f2fSMark Adams {
1201ce7c7f2fSMark Adams   PetscFunctionBegin;
1202ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1203cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));
1204ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1205ce7c7f2fSMark Adams }
1206ce7c7f2fSMark Adams 
1207ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1208ce7c7f2fSMark Adams {
1209ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1210ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1211ce7c7f2fSMark Adams 
1212ce7c7f2fSMark Adams   PetscFunctionBegin;
1213ce7c7f2fSMark Adams   pc_gamg->layout_type = flg;
1214ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1215ce7c7f2fSMark Adams }
1216ce7c7f2fSMark Adams 
1217ce7c7f2fSMark Adams /*@
12181cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
12194ef23d27SMark F. Adams 
12204ef23d27SMark F. Adams    Not collective on PC
12214ef23d27SMark F. Adams 
12224ef23d27SMark F. Adams    Input Parameters:
12231cc46a46SBarry Smith +  pc - the preconditioner
12241cc46a46SBarry Smith -  n - the maximum number of levels to use
12254ef23d27SMark F. Adams 
12264ef23d27SMark F. Adams    Options Database Key:
1227147403d9SBarry Smith .  -pc_mg_levels <n> - set the maximum number of levels to allow
12284ef23d27SMark F. Adams 
12294ef23d27SMark F. Adams    Level: intermediate
12304ef23d27SMark F. Adams 
12314ef23d27SMark F. Adams @*/
12324ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12334ef23d27SMark F. Adams {
12344ef23d27SMark F. Adams   PetscFunctionBegin;
12354ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1236cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
12374ef23d27SMark F. Adams   PetscFunctionReturn(0);
12384ef23d27SMark F. Adams }
12394ef23d27SMark F. Adams 
12401e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12414ef23d27SMark F. Adams {
12424ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
12434ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12444ef23d27SMark F. Adams 
12454ef23d27SMark F. Adams   PetscFunctionBegin;
12469d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12474ef23d27SMark F. Adams   PetscFunctionReturn(0);
12484ef23d27SMark F. Adams }
12494ef23d27SMark F. Adams 
12503542efc5SMark F. Adams /*@
12513542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12523542efc5SMark F. Adams 
12533542efc5SMark F. Adams    Not collective on PC
12543542efc5SMark F. Adams 
12553542efc5SMark F. Adams    Input Parameters:
12561cc46a46SBarry Smith +  pc - the preconditioner context
1257c9567895SMark .  v - array of threshold values for finest n levels; 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
1258055c8bd0SJed Brown -  n - number of threshold values provided in array
12593542efc5SMark F. Adams 
12603542efc5SMark F. Adams    Options Database Key:
1261147403d9SBarry Smith .  -pc_gamg_threshold <threshold> - the threshold to drop edges
12623542efc5SMark F. Adams 
126395452b02SPatrick Sanan    Notes:
1264af3c827dSMark Adams     Increasing the threshold decreases the rate of coarsening. Conversely reducing the threshold increases the rate of coarsening (aggressive coarsening) and thereby reduces the complexity of the coarse grids, and generally results in slower solver converge rates. Reducing coarse grid complexity reduced the complexity of Galerkin coarse grid construction considerably.
1265af3c827dSMark Adams     Before coarsening or aggregating the graph, GAMG removes small values from the graph with this threshold, and thus reducing the coupling in the graph and a different (perhaps better) coarser set of points.
1266cab9ed1eSBarry Smith 
1267055c8bd0SJed Brown     If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening.
1268055c8bd0SJed Brown     In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold().
1269055c8bd0SJed Brown     If n is greater than the total number of levels, the excess entries in threshold will not be used.
1270055c8bd0SJed Brown 
12713542efc5SMark F. Adams    Level: intermediate
12723542efc5SMark F. Adams 
1273db781477SPatrick Sanan .seealso: `PCGAMGFilterGraph()`, `PCGAMGSetSquareGraph()`
12743542efc5SMark F. Adams @*/
1275c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
12763542efc5SMark F. Adams {
12773542efc5SMark F. Adams   PetscFunctionBegin;
12783542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1279055c8bd0SJed Brown   if (n) PetscValidRealPointer(v,2);
1280cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));
12813542efc5SMark F. Adams   PetscFunctionReturn(0);
12823542efc5SMark F. Adams }
12833542efc5SMark F. Adams 
1284c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
12853542efc5SMark F. Adams {
1286c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1287c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1288c1eae691SMark Adams   PetscInt i;
1289c1eae691SMark Adams   PetscFunctionBegin;
1290055c8bd0SJed Brown   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1291055c8bd0SJed Brown   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;
1292c1eae691SMark Adams   PetscFunctionReturn(0);
1293c1eae691SMark Adams }
1294c1eae691SMark Adams 
1295c1eae691SMark Adams /*@
1296147403d9SBarry Smith    PCGAMGSetRankReductionFactors - Set manual schedule for process reduction on coarse grids
1297c9567895SMark 
1298c9567895SMark    Collective on PC
1299c9567895SMark 
1300c9567895SMark    Input Parameters:
1301c9567895SMark +  pc - the preconditioner context
1302c9567895SMark .  v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda
1303c9567895SMark -  n - number of values provided in array
1304c9567895SMark 
1305c9567895SMark    Options Database Key:
1306147403d9SBarry Smith .  -pc_gamg_rank_reduction_factors <factors> - provide the schedule
1307c9567895SMark 
1308c9567895SMark    Level: intermediate
1309c9567895SMark 
1310db781477SPatrick Sanan .seealso: `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`
1311c9567895SMark @*/
1312c9567895SMark PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1313c9567895SMark {
1314c9567895SMark   PetscFunctionBegin;
1315c9567895SMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1316c9567895SMark   if (n) PetscValidIntPointer(v,2);
1317cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n));
1318c9567895SMark   PetscFunctionReturn(0);
1319c9567895SMark }
1320c9567895SMark 
1321c9567895SMark static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1322c9567895SMark {
1323c9567895SMark   PC_MG   *mg      = (PC_MG*)pc->data;
1324c9567895SMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1325c9567895SMark   PetscInt i;
1326c9567895SMark   PetscFunctionBegin;
1327c9567895SMark   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1328c9567895SMark   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
1329c9567895SMark   PetscFunctionReturn(0);
1330c9567895SMark }
1331c9567895SMark 
1332c9567895SMark /*@
1333c1eae691SMark Adams    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1334c1eae691SMark Adams 
1335c1eae691SMark Adams    Not collective on PC
1336c1eae691SMark Adams 
1337c1eae691SMark Adams    Input Parameters:
1338c1eae691SMark Adams +  pc - the preconditioner context
1339147403d9SBarry Smith -  scale - the threshold value reduction, usually < 1.0
1340c1eae691SMark Adams 
1341c1eae691SMark Adams    Options Database Key:
1342147403d9SBarry Smith .  -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level
1343c1eae691SMark Adams 
1344055c8bd0SJed Brown    Notes:
1345055c8bd0SJed Brown    The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold().
1346055c8bd0SJed Brown    This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold().
1347055c8bd0SJed Brown 
1348c1eae691SMark Adams    Level: advanced
1349c1eae691SMark Adams 
1350db781477SPatrick Sanan .seealso: `PCGAMGSetThreshold()`
1351c1eae691SMark Adams @*/
1352c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1353c1eae691SMark Adams {
13543542efc5SMark F. Adams   PetscFunctionBegin;
1355c1eae691SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1356cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));
1357c1eae691SMark Adams   PetscFunctionReturn(0);
1358c1eae691SMark Adams }
1359c1eae691SMark Adams 
1360c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1361c1eae691SMark Adams {
1362c1eae691SMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1363c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1364c1eae691SMark Adams   PetscFunctionBegin;
1365c1eae691SMark Adams   pc_gamg->threshold_scale = v;
13663542efc5SMark F. Adams   PetscFunctionReturn(0);
13673542efc5SMark F. Adams }
13683542efc5SMark F. Adams 
1369e20c40e8SBarry Smith /*@C
1370c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1371676e1743SMark F. Adams 
1372676e1743SMark F. Adams    Collective on PC
1373676e1743SMark F. Adams 
1374676e1743SMark F. Adams    Input Parameters:
1375c60c7ad4SBarry Smith +  pc - the preconditioner context
1376c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1377676e1743SMark F. Adams 
1378676e1743SMark F. Adams    Options Database Key:
1379cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1380676e1743SMark F. Adams 
1381676e1743SMark F. Adams    Level: intermediate
1382676e1743SMark F. Adams 
1383db781477SPatrick Sanan .seealso: `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType`
1384676e1743SMark F. Adams @*/
138519fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1386676e1743SMark F. Adams {
1387676e1743SMark F. Adams   PetscFunctionBegin;
1388676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1389cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1390676e1743SMark F. Adams   PetscFunctionReturn(0);
1391676e1743SMark F. Adams }
1392676e1743SMark F. Adams 
1393e20c40e8SBarry Smith /*@C
1394c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1395c60c7ad4SBarry Smith 
1396c60c7ad4SBarry Smith    Collective on PC
1397c60c7ad4SBarry Smith 
1398c60c7ad4SBarry Smith    Input Parameter:
1399c60c7ad4SBarry Smith .  pc - the preconditioner context
1400c60c7ad4SBarry Smith 
1401c60c7ad4SBarry Smith    Output Parameter:
1402c60c7ad4SBarry Smith .  type - the type of algorithm used
1403c60c7ad4SBarry Smith 
1404c60c7ad4SBarry Smith    Level: intermediate
1405c60c7ad4SBarry Smith 
1406db781477SPatrick Sanan .seealso: `PCGAMGSetType()`, `PCGAMGType`
1407c60c7ad4SBarry Smith @*/
1408c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1409c60c7ad4SBarry Smith {
1410c60c7ad4SBarry Smith   PetscFunctionBegin;
1411c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1412cac4c232SBarry Smith   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1413c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1414c60c7ad4SBarry Smith }
1415c60c7ad4SBarry Smith 
1416c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1417c60c7ad4SBarry Smith {
1418c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1419c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1420c60c7ad4SBarry Smith 
1421c60c7ad4SBarry Smith   PetscFunctionBegin;
1422c60c7ad4SBarry Smith   *type = pc_gamg->type;
1423c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1424c60c7ad4SBarry Smith }
1425c60c7ad4SBarry Smith 
14261e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1427676e1743SMark F. Adams {
14281ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
14291ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
14305f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(PC);
1431676e1743SMark F. Adams 
1432676e1743SMark F. Adams   PetscFunctionBegin;
1433c60c7ad4SBarry Smith   pc_gamg->type = type;
14349566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(GAMGList,type,&r));
14355f80ce2aSJacob Faibussowitsch   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
14361ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
14379566063dSJacob Faibussowitsch     PetscCall((*pc_gamg->ops->destroy)(pc));
14389566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps)));
1439e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
14403ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
14413ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
14423ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
14433ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
14443ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
14459566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
14463ae0bb68SMark Adams     pc_gamg->data_sz = 0;
14471ab5ffc9SJed Brown   }
14489566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->gamg_type_name));
14499566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(type,&pc_gamg->gamg_type_name));
14509566063dSJacob Faibussowitsch   PetscCall((*r)(pc));
1451676e1743SMark F. Adams   PetscFunctionReturn(0);
1452676e1743SMark F. Adams }
1453676e1743SMark F. Adams 
14545adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
14555adeb434SBarry Smith {
14565adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
14575adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1458e7d4b4cbSMark Adams   PetscReal       gc=0, oc=0;
145990db8557SMark Adams 
14605adeb434SBarry Smith   PetscFunctionBegin;
14619566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n"));
14629566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level ="));
14639566063dSJacob Faibussowitsch   for (PetscInt i=0;i<mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]));
14649566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
14659566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale));
1466cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
14679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n"));
1468cab9ed1eSBarry Smith   }
1469171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
14709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n"));
1471171cca9aSMark Adams   }
1472ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1473ce7c7f2fSMark Adams   if (pc_gamg->cpu_pin_coarse_grids) {
14749566063dSJacob Faibussowitsch     /* PetscCall(PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n")); */
1475ce7c7f2fSMark Adams   }
1476ce7c7f2fSMark Adams #endif
1477ce7c7f2fSMark Adams   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
14789566063dSJacob Faibussowitsch   /*   PetscCall(PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n")); */
1479ce7c7f2fSMark Adams   /* } else { */
14809566063dSJacob Faibussowitsch   /*   PetscCall(PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n")); */
1481ce7c7f2fSMark Adams   /* } */
14825adeb434SBarry Smith   if (pc_gamg->ops->view) {
14839566063dSJacob Faibussowitsch     PetscCall((*pc_gamg->ops->view)(pc,viewer));
14845adeb434SBarry Smith   }
14859566063dSJacob Faibussowitsch   PetscCall(PCMGGetGridComplexity(pc,&gc,&oc));
148663a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g    operator = %g\n",(double)gc,(double)oc));
14875adeb434SBarry Smith   PetscFunctionReturn(0);
14885adeb434SBarry Smith }
14895adeb434SBarry Smith 
14904416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
14915b89ad90SMark F. Adams {
1492676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1493676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
14947e6512fdSJed Brown   PetscBool      flag;
14953b4367a7SBarry Smith   MPI_Comm       comm;
149618c3aa7eSMark   char           prefix[256],tname[32];
1497c1eae691SMark Adams   PetscInt       i,n;
149814a9496bSBarry Smith   const char     *pcpre;
14990a545947SLisandro Dalcin   static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL};
15005b89ad90SMark F. Adams   PetscFunctionBegin;
15019566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
1502d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"GAMG options");
15039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag));
1504bd94a7aaSJed Brown   if (flag) {
15059566063dSJacob Faibussowitsch     PetscCall(PCGAMGSetType(pc,tname));
15061ab5ffc9SJed Brown   }
15079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL));
15089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",pc_gamg->use_sa_esteig,&pc_gamg->use_sa_esteig,NULL));
15099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL));
15109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation aggregates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL));
15119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_use_parallel_coarse_grid_solver","Use parallel coarse grid solver (otherwise put last grid on one process)","PCGAMGSetUseParallelCoarseGridSolve",pc_gamg->use_parallel_coarse_grid_solver,&pc_gamg->use_parallel_coarse_grid_solver,NULL));
15129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_cpu_pin_coarse_grids","Pin coarse grids to the CPU","PCGAMGSetCpuPinCoarseGrids",pc_gamg->cpu_pin_coarse_grids,&pc_gamg->cpu_pin_coarse_grids,NULL));
15139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_gamg_coarse_grid_layout_type","compact: place reduced grids on processes in natural order; spread: distribute to whole machine for more memory bandwidth","PCGAMGSetCoarseGridLayoutType",LayoutTypes,(PetscEnum)pc_gamg->layout_type,(PetscEnum*)&pc_gamg->layout_type,NULL));
15149566063dSJacob Faibussowitsch   PetscCall(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));
15159566063dSJacob Faibussowitsch   PetscCall(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));
15169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL));
151718c3aa7eSMark   n = PETSC_MG_MAXLEVELS;
15189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag));
151918c3aa7eSMark   if (!flag || n < PETSC_MG_MAXLEVELS) {
1520efd3c5ceSMark Adams     if (!flag) n = 1;
1521c1eae691SMark Adams     i = n;
152218c3aa7eSMark     do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS);
1523c1eae691SMark Adams   }
1524c9567895SMark   n = PETSC_MG_MAXLEVELS;
15259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-pc_gamg_rank_reduction_factors","Manual schedule of coarse grid reduction factors that overrides internal heuristics (0 for first reduction puts one process/device)","PCGAMGSetRankReductionFactors",pc_gamg->level_reduction_factors,&n,&flag));
1526c9567895SMark   if (!flag) i = 0;
1527c9567895SMark   else i = n;
1528c9567895SMark   do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS);
15299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL));
153018c3aa7eSMark   {
153118c3aa7eSMark     PetscReal eminmax[2] = {0., 0.};
153218c3aa7eSMark     n = 2;
15339566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag));
153418c3aa7eSMark     if (flag) {
153508401ef6SPierre Jolivet       PetscCheck(n == 2,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
15369566063dSJacob Faibussowitsch       PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]));
153718c3aa7eSMark     }
153818c3aa7eSMark   }
1539b7cbab4eSMark Adams   /* set options for subtype */
15409566063dSJacob Faibussowitsch   if (pc_gamg->ops->setfromoptions) PetscCall((*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc));
154118c3aa7eSMark 
15429566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &pcpre));
15439566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : ""));
1544d0609cedSBarry Smith   PetscOptionsHeadEnd();
15455b89ad90SMark F. Adams   PetscFunctionReturn(0);
15465b89ad90SMark F. Adams }
15475b89ad90SMark F. Adams 
15485b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
15495b89ad90SMark F. Adams /*MC
15501cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
15515b89ad90SMark F. Adams 
1552280d9858SJed Brown    Options Database Keys:
1553cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1554cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1555cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1556cab9ed1eSBarry Smith .   -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother
1557cab9ed1eSBarry Smith .   -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit>
1558cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1559cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
15606008e27bSRichard Tran Mills .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1561c1eae691SMark Adams -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1562cab9ed1eSBarry Smith 
1563cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1564cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1565cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1566cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1567cab9ed1eSBarry Smith 
1568db9745e2SBarry Smith    Multigrid options:
1569db9745e2SBarry Smith +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1570db9745e2SBarry Smith .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1571db9745e2SBarry Smith .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1572cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
15735b89ad90SMark F. Adams 
157495452b02SPatrick Sanan   Notes:
157595452b02SPatrick Sanan     In order to obtain good performance for PCGAMG for vector valued problems you must
1576db9745e2SBarry Smith        Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1577db9745e2SBarry Smith        Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1578db9745e2SBarry Smith        See the Users Manual Chapter 4 for more details
15791cc46a46SBarry Smith 
15805b89ad90SMark F. Adams   Level: intermediate
1581280d9858SJed Brown 
1582db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`,
1583db781477SPatrick Sanan           `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGSetUseSAEstEig()`
15845b89ad90SMark F. Adams M*/
1585b2573a8aSBarry Smith 
15868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
15875b89ad90SMark F. Adams {
15885b89ad90SMark F. Adams   PC_GAMG *pc_gamg;
15895b89ad90SMark F. Adams   PC_MG   *mg;
15905b89ad90SMark F. Adams 
15915b89ad90SMark F. Adams   PetscFunctionBegin;
15921c1aac46SBarry Smith    /* register AMG type */
15939566063dSJacob Faibussowitsch   PetscCall(PCGAMGInitializePackage());
15941c1aac46SBarry Smith 
15955b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
15969566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCMG));
15979566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG));
15985b89ad90SMark F. Adams 
15995b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
16009566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&pc_gamg));
16019566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL));
16025b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
16035b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
16045b89ad90SMark F. Adams 
16059566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&pc_gamg->ops));
16061ab5ffc9SJed Brown 
16079d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
16089d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
16090a545947SLisandro Dalcin   pc_gamg->data    = NULL;
16105b89ad90SMark F. Adams 
16119d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
16125b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
16135b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
16145b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
16155b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
16165adeb434SBarry Smith   mg->view                = PCView_GAMG;
16175b89ad90SMark F. Adams 
16189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG));
16199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG));
16209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG));
16219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG));
16229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG));
16239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG));
16249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG));
16259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG));
16269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG));
16279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG));
16289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG));
16299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG));
16309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG));
16319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG));
16329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG));
16339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG));
16349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG));
16359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG));
16369d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1637d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
16380c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1639171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1640a0095786SMark   pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1641a0095786SMark   pc_gamg->layout_type      = PCGAMG_LAYOUT_SPREAD;
1642038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
164325a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
16449566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pc_gamg->threshold,PETSC_MG_MAXLEVELS));
1645c1eae691SMark Adams   pc_gamg->threshold_scale = 1.;
164618c3aa7eSMark   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
16479ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
16487e6512fdSJed Brown   pc_gamg->use_sa_esteig    = PETSC_TRUE;
164918c3aa7eSMark   pc_gamg->emin             = 0;
165018c3aa7eSMark   pc_gamg->emax             = 0;
165118c3aa7eSMark 
1652c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
16539d5b6da9SMark F. Adams 
1654bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
16559566063dSJacob Faibussowitsch   PetscCall(PCGAMGSetType(pc,PCGAMGAGG));
16565b89ad90SMark F. Adams   PetscFunctionReturn(0);
16575b89ad90SMark F. Adams }
16583e3471ccSMark Adams 
16593e3471ccSMark Adams /*@C
16603e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
16618a690491SBarry Smith     from PCInitializePackage().
16623e3471ccSMark Adams 
16633e3471ccSMark Adams  Level: developer
16643e3471ccSMark Adams 
1665db781477SPatrick Sanan  .seealso: `PetscInitialize()`
16663e3471ccSMark Adams @*/
16673e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
16683e3471ccSMark Adams {
16694555aa8cSStefano Zampini   PetscInt       l;
16703e3471ccSMark Adams 
16713e3471ccSMark Adams   PetscFunctionBegin;
16723e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
16733e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
16749566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO));
16759566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG));
16769566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical));
16779566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage));
1678c1c463dbSMark Adams 
1679c1c463dbSMark Adams   /* general events */
1680849bee69SMark Adams   PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP]));
1681849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH]));
1682849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGFilter", PC_CLASSID, &petsc_gamg_setup_events[GAMG_FILTER]));
1683849bee69SMark Adams   PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN]));
1684849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG Square G", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SQUARE]));
1685849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS]));
1686849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL]));
1687849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA]));
1688849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB]));
1689849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT]));
1690849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM]));
1691849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL]));
1692849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG PtAP",PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP]));
1693849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG Reduce",PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE]));
1694849bee69SMark Adams   PetscCall(PetscLogEventRegister("   GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART]));
1695849bee69SMark Adams   /* PetscCall(PetscLogEventRegister("   GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */
1696849bee69SMark Adams   /* PetscCall(PetscLogEventRegister("   GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */
1697849bee69SMark Adams   /* PetscCall(PetscLogEventRegister("   GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */
16984555aa8cSStefano Zampini   for (l=0;l<PETSC_MG_MAXLEVELS;l++) {
16994555aa8cSStefano Zampini     char ename[32];
17005b89ad90SMark F. Adams 
170163a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02" PetscInt_FMT,l));
17029566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]));
170363a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02" PetscInt_FMT,l));
17049566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]));
170563a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02" PetscInt_FMT,l));
17069566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]));
17074555aa8cSStefano Zampini   }
17084555aa8cSStefano Zampini #if defined(GAMG_STAGES)
1709849bee69SMark Adams   {  /* create timer stages */
17105b89ad90SMark F. Adams     char     str[32];
1711849bee69SMark Adams     sprintf(str,"GAMG Level %d",0);
17129566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister(str, &gamg_stages[0]));
17135b89ad90SMark F. Adams   }
17145b89ad90SMark F. Adams #endif
17153e3471ccSMark Adams   PetscFunctionReturn(0);
17163e3471ccSMark Adams }
17173e3471ccSMark Adams 
17183e3471ccSMark Adams /*@C
17191c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
17201c1aac46SBarry Smith     called from PetscFinalize() automatically.
17213e3471ccSMark Adams 
17223e3471ccSMark Adams  Level: developer
17233e3471ccSMark Adams 
1724db781477SPatrick Sanan  .seealso: `PetscFinalize()`
17253e3471ccSMark Adams @*/
17263e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
17273e3471ccSMark Adams {
17283e3471ccSMark Adams   PetscFunctionBegin;
17293e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
17309566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&GAMGList));
17313e3471ccSMark Adams   PetscFunctionReturn(0);
17323e3471ccSMark Adams }
1733a36cf38bSToby Isaac 
1734a36cf38bSToby Isaac /*@C
1735a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1736a36cf38bSToby Isaac 
1737a36cf38bSToby Isaac  Input Parameters:
1738a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1739a36cf38bSToby Isaac  - create - function for creating the gamg context.
1740a36cf38bSToby Isaac 
1741a36cf38bSToby Isaac   Level: advanced
1742a36cf38bSToby Isaac 
1743db781477SPatrick Sanan  .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
1744a36cf38bSToby Isaac @*/
1745a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1746a36cf38bSToby Isaac {
1747a36cf38bSToby Isaac   PetscFunctionBegin;
17489566063dSJacob Faibussowitsch   PetscCall(PCGAMGInitializePackage());
17499566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList,type,create));
1750a36cf38bSToby Isaac   PetscFunctionReturn(0);
1751a36cf38bSToby Isaac }
1752