xref: /petsc/src/sys/error/checkptr.c (revision c2a741ee9832612cf60ba3050316232832f029eb)
1af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
2022afb99SBarry Smith #include <petscvalgrind.h>
3*c2a741eeSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
4*c2a741eeSJunchao Zhang #include <cuda_runtime.h>
5*c2a741eeSJunchao Zhang #include <petsccublas.h>
6*c2a741eeSJunchao Zhang #endif
7d96cc911SJed Brown 
828559dc8SJed Brown static PetscInt petsc_checkpointer_intensity = 1;
928559dc8SJed Brown 
1028559dc8SJed Brown /*@
1128559dc8SJed Brown    PetscCheckPointerSetIntensity - An intense pointer check registers a signal handler and attempts to dereference to
1228559dc8SJed Brown    confirm whether the address is valid.  An intensity of 0 never uses signal handlers, 1 uses them when not in a "hot"
1328559dc8SJed Brown    function, and intensity of 2 always uses a signal handler.
1428559dc8SJed Brown 
1528559dc8SJed Brown    Not Collective
1628559dc8SJed Brown 
1728559dc8SJed Brown    Input Arguments:
1828559dc8SJed Brown .  intensity - how much to check pointers for validity
1928559dc8SJed Brown 
20c2f74817SBarry Smith    Options Database:
215789d1f5SJed Brown .  -check_pointer_intensity - intensity (0, 1, or 2)
22c2f74817SBarry Smith 
2328559dc8SJed Brown    Level: advanced
2428559dc8SJed Brown 
255789d1f5SJed Brown .seealso: PetscCheckPointer(), PetscFunctionBeginHot()
2628559dc8SJed Brown @*/
2728559dc8SJed Brown PetscErrorCode PetscCheckPointerSetIntensity(PetscInt intensity)
2828559dc8SJed Brown {
2928559dc8SJed Brown 
3028559dc8SJed Brown   PetscFunctionBegin;
3128559dc8SJed Brown   switch (intensity) {
3228559dc8SJed Brown   case 0:
3328559dc8SJed Brown   case 1:
3428559dc8SJed Brown   case 2:
3528559dc8SJed Brown     petsc_checkpointer_intensity = intensity;
3628559dc8SJed Brown     break;
3728559dc8SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Intensity %D not in 0,1,2",intensity);
3828559dc8SJed Brown   }
3928559dc8SJed Brown   PetscFunctionReturn(0);
4028559dc8SJed Brown }
4128559dc8SJed Brown 
42d96cc911SJed Brown /* ---------------------------------------------------------------------------------------*/
43718fc407SJed Brown 
44718fc407SJed Brown #if defined(PETSC_HAVE_SETJMP_H)
45d96cc911SJed Brown #include <setjmp.h>
46f8a67e6dSJed Brown static jmp_buf PetscSegvJumpBuf;
47f8a67e6dSJed Brown static PetscBool PetscSegvJumpBuf_set;
48f8a67e6dSJed Brown 
49f8a67e6dSJed Brown /*@C
50*c2a741eeSJunchao Zhang    PetscSignalSegvCheckPointerOrMpi - To be called from a signal handler for SIGSEGV.  If the signal was received
51*c2a741eeSJunchao Zhang    while executing PetscCheckPointer()/PetscCheckMpiGpuAwareness(), this function longjmps back there, otherwise returns
52*c2a741eeSJunchao Zhang    with no effect. This function is called automatically by PetscSignalHandlerDefault().
53f8a67e6dSJed Brown 
54f8a67e6dSJed Brown    Not Collective
55f8a67e6dSJed Brown 
56f8a67e6dSJed Brown    Level: developer
57f8a67e6dSJed Brown 
58f8a67e6dSJed Brown .seealso: PetscPushSignalHandler()
59f8a67e6dSJed Brown @*/
60*c2a741eeSJunchao Zhang void PetscSignalSegvCheckPointerOrMpi(void) {
61f8a67e6dSJed Brown   if (PetscSegvJumpBuf_set) longjmp(PetscSegvJumpBuf,1);
62f8a67e6dSJed Brown }
63d076d156SJed Brown 
64d96cc911SJed Brown /*@C
65d96cc911SJed Brown      PetscCheckPointer - Returns PETSC_TRUE if a pointer points to accessible data
66d96cc911SJed Brown 
67d96cc911SJed Brown    Not Collective
68d96cc911SJed Brown 
69d96cc911SJed Brown    Input Parameters:
70d96cc911SJed Brown +     ptr - the pointer
71d96cc911SJed Brown -     dtype - the type of data the pointer is suppose to point to
72d96cc911SJed Brown 
73d96cc911SJed Brown    Level: developer
74d96cc911SJed Brown 
755789d1f5SJed Brown .seealso: PetscCheckPointerSetIntensity()
76d96cc911SJed Brown @*/
77d96cc911SJed Brown PetscBool PetscCheckPointer(const void *ptr,PetscDataType dtype)
78d96cc911SJed Brown {
79d96cc911SJed Brown 
80d96cc911SJed Brown   if (PETSC_RUNNING_ON_VALGRIND) return PETSC_TRUE;
81d96cc911SJed Brown   if (!ptr) return PETSC_FALSE;
8228559dc8SJed Brown   if (petsc_checkpointer_intensity < 1) return PETSC_TRUE;
83d96cc911SJed Brown 
84a2f94806SJed Brown   /* Skip the verbose check if we are inside a hot function. */
855c25fcd7SBarry Smith   if (petscstack && petscstack->hotdepth > 0 && petsc_checkpointer_intensity < 2) return PETSC_TRUE;
86a2f94806SJed Brown 
87718fc407SJed Brown   PetscSegvJumpBuf_set = PETSC_TRUE;
88d96cc911SJed Brown 
89d96cc911SJed Brown   if (setjmp(PetscSegvJumpBuf)) {
90d96cc911SJed Brown     /* A segv was triggered in the code below hence we return with an error code */
91718fc407SJed Brown     PetscSegvJumpBuf_set = PETSC_FALSE;
92d96cc911SJed Brown     return PETSC_FALSE;
93d96cc911SJed Brown   } else {
94d96cc911SJed Brown     switch (dtype) {
95d96cc911SJed Brown     case PETSC_INT:{
96d96cc911SJed Brown       PETSC_UNUSED PetscInt x = (PetscInt)*(volatile PetscInt*)ptr;
97d96cc911SJed Brown       break;
98d96cc911SJed Brown     }
99d96cc911SJed Brown #if defined(PETSC_USE_COMPLEX)
100d96cc911SJed Brown     case PETSC_SCALAR:{         /* C++ is seriously dysfunctional with volatile std::complex. */
10196d2aba5SSatish Balay #if defined(PETSC_USE_CXXCOMPLEX)
102d96cc911SJed Brown       PetscReal xreal = ((volatile PetscReal*)ptr)[0],ximag = ((volatile PetscReal*)ptr)[1];
103d96cc911SJed Brown       PETSC_UNUSED volatile PetscScalar x = xreal + PETSC_i*ximag;
10496d2aba5SSatish Balay #else
10596d2aba5SSatish Balay       PETSC_UNUSED PetscScalar x = *(volatile PetscScalar*)ptr;
10696d2aba5SSatish Balay #endif
107d96cc911SJed Brown       break;
108d96cc911SJed Brown     }
109d96cc911SJed Brown #endif
110d96cc911SJed Brown     case PETSC_REAL:{
111d96cc911SJed Brown       PETSC_UNUSED PetscReal x = *(volatile PetscReal*)ptr;
112d96cc911SJed Brown       break;
113d96cc911SJed Brown     }
114d96cc911SJed Brown     case PETSC_BOOL:{
115d96cc911SJed Brown       PETSC_UNUSED PetscBool x = *(volatile PetscBool*)ptr;
116d96cc911SJed Brown       break;
117d96cc911SJed Brown     }
118d96cc911SJed Brown     case PETSC_ENUM:{
119d96cc911SJed Brown       PETSC_UNUSED PetscEnum x = *(volatile PetscEnum*)ptr;
120d96cc911SJed Brown       break;
121d96cc911SJed Brown     }
122d96cc911SJed Brown     case PETSC_CHAR:{
123f4e06bcbSJed Brown       PETSC_UNUSED char x = *(volatile char*)ptr;
124d96cc911SJed Brown       break;
125d96cc911SJed Brown     }
126d96cc911SJed Brown     case PETSC_OBJECT:{
127d96cc911SJed Brown       PETSC_UNUSED volatile PetscClassId classid = ((PetscObject)ptr)->classid;
128d96cc911SJed Brown       break;
129d96cc911SJed Brown     }
130d96cc911SJed Brown     default:;
131d96cc911SJed Brown     }
132d96cc911SJed Brown   }
133718fc407SJed Brown   PetscSegvJumpBuf_set = PETSC_FALSE;
134d96cc911SJed Brown   return PETSC_TRUE;
135d96cc911SJed Brown }
136*c2a741eeSJunchao Zhang 
137*c2a741eeSJunchao Zhang #if defined (PETSC_HAVE_CUDA)
138*c2a741eeSJunchao Zhang PetscBool PetscCheckMpiGpuAwareness(void)
139*c2a741eeSJunchao Zhang {
140*c2a741eeSJunchao Zhang   cudaError_t cerr=cudaSuccess;
141*c2a741eeSJunchao Zhang   int         ierr,hbuf[2]={1,0},*dbuf=NULL;
142*c2a741eeSJunchao Zhang   PetscBool   awareness=PETSC_FALSE;
143*c2a741eeSJunchao Zhang 
144*c2a741eeSJunchao Zhang   cerr = cudaMalloc((void**)&dbuf,sizeof(int)*2);if (cerr != cudaSuccess) return PETSC_FALSE;
145*c2a741eeSJunchao Zhang   cerr = cudaMemcpy(dbuf,hbuf,sizeof(int)*2,cudaMemcpyHostToDevice);if (cerr != cudaSuccess) return PETSC_FALSE;
146*c2a741eeSJunchao Zhang 
147*c2a741eeSJunchao Zhang   PetscSegvJumpBuf_set = PETSC_TRUE;
148*c2a741eeSJunchao Zhang   if (setjmp(PetscSegvJumpBuf)) {
149*c2a741eeSJunchao Zhang     /* If a segv was triggered in the MPI_Allreduce below, it is very likely due to the MPI is not GPU-aware */
150*c2a741eeSJunchao Zhang     awareness = PETSC_FALSE;
151*c2a741eeSJunchao Zhang   } else {
152*c2a741eeSJunchao Zhang     ierr = MPI_Allreduce(dbuf,dbuf+1,1,MPI_INT,MPI_SUM,PETSC_COMM_WORLD);
153*c2a741eeSJunchao Zhang     if (!ierr) awareness = PETSC_TRUE;
154*c2a741eeSJunchao Zhang   }
155*c2a741eeSJunchao Zhang   PetscSegvJumpBuf_set = PETSC_FALSE;
156*c2a741eeSJunchao Zhang   cerr = cudaFree(dbuf);if (cerr != cudaSuccess) return PETSC_FALSE;
157*c2a741eeSJunchao Zhang   return awareness;
158*c2a741eeSJunchao Zhang }
159*c2a741eeSJunchao Zhang #endif
160d96cc911SJed Brown #else
161*c2a741eeSJunchao Zhang void PetscSignalSegvCheckPointerOrMpi(void) {
162f8a67e6dSJed Brown   return;
163f8a67e6dSJed Brown }
164f8a67e6dSJed Brown 
165d96cc911SJed Brown PetscBool PetscCheckPointer(const void *ptr,PETSC_UNUSED PetscDataType dtype)
166d96cc911SJed Brown {
167d96cc911SJed Brown   if (!ptr) return PETSC_FALSE;
168d96cc911SJed Brown   return PETSC_TRUE;
169d96cc911SJed Brown }
170*c2a741eeSJunchao Zhang 
171*c2a741eeSJunchao Zhang PetscBool PetscCheckMpiGpuAwareness(void)
172*c2a741eeSJunchao Zhang {
173*c2a741eeSJunchao Zhang   /* If no setjmp (rare), return true and let users code run (and segfault if they should) */
174*c2a741eeSJunchao Zhang   return PETSC_TRUE;
175*c2a741eeSJunchao Zhang }
176d96cc911SJed Brown #endif
177