Actual source code: checkptr.c
petsc-3.4.5 2014-06-29
1: #include <petsc-private/petscimpl.h>
3: /* ---------------------------------------------------------------------------------------*/
4: #if defined(PETSC_HAVE_SETJMP_H) && defined(PETSC_HAVE_SIGINFO_T)
5: #include <signal.h>
6: #include <setjmp.h>
7: PETSC_INTERN jmp_buf PetscSegvJumpBuf;
8: PETSC_INTERN void PetscSegv_sigaction(int, siginfo_t*, void *);
10: /*@C
13: Not Collective
15: Input Parameters:
16: + ptr - the pointer
17: - dtype - the type of data the pointer is suppose to point to
19: Level: developer
21: @*/
23: {
24: struct sigaction sa,oldsa;
26: if (PETSC_RUNNING_ON_VALGRIND) return PETSC_TRUE;
27: if (!ptr) return PETSC_FALSE;
29: sigemptyset(&sa.sa_mask);
30: sa.sa_sigaction = PetscSegv_sigaction;
31: sa.sa_flags = SA_SIGINFO;
32: sigaction(SIGSEGV, &sa, &oldsa);
34: if (setjmp(PetscSegvJumpBuf)) {
35: /* A segv was triggered in the code below hence we return with an error code */
36: sigaction(SIGSEGV, &oldsa, NULL);/* reset old signal hanlder */
37: return PETSC_FALSE;
38: } else {
39: switch (dtype) {
40: case PETSC_INT:{
41: PETSC_UNUSED PetscInt x = (PetscInt)*(volatile PetscInt*)ptr;
42: break;
43: }
44: #if defined(PETSC_USE_COMPLEX)
45: case PETSC_SCALAR:{ /* C++ is seriously dysfunctional with volatile std::complex. */
46: PetscReal xreal = ((volatile PetscReal*)ptr)[0],ximag = ((volatile PetscReal*)ptr)[1];
47: PETSC_UNUSED volatile PetscScalar x = xreal + PETSC_i*ximag;
48: break;
49: }
50: #endif
51: case PETSC_REAL:{
52: PETSC_UNUSED PetscReal x = *(volatile PetscReal*)ptr;
53: break;
54: }
55: case PETSC_BOOL:{
56: PETSC_UNUSED PetscBool x = *(volatile PetscBool*)ptr;
57: break;
58: }
59: case PETSC_ENUM:{
60: PETSC_UNUSED PetscEnum x = *(volatile PetscEnum*)ptr;
61: break;
62: }
63: case PETSC_CHAR:{
64: PETSC_UNUSED char *x = *(char*volatile*)ptr;
65: break;
66: }
67: case PETSC_OBJECT:{
68: PETSC_UNUSED volatile PetscClassId classid = ((PetscObject)ptr)->classid;
69: break;
70: }
71: default:;
72: }
73: }
74: sigaction(SIGSEGV, &oldsa, NULL); /* reset old signal hanlder */
75: return PETSC_TRUE;
76: }
77: #else
79: {
80: if (!ptr) return PETSC_FALSE;
81: return PETSC_TRUE;
82: }
83: #endif