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 27d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset_GAMG(PC pc) 28d71ae5a4SJacob Faibussowitsch { 29d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 30d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 31d3d6bff4SMark F. Adams 32d3d6bff4SMark F. Adams PetscFunctionBegin; 339566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 341c1aac46SBarry Smith pc_gamg->data_sz = 0; 359566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->orig_data)); 365f80ce2aSJacob Faibussowitsch for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS; level++) { 3718c3aa7eSMark mg->min_eigen_DinvA[level] = 0; 3818c3aa7eSMark mg->max_eigen_DinvA[level] = 0; 3918c3aa7eSMark } 4018c3aa7eSMark pc_gamg->emin = 0; 4118c3aa7eSMark pc_gamg->emax = 0; 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43a2f3521dSMark F. Adams } 44a2f3521dSMark F. Adams 455b89ad90SMark F. Adams /* 46c238b0ebSToby Isaac PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 47a147abb0SMark F. Adams of active processors. 485b89ad90SMark F. Adams 495b89ad90SMark F. Adams Input Parameter: 50a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 51a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 529d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 53c5bfad50SMark F. Adams . cr_bs - coarse block size 543530afc2SMark F. Adams In/Output Parameter: 55a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 56afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 5711e60469SMark F. Adams Output Parameter: 583530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 595b89ad90SMark F. Adams */ 60d71ae5a4SJacob Faibussowitsch 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) 61d71ae5a4SJacob Faibussowitsch { 629d5b6da9SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 63486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 64a2f3521dSMark F. Adams Mat Cmat, Pold = *a_P_inout; 653b4367a7SBarry Smith MPI_Comm comm; 66c5df96a5SBarry Smith PetscMPIInt rank, size, new_size, nactive = *a_nactive_proc; 673ae0bb68SMark Adams PetscInt ncrs_eq, ncrs, f_bs; 685b89ad90SMark F. Adams 695b89ad90SMark F. Adams PetscFunctionBegin; 709566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat_fine, &comm)); 719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 739566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat_fine, &f_bs)); 74849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 759566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 769566063dSJacob Faibussowitsch PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat)); 779566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 78849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 79038e3b61SMark F. Adams 80ce7c7f2fSMark Adams if (Pcolumnperm) *Pcolumnperm = NULL; 81ce7c7f2fSMark Adams 823ae0bb68SMark Adams /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 839566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL)); 843ae0bb68SMark Adams if (pc_gamg->data_cell_rows > 0) { 853ae0bb68SMark Adams ncrs = pc_gamg->data_sz / pc_gamg->data_cell_cols / pc_gamg->data_cell_rows; 8673911c69SBarry Smith } else { 873ae0bb68SMark Adams PetscInt bs; 889566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Cmat, &bs)); 893ae0bb68SMark Adams ncrs = ncrs_eq / bs; 903ae0bb68SMark Adams } 91c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 92c9567895SMark 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. */ 93c9567895SMark #if defined(PETSC_HAVE_CUDA) 94c9567895SMark PetscShmComm pshmcomm; 95c9567895SMark PetscMPIInt locrank; 96c9567895SMark MPI_Comm loccomm; 97c9567895SMark PetscInt s_nnodes, r_nnodes, new_new_size; 98c9567895SMark cudaError_t cerr; 99c9567895SMark int devCount; 1009566063dSJacob Faibussowitsch PetscCall(PetscShmCommGet(comm, &pshmcomm)); 1019566063dSJacob Faibussowitsch PetscCall(PetscShmCommGetMpiShmComm(pshmcomm, &loccomm)); 1029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(loccomm, &locrank)); 103c9567895SMark s_nnodes = !locrank; 104712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&s_nnodes, &r_nnodes, 1, MPIU_INT, MPI_SUM, comm)); 10563a3b9bcSJacob Faibussowitsch PetscCheck((size % r_nnodes) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "odd number of nodes np=%d nnodes%" PetscInt_FMT, size, r_nnodes); 106c9567895SMark devCount = 0; 107c9567895SMark cerr = cudaGetDeviceCount(&devCount); 108c9567895SMark cudaGetLastError(); /* Reset the last error */ 109c9567895SMark if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */ 110c9567895SMark new_new_size = r_nnodes * devCount; 111c9567895SMark new_size = new_new_size; 11263a3b9bcSJacob 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)); 113c9567895SMark } else { 11463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: With Cuda but no device. Use heuristics.", ((PetscObject)pc)->prefix)); 115c9567895SMark goto HEURISTIC; 116c9567895SMark } 117c9567895SMark #else 118c9567895SMark SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "should not be here"); 119c9567895SMark #endif 120c9567895SMark } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) { 12163a3b9bcSJacob 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]); 122c9567895SMark new_size = nactive / pc_gamg->level_reduction_factors[pc_gamg->current_level]; 12363a3b9bcSJacob 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])); 124c9567895SMark } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) { 125c9567895SMark new_size = 1; 1269566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Force coarsest grid reduction to %d active processes\n", ((PetscObject)pc)->prefix, new_size)); 127c9567895SMark } else { 128472110cdSMark F. Adams PetscInt ncrs_eq_glob; 129c9567895SMark #if defined(PETSC_HAVE_CUDA) 130c9567895SMark HEURISTIC: 131c9567895SMark #endif 1329566063dSJacob Faibussowitsch PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL)); 133a90e85d9SMark Adams new_size = (PetscMPIInt)((float)ncrs_eq_glob / (float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 134da81f932SPierre Jolivet if (!new_size) new_size = 1; /* not likely, possible? */ 135c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 1369566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Coarse grid reduction from %d to %d active processes\n", ((PetscObject)pc)->prefix, nactive, new_size)); 137a2f3521dSMark F. Adams } 1382e3501ffSMark Adams if (new_size == nactive) { 139ef3f0257SMark Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 140ce7c7f2fSMark Adams if (new_size < size) { 141ce7c7f2fSMark Adams /* odd case where multiple coarse grids are on one processor or no coarsening ... */ 1429566063dSJacob 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)); 143ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 1449566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE)); 1459566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE)); 146ce7c7f2fSMark Adams } 147ce7c7f2fSMark Adams } 148ef3f0257SMark Adams /* we know that the grid structure can be reused in MatPtAP */ 1492e3501ffSMark Adams } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */ 150192c0e8bSMark Adams PetscInt *counts, *newproc_idx, ii, jj, kk, strideNew, *tidx, ncrs_new, ncrs_eq_new, nloc_old, expand_factor = 1, rfactor = 1; 151885364a3SMark Adams IS is_eq_newproc, is_eq_num, is_eq_num_prim, new_eq_indices; 152849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 15371959b99SBarry Smith nloc_old = ncrs_eq / cr_bs; 15463a3b9bcSJacob 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); 155ce7c7f2fSMark Adams /* get new_size and rfactor */ 156ce7c7f2fSMark Adams if (pc_gamg->layout_type == PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) { 157ce7c7f2fSMark Adams /* find factor */ 158ce7c7f2fSMark Adams if (new_size == 1) rfactor = size; /* don't modify */ 159ce7c7f2fSMark Adams else { 160ce7c7f2fSMark Adams PetscReal best_fact = 0.; 161ce7c7f2fSMark Adams jj = -1; 162ce7c7f2fSMark Adams for (kk = 1; kk <= size; kk++) { 163ce7c7f2fSMark Adams if (!(size % kk)) { /* a candidate */ 164ce7c7f2fSMark Adams PetscReal nactpe = (PetscReal)size / (PetscReal)kk, fact = nactpe / (PetscReal)new_size; 165ce7c7f2fSMark Adams if (fact > 1.0) fact = 1. / fact; /* keep fact < 1 */ 166ce7c7f2fSMark Adams if (fact > best_fact) { 1679371c9d4SSatish Balay best_fact = fact; 1689371c9d4SSatish Balay jj = kk; 169ce7c7f2fSMark Adams } 170ce7c7f2fSMark Adams } 171ce7c7f2fSMark Adams } 172ce7c7f2fSMark Adams if (jj != -1) rfactor = jj; 173ce7c7f2fSMark Adams else rfactor = 1; /* a prime */ 174ce7c7f2fSMark Adams if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1; 175ce7c7f2fSMark Adams else expand_factor = rfactor; 176ce7c7f2fSMark Adams } 177ce7c7f2fSMark Adams new_size = size / rfactor; /* make new size one that is factor */ 1784cdfd227SMark if (new_size == nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */ 1794cdfd227SMark *a_Amat_crs = Cmat; 18063a3b9bcSJacob 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)); 181849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 183ce7c7f2fSMark Adams } 184ce7c7f2fSMark Adams } 185a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 1869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &counts)); 187849bee69SMark Adams if (pc_gamg->repart) { /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */ 1885a9b9e01SMark F. Adams Mat adj; 189849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0)); 19063a3b9bcSJacob 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")); 191a2f3521dSMark F. Adams /* get 'adj' */ 192c5bfad50SMark F. Adams if (cr_bs == 1) { 1939566063dSJacob Faibussowitsch PetscCall(MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj)); 194806fa848SBarry Smith } else { 195a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 196eb07cef2SMark F. Adams Mat tMat; 197a2f3521dSMark F. Adams PetscInt Istart_crs, Iend_crs, ncols, jj, Ii; 198b4fbaa2aSMark F. Adams const PetscScalar *vals; 199b4fbaa2aSMark F. Adams const PetscInt *idx; 200a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 20139d09545SMark Adams static PetscInt llev = 0; /* ugly but just used for debugging */ 202d9558ea9SBarry Smith MatType mtype; 203b4fbaa2aSMark F. Adams 2049566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(ncrs, &d_nnz, ncrs, &o_nnz)); 2059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs)); 2069566063dSJacob Faibussowitsch PetscCall(MatGetSize(Cmat, &M, &N)); 207c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 2089566063dSJacob Faibussowitsch PetscCall(MatGetRow(Cmat, Ii, &ncols, NULL, NULL)); 209c5bfad50SMark F. Adams d_nnz[jj] = ncols / cr_bs; 210c5bfad50SMark F. Adams o_nnz[jj] = ncols / cr_bs; 2119566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Cmat, Ii, &ncols, NULL, NULL)); 2123ae0bb68SMark Adams if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 2133ae0bb68SMark Adams if (o_nnz[jj] > (M / cr_bs - ncrs)) o_nnz[jj] = M / cr_bs - ncrs; 21458471d46SMark F. Adams } 2156876a03eSMark F. Adams 2169566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat_fine, &mtype)); 2179566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &tMat)); 2189566063dSJacob Faibussowitsch PetscCall(MatSetSizes(tMat, ncrs, ncrs, PETSC_DETERMINE, PETSC_DETERMINE)); 2199566063dSJacob Faibussowitsch PetscCall(MatSetType(tMat, mtype)); 2209566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(tMat, 0, d_nnz)); 2219566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(tMat, 0, d_nnz, 0, o_nnz)); 2229566063dSJacob Faibussowitsch PetscCall(PetscFree2(d_nnz, o_nnz)); 223eb07cef2SMark F. Adams 224a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 225c5bfad50SMark F. Adams PetscInt dest_row = ii / cr_bs; 2269566063dSJacob Faibussowitsch PetscCall(MatGetRow(Cmat, ii, &ncols, &idx, &vals)); 227eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 228c5bfad50SMark F. Adams PetscInt dest_col = idx[jj] / cr_bs; 229eb07cef2SMark F. Adams PetscScalar v = 1.0; 2309566063dSJacob Faibussowitsch PetscCall(MatSetValues(tMat, 1, &dest_row, 1, &dest_col, &v, ADD_VALUES)); 231eb07cef2SMark F. Adams } 2329566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Cmat, ii, &ncols, &idx, &vals)); 233eb07cef2SMark F. Adams } 2349566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tMat, MAT_FINAL_ASSEMBLY)); 2359566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tMat, MAT_FINAL_ASSEMBLY)); 236eb07cef2SMark F. Adams 237b4fbaa2aSMark F. Adams if (llev++ == -1) { 2389371c9d4SSatish Balay PetscViewer viewer; 2399371c9d4SSatish Balay char fname[32]; 24063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fname, sizeof(fname), "part_mat_%" PetscInt_FMT ".mat", llev)); 2413ba16761SJacob Faibussowitsch PetscCall(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++) { 2709371c9d4SSatish Balay for (ii = 0; ii < cr_bs; ii++, jj++) { newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ } 2715ef31b24SMark F. Adams } 2729566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(proc_is, &is_idx)); 2739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&proc_is)); 2745ef31b24SMark F. Adams } 2759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 2765a9b9e01SMark F. Adams 2779566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFree(newproc_idx)); 279849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0)); 28031cb4603SMark Adams } else { /* simple aggregation of parts -- 'is_eq_newproc' */ 281ce7c7f2fSMark Adams PetscInt targetPE; 28208401ef6SPierre Jolivet PetscCheck(new_size != nactive, PETSC_COMM_SELF, PETSC_ERR_PLIB, "new_size==nactive. Should not happen"); 28363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Number of equations (loc) %" PetscInt_FMT " with simple aggregation\n", ((PetscObject)pc)->prefix, ncrs_eq)); 284ce7c7f2fSMark Adams targetPE = (rank / rfactor) * expand_factor; 2859566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc)); 286a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 287e33ef3b1SMark F. Adams 28811e60469SMark F. Adams /* 289a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 29011e60469SMark F. Adams */ 2919566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num)); 2927700e67bSMark Adams is_eq_num_prim = is_eq_num; 29311e60469SMark F. Adams /* 294a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 29511e60469SMark F. Adams */ 2969566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(is_eq_newproc, size, counts)); 297c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 2989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_eq_newproc)); 299ce7c7f2fSMark Adams ncrs_new = ncrs_eq_new / cr_bs; 300a2f3521dSMark F. Adams 3019566063dSJacob Faibussowitsch PetscCall(PetscFree(counts)); 3026aad120cSJose 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 */ 303885364a3SMark Adams { 304885364a3SMark Adams Vec src_crd, dest_crd; 305885364a3SMark 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; 306885364a3SMark Adams VecScatter vecscat; 307885364a3SMark Adams PetscScalar *array; 308885364a3SMark Adams IS isscat; 309a2f3521dSMark F. Adams /* move data (for primal equations only) */ 31022063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 3119566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &dest_crd)); 3129566063dSJacob Faibussowitsch PetscCall(VecSetSizes(dest_crd, node_data_sz * ncrs_new, PETSC_DECIDE)); 3139566063dSJacob Faibussowitsch PetscCall(VecSetType(dest_crd, VECSTANDARD)); /* this is needed! */ 31411e60469SMark F. Adams /* 3159d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 316c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 31711e60469SMark F. Adams */ 3189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncrs * node_data_sz, &tidx)); 3199566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is_eq_num_prim, &idx)); 3203ae0bb68SMark Adams for (ii = 0, jj = 0; ii < ncrs; ii++) { 321c5bfad50SMark F. Adams PetscInt id = idx[ii * cr_bs] / cr_bs; /* get node back */ 322a2f3521dSMark F. Adams for (kk = 0; kk < node_data_sz; kk++, jj++) tidx[jj] = id * node_data_sz + kk; 32311e60469SMark F. Adams } 3249566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is_eq_num_prim, &idx)); 3259566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, node_data_sz * ncrs, tidx, PETSC_COPY_VALUES, &isscat)); 3269566063dSJacob Faibussowitsch PetscCall(PetscFree(tidx)); 32711e60469SMark F. Adams /* 32811e60469SMark F. Adams Create a vector to contain the original vertex information for each element 32911e60469SMark F. Adams */ 3309566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, node_data_sz * ncrs, &src_crd)); 3319d5b6da9SMark F. Adams for (jj = 0; jj < ndata_cols; jj++) { 3323ae0bb68SMark Adams const PetscInt stride0 = ncrs * pc_gamg->data_cell_rows; 3333ae0bb68SMark Adams for (ii = 0; ii < ncrs; ii++) { 3349d5b6da9SMark F. Adams for (kk = 0; kk < ndata_rows; kk++) { 335a2f3521dSMark F. Adams PetscInt ix = ii * ndata_rows + kk + jj * stride0, jx = ii * node_data_sz + kk * ndata_cols + jj; 336c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 3379566063dSJacob Faibussowitsch PetscCall(VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES)); 338d3d6bff4SMark F. Adams } 339038e3b61SMark F. Adams } 340eb07cef2SMark F. Adams } 3419566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(src_crd)); 3429566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(src_crd)); 34311e60469SMark F. Adams /* 34411e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 34511e60469SMark F. Adams to the correct processor 34611e60469SMark F. Adams */ 3479566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat)); 3489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isscat)); 3499566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD)); 3509566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD)); 3519566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vecscat)); 3529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&src_crd)); 35311e60469SMark F. Adams /* 35411e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 35511e60469SMark F. Adams */ 3569566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 3579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(node_data_sz * ncrs_new, &pc_gamg->data)); 3582fa5cd67SKarl Rupp 3593ae0bb68SMark Adams pc_gamg->data_sz = node_data_sz * ncrs_new; 3603ae0bb68SMark Adams strideNew = ncrs_new * ndata_rows; 3612fa5cd67SKarl Rupp 3629566063dSJacob Faibussowitsch PetscCall(VecGetArray(dest_crd, &array)); 3639d5b6da9SMark F. Adams for (jj = 0; jj < ndata_cols; jj++) { 3643ae0bb68SMark Adams for (ii = 0; ii < ncrs_new; ii++) { 3659d5b6da9SMark F. Adams for (kk = 0; kk < ndata_rows; kk++) { 366a2f3521dSMark F. Adams PetscInt ix = ii * ndata_rows + kk + jj * strideNew, jx = ii * node_data_sz + kk * ndata_cols + jj; 367c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 368d3d6bff4SMark F. Adams } 369038e3b61SMark F. Adams } 370038e3b61SMark F. Adams } 3719566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(dest_crd, &array)); 3729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dest_crd)); 373885364a3SMark Adams } 374a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 375849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0)); */ 37611e60469SMark F. Adams /* 3777dae84e0SHong Zhang Invert for MatCreateSubMatrix 37811e60469SMark F. Adams */ 3799566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices)); 3809566063dSJacob Faibussowitsch PetscCall(ISSort(new_eq_indices)); /* is this needed? */ 3819566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(new_eq_indices, cr_bs)); 3829371c9d4SSatish Balay if (is_eq_num != is_eq_num_prim) { PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */ } 3833cb8563fSToby Isaac if (Pcolumnperm) { 3849566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)new_eq_indices)); 3853cb8563fSToby Isaac *Pcolumnperm = new_eq_indices; 3863cb8563fSToby Isaac } 3879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_eq_num)); 388849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0)); */ 389849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 390849bee69SMark Adams 391a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 392a2f3521dSMark F. Adams { 393a2f3521dSMark F. Adams Mat mat; 394b94d7dedSBarry Smith PetscBool isset, isspd, isher; 39590db8557SMark Adams #if !defined(PETSC_USE_COMPLEX) 396b94d7dedSBarry Smith PetscBool issym; 397b94d7dedSBarry Smith #endif 398b94d7dedSBarry Smith 399b94d7dedSBarry Smith PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat)); 400b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(Cmat, &isset, &isspd)); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ? 401b94d7dedSBarry Smith if (isset) PetscCall(MatSetOption(mat, MAT_SPD, isspd)); 402b94d7dedSBarry Smith else { 403b94d7dedSBarry Smith PetscCall(MatIsHermitianKnown(Cmat, &isset, &isher)); 404b94d7dedSBarry Smith if (isset) PetscCall(MatSetOption(mat, MAT_HERMITIAN, isher)); 405b94d7dedSBarry Smith else { 406b94d7dedSBarry Smith #if !defined(PETSC_USE_COMPLEX) 407b94d7dedSBarry Smith PetscCall(MatIsSymmetricKnown(Cmat, &isset, &issym)); 408b94d7dedSBarry Smith if (isset) PetscCall(MatSetOption(mat, MAT_SYMMETRIC, issym)); 40990db8557SMark Adams #endif 41090db8557SMark Adams } 41190db8557SMark Adams } 412a2f3521dSMark F. Adams *a_Amat_crs = mat; 413a2f3521dSMark F. Adams } 4149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cmat)); 415a2f3521dSMark F. Adams 416849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 41711e60469SMark F. Adams /* prolongator */ 41811e60469SMark F. Adams { 41911e60469SMark F. Adams IS findices; 420a2f3521dSMark F. Adams PetscInt Istart, Iend; 421a2f3521dSMark F. Adams Mat Pnew; 42262294041SBarry Smith 4239566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend)); 424849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 4259566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, Iend - Istart, Istart, 1, &findices)); 4269566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(findices, f_bs)); 4279566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew)); 4289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&findices)); 4299566063dSJacob Faibussowitsch PetscCall(MatSetOption(Pnew, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 430c5bfad50SMark F. Adams 431849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 4329566063dSJacob Faibussowitsch PetscCall(MatDestroy(a_P_inout)); 433a2f3521dSMark F. Adams 434a2f3521dSMark F. Adams /* output - repartitioned */ 435a2f3521dSMark F. Adams *a_P_inout = Pnew; 436e33ef3b1SMark F. Adams } 4379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&new_eq_indices)); 4385b89ad90SMark F. Adams 439c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 440ce7c7f2fSMark Adams 441ce7c7f2fSMark Adams /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */ 442ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 443ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 4448bca76a6SMark Adams static PetscInt llev = 2; 44563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Pinning level %" PetscInt_FMT " to the CPU\n", ((PetscObject)pc)->prefix, llev++)); 446ce7c7f2fSMark Adams #endif 4479566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE)); 4489566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE)); 449adf5291fSStefano Zampini if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */ 450ce7c7f2fSMark Adams Mat A = *a_Amat_crs, P = *a_P_inout; 451ce7c7f2fSMark Adams PetscMPIInt size; 4529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 453ce7c7f2fSMark Adams if (size > 1) { 454ce7c7f2fSMark Adams Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data, *p = (Mat_MPIAIJ *)P->data; 4559566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(a->lvec, PETSC_TRUE)); 4569566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(p->lvec, PETSC_TRUE)); 457ce7c7f2fSMark Adams } 458ce7c7f2fSMark Adams } 459ce7c7f2fSMark Adams } 460849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 461849bee69SMark Adams } // processor reduce 4623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4635b89ad90SMark F. Adams } 4645b89ad90SMark F. Adams 465bae903cbSmarkadams4 // used in GEO 466d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat *Gmat2) 467d71ae5a4SJacob Faibussowitsch { 4684b1575e2SStefano Zampini const char *prefix; 4694b1575e2SStefano Zampini char addp[32]; 4704b1575e2SStefano Zampini PC_MG *mg = (PC_MG *)a_pc->data; 4714b1575e2SStefano Zampini PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 4724b1575e2SStefano Zampini 4734b1575e2SStefano Zampini PetscFunctionBegin; 4749566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(a_pc, &prefix)); 47563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(a_pc, "%s: Square Graph on level %" PetscInt_FMT "\n", ((PetscObject)a_pc)->prefix, pc_gamg->current_level + 1)); 4769566063dSJacob Faibussowitsch PetscCall(MatProductCreate(Gmat1, Gmat1, NULL, Gmat2)); 4779566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(*Gmat2, prefix)); 47863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_square_%" PetscInt_FMT "_", pc_gamg->current_level)); 4799566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(*Gmat2, addp)); 480b94d7dedSBarry Smith if ((*Gmat2)->structurally_symmetric == PETSC_BOOL3_TRUE) { 4819566063dSJacob Faibussowitsch PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AB)); 482b4da3a1bSStefano Zampini } else { 4839566063dSJacob Faibussowitsch PetscCall(MatSetOption(Gmat1, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 4849566063dSJacob Faibussowitsch PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AtB)); 485b4da3a1bSStefano Zampini } 4869566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(*Gmat2)); 4879566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0)); 4889566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(*Gmat2)); 4899566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0)); 4909566063dSJacob Faibussowitsch PetscCall(MatProductClear(*Gmat2)); 4914b1575e2SStefano Zampini /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */ 4924b1575e2SStefano Zampini (*Gmat2)->assembled = PETSC_TRUE; 4933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4944b1575e2SStefano Zampini } 4954b1575e2SStefano Zampini 4965b89ad90SMark F. Adams /* 4975b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4985b89ad90SMark F. Adams by setting data structures and options. 4995b89ad90SMark F. Adams 5005b89ad90SMark F. Adams Input Parameter: 5015b89ad90SMark F. Adams . pc - the preconditioner context 5025b89ad90SMark F. Adams 5035b89ad90SMark F. Adams */ 504d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp_GAMG(PC pc) 505d71ae5a4SJacob Faibussowitsch { 5069d5b6da9SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 5075b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 5082adcac29SMark F. Adams Mat Pmat = pc->pmat; 50918c3aa7eSMark PetscInt fine_level, level, level1, bs, M, N, qq, lidx, nASMBlocksArr[PETSC_MG_MAXLEVELS]; 5103b4367a7SBarry Smith MPI_Comm comm; 511c5df96a5SBarry Smith PetscMPIInt rank, size, nactivepe; 51218c3aa7eSMark Mat Aarr[PETSC_MG_MAXLEVELS], Parr[PETSC_MG_MAXLEVELS]; 51318c3aa7eSMark IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 5144279555eSSatish Balay PetscBool is_last = PETSC_FALSE; 5154279555eSSatish Balay #if defined(PETSC_USE_INFO) 516a2f3521dSMark F. Adams PetscLogDouble nnz0 = 0., nnztot = 0.; 517569f4572SMark Adams MatInfo info; 5184279555eSSatish Balay #endif 5195ef31b24SMark F. Adams 5205b89ad90SMark F. Adams PetscFunctionBegin; 5219566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 5229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 524849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 5258abdc6daSStefano Zampini if (pc->setupcalled) { 5268abdc6daSStefano Zampini if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) { 527878e152fSMark F. Adams /* reset everything */ 5289566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 529878e152fSMark F. Adams pc->setupcalled = 0; 530806fa848SBarry Smith } else { 53184d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 53203a628feSMark F. Adams /* just do Galerkin grids */ 53358471d46SMark F. Adams Mat B, dA, dB; 5349d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 5354555aa8cSStefano Zampini PetscInt gl; 53658471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 5379566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, &dA, &dB)); 53858471d46SMark F. Adams /* (re)set to get dirty flag */ 5399566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, dA, dB)); 54058471d46SMark F. Adams 5414555aa8cSStefano Zampini for (level = pc_gamg->Nlevels - 2, gl = 0; level >= 0; level--, gl++) { 5428abdc6daSStefano Zampini MatReuse reuse = MAT_INITIAL_MATRIX; 543849bee69SMark Adams #if defined(GAMG_STAGES) 544849bee69SMark Adams PetscCall(PetscLogStagePush(gamg_stages[gl])); 545849bee69SMark Adams #endif 5468abdc6daSStefano Zampini /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 5479566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[level]->smoothd, NULL, &B)); 5488abdc6daSStefano Zampini if (B->product) { 549ad540459SPierre Jolivet if (B->product->A == dB && B->product->B == mglevels[level + 1]->interpolate) reuse = MAT_REUSE_MATRIX; 5508abdc6daSStefano Zampini } 5519566063dSJacob Faibussowitsch if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A)); 5528abdc6daSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 55363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: RAP after first solve, reuse matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 5548abdc6daSStefano Zampini } else { 55563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: RAP after first solve, new matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 5568abdc6daSStefano Zampini } 5579566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0)); 5589566063dSJacob Faibussowitsch PetscCall(MatPtAP(dB, mglevels[level + 1]->interpolate, reuse, PETSC_DEFAULT, &B)); 5599566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0)); 56063b77682SMark Adams if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B; 5619566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[level]->smoothd, B, B)); 56258471d46SMark F. Adams dB = B; 563849bee69SMark Adams #if defined(GAMG_STAGES) 564849bee69SMark Adams PetscCall(PetscLogStagePop()); 565849bee69SMark Adams #endif 56658471d46SMark F. Adams } 5675f8cf99dSMark F. Adams } 568d5280255SMark F. Adams 5699566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 570849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 5713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 572eb07cef2SMark F. Adams } 573878e152fSMark F. Adams } 574f6536408SMark F. Adams 575878e152fSMark F. Adams if (!pc_gamg->data) { 576878e152fSMark F. Adams if (pc_gamg->orig_data) { 5779566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Pmat, &bs)); 5789566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Pmat, &qq, NULL)); 5792fa5cd67SKarl Rupp 580878e152fSMark F. Adams pc_gamg->data_sz = (qq / bs) * pc_gamg->orig_data_cell_rows * pc_gamg->orig_data_cell_cols; 581878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 582878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 5832fa5cd67SKarl Rupp 5849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data)); 585878e152fSMark F. Adams for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 586806fa848SBarry Smith } else { 5875f80ce2aSJacob Faibussowitsch PetscCheck(pc_gamg->ops->createdefaultdata, comm, PETSC_ERR_PLIB, "'createdefaultdata' not set(?) need to support NULL data"); 5889566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->createdefaultdata(pc, Pmat)); 5899d5b6da9SMark F. Adams } 590878e152fSMark F. Adams } 591878e152fSMark F. Adams 592878e152fSMark F. Adams /* cache original data for reuse */ 5931c1aac46SBarry Smith if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 5949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data)); 595878e152fSMark F. Adams for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 596878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 597878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 598878e152fSMark F. Adams } 599038e3b61SMark F. Adams 600302f38e8SMark F. Adams /* get basic dims */ 6019566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Pmat, &bs)); 6029566063dSJacob Faibussowitsch PetscCall(MatGetSize(Pmat, &M, &N)); 60384d3f75bSMark F. Adams 6044279555eSSatish Balay #if defined(PETSC_USE_INFO) 6059566063dSJacob Faibussowitsch PetscCall(MatGetInfo(Pmat, MAT_GLOBAL_SUM, &info)); /* global reduction */ 606569f4572SMark Adams nnz0 = info.nz_used; 607569f4572SMark Adams nnztot = info.nz_used; 6084279555eSSatish Balay #endif 609bae903cbSmarkadams4 PetscCall(PetscInfo(pc, "%s: level %d) N=%" PetscInt_FMT ", n data rows=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", np=%d\n", ((PetscObject)pc)->prefix, 0, M, pc_gamg->data_cell_rows, pc_gamg->data_cell_cols, (PetscInt)(nnz0 / (PetscReal)M + 0.5), size)); 610569f4572SMark Adams 611a2f3521dSMark F. Adams /* Get A_i and R_i */ 61262294041SBarry Smith for (level = 0, Aarr[0] = Pmat, nactivepe = size; level < (pc_gamg->Nlevels - 1) && (!level || M > pc_gamg->coarse_eq_limit); level++) { 6139ab59c8bSMark Adams pc_gamg->current_level = level; 61463a3b9bcSJacob Faibussowitsch PetscCheck(level < PETSC_MG_MAXLEVELS, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Too many levels %" PetscInt_FMT, level); 6155b89ad90SMark F. Adams level1 = level + 1; 6164555aa8cSStefano Zampini #if defined(GAMG_STAGES) 617849bee69SMark Adams if (!gamg_stages[level]) { 618849bee69SMark Adams char str[32]; 619a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", (int)level)); 620849bee69SMark Adams PetscCall(PetscLogStageRegister(str, &gamg_stages[level])); 621849bee69SMark Adams } 6229566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(gamg_stages[level])); 623b4fbaa2aSMark F. Adams #endif 624c8b0795cSMark F. Adams { /* construct prolongator */ 625725b86d8SJed Brown Mat Gmat; 6260cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 6277700e67bSMark Adams Mat Prol11; 628c8b0795cSMark F. Adams 6292d776b49SBarry Smith PetscCall(PCGAMGCreateGraph(pc, Aarr[level], &Gmat)); 6309566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists)); 6319566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->prolongator(pc, Aarr[level], Gmat, agg_lists, &Prol11)); 632c8b0795cSMark F. Adams 633a2f3521dSMark F. Adams /* could have failed to create new level */ 634a2f3521dSMark F. Adams if (Prol11) { 635f7df55f0SStefano Zampini const char *prefix; 636f7df55f0SStefano Zampini char addp[32]; 637f7df55f0SStefano Zampini 6389d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 6399566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(Prol11, NULL, &bs)); 640a2f3521dSMark F. Adams 641fd1112cbSBarry Smith if (pc_gamg->ops->optprolongator) { 642c8b0795cSMark F. Adams /* smooth */ 6439566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11)); 644c8b0795cSMark F. Adams } 645c8b0795cSMark F. Adams 6460c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 6471b18a24aSMark Adams PetscInt bs; 648849bee69SMark Adams PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // not timed directly, ugly, could remove, but good ASM method 6499566063dSJacob Faibussowitsch PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level])); 650ffc955d6SMark F. Adams } 651ffc955d6SMark F. Adams 6529566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 6539566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(Prol11, prefix)); 6549566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_prolongator_%d_", (int)level)); 6559566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(Prol11, addp)); 65691f31d3dSStefano Zampini /* Always generate the transpose with CUDA 657f7df55f0SStefano Zampini Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */ 6589566063dSJacob Faibussowitsch PetscCall(MatSetOption(Prol11, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 6599566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Prol11)); 6604bde40a0SMark Adams Parr[level1] = Prol11; 6614bde40a0SMark Adams } else Parr[level1] = NULL; /* failed to coarsen */ 6624bde40a0SMark Adams 6639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Gmat)); 6649566063dSJacob Faibussowitsch PetscCall(PetscCDDestroy(agg_lists)); 665a2f3521dSMark F. Adams } /* construct prolongator scope */ 6667f66b68fSMark Adams if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 667171cca9aSMark Adams if (!Parr[level1]) { /* failed to coarsen */ 66863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Stop gridding, level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 6694555aa8cSStefano Zampini #if defined(GAMG_STAGES) 6709566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 671a90e85d9SMark Adams #endif 672c8b0795cSMark F. Adams break; 673c8b0795cSMark F. Adams } 6749566063dSJacob Faibussowitsch PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */ 6752472a847SBarry Smith PetscCheck(!is_last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Is last ?"); 676171cca9aSMark Adams if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 6770e2909e1SMark Adams if (level1 == pc_gamg->Nlevels - 1) is_last = PETSC_TRUE; 678849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0)); 6799566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last)); 680849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0)); 681a2f3521dSMark F. Adams 6829566063dSJacob Faibussowitsch PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */ 6834279555eSSatish Balay #if defined(PETSC_USE_INFO) 6849566063dSJacob Faibussowitsch PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info)); 685569f4572SMark Adams nnztot += info.nz_used; 6864279555eSSatish Balay #endif 687bae903cbSmarkadams4 PetscCall(PetscInfo(pc, "%s: %d) N=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", %d active pes\n", ((PetscObject)pc)->prefix, (int)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 701ba4ceb06Smarkadams4 PetscCall(PetscInfo(pc, "%s: %" PetscInt_FMT " levels, operator 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; 7412eab5cd7Smarkadams4 PetscCall(PCASMSetLocalSubdomains(subpc, sz, iss, NULL)); 74248a46eb9SPierre Jolivet for (kk = 0; kk < sz; kk++) PetscCall(ISDestroy(&iss[kk])); 7439566063dSJacob Faibussowitsch PetscCall(PetscFree(iss)); 744ffc955d6SMark F. Adams } 7450298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 746ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 747806fa848SBarry Smith } else { 7489566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCJACOBI)); 749ffc955d6SMark F. Adams } 750d5280255SMark F. Adams } 751d5280255SMark F. Adams { 752d5280255SMark F. Adams /* coarse grid */ 7539371c9d4SSatish Balay KSP smoother, *k2; 7549371c9d4SSatish Balay PC subpc, pc2; 7559371c9d4SSatish Balay PetscInt ii, first; 7569371c9d4SSatish Balay Mat Lmat = Aarr[(level = pc_gamg->Nlevels - 1)]; 7579371c9d4SSatish Balay lidx = 0; 7580ed2132dSStefano Zampini 7599566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7609566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Lmat, Lmat)); 761cf8ae1d3SMark Adams if (!pc_gamg->use_parallel_coarse_grid_solver) { 7629566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 7639566063dSJacob Faibussowitsch PetscCall(KSPGetPC(smoother, &subpc)); 7649566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCBJACOBI)); 7659566063dSJacob Faibussowitsch PetscCall(PCSetUp(subpc)); 7669566063dSJacob Faibussowitsch PetscCall(PCBJacobiGetSubKSP(subpc, &ii, &first, &k2)); 76763a3b9bcSJacob Faibussowitsch PetscCheck(ii == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " is not one", ii); 7689566063dSJacob Faibussowitsch PetscCall(KSPGetPC(k2[0], &pc2)); 7699566063dSJacob Faibussowitsch PetscCall(PCSetType(pc2, PCLU)); 7709566063dSJacob Faibussowitsch PetscCall(PCFactorSetShiftType(pc2, MAT_SHIFT_INBLOCKS)); 7719566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(k2[0], PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1)); 7729566063dSJacob Faibussowitsch PetscCall(KSPSetType(k2[0], KSPPREONLY)); 773d5280255SMark F. Adams } 774cf8ae1d3SMark Adams } 775d5280255SMark F. Adams 776d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 777d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)pc); 778dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); 779d0609cedSBarry Smith PetscOptionsEnd(); 7809566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 781d5280255SMark F. Adams 782f1580f4eSBarry Smith /* set cheby eigen estimates from SA to use in the solver */ 7837e6512fdSJed Brown if (pc_gamg->use_sa_esteig) { 78418c3aa7eSMark for (lidx = 1, level = pc_gamg->Nlevels - 2; level >= 0; lidx++, level--) { 78518c3aa7eSMark KSP smoother; 78618c3aa7eSMark PetscBool ischeb; 7870ed2132dSStefano Zampini 7889566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7899566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb)); 79018c3aa7eSMark if (ischeb) { 79118c3aa7eSMark KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data; 7920ed2132dSStefano Zampini 7932de708cbSJed Brown // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG 794efe053fcSJed Brown if (mg->max_eigen_DinvA[level] > 0) { 7957e6512fdSJed Brown // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it. 7967e6512fdSJed Brown // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.) 79718c3aa7eSMark PetscReal emax, emin; 7980ed2132dSStefano Zampini 79918c3aa7eSMark emin = mg->min_eigen_DinvA[level]; 80018c3aa7eSMark emax = mg->max_eigen_DinvA[level]; 80163a3b9bcSJacob 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)); 8020a94a983SJed Brown cheb->emin_provided = emin; 8030a94a983SJed Brown cheb->emax_provided = emax; 80418c3aa7eSMark } 80518c3aa7eSMark } 80618c3aa7eSMark } 80718c3aa7eSMark } 8080ed2132dSStefano Zampini 8099566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 8100ed2132dSStefano Zampini 811d5280255SMark F. Adams /* clean up */ 812d5280255SMark F. Adams for (level = 1; level < pc_gamg->Nlevels; level++) { 8139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Parr[level])); 8149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aarr[level])); 8155b89ad90SMark F. Adams } 816806fa848SBarry Smith } else { 8175f8cf99dSMark F. Adams KSP smoother; 8180ed2132dSStefano Zampini 81963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: One level solver used (system is seen as DD). Using default solver.\n", ((PetscObject)pc)->prefix)); 8209566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, 0, &smoother)); 8219566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0])); 8229566063dSJacob Faibussowitsch PetscCall(KSPSetType(smoother, KSPPREONLY)); 8239566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 8245f8cf99dSMark F. Adams } 825849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 8263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8275b89ad90SMark F. Adams } 8285b89ad90SMark 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 */ 838d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy_GAMG(PC pc) 839d71ae5a4SJacob Faibussowitsch { 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)); 8451baa6e33SBarry Smith if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc)); 8469566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->ops)); 8479566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->gamg_type_name)); 8489566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg)); 8492e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", NULL)); 8502e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", NULL)); 8512e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", NULL)); 8522e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", NULL)); 8532e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", NULL)); 8542e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", NULL)); 8552e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", NULL)); 8562e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", NULL)); 8572e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", NULL)); 8582e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", NULL)); 8592e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", NULL)); 8602e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", NULL)); 8612e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", NULL)); 8622e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", NULL)); 8632e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", NULL)); 8642e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", NULL)); 8659566063dSJacob Faibussowitsch PetscCall(PCDestroy_MG(pc)); 8663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8675b89ad90SMark F. Adams } 8685b89ad90SMark F. Adams 869676e1743SMark F. Adams /*@ 870f1580f4eSBarry Smith PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction in `PCGAMG` 871676e1743SMark F. Adams 872c3339decSBarry Smith Logically Collective 873676e1743SMark F. Adams 874676e1743SMark F. Adams Input Parameters: 8751cc46a46SBarry Smith + pc - the preconditioner context 8761cc46a46SBarry Smith - n - the number of equations 877676e1743SMark F. Adams 878676e1743SMark F. Adams Options Database Key: 879147403d9SBarry Smith . -pc_gamg_process_eq_limit <limit> - set the limit 880676e1743SMark F. Adams 88120f4b53cSBarry Smith Level: intermediate 88220f4b53cSBarry Smith 883f1580f4eSBarry Smith Note: 884f1580f4eSBarry Smith `PCGAMG` will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process 885cab9ed1eSBarry Smith that has degrees of freedom 886cab9ed1eSBarry Smith 887f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()` 888676e1743SMark F. Adams @*/ 889d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 890d71ae5a4SJacob Faibussowitsch { 891676e1743SMark F. Adams PetscFunctionBegin; 892676e1743SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 893cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetProcEqLim_C", (PC, PetscInt), (pc, n)); 8943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 895676e1743SMark F. Adams } 896676e1743SMark F. Adams 897d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 898d71ae5a4SJacob Faibussowitsch { 899c20e4228SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 900c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 901676e1743SMark F. Adams 902676e1743SMark F. Adams PetscFunctionBegin; 9039d5b6da9SMark F. Adams if (n > 0) pc_gamg->min_eq_proc = n; 9043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 905676e1743SMark F. Adams } 906676e1743SMark F. Adams 907389730f3SMark F. Adams /*@ 908f1580f4eSBarry Smith PCGAMGSetCoarseEqLim - Set maximum number of equations on the coarsest grid of `PCGAMG` 909389730f3SMark F. Adams 910c3339decSBarry Smith Collective 911389730f3SMark F. Adams 912389730f3SMark F. Adams Input Parameters: 9131cc46a46SBarry Smith + pc - the preconditioner context 9141cc46a46SBarry Smith - n - maximum number of equations to aim for 915389730f3SMark F. Adams 916389730f3SMark F. Adams Options Database Key: 917147403d9SBarry Smith . -pc_gamg_coarse_eq_limit <limit> - set the limit 918389730f3SMark F. Adams 91920f4b53cSBarry Smith Level: intermediate 92020f4b53cSBarry Smith 921f1580f4eSBarry Smith Note: 922f1580f4eSBarry Smith For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 92374329af1SBarry Smith has less than 1000 unknowns. 92474329af1SBarry Smith 925f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()` 926389730f3SMark F. Adams @*/ 927d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 928d71ae5a4SJacob Faibussowitsch { 929389730f3SMark F. Adams PetscFunctionBegin; 930389730f3SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 931cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetCoarseEqLim_C", (PC, PetscInt), (pc, n)); 9323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 933389730f3SMark F. Adams } 934389730f3SMark F. Adams 935d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 936d71ae5a4SJacob Faibussowitsch { 937389730f3SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 938389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 939389730f3SMark F. Adams 940389730f3SMark F. Adams PetscFunctionBegin; 9419d5b6da9SMark F. Adams if (n > 0) pc_gamg->coarse_eq_limit = n; 9423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 943389730f3SMark F. Adams } 944389730f3SMark F. Adams 945676e1743SMark F. Adams /*@ 946f1580f4eSBarry Smith PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids when reducing the number of MPI ranks to use 947676e1743SMark F. Adams 948c3339decSBarry Smith Collective 949676e1743SMark F. Adams 950676e1743SMark F. Adams Input Parameters: 9511cc46a46SBarry Smith + pc - the preconditioner context 952f1580f4eSBarry Smith - n - `PETSC_TRUE` or `PETSC_FALSE` 953676e1743SMark F. Adams 954676e1743SMark F. Adams Options Database Key: 955147403d9SBarry Smith . -pc_gamg_repartition <true,false> - turn on the repartitioning 956676e1743SMark F. Adams 95720f4b53cSBarry Smith Level: intermediate 95820f4b53cSBarry Smith 959f1580f4eSBarry Smith Note: 960f1580f4eSBarry Smith This will generally improve the loading balancing of the work on each level 961cab9ed1eSBarry Smith 962f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()` 963676e1743SMark F. Adams @*/ 964d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) 965d71ae5a4SJacob Faibussowitsch { 966676e1743SMark F. Adams PetscFunctionBegin; 967676e1743SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 968cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetRepartition_C", (PC, PetscBool), (pc, n)); 9693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 970676e1743SMark F. Adams } 971676e1743SMark F. Adams 972d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) 973d71ae5a4SJacob Faibussowitsch { 974c20e4228SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 975c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 976676e1743SMark F. Adams 977676e1743SMark F. Adams PetscFunctionBegin; 9789d5b6da9SMark F. Adams pc_gamg->repart = n; 9793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 980676e1743SMark F. Adams } 981676e1743SMark F. Adams 982dfd5c07aSMark F. Adams /*@ 983f1580f4eSBarry Smith PCGAMGSetUseSAEstEig - Use the eigen estimate from smoothed aggregation for the Chebyshev smoother during the solution process 98418c3aa7eSMark 985c3339decSBarry Smith Collective 98618c3aa7eSMark 98718c3aa7eSMark Input Parameters: 98818c3aa7eSMark + pc - the preconditioner context 98918c3aa7eSMark - n - number of its 99018c3aa7eSMark 99118c3aa7eSMark Options Database Key: 992147403d9SBarry Smith . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate 99318c3aa7eSMark 99420f4b53cSBarry Smith Level: advanced 99520f4b53cSBarry Smith 99618c3aa7eSMark Notes: 9977e6512fdSJed 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$. 998f1580f4eSBarry Smith Eigenvalue estimates (based on a few `PCCG` or `PCGMRES` iterations) are computed to choose $\omega$ so that this is a stable smoothing operation. 999f1580f4eSBarry Smith If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused during the solution process 1000f1580f4eSBarry Smith This option is only used when the smoother uses Jacobi, and should be turned off if a different `PCJacobiType` is used. 1001efe053fcSJed Brown It became default in PETSc 3.17. 100218c3aa7eSMark 1003f1580f4eSBarry Smith .seealso: `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()` 100418c3aa7eSMark @*/ 1005d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n) 1006d71ae5a4SJacob Faibussowitsch { 100718c3aa7eSMark PetscFunctionBegin; 100818c3aa7eSMark PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1009cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetUseSAEstEig_C", (PC, PetscBool), (pc, n)); 10103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 101118c3aa7eSMark } 101218c3aa7eSMark 1013d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n) 1014d71ae5a4SJacob Faibussowitsch { 101518c3aa7eSMark PC_MG *mg = (PC_MG *)pc->data; 101618c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 101718c3aa7eSMark 101818c3aa7eSMark PetscFunctionBegin; 10197e6512fdSJed Brown pc_gamg->use_sa_esteig = n; 10203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 102118c3aa7eSMark } 102218c3aa7eSMark 102318c3aa7eSMark /*@ 1024f1580f4eSBarry Smith PCGAMGSetEigenvalues - Set WHAT eigenvalues WHY? 102518c3aa7eSMark 1026c3339decSBarry Smith Collective 102718c3aa7eSMark 102818c3aa7eSMark Input Parameters: 102918c3aa7eSMark + pc - the preconditioner context 103018c3aa7eSMark - emax - max eigenvalue 103118c3aa7eSMark . emin - min eigenvalue 103218c3aa7eSMark 103318c3aa7eSMark Options Database Key: 1034147403d9SBarry Smith . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues 103518c3aa7eSMark 103618c3aa7eSMark Level: intermediate 103718c3aa7eSMark 1038f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetUseSAEstEig()` 103918c3aa7eSMark @*/ 1040d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax, PetscReal emin) 1041d71ae5a4SJacob Faibussowitsch { 104218c3aa7eSMark PetscFunctionBegin; 104318c3aa7eSMark PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1044cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetEigenvalues_C", (PC, PetscReal, PetscReal), (pc, emax, emin)); 10453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 104618c3aa7eSMark } 104741ffd417SStefano Zampini 1048d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc, PetscReal emax, PetscReal emin) 1049d71ae5a4SJacob Faibussowitsch { 105018c3aa7eSMark PC_MG *mg = (PC_MG *)pc->data; 105118c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 105218c3aa7eSMark 105318c3aa7eSMark PetscFunctionBegin; 10545f80ce2aSJacob 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); 10555f80ce2aSJacob 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); 105618c3aa7eSMark pc_gamg->emax = emax; 105718c3aa7eSMark pc_gamg->emin = emin; 10583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 105918c3aa7eSMark } 106018c3aa7eSMark 106118c3aa7eSMark /*@ 1062f1580f4eSBarry Smith PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding a `PCGAMG` algebraic multigrid preconditioner 1063dfd5c07aSMark F. Adams 1064c3339decSBarry Smith Collective 1065dfd5c07aSMark F. Adams 1066dfd5c07aSMark F. Adams Input Parameters: 10671cc46a46SBarry Smith + pc - the preconditioner context 1068f1580f4eSBarry Smith - n - `PETSC_TRUE` or `PETSC_FALSE` 1069dfd5c07aSMark F. Adams 1070dfd5c07aSMark F. Adams Options Database Key: 1071147403d9SBarry Smith . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation 1072dfd5c07aSMark F. Adams 1073dfd5c07aSMark F. Adams Level: intermediate 1074dfd5c07aSMark F. Adams 1075f1580f4eSBarry Smith Note: 1076147403d9SBarry Smith May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1077cab9ed1eSBarry Smith rebuilding the preconditioner quicker. 1078cab9ed1eSBarry Smith 1079f1580f4eSBarry Smith .seealso: `PCGAMG` 1080dfd5c07aSMark F. Adams @*/ 1081d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 1082d71ae5a4SJacob Faibussowitsch { 1083dfd5c07aSMark F. Adams PetscFunctionBegin; 1084dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1085cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetReuseInterpolation_C", (PC, PetscBool), (pc, n)); 10863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1087dfd5c07aSMark F. Adams } 1088dfd5c07aSMark F. Adams 1089d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 1090d71ae5a4SJacob Faibussowitsch { 1091dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1092dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1093dfd5c07aSMark F. Adams 1094dfd5c07aSMark F. Adams PetscFunctionBegin; 1095dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 10963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1097dfd5c07aSMark F. Adams } 1098dfd5c07aSMark F. Adams 1099ffc955d6SMark F. Adams /*@ 1100f1580f4eSBarry 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 1101f1580f4eSBarry Smith used as the smoother 1102ffc955d6SMark F. Adams 1103c3339decSBarry Smith Collective 1104ffc955d6SMark F. Adams 1105ffc955d6SMark F. Adams Input Parameters: 1106cab9ed1eSBarry Smith + pc - the preconditioner context 1107f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use aggregates, `PETSC_FALSE` to not 1108ffc955d6SMark F. Adams 1109ffc955d6SMark F. Adams Options Database Key: 1110147403d9SBarry Smith . -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains 1111ffc955d6SMark F. Adams 1112ffc955d6SMark F. Adams Level: intermediate 1113ffc955d6SMark F. Adams 1114f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCASM`, `PCSetType` 1115ffc955d6SMark F. Adams @*/ 1116d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) 1117d71ae5a4SJacob Faibussowitsch { 1118ffc955d6SMark F. Adams PetscFunctionBegin; 1119ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1120cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGASMSetUseAggs_C", (PC, PetscBool), (pc, flg)); 11213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1122ffc955d6SMark F. Adams } 1123ffc955d6SMark F. Adams 1124d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) 1125d71ae5a4SJacob Faibussowitsch { 1126ffc955d6SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1127ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1128ffc955d6SMark F. Adams 1129ffc955d6SMark F. Adams PetscFunctionBegin; 1130cab9ed1eSBarry Smith pc_gamg->use_aggs_in_asm = flg; 11313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1132ffc955d6SMark F. Adams } 1133ffc955d6SMark F. Adams 1134171cca9aSMark Adams /*@ 1135cf8ae1d3SMark Adams PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver 1136171cca9aSMark Adams 1137c3339decSBarry Smith Collective 1138171cca9aSMark Adams 1139171cca9aSMark Adams Input Parameters: 1140171cca9aSMark Adams + pc - the preconditioner context 1141f1580f4eSBarry Smith - flg - `PETSC_TRUE` to not force coarse grid onto one processor 1142171cca9aSMark Adams 1143171cca9aSMark Adams Options Database Key: 1144147403d9SBarry Smith . -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver 1145171cca9aSMark Adams 1146171cca9aSMark Adams Level: intermediate 1147171cca9aSMark Adams 1148f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()` 1149171cca9aSMark Adams @*/ 1150d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg) 1151d71ae5a4SJacob Faibussowitsch { 1152171cca9aSMark Adams PetscFunctionBegin; 1153171cca9aSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1154cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetUseParallelCoarseGridSolve_C", (PC, PetscBool), (pc, flg)); 11553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1156171cca9aSMark Adams } 1157171cca9aSMark Adams 1158d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) 1159d71ae5a4SJacob Faibussowitsch { 1160171cca9aSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1161171cca9aSMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1162171cca9aSMark Adams 1163171cca9aSMark Adams PetscFunctionBegin; 1164171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = flg; 11653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1166ffc955d6SMark F. Adams } 1167ffc955d6SMark F. Adams 11684ef23d27SMark F. Adams /*@ 1169f1580f4eSBarry Smith PCGAMGSetCpuPinCoarseGrids - pin the coarse grids created in `PCGAMG` to run only on the CPU since the problems may be too small to run efficiently on the GPUs 1170ce7c7f2fSMark Adams 1171c3339decSBarry Smith Collective 1172ce7c7f2fSMark Adams 1173ce7c7f2fSMark Adams Input Parameters: 1174ce7c7f2fSMark Adams + pc - the preconditioner context 1175f1580f4eSBarry Smith - flg - `PETSC_TRUE` to pin coarse grids to the CPU 1176ce7c7f2fSMark Adams 1177ce7c7f2fSMark Adams Options Database Key: 1178147403d9SBarry Smith . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU 1179ce7c7f2fSMark Adams 1180ce7c7f2fSMark Adams Level: intermediate 1181ce7c7f2fSMark Adams 1182f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetUseParallelCoarseGridSolve()` 1183ce7c7f2fSMark Adams @*/ 1184d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) 1185d71ae5a4SJacob Faibussowitsch { 1186ce7c7f2fSMark Adams PetscFunctionBegin; 1187ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1188cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetCpuPinCoarseGrids_C", (PC, PetscBool), (pc, flg)); 11893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1190ce7c7f2fSMark Adams } 1191ce7c7f2fSMark Adams 1192d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) 1193d71ae5a4SJacob Faibussowitsch { 1194ce7c7f2fSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1195ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1196ce7c7f2fSMark Adams 1197ce7c7f2fSMark Adams PetscFunctionBegin; 1198ce7c7f2fSMark Adams pc_gamg->cpu_pin_coarse_grids = flg; 11993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1200ce7c7f2fSMark Adams } 1201ce7c7f2fSMark Adams 1202ce7c7f2fSMark Adams /*@ 1203147403d9SBarry Smith PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type) 1204ce7c7f2fSMark Adams 1205c3339decSBarry Smith Collective 1206ce7c7f2fSMark Adams 1207ce7c7f2fSMark Adams Input Parameters: 1208ce7c7f2fSMark Adams + pc - the preconditioner context 1209f1580f4eSBarry Smith - flg - `PCGAMGLayoutType` type, either `PCGAMG_LAYOUT_COMPACT` or `PCGAMG_LAYOUT_SPREAD` 1210ce7c7f2fSMark Adams 1211ce7c7f2fSMark Adams Options Database Key: 1212147403d9SBarry Smith . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering 1213ce7c7f2fSMark Adams 1214ce7c7f2fSMark Adams Level: intermediate 1215ce7c7f2fSMark Adams 1216f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGLayoutType`, `PCGAMG_LAYOUT_COMPACT`, `PCGAMG_LAYOUT_SPREAD` 1217ce7c7f2fSMark Adams @*/ 1218d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) 1219d71ae5a4SJacob Faibussowitsch { 1220ce7c7f2fSMark Adams PetscFunctionBegin; 1221ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1222cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetCoarseGridLayoutType_C", (PC, PCGAMGLayoutType), (pc, flg)); 12233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1224ce7c7f2fSMark Adams } 1225ce7c7f2fSMark Adams 1226d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) 1227d71ae5a4SJacob Faibussowitsch { 1228ce7c7f2fSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1229ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1230ce7c7f2fSMark Adams 1231ce7c7f2fSMark Adams PetscFunctionBegin; 1232ce7c7f2fSMark Adams pc_gamg->layout_type = flg; 12333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1234ce7c7f2fSMark Adams } 1235ce7c7f2fSMark Adams 1236ce7c7f2fSMark Adams /*@ 1237f1580f4eSBarry Smith PCGAMGSetNlevels - Sets the maximum number of levels `PCGAMG` will use 12384ef23d27SMark F. Adams 123920f4b53cSBarry Smith Not Collective 12404ef23d27SMark F. Adams 12414ef23d27SMark F. Adams Input Parameters: 12421cc46a46SBarry Smith + pc - the preconditioner 12431cc46a46SBarry Smith - n - the maximum number of levels to use 12444ef23d27SMark F. Adams 12454ef23d27SMark F. Adams Options Database Key: 1246147403d9SBarry Smith . -pc_mg_levels <n> - set the maximum number of levels to allow 12474ef23d27SMark F. Adams 12484ef23d27SMark F. Adams Level: intermediate 12494ef23d27SMark F. Adams 1250f1580f4eSBarry Smith Developer Note: 1251f1580f4eSBarry Smith Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG` 1252f1580f4eSBarry Smith 1253f1580f4eSBarry Smith .seealso: `PCGAMG` 12544ef23d27SMark F. Adams @*/ 1255d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 1256d71ae5a4SJacob Faibussowitsch { 12574ef23d27SMark F. Adams PetscFunctionBegin; 12584ef23d27SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1259cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetNlevels_C", (PC, PetscInt), (pc, n)); 12603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12614ef23d27SMark F. Adams } 12624ef23d27SMark F. Adams 1263d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 1264d71ae5a4SJacob Faibussowitsch { 12654ef23d27SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 12664ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 12674ef23d27SMark F. Adams 12684ef23d27SMark F. Adams PetscFunctionBegin; 12699d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 12703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12714ef23d27SMark F. Adams } 12724ef23d27SMark F. Adams 12733542efc5SMark F. Adams /*@ 12743542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 12753542efc5SMark F. Adams 127620f4b53cSBarry Smith Not Collective 12773542efc5SMark F. Adams 12783542efc5SMark F. Adams Input Parameters: 12791cc46a46SBarry Smith + pc - the preconditioner context 1280c9567895SMark . 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 1281055c8bd0SJed Brown - n - number of threshold values provided in array 12823542efc5SMark F. Adams 12833542efc5SMark F. Adams Options Database Key: 1284147403d9SBarry Smith . -pc_gamg_threshold <threshold> - the threshold to drop edges 12853542efc5SMark F. Adams 128620f4b53cSBarry Smith Level: intermediate 128720f4b53cSBarry Smith 128895452b02SPatrick Sanan Notes: 1289af3c827dSMark 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. 1290f1580f4eSBarry Smith Before coarsening or aggregating the graph, `PCGAMG` 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. 1291cab9ed1eSBarry Smith 129220f4b53cSBarry Smith If `n` is less than the total number of coarsenings (see `PCGAMGSetNlevels()`), then threshold scaling (see `PCGAMGSetThresholdScale()`) is used for each successive coarsening. 1293f1580f4eSBarry Smith In this case, `PCGAMGSetThresholdScale()` must be called before `PCGAMGSetThreshold()`. 129420f4b53cSBarry Smith If `n` is greater than the total number of levels, the excess entries in threshold will not be used. 12953542efc5SMark F. Adams 12962eab5cd7Smarkadams4 .seealso: `PCGAMG`, `PCGAMGSetAggressiveLevels()`, `PCGAMGSetThresholdScale()` 12973542efc5SMark F. Adams @*/ 1298d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) 1299d71ae5a4SJacob Faibussowitsch { 13003542efc5SMark F. Adams PetscFunctionBegin; 13013542efc5SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1302055c8bd0SJed Brown if (n) PetscValidRealPointer(v, 2); 1303cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetThreshold_C", (PC, PetscReal[], PetscInt), (pc, v, n)); 13043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13053542efc5SMark F. Adams } 13063542efc5SMark F. Adams 1307d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) 1308d71ae5a4SJacob Faibussowitsch { 1309c20e4228SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1310c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1311c1eae691SMark Adams PetscInt i; 1312c1eae691SMark Adams PetscFunctionBegin; 1313055c8bd0SJed Brown for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i]; 1314055c8bd0SJed Brown for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale; 13153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1316c1eae691SMark Adams } 1317c1eae691SMark Adams 1318c1eae691SMark Adams /*@ 1319f1580f4eSBarry Smith PCGAMGSetRankReductionFactors - Set a manual schedule for MPI rank reduction on coarse grids 1320c9567895SMark 1321c3339decSBarry Smith Collective 1322c9567895SMark 1323c9567895SMark Input Parameters: 1324c9567895SMark + pc - the preconditioner context 1325f1580f4eSBarry Smith . v - array of reduction factors. 0 for first value forces a reduction to one process/device on first level in CUDA 1326c9567895SMark - n - number of values provided in array 1327c9567895SMark 1328c9567895SMark Options Database Key: 1329147403d9SBarry Smith . -pc_gamg_rank_reduction_factors <factors> - provide the schedule 1330c9567895SMark 1331c9567895SMark Level: intermediate 1332c9567895SMark 1333f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()` 1334c9567895SMark @*/ 1335d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n) 1336d71ae5a4SJacob Faibussowitsch { 1337c9567895SMark PetscFunctionBegin; 1338c9567895SMark PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1339c9567895SMark if (n) PetscValidIntPointer(v, 2); 1340cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetRankReductionFactors_C", (PC, PetscInt[], PetscInt), (pc, v, n)); 13413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1342c9567895SMark } 1343c9567895SMark 1344d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n) 1345d71ae5a4SJacob Faibussowitsch { 1346c9567895SMark PC_MG *mg = (PC_MG *)pc->data; 1347c9567895SMark PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1348c9567895SMark PetscInt i; 1349c9567895SMark PetscFunctionBegin; 1350c9567895SMark for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i]; 1351c9567895SMark for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */ 13523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1353c9567895SMark } 1354c9567895SMark 1355c9567895SMark /*@ 1356c1eae691SMark Adams PCGAMGSetThresholdScale - Relative threshold reduction at each level 1357c1eae691SMark Adams 135820f4b53cSBarry Smith Not Collective 1359c1eae691SMark Adams 1360c1eae691SMark Adams Input Parameters: 1361c1eae691SMark Adams + pc - the preconditioner context 1362147403d9SBarry Smith - scale - the threshold value reduction, usually < 1.0 1363c1eae691SMark Adams 1364c1eae691SMark Adams Options Database Key: 1365147403d9SBarry Smith . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level 1366c1eae691SMark Adams 136720f4b53cSBarry Smith Level: advanced 136820f4b53cSBarry Smith 1369f1580f4eSBarry Smith Note: 1370f1580f4eSBarry Smith The initial threshold (for an arbitrary number of levels starting from the finest) can be set with `PCGAMGSetThreshold()`. 1371f1580f4eSBarry Smith This scaling is used for each subsequent coarsening, but must be called before `PCGAMGSetThreshold()`. 1372055c8bd0SJed Brown 1373f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetThreshold()` 1374c1eae691SMark Adams @*/ 1375d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) 1376d71ae5a4SJacob Faibussowitsch { 13773542efc5SMark F. Adams PetscFunctionBegin; 1378c1eae691SMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1379cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetThresholdScale_C", (PC, PetscReal), (pc, v)); 13803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1381c1eae691SMark Adams } 1382c1eae691SMark Adams 1383d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) 1384d71ae5a4SJacob Faibussowitsch { 1385c1eae691SMark Adams PC_MG *mg = (PC_MG *)pc->data; 1386c1eae691SMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1387c1eae691SMark Adams PetscFunctionBegin; 1388c1eae691SMark Adams pc_gamg->threshold_scale = v; 13893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13903542efc5SMark F. Adams } 13913542efc5SMark F. Adams 1392e20c40e8SBarry Smith /*@C 1393f1580f4eSBarry Smith PCGAMGSetType - Set the type of algorithm `PCGAMG` should use 1394676e1743SMark F. Adams 1395c3339decSBarry Smith Collective 1396676e1743SMark F. Adams 1397676e1743SMark F. Adams Input Parameters: 1398c60c7ad4SBarry Smith + pc - the preconditioner context 1399f1580f4eSBarry Smith - type - `PCGAMGAGG`, `PCGAMGGEO`, or `PCGAMGCLASSICAL` 1400676e1743SMark F. Adams 1401676e1743SMark F. Adams Options Database Key: 1402*5e4ac8c8Smarkadams4 . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply - only agg supported 1403676e1743SMark F. Adams 1404676e1743SMark F. Adams Level: intermediate 1405676e1743SMark F. Adams 1406db781477SPatrick Sanan .seealso: `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType` 1407676e1743SMark F. Adams @*/ 1408d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1409d71ae5a4SJacob Faibussowitsch { 1410676e1743SMark F. Adams PetscFunctionBegin; 1411676e1743SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1412cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetType_C", (PC, PCGAMGType), (pc, type)); 14133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1414676e1743SMark F. Adams } 1415676e1743SMark F. Adams 1416e20c40e8SBarry Smith /*@C 1417f1580f4eSBarry Smith PCGAMGGetType - Get the type of algorithm `PCGAMG` will use 1418c60c7ad4SBarry Smith 1419c3339decSBarry Smith Collective 1420c60c7ad4SBarry Smith 1421c60c7ad4SBarry Smith Input Parameter: 1422c60c7ad4SBarry Smith . pc - the preconditioner context 1423c60c7ad4SBarry Smith 1424c60c7ad4SBarry Smith Output Parameter: 1425c60c7ad4SBarry Smith . type - the type of algorithm used 1426c60c7ad4SBarry Smith 1427c60c7ad4SBarry Smith Level: intermediate 1428c60c7ad4SBarry Smith 1429f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetType()`, `PCGAMGType` 1430c60c7ad4SBarry Smith @*/ 1431d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1432d71ae5a4SJacob Faibussowitsch { 1433c60c7ad4SBarry Smith PetscFunctionBegin; 1434c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1435cac4c232SBarry Smith PetscUseMethod(pc, "PCGAMGGetType_C", (PC, PCGAMGType *), (pc, type)); 14363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1437c60c7ad4SBarry Smith } 1438c60c7ad4SBarry Smith 1439d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1440d71ae5a4SJacob Faibussowitsch { 1441c60c7ad4SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1442c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1443c60c7ad4SBarry Smith 1444c60c7ad4SBarry Smith PetscFunctionBegin; 1445c60c7ad4SBarry Smith *type = pc_gamg->type; 14463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1447c60c7ad4SBarry Smith } 1448c60c7ad4SBarry Smith 1449d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1450d71ae5a4SJacob Faibussowitsch { 14511ab5ffc9SJed Brown PC_MG *mg = (PC_MG *)pc->data; 14521ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 14535f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(PC); 1454676e1743SMark F. Adams 1455676e1743SMark F. Adams PetscFunctionBegin; 1456c60c7ad4SBarry Smith pc_gamg->type = type; 14579566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(GAMGList, type, &r)); 14586adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown GAMG type %s given", type); 14591ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 14609566063dSJacob Faibussowitsch PetscCall((*pc_gamg->ops->destroy)(pc)); 14619566063dSJacob Faibussowitsch PetscCall(PetscMemzero(pc_gamg->ops, sizeof(struct _PCGAMGOps))); 1462e616c208SToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1463da81f932SPierre Jolivet /* cleaning up common data in pc_gamg - this should disappear someday */ 14643ae0bb68SMark Adams pc_gamg->data_cell_cols = 0; 14653ae0bb68SMark Adams pc_gamg->data_cell_rows = 0; 14663ae0bb68SMark Adams pc_gamg->orig_data_cell_cols = 0; 14673ae0bb68SMark Adams pc_gamg->orig_data_cell_rows = 0; 14689566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 14693ae0bb68SMark Adams pc_gamg->data_sz = 0; 14701ab5ffc9SJed Brown } 14719566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->gamg_type_name)); 14729566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(type, &pc_gamg->gamg_type_name)); 14739566063dSJacob Faibussowitsch PetscCall((*r)(pc)); 14743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1475676e1743SMark F. Adams } 1476676e1743SMark F. Adams 1477d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG(PC pc, PetscViewer viewer) 1478d71ae5a4SJacob Faibussowitsch { 14795adeb434SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 14805adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1481e7d4b4cbSMark Adams PetscReal gc = 0, oc = 0; 148290db8557SMark Adams 14835adeb434SBarry Smith PetscFunctionBegin; 14849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GAMG specific options\n")); 14859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Threshold for dropping small values in graph on each level =")); 14869566063dSJacob Faibussowitsch for (PetscInt i = 0; i < mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)pc_gamg->threshold[i])); 14879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 14889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Threshold scaling factor for each level not specified = %g\n", (double)pc_gamg->threshold_scale)); 148948a46eb9SPierre Jolivet if (pc_gamg->use_aggs_in_asm) PetscCall(PetscViewerASCIIPrintf(viewer, " Using aggregates from coarsening process to define subdomains for PCASM\n")); 149048a46eb9SPierre Jolivet if (pc_gamg->use_parallel_coarse_grid_solver) PetscCall(PetscViewerASCIIPrintf(viewer, " Using parallel coarse grid solver (all coarse grid equations not put on one process)\n")); 14911baa6e33SBarry Smith if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc, viewer)); 14929566063dSJacob Faibussowitsch PetscCall(PCMGGetGridComplexity(pc, &gc, &oc)); 149363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Complexity: grid = %g operator = %g\n", (double)gc, (double)oc)); 14943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14955adeb434SBarry Smith } 14965adeb434SBarry Smith 1497d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFromOptions_GAMG(PC pc, PetscOptionItems *PetscOptionsObject) 1498d71ae5a4SJacob Faibussowitsch { 1499676e1743SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1500676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 15017e6512fdSJed Brown PetscBool flag; 15023b4367a7SBarry Smith MPI_Comm comm; 150318c3aa7eSMark char prefix[256], tname[32]; 1504c1eae691SMark Adams PetscInt i, n; 150514a9496bSBarry Smith const char *pcpre; 15060a545947SLisandro Dalcin static const char *LayoutTypes[] = {"compact", "spread", "PCGAMGLayoutType", "PC_GAMG_LAYOUT", NULL}; 15075b89ad90SMark F. Adams PetscFunctionBegin; 15089566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 1509d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG options"); 1510*5e4ac8c8Smarkadams4 PetscCall(PetscOptionsFList("-pc_gamg_type", "Type of AMG method (only 'agg' supported and useful)", "PCGAMGSetType", GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag)); 15111baa6e33SBarry Smith if (flag) PetscCall(PCGAMGSetType(pc, tname)); 15129566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_repartition", "Repartion coarse grids", "PCGAMGSetRepartition", pc_gamg->repart, &pc_gamg->repart, NULL)); 1513f1580f4eSBarry Smith 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)); 15149566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation", "Reuse prolongation operator", "PCGAMGReuseInterpolation", pc_gamg->reuse_prol, &pc_gamg->reuse_prol, NULL)); 15159566063dSJacob 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)); 15169566063dSJacob 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)); 15179566063dSJacob 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)); 15189371c9d4SSatish Balay 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, 15199371c9d4SSatish Balay (PetscEnum)pc_gamg->layout_type, (PetscEnum *)&pc_gamg->layout_type, NULL)); 15209566063dSJacob 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)); 15219566063dSJacob 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)); 15229566063dSJacob 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)); 152318c3aa7eSMark n = PETSC_MG_MAXLEVELS; 15249566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", pc_gamg->threshold, &n, &flag)); 152518c3aa7eSMark if (!flag || n < PETSC_MG_MAXLEVELS) { 1526efd3c5ceSMark Adams if (!flag) n = 1; 1527c1eae691SMark Adams i = n; 1528d71ae5a4SJacob Faibussowitsch do { 1529d71ae5a4SJacob Faibussowitsch pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale; 1530d71ae5a4SJacob Faibussowitsch } while (++i < PETSC_MG_MAXLEVELS); 1531c1eae691SMark Adams } 1532c9567895SMark n = PETSC_MG_MAXLEVELS; 15339566063dSJacob 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)); 1534c9567895SMark if (!flag) i = 0; 1535c9567895SMark else i = n; 1536d71ae5a4SJacob Faibussowitsch do { 1537d71ae5a4SJacob Faibussowitsch pc_gamg->level_reduction_factors[i] = -1; 1538d71ae5a4SJacob Faibussowitsch } while (++i < PETSC_MG_MAXLEVELS); 15399566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_levels", "Set number of MG levels", "PCGAMGSetNlevels", pc_gamg->Nlevels, &pc_gamg->Nlevels, NULL)); 154018c3aa7eSMark { 154118c3aa7eSMark PetscReal eminmax[2] = {0., 0.}; 154218c3aa7eSMark n = 2; 15439566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues", "extreme eigenvalues for smoothed aggregation", "PCGAMGSetEigenvalues", eminmax, &n, &flag)); 154418c3aa7eSMark if (flag) { 154508401ef6SPierre Jolivet PetscCheck(n == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 15469566063dSJacob Faibussowitsch PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0])); 154718c3aa7eSMark } 154818c3aa7eSMark } 1549b7cbab4eSMark Adams /* set options for subtype */ 1550dbbe0bcdSBarry Smith PetscCall((*pc_gamg->ops->setfromoptions)(pc, PetscOptionsObject)); 155118c3aa7eSMark 15529566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 15539566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : "")); 1554d0609cedSBarry Smith PetscOptionsHeadEnd(); 15553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15565b89ad90SMark F. Adams } 15575b89ad90SMark F. Adams 15585b89ad90SMark F. Adams /*MC 15591cc46a46SBarry Smith PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 15605b89ad90SMark F. Adams 1561280d9858SJed Brown Options Database Keys: 1562*5e4ac8c8Smarkadams4 + -pc_gamg_type <type,default=agg> - one of agg, geo, or classical (only smoothed aggregation, agg, supported) 1563da81f932SPierre Jolivet . -pc_gamg_repartition <bool,default=false> - repartition the degrees of freedom across the coarse grids as they are determined 156421d928e4Smarkadams4 . -pc_gamg_asm_use_agg <bool,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother 1565da81f932SPierre Jolivet . -pc_gamg_process_eq_limit <limit, default=50> - `PCGAMG` will reduce the number of MPI ranks used directly on the coarse grids so that there are around <limit> 1566cab9ed1eSBarry Smith equations on each process that has degrees of freedom 15672d776b49SBarry Smith . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 156821d928e4Smarkadams4 . -pc_gamg_reuse_interpolation <bool,default=true> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations (should always be true) 156921d928e4Smarkadams4 . -pc_gamg_threshold[] <thresh,default=[-1,...]> - Before aggregating the graph `PCGAMG` will remove small values from the graph on each level (< 0 does no filtering) 15702d776b49SBarry Smith - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1571cab9ed1eSBarry Smith 1572f1580f4eSBarry Smith Options Database Keys for Aggregation: 1573cab9ed1eSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 15742d776b49SBarry Smith . -pc_gamg_square_graph <n,default=1> - alias for pc_gamg_aggressive_coarsening (deprecated) 1575bae903cbSmarkadams4 - -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest. 1576cab9ed1eSBarry Smith 1577f1580f4eSBarry Smith Options Database Keys for Multigrid: 1578a9f5add0SYANG Zongze + -pc_mg_cycle_type <v> - v or w, see `PCMGSetCycleType()` 1579db9745e2SBarry Smith . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1580db9745e2SBarry Smith . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 158121d928e4Smarkadams4 - -pc_mg_levels <levels> - Number of levels of multigrid to use. GAMG has a heuristic so pc_mg_levels is not usually used with GAMG 15825b89ad90SMark F. Adams 158320f4b53cSBarry Smith Level: intermediate 158420f4b53cSBarry Smith 158595452b02SPatrick Sanan Notes: 1586f1580f4eSBarry Smith To obtain good performance for `PCGAMG` for vector valued problems you must 1587f1580f4eSBarry Smith call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point 1588f1580f4eSBarry Smith call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator 1589f1580f4eSBarry Smith 1590f1580f4eSBarry Smith See [the Users Manual section on PCGAMG](sec_amg) for more details. 15911cc46a46SBarry Smith 1592db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, 1593db781477SPatrick Sanan `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGSetUseSAEstEig()` 15945b89ad90SMark F. Adams M*/ 1595b2573a8aSBarry Smith 1596d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 1597d71ae5a4SJacob Faibussowitsch { 15985b89ad90SMark F. Adams PC_GAMG *pc_gamg; 15995b89ad90SMark F. Adams PC_MG *mg; 16005b89ad90SMark F. Adams 16015b89ad90SMark F. Adams PetscFunctionBegin; 16021c1aac46SBarry Smith /* register AMG type */ 16039566063dSJacob Faibussowitsch PetscCall(PCGAMGInitializePackage()); 16041c1aac46SBarry Smith 16055b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 16069566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCMG)); 16079566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG)); 16085b89ad90SMark F. Adams 16095b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 16104dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg)); 16119566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 16125b89ad90SMark F. Adams mg = (PC_MG *)pc->data; 16135b89ad90SMark F. Adams mg->innerctx = pc_gamg; 16145b89ad90SMark F. Adams 16154dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg->ops)); 16161ab5ffc9SJed Brown 16179d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 16189d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 16190a545947SLisandro Dalcin pc_gamg->data = NULL; 16205b89ad90SMark F. Adams 16219d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 16225b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 16235b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 16245b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 16255b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 16265adeb434SBarry Smith mg->view = PCView_GAMG; 16275b89ad90SMark F. Adams 16289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", PCGAMGSetProcEqLim_GAMG)); 16299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", PCGAMGSetCoarseEqLim_GAMG)); 16309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", PCGAMGSetRepartition_GAMG)); 16319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", PCGAMGSetEigenvalues_GAMG)); 16329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", PCGAMGSetUseSAEstEig_GAMG)); 16339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", PCGAMGSetReuseInterpolation_GAMG)); 16349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", PCGAMGASMSetUseAggs_GAMG)); 16359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", PCGAMGSetUseParallelCoarseGridSolve_GAMG)); 16369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", PCGAMGSetCpuPinCoarseGrids_GAMG)); 16379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", PCGAMGSetCoarseGridLayoutType_GAMG)); 16389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", PCGAMGSetThreshold_GAMG)); 16399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", PCGAMGSetRankReductionFactors_GAMG)); 16409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", PCGAMGSetThresholdScale_GAMG)); 16419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", PCGAMGSetType_GAMG)); 16429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", PCGAMGGetType_GAMG)); 16439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", PCGAMGSetNlevels_GAMG)); 16449d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 164521d928e4Smarkadams4 pc_gamg->reuse_prol = PETSC_TRUE; 16460c3bc534SBarry Smith pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1647171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1648a0095786SMark pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1649a0095786SMark pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1650038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 165125a145a7SMark Adams pc_gamg->coarse_eq_limit = 50; 165253134ebeSMark Adams for (int i = 0; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = -1; 1653c1eae691SMark Adams pc_gamg->threshold_scale = 1.; 165418c3aa7eSMark pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 16559ab59c8bSMark Adams pc_gamg->current_level = 0; /* don't need to init really */ 16567e6512fdSJed Brown pc_gamg->use_sa_esteig = PETSC_TRUE; 165718c3aa7eSMark pc_gamg->emin = 0; 165818c3aa7eSMark pc_gamg->emax = 0; 165918c3aa7eSMark 1660c238b0ebSToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 16619d5b6da9SMark F. Adams 1662bd94a7aaSJed Brown /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 16639566063dSJacob Faibussowitsch PetscCall(PCGAMGSetType(pc, PCGAMGAGG)); 16643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16655b89ad90SMark F. Adams } 16663e3471ccSMark Adams 16673e3471ccSMark Adams /*@C 1668f1580f4eSBarry Smith PCGAMGInitializePackage - This function initializes everything in the `PCGAMG` package. It is called 1669f1580f4eSBarry Smith from `PCInitializePackage()`. 16703e3471ccSMark Adams 16713e3471ccSMark Adams Level: developer 16723e3471ccSMark Adams 1673db781477SPatrick Sanan .seealso: `PetscInitialize()` 16743e3471ccSMark Adams @*/ 1675d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGInitializePackage(void) 1676d71ae5a4SJacob Faibussowitsch { 16774555aa8cSStefano Zampini PetscInt l; 16783e3471ccSMark Adams 16793e3471ccSMark Adams PetscFunctionBegin; 16803ba16761SJacob Faibussowitsch if (PCGAMGPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 16813e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 16829566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGGEO, PCCreateGAMG_GEO)); 16839566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGAGG, PCCreateGAMG_AGG)); 16849566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGCLASSICAL, PCCreateGAMG_Classical)); 16859566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage)); 1686c1c463dbSMark Adams 1687c1c463dbSMark Adams /* general events */ 1688849bee69SMark Adams PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP])); 1689849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH])); 1690849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN])); 1691849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS])); 1692849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL])); 1693849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA])); 1694849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB])); 1695849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT])); 1696849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM])); 1697849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL])); 1698849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG PtAP", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP])); 1699849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Reduce", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE])); 1700849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART])); 1701849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */ 1702849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */ 1703849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */ 17044555aa8cSStefano Zampini for (l = 0; l < PETSC_MG_MAXLEVELS; l++) { 17054555aa8cSStefano Zampini char ename[32]; 17065b89ad90SMark F. Adams 170763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Squ l%02" PetscInt_FMT, l)); 17089566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0])); 170963a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Gal l%02" PetscInt_FMT, l)); 17109566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1])); 171163a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Opt l%02" PetscInt_FMT, l)); 17129566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2])); 17134555aa8cSStefano Zampini } 17144555aa8cSStefano Zampini #if defined(GAMG_STAGES) 1715849bee69SMark Adams { /* create timer stages */ 17165b89ad90SMark F. Adams char str[32]; 1717a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", 0)); 17189566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(str, &gamg_stages[0])); 17195b89ad90SMark F. Adams } 17205b89ad90SMark F. Adams #endif 17213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17223e3471ccSMark Adams } 17233e3471ccSMark Adams 17243e3471ccSMark Adams /*@C 1725f1580f4eSBarry Smith PCGAMGFinalizePackage - This function frees everything from the `PCGAMG` package. It is 1726f1580f4eSBarry Smith called from `PetscFinalize()` automatically. 17273e3471ccSMark Adams 17283e3471ccSMark Adams Level: developer 17293e3471ccSMark Adams 1730db781477SPatrick Sanan .seealso: `PetscFinalize()` 17313e3471ccSMark Adams @*/ 1732d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGFinalizePackage(void) 1733d71ae5a4SJacob Faibussowitsch { 17343e3471ccSMark Adams PetscFunctionBegin; 17353e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 17369566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&GAMGList)); 17373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17383e3471ccSMark Adams } 1739a36cf38bSToby Isaac 1740a36cf38bSToby Isaac /*@C 1741f1580f4eSBarry Smith PCGAMGRegister - Register a `PCGAMG` implementation. 1742a36cf38bSToby Isaac 1743a36cf38bSToby Isaac Input Parameters: 1744f1580f4eSBarry Smith + type - string that will be used as the name of the `PCGAMG` type. 1745a36cf38bSToby Isaac - create - function for creating the gamg context. 1746a36cf38bSToby Isaac 1747f1580f4eSBarry Smith Level: developer 1748a36cf38bSToby Isaac 1749db781477SPatrick Sanan .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()` 1750a36cf38bSToby Isaac @*/ 1751d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1752d71ae5a4SJacob Faibussowitsch { 1753a36cf38bSToby Isaac PetscFunctionBegin; 17549566063dSJacob Faibussowitsch PetscCall(PCGAMGInitializePackage()); 17559566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, type, create)); 17563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1757a36cf38bSToby Isaac } 17582d776b49SBarry Smith 17592d776b49SBarry Smith /*@C 17602d776b49SBarry Smith PCGAMGCreateGraph - Creates a graph that is used by the ``PCGAMGType`` in the coarsening process 17612d776b49SBarry Smith 17622d776b49SBarry Smith Input Parameters: 17632d776b49SBarry Smith + pc - the `PCGAMG` 17642d776b49SBarry Smith - A - the matrix, for any level 17652d776b49SBarry Smith 17662d776b49SBarry Smith Output Parameter: 17672d776b49SBarry Smith . G - the graph 17682d776b49SBarry Smith 17692d776b49SBarry Smith Level: advanced 17702d776b49SBarry Smith 17712d776b49SBarry Smith .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()` 17722d776b49SBarry Smith @*/ 17732d776b49SBarry Smith PetscErrorCode PCGAMGCreateGraph(PC pc, Mat A, Mat *G) 17742d776b49SBarry Smith { 17752d776b49SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 17762d776b49SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 17772d776b49SBarry Smith 17782d776b49SBarry Smith PetscFunctionBegin; 17792d776b49SBarry Smith PetscCall(pc_gamg->ops->creategraph(pc, A, G)); 17803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17812d776b49SBarry Smith } 1782