:orphan: # PCMGSetLevels Sets the number of levels to use with `PCMG`. Must be called before any other `PCMG` routine. ## Synopsis ``` #include "petscksp.h" PetscErrorCode PCMGSetLevels(PC pc, PetscInt levels, MPI_Comm *comms) ``` Logically Collective ## Input Parameters - ***pc -*** the preconditioner context - ***levels -*** the number of levels - ***comms -*** optional communicators for each level; this is to allow solving the coarser problems on smaller sets of processes. For processes that are not included in the computation you must pass `MPI_COMM_NULL`. Use comms = `NULL` to specify that all processes should participate in each level of problem. ## Notes If the number of levels is one then the multigrid uses the `-mg_levels` prefix for setting the level options rather than the `-mg_coarse` prefix. You can free the information in comms after this routine is called. The array of MPI communicators must contain `MPI_COMM_NULL` for those ranks that at each level are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on the two levels, rank 0 in the original communicator will pass in an array of 2 communicators of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators the first of size 2 and the second of value `MPI_COMM_NULL` since the rank 1 does not participate in the coarse grid solve. Since each coarser level may have a new `MPI_Comm` with fewer ranks than the previous, one must take special care in providing the restriction and interpolation operation. We recommend providing these as two step operations; first perform a standard restriction or interpolation on the full number of ranks for that level and then use an MPI call to copy the resulting vector array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, note in both cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and receives or `MPI_AlltoAllv()` could be used to do the reshuffling of the vector entries. ## Fortran Note Use comms = `PETSC_NULL_MPI_COMM` as the equivalent of `NULL` in the C interface. Note `PETSC_NULL_MPI_COMM` is not `MPI_COMM_NULL`. It is more like `PETSC_NULL_INTEGER`, `PETSC_NULL_REAL` etc. ## See Also `PCMGSetType()`, `PCMGGetLevels()` ## Level intermediate ## Location src/ksp/pc/impls/mg/mg.c ## Examples src/dm/impls/stag/tutorials/ex4.c
src/ksp/ksp/tutorials/ex35.cxx
src/ksp/ksp/tutorials/ex36.cxx
src/ksp/ksp/tutorials/ex42.c
## Implementations PCMGSetLevels_MG in src/ksp/pc/impls/mg/mg.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/pc/impls/mg/mg.c) [Index of all PC routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)