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