1 /* 2 GAMG geometric-algebric multigrid PC - Mark Adams 2011 3 */ 4 #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 5 #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h> /*I "petscksp.h" I*/ 6 7 #if defined(PETSC_HAVE_CUDA) 8 #include <cuda_runtime.h> 9 #endif 10 11 #if defined(PETSC_HAVE_HIP) 12 #include <hip/hip_runtime.h> 13 #endif 14 15 PetscLogEvent petsc_gamg_setup_events[GAMG_NUM_SET]; 16 PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3]; 17 18 // #define GAMG_STAGES 19 #if defined(GAMG_STAGES) 20 static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS]; 21 #endif 22 23 static PetscFunctionList GAMGList = NULL; 24 static PetscBool PCGAMGPackageInitialized; 25 26 static PetscErrorCode PCReset_GAMG(PC pc) 27 { 28 PC_MG *mg = (PC_MG *)pc->data; 29 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 30 31 PetscFunctionBegin; 32 PetscCall(PetscFree(pc_gamg->data)); 33 pc_gamg->data_sz = 0; 34 PetscCall(PetscFree(pc_gamg->orig_data)); 35 for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS; level++) { 36 mg->min_eigen_DinvA[level] = 0; 37 mg->max_eigen_DinvA[level] = 0; 38 } 39 pc_gamg->emin = 0; 40 pc_gamg->emax = 0; 41 PetscCall(PCReset_MG(pc)); 42 PetscCall(MatCoarsenDestroy(&pc_gamg->asm_crs)); 43 PetscFunctionReturn(PETSC_SUCCESS); 44 } 45 46 /* 47 PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 48 of active processors. 49 50 Input Parameter: 51 . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 52 'pc_gamg->data_sz' are changed via repartitioning/reduction. 53 . Amat_fine - matrix on this fine (k) level 54 . cr_bs - coarse block size 55 In/Output Parameter: 56 . a_P_inout - prolongation operator to the next level (k-->k-1) 57 . a_nactive_proc - number of active procs 58 Output Parameter: 59 . a_Amat_crs - coarse matrix that is created (k-1) 60 */ 61 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) 62 { 63 PC_MG *mg = (PC_MG *)pc->data; 64 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 65 Mat Cmat = NULL, Pold = *a_P_inout; 66 MPI_Comm comm; 67 PetscMPIInt rank, size, new_size, nactive = *a_nactive_proc; 68 PetscInt ncrs_eq, ncrs, f_bs; 69 70 PetscFunctionBegin; 71 PetscCall(PetscObjectGetComm((PetscObject)Amat_fine, &comm)); 72 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 73 PetscCallMPI(MPI_Comm_size(comm, &size)); 74 PetscCall(MatGetBlockSize(Amat_fine, &f_bs)); 75 76 if (Pcolumnperm) *Pcolumnperm = NULL; 77 78 /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 79 PetscCall(MatGetLocalSize(Pold, NULL, &ncrs_eq)); 80 if (pc_gamg->data_cell_rows > 0) { 81 ncrs = pc_gamg->data_sz / pc_gamg->data_cell_cols / pc_gamg->data_cell_rows; 82 } else { 83 PetscInt bs; 84 PetscCall(MatGetBlockSizes(Pold, NULL, &bs)); 85 ncrs = ncrs_eq / bs; 86 } 87 /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 88 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. */ 89 #if defined(PETSC_HAVE_CUDA) 90 PetscShmComm pshmcomm; 91 PetscMPIInt locrank; 92 MPI_Comm loccomm; 93 PetscInt s_nnodes, r_nnodes, new_new_size; 94 cudaError_t cerr; 95 int devCount; 96 PetscCall(PetscShmCommGet(comm, &pshmcomm)); 97 PetscCall(PetscShmCommGetMpiShmComm(pshmcomm, &loccomm)); 98 PetscCallMPI(MPI_Comm_rank(loccomm, &locrank)); 99 s_nnodes = !locrank; 100 PetscCall(MPIU_Allreduce(&s_nnodes, &r_nnodes, 1, MPIU_INT, MPI_SUM, comm)); 101 PetscCheck((size % r_nnodes) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "odd number of nodes np=%d nnodes%" PetscInt_FMT, size, r_nnodes); 102 devCount = 0; 103 cerr = cudaGetDeviceCount(&devCount); 104 cudaGetLastError(); /* Reset the last error */ 105 if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */ 106 new_new_size = r_nnodes * devCount; 107 new_size = new_new_size; 108 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)); 109 } else { 110 PetscCall(PetscInfo(pc, "%s: With Cuda but no device. Use heuristics.\n", ((PetscObject)pc)->prefix)); 111 goto HEURISTIC; 112 } 113 #else 114 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "should not be here"); 115 #endif 116 } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) { 117 if (nactive < pc_gamg->level_reduction_factors[pc_gamg->current_level]) { 118 new_size = 1; 119 PetscCall(PetscInfo(pc, "%s: reduction factor too small for %d active processes: reduce to one process\n", ((PetscObject)pc)->prefix, new_size)); 120 } else { 121 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]); 122 new_size = nactive / pc_gamg->level_reduction_factors[pc_gamg->current_level]; 123 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])); 124 } 125 } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) { 126 new_size = 1; 127 PetscCall(PetscInfo(pc, "%s: Force coarsest grid reduction to %d active processes\n", ((PetscObject)pc)->prefix, new_size)); 128 } else { 129 PetscInt ncrs_eq_glob; 130 #if defined(PETSC_HAVE_CUDA) 131 HEURISTIC: 132 #endif 133 PetscCall(MatGetSize(Pold, NULL, &ncrs_eq_glob)); 134 new_size = (PetscMPIInt)((float)ncrs_eq_glob / (float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 135 if (!new_size) new_size = 1; /* not likely, possible? */ 136 else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 137 PetscCall(PetscInfo(pc, "%s: Coarse grid reduction from %d to %d active processes\n", ((PetscObject)pc)->prefix, nactive, new_size)); 138 } 139 if (new_size == nactive) { 140 /* output - no repartitioning or reduction - could bail here 141 we know that the grid structure can be reused in MatPtAP */ 142 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 143 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 144 PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, a_Amat_crs)); 145 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 146 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 147 if (new_size < size) { 148 /* odd case where multiple coarse grids are on one processor or no coarsening ... */ 149 PetscCall(PetscInfo(pc, "%s: reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n", ((PetscObject)pc)->prefix, nactive)); 150 if (pc_gamg->cpu_pin_coarse_grids) { 151 PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE)); 152 PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE)); 153 } 154 } 155 } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */ 156 PetscInt *counts, *newproc_idx, ii, jj, kk, strideNew, *tidx, ncrs_new, ncrs_eq_new, nloc_old, expand_factor = 1, rfactor = 1; 157 IS is_eq_newproc, is_eq_num, new_eq_indices; 158 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 159 nloc_old = ncrs_eq / cr_bs; 160 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); 161 /* get new_size and rfactor */ 162 if (pc_gamg->layout_type == PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) { 163 /* find factor */ 164 if (new_size == 1) rfactor = size; /* don't modify */ 165 else { 166 PetscReal best_fact = 0.; 167 jj = -1; 168 for (kk = 1; kk <= size; kk++) { 169 if (!(size % kk)) { /* a candidate */ 170 PetscReal nactpe = (PetscReal)size / (PetscReal)kk, fact = nactpe / (PetscReal)new_size; 171 if (fact > 1.0) fact = 1. / fact; /* keep fact < 1 */ 172 if (fact > best_fact) { 173 best_fact = fact; 174 jj = kk; 175 } 176 } 177 } 178 if (jj != -1) rfactor = jj; 179 else rfactor = 1; /* a prime */ 180 if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1; 181 else expand_factor = rfactor; 182 } 183 new_size = size / rfactor; /* make new size one that is factor */ 184 if (new_size == nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */ 185 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)); 186 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 187 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 188 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 189 PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, a_Amat_crs)); 190 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 191 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 192 PetscFunctionReturn(PETSC_SUCCESS); 193 } 194 } 195 /* make 'is_eq_newproc' */ 196 if (pc_gamg->repart) { /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */ 197 Mat adj; 198 199 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 200 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 201 PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat)); 202 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 203 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 204 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0)); 205 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")); 206 /* get 'adj' */ 207 if (cr_bs == 1) { 208 PetscCall(MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj)); 209 } else { 210 /* make a scalar matrix to partition (no Stokes here) */ 211 Mat tMat; 212 PetscInt Istart_crs, Iend_crs, ncols, jj, Ii; 213 const PetscScalar *vals; 214 const PetscInt *idx; 215 PetscInt *d_nnz, *o_nnz, M, N; 216 static PetscInt llev = 0; /* ugly but just used for debugging */ 217 MatType mtype; 218 219 PetscCall(PetscMalloc2(ncrs, &d_nnz, ncrs, &o_nnz)); 220 PetscCall(MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs)); 221 PetscCall(MatGetSize(Cmat, &M, &N)); 222 for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 223 PetscCall(MatGetRow(Cmat, Ii, &ncols, NULL, NULL)); 224 d_nnz[jj] = ncols / cr_bs; 225 o_nnz[jj] = ncols / cr_bs; 226 PetscCall(MatRestoreRow(Cmat, Ii, &ncols, NULL, NULL)); 227 if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 228 if (o_nnz[jj] > (M / cr_bs - ncrs)) o_nnz[jj] = M / cr_bs - ncrs; 229 } 230 231 PetscCall(MatGetType(Amat_fine, &mtype)); 232 PetscCall(MatCreate(comm, &tMat)); 233 PetscCall(MatSetSizes(tMat, ncrs, ncrs, PETSC_DETERMINE, PETSC_DETERMINE)); 234 PetscCall(MatSetType(tMat, mtype)); 235 PetscCall(MatSeqAIJSetPreallocation(tMat, 0, d_nnz)); 236 PetscCall(MatMPIAIJSetPreallocation(tMat, 0, d_nnz, 0, o_nnz)); 237 PetscCall(PetscFree2(d_nnz, o_nnz)); 238 239 for (ii = Istart_crs; ii < Iend_crs; ii++) { 240 PetscInt dest_row = ii / cr_bs; 241 PetscCall(MatGetRow(Cmat, ii, &ncols, &idx, &vals)); 242 for (jj = 0; jj < ncols; jj++) { 243 PetscInt dest_col = idx[jj] / cr_bs; 244 PetscScalar v = 1.0; 245 PetscCall(MatSetValues(tMat, 1, &dest_row, 1, &dest_col, &v, ADD_VALUES)); 246 } 247 PetscCall(MatRestoreRow(Cmat, ii, &ncols, &idx, &vals)); 248 } 249 PetscCall(MatAssemblyBegin(tMat, MAT_FINAL_ASSEMBLY)); 250 PetscCall(MatAssemblyEnd(tMat, MAT_FINAL_ASSEMBLY)); 251 252 if (llev++ == -1) { 253 PetscViewer viewer; 254 char fname[32]; 255 PetscCall(PetscSNPrintf(fname, sizeof(fname), "part_mat_%" PetscInt_FMT ".mat", llev)); 256 PetscCall(PetscViewerBinaryOpen(comm, fname, FILE_MODE_WRITE, &viewer)); 257 PetscCall(MatView(tMat, viewer)); 258 PetscCall(PetscViewerDestroy(&viewer)); 259 } 260 PetscCall(MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj)); 261 PetscCall(MatDestroy(&tMat)); 262 } /* create 'adj' */ 263 264 { /* partition: get newproc_idx */ 265 char prefix[256]; 266 const char *pcpre; 267 const PetscInt *is_idx; 268 MatPartitioning mpart; 269 IS proc_is; 270 271 PetscCall(MatPartitioningCreate(comm, &mpart)); 272 PetscCall(MatPartitioningSetAdjacency(mpart, adj)); 273 PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 274 PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : "")); 275 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix)); 276 PetscCall(MatPartitioningSetFromOptions(mpart)); 277 PetscCall(MatPartitioningSetNParts(mpart, new_size)); 278 PetscCall(MatPartitioningApply(mpart, &proc_is)); 279 PetscCall(MatPartitioningDestroy(&mpart)); 280 281 /* collect IS info */ 282 PetscCall(PetscMalloc1(ncrs_eq, &newproc_idx)); 283 PetscCall(ISGetIndices(proc_is, &is_idx)); 284 for (kk = jj = 0; kk < nloc_old; kk++) { 285 for (ii = 0; ii < cr_bs; ii++, jj++) { newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ } 286 } 287 PetscCall(ISRestoreIndices(proc_is, &is_idx)); 288 PetscCall(ISDestroy(&proc_is)); 289 } 290 PetscCall(MatDestroy(&adj)); 291 292 PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_OWN_POINTER, &is_eq_newproc)); 293 /* 294 Create an index set from the is_eq_newproc index set to indicate the mapping TO 295 */ 296 PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num)); 297 /* 298 Determine how many equations/vertices are assigned to each processor 299 */ 300 PetscCall(PetscMalloc1(size, &counts)); 301 PetscCall(ISPartitioningCount(is_eq_newproc, size, counts)); 302 ncrs_eq_new = counts[rank]; 303 PetscCall(ISDestroy(&is_eq_newproc)); 304 PetscCall(PetscFree(counts)); 305 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0)); 306 } else { /* simple aggregation of parts -- 'is_eq_newproc' */ 307 const PetscInt *ranges; 308 PetscInt newstart = 0; 309 PetscLayout ilay; 310 311 PetscCheck(new_size != nactive, PETSC_COMM_SELF, PETSC_ERR_PLIB, "new_size==nactive. Should not happen"); 312 PetscCall(PetscInfo(pc, "%s: Number of equations (loc) %" PetscInt_FMT " with simple aggregation\n", ((PetscObject)pc)->prefix, ncrs_eq)); 313 PetscCallMPI(MPI_Exscan(&ncrs_eq, &newstart, 1, MPIU_INT, MPI_SUM, comm)); 314 PetscCall(ISCreateStride(comm, ncrs_eq, newstart, 1, &is_eq_num)); 315 PetscCall(ISSetPermutation(is_eq_num)); 316 PetscCall(ISGetLayout(is_eq_num, &ilay)); 317 PetscCall(PetscLayoutGetRanges(ilay, &ranges)); 318 ncrs_eq_new = 0; 319 for (PetscInt r = 0; r < size; r++) 320 if (rank == (r / rfactor) * expand_factor) ncrs_eq_new += ranges[r + 1] - ranges[r]; 321 //targetPE = (rank / rfactor) * expand_factor; 322 //PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc)); 323 //PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num)); 324 //PetscCall(PetscMalloc1(size, &counts)); 325 //PetscCall(ISPartitioningCount(is_eq_newproc, size, counts)); 326 //ncrs_eq_new = counts[rank]; 327 //PetscCall(ISDestroy(&is_eq_newproc)); 328 //PetscCall(PetscFree(counts)); 329 } /* end simple 'is_eq_newproc' */ 330 331 ncrs_new = ncrs_eq_new / cr_bs; 332 333 /* 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 */ 334 { 335 Vec src_crd, dest_crd; 336 const PetscInt *idx, ndata_rows = pc_gamg->data_cell_rows, ndata_cols = pc_gamg->data_cell_cols, node_data_sz = ndata_rows * ndata_cols; 337 VecScatter vecscat; 338 PetscScalar *array; 339 IS isscat; 340 /* move data (for primal equations only) */ 341 /* Create a vector to contain the newly ordered element information */ 342 PetscCall(VecCreate(comm, &dest_crd)); 343 PetscCall(VecSetSizes(dest_crd, node_data_sz * ncrs_new, PETSC_DECIDE)); 344 PetscCall(VecSetType(dest_crd, VECSTANDARD)); /* this is needed! */ 345 /* 346 There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 347 a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 348 */ 349 PetscCall(PetscMalloc1(ncrs * node_data_sz, &tidx)); 350 PetscCall(ISGetIndices(is_eq_num, &idx)); 351 for (ii = 0, jj = 0; ii < ncrs; ii++) { 352 PetscInt id = idx[ii * cr_bs] / cr_bs; /* get node back */ 353 for (kk = 0; kk < node_data_sz; kk++, jj++) tidx[jj] = id * node_data_sz + kk; 354 } 355 PetscCall(ISRestoreIndices(is_eq_num, &idx)); 356 PetscCall(ISCreateGeneral(comm, node_data_sz * ncrs, tidx, PETSC_COPY_VALUES, &isscat)); 357 PetscCall(PetscFree(tidx)); 358 /* 359 Create a vector to contain the original vertex information for each element 360 */ 361 PetscCall(VecCreateSeq(PETSC_COMM_SELF, node_data_sz * ncrs, &src_crd)); 362 for (jj = 0; jj < ndata_cols; jj++) { 363 const PetscInt stride0 = ncrs * pc_gamg->data_cell_rows; 364 for (ii = 0; ii < ncrs; ii++) { 365 for (kk = 0; kk < ndata_rows; kk++) { 366 PetscInt ix = ii * ndata_rows + kk + jj * stride0, jx = ii * node_data_sz + kk * ndata_cols + jj; 367 PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 368 PetscCall(VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES)); 369 } 370 } 371 } 372 PetscCall(VecAssemblyBegin(src_crd)); 373 PetscCall(VecAssemblyEnd(src_crd)); 374 /* 375 Scatter the element vertex information (still in the original vertex ordering) 376 to the correct processor 377 */ 378 PetscCall(VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat)); 379 PetscCall(ISDestroy(&isscat)); 380 PetscCall(VecScatterBegin(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD)); 381 PetscCall(VecScatterEnd(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD)); 382 PetscCall(VecScatterDestroy(&vecscat)); 383 PetscCall(VecDestroy(&src_crd)); 384 /* 385 Put the element vertex data into a new allocation of the gdata->ele 386 */ 387 PetscCall(PetscFree(pc_gamg->data)); 388 PetscCall(PetscMalloc1(node_data_sz * ncrs_new, &pc_gamg->data)); 389 390 pc_gamg->data_sz = node_data_sz * ncrs_new; 391 strideNew = ncrs_new * ndata_rows; 392 393 PetscCall(VecGetArray(dest_crd, &array)); 394 for (jj = 0; jj < ndata_cols; jj++) { 395 for (ii = 0; ii < ncrs_new; ii++) { 396 for (kk = 0; kk < ndata_rows; kk++) { 397 PetscInt ix = ii * ndata_rows + kk + jj * strideNew, jx = ii * node_data_sz + kk * ndata_cols + jj; 398 pc_gamg->data[ix] = PetscRealPart(array[jx]); 399 } 400 } 401 } 402 PetscCall(VecRestoreArray(dest_crd, &array)); 403 PetscCall(VecDestroy(&dest_crd)); 404 } 405 /* move A and P (columns) with new layout */ 406 /* 407 Invert for MatCreateSubMatrix 408 */ 409 PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices)); 410 PetscCall(ISSort(new_eq_indices)); 411 PetscCall(ISSetBlockSize(new_eq_indices, cr_bs)); 412 if (Pcolumnperm) { 413 PetscCall(PetscObjectReference((PetscObject)new_eq_indices)); 414 *Pcolumnperm = new_eq_indices; 415 } 416 PetscCall(ISDestroy(&is_eq_num)); 417 418 /* 'a_Amat_crs' output */ 419 if (Cmat) { /* repartitioning from Cmat adjacency case */ 420 Mat mat; 421 PetscBool isset, isspd, isher; 422 #if !defined(PETSC_USE_COMPLEX) 423 PetscBool issym; 424 #endif 425 426 PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat)); 427 PetscCall(MatIsSPDKnown(Cmat, &isset, &isspd)); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ? 428 if (isset) PetscCall(MatSetOption(mat, MAT_SPD, isspd)); 429 else { 430 PetscCall(MatIsHermitianKnown(Cmat, &isset, &isher)); 431 if (isset) PetscCall(MatSetOption(mat, MAT_HERMITIAN, isher)); 432 else { 433 #if !defined(PETSC_USE_COMPLEX) 434 PetscCall(MatIsSymmetricKnown(Cmat, &isset, &issym)); 435 if (isset) PetscCall(MatSetOption(mat, MAT_SYMMETRIC, issym)); 436 #endif 437 } 438 } 439 *a_Amat_crs = mat; 440 } 441 442 /* prolongator */ 443 { 444 IS findices; 445 PetscInt Istart, Iend; 446 Mat Pnew; 447 448 PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend)); 449 PetscCall(ISCreateStride(comm, Iend - Istart, Istart, 1, &findices)); 450 PetscCall(ISSetBlockSize(findices, f_bs)); 451 PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew)); 452 PetscCall(ISDestroy(&findices)); 453 PetscCall(MatSetOption(Pnew, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 454 455 PetscCall(MatDestroy(a_P_inout)); 456 457 /* output - repartitioned */ 458 *a_P_inout = Pnew; 459 } 460 461 if (!Cmat) { /* simple repartitioning case */ 462 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 463 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 464 PetscCall(MatPtAP(Amat_fine, *a_P_inout, MAT_INITIAL_MATRIX, 2.0, a_Amat_crs)); 465 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0)); 466 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0)); 467 } 468 PetscCall(MatDestroy(&Cmat)); 469 PetscCall(ISDestroy(&new_eq_indices)); 470 471 *a_nactive_proc = new_size; /* output */ 472 473 /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */ 474 if (pc_gamg->cpu_pin_coarse_grids) { 475 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 476 static PetscInt llev = 2; 477 PetscCall(PetscInfo(pc, "%s: Pinning level %" PetscInt_FMT " to the CPU\n", ((PetscObject)pc)->prefix, llev++)); 478 #endif 479 PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE)); 480 PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE)); 481 if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */ 482 Mat A = *a_Amat_crs, P = *a_P_inout; 483 PetscMPIInt size; 484 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 485 if (size > 1) { 486 Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data, *p = (Mat_MPIAIJ *)P->data; 487 PetscCall(VecBindToCPU(a->lvec, PETSC_TRUE)); 488 PetscCall(VecBindToCPU(p->lvec, PETSC_TRUE)); 489 } 490 } 491 } 492 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0)); 493 } // processor reduce 494 PetscFunctionReturn(PETSC_SUCCESS); 495 } 496 497 // used in GEO 498 PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat *Gmat2) 499 { 500 const char *prefix; 501 char addp[32]; 502 PC_MG *mg = (PC_MG *)a_pc->data; 503 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 504 505 PetscFunctionBegin; 506 PetscCall(PCGetOptionsPrefix(a_pc, &prefix)); 507 PetscCall(PetscInfo(a_pc, "%s: Square Graph on level %" PetscInt_FMT "\n", ((PetscObject)a_pc)->prefix, pc_gamg->current_level + 1)); 508 PetscCall(MatProductCreate(Gmat1, Gmat1, NULL, Gmat2)); 509 PetscCall(MatSetOptionsPrefix(*Gmat2, prefix)); 510 PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_square_%" PetscInt_FMT "_", pc_gamg->current_level)); 511 PetscCall(MatAppendOptionsPrefix(*Gmat2, addp)); 512 if ((*Gmat2)->structurally_symmetric == PETSC_BOOL3_TRUE) { 513 PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AB)); 514 } else { 515 PetscCall(MatSetOption(Gmat1, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 516 PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AtB)); 517 } 518 PetscCall(MatProductSetFromOptions(*Gmat2)); 519 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0)); 520 PetscCall(MatProductSymbolic(*Gmat2)); 521 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0)); 522 PetscCall(MatProductClear(*Gmat2)); 523 /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */ 524 (*Gmat2)->assembled = PETSC_TRUE; 525 PetscFunctionReturn(PETSC_SUCCESS); 526 } 527 528 /* 529 PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 530 by setting data structures and options. 531 532 Input Parameter: 533 . pc - the preconditioner context 534 535 */ 536 static PetscErrorCode PCSetUp_GAMG(PC pc) 537 { 538 PC_MG *mg = (PC_MG *)pc->data; 539 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 540 Mat Pmat = pc->pmat; 541 PetscInt fine_level, level, level1, bs, M, N, qq, lidx, nASMBlocksArr[PETSC_MG_MAXLEVELS], cr_bs; 542 MPI_Comm comm; 543 PetscMPIInt rank, size, nactivepe; 544 Mat Aarr[PETSC_MG_MAXLEVELS], Parr[PETSC_MG_MAXLEVELS]; 545 IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 546 PetscBool is_last = PETSC_FALSE; 547 #if defined(PETSC_USE_INFO) 548 PetscLogDouble nnz0 = 0., nnztot = 0.; 549 MatInfo info; 550 #endif 551 552 PetscFunctionBegin; 553 PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 554 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 555 PetscCallMPI(MPI_Comm_size(comm, &size)); 556 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 557 if (pc->setupcalled) { 558 if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) { 559 /* reset everything */ 560 PetscCall(PCReset_MG(pc)); 561 pc->setupcalled = 0; 562 } else { 563 PC_MG_Levels **mglevels = mg->levels; 564 /* just do Galerkin grids */ 565 Mat B, dA, dB; 566 if (pc_gamg->Nlevels > 1) { 567 PetscInt gl; 568 /* currently only handle case where mat and pmat are the same on coarser levels */ 569 PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, &dA, &dB)); 570 /* (re)set to get dirty flag */ 571 PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, dA, dB)); 572 573 for (level = pc_gamg->Nlevels - 2, gl = 0; level >= 0; level--, gl++) { 574 MatReuse reuse = MAT_INITIAL_MATRIX; 575 #if defined(GAMG_STAGES) 576 PetscCall(PetscLogStagePush(gamg_stages[gl])); 577 #endif 578 /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 579 PetscCall(KSPGetOperators(mglevels[level]->smoothd, NULL, &B)); 580 if (B->product) { 581 if (B->product->A == dB && B->product->B == mglevels[level + 1]->interpolate) reuse = MAT_REUSE_MATRIX; 582 } 583 if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A)); 584 if (reuse == MAT_REUSE_MATRIX) { 585 PetscCall(PetscInfo(pc, "%s: RAP after initial setup, reuse matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 586 } else { 587 PetscCall(PetscInfo(pc, "%s: RAP after initial setup, with repartitioning (new matrix) level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 588 } 589 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0)); 590 PetscCall(MatPtAP(dB, mglevels[level + 1]->interpolate, reuse, PETSC_DETERMINE, &B)); 591 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0)); 592 if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B; 593 PetscCall(KSPSetOperators(mglevels[level]->smoothd, B, B)); 594 // check for redoing eigen estimates 595 if (pc_gamg->recompute_esteig) { 596 PetscBool ischeb; 597 KSP smoother; 598 PetscCall(PCMGGetSmoother(pc, level + 1, &smoother)); 599 PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb)); 600 if (ischeb) { 601 KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data; 602 cheb->emin_provided = 0; 603 cheb->emax_provided = 0; 604 } 605 /* we could call PetscCall(KSPChebyshevSetEigenvalues(smoother, 0, 0)); but the logic does not work properly */ 606 } 607 // inc 608 dB = B; 609 #if defined(GAMG_STAGES) 610 PetscCall(PetscLogStagePop()); 611 #endif 612 } 613 } 614 615 PetscCall(PCSetUp_MG(pc)); 616 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 617 PetscFunctionReturn(PETSC_SUCCESS); 618 } 619 } 620 621 if (!pc_gamg->data) { 622 if (pc_gamg->orig_data) { 623 PetscCall(MatGetBlockSize(Pmat, &bs)); 624 PetscCall(MatGetLocalSize(Pmat, &qq, NULL)); 625 626 pc_gamg->data_sz = (qq / bs) * pc_gamg->orig_data_cell_rows * pc_gamg->orig_data_cell_cols; 627 pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 628 pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 629 630 PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data)); 631 for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 632 } else { 633 PetscCheck(pc_gamg->ops->createdefaultdata, comm, PETSC_ERR_PLIB, "'createdefaultdata' not set(?) need to support NULL data"); 634 PetscCall(pc_gamg->ops->createdefaultdata(pc, Pmat)); 635 } 636 } 637 638 /* cache original data for reuse */ 639 if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 640 PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data)); 641 for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 642 pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 643 pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 644 } 645 646 /* get basic dims */ 647 PetscCall(MatGetBlockSize(Pmat, &bs)); 648 PetscCall(MatGetSize(Pmat, &M, NULL)); 649 650 #if defined(PETSC_USE_INFO) 651 PetscCall(MatGetInfo(Pmat, MAT_GLOBAL_SUM, &info)); /* global reduction */ 652 nnz0 = info.nz_used; 653 nnztot = info.nz_used; 654 #endif 655 PetscCall(PetscInfo(pc, "%s: level %d) N=%" PetscInt_FMT ", n data rows=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", block size %d, np=%d\n", ((PetscObject)pc)->prefix, 0, M, pc_gamg->data_cell_rows, pc_gamg->data_cell_cols, (PetscInt)(nnz0 / (PetscReal)M + 0.5), (int)bs, size)); 656 657 /* Get A_i and R_i */ 658 for (level = 0, Aarr[0] = Pmat, nactivepe = size; level < (pc_gamg->Nlevels - 1) && (level == 0 || M > pc_gamg->coarse_eq_limit); level++) { 659 pc_gamg->current_level = level; 660 level1 = level + 1; 661 #if defined(GAMG_STAGES) 662 if (!gamg_stages[level]) { 663 char str[32]; 664 PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", (int)level)); 665 PetscCall(PetscLogStageRegister(str, &gamg_stages[level])); 666 } 667 PetscCall(PetscLogStagePush(gamg_stages[level])); 668 #endif 669 /* construct prolongator - Parr[level1] */ 670 if (level == 0 && pc_gamg->injection_index_size > 0) { 671 Mat Prol; 672 MatType mtype; 673 PetscInt prol_m, prol_n, Prol_N = (M / bs) * pc_gamg->injection_index_size, Istart, Iend, nn, row; 674 PetscCall(PetscInfo(pc, "Create fine grid injection space prolongation %" PetscInt_FMT " x %" PetscInt_FMT ". %s\n", M, Prol_N, pc_gamg->data ? "delete null space data" : "")); 675 PetscCall(MatGetOwnershipRange(Pmat, &Istart, &Iend)); 676 PetscCall(MatGetLocalSize(Pmat, &prol_m, NULL)); // rows m x n 677 prol_n = (prol_m / bs) * pc_gamg->injection_index_size; 678 PetscCheck(pc_gamg->injection_index_size < bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Injection size %d must be less that block size %d", (int)pc_gamg->injection_index_size, (int)bs); 679 PetscCall(MatGetType(Pmat, &mtype)); 680 PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &Prol)); 681 PetscCall(MatSetBlockSizes(Prol, bs, pc_gamg->injection_index_size)); 682 PetscCall(MatSetSizes(Prol, prol_m, prol_n, M, Prol_N)); 683 PetscCall(MatSetType(Prol, mtype)); 684 #if PetscDefined(HAVE_DEVICE) 685 PetscBool flg; 686 PetscCall(MatBoundToCPU(Pmat, &flg)); 687 PetscCall(MatBindToCPU(Prol, flg)); 688 if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE)); 689 #endif 690 PetscCall(MatSeqAIJSetPreallocation(Prol, 1, NULL)); 691 PetscCall(MatMPIAIJSetPreallocation(Prol, 1, NULL, 0, NULL)); 692 // set I \kron [1, 1, ... ]^T 693 for (PetscInt ii = Istart, col = (Istart / bs) * pc_gamg->injection_index_size; ii < Iend; ii += bs) { 694 const PetscScalar one = 1; 695 for (PetscInt jj = 0; jj < pc_gamg->injection_index_size; jj++, col++) { 696 PetscInt row = ii + pc_gamg->injection_index[jj]; 697 PetscCall(MatSetValues(Prol, 1, &row, 1, &col, &one, INSERT_VALUES)); 698 } 699 } 700 PetscCall(MatAssemblyBegin(Prol, MAT_FINAL_ASSEMBLY)); 701 PetscCall(MatAssemblyEnd(Prol, MAT_FINAL_ASSEMBLY)); 702 PetscCall(MatViewFromOptions(Prol, NULL, "-mat_view_injection")); 703 PetscCall(MatGetBlockSizes(Prol, NULL, &cr_bs)); // column size 704 Parr[level1] = Prol; 705 // can not deal with null space -- with array of 'injection cols' we could take 'injection rows and 'injection cols' to 'data' 706 if (pc_gamg->data) { 707 pc_gamg->data_cell_cols = pc_gamg->injection_index_size; 708 pc_gamg->data_cell_rows = pc_gamg->injection_index_size; 709 pc_gamg->orig_data_cell_cols = 0; 710 pc_gamg->orig_data_cell_rows = 0; 711 PetscCall(PetscFree(pc_gamg->data)); 712 pc_gamg->data_sz = pc_gamg->injection_index_size * prol_n; 713 PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data)); 714 for (row = nn = 0; row < prol_n; row += pc_gamg->injection_index_size) { 715 for (int jj = 0; jj < pc_gamg->injection_index_size; jj++) { 716 int idx = row * pc_gamg->injection_index_size + jj * pc_gamg->injection_index_size; 717 for (int kk = 0; kk < pc_gamg->injection_index_size; kk++, nn++) { pc_gamg->data[idx + kk] = (jj == kk) ? 1 : 0; } 718 } 719 } 720 PetscCheck(nn == pc_gamg->data_sz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nn != pc_gamg->data_sz %" PetscInt_FMT " %" PetscInt_FMT, pc_gamg->data_sz, nn); 721 } 722 } else { 723 Mat Gmat, mat; 724 PetscCoarsenData *agg_lists; 725 Mat Prol11; 726 727 PetscCall(PCGAMGCreateGraph(pc, Aarr[level], &Gmat)); 728 PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists)); // Gmat may have ghosts for QR aggregates not in matrix 729 PetscCall(PetscCDGetMat(agg_lists, &mat)); 730 if (!mat) PetscCall(PetscCDSetMat(agg_lists, Gmat)); 731 PetscCall(pc_gamg->ops->prolongator(pc, Aarr[level], agg_lists, &Prol11)); 732 /* could have failed to create new level */ 733 if (Prol11) { 734 const char *prefix; 735 char addp[32]; 736 737 /* get new block size of coarse matrices */ 738 PetscCall(MatGetBlockSizes(Prol11, NULL, &cr_bs)); // column size 739 740 if (pc_gamg->ops->optprolongator) { 741 /* smooth */ 742 PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11)); 743 } 744 745 if (pc_gamg->use_aggs_in_asm) { 746 PetscInt bs; 747 PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // row block size 748 PetscCall(PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level])); 749 PetscCall(PetscInfo(pc, "%d: %" PetscInt_FMT " ASM local domains, bs = %d\n", (int)level, nASMBlocksArr[level], (int)bs)); 750 } else if (pc_gamg->asm_hem_aggs) { 751 const char *prefix; 752 PetscInt bs; 753 754 /* 755 Do not use aggs created for defining coarser problems, instead create aggs specifically to use 756 to define PCASM blocks. 757 */ 758 PetscCall(PetscCDGetMat(agg_lists, &mat)); 759 if (mat == Gmat) PetscCall(PetscCDClearMat(agg_lists)); // take the Mat away from the list (yuck) 760 PetscCall(PetscCDDestroy(agg_lists)); 761 PetscCall(PetscInfo(pc, "HEM ASM passes = %d\n", (int)pc_gamg->asm_hem_aggs)); 762 PetscCall(MatCoarsenDestroy(&pc_gamg->asm_crs)); 763 PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg->asm_crs)); 764 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 765 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg->asm_crs, prefix)); 766 PetscCall(MatCoarsenSetFromOptions(pc_gamg->asm_crs)); // get strength args 767 PetscCall(MatCoarsenSetType(pc_gamg->asm_crs, MATCOARSENHEM)); 768 PetscCall(MatCoarsenSetMaximumIterations(pc_gamg->asm_crs, pc_gamg->asm_hem_aggs)); 769 PetscCall(MatCoarsenSetAdjacency(pc_gamg->asm_crs, Gmat)); 770 PetscCall(MatCoarsenSetStrictAggs(pc_gamg->asm_crs, PETSC_TRUE)); 771 PetscCall(MatCoarsenApply(pc_gamg->asm_crs)); 772 PetscCall(MatCoarsenGetData(pc_gamg->asm_crs, &agg_lists)); /* output */ 773 // create aggregates 774 PetscCall(MatGetBlockSizes(Aarr[level], &bs, NULL)); // row block size 775 PetscCall(PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level])); 776 } 777 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 778 PetscCall(MatSetOptionsPrefix(Prol11, prefix)); 779 PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_prolongator_%d_", (int)level)); 780 PetscCall(MatAppendOptionsPrefix(Prol11, addp)); 781 /* Always generate the transpose with CUDA 782 Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */ 783 PetscCall(MatSetOption(Prol11, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 784 PetscCall(MatSetFromOptions(Prol11)); 785 Parr[level1] = Prol11; 786 } else Parr[level1] = NULL; /* failed to coarsen */ 787 PetscCall(PetscCDGetMat(agg_lists, &mat)); 788 if (mat == Gmat) PetscCall(PetscCDClearMat(agg_lists)); // take the Mat away from the list (yuck) 789 PetscCall(MatDestroy(&Gmat)); 790 PetscCall(PetscCDDestroy(agg_lists)); 791 } /* construct prolongator scope */ 792 if (level == 0) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 793 if (!Parr[level1]) { /* failed to coarsen */ 794 PetscCall(PetscInfo(pc, "%s: Stop gridding, level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 795 #if defined(GAMG_STAGES) 796 PetscCall(PetscLogStagePop()); 797 #endif 798 break; 799 } 800 PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */ 801 PetscCheck(!is_last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Is last ?"); 802 if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 803 if (level1 == pc_gamg->Nlevels - 1) is_last = PETSC_TRUE; 804 if (level == PETSC_MG_MAXLEVELS - 2) is_last = PETSC_TRUE; 805 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0)); 806 PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], cr_bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last)); 807 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0)); 808 809 PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */ 810 #if defined(PETSC_USE_INFO) 811 PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info)); 812 nnztot += info.nz_used; 813 #endif 814 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)); 815 816 #if defined(GAMG_STAGES) 817 PetscCall(PetscLogStagePop()); 818 #endif 819 /* stop if one node or one proc -- could pull back for singular problems */ 820 if ((pc_gamg->data_cell_cols && M / pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M / bs < 2)) { 821 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)); 822 level++; 823 break; 824 } else if (level == PETSC_MG_MAXLEVELS - 2) { /* stop if we are limited by PC_MG_MAXLEVELS */ 825 PetscCall(PetscInfo(pc, "%s: HARD stop of coarsening on level %" PetscInt_FMT ". PC_MG_MAXLEVELS reached\n", ((PetscObject)pc)->prefix, level)); 826 level++; 827 break; 828 } 829 } /* levels */ 830 PetscCall(PetscFree(pc_gamg->data)); 831 832 PetscCall(PetscInfo(pc, "%s: %" PetscInt_FMT " levels, operator complexity = %g\n", ((PetscObject)pc)->prefix, level + 1, nnztot / nnz0)); 833 pc_gamg->Nlevels = level + 1; 834 fine_level = level; 835 PetscCall(PCMGSetLevels(pc, pc_gamg->Nlevels, NULL)); 836 837 if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 838 839 /* set default smoothers & set operators */ 840 for (lidx = 1, level = pc_gamg->Nlevels - 2; lidx <= fine_level; lidx++, level--) { 841 KSP smoother; 842 PC subpc; 843 844 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 845 PetscCall(KSPGetPC(smoother, &subpc)); 846 847 PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 848 /* set ops */ 849 PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level])); 850 PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level + 1])); 851 852 /* set defaults */ 853 PetscCall(KSPSetType(smoother, KSPCHEBYSHEV)); 854 855 /* set blocks for ASM smoother that uses the 'aggregates' */ 856 if (pc_gamg->use_aggs_in_asm || pc_gamg->asm_hem_aggs) { 857 PetscInt sz; 858 IS *iss; 859 860 sz = nASMBlocksArr[level]; 861 iss = ASMLocalIDsArr[level]; 862 PetscCall(PCSetType(subpc, PCASM)); 863 PetscCall(PCASMSetOverlap(subpc, 0)); 864 PetscCall(PCASMSetType(subpc, PC_ASM_BASIC)); 865 if (!sz) { 866 IS is; 867 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is)); 868 PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is)); 869 PetscCall(ISDestroy(&is)); 870 } else { 871 PetscInt kk; 872 PetscCall(PCASMSetLocalSubdomains(subpc, sz, iss, NULL)); 873 for (kk = 0; kk < sz; kk++) PetscCall(ISDestroy(&iss[kk])); 874 PetscCall(PetscFree(iss)); 875 } 876 ASMLocalIDsArr[level] = NULL; 877 nASMBlocksArr[level] = 0; 878 } else { 879 PetscCall(PCSetType(subpc, PCJACOBI)); 880 } 881 } 882 { 883 /* coarse grid */ 884 KSP smoother, *k2; 885 PC subpc, pc2; 886 PetscInt ii, first; 887 Mat Lmat = Aarr[(level = pc_gamg->Nlevels - 1)]; 888 lidx = 0; 889 890 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 891 PetscCall(KSPSetOperators(smoother, Lmat, Lmat)); 892 if (!pc_gamg->use_parallel_coarse_grid_solver) { 893 PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 894 PetscCall(KSPGetPC(smoother, &subpc)); 895 PetscCall(PCSetType(subpc, PCBJACOBI)); 896 PetscCall(PCSetUp(subpc)); 897 PetscCall(PCBJacobiGetSubKSP(subpc, &ii, &first, &k2)); 898 PetscCheck(ii == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " is not one", ii); 899 PetscCall(KSPGetPC(k2[0], &pc2)); 900 PetscCall(PCSetType(pc2, PCLU)); 901 PetscCall(PCFactorSetShiftType(pc2, MAT_SHIFT_INBLOCKS)); 902 PetscCall(KSPSetTolerances(k2[0], PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, 1)); 903 PetscCall(KSPSetType(k2[0], KSPPREONLY)); 904 } 905 } 906 907 /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 908 PetscObjectOptionsBegin((PetscObject)pc); 909 PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); 910 PetscOptionsEnd(); 911 PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 912 913 /* set cheby eigen estimates from SA to use in the solver */ 914 if (pc_gamg->use_sa_esteig) { 915 for (lidx = 1, level = pc_gamg->Nlevels - 2; level >= 0; lidx++, level--) { 916 KSP smoother; 917 PetscBool ischeb; 918 919 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 920 PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb)); 921 if (ischeb) { 922 KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data; 923 924 // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG 925 if (mg->max_eigen_DinvA[level] > 0) { 926 // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it. 927 // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.) 928 PetscReal emax, emin; 929 930 emin = mg->min_eigen_DinvA[level]; 931 emax = mg->max_eigen_DinvA[level]; 932 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)); 933 cheb->emin_provided = emin; 934 cheb->emax_provided = emax; 935 } 936 } 937 } 938 } 939 940 PetscCall(PCSetUp_MG(pc)); 941 942 /* clean up */ 943 for (level = 1; level < pc_gamg->Nlevels; level++) { 944 PetscCall(MatDestroy(&Parr[level])); 945 PetscCall(MatDestroy(&Aarr[level])); 946 } 947 } else { 948 KSP smoother; 949 950 PetscCall(PetscInfo(pc, "%s: One level solver used (system is seen as DD). Using default solver.\n", ((PetscObject)pc)->prefix)); 951 PetscCall(PCMGGetSmoother(pc, 0, &smoother)); 952 PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0])); 953 PetscCall(KSPSetType(smoother, KSPPREONLY)); 954 PetscCall(PCSetUp_MG(pc)); 955 } 956 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 957 PetscFunctionReturn(PETSC_SUCCESS); 958 } 959 960 /* 961 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 962 that was created with PCCreate_GAMG(). 963 964 Input Parameter: 965 . pc - the preconditioner context 966 967 Application Interface Routine: PCDestroy() 968 */ 969 PetscErrorCode PCDestroy_GAMG(PC pc) 970 { 971 PC_MG *mg = (PC_MG *)pc->data; 972 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 973 974 PetscFunctionBegin; 975 PetscCall(PCReset_GAMG(pc)); 976 if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc)); 977 PetscCall(PetscFree(pc_gamg->ops)); 978 PetscCall(PetscFree(pc_gamg->gamg_type_name)); 979 PetscCall(PetscFree(pc_gamg)); 980 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", NULL)); 981 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", NULL)); 982 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", NULL)); 983 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", NULL)); 984 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", NULL)); 985 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRecomputeEstEig_C", NULL)); 986 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", NULL)); 987 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", NULL)); 988 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetParallelCoarseGridSolve_C", NULL)); 989 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", NULL)); 990 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", NULL)); 991 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", NULL)); 992 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", NULL)); 993 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", NULL)); 994 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", NULL)); 995 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", NULL)); 996 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", NULL)); 997 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetHEM_C", NULL)); 998 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetInjectionIndices_C", NULL)); 999 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetInjectionIndex_C", NULL)); 1000 PetscCall(PCDestroy_MG(pc)); 1001 PetscFunctionReturn(PETSC_SUCCESS); 1002 } 1003 1004 /*@ 1005 PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction in `PCGAMG` 1006 1007 Logically Collective 1008 1009 Input Parameters: 1010 + pc - the preconditioner context 1011 - n - the number of equations 1012 1013 Options Database Key: 1014 . -pc_gamg_process_eq_limit <limit> - set the limit 1015 1016 Level: intermediate 1017 1018 Note: 1019 `PCGAMG` will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process 1020 that has degrees of freedom 1021 1022 Developer Note: 1023 Should be named `PCGAMGSetProcessEquationLimit()`. 1024 1025 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()` 1026 @*/ 1027 PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 1028 { 1029 PetscFunctionBegin; 1030 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1031 PetscTryMethod(pc, "PCGAMGSetProcEqLim_C", (PC, PetscInt), (pc, n)); 1032 PetscFunctionReturn(PETSC_SUCCESS); 1033 } 1034 1035 static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 1036 { 1037 PC_MG *mg = (PC_MG *)pc->data; 1038 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1039 1040 PetscFunctionBegin; 1041 if (n > 0) pc_gamg->min_eq_proc = n; 1042 PetscFunctionReturn(PETSC_SUCCESS); 1043 } 1044 1045 /*@ 1046 PCGAMGSetCoarseEqLim - Set maximum number of equations on the coarsest grid of `PCGAMG` 1047 1048 Collective 1049 1050 Input Parameters: 1051 + pc - the preconditioner context 1052 - n - maximum number of equations to aim for 1053 1054 Options Database Key: 1055 . -pc_gamg_coarse_eq_limit <limit> - set the limit 1056 1057 Level: intermediate 1058 1059 Note: 1060 For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 1061 has less than 1000 unknowns. 1062 1063 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`, 1064 `PCGAMGSetParallelCoarseGridSolve()` 1065 @*/ 1066 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 1067 { 1068 PetscFunctionBegin; 1069 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1070 PetscTryMethod(pc, "PCGAMGSetCoarseEqLim_C", (PC, PetscInt), (pc, n)); 1071 PetscFunctionReturn(PETSC_SUCCESS); 1072 } 1073 1074 static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 1075 { 1076 PC_MG *mg = (PC_MG *)pc->data; 1077 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1078 1079 PetscFunctionBegin; 1080 if (n > 0) pc_gamg->coarse_eq_limit = n; 1081 PetscFunctionReturn(PETSC_SUCCESS); 1082 } 1083 1084 /*@ 1085 PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids when reducing the number of MPI ranks to use 1086 1087 Collective 1088 1089 Input Parameters: 1090 + pc - the preconditioner context 1091 - n - `PETSC_TRUE` or `PETSC_FALSE` 1092 1093 Options Database Key: 1094 . -pc_gamg_repartition <true,false> - turn on the repartitioning 1095 1096 Level: intermediate 1097 1098 Note: 1099 This will generally improve the loading balancing of the work on each level so the solves will be faster but it adds to the 1100 preconditioner setup time. 1101 1102 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()` 1103 @*/ 1104 PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) 1105 { 1106 PetscFunctionBegin; 1107 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1108 PetscTryMethod(pc, "PCGAMGSetRepartition_C", (PC, PetscBool), (pc, n)); 1109 PetscFunctionReturn(PETSC_SUCCESS); 1110 } 1111 1112 static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) 1113 { 1114 PC_MG *mg = (PC_MG *)pc->data; 1115 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1116 1117 PetscFunctionBegin; 1118 pc_gamg->repart = n; 1119 PetscFunctionReturn(PETSC_SUCCESS); 1120 } 1121 1122 /*@ 1123 PCGAMGSetUseSAEstEig - Use the eigen estimate from smoothed aggregation for the Chebyshev smoother during the solution process 1124 1125 Collective 1126 1127 Input Parameters: 1128 + pc - the preconditioner context 1129 - b - flag 1130 1131 Options Database Key: 1132 . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate 1133 1134 Level: advanced 1135 1136 Notes: 1137 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$. 1138 Eigenvalue estimates (based on a few `PCCG` or `PCGMRES` iterations) are computed to choose $\omega$ so that this is a stable smoothing operation. 1139 If `KSPCHEBYSHEV` with `PCJACOBI` (diagonal) preconditioning is used for smoothing on the finest level, then the eigenvalue estimates 1140 can be reused during the solution process. 1141 This option is only used when the smoother uses `PCJACOBI`, and should be turned off when a different `PCJacobiType` is used. 1142 It became default in PETSc 3.17. 1143 1144 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`, `PCGAMGSetRecomputeEstEig()` 1145 @*/ 1146 PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool b) 1147 { 1148 PetscFunctionBegin; 1149 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1150 PetscTryMethod(pc, "PCGAMGSetUseSAEstEig_C", (PC, PetscBool), (pc, b)); 1151 PetscFunctionReturn(PETSC_SUCCESS); 1152 } 1153 1154 static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool b) 1155 { 1156 PC_MG *mg = (PC_MG *)pc->data; 1157 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1158 1159 PetscFunctionBegin; 1160 pc_gamg->use_sa_esteig = b; 1161 PetscFunctionReturn(PETSC_SUCCESS); 1162 } 1163 1164 /*@ 1165 PCGAMGSetRecomputeEstEig - Set flag for Chebyshev smoothers to recompute the eigen estimates when a new matrix is used 1166 1167 Collective 1168 1169 Input Parameters: 1170 + pc - the preconditioner context 1171 - b - flag, default is `PETSC_TRUE` 1172 1173 Options Database Key: 1174 . -pc_gamg_recompute_esteig <true> - use the eigen estimate 1175 1176 Level: advanced 1177 1178 Note: 1179 If the matrix changes only slightly in a new solve using ``PETSC_FALSE`` will save time in the setting up of the preconditioner 1180 and may not affect the solution time much. 1181 1182 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()` 1183 @*/ 1184 PetscErrorCode PCGAMGSetRecomputeEstEig(PC pc, PetscBool b) 1185 { 1186 PetscFunctionBegin; 1187 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1188 PetscTryMethod(pc, "PCGAMGSetRecomputeEstEig_C", (PC, PetscBool), (pc, b)); 1189 PetscFunctionReturn(PETSC_SUCCESS); 1190 } 1191 1192 static PetscErrorCode PCGAMGSetRecomputeEstEig_GAMG(PC pc, PetscBool b) 1193 { 1194 PC_MG *mg = (PC_MG *)pc->data; 1195 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1196 1197 PetscFunctionBegin; 1198 pc_gamg->recompute_esteig = b; 1199 PetscFunctionReturn(PETSC_SUCCESS); 1200 } 1201 1202 /*@ 1203 PCGAMGSetEigenvalues - Set WHAT eigenvalues WHY? 1204 1205 Collective 1206 1207 Input Parameters: 1208 + pc - the preconditioner context 1209 . emax - max eigenvalue 1210 - emin - min eigenvalue 1211 1212 Options Database Key: 1213 . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues 1214 1215 Level: intermediate 1216 1217 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetUseSAEstEig()` 1218 @*/ 1219 PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax, PetscReal emin) 1220 { 1221 PetscFunctionBegin; 1222 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1223 PetscTryMethod(pc, "PCGAMGSetEigenvalues_C", (PC, PetscReal, PetscReal), (pc, emax, emin)); 1224 PetscFunctionReturn(PETSC_SUCCESS); 1225 } 1226 1227 static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc, PetscReal emax, PetscReal emin) 1228 { 1229 PC_MG *mg = (PC_MG *)pc->data; 1230 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1231 1232 PetscFunctionBegin; 1233 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); 1234 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); 1235 pc_gamg->emax = emax; 1236 pc_gamg->emin = emin; 1237 PetscFunctionReturn(PETSC_SUCCESS); 1238 } 1239 1240 /*@ 1241 PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding a `PCGAMG` algebraic multigrid preconditioner 1242 1243 Collective 1244 1245 Input Parameters: 1246 + pc - the preconditioner context 1247 - n - `PETSC_TRUE` or `PETSC_FALSE` 1248 1249 Options Database Key: 1250 . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation 1251 1252 Level: intermediate 1253 1254 Note: 1255 May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1256 rebuilding the preconditioner quicker. 1257 1258 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG` 1259 @*/ 1260 PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 1261 { 1262 PetscFunctionBegin; 1263 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1264 PetscTryMethod(pc, "PCGAMGSetReuseInterpolation_C", (PC, PetscBool), (pc, n)); 1265 PetscFunctionReturn(PETSC_SUCCESS); 1266 } 1267 1268 static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 1269 { 1270 PC_MG *mg = (PC_MG *)pc->data; 1271 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1272 1273 PetscFunctionBegin; 1274 pc_gamg->reuse_prol = n; 1275 PetscFunctionReturn(PETSC_SUCCESS); 1276 } 1277 1278 /*@ 1279 PCGAMGASMSetUseAggs - Have the `PCGAMG` smoother on each level use `PCASM` where the aggregates defined by the coarsening process are 1280 the subdomains for the additive Schwarz preconditioner used as the smoother 1281 1282 Collective 1283 1284 Input Parameters: 1285 + pc - the preconditioner context 1286 - flg - `PETSC_TRUE` to use aggregates, `PETSC_FALSE` to not 1287 1288 Options Database Key: 1289 . -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains 1290 1291 Level: intermediate 1292 1293 Note: 1294 This option automatically sets the preconditioner on the levels to be `PCASM`. 1295 1296 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCASM`, `PCSetType` 1297 @*/ 1298 PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) 1299 { 1300 PetscFunctionBegin; 1301 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1302 PetscTryMethod(pc, "PCGAMGASMSetUseAggs_C", (PC, PetscBool), (pc, flg)); 1303 PetscFunctionReturn(PETSC_SUCCESS); 1304 } 1305 1306 static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) 1307 { 1308 PC_MG *mg = (PC_MG *)pc->data; 1309 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1310 1311 PetscFunctionBegin; 1312 pc_gamg->use_aggs_in_asm = flg; 1313 PetscFunctionReturn(PETSC_SUCCESS); 1314 } 1315 1316 /*@ 1317 PCGAMGSetParallelCoarseGridSolve - allow a parallel coarse grid solver 1318 1319 Collective 1320 1321 Input Parameters: 1322 + pc - the preconditioner context 1323 - flg - `PETSC_TRUE` to not force coarse grid onto one processor 1324 1325 Options Database Key: 1326 . -pc_gamg_parallel_coarse_grid_solver - use a parallel coarse grid solver 1327 1328 Level: intermediate 1329 1330 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGSetRankReductionFactors()` 1331 @*/ 1332 PetscErrorCode PCGAMGSetParallelCoarseGridSolve(PC pc, PetscBool flg) 1333 { 1334 PetscFunctionBegin; 1335 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1336 PetscTryMethod(pc, "PCGAMGSetParallelCoarseGridSolve_C", (PC, PetscBool), (pc, flg)); 1337 PetscFunctionReturn(PETSC_SUCCESS); 1338 } 1339 1340 static PetscErrorCode PCGAMGSetParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) 1341 { 1342 PC_MG *mg = (PC_MG *)pc->data; 1343 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1344 1345 PetscFunctionBegin; 1346 pc_gamg->use_parallel_coarse_grid_solver = flg; 1347 PetscFunctionReturn(PETSC_SUCCESS); 1348 } 1349 1350 /*@ 1351 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 1352 1353 Collective 1354 1355 Input Parameters: 1356 + pc - the preconditioner context 1357 - flg - `PETSC_TRUE` to pin coarse grids to the CPU 1358 1359 Options Database Key: 1360 . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU 1361 1362 Level: intermediate 1363 1364 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetParallelCoarseGridSolve()` 1365 @*/ 1366 PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) 1367 { 1368 PetscFunctionBegin; 1369 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1370 PetscTryMethod(pc, "PCGAMGSetCpuPinCoarseGrids_C", (PC, PetscBool), (pc, flg)); 1371 PetscFunctionReturn(PETSC_SUCCESS); 1372 } 1373 1374 static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) 1375 { 1376 PC_MG *mg = (PC_MG *)pc->data; 1377 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1378 1379 PetscFunctionBegin; 1380 pc_gamg->cpu_pin_coarse_grids = flg; 1381 PetscFunctionReturn(PETSC_SUCCESS); 1382 } 1383 1384 /*@ 1385 PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type) 1386 1387 Collective 1388 1389 Input Parameters: 1390 + pc - the preconditioner context 1391 - flg - `PCGAMGLayoutType` type, either `PCGAMG_LAYOUT_COMPACT` or `PCGAMG_LAYOUT_SPREAD` 1392 1393 Options Database Key: 1394 . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering 1395 1396 Level: intermediate 1397 1398 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGLayoutType`, `PCGAMG_LAYOUT_COMPACT`, `PCGAMG_LAYOUT_SPREAD` 1399 @*/ 1400 PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) 1401 { 1402 PetscFunctionBegin; 1403 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1404 PetscTryMethod(pc, "PCGAMGSetCoarseGridLayoutType_C", (PC, PCGAMGLayoutType), (pc, flg)); 1405 PetscFunctionReturn(PETSC_SUCCESS); 1406 } 1407 1408 static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) 1409 { 1410 PC_MG *mg = (PC_MG *)pc->data; 1411 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1412 1413 PetscFunctionBegin; 1414 pc_gamg->layout_type = flg; 1415 PetscFunctionReturn(PETSC_SUCCESS); 1416 } 1417 1418 /*@ 1419 PCGAMGSetNlevels - Sets the maximum number of levels `PCGAMG` will use 1420 1421 Collective 1422 1423 Input Parameters: 1424 + pc - the preconditioner 1425 - n - the maximum number of levels to use 1426 1427 Options Database Key: 1428 . -pc_mg_levels <n> - set the maximum number of levels to allow 1429 1430 Level: intermediate 1431 1432 Developer Notes: 1433 Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG` 1434 1435 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG` 1436 @*/ 1437 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 1438 { 1439 PetscFunctionBegin; 1440 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1441 PetscTryMethod(pc, "PCGAMGSetNlevels_C", (PC, PetscInt), (pc, n)); 1442 PetscFunctionReturn(PETSC_SUCCESS); 1443 } 1444 1445 static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 1446 { 1447 PC_MG *mg = (PC_MG *)pc->data; 1448 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1449 1450 PetscFunctionBegin; 1451 pc_gamg->Nlevels = n; 1452 PetscFunctionReturn(PETSC_SUCCESS); 1453 } 1454 1455 /*@ 1456 PCGAMGASMSetHEM - Sets the number of HEM matching passed 1457 1458 Collective 1459 1460 Input Parameters: 1461 + pc - the preconditioner 1462 - n - number of HEM matching passed to construct ASM subdomains 1463 1464 Options Database Key: 1465 . -pc_gamg_asm_hem <n> - set the number of HEM matching passed 1466 1467 Level: intermediate 1468 1469 Developer Notes: 1470 Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG` 1471 1472 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG` 1473 @*/ 1474 PetscErrorCode PCGAMGASMSetHEM(PC pc, PetscInt n) 1475 { 1476 PetscFunctionBegin; 1477 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1478 PetscTryMethod(pc, "PCGAMGASMSetHEM_C", (PC, PetscInt), (pc, n)); 1479 PetscFunctionReturn(PETSC_SUCCESS); 1480 } 1481 1482 static PetscErrorCode PCGAMGASMSetHEM_GAMG(PC pc, PetscInt n) 1483 { 1484 PC_MG *mg = (PC_MG *)pc->data; 1485 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1486 1487 PetscFunctionBegin; 1488 pc_gamg->asm_hem_aggs = n; 1489 PetscFunctionReturn(PETSC_SUCCESS); 1490 } 1491 1492 /*@ 1493 PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 1494 1495 Not Collective 1496 1497 Input Parameters: 1498 + pc - the preconditioner context 1499 . 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 1500 - n - number of threshold values provided in array 1501 1502 Options Database Key: 1503 . -pc_gamg_threshold <threshold> - the threshold to drop edges 1504 1505 Level: intermediate 1506 1507 Notes: 1508 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. 1509 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. 1510 1511 If `n` is less than the total number of coarsenings (see `PCGAMGSetNlevels()`), then threshold scaling (see `PCGAMGSetThresholdScale()`) is used for each successive coarsening. 1512 In this case, `PCGAMGSetThresholdScale()` must be called before `PCGAMGSetThreshold()`. 1513 If `n` is greater than the total number of levels, the excess entries in threshold will not be used. 1514 1515 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetMinDegreeOrderingMISk()`, `PCGAMGSetThresholdScale()` 1516 @*/ 1517 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) 1518 { 1519 PetscFunctionBegin; 1520 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1521 if (n) PetscAssertPointer(v, 2); 1522 PetscTryMethod(pc, "PCGAMGSetThreshold_C", (PC, PetscReal[], PetscInt), (pc, v, n)); 1523 PetscFunctionReturn(PETSC_SUCCESS); 1524 } 1525 1526 static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) 1527 { 1528 PC_MG *mg = (PC_MG *)pc->data; 1529 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1530 PetscInt i; 1531 1532 PetscFunctionBegin; 1533 for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i]; 1534 for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale; 1535 PetscFunctionReturn(PETSC_SUCCESS); 1536 } 1537 1538 /*@ 1539 PCGAMGSetRankReductionFactors - Set a manual schedule for MPI rank reduction on coarse grids 1540 1541 Collective 1542 1543 Input Parameters: 1544 + pc - the preconditioner context 1545 . v - array of reduction factors. 0 for first value forces a reduction to one process/device on first level in CUDA 1546 - n - number of values provided in array 1547 1548 Options Database Key: 1549 . -pc_gamg_rank_reduction_factors <factors> - provide the schedule 1550 1551 Level: intermediate 1552 1553 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetParallelCoarseGridSolve()` 1554 @*/ 1555 PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n) 1556 { 1557 PetscFunctionBegin; 1558 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1559 if (n) PetscAssertPointer(v, 2); 1560 PetscTryMethod(pc, "PCGAMGSetRankReductionFactors_C", (PC, PetscInt[], PetscInt), (pc, v, n)); 1561 PetscFunctionReturn(PETSC_SUCCESS); 1562 } 1563 1564 static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n) 1565 { 1566 PC_MG *mg = (PC_MG *)pc->data; 1567 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1568 PetscInt i; 1569 1570 PetscFunctionBegin; 1571 for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i]; 1572 for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */ 1573 PetscFunctionReturn(PETSC_SUCCESS); 1574 } 1575 1576 /*@ 1577 PCGAMGSetThresholdScale - Relative threshold reduction at each level 1578 1579 Not Collective 1580 1581 Input Parameters: 1582 + pc - the preconditioner context 1583 - v - the threshold value reduction, usually < 1.0 1584 1585 Options Database Key: 1586 . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level 1587 1588 Level: advanced 1589 1590 Note: 1591 The initial threshold (for an arbitrary number of levels starting from the finest) can be set with `PCGAMGSetThreshold()`. 1592 This scaling is used for each subsequent coarsening, but must be called before `PCGAMGSetThreshold()`. 1593 1594 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()` 1595 @*/ 1596 PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) 1597 { 1598 PetscFunctionBegin; 1599 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1600 PetscTryMethod(pc, "PCGAMGSetThresholdScale_C", (PC, PetscReal), (pc, v)); 1601 PetscFunctionReturn(PETSC_SUCCESS); 1602 } 1603 1604 static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) 1605 { 1606 PC_MG *mg = (PC_MG *)pc->data; 1607 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1608 1609 PetscFunctionBegin; 1610 pc_gamg->threshold_scale = v; 1611 PetscFunctionReturn(PETSC_SUCCESS); 1612 } 1613 1614 /*@ 1615 PCGAMGSetType - Set the type of algorithm `PCGAMG` should use 1616 1617 Collective 1618 1619 Input Parameters: 1620 + pc - the preconditioner context 1621 - type - `PCGAMGAGG`, `PCGAMGGEO`, or `PCGAMGCLASSICAL` 1622 1623 Options Database Key: 1624 . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply - only agg is supported 1625 1626 Level: intermediate 1627 1628 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType` 1629 @*/ 1630 PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1631 { 1632 PetscFunctionBegin; 1633 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1634 PetscTryMethod(pc, "PCGAMGSetType_C", (PC, PCGAMGType), (pc, type)); 1635 PetscFunctionReturn(PETSC_SUCCESS); 1636 } 1637 1638 /*@ 1639 PCGAMGGetType - Get the type of algorithm `PCGAMG` will use 1640 1641 Collective 1642 1643 Input Parameter: 1644 . pc - the preconditioner context 1645 1646 Output Parameter: 1647 . type - the type of algorithm used 1648 1649 Level: intermediate 1650 1651 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetType()`, `PCGAMGType` 1652 @*/ 1653 PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1654 { 1655 PetscFunctionBegin; 1656 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1657 PetscUseMethod(pc, "PCGAMGGetType_C", (PC, PCGAMGType *), (pc, type)); 1658 PetscFunctionReturn(PETSC_SUCCESS); 1659 } 1660 1661 static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1662 { 1663 PC_MG *mg = (PC_MG *)pc->data; 1664 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1665 1666 PetscFunctionBegin; 1667 *type = pc_gamg->type; 1668 PetscFunctionReturn(PETSC_SUCCESS); 1669 } 1670 1671 static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1672 { 1673 PC_MG *mg = (PC_MG *)pc->data; 1674 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1675 PetscErrorCode (*r)(PC); 1676 1677 PetscFunctionBegin; 1678 pc_gamg->type = type; 1679 PetscCall(PetscFunctionListFind(GAMGList, type, &r)); 1680 PetscCheck(r, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown GAMG type %s given", type); 1681 if (pc_gamg->ops->destroy) { 1682 PetscCall((*pc_gamg->ops->destroy)(pc)); 1683 PetscCall(PetscMemzero(pc_gamg->ops, sizeof(struct _PCGAMGOps))); 1684 pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1685 /* cleaning up common data in pc_gamg - this should disappear someday */ 1686 pc_gamg->data_cell_cols = 0; 1687 pc_gamg->data_cell_rows = 0; 1688 pc_gamg->orig_data_cell_cols = 0; 1689 pc_gamg->orig_data_cell_rows = 0; 1690 PetscCall(PetscFree(pc_gamg->data)); 1691 pc_gamg->data_sz = 0; 1692 } 1693 PetscCall(PetscFree(pc_gamg->gamg_type_name)); 1694 PetscCall(PetscStrallocpy(type, &pc_gamg->gamg_type_name)); 1695 PetscCall((*r)(pc)); 1696 PetscFunctionReturn(PETSC_SUCCESS); 1697 } 1698 1699 static PetscErrorCode PCView_GAMG(PC pc, PetscViewer viewer) 1700 { 1701 PC_MG *mg = (PC_MG *)pc->data; 1702 PC_MG_Levels **mglevels = mg->levels; 1703 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1704 PetscReal gc, oc; 1705 1706 PetscFunctionBegin; 1707 PetscCall(PetscViewerASCIIPrintf(viewer, " GAMG specific options\n")); 1708 PetscCall(PetscViewerASCIIPrintf(viewer, " Threshold for dropping small values in graph on each level =")); 1709 for (PetscInt i = 0; i < mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)pc_gamg->threshold[i])); 1710 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 1711 PetscCall(PetscViewerASCIIPrintf(viewer, " Threshold scaling factor for each level not specified = %g\n", (double)pc_gamg->threshold_scale)); 1712 if (pc_gamg->use_aggs_in_asm) PetscCall(PetscViewerASCIIPrintf(viewer, " Using aggregates from GAMG coarsening process to define subdomains for PCASM\n")); // this take presedence 1713 else if (pc_gamg->asm_hem_aggs) { 1714 PetscCall(PetscViewerASCIIPrintf(viewer, " Using aggregates made with %d applications of heavy edge matching (HEM) to define subdomains for PCASM\n", (int)pc_gamg->asm_hem_aggs)); 1715 PetscCall(PetscViewerASCIIPushTab(viewer)); 1716 PetscCall(PetscViewerASCIIPushTab(viewer)); 1717 PetscCall(PetscViewerASCIIPushTab(viewer)); 1718 PetscCall(PetscViewerASCIIPushTab(viewer)); 1719 PetscCall(MatCoarsenView(pc_gamg->asm_crs, viewer)); 1720 PetscCall(PetscViewerASCIIPopTab(viewer)); 1721 PetscCall(PetscViewerASCIIPopTab(viewer)); 1722 PetscCall(PetscViewerASCIIPopTab(viewer)); 1723 } 1724 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")); 1725 if (pc_gamg->injection_index_size) { 1726 PetscCall(PetscViewerASCIIPrintf(viewer, " Using injection restriction/prolongation on first level, dofs:")); 1727 for (int i = 0; i < pc_gamg->injection_index_size; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %d", (int)pc_gamg->injection_index[i])); 1728 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 1729 } 1730 if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc, viewer)); 1731 gc = oc = 0; 1732 PetscCall(PCMGGetGridComplexity(pc, &gc, &oc)); 1733 PetscCall(PetscViewerASCIIPrintf(viewer, " Complexity: grid = %g operator = %g\n", (double)gc, (double)oc)); 1734 PetscCall(PetscViewerASCIIPrintf(viewer, " Per-level complexity: op = operator, int = interpolation\n")); 1735 PetscCall(PetscViewerASCIIPrintf(viewer, " #equations | #active PEs | avg nnz/row op | avg nnz/row int\n")); 1736 for (PetscInt i = 0; i < mg->nlevels; i++) { 1737 MatInfo info; 1738 Mat A; 1739 PetscReal rd[3]; 1740 PetscInt rst, ren, N; 1741 1742 PetscCall(KSPGetOperators(mglevels[i]->smoothd, NULL, &A)); 1743 PetscCall(MatGetOwnershipRange(A, &rst, &ren)); 1744 PetscCall(MatGetSize(A, &N, NULL)); 1745 PetscCall(MatGetInfo(A, MAT_LOCAL, &info)); 1746 rd[0] = (ren - rst > 0) ? 1 : 0; 1747 rd[1] = info.nz_used; 1748 rd[2] = 0; 1749 if (i) { 1750 Mat P; 1751 PetscCall(PCMGGetInterpolation(pc, i, &P)); 1752 PetscCall(MatGetInfo(P, MAT_LOCAL, &info)); 1753 rd[2] = info.nz_used; 1754 } 1755 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, rd, 3, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)pc))); 1756 PetscCall(PetscViewerASCIIPrintf(viewer, " %12" PetscInt_FMT " %12" PetscInt_FMT " %12" PetscInt_FMT " %12" PetscInt_FMT "\n", N, (PetscInt)rd[0], (PetscInt)PetscCeilReal(rd[1] / N), (PetscInt)PetscCeilReal(rd[2] / N))); 1757 } 1758 PetscFunctionReturn(PETSC_SUCCESS); 1759 } 1760 1761 /*@ 1762 PCGAMGSetInjectionIndex - Array of subset of variables per vertex to inject into coarse grid space 1763 1764 Logically Collective 1765 1766 Input Parameters: 1767 + pc - the coarsen context 1768 . n - number of indices 1769 - idx - array of indices 1770 1771 Options Database Key: 1772 . -pc_gamg_injection_index - array of subset of variables per vertex to use for injection coarse grid space 1773 1774 Level: intermediate 1775 1776 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG` 1777 @*/ 1778 PetscErrorCode PCGAMGSetInjectionIndex(PC pc, PetscInt n, PetscInt idx[]) 1779 { 1780 PetscFunctionBegin; 1781 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1782 PetscValidLogicalCollectiveInt(pc, n, 2); 1783 PetscTryMethod(pc, "PCGAMGSetInjectionIndex_C", (PC, PetscInt, PetscInt[]), (pc, n, idx)); 1784 PetscFunctionReturn(PETSC_SUCCESS); 1785 } 1786 1787 static PetscErrorCode PCGAMGSetInjectionIndex_GAMG(PC pc, PetscInt n, PetscInt idx[]) 1788 { 1789 PC_MG *mg = (PC_MG *)pc->data; 1790 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1791 1792 PetscFunctionBegin; 1793 pc_gamg->injection_index_size = n; 1794 PetscCheck(n < MAT_COARSEN_STRENGTH_INDEX_SIZE, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "array size %d larger than max %d", (int)n, MAT_COARSEN_STRENGTH_INDEX_SIZE); 1795 for (PetscInt i = 0; i < n; i++) pc_gamg->injection_index[i] = idx[i]; 1796 PetscFunctionReturn(PETSC_SUCCESS); 1797 } 1798 1799 static PetscErrorCode PCSetFromOptions_GAMG(PC pc, PetscOptionItems *PetscOptionsObject) 1800 { 1801 PC_MG *mg = (PC_MG *)pc->data; 1802 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1803 PetscBool flag; 1804 MPI_Comm comm; 1805 char prefix[256], tname[32]; 1806 PetscInt i, n; 1807 const char *pcpre; 1808 static const char *LayoutTypes[] = {"compact", "spread", "PCGAMGLayoutType", "PC_GAMG_LAYOUT", NULL}; 1809 1810 PetscFunctionBegin; 1811 PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 1812 PetscOptionsHeadBegin(PetscOptionsObject, "GAMG options"); 1813 PetscCall(PetscOptionsInt("-pc_mg_levels", "Number of Levels", "PCMGSetLevels", -1, &n, &flag)); 1814 PetscCheck(!flag, comm, PETSC_ERR_ARG_WRONG, "Invalid flag -pc_mg_levels. GAMG does not allow the number of levels to be set."); 1815 PetscCall(PetscOptionsFList("-pc_gamg_type", "Type of AMG method (only 'agg' supported and useful)", "PCGAMGSetType", GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag)); 1816 if (flag) PetscCall(PCGAMGSetType(pc, tname)); 1817 PetscCall(PetscOptionsBool("-pc_gamg_repartition", "Repartion coarse grids", "PCGAMGSetRepartition", pc_gamg->repart, &pc_gamg->repart, NULL)); 1818 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)); 1819 PetscCall(PetscOptionsBool("-pc_gamg_recompute_esteig", "Set flag to recompute eigen estimates for Chebyshev when matrix changes", "PCGAMGSetRecomputeEstEig", pc_gamg->recompute_esteig, &pc_gamg->recompute_esteig, NULL)); 1820 PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation", "Reuse prolongation operator", "PCGAMGReuseInterpolation", pc_gamg->reuse_prol, &pc_gamg->reuse_prol, NULL)); 1821 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)); 1822 PetscCall(PetscOptionsBool("-pc_gamg_parallel_coarse_grid_solver", "Use parallel coarse grid solver (otherwise put last grid on one process)", "PCGAMGSetParallelCoarseGridSolve", pc_gamg->use_parallel_coarse_grid_solver, &pc_gamg->use_parallel_coarse_grid_solver, NULL)); 1823 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)); 1824 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, 1825 (PetscEnum)pc_gamg->layout_type, (PetscEnum *)&pc_gamg->layout_type, NULL)); 1826 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)); 1827 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)); 1828 PetscCall(PetscOptionsInt("-pc_gamg_asm_hem_aggs", "Number of HEM matching passed in aggregates for ASM smoother", "PCGAMGASMSetHEM", pc_gamg->asm_hem_aggs, &pc_gamg->asm_hem_aggs, NULL)); 1829 PetscCall(PetscOptionsReal("-pc_gamg_threshold_scale", "Scaling of threshold for each level not specified", "PCGAMGSetThresholdScale", pc_gamg->threshold_scale, &pc_gamg->threshold_scale, NULL)); 1830 n = PETSC_MG_MAXLEVELS; 1831 PetscCall(PetscOptionsRealArray("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", pc_gamg->threshold, &n, &flag)); 1832 if (!flag || n < PETSC_MG_MAXLEVELS) { 1833 if (!flag) n = 1; 1834 i = n; 1835 do { 1836 pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale; 1837 } while (++i < PETSC_MG_MAXLEVELS); 1838 } 1839 PetscCall(PetscOptionsInt("-pc_mg_levels", "Set number of MG levels (should get from base class)", "PCGAMGSetNlevels", pc_gamg->Nlevels, &pc_gamg->Nlevels, NULL)); 1840 PetscCheck(pc_gamg->Nlevels <= PETSC_MG_MAXLEVELS, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_mg_levels (%d) >= PETSC_MG_MAXLEVELS (%d)", (int)pc_gamg->Nlevels, (int)PETSC_MG_MAXLEVELS); 1841 n = PETSC_MG_MAXLEVELS; 1842 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)); 1843 if (!flag) i = 0; 1844 else i = n; 1845 do { 1846 pc_gamg->level_reduction_factors[i] = -1; 1847 } while (++i < PETSC_MG_MAXLEVELS); 1848 { 1849 PetscReal eminmax[2] = {0., 0.}; 1850 n = 2; 1851 PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues", "extreme eigenvalues for smoothed aggregation", "PCGAMGSetEigenvalues", eminmax, &n, &flag)); 1852 if (flag) { 1853 PetscCheck(n == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 1854 PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0])); 1855 } 1856 } 1857 pc_gamg->injection_index_size = MAT_COARSEN_STRENGTH_INDEX_SIZE; 1858 PetscCall(PetscOptionsIntArray("-pc_gamg_injection_index", "Array of indices to use to use injection coarse grid space", "PCGAMGSetInjectionIndex", pc_gamg->injection_index, &pc_gamg->injection_index_size, NULL)); 1859 /* set options for subtype */ 1860 PetscCall((*pc_gamg->ops->setfromoptions)(pc, PetscOptionsObject)); 1861 1862 PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 1863 PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : "")); 1864 PetscOptionsHeadEnd(); 1865 PetscFunctionReturn(PETSC_SUCCESS); 1866 } 1867 1868 /*MC 1869 PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 1870 1871 Options Database Keys: 1872 + -pc_gamg_type <type,default=agg> - one of agg, geo, or classical (only smoothed aggregation, agg, supported) 1873 . -pc_gamg_repartition <bool,default=false> - repartition the degrees of freedom across the coarse grids as they are determined 1874 . -pc_gamg_asm_use_agg <bool,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother 1875 . -pc_gamg_process_eq_limit <limit, default=50> - `PCGAMG` will reduce the number of MPI ranks used directly on the coarse grids so that there are around <limit> 1876 equations on each process that has degrees of freedom 1877 . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 1878 . -pc_gamg_reuse_interpolation <bool,default=true> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations (should always be true) 1879 . -pc_gamg_threshold[] <thresh,default=[-1,...]> - Before aggregating the graph `PCGAMG` will remove small values from the graph on each level (< 0 does no filtering) 1880 - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1881 1882 Options Database Keys for Aggregation: 1883 + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation to construct prolongation 1884 . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest. 1885 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening 1886 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm 1887 . -pc_gamg_pc_gamg_asm_hem_aggs <n,default=0> - Number of HEM aggregation steps for ASM smoother 1888 - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive') 1889 1890 Options Database Keys for Multigrid: 1891 + -pc_mg_cycle_type <v> - v or w, see `PCMGSetCycleType()` 1892 . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1893 . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1894 - -pc_mg_levels <levels> - Number of levels of multigrid to use. GAMG has a heuristic so pc_mg_levels is not usually used with GAMG 1895 1896 Level: intermediate 1897 1898 Notes: 1899 To obtain good performance for `PCGAMG` for vector valued problems you must 1900 call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point 1901 call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator 1902 1903 The many options for `PCMG` also work directly for `PCGAMG` such as controlling the smoothers on each level etc. 1904 1905 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCCreate()`, `PCSetType()`, 1906 `MatSetBlockSize()`, 1907 `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, 1908 `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, 1909 `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()` 1910 M*/ 1911 PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 1912 { 1913 PC_GAMG *pc_gamg; 1914 PC_MG *mg; 1915 1916 PetscFunctionBegin; 1917 /* register AMG type */ 1918 PetscCall(PCGAMGInitializePackage()); 1919 1920 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 1921 PetscCall(PCSetType(pc, PCMG)); 1922 PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG)); 1923 1924 /* create a supporting struct and attach it to pc */ 1925 PetscCall(PetscNew(&pc_gamg)); 1926 PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 1927 mg = (PC_MG *)pc->data; 1928 mg->innerctx = pc_gamg; 1929 1930 PetscCall(PetscNew(&pc_gamg->ops)); 1931 1932 /* these should be in subctx but repartitioning needs simple arrays */ 1933 pc_gamg->data_sz = 0; 1934 pc_gamg->data = NULL; 1935 1936 /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 1937 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 1938 pc->ops->setup = PCSetUp_GAMG; 1939 pc->ops->reset = PCReset_GAMG; 1940 pc->ops->destroy = PCDestroy_GAMG; 1941 mg->view = PCView_GAMG; 1942 1943 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", PCGAMGSetProcEqLim_GAMG)); 1944 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", PCGAMGSetCoarseEqLim_GAMG)); 1945 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", PCGAMGSetRepartition_GAMG)); 1946 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", PCGAMGSetEigenvalues_GAMG)); 1947 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", PCGAMGSetUseSAEstEig_GAMG)); 1948 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRecomputeEstEig_C", PCGAMGSetRecomputeEstEig_GAMG)); 1949 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", PCGAMGSetReuseInterpolation_GAMG)); 1950 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", PCGAMGASMSetUseAggs_GAMG)); 1951 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetParallelCoarseGridSolve_C", PCGAMGSetParallelCoarseGridSolve_GAMG)); 1952 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", PCGAMGSetCpuPinCoarseGrids_GAMG)); 1953 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", PCGAMGSetCoarseGridLayoutType_GAMG)); 1954 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", PCGAMGSetThreshold_GAMG)); 1955 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", PCGAMGSetRankReductionFactors_GAMG)); 1956 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", PCGAMGSetThresholdScale_GAMG)); 1957 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", PCGAMGSetType_GAMG)); 1958 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", PCGAMGGetType_GAMG)); 1959 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", PCGAMGSetNlevels_GAMG)); 1960 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetHEM_C", PCGAMGASMSetHEM_GAMG)); 1961 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetInjectionIndex_C", PCGAMGSetInjectionIndex_GAMG)); 1962 pc_gamg->repart = PETSC_FALSE; 1963 pc_gamg->reuse_prol = PETSC_TRUE; 1964 pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1965 pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1966 pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1967 pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1968 pc_gamg->min_eq_proc = 50; 1969 pc_gamg->asm_hem_aggs = 0; 1970 pc_gamg->coarse_eq_limit = 50; 1971 for (int i = 0; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = -1; 1972 pc_gamg->threshold_scale = 1.; 1973 pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 1974 pc_gamg->current_level = 0; /* don't need to init really */ 1975 pc_gamg->use_sa_esteig = PETSC_TRUE; 1976 pc_gamg->recompute_esteig = PETSC_TRUE; 1977 pc_gamg->emin = 0; 1978 pc_gamg->emax = 0; 1979 1980 pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1981 1982 /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1983 PetscCall(PCGAMGSetType(pc, PCGAMGAGG)); 1984 PetscFunctionReturn(PETSC_SUCCESS); 1985 } 1986 1987 /*@C 1988 PCGAMGInitializePackage - This function initializes everything in the `PCGAMG` package. It is called 1989 from `PCInitializePackage()`. 1990 1991 Level: developer 1992 1993 .seealso: [](ch_ksp), `PetscInitialize()` 1994 @*/ 1995 PetscErrorCode PCGAMGInitializePackage(void) 1996 { 1997 PetscInt l; 1998 1999 PetscFunctionBegin; 2000 if (PCGAMGPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 2001 PCGAMGPackageInitialized = PETSC_TRUE; 2002 PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGGEO, PCCreateGAMG_GEO)); 2003 PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGAGG, PCCreateGAMG_AGG)); 2004 PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGCLASSICAL, PCCreateGAMG_Classical)); 2005 PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage)); 2006 2007 /* general events */ 2008 PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP])); 2009 PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH])); 2010 PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN])); 2011 PetscCall(PetscLogEventRegister(" GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS])); 2012 PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL])); 2013 PetscCall(PetscLogEventRegister(" GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA])); 2014 PetscCall(PetscLogEventRegister(" GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB])); 2015 PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT])); 2016 PetscCall(PetscLogEventRegister(" GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM])); 2017 PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL])); 2018 PetscCall(PetscLogEventRegister(" GAMG PtAP", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP])); 2019 PetscCall(PetscLogEventRegister(" GAMG Reduce", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE])); 2020 PetscCall(PetscLogEventRegister(" GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART])); 2021 /* PetscCall(PetscLogEventRegister(" GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */ 2022 /* PetscCall(PetscLogEventRegister(" GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */ 2023 /* PetscCall(PetscLogEventRegister(" GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */ 2024 for (l = 0; l < PETSC_MG_MAXLEVELS; l++) { 2025 char ename[32]; 2026 2027 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Squ l%02" PetscInt_FMT, l)); 2028 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0])); 2029 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Gal l%02" PetscInt_FMT, l)); 2030 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1])); 2031 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Opt l%02" PetscInt_FMT, l)); 2032 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2])); 2033 } 2034 #if defined(GAMG_STAGES) 2035 { /* create timer stages */ 2036 char str[32]; 2037 PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", 0)); 2038 PetscCall(PetscLogStageRegister(str, &gamg_stages[0])); 2039 } 2040 #endif 2041 PetscFunctionReturn(PETSC_SUCCESS); 2042 } 2043 2044 /*@C 2045 PCGAMGFinalizePackage - This function frees everything from the `PCGAMG` package. It is 2046 called from `PetscFinalize()` automatically. 2047 2048 Level: developer 2049 2050 .seealso: [](ch_ksp), `PetscFinalize()` 2051 @*/ 2052 PetscErrorCode PCGAMGFinalizePackage(void) 2053 { 2054 PetscFunctionBegin; 2055 PCGAMGPackageInitialized = PETSC_FALSE; 2056 PetscCall(PetscFunctionListDestroy(&GAMGList)); 2057 PetscFunctionReturn(PETSC_SUCCESS); 2058 } 2059 2060 /*@C 2061 PCGAMGRegister - Register a `PCGAMG` implementation. 2062 2063 Input Parameters: 2064 + type - string that will be used as the name of the `PCGAMG` type. 2065 - create - function for creating the gamg context. 2066 2067 Level: developer 2068 2069 .seealso: [](ch_ksp), `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()` 2070 @*/ 2071 PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 2072 { 2073 PetscFunctionBegin; 2074 PetscCall(PCGAMGInitializePackage()); 2075 PetscCall(PetscFunctionListAdd(&GAMGList, type, create)); 2076 PetscFunctionReturn(PETSC_SUCCESS); 2077 } 2078 2079 /*@ 2080 PCGAMGCreateGraph - Creates a graph that is used by the ``PCGAMGType`` in the coarsening process 2081 2082 Input Parameters: 2083 + pc - the `PCGAMG` 2084 - A - the matrix, for any level 2085 2086 Output Parameter: 2087 . G - the graph 2088 2089 Level: advanced 2090 2091 .seealso: [](ch_ksp), `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()` 2092 @*/ 2093 PetscErrorCode PCGAMGCreateGraph(PC pc, Mat A, Mat *G) 2094 { 2095 PC_MG *mg = (PC_MG *)pc->data; 2096 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 2097 2098 PetscFunctionBegin; 2099 PetscCall(pc_gamg->ops->creategraph(pc, A, G)); 2100 PetscFunctionReturn(PETSC_SUCCESS); 2101 } 2102