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 /* ----------------------------------------------------------------------------- */ 28 PetscErrorCode PCReset_GAMG(PC pc) { 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 /* 47 PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 48 of active processors. 49 50 Input Parameter: 51 . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 52 'pc_gamg->data_sz' are changed via repartitioning/reduction. 53 . Amat_fine - matrix on this fine (k) level 54 . cr_bs - coarse block size 55 In/Output Parameter: 56 . a_P_inout - prolongation operator to the next level (k-->k-1) 57 . a_nactive_proc - number of active procs 58 Output Parameter: 59 . a_Amat_crs - coarse matrix that is created (k-1) 60 */ 61 static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc, Mat Amat_fine, PetscInt cr_bs, Mat *a_P_inout, Mat *a_Amat_crs, PetscMPIInt *a_nactive_proc, IS *Pcolumnperm, PetscBool is_last) { 62 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 const char *prefix; 468 char addp[32]; 469 PC_MG *mg = (PC_MG *)a_pc->data; 470 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 471 472 PetscFunctionBegin; 473 PetscCall(PCGetOptionsPrefix(a_pc, &prefix)); 474 PetscCall(PetscInfo(a_pc, "%s: Square Graph on level %" PetscInt_FMT "\n", ((PetscObject)a_pc)->prefix, pc_gamg->current_level + 1)); 475 PetscCall(MatProductCreate(Gmat1, Gmat1, NULL, Gmat2)); 476 PetscCall(MatSetOptionsPrefix(*Gmat2, prefix)); 477 PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_square_%" PetscInt_FMT "_", pc_gamg->current_level)); 478 PetscCall(MatAppendOptionsPrefix(*Gmat2, addp)); 479 if ((*Gmat2)->structurally_symmetric == PETSC_BOOL3_TRUE) { 480 PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AB)); 481 } else { 482 PetscCall(MatSetOption(Gmat1, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 483 PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AtB)); 484 } 485 PetscCall(MatProductSetFromOptions(*Gmat2)); 486 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0)); 487 PetscCall(MatProductSymbolic(*Gmat2)); 488 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0)); 489 PetscCall(MatProductClear(*Gmat2)); 490 /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */ 491 (*Gmat2)->assembled = PETSC_TRUE; 492 PetscFunctionReturn(0); 493 } 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 PC_MG *mg = (PC_MG *)pc->data; 506 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 507 Mat Pmat = pc->pmat; 508 PetscInt fine_level, level, level1, bs, M, N, qq, lidx, nASMBlocksArr[PETSC_MG_MAXLEVELS]; 509 MPI_Comm comm; 510 PetscMPIInt rank, size, nactivepe; 511 Mat Aarr[PETSC_MG_MAXLEVELS], Parr[PETSC_MG_MAXLEVELS]; 512 IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 513 PetscLogDouble nnz0 = 0., nnztot = 0.; 514 MatInfo info; 515 PetscBool is_last = PETSC_FALSE; 516 517 PetscFunctionBegin; 518 PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 519 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 520 PetscCallMPI(MPI_Comm_size(comm, &size)); 521 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 522 if (pc->setupcalled) { 523 if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) { 524 /* reset everything */ 525 PetscCall(PCReset_MG(pc)); 526 pc->setupcalled = 0; 527 } else { 528 PC_MG_Levels **mglevels = mg->levels; 529 /* just do Galerkin grids */ 530 Mat B, dA, dB; 531 if (pc_gamg->Nlevels > 1) { 532 PetscInt gl; 533 /* currently only handle case where mat and pmat are the same on coarser levels */ 534 PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, &dA, &dB)); 535 /* (re)set to get dirty flag */ 536 PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, dA, dB)); 537 538 for (level = pc_gamg->Nlevels - 2, gl = 0; level >= 0; level--, gl++) { 539 MatReuse reuse = MAT_INITIAL_MATRIX; 540 #if defined(GAMG_STAGES) 541 PetscCall(PetscLogStagePush(gamg_stages[gl])); 542 #endif 543 /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 544 PetscCall(KSPGetOperators(mglevels[level]->smoothd, NULL, &B)); 545 if (B->product) { 546 if (B->product->A == dB && B->product->B == mglevels[level + 1]->interpolate) { reuse = MAT_REUSE_MATRIX; } 547 } 548 if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A)); 549 if (reuse == MAT_REUSE_MATRIX) { 550 PetscCall(PetscInfo(pc, "%s: RAP after first solve, reuse matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 551 } else { 552 PetscCall(PetscInfo(pc, "%s: RAP after first solve, new matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 553 } 554 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0)); 555 PetscCall(MatPtAP(dB, mglevels[level + 1]->interpolate, reuse, PETSC_DEFAULT, &B)); 556 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0)); 557 if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B; 558 PetscCall(KSPSetOperators(mglevels[level]->smoothd, B, B)); 559 dB = B; 560 #if defined(GAMG_STAGES) 561 PetscCall(PetscLogStagePop()); 562 #endif 563 } 564 } 565 566 PetscCall(PCSetUp_MG(pc)); 567 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 568 PetscFunctionReturn(0); 569 } 570 } 571 572 if (!pc_gamg->data) { 573 if (pc_gamg->orig_data) { 574 PetscCall(MatGetBlockSize(Pmat, &bs)); 575 PetscCall(MatGetLocalSize(Pmat, &qq, NULL)); 576 577 pc_gamg->data_sz = (qq / bs) * pc_gamg->orig_data_cell_rows * pc_gamg->orig_data_cell_cols; 578 pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 579 pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 580 581 PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data)); 582 for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 583 } else { 584 PetscCheck(pc_gamg->ops->createdefaultdata, comm, PETSC_ERR_PLIB, "'createdefaultdata' not set(?) need to support NULL data"); 585 PetscCall(pc_gamg->ops->createdefaultdata(pc, Pmat)); 586 } 587 } 588 589 /* cache original data for reuse */ 590 if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 591 PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data)); 592 for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 593 pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 594 pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 595 } 596 597 /* get basic dims */ 598 PetscCall(MatGetBlockSize(Pmat, &bs)); 599 PetscCall(MatGetSize(Pmat, &M, &N)); 600 601 PetscCall(MatGetInfo(Pmat, MAT_GLOBAL_SUM, &info)); /* global reduction */ 602 nnz0 = info.nz_used; 603 nnztot = info.nz_used; 604 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)); 605 606 /* Get A_i and R_i */ 607 for (level = 0, Aarr[0] = Pmat, nactivepe = size; level < (pc_gamg->Nlevels - 1) && (!level || M > pc_gamg->coarse_eq_limit); level++) { 608 pc_gamg->current_level = level; 609 PetscCheck(level < PETSC_MG_MAXLEVELS, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Too many levels %" PetscInt_FMT, level); 610 level1 = level + 1; 611 #if defined(GAMG_STAGES) 612 if (!gamg_stages[level]) { 613 char str[32]; 614 sprintf(str, "GAMG Level %d", (int)level); 615 PetscCall(PetscLogStageRegister(str, &gamg_stages[level])); 616 } 617 PetscCall(PetscLogStagePush(gamg_stages[level])); 618 #endif 619 { /* construct prolongator */ 620 Mat Gmat; 621 PetscCoarsenData *agg_lists; 622 Mat Prol11; 623 624 PetscCall(pc_gamg->ops->graph(pc, Aarr[level], &Gmat)); 625 PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists)); 626 PetscCall(pc_gamg->ops->prolongator(pc, Aarr[level], Gmat, agg_lists, &Prol11)); 627 628 /* could have failed to create new level */ 629 if (Prol11) { 630 const char *prefix; 631 char addp[32]; 632 633 /* get new block size of coarse matrices */ 634 PetscCall(MatGetBlockSizes(Prol11, NULL, &bs)); 635 636 if (pc_gamg->ops->optprolongator) { 637 /* smooth */ 638 PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11)); 639 } 640 641 if (pc_gamg->use_aggs_in_asm) { 642 PetscInt bs; 643 PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // not timed directly, ugly, could remove, but good ASM method 644 PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level])); 645 } 646 647 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 648 PetscCall(MatSetOptionsPrefix(Prol11, prefix)); 649 PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_prolongator_%d_", (int)level)); 650 PetscCall(MatAppendOptionsPrefix(Prol11, addp)); 651 /* Always generate the transpose with CUDA 652 Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */ 653 PetscCall(MatSetOption(Prol11, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 654 PetscCall(MatSetFromOptions(Prol11)); 655 Parr[level1] = Prol11; 656 } else Parr[level1] = NULL; /* failed to coarsen */ 657 658 PetscCall(MatDestroy(&Gmat)); 659 PetscCall(PetscCDDestroy(agg_lists)); 660 } /* construct prolongator scope */ 661 if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 662 if (!Parr[level1]) { /* failed to coarsen */ 663 PetscCall(PetscInfo(pc, "%s: Stop gridding, level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level)); 664 #if defined(GAMG_STAGES) 665 PetscCall(PetscLogStagePop()); 666 #endif 667 break; 668 } 669 PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */ 670 PetscCheck(!is_last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Is last ?"); 671 if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 672 if (level1 == pc_gamg->Nlevels - 1) is_last = PETSC_TRUE; 673 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0)); 674 PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last)); 675 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0)); 676 677 PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */ 678 PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info)); 679 nnztot += info.nz_used; 680 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)); 681 682 #if defined(GAMG_STAGES) 683 PetscCall(PetscLogStagePop()); 684 #endif 685 /* stop if one node or one proc -- could pull back for singular problems */ 686 if ((pc_gamg->data_cell_cols && M / pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M / bs < 2)) { 687 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)); 688 level++; 689 break; 690 } 691 } /* levels */ 692 PetscCall(PetscFree(pc_gamg->data)); 693 694 PetscCall(PetscInfo(pc, "%s: %" PetscInt_FMT " levels, grid complexity = %g\n", ((PetscObject)pc)->prefix, level + 1, nnztot / nnz0)); 695 pc_gamg->Nlevels = level + 1; 696 fine_level = level; 697 PetscCall(PCMGSetLevels(pc, pc_gamg->Nlevels, NULL)); 698 699 if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 700 701 /* set default smoothers & set operators */ 702 for (lidx = 1, level = pc_gamg->Nlevels - 2; lidx <= fine_level; lidx++, level--) { 703 KSP smoother; 704 PC subpc; 705 706 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 707 PetscCall(KSPGetPC(smoother, &subpc)); 708 709 PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 710 /* set ops */ 711 PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level])); 712 PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level + 1])); 713 714 /* set defaults */ 715 PetscCall(KSPSetType(smoother, KSPCHEBYSHEV)); 716 717 /* set blocks for ASM smoother that uses the 'aggregates' */ 718 if (pc_gamg->use_aggs_in_asm) { 719 PetscInt sz; 720 IS *iss; 721 722 sz = nASMBlocksArr[level]; 723 iss = ASMLocalIDsArr[level]; 724 PetscCall(PCSetType(subpc, PCASM)); 725 PetscCall(PCASMSetOverlap(subpc, 0)); 726 PetscCall(PCASMSetType(subpc, PC_ASM_BASIC)); 727 if (!sz) { 728 IS is; 729 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is)); 730 PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is)); 731 PetscCall(ISDestroy(&is)); 732 } else { 733 PetscInt kk; 734 PetscCall(PCASMSetLocalSubdomains(subpc, sz, NULL, iss)); 735 for (kk = 0; kk < sz; kk++) { PetscCall(ISDestroy(&iss[kk])); } 736 PetscCall(PetscFree(iss)); 737 } 738 ASMLocalIDsArr[level] = NULL; 739 nASMBlocksArr[level] = 0; 740 } else { 741 PetscCall(PCSetType(subpc, PCJACOBI)); 742 } 743 } 744 { 745 /* coarse grid */ 746 KSP smoother, *k2; 747 PC subpc, pc2; 748 PetscInt ii, first; 749 Mat Lmat = Aarr[(level = pc_gamg->Nlevels - 1)]; 750 lidx = 0; 751 752 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 753 PetscCall(KSPSetOperators(smoother, Lmat, Lmat)); 754 if (!pc_gamg->use_parallel_coarse_grid_solver) { 755 PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 756 PetscCall(KSPGetPC(smoother, &subpc)); 757 PetscCall(PCSetType(subpc, PCBJACOBI)); 758 PetscCall(PCSetUp(subpc)); 759 PetscCall(PCBJacobiGetSubKSP(subpc, &ii, &first, &k2)); 760 PetscCheck(ii == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " is not one", ii); 761 PetscCall(KSPGetPC(k2[0], &pc2)); 762 PetscCall(PCSetType(pc2, PCLU)); 763 PetscCall(PCFactorSetShiftType(pc2, MAT_SHIFT_INBLOCKS)); 764 PetscCall(KSPSetTolerances(k2[0], PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1)); 765 PetscCall(KSPSetType(k2[0], KSPPREONLY)); 766 } 767 } 768 769 /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 770 PetscObjectOptionsBegin((PetscObject)pc); 771 PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); 772 PetscOptionsEnd(); 773 PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 774 775 /* setup cheby eigen estimates from SA */ 776 if (pc_gamg->use_sa_esteig) { 777 for (lidx = 1, level = pc_gamg->Nlevels - 2; level >= 0; lidx++, level--) { 778 KSP smoother; 779 PetscBool ischeb; 780 781 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 782 PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb)); 783 if (ischeb) { 784 KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data; 785 786 // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG 787 if (mg->max_eigen_DinvA[level] > 0) { 788 // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it. 789 // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.) 790 PetscReal emax, emin; 791 792 emin = mg->min_eigen_DinvA[level]; 793 emax = mg->max_eigen_DinvA[level]; 794 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)); 795 cheb->emin_provided = emin; 796 cheb->emax_provided = emax; 797 } 798 } 799 } 800 } 801 802 PetscCall(PCSetUp_MG(pc)); 803 804 /* clean up */ 805 for (level = 1; level < pc_gamg->Nlevels; level++) { 806 PetscCall(MatDestroy(&Parr[level])); 807 PetscCall(MatDestroy(&Aarr[level])); 808 } 809 } else { 810 KSP smoother; 811 812 PetscCall(PetscInfo(pc, "%s: One level solver used (system is seen as DD). Using default solver.\n", ((PetscObject)pc)->prefix)); 813 PetscCall(PCMGGetSmoother(pc, 0, &smoother)); 814 PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0])); 815 PetscCall(KSPSetType(smoother, KSPPREONLY)); 816 PetscCall(PCSetUp_MG(pc)); 817 } 818 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0)); 819 PetscFunctionReturn(0); 820 } 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 PC_MG *mg = (PC_MG *)pc->data; 834 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 835 836 PetscFunctionBegin; 837 PetscCall(PCReset_GAMG(pc)); 838 if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc)); 839 PetscCall(PetscFree(pc_gamg->ops)); 840 PetscCall(PetscFree(pc_gamg->gamg_type_name)); 841 PetscCall(PetscFree(pc_gamg)); 842 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", NULL)); 843 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", NULL)); 844 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", NULL)); 845 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", NULL)); 846 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", NULL)); 847 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", NULL)); 848 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", NULL)); 849 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", NULL)); 850 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", NULL)); 851 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", NULL)); 852 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", NULL)); 853 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", NULL)); 854 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", NULL)); 855 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", NULL)); 856 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", NULL)); 857 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", NULL)); 858 PetscCall(PCDestroy_MG(pc)); 859 PetscFunctionReturn(0); 860 } 861 862 /*@ 863 PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction. 864 865 Logically Collective on PC 866 867 Input Parameters: 868 + pc - the preconditioner context 869 - n - the number of equations 870 871 Options Database Key: 872 . -pc_gamg_process_eq_limit <limit> - set the limit 873 874 Notes: 875 GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process 876 that has degrees of freedom 877 878 Level: intermediate 879 880 .seealso: `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()` 881 @*/ 882 PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) { 883 PetscFunctionBegin; 884 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 885 PetscTryMethod(pc, "PCGAMGSetProcEqLim_C", (PC, PetscInt), (pc, n)); 886 PetscFunctionReturn(0); 887 } 888 889 static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) { 890 PC_MG *mg = (PC_MG *)pc->data; 891 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 892 893 PetscFunctionBegin; 894 if (n > 0) pc_gamg->min_eq_proc = n; 895 PetscFunctionReturn(0); 896 } 897 898 /*@ 899 PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid. 900 901 Collective on PC 902 903 Input Parameters: 904 + pc - the preconditioner context 905 - n - maximum number of equations to aim for 906 907 Options Database Key: 908 . -pc_gamg_coarse_eq_limit <limit> - set the limit 909 910 Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 911 has less than 1000 unknowns. 912 913 Level: intermediate 914 915 .seealso: `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()` 916 @*/ 917 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) { 918 PetscFunctionBegin; 919 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 920 PetscTryMethod(pc, "PCGAMGSetCoarseEqLim_C", (PC, PetscInt), (pc, n)); 921 PetscFunctionReturn(0); 922 } 923 924 static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) { 925 PC_MG *mg = (PC_MG *)pc->data; 926 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 927 928 PetscFunctionBegin; 929 if (n > 0) pc_gamg->coarse_eq_limit = n; 930 PetscFunctionReturn(0); 931 } 932 933 /*@ 934 PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids 935 936 Collective on PC 937 938 Input Parameters: 939 + pc - the preconditioner context 940 - n - PETSC_TRUE or PETSC_FALSE 941 942 Options Database Key: 943 . -pc_gamg_repartition <true,false> - turn on the repartitioning 944 945 Notes: 946 this will generally improve the loading balancing of the work on each level 947 948 Level: intermediate 949 950 @*/ 951 PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) { 952 PetscFunctionBegin; 953 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 954 PetscTryMethod(pc, "PCGAMGSetRepartition_C", (PC, PetscBool), (pc, n)); 955 PetscFunctionReturn(0); 956 } 957 958 static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) { 959 PC_MG *mg = (PC_MG *)pc->data; 960 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 961 962 PetscFunctionBegin; 963 pc_gamg->repart = n; 964 PetscFunctionReturn(0); 965 } 966 967 /*@ 968 PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Chebyshev smoother 969 970 Collective on PC 971 972 Input Parameters: 973 + pc - the preconditioner context 974 - n - number of its 975 976 Options Database Key: 977 . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate 978 979 Notes: 980 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$. 981 Eigenvalue estimates (based on a few CG or GMRES iterations) are computed to choose $\omega$ so that this is a stable smoothing operation. 982 If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused. 983 This option is only used when the smoother uses Jacobi, and should be turned off if a different PCJacobiType is used. 984 It became default in PETSc 3.17. 985 986 Level: advanced 987 988 .seealso: `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()` 989 @*/ 990 PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n) { 991 PetscFunctionBegin; 992 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 993 PetscTryMethod(pc, "PCGAMGSetUseSAEstEig_C", (PC, PetscBool), (pc, n)); 994 PetscFunctionReturn(0); 995 } 996 997 static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n) { 998 PC_MG *mg = (PC_MG *)pc->data; 999 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1000 1001 PetscFunctionBegin; 1002 pc_gamg->use_sa_esteig = n; 1003 PetscFunctionReturn(0); 1004 } 1005 1006 /*@ 1007 PCGAMGSetEigenvalues - Set eigenvalues 1008 1009 Collective on PC 1010 1011 Input Parameters: 1012 + pc - the preconditioner context 1013 - emax - max eigenvalue 1014 . emin - min eigenvalue 1015 1016 Options Database Key: 1017 . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues 1018 1019 Level: intermediate 1020 1021 .seealso: `PCGAMGSetUseSAEstEig()` 1022 @*/ 1023 PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax, PetscReal emin) { 1024 PetscFunctionBegin; 1025 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1026 PetscTryMethod(pc, "PCGAMGSetEigenvalues_C", (PC, PetscReal, PetscReal), (pc, emax, emin)); 1027 PetscFunctionReturn(0); 1028 } 1029 1030 static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc, PetscReal emax, PetscReal emin) { 1031 PC_MG *mg = (PC_MG *)pc->data; 1032 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1033 1034 PetscFunctionBegin; 1035 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); 1036 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); 1037 pc_gamg->emax = emax; 1038 pc_gamg->emin = emin; 1039 PetscFunctionReturn(0); 1040 } 1041 1042 /*@ 1043 PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner 1044 1045 Collective on PC 1046 1047 Input Parameters: 1048 + pc - the preconditioner context 1049 - n - PETSC_TRUE or PETSC_FALSE 1050 1051 Options Database Key: 1052 . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation 1053 1054 Level: intermediate 1055 1056 Notes: 1057 May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1058 rebuilding the preconditioner quicker. 1059 1060 @*/ 1061 PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) { 1062 PetscFunctionBegin; 1063 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1064 PetscTryMethod(pc, "PCGAMGSetReuseInterpolation_C", (PC, PetscBool), (pc, n)); 1065 PetscFunctionReturn(0); 1066 } 1067 1068 static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) { 1069 PC_MG *mg = (PC_MG *)pc->data; 1070 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1071 1072 PetscFunctionBegin; 1073 pc_gamg->reuse_prol = n; 1074 PetscFunctionReturn(0); 1075 } 1076 1077 /*@ 1078 PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner. 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 @*/ 1092 PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) { 1093 PetscFunctionBegin; 1094 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1095 PetscTryMethod(pc, "PCGAMGASMSetUseAggs_C", (PC, PetscBool), (pc, flg)); 1096 PetscFunctionReturn(0); 1097 } 1098 1099 static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) { 1100 PC_MG *mg = (PC_MG *)pc->data; 1101 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1102 1103 PetscFunctionBegin; 1104 pc_gamg->use_aggs_in_asm = flg; 1105 PetscFunctionReturn(0); 1106 } 1107 1108 /*@ 1109 PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver 1110 1111 Collective on PC 1112 1113 Input Parameters: 1114 + pc - the preconditioner context 1115 - flg - PETSC_TRUE to not force coarse grid onto one processor 1116 1117 Options Database Key: 1118 . -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver 1119 1120 Level: intermediate 1121 1122 .seealso: `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()` 1123 @*/ 1124 PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg) { 1125 PetscFunctionBegin; 1126 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1127 PetscTryMethod(pc, "PCGAMGSetUseParallelCoarseGridSolve_C", (PC, PetscBool), (pc, flg)); 1128 PetscFunctionReturn(0); 1129 } 1130 1131 static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) { 1132 PC_MG *mg = (PC_MG *)pc->data; 1133 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1134 1135 PetscFunctionBegin; 1136 pc_gamg->use_parallel_coarse_grid_solver = flg; 1137 PetscFunctionReturn(0); 1138 } 1139 1140 /*@ 1141 PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU 1142 1143 Collective on PC 1144 1145 Input Parameters: 1146 + pc - the preconditioner context 1147 - flg - PETSC_TRUE to pin coarse grids to CPU 1148 1149 Options Database Key: 1150 . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU 1151 1152 Level: intermediate 1153 1154 .seealso: `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetUseParallelCoarseGridSolve()` 1155 @*/ 1156 PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) { 1157 PetscFunctionBegin; 1158 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1159 PetscTryMethod(pc, "PCGAMGSetCpuPinCoarseGrids_C", (PC, PetscBool), (pc, flg)); 1160 PetscFunctionReturn(0); 1161 } 1162 1163 static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) { 1164 PC_MG *mg = (PC_MG *)pc->data; 1165 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1166 1167 PetscFunctionBegin; 1168 pc_gamg->cpu_pin_coarse_grids = flg; 1169 PetscFunctionReturn(0); 1170 } 1171 1172 /*@ 1173 PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type) 1174 1175 Collective on PC 1176 1177 Input Parameters: 1178 + pc - the preconditioner context 1179 - flg - Layout type 1180 1181 Options Database Key: 1182 . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering 1183 1184 Level: intermediate 1185 1186 .seealso: `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()` 1187 @*/ 1188 PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) { 1189 PetscFunctionBegin; 1190 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1191 PetscTryMethod(pc, "PCGAMGSetCoarseGridLayoutType_C", (PC, PCGAMGLayoutType), (pc, flg)); 1192 PetscFunctionReturn(0); 1193 } 1194 1195 static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) { 1196 PC_MG *mg = (PC_MG *)pc->data; 1197 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1198 1199 PetscFunctionBegin; 1200 pc_gamg->layout_type = flg; 1201 PetscFunctionReturn(0); 1202 } 1203 1204 /*@ 1205 PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use 1206 1207 Not collective on PC 1208 1209 Input Parameters: 1210 + pc - the preconditioner 1211 - n - the maximum number of levels to use 1212 1213 Options Database Key: 1214 . -pc_mg_levels <n> - set the maximum number of levels to allow 1215 1216 Level: intermediate 1217 1218 @*/ 1219 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) { 1220 PetscFunctionBegin; 1221 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1222 PetscTryMethod(pc, "PCGAMGSetNlevels_C", (PC, PetscInt), (pc, n)); 1223 PetscFunctionReturn(0); 1224 } 1225 1226 static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) { 1227 PC_MG *mg = (PC_MG *)pc->data; 1228 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1229 1230 PetscFunctionBegin; 1231 pc_gamg->Nlevels = n; 1232 PetscFunctionReturn(0); 1233 } 1234 1235 /*@ 1236 PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 1237 1238 Not collective on PC 1239 1240 Input Parameters: 1241 + pc - the preconditioner context 1242 . 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 1243 - n - number of threshold values provided in array 1244 1245 Options Database Key: 1246 . -pc_gamg_threshold <threshold> - the threshold to drop edges 1247 1248 Notes: 1249 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. 1250 Before coarsening or aggregating the graph, GAMG 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. 1251 1252 If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening. 1253 In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold(). 1254 If n is greater than the total number of levels, the excess entries in threshold will not be used. 1255 1256 Level: intermediate 1257 1258 .seealso: `PCGAMGFilterGraph()`, `PCGAMGSetAggressiveLevels()` 1259 @*/ 1260 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) { 1261 PetscFunctionBegin; 1262 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1263 if (n) PetscValidRealPointer(v, 2); 1264 PetscTryMethod(pc, "PCGAMGSetThreshold_C", (PC, PetscReal[], PetscInt), (pc, v, n)); 1265 PetscFunctionReturn(0); 1266 } 1267 1268 static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) { 1269 PC_MG *mg = (PC_MG *)pc->data; 1270 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1271 PetscInt i; 1272 PetscFunctionBegin; 1273 for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i]; 1274 for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale; 1275 PetscFunctionReturn(0); 1276 } 1277 1278 /*@ 1279 PCGAMGSetRankReductionFactors - Set manual schedule for process reduction on coarse grids 1280 1281 Collective on PC 1282 1283 Input Parameters: 1284 + pc - the preconditioner context 1285 . v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda 1286 - n - number of values provided in array 1287 1288 Options Database Key: 1289 . -pc_gamg_rank_reduction_factors <factors> - provide the schedule 1290 1291 Level: intermediate 1292 1293 .seealso: `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()` 1294 @*/ 1295 PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n) { 1296 PetscFunctionBegin; 1297 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1298 if (n) PetscValidIntPointer(v, 2); 1299 PetscTryMethod(pc, "PCGAMGSetRankReductionFactors_C", (PC, PetscInt[], PetscInt), (pc, v, n)); 1300 PetscFunctionReturn(0); 1301 } 1302 1303 static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n) { 1304 PC_MG *mg = (PC_MG *)pc->data; 1305 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1306 PetscInt i; 1307 PetscFunctionBegin; 1308 for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i]; 1309 for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */ 1310 PetscFunctionReturn(0); 1311 } 1312 1313 /*@ 1314 PCGAMGSetThresholdScale - Relative threshold reduction at each level 1315 1316 Not collective on PC 1317 1318 Input Parameters: 1319 + pc - the preconditioner context 1320 - scale - the threshold value reduction, usually < 1.0 1321 1322 Options Database Key: 1323 . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level 1324 1325 Notes: 1326 The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold(). 1327 This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold(). 1328 1329 Level: advanced 1330 1331 .seealso: `PCGAMGSetThreshold()` 1332 @*/ 1333 PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) { 1334 PetscFunctionBegin; 1335 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1336 PetscTryMethod(pc, "PCGAMGSetThresholdScale_C", (PC, PetscReal), (pc, v)); 1337 PetscFunctionReturn(0); 1338 } 1339 1340 static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) { 1341 PC_MG *mg = (PC_MG *)pc->data; 1342 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1343 PetscFunctionBegin; 1344 pc_gamg->threshold_scale = v; 1345 PetscFunctionReturn(0); 1346 } 1347 1348 /*@C 1349 PCGAMGSetType - Set solution method 1350 1351 Collective on PC 1352 1353 Input Parameters: 1354 + pc - the preconditioner context 1355 - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1356 1357 Options Database Key: 1358 . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply 1359 1360 Level: intermediate 1361 1362 .seealso: `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType` 1363 @*/ 1364 PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) { 1365 PetscFunctionBegin; 1366 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1367 PetscTryMethod(pc, "PCGAMGSetType_C", (PC, PCGAMGType), (pc, type)); 1368 PetscFunctionReturn(0); 1369 } 1370 1371 /*@C 1372 PCGAMGGetType - Get solution method 1373 1374 Collective on PC 1375 1376 Input Parameter: 1377 . pc - the preconditioner context 1378 1379 Output Parameter: 1380 . type - the type of algorithm used 1381 1382 Level: intermediate 1383 1384 .seealso: `PCGAMGSetType()`, `PCGAMGType` 1385 @*/ 1386 PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) { 1387 PetscFunctionBegin; 1388 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1389 PetscUseMethod(pc, "PCGAMGGetType_C", (PC, PCGAMGType *), (pc, type)); 1390 PetscFunctionReturn(0); 1391 } 1392 1393 static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) { 1394 PC_MG *mg = (PC_MG *)pc->data; 1395 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1396 1397 PetscFunctionBegin; 1398 *type = pc_gamg->type; 1399 PetscFunctionReturn(0); 1400 } 1401 1402 static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) { 1403 PC_MG *mg = (PC_MG *)pc->data; 1404 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1405 PetscErrorCode (*r)(PC); 1406 1407 PetscFunctionBegin; 1408 pc_gamg->type = type; 1409 PetscCall(PetscFunctionListFind(GAMGList, type, &r)); 1410 PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown GAMG type %s given", type); 1411 if (pc_gamg->ops->destroy) { 1412 PetscCall((*pc_gamg->ops->destroy)(pc)); 1413 PetscCall(PetscMemzero(pc_gamg->ops, sizeof(struct _PCGAMGOps))); 1414 pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1415 /* cleaning up common data in pc_gamg - this should disapear someday */ 1416 pc_gamg->data_cell_cols = 0; 1417 pc_gamg->data_cell_rows = 0; 1418 pc_gamg->orig_data_cell_cols = 0; 1419 pc_gamg->orig_data_cell_rows = 0; 1420 PetscCall(PetscFree(pc_gamg->data)); 1421 pc_gamg->data_sz = 0; 1422 } 1423 PetscCall(PetscFree(pc_gamg->gamg_type_name)); 1424 PetscCall(PetscStrallocpy(type, &pc_gamg->gamg_type_name)); 1425 PetscCall((*r)(pc)); 1426 PetscFunctionReturn(0); 1427 } 1428 1429 static PetscErrorCode PCView_GAMG(PC pc, PetscViewer viewer) { 1430 PC_MG *mg = (PC_MG *)pc->data; 1431 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1432 PetscReal gc = 0, oc = 0; 1433 1434 PetscFunctionBegin; 1435 PetscCall(PetscViewerASCIIPrintf(viewer, " GAMG specific options\n")); 1436 PetscCall(PetscViewerASCIIPrintf(viewer, " Threshold for dropping small values in graph on each level =")); 1437 for (PetscInt i = 0; i < mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)pc_gamg->threshold[i])); 1438 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 1439 PetscCall(PetscViewerASCIIPrintf(viewer, " Threshold scaling factor for each level not specified = %g\n", (double)pc_gamg->threshold_scale)); 1440 if (pc_gamg->use_aggs_in_asm) { PetscCall(PetscViewerASCIIPrintf(viewer, " Using aggregates from coarsening process to define subdomains for PCASM\n")); } 1441 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")); } 1442 if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc, viewer)); 1443 PetscCall(PCMGGetGridComplexity(pc, &gc, &oc)); 1444 PetscCall(PetscViewerASCIIPrintf(viewer, " Complexity: grid = %g operator = %g\n", (double)gc, (double)oc)); 1445 PetscFunctionReturn(0); 1446 } 1447 1448 PetscErrorCode PCSetFromOptions_GAMG(PC pc, PetscOptionItems *PetscOptionsObject) { 1449 PC_MG *mg = (PC_MG *)pc->data; 1450 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1451 PetscBool flag; 1452 MPI_Comm comm; 1453 char prefix[256], tname[32]; 1454 PetscInt i, n; 1455 const char *pcpre; 1456 static const char *LayoutTypes[] = {"compact", "spread", "PCGAMGLayoutType", "PC_GAMG_LAYOUT", NULL}; 1457 PetscFunctionBegin; 1458 PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 1459 PetscOptionsHeadBegin(PetscOptionsObject, "GAMG options"); 1460 PetscCall(PetscOptionsFList("-pc_gamg_type", "Type of AMG method", "PCGAMGSetType", GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag)); 1461 if (flag) PetscCall(PCGAMGSetType(pc, tname)); 1462 PetscCall(PetscOptionsBool("-pc_gamg_repartition", "Repartion coarse grids", "PCGAMGSetRepartition", pc_gamg->repart, &pc_gamg->repart, NULL)); 1463 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)); 1464 PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation", "Reuse prolongation operator", "PCGAMGReuseInterpolation", pc_gamg->reuse_prol, &pc_gamg->reuse_prol, NULL)); 1465 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)); 1466 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)); 1467 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)); 1468 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, 1469 (PetscEnum)pc_gamg->layout_type, (PetscEnum *)&pc_gamg->layout_type, NULL)); 1470 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)); 1471 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)); 1472 PetscCall(PetscOptionsReal("-pc_gamg_threshold_scale", "Scaling of threshold for each level not specified", "PCGAMGSetThresholdScale", pc_gamg->threshold_scale, &pc_gamg->threshold_scale, NULL)); 1473 n = PETSC_MG_MAXLEVELS; 1474 PetscCall(PetscOptionsRealArray("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", pc_gamg->threshold, &n, &flag)); 1475 if (!flag || n < PETSC_MG_MAXLEVELS) { 1476 if (!flag) n = 1; 1477 i = n; 1478 do { pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale; } while (++i < PETSC_MG_MAXLEVELS); 1479 } 1480 n = PETSC_MG_MAXLEVELS; 1481 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)); 1482 if (!flag) i = 0; 1483 else i = n; 1484 do { pc_gamg->level_reduction_factors[i] = -1; } while (++i < PETSC_MG_MAXLEVELS); 1485 PetscCall(PetscOptionsInt("-pc_mg_levels", "Set number of MG levels", "PCGAMGSetNlevels", pc_gamg->Nlevels, &pc_gamg->Nlevels, NULL)); 1486 { 1487 PetscReal eminmax[2] = {0., 0.}; 1488 n = 2; 1489 PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues", "extreme eigenvalues for smoothed aggregation", "PCGAMGSetEigenvalues", eminmax, &n, &flag)); 1490 if (flag) { 1491 PetscCheck(n == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 1492 PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0])); 1493 } 1494 } 1495 /* set options for subtype */ 1496 PetscCall((*pc_gamg->ops->setfromoptions)(pc, PetscOptionsObject)); 1497 1498 PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 1499 PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : "")); 1500 PetscOptionsHeadEnd(); 1501 PetscFunctionReturn(0); 1502 } 1503 1504 /* -------------------------------------------------------------------------- */ 1505 /*MC 1506 PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 1507 1508 Options Database Keys: 1509 + -pc_gamg_type <type> - one of agg, geo, or classical 1510 . -pc_gamg_repartition <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined 1511 - -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations 1512 + -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 1513 . -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> 1514 equations on each process that has degrees of freedom 1515 - -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 1516 + -pc_gamg_threshold[] <thresh,default=-1> - Before aggregating the graph GAMG will remove small values from the graph on each level (< 0 is no filtering) 1517 . -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1518 1519 Options Database Keys for default Aggregation: 1520 + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 1521 . -pc_gamg_symmetrize_graph <true,default=false> - symmetrize the graph before computing the aggregation 1522 - -pc_gamg_square_graph <n,default=1> - alias for pc_gamg_aggressive_coarsening (deprecated) 1523 - -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest. 1524 1525 Multigrid options: 1526 + -pc_mg_cycles <v> - v or w, see PCMGSetCycleType() 1527 . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1528 . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1529 - -pc_mg_levels <levels> - Number of levels of multigrid to use. 1530 1531 Notes: 1532 In order to obtain good performance for PCGAMG for vector valued problems you must 1533 Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point 1534 Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator 1535 See the Users Manual Chapter 4 for more details 1536 1537 Level: intermediate 1538 1539 .seealso: `PCCreate()`, `PCSetType()`, `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, 1540 `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGSetUseSAEstEig()` 1541 M*/ 1542 1543 PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) { 1544 PC_GAMG *pc_gamg; 1545 PC_MG *mg; 1546 1547 PetscFunctionBegin; 1548 /* register AMG type */ 1549 PetscCall(PCGAMGInitializePackage()); 1550 1551 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 1552 PetscCall(PCSetType(pc, PCMG)); 1553 PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG)); 1554 1555 /* create a supporting struct and attach it to pc */ 1556 PetscCall(PetscNewLog(pc, &pc_gamg)); 1557 PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL)); 1558 mg = (PC_MG *)pc->data; 1559 mg->innerctx = pc_gamg; 1560 1561 PetscCall(PetscNewLog(pc, &pc_gamg->ops)); 1562 1563 /* these should be in subctx but repartitioning needs simple arrays */ 1564 pc_gamg->data_sz = 0; 1565 pc_gamg->data = NULL; 1566 1567 /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 1568 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 1569 pc->ops->setup = PCSetUp_GAMG; 1570 pc->ops->reset = PCReset_GAMG; 1571 pc->ops->destroy = PCDestroy_GAMG; 1572 mg->view = PCView_GAMG; 1573 1574 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", PCGAMGSetProcEqLim_GAMG)); 1575 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", PCGAMGSetCoarseEqLim_GAMG)); 1576 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", PCGAMGSetRepartition_GAMG)); 1577 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", PCGAMGSetEigenvalues_GAMG)); 1578 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", PCGAMGSetUseSAEstEig_GAMG)); 1579 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", PCGAMGSetReuseInterpolation_GAMG)); 1580 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", PCGAMGASMSetUseAggs_GAMG)); 1581 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", PCGAMGSetUseParallelCoarseGridSolve_GAMG)); 1582 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", PCGAMGSetCpuPinCoarseGrids_GAMG)); 1583 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", PCGAMGSetCoarseGridLayoutType_GAMG)); 1584 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", PCGAMGSetThreshold_GAMG)); 1585 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", PCGAMGSetRankReductionFactors_GAMG)); 1586 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", PCGAMGSetThresholdScale_GAMG)); 1587 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", PCGAMGSetType_GAMG)); 1588 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", PCGAMGGetType_GAMG)); 1589 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", PCGAMGSetNlevels_GAMG)); 1590 pc_gamg->repart = PETSC_FALSE; 1591 pc_gamg->reuse_prol = PETSC_FALSE; 1592 pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1593 pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1594 pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1595 pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1596 pc_gamg->min_eq_proc = 50; 1597 pc_gamg->coarse_eq_limit = 50; 1598 for (int i = 0; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = -1; 1599 pc_gamg->threshold_scale = 1.; 1600 pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 1601 pc_gamg->current_level = 0; /* don't need to init really */ 1602 pc_gamg->use_sa_esteig = PETSC_TRUE; 1603 pc_gamg->emin = 0; 1604 pc_gamg->emax = 0; 1605 1606 pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1607 1608 /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1609 PetscCall(PCGAMGSetType(pc, PCGAMGAGG)); 1610 PetscFunctionReturn(0); 1611 } 1612 1613 /*@C 1614 PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 1615 from PCInitializePackage(). 1616 1617 Level: developer 1618 1619 .seealso: `PetscInitialize()` 1620 @*/ 1621 PetscErrorCode PCGAMGInitializePackage(void) { 1622 PetscInt l; 1623 1624 PetscFunctionBegin; 1625 if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 1626 PCGAMGPackageInitialized = PETSC_TRUE; 1627 PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGGEO, PCCreateGAMG_GEO)); 1628 PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGAGG, PCCreateGAMG_AGG)); 1629 PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGCLASSICAL, PCCreateGAMG_Classical)); 1630 PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage)); 1631 1632 /* general events */ 1633 PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP])); 1634 PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH])); 1635 PetscCall(PetscLogEventRegister(" PCGAMGFilter", PC_CLASSID, &petsc_gamg_setup_events[GAMG_FILTER])); 1636 PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN])); 1637 PetscCall(PetscLogEventRegister(" GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS])); 1638 PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL])); 1639 PetscCall(PetscLogEventRegister(" GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA])); 1640 PetscCall(PetscLogEventRegister(" GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB])); 1641 PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT])); 1642 PetscCall(PetscLogEventRegister(" GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM])); 1643 PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL])); 1644 PetscCall(PetscLogEventRegister(" GAMG PtAP", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP])); 1645 PetscCall(PetscLogEventRegister(" GAMG Reduce", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE])); 1646 PetscCall(PetscLogEventRegister(" GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART])); 1647 /* PetscCall(PetscLogEventRegister(" GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */ 1648 /* PetscCall(PetscLogEventRegister(" GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */ 1649 /* PetscCall(PetscLogEventRegister(" GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */ 1650 for (l = 0; l < PETSC_MG_MAXLEVELS; l++) { 1651 char ename[32]; 1652 1653 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Squ l%02" PetscInt_FMT, l)); 1654 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0])); 1655 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Gal l%02" PetscInt_FMT, l)); 1656 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1])); 1657 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Opt l%02" PetscInt_FMT, l)); 1658 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2])); 1659 } 1660 #if defined(GAMG_STAGES) 1661 { /* create timer stages */ 1662 char str[32]; 1663 sprintf(str, "GAMG Level %d", 0); 1664 PetscCall(PetscLogStageRegister(str, &gamg_stages[0])); 1665 } 1666 #endif 1667 PetscFunctionReturn(0); 1668 } 1669 1670 /*@C 1671 PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is 1672 called from PetscFinalize() automatically. 1673 1674 Level: developer 1675 1676 .seealso: `PetscFinalize()` 1677 @*/ 1678 PetscErrorCode PCGAMGFinalizePackage(void) { 1679 PetscFunctionBegin; 1680 PCGAMGPackageInitialized = PETSC_FALSE; 1681 PetscCall(PetscFunctionListDestroy(&GAMGList)); 1682 PetscFunctionReturn(0); 1683 } 1684 1685 /*@C 1686 PCGAMGRegister - Register a PCGAMG implementation. 1687 1688 Input Parameters: 1689 + type - string that will be used as the name of the GAMG type. 1690 - create - function for creating the gamg context. 1691 1692 Level: advanced 1693 1694 .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()` 1695 @*/ 1696 PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) { 1697 PetscFunctionBegin; 1698 PetscCall(PCGAMGInitializePackage()); 1699 PetscCall(PetscFunctionListAdd(&GAMGList, type, create)); 1700 PetscFunctionReturn(0); 1701 } 1702