Actual source code: checkptr.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1:  #include <petsc/private/petscimpl.h>
  2:  #include <petscvalgrind.h>
  3: #if defined(PETSC_HAVE_CUDA)
  4: #include <cuda_runtime.h>
  5:  #include <petsccublas.h>
  6: #endif

  8: static PetscInt petsc_checkpointer_intensity = 1;

 10: /*@
 12:    confirm whether the address is valid.  An intensity of 0 never uses signal handlers, 1 uses them when not in a "hot"
 13:    function, and intensity of 2 always uses a signal handler.

 15:    Not Collective

 17:    Input Arguments:
 18: .  intensity - how much to check pointers for validity

 20:    Options Database:
 21: .  -check_pointer_intensity - intensity (0, 1, or 2)

 23:    Level: advanced

 26: @*/
 28: {

 31:   switch (intensity) {
 32:   case 0:
 33:   case 1:
 34:   case 2:
 35:     petsc_checkpointer_intensity = intensity;
 36:     break;
 37:   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Intensity %D not in 0,1,2",intensity);
 38:   }
 39:   return(0);
 40: }

 42: /* ---------------------------------------------------------------------------------------*/

 44: #if defined(PETSC_HAVE_SETJMP_H)
 45: #include <setjmp.h>
 46: static jmp_buf PetscSegvJumpBuf;
 47: static PetscBool PetscSegvJumpBuf_set;

 49: /*@C
 50:    PetscSignalSegvCheckPointerOrMpi - To be called from a signal handler for SIGSEGV.  If the signal was received
 52:    with no effect. This function is called automatically by PetscSignalHandlerDefault().

 54:    Not Collective

 56:    Level: developer

 58: .seealso: PetscPushSignalHandler()
 59: @*/
 60: void PetscSignalSegvCheckPointerOrMpi(void) {
 61:   if (PetscSegvJumpBuf_set) longjmp(PetscSegvJumpBuf,1);
 62: }

 64: /*@C

 67:    Not Collective

 69:    Input Parameters:
 70: +     ptr - the pointer
 71: -     dtype - the type of data the pointer is suppose to point to

 73:    Level: developer

 76: @*/
 78: {

 80:   if (PETSC_RUNNING_ON_VALGRIND) return PETSC_TRUE;
 81:   if (!ptr) return PETSC_FALSE;
 82:   if (petsc_checkpointer_intensity < 1) return PETSC_TRUE;

 84:   /* Skip the verbose check if we are inside a hot function. */
 85:   if (petscstack && petscstack->hotdepth > 0 && petsc_checkpointer_intensity < 2) return PETSC_TRUE;

 87:   PetscSegvJumpBuf_set = PETSC_TRUE;

 89:   if (setjmp(PetscSegvJumpBuf)) {
 90:     /* A segv was triggered in the code below hence we return with an error code */
 91:     PetscSegvJumpBuf_set = PETSC_FALSE;
 92:     return PETSC_FALSE;
 93:   } else {
 94:     switch (dtype) {
 95:     case PETSC_INT:{
 96:       PETSC_UNUSED PetscInt x = (PetscInt)*(volatile PetscInt*)ptr;
 97:       break;
 98:     }
 99: #if defined(PETSC_USE_COMPLEX)
100:     case PETSC_SCALAR:{         /* C++ is seriously dysfunctional with volatile std::complex. */
101: #if defined(PETSC_USE_CXXCOMPLEX)
102:       PetscReal xreal = ((volatile PetscReal*)ptr)[0],ximag = ((volatile PetscReal*)ptr)[1];
103:       PETSC_UNUSED volatile PetscScalar x = xreal + PETSC_i*ximag;
104: #else
105:       PETSC_UNUSED PetscScalar x = *(volatile PetscScalar*)ptr;
106: #endif
107:       break;
108:     }
109: #endif
110:     case PETSC_REAL:{
111:       PETSC_UNUSED PetscReal x = *(volatile PetscReal*)ptr;
112:       break;
113:     }
114:     case PETSC_BOOL:{
115:       PETSC_UNUSED PetscBool x = *(volatile PetscBool*)ptr;
116:       break;
117:     }
118:     case PETSC_ENUM:{
119:       PETSC_UNUSED PetscEnum x = *(volatile PetscEnum*)ptr;
120:       break;
121:     }
122:     case PETSC_CHAR:{
123:       PETSC_UNUSED char x = *(volatile char*)ptr;
124:       break;
125:     }
126:     case PETSC_OBJECT:{
127:       PETSC_UNUSED volatile PetscClassId classid = ((PetscObject)ptr)->classid;
128:       break;
129:     }
130:     default:;
131:     }
132:   }
133:   PetscSegvJumpBuf_set = PETSC_FALSE;
134:   return PETSC_TRUE;
135: }

137: #if defined (PETSC_HAVE_CUDA)
139: {
140:   cudaError_t cerr=cudaSuccess;
141:   int         ierr,hbuf[2]={1,0},*dbuf=NULL;
142:   PetscBool   awareness=PETSC_FALSE;

144:   cerr = cudaMalloc((void**)&dbuf,sizeof(int)*2);if (cerr != cudaSuccess) return PETSC_FALSE;
145:   cerr = cudaMemcpy(dbuf,hbuf,sizeof(int)*2,cudaMemcpyHostToDevice);if (cerr != cudaSuccess) return PETSC_FALSE;

147:   PetscSegvJumpBuf_set = PETSC_TRUE;
148:   if (setjmp(PetscSegvJumpBuf)) {
149:     /* If a segv was triggered in the MPI_Allreduce below, it is very likely due to the MPI is not GPU-aware */
150:     awareness = PETSC_FALSE;
151:   } else {
152:     MPI_Allreduce(dbuf,dbuf+1,1,MPI_INT,MPI_SUM,PETSC_COMM_WORLD);
153:     if (!ierr) awareness = PETSC_TRUE;
154:   }
155:   PetscSegvJumpBuf_set = PETSC_FALSE;
156:   cerr = cudaFree(dbuf);if (cerr != cudaSuccess) return PETSC_FALSE;
157:   return awareness;
158: }
159: #endif
160: #else
161: void PetscSignalSegvCheckPointerOrMpi(void) {
162:   return;
163: }

166: {
167:   if (!ptr) return PETSC_FALSE;
168:   return PETSC_TRUE;
169: }

172: {
173:   /* If no setjmp (rare), return true and let users code run (and segfault if they should) */
174:   return PETSC_TRUE;
175: }
176: #endif