xref: /petsc/src/sys/error/checkptr.c (revision 042217e806e4e3ebc5975c28dbf27a273923695f)
1 #include <petsc/private/petscimpl.h>
2 #include <petscvalgrind.h>
3 
4 #if defined(PETSC_HAVE_CUDA)
5   #include <cuda_runtime.h>
6 #endif
7 
8 #if defined(PETSC_HAVE_HIP)
9   #include <hip/hip_runtime.h>
10 #endif
11 
12 static PetscInt petsc_checkpointer_intensity = 1;
13 
14 /*@
15    PetscCheckPointerSetIntensity - An intense pointer check registers a signal handler and attempts to dereference to
16    confirm whether the address is valid.  An intensity of 0 never uses signal handlers, 1 uses them when not in a "hot"
17    function, and intensity of 2 always uses a signal handler.
18 
19    Not Collective
20 
21    Input Arguments:
22 .  intensity - how much to check pointers for validity
23 
24    Options Database:
25 .  -check_pointer_intensity - intensity (0, 1, or 2)
26 
27    Level: advanced
28 
29 .seealso: PetscCheckPointer(), PetscFunctionBeginHot()
30 @*/
31 PetscErrorCode PetscCheckPointerSetIntensity(PetscInt intensity)
32 {
33 
34   PetscFunctionBegin;
35   switch (intensity) {
36   case 0:
37   case 1:
38   case 2:
39     petsc_checkpointer_intensity = intensity;
40     break;
41   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Intensity %D not in 0,1,2",intensity);
42   }
43   PetscFunctionReturn(0);
44 }
45 
46 /* ---------------------------------------------------------------------------------------*/
47 
48 #if defined(PETSC_HAVE_SETJMP_H)
49 #include <setjmp.h>
50 static jmp_buf PetscSegvJumpBuf;
51 static PetscBool PetscSegvJumpBuf_set;
52 
53 /*@C
54    PetscSignalSegvCheckPointerOrMpi - To be called from a signal handler for SIGSEGV.  If the signal was received
55    while executing PetscCheckPointer()/PetscCheckMpiXxxAwareness(), this function longjmps back there, otherwise returns
56    with no effect. This function is called automatically by PetscSignalHandlerDefault().
57 
58    Not Collective
59 
60    Level: developer
61 
62 .seealso: PetscPushSignalHandler()
63 @*/
64 void PetscSignalSegvCheckPointerOrMpi(void)
65 {
66   if (PetscSegvJumpBuf_set) longjmp(PetscSegvJumpBuf,1);
67 }
68 
69 /*@C
70      PetscCheckPointer - Returns PETSC_TRUE if a pointer points to accessible data
71 
72    Not Collective
73 
74    Input Parameters:
75 +     ptr - the pointer
76 -     dtype - the type of data the pointer is suppose to point to
77 
78    Level: developer
79 
80 .seealso: PetscCheckPointerSetIntensity()
81 @*/
82 PetscBool PetscCheckPointer(const void *ptr,PetscDataType dtype)
83 {
84 
85   if (PETSC_RUNNING_ON_VALGRIND) return PETSC_TRUE;
86   if (!ptr) return PETSC_FALSE;
87   if (petsc_checkpointer_intensity < 1) return PETSC_TRUE;
88 
89   /* Skip the verbose check if we are inside a hot function. */
90   if (petscstack && petscstack->hotdepth > 0 && petsc_checkpointer_intensity < 2) return PETSC_TRUE;
91 
92   PetscSegvJumpBuf_set = PETSC_TRUE;
93 
94   if (setjmp(PetscSegvJumpBuf)) {
95     /* A segv was triggered in the code below hence we return with an error code */
96     PetscSegvJumpBuf_set = PETSC_FALSE;
97     return PETSC_FALSE;
98   } else {
99     switch (dtype) {
100     case PETSC_INT:{
101       PETSC_UNUSED PetscInt x = (PetscInt)*(volatile PetscInt*)ptr;
102       break;
103     }
104 #if defined(PETSC_USE_COMPLEX)
105     case PETSC_SCALAR:{         /* C++ is seriously dysfunctional with volatile std::complex. */
106 #if defined(PETSC_USE_CXXCOMPLEX)
107       PetscReal xreal = ((volatile PetscReal*)ptr)[0],ximag = ((volatile PetscReal*)ptr)[1];
108       PETSC_UNUSED volatile PetscScalar x = xreal + PETSC_i*ximag;
109 #else
110       PETSC_UNUSED PetscScalar x = *(volatile PetscScalar*)ptr;
111 #endif
112       break;
113     }
114 #endif
115     case PETSC_REAL:{
116       PETSC_UNUSED PetscReal x = *(volatile PetscReal*)ptr;
117       break;
118     }
119     case PETSC_BOOL:{
120       PETSC_UNUSED PetscBool x = *(volatile PetscBool*)ptr;
121       break;
122     }
123     case PETSC_ENUM:{
124       PETSC_UNUSED PetscEnum x = *(volatile PetscEnum*)ptr;
125       break;
126     }
127     case PETSC_CHAR:{
128       PETSC_UNUSED char x = *(volatile char*)ptr;
129       break;
130     }
131     case PETSC_OBJECT:{
132       PETSC_UNUSED volatile PetscClassId classid = ((PetscObject)ptr)->classid;
133       break;
134     }
135     default:;
136     }
137   }
138   PetscSegvJumpBuf_set = PETSC_FALSE;
139   return PETSC_TRUE;
140 }
141 
142 #define PetscMPICUPMAwarnessCheckFunction  \
143 PetscBool PetscMPICUPMAwarenessCheck(void)  \
144 { \
145   cupmError_t cerr=cupmSuccess; \
146   int         ierr,hbuf[2]={1,0},*dbuf=NULL; \
147   PetscBool   awareness=PETSC_FALSE; \
148   cerr = cupmMalloc((void**)&dbuf,sizeof(int)*2);if (cerr != cupmSuccess) return PETSC_FALSE; \
149   cerr = cupmMemcpy(dbuf,hbuf,sizeof(int)*2,cupmMemcpyHostToDevice);if (cerr != cupmSuccess) return PETSC_FALSE; \
150   PetscSegvJumpBuf_set = PETSC_TRUE; \
151   if (setjmp(PetscSegvJumpBuf)) { \
152     /* If a segv was triggered in the MPI_Allreduce below, it is very likely due to the MPI is not GPU-aware */ \
153     awareness = PETSC_FALSE; \
154   } else { \
155     ierr = MPI_Allreduce(dbuf,dbuf+1,1,MPI_INT,MPI_SUM,PETSC_COMM_SELF); \
156     if (!ierr) awareness = PETSC_TRUE; \
157   } \
158   PetscSegvJumpBuf_set = PETSC_FALSE; \
159   cerr = cupmFree(dbuf);if (cerr != cupmSuccess) return PETSC_FALSE; \
160   return awareness; \
161 }
162 
163 #if defined(PETSC_HAVE_CUDA)
164   #define cupmError_t                   cudaError_t
165   #define cupmMalloc                    cudaMalloc
166   #define cupmMemcpy                    cudaMemcpy
167   #define cupmFree                      cudaFree
168   #define cupmSuccess                   cudaSuccess
169   #define cupmMemcpyHostToDevice        cudaMemcpyHostToDevice
170   #define PetscMPICUPMAwarenessCheck    PetscMPICUDAAwarenessCheck
171   PetscMPICUPMAwarnessCheckFunction
172 #endif
173 
174 #if defined(PETSC_HAVE_HIP)
175   #define cupmError_t                   hipError_t
176   #define cupmMalloc                    hipMalloc
177   #define cupmMemcpy                    hipMemcpy
178   #define cupmFree                      hipFree
179   #define cupmSuccess                   hipSuccess
180   #define cupmMemcpyHostToDevice        hipMemcpyHostToDevice
181   #define PetscMPICUPMAwarenessCheck    PetscMPIHIPAwarenessCheck
182   PetscMPICUPMAwarnessCheckFunction
183 #endif
184 
185 #else
186 void PetscSignalSegvCheckPointerOrMpi(void)
187 {
188   return;
189 }
190 
191 PetscBool PetscCheckPointer(const void *ptr,PETSC_UNUSED PetscDataType dtype)
192 {
193   if (!ptr) return PETSC_FALSE;
194   return PETSC_TRUE;
195 }
196 
197 #if defined (PETSC_HAVE_CUDA)
198 PetscBool PetscMPICUDAAwarenessCheck(void)
199 {
200   /* If no setjmp (rare), return true and let users code run (and segfault if they should) */
201   return PETSC_TRUE;
202 }
203 #endif
204 
205 #if defined (PETSC_HAVE_HIP)
206 PetscBool PetscMPIHIPAwarenessCheck(void)
207 {
208   /* If no setjmp (rare), return true and let users code run (and segfault if they should) */
209   return PETSC_TRUE;
210 }
211 #endif
212 
213 #endif
214