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 279371c9d4SSatish Balay PetscErrorCode PCReset_GAMG(PC pc) { 28d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 29d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 30d3d6bff4SMark F. Adams 31d3d6bff4SMark F. Adams PetscFunctionBegin; 329566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 331c1aac46SBarry Smith pc_gamg->data_sz = 0; 349566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->orig_data)); 355f80ce2aSJacob Faibussowitsch for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS; level++) { 3618c3aa7eSMark mg->min_eigen_DinvA[level] = 0; 3718c3aa7eSMark mg->max_eigen_DinvA[level] = 0; 3818c3aa7eSMark } 3918c3aa7eSMark pc_gamg->emin = 0; 4018c3aa7eSMark pc_gamg->emax = 0; 41a2f3521dSMark F. Adams PetscFunctionReturn(0); 42a2f3521dSMark F. Adams } 43a2f3521dSMark F. Adams 445b89ad90SMark F. Adams /* 45c238b0ebSToby Isaac PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 46a147abb0SMark F. Adams of active processors. 475b89ad90SMark F. Adams 485b89ad90SMark F. Adams Input Parameter: 49a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 50a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 519d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 52c5bfad50SMark F. Adams . cr_bs - coarse block size 533530afc2SMark F. Adams In/Output Parameter: 54a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 55afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 5611e60469SMark F. Adams Output Parameter: 573530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 585b89ad90SMark F. Adams */ 599371c9d4SSatish Balay 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) { 609d5b6da9SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 61486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 62a2f3521dSMark F. Adams Mat Cmat, Pold = *a_P_inout; 633b4367a7SBarry Smith MPI_Comm comm; 64c5df96a5SBarry Smith PetscMPIInt rank, size, new_size, nactive = *a_nactive_proc; 653ae0bb68SMark Adams PetscInt ncrs_eq, ncrs, f_bs; 665b89ad90SMark F. Adams 675b89ad90SMark F. Adams PetscFunctionBegin; 689566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat_fine, &comm)); 699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 719566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat_fine, &f_bs)); 72849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 739566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 749566063dSJacob Faibussowitsch PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat)); 759566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 76849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 77038e3b61SMark F. Adams 78ce7c7f2fSMark Adams if (Pcolumnperm) *Pcolumnperm = NULL; 79ce7c7f2fSMark Adams 803ae0bb68SMark Adams /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 819566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL)); 823ae0bb68SMark Adams if (pc_gamg->data_cell_rows > 0) { 833ae0bb68SMark Adams ncrs = pc_gamg->data_sz / pc_gamg->data_cell_cols / pc_gamg->data_cell_rows; 8473911c69SBarry Smith } else { 853ae0bb68SMark Adams PetscInt bs; 869566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Cmat, &bs)); 873ae0bb68SMark Adams ncrs = ncrs_eq / bs; 883ae0bb68SMark Adams } 89c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 90c9567895SMark 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. */ 91c9567895SMark #if defined(PETSC_HAVE_CUDA) 92c9567895SMark PetscShmComm pshmcomm; 93c9567895SMark PetscMPIInt locrank; 94c9567895SMark MPI_Comm loccomm; 95c9567895SMark PetscInt s_nnodes, r_nnodes, new_new_size; 96c9567895SMark cudaError_t cerr; 97c9567895SMark int devCount; 989566063dSJacob Faibussowitsch PetscCall(PetscShmCommGet(comm, &pshmcomm)); 999566063dSJacob Faibussowitsch PetscCall(PetscShmCommGetMpiShmComm(pshmcomm, &loccomm)); 1009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(loccomm, &locrank)); 101c9567895SMark s_nnodes = !locrank; 1029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&s_nnodes, &r_nnodes, 1, MPIU_INT, MPI_SUM, comm)); 10363a3b9bcSJacob Faibussowitsch PetscCheck((size % r_nnodes) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "odd number of nodes np=%d nnodes%" PetscInt_FMT, size, r_nnodes); 104c9567895SMark devCount = 0; 105c9567895SMark cerr = cudaGetDeviceCount(&devCount); 106c9567895SMark cudaGetLastError(); /* Reset the last error */ 107c9567895SMark if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */ 108c9567895SMark new_new_size = r_nnodes * devCount; 109c9567895SMark new_size = new_new_size; 11063a3b9bcSJacob 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)); 111c9567895SMark } else { 11263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: With Cuda but no device. Use heuristics.", ((PetscObject)pc)->prefix)); 113c9567895SMark goto HEURISTIC; 114c9567895SMark } 115c9567895SMark #else 116c9567895SMark SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "should not be here"); 117c9567895SMark #endif 118c9567895SMark } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) { 11963a3b9bcSJacob 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]); 120c9567895SMark new_size = nactive / pc_gamg->level_reduction_factors[pc_gamg->current_level]; 12163a3b9bcSJacob 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])); 122c9567895SMark } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) { 123c9567895SMark new_size = 1; 1249566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Force coarsest grid reduction to %d active processes\n", ((PetscObject)pc)->prefix, new_size)); 125c9567895SMark } else { 126472110cdSMark F. Adams PetscInt ncrs_eq_glob; 127c9567895SMark #if defined(PETSC_HAVE_CUDA) 128c9567895SMark HEURISTIC: 129c9567895SMark #endif 1309566063dSJacob Faibussowitsch PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL)); 131a90e85d9SMark Adams new_size = (PetscMPIInt)((float)ncrs_eq_glob / (float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 1327f66b68fSMark Adams if (!new_size) new_size = 1; /* not likely, posible? */ 133c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 1349566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Coarse grid reduction from %d to %d active processes\n", ((PetscObject)pc)->prefix, nactive, new_size)); 135a2f3521dSMark F. Adams } 1362e3501ffSMark Adams if (new_size == nactive) { 137ef3f0257SMark Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 138ce7c7f2fSMark Adams if (new_size < size) { 139ce7c7f2fSMark Adams /* odd case where multiple coarse grids are on one processor or no coarsening ... */ 1409566063dSJacob 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)); 141ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 1429566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE)); 1439566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE)); 144ce7c7f2fSMark Adams } 145ce7c7f2fSMark Adams } 146ef3f0257SMark Adams /* we know that the grid structure can be reused in MatPtAP */ 1472e3501ffSMark Adams } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */ 148192c0e8bSMark Adams PetscInt *counts, *newproc_idx, ii, jj, kk, strideNew, *tidx, ncrs_new, ncrs_eq_new, nloc_old, expand_factor = 1, rfactor = 1; 149885364a3SMark Adams IS is_eq_newproc, is_eq_num, is_eq_num_prim, new_eq_indices; 150849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 15171959b99SBarry Smith nloc_old = ncrs_eq / cr_bs; 15263a3b9bcSJacob 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); 153ce7c7f2fSMark Adams /* get new_size and rfactor */ 154ce7c7f2fSMark Adams if (pc_gamg->layout_type == PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) { 155ce7c7f2fSMark Adams /* find factor */ 156ce7c7f2fSMark Adams if (new_size == 1) rfactor = size; /* don't modify */ 157ce7c7f2fSMark Adams else { 158ce7c7f2fSMark Adams PetscReal best_fact = 0.; 159ce7c7f2fSMark Adams jj = -1; 160ce7c7f2fSMark Adams for (kk = 1; kk <= size; kk++) { 161ce7c7f2fSMark Adams if (!(size % kk)) { /* a candidate */ 162ce7c7f2fSMark Adams PetscReal nactpe = (PetscReal)size / (PetscReal)kk, fact = nactpe / (PetscReal)new_size; 163ce7c7f2fSMark Adams if (fact > 1.0) fact = 1. / fact; /* keep fact < 1 */ 164ce7c7f2fSMark Adams if (fact > best_fact) { 1659371c9d4SSatish Balay best_fact = fact; 1669371c9d4SSatish Balay jj = kk; 167ce7c7f2fSMark Adams } 168ce7c7f2fSMark Adams } 169ce7c7f2fSMark Adams } 170ce7c7f2fSMark Adams if (jj != -1) rfactor = jj; 171ce7c7f2fSMark Adams else rfactor = 1; /* a prime */ 172ce7c7f2fSMark Adams if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1; 173ce7c7f2fSMark Adams else expand_factor = rfactor; 174ce7c7f2fSMark Adams } 175ce7c7f2fSMark Adams new_size = size / rfactor; /* make new size one that is factor */ 1764cdfd227SMark if (new_size == nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */ 1774cdfd227SMark *a_Amat_crs = Cmat; 17863a3b9bcSJacob 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)); 179849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 180ce7c7f2fSMark Adams PetscFunctionReturn(0); 181ce7c7f2fSMark Adams } 182ce7c7f2fSMark Adams } 183a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 1849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &counts)); 185849bee69SMark Adams if (pc_gamg->repart) { /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */ 1865a9b9e01SMark F. Adams Mat adj; 187849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0)); 18863a3b9bcSJacob 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")); 189a2f3521dSMark F. Adams /* get 'adj' */ 190c5bfad50SMark F. Adams if (cr_bs == 1) { 1919566063dSJacob Faibussowitsch PetscCall(MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj)); 192806fa848SBarry Smith } else { 193a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 194eb07cef2SMark F. Adams Mat tMat; 195a2f3521dSMark F. Adams PetscInt Istart_crs, Iend_crs, ncols, jj, Ii; 196b4fbaa2aSMark F. Adams const PetscScalar *vals; 197b4fbaa2aSMark F. Adams const PetscInt *idx; 198a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 19939d09545SMark Adams static PetscInt llev = 0; /* ugly but just used for debugging */ 200d9558ea9SBarry Smith MatType mtype; 201b4fbaa2aSMark F. Adams 2029566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(ncrs, &d_nnz, ncrs, &o_nnz)); 2039566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs)); 2049566063dSJacob Faibussowitsch PetscCall(MatGetSize(Cmat, &M, &N)); 205c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 2069566063dSJacob Faibussowitsch PetscCall(MatGetRow(Cmat, Ii, &ncols, NULL, NULL)); 207c5bfad50SMark F. Adams d_nnz[jj] = ncols / cr_bs; 208c5bfad50SMark F. Adams o_nnz[jj] = ncols / cr_bs; 2099566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Cmat, Ii, &ncols, NULL, NULL)); 2103ae0bb68SMark Adams if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 2113ae0bb68SMark Adams if (o_nnz[jj] > (M / cr_bs - ncrs)) o_nnz[jj] = M / cr_bs - ncrs; 21258471d46SMark F. Adams } 2136876a03eSMark F. Adams 2149566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat_fine, &mtype)); 2159566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &tMat)); 2169566063dSJacob Faibussowitsch PetscCall(MatSetSizes(tMat, ncrs, ncrs, PETSC_DETERMINE, PETSC_DETERMINE)); 2179566063dSJacob Faibussowitsch PetscCall(MatSetType(tMat, mtype)); 2189566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(tMat, 0, d_nnz)); 2199566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(tMat, 0, d_nnz, 0, o_nnz)); 2209566063dSJacob Faibussowitsch PetscCall(PetscFree2(d_nnz, o_nnz)); 221eb07cef2SMark F. Adams 222a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 223c5bfad50SMark F. Adams PetscInt dest_row = ii / cr_bs; 2249566063dSJacob Faibussowitsch PetscCall(MatGetRow(Cmat, ii, &ncols, &idx, &vals)); 225eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 226c5bfad50SMark F. Adams PetscInt dest_col = idx[jj] / cr_bs; 227eb07cef2SMark F. Adams PetscScalar v = 1.0; 2289566063dSJacob Faibussowitsch PetscCall(MatSetValues(tMat, 1, &dest_row, 1, &dest_col, &v, ADD_VALUES)); 229eb07cef2SMark F. Adams } 2309566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Cmat, ii, &ncols, &idx, &vals)); 231eb07cef2SMark F. Adams } 2329566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tMat, MAT_FINAL_ASSEMBLY)); 2339566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tMat, MAT_FINAL_ASSEMBLY)); 234eb07cef2SMark F. Adams 235b4fbaa2aSMark F. Adams if (llev++ == -1) { 2369371c9d4SSatish Balay PetscViewer viewer; 2379371c9d4SSatish Balay char fname[32]; 23863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fname, sizeof(fname), "part_mat_%" PetscInt_FMT ".mat", llev)); 2393b4367a7SBarry Smith PetscViewerBinaryOpen(comm, fname, FILE_MODE_WRITE, &viewer); 2409566063dSJacob Faibussowitsch PetscCall(MatView(tMat, viewer)); 2419566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 242b4fbaa2aSMark F. Adams } 2439566063dSJacob Faibussowitsch PetscCall(MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj)); 2449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tMat)); 245a2f3521dSMark F. Adams } /* create 'adj' */ 246f150b916SMark F. Adams 247a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 2485a9b9e01SMark F. Adams char prefix[256]; 2495a9b9e01SMark F. Adams const char *pcpre; 250b4fbaa2aSMark F. Adams const PetscInt *is_idx; 251b4fbaa2aSMark F. Adams MatPartitioning mpart; 252a4b7d37bSMark F. Adams IS proc_is; 2532f03bc48SMark F. Adams 2549566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(comm, &mpart)); 2559566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(mpart, adj)); 2569566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 2579566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : "")); 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix)); 2599566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(mpart)); 2609566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(mpart, new_size)); 2619566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(mpart, &proc_is)); 2629566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&mpart)); 2635a9b9e01SMark F. Adams 2645ef31b24SMark F. Adams /* collect IS info */ 2659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncrs_eq, &newproc_idx)); 2669566063dSJacob Faibussowitsch PetscCall(ISGetIndices(proc_is, &is_idx)); 267a2f3521dSMark F. Adams for (kk = jj = 0; kk < nloc_old; kk++) { 2689371c9d4SSatish Balay for (ii = 0; ii < cr_bs; ii++, jj++) { newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ } 2695ef31b24SMark F. Adams } 2709566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(proc_is, &is_idx)); 2719566063dSJacob Faibussowitsch PetscCall(ISDestroy(&proc_is)); 2725ef31b24SMark F. Adams } 2739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 2745a9b9e01SMark F. Adams 2759566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFree(newproc_idx)); 277849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0)); 27831cb4603SMark Adams } else { /* simple aggregation of parts -- 'is_eq_newproc' */ 279ce7c7f2fSMark Adams PetscInt targetPE; 28008401ef6SPierre Jolivet PetscCheck(new_size != nactive, PETSC_COMM_SELF, PETSC_ERR_PLIB, "new_size==nactive. Should not happen"); 28163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Number of equations (loc) %" PetscInt_FMT " with simple aggregation\n", ((PetscObject)pc)->prefix, ncrs_eq)); 282ce7c7f2fSMark Adams targetPE = (rank / rfactor) * expand_factor; 2839566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc)); 284a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 285e33ef3b1SMark F. Adams 28611e60469SMark F. Adams /* 287a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 28811e60469SMark F. Adams */ 2899566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num)); 2907700e67bSMark Adams is_eq_num_prim = is_eq_num; 29111e60469SMark F. Adams /* 292a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 29311e60469SMark F. Adams */ 2949566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(is_eq_newproc, size, counts)); 295c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 2969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_eq_newproc)); 297ce7c7f2fSMark Adams ncrs_new = ncrs_eq_new / cr_bs; 298a2f3521dSMark F. Adams 2999566063dSJacob Faibussowitsch PetscCall(PetscFree(counts)); 3006aad120cSJose 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 */ 301885364a3SMark Adams { 302885364a3SMark Adams Vec src_crd, dest_crd; 303885364a3SMark 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; 304885364a3SMark Adams VecScatter vecscat; 305885364a3SMark Adams PetscScalar *array; 306885364a3SMark Adams IS isscat; 307a2f3521dSMark F. Adams /* move data (for primal equations only) */ 30822063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 3099566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &dest_crd)); 3109566063dSJacob Faibussowitsch PetscCall(VecSetSizes(dest_crd, node_data_sz * ncrs_new, PETSC_DECIDE)); 3119566063dSJacob Faibussowitsch PetscCall(VecSetType(dest_crd, VECSTANDARD)); /* this is needed! */ 31211e60469SMark F. Adams /* 3139d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 314c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 31511e60469SMark F. Adams */ 3169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncrs * node_data_sz, &tidx)); 3179566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is_eq_num_prim, &idx)); 3183ae0bb68SMark Adams for (ii = 0, jj = 0; ii < ncrs; ii++) { 319c5bfad50SMark F. Adams PetscInt id = idx[ii * cr_bs] / cr_bs; /* get node back */ 320a2f3521dSMark F. Adams for (kk = 0; kk < node_data_sz; kk++, jj++) tidx[jj] = id * node_data_sz + kk; 32111e60469SMark F. Adams } 3229566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is_eq_num_prim, &idx)); 3239566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, node_data_sz * ncrs, tidx, PETSC_COPY_VALUES, &isscat)); 3249566063dSJacob Faibussowitsch PetscCall(PetscFree(tidx)); 32511e60469SMark F. Adams /* 32611e60469SMark F. Adams Create a vector to contain the original vertex information for each element 32711e60469SMark F. Adams */ 3289566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, node_data_sz * ncrs, &src_crd)); 3299d5b6da9SMark F. Adams for (jj = 0; jj < ndata_cols; jj++) { 3303ae0bb68SMark Adams const PetscInt stride0 = ncrs * pc_gamg->data_cell_rows; 3313ae0bb68SMark Adams for (ii = 0; ii < ncrs; ii++) { 3329d5b6da9SMark F. Adams for (kk = 0; kk < ndata_rows; kk++) { 333a2f3521dSMark F. Adams PetscInt ix = ii * ndata_rows + kk + jj * stride0, jx = ii * node_data_sz + kk * ndata_cols + jj; 334c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 3359566063dSJacob Faibussowitsch PetscCall(VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES)); 336d3d6bff4SMark F. Adams } 337038e3b61SMark F. Adams } 338eb07cef2SMark F. Adams } 3399566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(src_crd)); 3409566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(src_crd)); 34111e60469SMark F. Adams /* 34211e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 34311e60469SMark F. Adams to the correct processor 34411e60469SMark F. Adams */ 3459566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat)); 3469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isscat)); 3479566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD)); 3489566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD)); 3499566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vecscat)); 3509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&src_crd)); 35111e60469SMark F. Adams /* 35211e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 35311e60469SMark F. Adams */ 3549566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 3559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(node_data_sz * ncrs_new, &pc_gamg->data)); 3562fa5cd67SKarl Rupp 3573ae0bb68SMark Adams pc_gamg->data_sz = node_data_sz * ncrs_new; 3583ae0bb68SMark Adams strideNew = ncrs_new * ndata_rows; 3592fa5cd67SKarl Rupp 3609566063dSJacob Faibussowitsch PetscCall(VecGetArray(dest_crd, &array)); 3619d5b6da9SMark F. Adams for (jj = 0; jj < ndata_cols; jj++) { 3623ae0bb68SMark Adams for (ii = 0; ii < ncrs_new; ii++) { 3639d5b6da9SMark F. Adams for (kk = 0; kk < ndata_rows; kk++) { 364a2f3521dSMark F. Adams PetscInt ix = ii * ndata_rows + kk + jj * strideNew, jx = ii * node_data_sz + kk * ndata_cols + jj; 365c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 366d3d6bff4SMark F. Adams } 367038e3b61SMark F. Adams } 368038e3b61SMark F. Adams } 3699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(dest_crd, &array)); 3709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dest_crd)); 371885364a3SMark Adams } 372a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 373849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0)); */ 37411e60469SMark F. Adams /* 3757dae84e0SHong Zhang Invert for MatCreateSubMatrix 37611e60469SMark F. Adams */ 3779566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices)); 3789566063dSJacob Faibussowitsch PetscCall(ISSort(new_eq_indices)); /* is this needed? */ 3799566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(new_eq_indices, cr_bs)); 3809371c9d4SSatish Balay if (is_eq_num != is_eq_num_prim) { PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */ } 3813cb8563fSToby Isaac if (Pcolumnperm) { 3829566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)new_eq_indices)); 3833cb8563fSToby Isaac *Pcolumnperm = new_eq_indices; 3843cb8563fSToby Isaac } 3859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_eq_num)); 386849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0)); */ 387849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 388849bee69SMark Adams 389a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 390a2f3521dSMark F. Adams { 391a2f3521dSMark F. Adams Mat mat; 392b94d7dedSBarry Smith PetscBool isset, isspd, isher; 39390db8557SMark Adams #if !defined(PETSC_USE_COMPLEX) 394b94d7dedSBarry Smith PetscBool issym; 395b94d7dedSBarry Smith #endif 396b94d7dedSBarry Smith 397b94d7dedSBarry Smith PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat)); 398b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(Cmat, &isset, &isspd)); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ? 399b94d7dedSBarry Smith if (isset) PetscCall(MatSetOption(mat, MAT_SPD, isspd)); 400b94d7dedSBarry Smith else { 401b94d7dedSBarry Smith PetscCall(MatIsHermitianKnown(Cmat, &isset, &isher)); 402b94d7dedSBarry Smith if (isset) PetscCall(MatSetOption(mat, MAT_HERMITIAN, isher)); 403b94d7dedSBarry Smith else { 404b94d7dedSBarry Smith #if !defined(PETSC_USE_COMPLEX) 405b94d7dedSBarry Smith PetscCall(MatIsSymmetricKnown(Cmat, &isset, &issym)); 406b94d7dedSBarry Smith if (isset) PetscCall(MatSetOption(mat, MAT_SYMMETRIC, issym)); 40790db8557SMark Adams #endif 40890db8557SMark Adams } 40990db8557SMark Adams } 410a2f3521dSMark F. Adams *a_Amat_crs = mat; 411a2f3521dSMark F. Adams } 4129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cmat)); 413a2f3521dSMark F. Adams 414849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 41511e60469SMark F. Adams /* prolongator */ 41611e60469SMark F. Adams { 41711e60469SMark F. Adams IS findices; 418a2f3521dSMark F. Adams PetscInt Istart, Iend; 419a2f3521dSMark F. Adams Mat Pnew; 42062294041SBarry Smith 4219566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend)); 422849bee69SMark Adams /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 4239566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, Iend - Istart, Istart, 1, &findices)); 4249566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(findices, f_bs)); 4259566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew)); 4269566063dSJacob Faibussowitsch PetscCall(ISDestroy(&findices)); 4279566063dSJacob Faibussowitsch PetscCall(MatSetOption(Pnew, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 428c5bfad50SMark F. Adams 429849bee69SMark Adams /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 4309566063dSJacob Faibussowitsch PetscCall(MatDestroy(a_P_inout)); 431a2f3521dSMark F. Adams 432a2f3521dSMark F. Adams /* output - repartitioned */ 433a2f3521dSMark F. Adams *a_P_inout = Pnew; 434e33ef3b1SMark F. Adams } 4359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&new_eq_indices)); 4365b89ad90SMark F. Adams 437c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 438ce7c7f2fSMark Adams 439ce7c7f2fSMark Adams /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */ 440ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 441ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 4428bca76a6SMark Adams static PetscInt llev = 2; 44363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Pinning level %" PetscInt_FMT " to the CPU\n", ((PetscObject)pc)->prefix, llev++)); 444ce7c7f2fSMark Adams #endif 4459566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE)); 4469566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE)); 447adf5291fSStefano Zampini if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */ 448ce7c7f2fSMark Adams Mat A = *a_Amat_crs, P = *a_P_inout; 449ce7c7f2fSMark Adams PetscMPIInt size; 4509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 451ce7c7f2fSMark Adams if (size > 1) { 452ce7c7f2fSMark Adams Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data, *p = (Mat_MPIAIJ *)P->data; 4539566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(a->lvec, PETSC_TRUE)); 4549566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(p->lvec, PETSC_TRUE)); 455ce7c7f2fSMark Adams } 456ce7c7f2fSMark Adams } 457ce7c7f2fSMark Adams } 458849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 459849bee69SMark Adams } // processor reduce 4605b89ad90SMark F. Adams PetscFunctionReturn(0); 4615b89ad90SMark F. Adams } 4625b89ad90SMark F. Adams 463bae903cbSmarkadams4 // used in GEO 4649371c9d4SSatish Balay PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat *Gmat2) { 4654b1575e2SStefano Zampini const char *prefix; 4664b1575e2SStefano Zampini char addp[32]; 4674b1575e2SStefano Zampini PC_MG *mg = (PC_MG *)a_pc->data; 4684b1575e2SStefano Zampini PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 4694b1575e2SStefano Zampini 4704b1575e2SStefano Zampini PetscFunctionBegin; 4719566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(a_pc, &prefix)); 47263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(a_pc, "%s: Square Graph on level %" PetscInt_FMT "\n", ((PetscObject)a_pc)->prefix, pc_gamg->current_level + 1)); 4739566063dSJacob Faibussowitsch PetscCall(MatProductCreate(Gmat1, Gmat1, NULL, Gmat2)); 4749566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(*Gmat2, prefix)); 47563a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_square_%" PetscInt_FMT "_", pc_gamg->current_level)); 4769566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(*Gmat2, addp)); 477b94d7dedSBarry Smith if ((*Gmat2)->structurally_symmetric == PETSC_BOOL3_TRUE) { 4789566063dSJacob Faibussowitsch PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AB)); 479b4da3a1bSStefano Zampini } else { 4809566063dSJacob Faibussowitsch PetscCall(MatSetOption(Gmat1, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 4819566063dSJacob Faibussowitsch PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AtB)); 482b4da3a1bSStefano Zampini } 4839566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(*Gmat2)); 4849566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0)); 4859566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(*Gmat2)); 4869566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0)); 4879566063dSJacob Faibussowitsch PetscCall(MatProductClear(*Gmat2)); 4884b1575e2SStefano Zampini /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */ 4894b1575e2SStefano Zampini (*Gmat2)->assembled = PETSC_TRUE; 4904b1575e2SStefano Zampini PetscFunctionReturn(0); 4914b1575e2SStefano Zampini } 4924b1575e2SStefano Zampini 4935b89ad90SMark F. Adams /* 4945b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4955b89ad90SMark F. Adams by setting data structures and options. 4965b89ad90SMark F. Adams 4975b89ad90SMark F. Adams Input Parameter: 4985b89ad90SMark F. Adams . pc - the preconditioner context 4995b89ad90SMark F. Adams 5005b89ad90SMark F. Adams */ 5019371c9d4SSatish Balay PetscErrorCode PCSetUp_GAMG(PC pc) { 5029d5b6da9SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 5035b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 5042adcac29SMark F. Adams Mat Pmat = pc->pmat; 50518c3aa7eSMark PetscInt fine_level, level, level1, bs, M, N, qq, lidx, nASMBlocksArr[PETSC_MG_MAXLEVELS]; 5063b4367a7SBarry Smith MPI_Comm comm; 507c5df96a5SBarry Smith PetscMPIInt rank, size, nactivepe; 50818c3aa7eSMark Mat Aarr[PETSC_MG_MAXLEVELS], Parr[PETSC_MG_MAXLEVELS]; 50918c3aa7eSMark IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 510a2f3521dSMark F. Adams PetscLogDouble nnz0 = 0., nnztot = 0.; 511569f4572SMark Adams MatInfo info; 512171cca9aSMark Adams PetscBool is_last = PETSC_FALSE; 5135ef31b24SMark F. Adams 5145b89ad90SMark F. Adams PetscFunctionBegin; 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 5169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 518849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 5198abdc6daSStefano Zampini if (pc->setupcalled) { 5208abdc6daSStefano Zampini if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) { 521878e152fSMark F. Adams /* reset everything */ 5229566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 523878e152fSMark F. Adams pc->setupcalled = 0; 524806fa848SBarry Smith } else { 52584d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 52603a628feSMark F. Adams /* just do Galerkin grids */ 52758471d46SMark F. Adams Mat B, dA, dB; 5289d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 5294555aa8cSStefano Zampini PetscInt gl; 53058471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 5319566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, &dA, &dB)); 53258471d46SMark F. Adams /* (re)set to get dirty flag */ 5339566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, dA, dB)); 53458471d46SMark F. Adams 5354555aa8cSStefano Zampini for (level = pc_gamg->Nlevels - 2, gl = 0; level >= 0; level--, gl++) { 5368abdc6daSStefano Zampini MatReuse reuse = MAT_INITIAL_MATRIX; 537849bee69SMark Adams #if defined(GAMG_STAGES) 538849bee69SMark Adams PetscCall(PetscLogStagePush(gamg_stages[gl])); 539849bee69SMark Adams #endif 5408abdc6daSStefano Zampini /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 5419566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[level]->smoothd, NULL, &B)); 5428abdc6daSStefano Zampini if (B->product) { 543ad540459SPierre Jolivet if (B->product->A == dB && B->product->B == mglevels[level + 1]->interpolate) reuse = MAT_REUSE_MATRIX; 5448abdc6daSStefano Zampini } 5459566063dSJacob Faibussowitsch if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A)); 5468abdc6daSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 54763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: RAP after first solve, reuse matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 5488abdc6daSStefano Zampini } else { 54963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: RAP after first solve, new matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 5508abdc6daSStefano Zampini } 5519566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0)); 5529566063dSJacob Faibussowitsch PetscCall(MatPtAP(dB, mglevels[level + 1]->interpolate, reuse, PETSC_DEFAULT, &B)); 5539566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0)); 55463b77682SMark Adams if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B; 5559566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[level]->smoothd, B, B)); 55658471d46SMark F. Adams dB = B; 557849bee69SMark Adams #if defined(GAMG_STAGES) 558849bee69SMark Adams PetscCall(PetscLogStagePop()); 559849bee69SMark Adams #endif 56058471d46SMark F. Adams } 5615f8cf99dSMark F. Adams } 562d5280255SMark F. Adams 5639566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 564849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 56558471d46SMark F. Adams PetscFunctionReturn(0); 566eb07cef2SMark F. Adams } 567878e152fSMark F. Adams } 568f6536408SMark F. Adams 569878e152fSMark F. Adams if (!pc_gamg->data) { 570878e152fSMark F. Adams if (pc_gamg->orig_data) { 5719566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Pmat, &bs)); 5729566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Pmat, &qq, NULL)); 5732fa5cd67SKarl Rupp 574878e152fSMark F. Adams pc_gamg->data_sz = (qq / bs) * pc_gamg->orig_data_cell_rows * pc_gamg->orig_data_cell_cols; 575878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 576878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 5772fa5cd67SKarl Rupp 5789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data)); 579878e152fSMark F. Adams for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 580806fa848SBarry Smith } else { 5815f80ce2aSJacob Faibussowitsch PetscCheck(pc_gamg->ops->createdefaultdata, comm, PETSC_ERR_PLIB, "'createdefaultdata' not set(?) need to support NULL data"); 5829566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->createdefaultdata(pc, Pmat)); 5839d5b6da9SMark F. Adams } 584878e152fSMark F. Adams } 585878e152fSMark F. Adams 586878e152fSMark F. Adams /* cache original data for reuse */ 5871c1aac46SBarry Smith if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 5889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data)); 589878e152fSMark F. Adams for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 590878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 591878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 592878e152fSMark F. Adams } 593038e3b61SMark F. Adams 594302f38e8SMark F. Adams /* get basic dims */ 5959566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Pmat, &bs)); 5969566063dSJacob Faibussowitsch PetscCall(MatGetSize(Pmat, &M, &N)); 59784d3f75bSMark F. Adams 5989566063dSJacob Faibussowitsch PetscCall(MatGetInfo(Pmat, MAT_GLOBAL_SUM, &info)); /* global reduction */ 599569f4572SMark Adams nnz0 = info.nz_used; 600569f4572SMark Adams nnztot = info.nz_used; 601bae903cbSmarkadams4 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)); 602569f4572SMark Adams 603a2f3521dSMark F. Adams /* Get A_i and R_i */ 60462294041SBarry Smith for (level = 0, Aarr[0] = Pmat, nactivepe = size; level < (pc_gamg->Nlevels - 1) && (!level || M > pc_gamg->coarse_eq_limit); level++) { 6059ab59c8bSMark Adams pc_gamg->current_level = level; 60663a3b9bcSJacob Faibussowitsch PetscCheck(level < PETSC_MG_MAXLEVELS, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Too many levels %" PetscInt_FMT, level); 6075b89ad90SMark F. Adams level1 = level + 1; 6084555aa8cSStefano Zampini #if defined(GAMG_STAGES) 609849bee69SMark Adams if (!gamg_stages[level]) { 610849bee69SMark Adams char str[32]; 611849bee69SMark Adams sprintf(str, "GAMG Level %d", (int)level); 612849bee69SMark Adams PetscCall(PetscLogStageRegister(str, &gamg_stages[level])); 613849bee69SMark Adams } 6149566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(gamg_stages[level])); 615b4fbaa2aSMark F. Adams #endif 616c8b0795cSMark F. Adams { /* construct prolongator */ 617725b86d8SJed Brown Mat Gmat; 6180cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 6197700e67bSMark Adams Mat Prol11; 620c8b0795cSMark F. Adams 6219566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->graph(pc, Aarr[level], &Gmat)); 6229566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists)); 6239566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->prolongator(pc, Aarr[level], Gmat, agg_lists, &Prol11)); 624c8b0795cSMark F. Adams 625a2f3521dSMark F. Adams /* could have failed to create new level */ 626a2f3521dSMark F. Adams if (Prol11) { 627f7df55f0SStefano Zampini const char *prefix; 628f7df55f0SStefano Zampini char addp[32]; 629f7df55f0SStefano Zampini 6309d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 6319566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(Prol11, NULL, &bs)); 632a2f3521dSMark F. Adams 633fd1112cbSBarry Smith if (pc_gamg->ops->optprolongator) { 634c8b0795cSMark F. Adams /* smooth */ 6359566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11)); 636c8b0795cSMark F. Adams } 637c8b0795cSMark F. Adams 6380c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 6391b18a24aSMark Adams PetscInt bs; 640849bee69SMark Adams PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // not timed directly, ugly, could remove, but good ASM method 6419566063dSJacob Faibussowitsch PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level])); 642ffc955d6SMark F. Adams } 643ffc955d6SMark F. Adams 6449566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 6459566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(Prol11, prefix)); 6469566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_prolongator_%d_", (int)level)); 6479566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(Prol11, addp)); 64891f31d3dSStefano Zampini /* Always generate the transpose with CUDA 649f7df55f0SStefano Zampini Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */ 6509566063dSJacob Faibussowitsch PetscCall(MatSetOption(Prol11, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 6519566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Prol11)); 6524bde40a0SMark Adams Parr[level1] = Prol11; 6534bde40a0SMark Adams } else Parr[level1] = NULL; /* failed to coarsen */ 6544bde40a0SMark Adams 6559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Gmat)); 6569566063dSJacob Faibussowitsch PetscCall(PetscCDDestroy(agg_lists)); 657a2f3521dSMark F. Adams } /* construct prolongator scope */ 6587f66b68fSMark Adams if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 659171cca9aSMark Adams if (!Parr[level1]) { /* failed to coarsen */ 66063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Stop gridding, level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 6614555aa8cSStefano Zampini #if defined(GAMG_STAGES) 6629566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 663a90e85d9SMark Adams #endif 664c8b0795cSMark F. Adams break; 665c8b0795cSMark F. Adams } 6669566063dSJacob Faibussowitsch PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */ 6672472a847SBarry Smith PetscCheck(!is_last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Is last ?"); 668171cca9aSMark Adams if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 6690e2909e1SMark Adams if (level1 == pc_gamg->Nlevels - 1) is_last = PETSC_TRUE; 670849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0)); 6719566063dSJacob Faibussowitsch PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last)); 672849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0)); 673a2f3521dSMark F. Adams 6749566063dSJacob Faibussowitsch PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */ 6759566063dSJacob Faibussowitsch PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info)); 676569f4572SMark Adams nnztot += info.nz_used; 677bae903cbSmarkadams4 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)); 678569f4572SMark Adams 6794555aa8cSStefano Zampini #if defined(GAMG_STAGES) 6809566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 681b4fbaa2aSMark F. Adams #endif 682a90e85d9SMark Adams /* stop if one node or one proc -- could pull back for singular problems */ 6839ab59c8bSMark Adams if ((pc_gamg->data_cell_cols && M / pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M / bs < 2)) { 68463a3b9bcSJacob 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)); 685a90e85d9SMark Adams level++; 686a90e85d9SMark Adams break; 687a90e85d9SMark Adams } 688c8b0795cSMark F. Adams } /* levels */ 6899566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 690c8b0795cSMark F. Adams 691ba4ceb06Smarkadams4 PetscCall(PetscInfo(pc, "%s: %" PetscInt_FMT " levels, operator complexity = %g\n", ((PetscObject)pc)->prefix, level + 1, nnztot / nnz0)); 6929d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 6935b89ad90SMark F. Adams fine_level = level; 6949566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, pc_gamg->Nlevels, NULL)); 6955b89ad90SMark F. Adams 69662294041SBarry Smith if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 6970ed2132dSStefano Zampini 698d5280255SMark F. Adams /* set default smoothers & set operators */ 69962294041SBarry Smith for (lidx = 1, level = pc_gamg->Nlevels - 2; lidx <= fine_level; lidx++, level--) { 700ffc955d6SMark F. Adams KSP smoother; 701ffc955d6SMark F. Adams PC subpc; 702a2f3521dSMark F. Adams 7039566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7049566063dSJacob Faibussowitsch PetscCall(KSPGetPC(smoother, &subpc)); 705ffc955d6SMark F. Adams 7069566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 707a2f3521dSMark F. Adams /* set ops */ 7089566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level])); 7099566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level + 1])); 710a2f3521dSMark F. Adams 711a2f3521dSMark F. Adams /* set defaults */ 7129566063dSJacob Faibussowitsch PetscCall(KSPSetType(smoother, KSPCHEBYSHEV)); 713a2f3521dSMark F. Adams 7140c3bc534SBarry Smith /* set blocks for ASM smoother that uses the 'aggregates' */ 7150c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 7162d3561bbSSatish Balay PetscInt sz; 7177a28f3e5SMark Adams IS *iss; 718a2f3521dSMark F. Adams 7192d3561bbSSatish Balay sz = nASMBlocksArr[level]; 7207a28f3e5SMark Adams iss = ASMLocalIDsArr[level]; 7219566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCASM)); 7229566063dSJacob Faibussowitsch PetscCall(PCASMSetOverlap(subpc, 0)); 7239566063dSJacob Faibussowitsch PetscCall(PCASMSetType(subpc, PC_ASM_BASIC)); 7247f66b68fSMark Adams if (!sz) { 725ffc955d6SMark F. Adams IS is; 7269566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is)); 7279566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is)); 7289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 729806fa848SBarry Smith } else { 730a94c3b12SMark F. Adams PetscInt kk; 7319566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(subpc, sz, NULL, iss)); 73248a46eb9SPierre Jolivet for (kk = 0; kk < sz; kk++) PetscCall(ISDestroy(&iss[kk])); 7339566063dSJacob Faibussowitsch PetscCall(PetscFree(iss)); 734ffc955d6SMark F. Adams } 7350298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 736ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 737806fa848SBarry Smith } else { 7389566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCJACOBI)); 739ffc955d6SMark F. Adams } 740d5280255SMark F. Adams } 741d5280255SMark F. Adams { 742d5280255SMark F. Adams /* coarse grid */ 7439371c9d4SSatish Balay KSP smoother, *k2; 7449371c9d4SSatish Balay PC subpc, pc2; 7459371c9d4SSatish Balay PetscInt ii, first; 7469371c9d4SSatish Balay Mat Lmat = Aarr[(level = pc_gamg->Nlevels - 1)]; 7479371c9d4SSatish Balay lidx = 0; 7480ed2132dSStefano Zampini 7499566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7509566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Lmat, Lmat)); 751cf8ae1d3SMark Adams if (!pc_gamg->use_parallel_coarse_grid_solver) { 7529566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 7539566063dSJacob Faibussowitsch PetscCall(KSPGetPC(smoother, &subpc)); 7549566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, PCBJACOBI)); 7559566063dSJacob Faibussowitsch PetscCall(PCSetUp(subpc)); 7569566063dSJacob Faibussowitsch PetscCall(PCBJacobiGetSubKSP(subpc, &ii, &first, &k2)); 75763a3b9bcSJacob Faibussowitsch PetscCheck(ii == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " is not one", ii); 7589566063dSJacob Faibussowitsch PetscCall(KSPGetPC(k2[0], &pc2)); 7599566063dSJacob Faibussowitsch PetscCall(PCSetType(pc2, PCLU)); 7609566063dSJacob Faibussowitsch PetscCall(PCFactorSetShiftType(pc2, MAT_SHIFT_INBLOCKS)); 7619566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(k2[0], PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1)); 7629566063dSJacob Faibussowitsch PetscCall(KSPSetType(k2[0], KSPPREONLY)); 763d5280255SMark F. Adams } 764cf8ae1d3SMark Adams } 765d5280255SMark F. Adams 766d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 767d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)pc); 768dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); 769d0609cedSBarry Smith PetscOptionsEnd(); 7709566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 771d5280255SMark F. Adams 772*f1580f4eSBarry Smith /* set cheby eigen estimates from SA to use in the solver */ 7737e6512fdSJed Brown if (pc_gamg->use_sa_esteig) { 77418c3aa7eSMark for (lidx = 1, level = pc_gamg->Nlevels - 2; level >= 0; lidx++, level--) { 77518c3aa7eSMark KSP smoother; 77618c3aa7eSMark PetscBool ischeb; 7770ed2132dSStefano Zampini 7789566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 7799566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb)); 78018c3aa7eSMark if (ischeb) { 78118c3aa7eSMark KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data; 7820ed2132dSStefano Zampini 7832de708cbSJed Brown // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG 784efe053fcSJed Brown if (mg->max_eigen_DinvA[level] > 0) { 7857e6512fdSJed Brown // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it. 7867e6512fdSJed Brown // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.) 78718c3aa7eSMark PetscReal emax, emin; 7880ed2132dSStefano Zampini 78918c3aa7eSMark emin = mg->min_eigen_DinvA[level]; 79018c3aa7eSMark emax = mg->max_eigen_DinvA[level]; 79163a3b9bcSJacob 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)); 7920a94a983SJed Brown cheb->emin_provided = emin; 7930a94a983SJed Brown cheb->emax_provided = emax; 79418c3aa7eSMark } 79518c3aa7eSMark } 79618c3aa7eSMark } 79718c3aa7eSMark } 7980ed2132dSStefano Zampini 7999566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 8000ed2132dSStefano Zampini 801d5280255SMark F. Adams /* clean up */ 802d5280255SMark F. Adams for (level = 1; level < pc_gamg->Nlevels; level++) { 8039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Parr[level])); 8049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aarr[level])); 8055b89ad90SMark F. Adams } 806806fa848SBarry Smith } else { 8075f8cf99dSMark F. Adams KSP smoother; 8080ed2132dSStefano Zampini 80963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: One level solver used (system is seen as DD). Using default solver.\n", ((PetscObject)pc)->prefix)); 8109566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, 0, &smoother)); 8119566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0])); 8129566063dSJacob Faibussowitsch PetscCall(KSPSetType(smoother, KSPPREONLY)); 8139566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 8145f8cf99dSMark F. Adams } 815849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 8165b89ad90SMark F. Adams PetscFunctionReturn(0); 8175b89ad90SMark F. Adams } 8185b89ad90SMark F. Adams 8195b89ad90SMark F. Adams /* 8205b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 8215b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 8225b89ad90SMark F. Adams 8235b89ad90SMark F. Adams Input Parameter: 8245b89ad90SMark F. Adams . pc - the preconditioner context 8255b89ad90SMark F. Adams 8265b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 8275b89ad90SMark F. Adams */ 8289371c9d4SSatish Balay PetscErrorCode PCDestroy_GAMG(PC pc) { 8295b89ad90SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 8305b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 8315b89ad90SMark F. Adams 8325b89ad90SMark F. Adams PetscFunctionBegin; 8339566063dSJacob Faibussowitsch PetscCall(PCReset_GAMG(pc)); 8341baa6e33SBarry Smith if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc)); 8359566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->ops)); 8369566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->gamg_type_name)); 8379566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg)); 8382e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", NULL)); 8392e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", NULL)); 8402e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", NULL)); 8412e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", NULL)); 8422e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", NULL)); 8432e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", NULL)); 8442e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", NULL)); 8452e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", NULL)); 8462e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", NULL)); 8472e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", NULL)); 8482e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", NULL)); 8492e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", NULL)); 8502e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", NULL)); 8512e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", NULL)); 8522e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", NULL)); 8532e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", NULL)); 8549566063dSJacob Faibussowitsch PetscCall(PCDestroy_MG(pc)); 8555b89ad90SMark F. Adams PetscFunctionReturn(0); 8565b89ad90SMark F. Adams } 8575b89ad90SMark F. Adams 858676e1743SMark F. Adams /*@ 859*f1580f4eSBarry Smith PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction in `PCGAMG` 860676e1743SMark F. Adams 861*f1580f4eSBarry Smith Logically Collective on pc 862676e1743SMark F. Adams 863676e1743SMark F. Adams Input Parameters: 8641cc46a46SBarry Smith + pc - the preconditioner context 8651cc46a46SBarry Smith - n - the number of equations 866676e1743SMark F. Adams 867676e1743SMark F. Adams Options Database Key: 868147403d9SBarry Smith . -pc_gamg_process_eq_limit <limit> - set the limit 869676e1743SMark F. Adams 870*f1580f4eSBarry Smith Note: 871*f1580f4eSBarry 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 872cab9ed1eSBarry Smith that has degrees of freedom 873cab9ed1eSBarry Smith 874676e1743SMark F. Adams Level: intermediate 875676e1743SMark F. Adams 876*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()` 877676e1743SMark F. Adams @*/ 8789371c9d4SSatish Balay PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) { 879676e1743SMark F. Adams PetscFunctionBegin; 880676e1743SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 881cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetProcEqLim_C", (PC, PetscInt), (pc, n)); 882676e1743SMark F. Adams PetscFunctionReturn(0); 883676e1743SMark F. Adams } 884676e1743SMark F. Adams 8859371c9d4SSatish Balay static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) { 886c20e4228SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 887c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 888676e1743SMark F. Adams 889676e1743SMark F. Adams PetscFunctionBegin; 8909d5b6da9SMark F. Adams if (n > 0) pc_gamg->min_eq_proc = n; 891676e1743SMark F. Adams PetscFunctionReturn(0); 892676e1743SMark F. Adams } 893676e1743SMark F. Adams 894389730f3SMark F. Adams /*@ 895*f1580f4eSBarry Smith PCGAMGSetCoarseEqLim - Set maximum number of equations on the coarsest grid of `PCGAMG` 896389730f3SMark F. Adams 897*f1580f4eSBarry Smith Collective on pc 898389730f3SMark F. Adams 899389730f3SMark F. Adams Input Parameters: 9001cc46a46SBarry Smith + pc - the preconditioner context 9011cc46a46SBarry Smith - n - maximum number of equations to aim for 902389730f3SMark F. Adams 903389730f3SMark F. Adams Options Database Key: 904147403d9SBarry Smith . -pc_gamg_coarse_eq_limit <limit> - set the limit 905389730f3SMark F. Adams 906*f1580f4eSBarry Smith Note: 907*f1580f4eSBarry Smith For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 90874329af1SBarry Smith has less than 1000 unknowns. 90974329af1SBarry Smith 910389730f3SMark F. Adams Level: intermediate 911389730f3SMark F. Adams 912*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()` 913389730f3SMark F. Adams @*/ 9149371c9d4SSatish Balay PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) { 915389730f3SMark F. Adams PetscFunctionBegin; 916389730f3SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 917cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetCoarseEqLim_C", (PC, PetscInt), (pc, n)); 918389730f3SMark F. Adams PetscFunctionReturn(0); 919389730f3SMark F. Adams } 920389730f3SMark F. Adams 9219371c9d4SSatish Balay static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) { 922389730f3SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 923389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 924389730f3SMark F. Adams 925389730f3SMark F. Adams PetscFunctionBegin; 9269d5b6da9SMark F. Adams if (n > 0) pc_gamg->coarse_eq_limit = n; 927389730f3SMark F. Adams PetscFunctionReturn(0); 928389730f3SMark F. Adams } 929389730f3SMark F. Adams 930676e1743SMark F. Adams /*@ 931*f1580f4eSBarry Smith PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids when reducing the number of MPI ranks to use 932676e1743SMark F. Adams 933*f1580f4eSBarry Smith Collective on pc 934676e1743SMark F. Adams 935676e1743SMark F. Adams Input Parameters: 9361cc46a46SBarry Smith + pc - the preconditioner context 937*f1580f4eSBarry Smith - n - `PETSC_TRUE` or `PETSC_FALSE` 938676e1743SMark F. Adams 939676e1743SMark F. Adams Options Database Key: 940147403d9SBarry Smith . -pc_gamg_repartition <true,false> - turn on the repartitioning 941676e1743SMark F. Adams 942*f1580f4eSBarry Smith Note: 943*f1580f4eSBarry Smith This will generally improve the loading balancing of the work on each level 944cab9ed1eSBarry Smith 945676e1743SMark F. Adams Level: intermediate 946676e1743SMark F. Adams 947*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()` 948676e1743SMark F. Adams @*/ 9499371c9d4SSatish Balay PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) { 950676e1743SMark F. Adams PetscFunctionBegin; 951676e1743SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 952cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetRepartition_C", (PC, PetscBool), (pc, n)); 953676e1743SMark F. Adams PetscFunctionReturn(0); 954676e1743SMark F. Adams } 955676e1743SMark F. Adams 9569371c9d4SSatish Balay static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) { 957c20e4228SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 958c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 959676e1743SMark F. Adams 960676e1743SMark F. Adams PetscFunctionBegin; 9619d5b6da9SMark F. Adams pc_gamg->repart = n; 962676e1743SMark F. Adams PetscFunctionReturn(0); 963676e1743SMark F. Adams } 964676e1743SMark F. Adams 965dfd5c07aSMark F. Adams /*@ 966*f1580f4eSBarry Smith PCGAMGSetUseSAEstEig - Use the eigen estimate from smoothed aggregation for the Chebyshev smoother during the solution process 96718c3aa7eSMark 968*f1580f4eSBarry Smith Collective on pc 96918c3aa7eSMark 97018c3aa7eSMark Input Parameters: 97118c3aa7eSMark + pc - the preconditioner context 97218c3aa7eSMark - n - number of its 97318c3aa7eSMark 97418c3aa7eSMark Options Database Key: 975147403d9SBarry Smith . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate 97618c3aa7eSMark 97718c3aa7eSMark Notes: 9787e6512fdSJed 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$. 979*f1580f4eSBarry Smith Eigenvalue estimates (based on a few `PCCG` or `PCGMRES` iterations) are computed to choose $\omega$ so that this is a stable smoothing operation. 980*f1580f4eSBarry Smith If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused during the solution process 981*f1580f4eSBarry Smith This option is only used when the smoother uses Jacobi, and should be turned off if a different `PCJacobiType` is used. 982efe053fcSJed Brown It became default in PETSc 3.17. 98318c3aa7eSMark 9847e6512fdSJed Brown Level: advanced 98518c3aa7eSMark 986*f1580f4eSBarry Smith .seealso: `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()` 98718c3aa7eSMark @*/ 9889371c9d4SSatish Balay PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n) { 98918c3aa7eSMark PetscFunctionBegin; 99018c3aa7eSMark PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 991cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetUseSAEstEig_C", (PC, PetscBool), (pc, n)); 99218c3aa7eSMark PetscFunctionReturn(0); 99318c3aa7eSMark } 99418c3aa7eSMark 9959371c9d4SSatish Balay static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n) { 99618c3aa7eSMark PC_MG *mg = (PC_MG *)pc->data; 99718c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 99818c3aa7eSMark 99918c3aa7eSMark PetscFunctionBegin; 10007e6512fdSJed Brown pc_gamg->use_sa_esteig = n; 100118c3aa7eSMark PetscFunctionReturn(0); 100218c3aa7eSMark } 100318c3aa7eSMark 100418c3aa7eSMark /*@ 1005*f1580f4eSBarry Smith PCGAMGSetEigenvalues - Set WHAT eigenvalues WHY? 100618c3aa7eSMark 1007*f1580f4eSBarry Smith Collective on pc 100818c3aa7eSMark 100918c3aa7eSMark Input Parameters: 101018c3aa7eSMark + pc - the preconditioner context 101118c3aa7eSMark - emax - max eigenvalue 101218c3aa7eSMark . emin - min eigenvalue 101318c3aa7eSMark 101418c3aa7eSMark Options Database Key: 1015147403d9SBarry Smith . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues 101618c3aa7eSMark 101718c3aa7eSMark Level: intermediate 101818c3aa7eSMark 1019*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetUseSAEstEig()` 102018c3aa7eSMark @*/ 10219371c9d4SSatish Balay PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax, PetscReal emin) { 102218c3aa7eSMark PetscFunctionBegin; 102318c3aa7eSMark PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1024cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetEigenvalues_C", (PC, PetscReal, PetscReal), (pc, emax, emin)); 102518c3aa7eSMark PetscFunctionReturn(0); 102618c3aa7eSMark } 102741ffd417SStefano Zampini 10289371c9d4SSatish Balay static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc, PetscReal emax, PetscReal emin) { 102918c3aa7eSMark PC_MG *mg = (PC_MG *)pc->data; 103018c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 103118c3aa7eSMark 103218c3aa7eSMark PetscFunctionBegin; 10335f80ce2aSJacob 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); 10345f80ce2aSJacob 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); 103518c3aa7eSMark pc_gamg->emax = emax; 103618c3aa7eSMark pc_gamg->emin = emin; 103718c3aa7eSMark PetscFunctionReturn(0); 103818c3aa7eSMark } 103918c3aa7eSMark 104018c3aa7eSMark /*@ 1041*f1580f4eSBarry Smith PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding a `PCGAMG` algebraic multigrid preconditioner 1042dfd5c07aSMark F. Adams 1043*f1580f4eSBarry Smith Collective on pc 1044dfd5c07aSMark F. Adams 1045dfd5c07aSMark F. Adams Input Parameters: 10461cc46a46SBarry Smith + pc - the preconditioner context 1047*f1580f4eSBarry Smith - n - `PETSC_TRUE` or `PETSC_FALSE` 1048dfd5c07aSMark F. Adams 1049dfd5c07aSMark F. Adams Options Database Key: 1050147403d9SBarry Smith . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation 1051dfd5c07aSMark F. Adams 1052dfd5c07aSMark F. Adams Level: intermediate 1053dfd5c07aSMark F. Adams 1054*f1580f4eSBarry Smith Note: 1055147403d9SBarry Smith May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1056cab9ed1eSBarry Smith rebuilding the preconditioner quicker. 1057cab9ed1eSBarry Smith 1058*f1580f4eSBarry Smith .seealso: `PCGAMG` 1059dfd5c07aSMark F. Adams @*/ 10609371c9d4SSatish Balay PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) { 1061dfd5c07aSMark F. Adams PetscFunctionBegin; 1062dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1063cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetReuseInterpolation_C", (PC, PetscBool), (pc, n)); 1064dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1065dfd5c07aSMark F. Adams } 1066dfd5c07aSMark F. Adams 10679371c9d4SSatish Balay static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) { 1068dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1069dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1070dfd5c07aSMark F. Adams 1071dfd5c07aSMark F. Adams PetscFunctionBegin; 1072dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 1073dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1074dfd5c07aSMark F. Adams } 1075dfd5c07aSMark F. Adams 1076ffc955d6SMark F. Adams /*@ 1077*f1580f4eSBarry 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 1078*f1580f4eSBarry Smith used as the smoother 1079ffc955d6SMark F. Adams 1080*f1580f4eSBarry Smith Collective on pc 1081ffc955d6SMark F. Adams 1082ffc955d6SMark F. Adams Input Parameters: 1083cab9ed1eSBarry Smith + pc - the preconditioner context 1084*f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use aggregates, `PETSC_FALSE` to not 1085ffc955d6SMark F. Adams 1086ffc955d6SMark F. Adams Options Database Key: 1087147403d9SBarry Smith . -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains 1088ffc955d6SMark F. Adams 1089ffc955d6SMark F. Adams Level: intermediate 1090ffc955d6SMark F. Adams 1091*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCASM`, `PCSetType` 1092ffc955d6SMark F. Adams @*/ 10939371c9d4SSatish Balay PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) { 1094ffc955d6SMark F. Adams PetscFunctionBegin; 1095ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1096cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGASMSetUseAggs_C", (PC, PetscBool), (pc, flg)); 1097ffc955d6SMark F. Adams PetscFunctionReturn(0); 1098ffc955d6SMark F. Adams } 1099ffc955d6SMark F. Adams 11009371c9d4SSatish Balay static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) { 1101ffc955d6SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1102ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1103ffc955d6SMark F. Adams 1104ffc955d6SMark F. Adams PetscFunctionBegin; 1105cab9ed1eSBarry Smith pc_gamg->use_aggs_in_asm = flg; 1106ffc955d6SMark F. Adams PetscFunctionReturn(0); 1107ffc955d6SMark F. Adams } 1108ffc955d6SMark F. Adams 1109171cca9aSMark Adams /*@ 1110cf8ae1d3SMark Adams PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver 1111171cca9aSMark Adams 1112*f1580f4eSBarry Smith Collective on pc 1113171cca9aSMark Adams 1114171cca9aSMark Adams Input Parameters: 1115171cca9aSMark Adams + pc - the preconditioner context 1116*f1580f4eSBarry Smith - flg - `PETSC_TRUE` to not force coarse grid onto one processor 1117171cca9aSMark Adams 1118171cca9aSMark Adams Options Database Key: 1119147403d9SBarry Smith . -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver 1120171cca9aSMark Adams 1121171cca9aSMark Adams Level: intermediate 1122171cca9aSMark Adams 1123*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()` 1124171cca9aSMark Adams @*/ 11259371c9d4SSatish Balay PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg) { 1126171cca9aSMark Adams PetscFunctionBegin; 1127171cca9aSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1128cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetUseParallelCoarseGridSolve_C", (PC, PetscBool), (pc, flg)); 1129171cca9aSMark Adams PetscFunctionReturn(0); 1130171cca9aSMark Adams } 1131171cca9aSMark Adams 11329371c9d4SSatish Balay static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) { 1133171cca9aSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1134171cca9aSMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1135171cca9aSMark Adams 1136171cca9aSMark Adams PetscFunctionBegin; 1137171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = flg; 1138ffc955d6SMark F. Adams PetscFunctionReturn(0); 1139ffc955d6SMark F. Adams } 1140ffc955d6SMark F. Adams 11414ef23d27SMark F. Adams /*@ 1142*f1580f4eSBarry 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 1143ce7c7f2fSMark Adams 1144*f1580f4eSBarry Smith Collective on pc 1145ce7c7f2fSMark Adams 1146ce7c7f2fSMark Adams Input Parameters: 1147ce7c7f2fSMark Adams + pc - the preconditioner context 1148*f1580f4eSBarry Smith - flg - `PETSC_TRUE` to pin coarse grids to the CPU 1149ce7c7f2fSMark Adams 1150ce7c7f2fSMark Adams Options Database Key: 1151147403d9SBarry Smith . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU 1152ce7c7f2fSMark Adams 1153ce7c7f2fSMark Adams Level: intermediate 1154ce7c7f2fSMark Adams 1155*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetUseParallelCoarseGridSolve()` 1156ce7c7f2fSMark Adams @*/ 11579371c9d4SSatish Balay PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) { 1158ce7c7f2fSMark Adams PetscFunctionBegin; 1159ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1160cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetCpuPinCoarseGrids_C", (PC, PetscBool), (pc, flg)); 1161ce7c7f2fSMark Adams PetscFunctionReturn(0); 1162ce7c7f2fSMark Adams } 1163ce7c7f2fSMark Adams 11649371c9d4SSatish Balay static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) { 1165ce7c7f2fSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1166ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1167ce7c7f2fSMark Adams 1168ce7c7f2fSMark Adams PetscFunctionBegin; 1169ce7c7f2fSMark Adams pc_gamg->cpu_pin_coarse_grids = flg; 1170ce7c7f2fSMark Adams PetscFunctionReturn(0); 1171ce7c7f2fSMark Adams } 1172ce7c7f2fSMark Adams 1173ce7c7f2fSMark Adams /*@ 1174147403d9SBarry Smith PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type) 1175ce7c7f2fSMark Adams 1176*f1580f4eSBarry Smith Collective on pc 1177ce7c7f2fSMark Adams 1178ce7c7f2fSMark Adams Input Parameters: 1179ce7c7f2fSMark Adams + pc - the preconditioner context 1180*f1580f4eSBarry Smith - flg - `PCGAMGLayoutType` type, either `PCGAMG_LAYOUT_COMPACT` or `PCGAMG_LAYOUT_SPREAD` 1181ce7c7f2fSMark Adams 1182ce7c7f2fSMark Adams Options Database Key: 1183147403d9SBarry Smith . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering 1184ce7c7f2fSMark Adams 1185ce7c7f2fSMark Adams Level: intermediate 1186ce7c7f2fSMark Adams 1187*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGLayoutType`, `PCGAMG_LAYOUT_COMPACT`, `PCGAMG_LAYOUT_SPREAD` 1188ce7c7f2fSMark Adams @*/ 11899371c9d4SSatish Balay PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) { 1190ce7c7f2fSMark Adams PetscFunctionBegin; 1191ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1192cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetCoarseGridLayoutType_C", (PC, PCGAMGLayoutType), (pc, flg)); 1193ce7c7f2fSMark Adams PetscFunctionReturn(0); 1194ce7c7f2fSMark Adams } 1195ce7c7f2fSMark Adams 11969371c9d4SSatish Balay static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) { 1197ce7c7f2fSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1198ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1199ce7c7f2fSMark Adams 1200ce7c7f2fSMark Adams PetscFunctionBegin; 1201ce7c7f2fSMark Adams pc_gamg->layout_type = flg; 1202ce7c7f2fSMark Adams PetscFunctionReturn(0); 1203ce7c7f2fSMark Adams } 1204ce7c7f2fSMark Adams 1205ce7c7f2fSMark Adams /*@ 1206*f1580f4eSBarry Smith PCGAMGSetNlevels - Sets the maximum number of levels `PCGAMG` will use 12074ef23d27SMark F. Adams 1208*f1580f4eSBarry Smith Not collective on pc 12094ef23d27SMark F. Adams 12104ef23d27SMark F. Adams Input Parameters: 12111cc46a46SBarry Smith + pc - the preconditioner 12121cc46a46SBarry Smith - n - the maximum number of levels to use 12134ef23d27SMark F. Adams 12144ef23d27SMark F. Adams Options Database Key: 1215147403d9SBarry Smith . -pc_mg_levels <n> - set the maximum number of levels to allow 12164ef23d27SMark F. Adams 12174ef23d27SMark F. Adams Level: intermediate 12184ef23d27SMark F. Adams 1219*f1580f4eSBarry Smith Developer Note: 1220*f1580f4eSBarry Smith Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG` 1221*f1580f4eSBarry Smith 1222*f1580f4eSBarry Smith .seealso: `PCGAMG` 12234ef23d27SMark F. Adams @*/ 12249371c9d4SSatish Balay PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) { 12254ef23d27SMark F. Adams PetscFunctionBegin; 12264ef23d27SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1227cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetNlevels_C", (PC, PetscInt), (pc, n)); 12284ef23d27SMark F. Adams PetscFunctionReturn(0); 12294ef23d27SMark F. Adams } 12304ef23d27SMark F. Adams 12319371c9d4SSatish Balay static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) { 12324ef23d27SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 12334ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 12344ef23d27SMark F. Adams 12354ef23d27SMark F. Adams PetscFunctionBegin; 12369d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 12374ef23d27SMark F. Adams PetscFunctionReturn(0); 12384ef23d27SMark F. Adams } 12394ef23d27SMark F. Adams 12403542efc5SMark F. Adams /*@ 12413542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 12423542efc5SMark F. Adams 1243*f1580f4eSBarry Smith Not collective on pc 12443542efc5SMark F. Adams 12453542efc5SMark F. Adams Input Parameters: 12461cc46a46SBarry Smith + pc - the preconditioner context 1247c9567895SMark . 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 1248055c8bd0SJed Brown - n - number of threshold values provided in array 12493542efc5SMark F. Adams 12503542efc5SMark F. Adams Options Database Key: 1251147403d9SBarry Smith . -pc_gamg_threshold <threshold> - the threshold to drop edges 12523542efc5SMark F. Adams 125395452b02SPatrick Sanan Notes: 1254af3c827dSMark 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. 1255*f1580f4eSBarry 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. 1256cab9ed1eSBarry Smith 1257*f1580f4eSBarry Smith If n is less than the total number of coarsenings (see `PCGAMGSetNlevels()`), then threshold scaling (see `PCGAMGSetThresholdScale()`) is used for each successive coarsening. 1258*f1580f4eSBarry Smith In this case, `PCGAMGSetThresholdScale()` must be called before `PCGAMGSetThreshold()`. 1259055c8bd0SJed Brown If n is greater than the total number of levels, the excess entries in threshold will not be used. 1260055c8bd0SJed Brown 12613542efc5SMark F. Adams Level: intermediate 12623542efc5SMark F. Adams 1263*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGFilterGraph()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGSetThresholdScale()` 12643542efc5SMark F. Adams @*/ 12659371c9d4SSatish Balay PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) { 12663542efc5SMark F. Adams PetscFunctionBegin; 12673542efc5SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1268055c8bd0SJed Brown if (n) PetscValidRealPointer(v, 2); 1269cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetThreshold_C", (PC, PetscReal[], PetscInt), (pc, v, n)); 12703542efc5SMark F. Adams PetscFunctionReturn(0); 12713542efc5SMark F. Adams } 12723542efc5SMark F. Adams 12739371c9d4SSatish Balay static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) { 1274c20e4228SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1275c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1276c1eae691SMark Adams PetscInt i; 1277c1eae691SMark Adams PetscFunctionBegin; 1278055c8bd0SJed Brown for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i]; 1279055c8bd0SJed Brown for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale; 1280c1eae691SMark Adams PetscFunctionReturn(0); 1281c1eae691SMark Adams } 1282c1eae691SMark Adams 1283c1eae691SMark Adams /*@ 1284*f1580f4eSBarry Smith PCGAMGSetRankReductionFactors - Set a manual schedule for MPI rank reduction on coarse grids 1285c9567895SMark 1286*f1580f4eSBarry Smith Collective on pc 1287c9567895SMark 1288c9567895SMark Input Parameters: 1289c9567895SMark + pc - the preconditioner context 1290*f1580f4eSBarry Smith . v - array of reduction factors. 0 for first value forces a reduction to one process/device on first level in CUDA 1291c9567895SMark - n - number of values provided in array 1292c9567895SMark 1293c9567895SMark Options Database Key: 1294147403d9SBarry Smith . -pc_gamg_rank_reduction_factors <factors> - provide the schedule 1295c9567895SMark 1296c9567895SMark Level: intermediate 1297c9567895SMark 1298*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()` 1299c9567895SMark @*/ 13009371c9d4SSatish Balay PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n) { 1301c9567895SMark PetscFunctionBegin; 1302c9567895SMark PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1303c9567895SMark if (n) PetscValidIntPointer(v, 2); 1304cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetRankReductionFactors_C", (PC, PetscInt[], PetscInt), (pc, v, n)); 1305c9567895SMark PetscFunctionReturn(0); 1306c9567895SMark } 1307c9567895SMark 13089371c9d4SSatish Balay static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n) { 1309c9567895SMark PC_MG *mg = (PC_MG *)pc->data; 1310c9567895SMark PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1311c9567895SMark PetscInt i; 1312c9567895SMark PetscFunctionBegin; 1313c9567895SMark for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i]; 1314c9567895SMark for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */ 1315c9567895SMark PetscFunctionReturn(0); 1316c9567895SMark } 1317c9567895SMark 1318c9567895SMark /*@ 1319c1eae691SMark Adams PCGAMGSetThresholdScale - Relative threshold reduction at each level 1320c1eae691SMark Adams 1321*f1580f4eSBarry Smith Not collective on pc 1322c1eae691SMark Adams 1323c1eae691SMark Adams Input Parameters: 1324c1eae691SMark Adams + pc - the preconditioner context 1325147403d9SBarry Smith - scale - the threshold value reduction, usually < 1.0 1326c1eae691SMark Adams 1327c1eae691SMark Adams Options Database Key: 1328147403d9SBarry Smith . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level 1329c1eae691SMark Adams 1330*f1580f4eSBarry Smith Note: 1331*f1580f4eSBarry Smith The initial threshold (for an arbitrary number of levels starting from the finest) can be set with `PCGAMGSetThreshold()`. 1332*f1580f4eSBarry Smith This scaling is used for each subsequent coarsening, but must be called before `PCGAMGSetThreshold()`. 1333055c8bd0SJed Brown 1334c1eae691SMark Adams Level: advanced 1335c1eae691SMark Adams 1336*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetThreshold()` 1337c1eae691SMark Adams @*/ 13389371c9d4SSatish Balay PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) { 13393542efc5SMark F. Adams PetscFunctionBegin; 1340c1eae691SMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1341cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetThresholdScale_C", (PC, PetscReal), (pc, v)); 1342c1eae691SMark Adams PetscFunctionReturn(0); 1343c1eae691SMark Adams } 1344c1eae691SMark Adams 13459371c9d4SSatish Balay static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) { 1346c1eae691SMark Adams PC_MG *mg = (PC_MG *)pc->data; 1347c1eae691SMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1348c1eae691SMark Adams PetscFunctionBegin; 1349c1eae691SMark Adams pc_gamg->threshold_scale = v; 13503542efc5SMark F. Adams PetscFunctionReturn(0); 13513542efc5SMark F. Adams } 13523542efc5SMark F. Adams 1353e20c40e8SBarry Smith /*@C 1354*f1580f4eSBarry Smith PCGAMGSetType - Set the type of algorithm `PCGAMG` should use 1355676e1743SMark F. Adams 1356*f1580f4eSBarry Smith Collective on pc 1357676e1743SMark F. Adams 1358676e1743SMark F. Adams Input Parameters: 1359c60c7ad4SBarry Smith + pc - the preconditioner context 1360*f1580f4eSBarry Smith - type - `PCGAMGAGG`, `PCGAMGGEO`, or `PCGAMGCLASSICAL` 1361676e1743SMark F. Adams 1362676e1743SMark F. Adams Options Database Key: 1363cab9ed1eSBarry Smith . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply 1364676e1743SMark F. Adams 1365676e1743SMark F. Adams Level: intermediate 1366676e1743SMark F. Adams 1367db781477SPatrick Sanan .seealso: `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType` 1368676e1743SMark F. Adams @*/ 13699371c9d4SSatish Balay PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) { 1370676e1743SMark F. Adams PetscFunctionBegin; 1371676e1743SMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1372cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetType_C", (PC, PCGAMGType), (pc, type)); 1373676e1743SMark F. Adams PetscFunctionReturn(0); 1374676e1743SMark F. Adams } 1375676e1743SMark F. Adams 1376e20c40e8SBarry Smith /*@C 1377*f1580f4eSBarry Smith PCGAMGGetType - Get the type of algorithm `PCGAMG` will use 1378c60c7ad4SBarry Smith 1379*f1580f4eSBarry Smith Collective on pc 1380c60c7ad4SBarry Smith 1381c60c7ad4SBarry Smith Input Parameter: 1382c60c7ad4SBarry Smith . pc - the preconditioner context 1383c60c7ad4SBarry Smith 1384c60c7ad4SBarry Smith Output Parameter: 1385c60c7ad4SBarry Smith . type - the type of algorithm used 1386c60c7ad4SBarry Smith 1387c60c7ad4SBarry Smith Level: intermediate 1388c60c7ad4SBarry Smith 1389*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetType()`, `PCGAMGType` 1390c60c7ad4SBarry Smith @*/ 13919371c9d4SSatish Balay PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) { 1392c60c7ad4SBarry Smith PetscFunctionBegin; 1393c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1394cac4c232SBarry Smith PetscUseMethod(pc, "PCGAMGGetType_C", (PC, PCGAMGType *), (pc, type)); 1395c60c7ad4SBarry Smith PetscFunctionReturn(0); 1396c60c7ad4SBarry Smith } 1397c60c7ad4SBarry Smith 13989371c9d4SSatish Balay static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) { 1399c60c7ad4SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1400c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1401c60c7ad4SBarry Smith 1402c60c7ad4SBarry Smith PetscFunctionBegin; 1403c60c7ad4SBarry Smith *type = pc_gamg->type; 1404c60c7ad4SBarry Smith PetscFunctionReturn(0); 1405c60c7ad4SBarry Smith } 1406c60c7ad4SBarry Smith 14079371c9d4SSatish Balay static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) { 14081ab5ffc9SJed Brown PC_MG *mg = (PC_MG *)pc->data; 14091ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 14105f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(PC); 1411676e1743SMark F. Adams 1412676e1743SMark F. Adams PetscFunctionBegin; 1413c60c7ad4SBarry Smith pc_gamg->type = type; 14149566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(GAMGList, type, &r)); 14155f80ce2aSJacob Faibussowitsch PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown GAMG type %s given", type); 14161ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 14179566063dSJacob Faibussowitsch PetscCall((*pc_gamg->ops->destroy)(pc)); 14189566063dSJacob Faibussowitsch PetscCall(PetscMemzero(pc_gamg->ops, sizeof(struct _PCGAMGOps))); 1419e616c208SToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 14203ae0bb68SMark Adams /* cleaning up common data in pc_gamg - this should disapear someday */ 14213ae0bb68SMark Adams pc_gamg->data_cell_cols = 0; 14223ae0bb68SMark Adams pc_gamg->data_cell_rows = 0; 14233ae0bb68SMark Adams pc_gamg->orig_data_cell_cols = 0; 14243ae0bb68SMark Adams pc_gamg->orig_data_cell_rows = 0; 14259566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 14263ae0bb68SMark Adams pc_gamg->data_sz = 0; 14271ab5ffc9SJed Brown } 14289566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->gamg_type_name)); 14299566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(type, &pc_gamg->gamg_type_name)); 14309566063dSJacob Faibussowitsch PetscCall((*r)(pc)); 1431676e1743SMark F. Adams PetscFunctionReturn(0); 1432676e1743SMark F. Adams } 1433676e1743SMark F. Adams 14349371c9d4SSatish Balay static PetscErrorCode PCView_GAMG(PC pc, PetscViewer viewer) { 14355adeb434SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 14365adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1437e7d4b4cbSMark Adams PetscReal gc = 0, oc = 0; 143890db8557SMark Adams 14395adeb434SBarry Smith PetscFunctionBegin; 14409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GAMG specific options\n")); 14419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Threshold for dropping small values in graph on each level =")); 14429566063dSJacob Faibussowitsch for (PetscInt i = 0; i < mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)pc_gamg->threshold[i])); 14439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 14449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Threshold scaling factor for each level not specified = %g\n", (double)pc_gamg->threshold_scale)); 144548a46eb9SPierre Jolivet if (pc_gamg->use_aggs_in_asm) PetscCall(PetscViewerASCIIPrintf(viewer, " Using aggregates from coarsening process to define subdomains for PCASM\n")); 144648a46eb9SPierre 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")); 14471baa6e33SBarry Smith if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc, viewer)); 14489566063dSJacob Faibussowitsch PetscCall(PCMGGetGridComplexity(pc, &gc, &oc)); 144963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Complexity: grid = %g operator = %g\n", (double)gc, (double)oc)); 14505adeb434SBarry Smith PetscFunctionReturn(0); 14515adeb434SBarry Smith } 14525adeb434SBarry Smith 14539371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_GAMG(PC pc, PetscOptionItems *PetscOptionsObject) { 1454676e1743SMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1455676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 14567e6512fdSJed Brown PetscBool flag; 14573b4367a7SBarry Smith MPI_Comm comm; 145818c3aa7eSMark char prefix[256], tname[32]; 1459c1eae691SMark Adams PetscInt i, n; 146014a9496bSBarry Smith const char *pcpre; 14610a545947SLisandro Dalcin static const char *LayoutTypes[] = {"compact", "spread", "PCGAMGLayoutType", "PC_GAMG_LAYOUT", NULL}; 14625b89ad90SMark F. Adams PetscFunctionBegin; 14639566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 1464d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG options"); 14659566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_gamg_type", "Type of AMG method", "PCGAMGSetType", GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag)); 14661baa6e33SBarry Smith if (flag) PetscCall(PCGAMGSetType(pc, tname)); 14679566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_repartition", "Repartion coarse grids", "PCGAMGSetRepartition", pc_gamg->repart, &pc_gamg->repart, NULL)); 1468*f1580f4eSBarry 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)); 14699566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation", "Reuse prolongation operator", "PCGAMGReuseInterpolation", pc_gamg->reuse_prol, &pc_gamg->reuse_prol, NULL)); 14709566063dSJacob 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)); 14719566063dSJacob 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)); 14729566063dSJacob 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)); 14739371c9d4SSatish 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, 14749371c9d4SSatish Balay (PetscEnum)pc_gamg->layout_type, (PetscEnum *)&pc_gamg->layout_type, NULL)); 14759566063dSJacob 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)); 14769566063dSJacob 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)); 14779566063dSJacob 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)); 147818c3aa7eSMark n = PETSC_MG_MAXLEVELS; 14799566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", pc_gamg->threshold, &n, &flag)); 148018c3aa7eSMark if (!flag || n < PETSC_MG_MAXLEVELS) { 1481efd3c5ceSMark Adams if (!flag) n = 1; 1482c1eae691SMark Adams i = n; 148318c3aa7eSMark do { pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale; } while (++i < PETSC_MG_MAXLEVELS); 1484c1eae691SMark Adams } 1485c9567895SMark n = PETSC_MG_MAXLEVELS; 14869566063dSJacob 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)); 1487c9567895SMark if (!flag) i = 0; 1488c9567895SMark else i = n; 1489c9567895SMark do { pc_gamg->level_reduction_factors[i] = -1; } while (++i < PETSC_MG_MAXLEVELS); 14909566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_levels", "Set number of MG levels", "PCGAMGSetNlevels", pc_gamg->Nlevels, &pc_gamg->Nlevels, NULL)); 149118c3aa7eSMark { 149218c3aa7eSMark PetscReal eminmax[2] = {0., 0.}; 149318c3aa7eSMark n = 2; 14949566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues", "extreme eigenvalues for smoothed aggregation", "PCGAMGSetEigenvalues", eminmax, &n, &flag)); 149518c3aa7eSMark if (flag) { 149608401ef6SPierre Jolivet PetscCheck(n == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 14979566063dSJacob Faibussowitsch PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0])); 149818c3aa7eSMark } 149918c3aa7eSMark } 1500b7cbab4eSMark Adams /* set options for subtype */ 1501dbbe0bcdSBarry Smith PetscCall((*pc_gamg->ops->setfromoptions)(pc, PetscOptionsObject)); 150218c3aa7eSMark 15039566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 15049566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : "")); 1505d0609cedSBarry Smith PetscOptionsHeadEnd(); 15065b89ad90SMark F. Adams PetscFunctionReturn(0); 15075b89ad90SMark F. Adams } 15085b89ad90SMark F. Adams 15095b89ad90SMark F. Adams /*MC 15101cc46a46SBarry Smith PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 15115b89ad90SMark F. Adams 1512280d9858SJed Brown Options Database Keys: 1513cab9ed1eSBarry Smith + -pc_gamg_type <type> - one of agg, geo, or classical 1514cab9ed1eSBarry Smith . -pc_gamg_repartition <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined 1515bae903cbSmarkadams4 - -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations 1516*f1580f4eSBarry Smith + -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the `PCASM` smoother 1517*f1580f4eSBarry Smith . -pc_gamg_process_eq_limit <limit, default=50> - `PCGAMG` will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> 1518cab9ed1eSBarry Smith equations on each process that has degrees of freedom 1519bae903cbSmarkadams4 - -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 1520*f1580f4eSBarry Smith + -pc_gamg_threshold[] <thresh,default=-1> - Before aggregating the graph `PCGAMG` will remove small values from the graph on each level (< 0 is no filtering) 1521bae903cbSmarkadams4 . -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1522cab9ed1eSBarry Smith 1523*f1580f4eSBarry Smith Options Database Keys for Aggregation: 1524cab9ed1eSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 1525bae903cbSmarkadams4 . -pc_gamg_symmetrize_graph <true,default=false> - symmetrize the graph before computing the aggregation 1526bae903cbSmarkadams4 - -pc_gamg_square_graph <n,default=1> - alias for pc_gamg_aggressive_coarsening (deprecated) 1527bae903cbSmarkadams4 - -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest. 1528cab9ed1eSBarry Smith 1529*f1580f4eSBarry Smith Options Database Keys for Multigrid: 1530*f1580f4eSBarry Smith + -pc_mg_cycles <v> - v or w, see `PCMGSetCycleType()` 1531db9745e2SBarry Smith . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1532db9745e2SBarry Smith . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1533cab9ed1eSBarry Smith - -pc_mg_levels <levels> - Number of levels of multigrid to use. 15345b89ad90SMark F. Adams 153595452b02SPatrick Sanan Notes: 1536*f1580f4eSBarry Smith To obtain good performance for `PCGAMG` for vector valued problems you must 1537*f1580f4eSBarry Smith call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point 1538*f1580f4eSBarry Smith call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator 1539*f1580f4eSBarry Smith 1540*f1580f4eSBarry Smith See [the Users Manual section on PCGAMG](sec_amg) for more details. 15411cc46a46SBarry Smith 15425b89ad90SMark F. Adams Level: intermediate 1543280d9858SJed Brown 1544db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, 1545db781477SPatrick Sanan `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGSetUseSAEstEig()` 15465b89ad90SMark F. Adams M*/ 1547b2573a8aSBarry Smith 15489371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) { 15495b89ad90SMark F. Adams PC_GAMG *pc_gamg; 15505b89ad90SMark F. Adams PC_MG *mg; 15515b89ad90SMark F. Adams 15525b89ad90SMark F. Adams PetscFunctionBegin; 15531c1aac46SBarry Smith /* register AMG type */ 15549566063dSJacob Faibussowitsch PetscCall(PCGAMGInitializePackage()); 15551c1aac46SBarry Smith 15565b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 15579566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCMG)); 15589566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG)); 15595b89ad90SMark F. Adams 15605b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 15619566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &pc_gamg)); 15629566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 15635b89ad90SMark F. Adams mg = (PC_MG *)pc->data; 15645b89ad90SMark F. Adams mg->innerctx = pc_gamg; 15655b89ad90SMark F. Adams 15669566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &pc_gamg->ops)); 15671ab5ffc9SJed Brown 15689d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 15699d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 15700a545947SLisandro Dalcin pc_gamg->data = NULL; 15715b89ad90SMark F. Adams 15729d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 15735b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 15745b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 15755b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 15765b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 15775adeb434SBarry Smith mg->view = PCView_GAMG; 15785b89ad90SMark F. Adams 15799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", PCGAMGSetProcEqLim_GAMG)); 15809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", PCGAMGSetCoarseEqLim_GAMG)); 15819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", PCGAMGSetRepartition_GAMG)); 15829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", PCGAMGSetEigenvalues_GAMG)); 15839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", PCGAMGSetUseSAEstEig_GAMG)); 15849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", PCGAMGSetReuseInterpolation_GAMG)); 15859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", PCGAMGASMSetUseAggs_GAMG)); 15869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", PCGAMGSetUseParallelCoarseGridSolve_GAMG)); 15879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", PCGAMGSetCpuPinCoarseGrids_GAMG)); 15889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", PCGAMGSetCoarseGridLayoutType_GAMG)); 15899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", PCGAMGSetThreshold_GAMG)); 15909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", PCGAMGSetRankReductionFactors_GAMG)); 15919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", PCGAMGSetThresholdScale_GAMG)); 15929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", PCGAMGSetType_GAMG)); 15939566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", PCGAMGGetType_GAMG)); 15949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", PCGAMGSetNlevels_GAMG)); 15959d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1596d3042614SMark Adams pc_gamg->reuse_prol = PETSC_FALSE; 15970c3bc534SBarry Smith pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1598171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1599a0095786SMark pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1600a0095786SMark pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1601038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 160225a145a7SMark Adams pc_gamg->coarse_eq_limit = 50; 160353134ebeSMark Adams for (int i = 0; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = -1; 1604c1eae691SMark Adams pc_gamg->threshold_scale = 1.; 160518c3aa7eSMark pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 16069ab59c8bSMark Adams pc_gamg->current_level = 0; /* don't need to init really */ 16077e6512fdSJed Brown pc_gamg->use_sa_esteig = PETSC_TRUE; 160818c3aa7eSMark pc_gamg->emin = 0; 160918c3aa7eSMark pc_gamg->emax = 0; 161018c3aa7eSMark 1611c238b0ebSToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 16129d5b6da9SMark F. Adams 1613bd94a7aaSJed Brown /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 16149566063dSJacob Faibussowitsch PetscCall(PCGAMGSetType(pc, PCGAMGAGG)); 16155b89ad90SMark F. Adams PetscFunctionReturn(0); 16165b89ad90SMark F. Adams } 16173e3471ccSMark Adams 16183e3471ccSMark Adams /*@C 1619*f1580f4eSBarry Smith PCGAMGInitializePackage - This function initializes everything in the `PCGAMG` package. It is called 1620*f1580f4eSBarry Smith from `PCInitializePackage()`. 16213e3471ccSMark Adams 16223e3471ccSMark Adams Level: developer 16233e3471ccSMark Adams 1624db781477SPatrick Sanan .seealso: `PetscInitialize()` 16253e3471ccSMark Adams @*/ 16269371c9d4SSatish Balay PetscErrorCode PCGAMGInitializePackage(void) { 16274555aa8cSStefano Zampini PetscInt l; 16283e3471ccSMark Adams 16293e3471ccSMark Adams PetscFunctionBegin; 16303e3471ccSMark Adams if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 16313e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 16329566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGGEO, PCCreateGAMG_GEO)); 16339566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGAGG, PCCreateGAMG_AGG)); 16349566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGCLASSICAL, PCCreateGAMG_Classical)); 16359566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage)); 1636c1c463dbSMark Adams 1637c1c463dbSMark Adams /* general events */ 1638849bee69SMark Adams PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP])); 1639849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH])); 1640849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGFilter", PC_CLASSID, &petsc_gamg_setup_events[GAMG_FILTER])); 1641849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN])); 1642849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS])); 1643849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL])); 1644849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA])); 1645849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB])); 1646849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT])); 1647849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM])); 1648849bee69SMark Adams PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL])); 1649849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG PtAP", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP])); 1650849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Reduce", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE])); 1651849bee69SMark Adams PetscCall(PetscLogEventRegister(" GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART])); 1652849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */ 1653849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */ 1654849bee69SMark Adams /* PetscCall(PetscLogEventRegister(" GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */ 16554555aa8cSStefano Zampini for (l = 0; l < PETSC_MG_MAXLEVELS; l++) { 16564555aa8cSStefano Zampini char ename[32]; 16575b89ad90SMark F. Adams 165863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Squ l%02" PetscInt_FMT, l)); 16599566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0])); 166063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Gal l%02" PetscInt_FMT, l)); 16619566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1])); 166263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Opt l%02" PetscInt_FMT, l)); 16639566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2])); 16644555aa8cSStefano Zampini } 16654555aa8cSStefano Zampini #if defined(GAMG_STAGES) 1666849bee69SMark Adams { /* create timer stages */ 16675b89ad90SMark F. Adams char str[32]; 1668849bee69SMark Adams sprintf(str, "GAMG Level %d", 0); 16699566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(str, &gamg_stages[0])); 16705b89ad90SMark F. Adams } 16715b89ad90SMark F. Adams #endif 16723e3471ccSMark Adams PetscFunctionReturn(0); 16733e3471ccSMark Adams } 16743e3471ccSMark Adams 16753e3471ccSMark Adams /*@C 1676*f1580f4eSBarry Smith PCGAMGFinalizePackage - This function frees everything from the `PCGAMG` package. It is 1677*f1580f4eSBarry Smith called from `PetscFinalize()` automatically. 16783e3471ccSMark Adams 16793e3471ccSMark Adams Level: developer 16803e3471ccSMark Adams 1681db781477SPatrick Sanan .seealso: `PetscFinalize()` 16823e3471ccSMark Adams @*/ 16839371c9d4SSatish Balay PetscErrorCode PCGAMGFinalizePackage(void) { 16843e3471ccSMark Adams PetscFunctionBegin; 16853e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 16869566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&GAMGList)); 16873e3471ccSMark Adams PetscFunctionReturn(0); 16883e3471ccSMark Adams } 1689a36cf38bSToby Isaac 1690a36cf38bSToby Isaac /*@C 1691*f1580f4eSBarry Smith PCGAMGRegister - Register a `PCGAMG` implementation. 1692a36cf38bSToby Isaac 1693a36cf38bSToby Isaac Input Parameters: 1694*f1580f4eSBarry Smith + type - string that will be used as the name of the `PCGAMG` type. 1695a36cf38bSToby Isaac - create - function for creating the gamg context. 1696a36cf38bSToby Isaac 1697*f1580f4eSBarry Smith Level: developer 1698a36cf38bSToby Isaac 1699db781477SPatrick Sanan .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()` 1700a36cf38bSToby Isaac @*/ 17019371c9d4SSatish Balay PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) { 1702a36cf38bSToby Isaac PetscFunctionBegin; 17039566063dSJacob Faibussowitsch PetscCall(PCGAMGInitializePackage()); 17049566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&GAMGList, type, create)); 1705a36cf38bSToby Isaac PetscFunctionReturn(0); 1706a36cf38bSToby Isaac } 1707