:orphan: # PetscPoisonMemoryRegion Poison a region in memory, making it undereferencable ## Synopsis ``` static inline PetscErrorCode PetscPoisonMemoryRegion(const void *ptr, size_t size) ``` Not Available From Fortran ## Input Parameters - ***ptr -*** The pointer to the start of the region - ***size -*** The size (in bytes) of the region to poison ## Notes `ptr` must not be `NULL`. It is OK to poison the same memory region repeatedly (it is a no-op). Any attempt to dereference the region after this routine returns results in an error being raised. The memory region may be un-poisoned using `PetscUnpoisonMemoryRegion()`, making it safe to dereference again. ## Example Usage ```none PetscInt *array; PetscMalloc1(15, &array); // OK, memory is normal array[0] = 10; array[1] = 15; PetscPoisonMemoryRegion(array, 15 * sizeof(*array)); // ERROR this region is poisoned! array[0] = 10; // ERROR reading is not allowed either! PetscInt v = array[15]; // OK can re-poison the region PetscPoisonMemoryRegion(array, 15 * sizeof(*array)); // OK can re-poison any subregion too PetscPoisonMemoryRegion(array + 5, 1 * sizeof(*array)); PetscUnpoisonMemoryRegion(array, 1 * sizeof(*array)); // OK the first entry has been unpoisoned array[0] = 10; // ERROR the rest of the region is still poisoned! array[1] = 12345; PetscUnpoisonMemoryRegion(array + 10, sizeof(*array)); // OK this region is unpoisoned (even though surrounding memory is still poisoned!) array[10] = 0; PetscInt stack_array[10]; // OK can poison stack memory as well PetscPoisonMemoryRegion(stack_array, 10 * sizeof(*stack_array)); // ERROR stack array is poisoned! stack_array[0] = 10; ``` ## See Also `PetscUnpoisonMemoryRegion()`, `PetscIsRegionPoisoned()` ## Level developer ## Location include/petsc/private/mempoison.h --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/include/petsc/private/mempoison.h) [Index of all Sys routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)