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