:orphan: # PetscSectionCreateGlobalSectionCensored Create a `PetscSection` describing the globallayout using a local (sequential) `PetscSection` on each MPI process and an `PetscSF` describing the section point overlap. ## Synopsis ``` #include "petscsection.h" PetscErrorCode PetscSectionCreateGlobalSectionCensored(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscInt numExcludes, const PetscInt excludes[], PetscSection *gsection) ``` ## Input Parameters - ***s -*** The `PetscSection` for the local field layout - ***sf -*** The `PetscSF` describing parallel layout of the section points - ***includeConstraints -*** By default this is `PETSC_FALSE`, meaning that the global vector will not possess constrained dofs - ***numExcludes -*** The number of exclusion ranges, this must have the same value on all MPI processes - ***excludes -*** An array [start_0, end_0, start_1, end_1, ...] where there are `numExcludes` pairs and must have the same values on all MPI processes ## Output Parameter - ***gsection -*** The `PetscSection` for the global field layout ## Notes On each MPI process `gsection` inherits the chart of the `s` on that process. This sets negative sizes and offsets to points not owned by this process as defined by `sf` but that are within the local value of the chart of `gsection`. In those locations the value of size is -(size+1) and the value of the offset on the remote process is -(off+1). This routine augments `PetscSectionCreateGlobalSection()` by allowing one to exclude certain ranges in the chart of the `PetscSection` ## Developer Note This is a terrible function name ## See Also [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`, `PetscSectionCreateGlobalSectionCensored()` ## Level advanced ## Location src/vec/is/section/interface/section.c --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/vec/is/section/interface/section.c) [Index of all PetscSection routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)