:orphan: # PCASMSetTotalSubdomains Sets the subdomains for all processors for the additive Schwarz preconditioner, `PCASM`. ## Synopsis ``` #include "petscpc.h" PetscErrorCode PCASMSetTotalSubdomains(PC pc, PetscInt N, IS is[], IS is_local[]) ``` Collective, all MPI ranks must pass in the same array of `IS` ## Input Parameters - ***pc -*** the preconditioner context - ***N -*** the number of subdomains for all processors - ***is -*** the index sets that define the subdomains for all processors (or `NULL` to ask PETSc to determine the subdomains) - ***is_local -*** the index sets that define the local part of the subdomains for this processor (or `NULL` to not provide this information) ## Options Database Key - ***-pc_asm_blocks -*** Sets total blocks ## Notes Currently you cannot use this to set the actual subdomains with the argument is or is_local. By default the `PCASM` preconditioner uses 1 block per processor. These index sets cannot be destroyed until after completion of the linear solves for which the `PCASM` preconditioner is being used. Use `PCASMSetLocalSubdomains()` to set local subdomains. The `IS` numbering is in the parallel, global numbering of the vector for both is and is_local ## See Also `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, `PCASMCreateSubdomains2D()`, `PCGASM` ## Level advanced ## Location src/ksp/pc/impls/asm/asm.c ## Examples src/ksp/ksp/tutorials/ex8.c
## Implementations PCASMSetTotalSubdomains_ASM in src/ksp/pc/impls/asm/asm.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/pc/impls/asm/asm.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)