17d0a6c19SBarry Smith 2e5c89e4eSSatish Balay /* 3e5c89e4eSSatish Balay Code that allows a user to dictate what malloc() PETSc uses. 4e5c89e4eSSatish Balay */ 5c6db04a5SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 6ba282f50SJed Brown #include <stdarg.h> 7e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H) 8e5c89e4eSSatish Balay #include <malloc.h> 9e5c89e4eSSatish Balay #endif 10de1d6c17SHong Zhang #if defined(PETSC_HAVE_MEMKIND) 11de1d6c17SHong Zhang #include <memkind.h> 12*e3acc61dSHong Zhang typedef enum {PETSC_MK_DEFAULT=0,PETSC_MK_HBW_PREFERRED=1} PetscMemkindType; 13*e3acc61dSHong Zhang PetscMemkindType currentmktype = PETSC_MK_HBW_PREFERRED; 14*e3acc61dSHong Zhang PetscMemkindType previousmktype = PETSC_MK_HBW_PREFERRED; 15de1d6c17SHong Zhang #endif 16e5c89e4eSSatish Balay /* 17e5c89e4eSSatish Balay We want to make sure that all mallocs of double or complex numbers are complex aligned. 18e5c89e4eSSatish Balay 1) on systems with memalign() we call that routine to get an aligned memory location 19e5c89e4eSSatish Balay 2) on systems without memalign() we 20e5c89e4eSSatish Balay - allocate one sizeof(PetscScalar) extra space 21e5c89e4eSSatish Balay - we shift the pointer up slightly if needed to get PetscScalar aligned 220700a824SBarry Smith - if shifted we store at ptr[-1] the amount of shift (plus a classid) 23e5c89e4eSSatish Balay */ 240700a824SBarry Smith #define SHIFT_CLASSID 456123 25e5c89e4eSSatish Balay 26efca3c55SSatish Balay PetscErrorCode PetscMallocAlign(size_t mem,int line,const char func[],const char file[],void **result) 27e5c89e4eSSatish Balay { 28f0ba7cfcSLisandro Dalcin if (!mem) { *result = NULL; return 0; } 29fc2a7144SHong Zhang #if defined(PETSC_HAVE_MEMKIND) 30fc2a7144SHong Zhang { 31fc2a7144SHong Zhang int ierr; 32fc2a7144SHong Zhang if (!currentmktype) ierr = memkind_posix_memalign(MEMKIND_DEFAULT,result,PETSC_MEMALIGN,mem); 33*e3acc61dSHong Zhang else ierr = memkind_posix_memalign(MEMKIND_HBW_PREFERRED,result,PETSC_MEMALIGN,mem); 34fc2a7144SHong Zhang if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MEM,"Memory requested with memkind %.0f",(PetscLogDouble)mem); 35fc2a7144SHong Zhang } 36fc2a7144SHong Zhang #else 37e5c89e4eSSatish Balay # if defined(PETSC_HAVE_DOUBLE_ALIGN_MALLOC) && (PETSC_MEMALIGN == 8) 38e5c89e4eSSatish Balay *result = malloc(mem); 39e5c89e4eSSatish Balay # elif defined(PETSC_HAVE_MEMALIGN) 40e5c89e4eSSatish Balay *result = memalign(PETSC_MEMALIGN,mem); 41e5c89e4eSSatish Balay # else 42e5c89e4eSSatish Balay { 43e5c89e4eSSatish Balay /* 44e5c89e4eSSatish Balay malloc space for two extra chunks and shift ptr 1 + enough to get it PetscScalar aligned 45e5c89e4eSSatish Balay */ 465b232624SHong Zhang int *ptr = (int*)malloc(mem + 2*PETSC_MEMALIGN); 47e5c89e4eSSatish Balay if (ptr) { 48f0ba7cfcSLisandro Dalcin int shift = (int)(((PETSC_UINTPTR_T) ptr) % PETSC_MEMALIGN); 49e5c89e4eSSatish Balay shift = (2*PETSC_MEMALIGN - shift)/sizeof(int); 500700a824SBarry Smith ptr[shift-1] = shift + SHIFT_CLASSID; 51e5c89e4eSSatish Balay ptr += shift; 52e5c89e4eSSatish Balay *result = (void*)ptr; 53e5c89e4eSSatish Balay } else { 54f0ba7cfcSLisandro Dalcin *result = NULL; 55e5c89e4eSSatish Balay } 56e5c89e4eSSatish Balay } 57e5c89e4eSSatish Balay # endif 58fc2a7144SHong Zhang #endif 59f0ba7cfcSLisandro Dalcin if (!*result) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_MEM,PETSC_ERROR_INITIAL,"Memory requested %.0f",(PetscLogDouble)mem); 60e5c89e4eSSatish Balay return 0; 61e5c89e4eSSatish Balay } 62e5c89e4eSSatish Balay 63efca3c55SSatish Balay PetscErrorCode PetscFreeAlign(void *ptr,int line,const char func[],const char file[]) 64e5c89e4eSSatish Balay { 65f0ba7cfcSLisandro Dalcin if (!ptr) return 0; 66fc2a7144SHong Zhang #if defined(PETSC_HAVE_MEMKIND) 67fc2a7144SHong Zhang memkind_free(0,ptr); /* specify the kind to 0 so that memkind will look up for the right type */ 68fc2a7144SHong Zhang #else 69e5c89e4eSSatish Balay # if (!(defined(PETSC_HAVE_DOUBLE_ALIGN_MALLOC) && (PETSC_MEMALIGN == 8)) && !defined(PETSC_HAVE_MEMALIGN)) 70f0ba7cfcSLisandro Dalcin { 71e5c89e4eSSatish Balay /* 72e5c89e4eSSatish Balay Previous int tells us how many ints the pointer has been shifted from 73e5c89e4eSSatish Balay the original address provided by the system malloc(). 74e5c89e4eSSatish Balay */ 75f0ba7cfcSLisandro Dalcin int shift = *(((int*)ptr)-1) - SHIFT_CLASSID; 76efca3c55SSatish Balay if (shift > PETSC_MEMALIGN-1) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"Likely memory corruption in heap"); 77efca3c55SSatish Balay if (shift < 0) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"Likely memory corruption in heap"); 78e5c89e4eSSatish Balay ptr = (void*)(((int*)ptr) - shift); 79e5c89e4eSSatish Balay } 80f0ba7cfcSLisandro Dalcin # endif 81e5c89e4eSSatish Balay 82e5c89e4eSSatish Balay # if defined(PETSC_HAVE_FREE_RETURN_INT) 83e5c89e4eSSatish Balay int err = free(ptr); 84efca3c55SSatish Balay if (err) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"System free returned error %d\n",err); 85e5c89e4eSSatish Balay # else 86e5c89e4eSSatish Balay free(ptr); 87e5c89e4eSSatish Balay # endif 88fc2a7144SHong Zhang #endif 89e5c89e4eSSatish Balay return 0; 90e5c89e4eSSatish Balay } 91e5c89e4eSSatish Balay 923221ece2SMatthew G. Knepley PetscErrorCode PetscReallocAlign(size_t mem, int line, const char func[], const char file[], void **result) 933221ece2SMatthew G. Knepley { 94c22f1541SToby Isaac PetscErrorCode ierr; 95c22f1541SToby Isaac 96c22f1541SToby Isaac if (!mem) { 97c22f1541SToby Isaac ierr = PetscFreeAlign(*result, line, func, file); 98c22f1541SToby Isaac if (ierr) return ierr; 99c22f1541SToby Isaac *result = NULL; 100c22f1541SToby Isaac return 0; 101c22f1541SToby Isaac } 102fc2a7144SHong Zhang #if defined(PETSC_HAVE_MEMKIND) 103fc2a7144SHong Zhang if (!currentmktype) *result = memkind_realloc(MEMKIND_DEFAULT,*result,mem); 104*e3acc61dSHong Zhang else *result = memkind_realloc(MEMKIND_HBW_PREFERRED,*result,mem); 105fc2a7144SHong Zhang #else 1063221ece2SMatthew G. Knepley # if (!(defined(PETSC_HAVE_DOUBLE_ALIGN_MALLOC) && (PETSC_MEMALIGN == 8)) && !defined(PETSC_HAVE_MEMALIGN)) 1073221ece2SMatthew G. Knepley { 1083221ece2SMatthew G. Knepley /* 1093221ece2SMatthew G. Knepley Previous int tells us how many ints the pointer has been shifted from 1103221ece2SMatthew G. Knepley the original address provided by the system malloc(). 1113221ece2SMatthew G. Knepley */ 1123221ece2SMatthew G. Knepley int shift = *(((int*)*result)-1) - SHIFT_CLASSID; 1133221ece2SMatthew G. Knepley if (shift > PETSC_MEMALIGN-1) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"Likely memory corruption in heap"); 1143221ece2SMatthew G. Knepley if (shift < 0) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"Likely memory corruption in heap"); 1153221ece2SMatthew G. Knepley *result = (void*)(((int*)*result) - shift); 1163221ece2SMatthew G. Knepley } 1173221ece2SMatthew G. Knepley # endif 1183221ece2SMatthew G. Knepley 119c22f1541SToby Isaac # if (defined(PETSC_HAVE_DOUBLE_ALIGN_MALLOC) && (PETSC_MEMALIGN == 8)) || defined(PETSC_HAVE_MEMALIGN) 1203221ece2SMatthew G. Knepley *result = realloc(*result, mem); 1213221ece2SMatthew G. Knepley # else 1223221ece2SMatthew G. Knepley { 1233221ece2SMatthew G. Knepley /* 1243221ece2SMatthew G. Knepley malloc space for two extra chunks and shift ptr 1 + enough to get it PetscScalar aligned 1253221ece2SMatthew G. Knepley */ 1263221ece2SMatthew G. Knepley int *ptr = (int *) realloc(*result, mem + 2*PETSC_MEMALIGN); 1273221ece2SMatthew G. Knepley if (ptr) { 1283221ece2SMatthew G. Knepley int shift = (int)(((PETSC_UINTPTR_T) ptr) % PETSC_MEMALIGN); 1293221ece2SMatthew G. Knepley shift = (2*PETSC_MEMALIGN - shift)/sizeof(int); 1303221ece2SMatthew G. Knepley ptr[shift-1] = shift + SHIFT_CLASSID; 1313221ece2SMatthew G. Knepley ptr += shift; 1323221ece2SMatthew G. Knepley *result = (void*)ptr; 1333221ece2SMatthew G. Knepley } else { 1343221ece2SMatthew G. Knepley *result = NULL; 1353221ece2SMatthew G. Knepley } 1363221ece2SMatthew G. Knepley } 1373221ece2SMatthew G. Knepley # endif 138fc2a7144SHong Zhang #endif 1393221ece2SMatthew G. Knepley if (!*result) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_MEM,PETSC_ERROR_INITIAL,"Memory requested %.0f",(PetscLogDouble)mem); 140c22f1541SToby Isaac #if defined(PETSC_HAVE_MEMALIGN) 141c22f1541SToby Isaac /* There are no standard guarantees that realloc() maintains the alignment of memalign(), so I think we have to 142c22f1541SToby Isaac * realloc and, if the alignment is wrong, malloc/copy/free. */ 143c22f1541SToby Isaac if (((size_t) (*result)) % PETSC_MEMALIGN) { 144c22f1541SToby Isaac void *newResult; 145fc2a7144SHong Zhang # if defined(PETSC_HAVE_MEMKIND) 146fc2a7144SHong Zhang { 147fc2a7144SHong Zhang int ierr; 148fc2a7144SHong Zhang if (!currentmktype) ierr = memkind_posix_memalign(MEMKIND_DEFAULT,&newResult,PETSC_MEMALIGN,mem); 149*e3acc61dSHong Zhang else ierr = memkind_posix_memalign(MEMKIND_HBW_PREFERRED,&newResult,PETSC_MEMALIGN,mem); 150fc2a7144SHong Zhang if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MEM,"Memory requested with memkind %.0f",(PetscLogDouble)mem); 151fc2a7144SHong Zhang } 152fc2a7144SHong Zhang # else 153c22f1541SToby Isaac newResult = memalign(PETSC_MEMALIGN,mem); 154fc2a7144SHong Zhang # endif 155c22f1541SToby Isaac if (!newResult) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_MEM,PETSC_ERROR_INITIAL,"Memory requested %.0f",(PetscLogDouble)mem); 156c22f1541SToby Isaac ierr = PetscMemcpy(newResult,*result,mem); 157c22f1541SToby Isaac if (ierr) return ierr; 158c22f1541SToby Isaac # if defined(PETSC_HAVE_FREE_RETURN_INT) 159c22f1541SToby Isaac { 160c22f1541SToby Isaac int err = free(*result); 161c22f1541SToby Isaac if (err) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"System free returned error %d\n",err); 162c22f1541SToby Isaac } 163c22f1541SToby Isaac # else 164de1d6c17SHong Zhang # if defined(PETSC_HAVE_MEMKIND) 165de1d6c17SHong Zhang memkind_free(0,*result); 166de1d6c17SHong Zhang # else 167c22f1541SToby Isaac free(*result); 168c22f1541SToby Isaac # endif 169de1d6c17SHong Zhang # endif 170c22f1541SToby Isaac *result = newResult; 171c22f1541SToby Isaac } 172c22f1541SToby Isaac #endif 1733221ece2SMatthew G. Knepley return 0; 1743221ece2SMatthew G. Knepley } 1753221ece2SMatthew G. Knepley 176efca3c55SSatish Balay PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**) = PetscMallocAlign; 177efca3c55SSatish Balay PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]) = PetscFreeAlign; 1783221ece2SMatthew G. Knepley PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**) = PetscReallocAlign; 179e5c89e4eSSatish Balay 180ace3abfcSBarry Smith PetscBool petscsetmallocvisited = PETSC_FALSE; 181e5c89e4eSSatish Balay 182e5c89e4eSSatish Balay /*@C 1831d1a0024SBarry Smith PetscMallocSet - Sets the routines used to do mallocs and frees. 184e5c89e4eSSatish Balay This routine MUST be called before PetscInitialize() and may be 185e5c89e4eSSatish Balay called only once. 186e5c89e4eSSatish Balay 187e5c89e4eSSatish Balay Not Collective 188e5c89e4eSSatish Balay 189e5c89e4eSSatish Balay Input Parameters: 190e5c89e4eSSatish Balay + malloc - the malloc routine 191e5c89e4eSSatish Balay - free - the free routine 192e5c89e4eSSatish Balay 193e5c89e4eSSatish Balay Level: developer 194e5c89e4eSSatish Balay 195e5c89e4eSSatish Balay Concepts: malloc 196e5c89e4eSSatish Balay Concepts: memory^allocation 197e5c89e4eSSatish Balay 198e5c89e4eSSatish Balay @*/ 199efca3c55SSatish Balay PetscErrorCode PetscMallocSet(PetscErrorCode (*imalloc)(size_t,int,const char[],const char[],void**), 200efca3c55SSatish Balay PetscErrorCode (*ifree)(void*,int,const char[],const char[])) 201e5c89e4eSSatish Balay { 202e5c89e4eSSatish Balay PetscFunctionBegin; 203e32f2f54SBarry Smith if (petscsetmallocvisited && (imalloc != PetscTrMalloc || ifree != PetscTrFree)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"cannot call multiple times"); 204e5c89e4eSSatish Balay PetscTrMalloc = imalloc; 205e5c89e4eSSatish Balay PetscTrFree = ifree; 206e5c89e4eSSatish Balay petscsetmallocvisited = PETSC_TRUE; 207e5c89e4eSSatish Balay PetscFunctionReturn(0); 208e5c89e4eSSatish Balay } 209e5c89e4eSSatish Balay 210e5c89e4eSSatish Balay /*@C 2111d1a0024SBarry Smith PetscMallocClear - Resets the routines used to do mallocs and frees to the 212e5c89e4eSSatish Balay defaults. 213e5c89e4eSSatish Balay 214e5c89e4eSSatish Balay Not Collective 215e5c89e4eSSatish Balay 216e5c89e4eSSatish Balay Level: developer 217e5c89e4eSSatish Balay 218e5c89e4eSSatish Balay Notes: 219e5c89e4eSSatish Balay In general one should never run a PETSc program with different malloc() and 220e5c89e4eSSatish Balay free() settings for different parts; this is because one NEVER wants to 221e5c89e4eSSatish Balay free() an address that was malloced by a different memory management system 222e5c89e4eSSatish Balay 223e5c89e4eSSatish Balay @*/ 2247087cfbeSBarry Smith PetscErrorCode PetscMallocClear(void) 225e5c89e4eSSatish Balay { 226e5c89e4eSSatish Balay PetscFunctionBegin; 227e5c89e4eSSatish Balay PetscTrMalloc = PetscMallocAlign; 228e5c89e4eSSatish Balay PetscTrFree = PetscFreeAlign; 229e5c89e4eSSatish Balay petscsetmallocvisited = PETSC_FALSE; 230e5c89e4eSSatish Balay PetscFunctionReturn(0); 231e5c89e4eSSatish Balay } 232b44d5720SBarry Smith 233b44d5720SBarry Smith PetscErrorCode PetscMemoryTrace(const char label[]) 234b44d5720SBarry Smith { 235b44d5720SBarry Smith PetscErrorCode ierr; 236b44d5720SBarry Smith PetscLogDouble mem,mal; 237b44d5720SBarry Smith static PetscLogDouble oldmem = 0,oldmal = 0; 238b44d5720SBarry Smith 239b44d5720SBarry Smith PetscFunctionBegin; 240b44d5720SBarry Smith ierr = PetscMemoryGetCurrentUsage(&mem);CHKERRQ(ierr); 241b44d5720SBarry Smith ierr = PetscMallocGetCurrentUsage(&mal);CHKERRQ(ierr); 242b44d5720SBarry Smith 243b44d5720SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"%s High water %8.3f MB increase %8.3f MB Current %8.3f MB increase %8.3f MB\n",label,mem*1e-6,(mem - oldmem)*1e-6,mal*1e-6,(mal - oldmal)*1e-6);CHKERRQ(ierr); 244b44d5720SBarry Smith oldmem = mem; 245b44d5720SBarry Smith oldmal = mal; 246b44d5720SBarry Smith PetscFunctionReturn(0); 247b44d5720SBarry Smith } 24813850c04SHong Zhang 24950a41461SHong Zhang static PetscErrorCode (*PetscTrMallocOld)(size_t,int,const char[],const char[],void**) = PetscMallocAlign; 25050a41461SHong Zhang static PetscErrorCode (*PetscTrFreeOld)(void*,int,const char[],const char[]) = PetscFreeAlign; 251de1d6c17SHong Zhang 252de1d6c17SHong Zhang /*@C 253de1d6c17SHong Zhang PetscMallocSetDRAM - Set PetscMalloc to use DRAM. 254de1d6c17SHong Zhang If memkind is available, change the memkind type. Otherwise, switch the 255de1d6c17SHong Zhang current malloc and free routines to the PetscMallocAlign and 256de1d6c17SHong Zhang PetscFreeAlign (PETSc default). 257de1d6c17SHong Zhang 258de1d6c17SHong Zhang Not Collective 259de1d6c17SHong Zhang 260de1d6c17SHong Zhang Level: developer 261de1d6c17SHong Zhang 262de1d6c17SHong Zhang Notes: 263de1d6c17SHong Zhang This provides a way to do the allocation on DRAM temporarily. One 264de1d6c17SHong Zhang can switch back to the previous choice by calling PetscMallocReset(). 265de1d6c17SHong Zhang 266de1d6c17SHong Zhang .seealso: PetscMallocReset() 267de1d6c17SHong Zhang @*/ 26813850c04SHong Zhang PetscErrorCode PetscMallocSetDRAM(void) 26913850c04SHong Zhang { 27013850c04SHong Zhang PetscFunctionBegin; 271de1d6c17SHong Zhang if (PetscTrMalloc == PetscMallocAlign) { 272de1d6c17SHong Zhang #if defined(PETSC_HAVE_MEMKIND) 273de1d6c17SHong Zhang previousmktype = currentmktype; 274de1d6c17SHong Zhang currentmktype = PETSC_MK_DEFAULT; 275de1d6c17SHong Zhang #endif 276de1d6c17SHong Zhang } else { 27713850c04SHong Zhang /* Save the previous choice */ 27813850c04SHong Zhang PetscTrMallocOld = PetscTrMalloc; 27913850c04SHong Zhang PetscTrFreeOld = PetscTrFree; 28013850c04SHong Zhang PetscTrMalloc = PetscMallocAlign; 28113850c04SHong Zhang PetscTrFree = PetscFreeAlign; 282de1d6c17SHong Zhang } 28313850c04SHong Zhang PetscFunctionReturn(0); 28413850c04SHong Zhang } 28513850c04SHong Zhang 286de1d6c17SHong Zhang /*@C 287de1d6c17SHong Zhang PetscMallocResetDRAM - Reset the changes made by PetscMallocSetDRAM 288de1d6c17SHong Zhang 289de1d6c17SHong Zhang Not Collective 290de1d6c17SHong Zhang 291de1d6c17SHong Zhang Level: developer 292de1d6c17SHong Zhang 293de1d6c17SHong Zhang .seealso: PetscMallocSetDRAM() 294de1d6c17SHong Zhang @*/ 29513850c04SHong Zhang PetscErrorCode PetscMallocResetDRAM(void) 29613850c04SHong Zhang { 29713850c04SHong Zhang PetscFunctionBegin; 298de1d6c17SHong Zhang if (PetscTrMalloc == PetscMallocAlign) { 299de1d6c17SHong Zhang #if defined(PETSC_HAVE_MEMKIND) 300de1d6c17SHong Zhang currentmktype = previousmktype; 301de1d6c17SHong Zhang #endif 302de1d6c17SHong Zhang } else { 30313850c04SHong Zhang /* Reset to the previous choice */ 30413850c04SHong Zhang PetscTrMalloc = PetscTrMallocOld; 30513850c04SHong Zhang PetscTrFree = PetscTrFreeOld; 306de1d6c17SHong Zhang } 30713850c04SHong Zhang PetscFunctionReturn(0); 30813850c04SHong Zhang } 309ba282f50SJed Brown 310ba282f50SJed Brown static PetscBool petscmalloccoalesce = 311ba282f50SJed Brown #if defined(PETSC_USE_MALLOC_COALESCED) 312ba282f50SJed Brown PETSC_TRUE; 313ba282f50SJed Brown #else 314ba282f50SJed Brown PETSC_FALSE; 315ba282f50SJed Brown #endif 316ba282f50SJed Brown 317ba282f50SJed Brown /*@C 318ba282f50SJed Brown PetscMallocSetCoalesce - Use coalesced malloc when allocating groups of objects 319ba282f50SJed Brown 320ba282f50SJed Brown Not Collective 321ba282f50SJed Brown 322ba282f50SJed Brown Input Parameters: 323ba282f50SJed Brown . coalesce - PETSC_TRUE to use coalesced malloc for multi-object allocation. 324ba282f50SJed Brown 325ba282f50SJed Brown Options Database Keys: 326ba282f50SJed Brown . -malloc_coalesce - turn coalesced malloc on or off 327ba282f50SJed Brown 328ba282f50SJed Brown Note: 329ba282f50SJed Brown PETSc uses coalesced malloc by default for optimized builds and not for debugging builds. This default can be changed via the command-line option -malloc_coalesce or by calling this function. 330ba282f50SJed Brown 331ba282f50SJed Brown Level: developer 332ba282f50SJed Brown 333ba282f50SJed Brown .seealso: PetscMallocA() 334ba282f50SJed Brown @*/ 335ba282f50SJed Brown PetscErrorCode PetscMallocSetCoalesce(PetscBool coalesce) 336ba282f50SJed Brown { 337ba282f50SJed Brown PetscFunctionBegin; 338ba282f50SJed Brown petscmalloccoalesce = coalesce; 339ba282f50SJed Brown PetscFunctionReturn(0); 340ba282f50SJed Brown } 341ba282f50SJed Brown 342ba282f50SJed Brown /*@C 343ba282f50SJed Brown PetscMallocA - Allocate and optionally clear one or more objects, possibly using coalesced malloc 344ba282f50SJed Brown 345ba282f50SJed Brown Not Collective 346ba282f50SJed Brown 347ba282f50SJed Brown Input Parameters: 348ba282f50SJed Brown + n - number of objects to allocate (at least 1) 349ba282f50SJed Brown . lineno - line number to attribute allocation (typically __LINE__) 350ba282f50SJed Brown . function - function to attribute allocation (typically PETSC_FUNCTION_NAME) 351ba282f50SJed Brown . filename - file name to attribute allocation (typically __FILE__) 352ba282f50SJed Brown - bytes0 - first of n object sizes 353ba282f50SJed Brown 354ba282f50SJed Brown Output Parameters: 355ba282f50SJed Brown . ptr0 - first of n pointers to allocate 356ba282f50SJed Brown 357ba282f50SJed Brown Notes: 358ba282f50SJed Brown This function is not normally called directly by users, but rather via the macros PetscMalloc1(), PetscMalloc2(), or PetscCalloc1(), etc. 359ba282f50SJed Brown 360ba282f50SJed Brown Level: developer 361ba282f50SJed Brown 362ba282f50SJed Brown .seealso: PetscMallocAlign(), PetscMallocSet(), PetscMalloc1(), PetscMalloc2(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(), PetscMalloc7(), PetscCalloc1(), PetscCalloc2(), PetscCalloc3(), PetscCalloc4(), PetscCalloc5(), PetscCalloc6(), PetscCalloc7(), PetscFreeA() 363ba282f50SJed Brown @*/ 364ba282f50SJed Brown PetscErrorCode PetscMallocA(int n,PetscBool clear,int lineno,const char *function,const char *filename,size_t bytes0,void *ptr0,...) 365ba282f50SJed Brown { 366ba282f50SJed Brown PetscErrorCode ierr; 367ba282f50SJed Brown va_list Argp; 368ba282f50SJed Brown size_t bytes[8],sumbytes; 369ba282f50SJed Brown void **ptr[8]; 370ba282f50SJed Brown int i; 371ba282f50SJed Brown 372ba282f50SJed Brown PetscFunctionBegin; 373ba282f50SJed Brown if (n > 8) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Attempt to allocate %d objects but only 8 supported",n); 374ba282f50SJed Brown bytes[0] = bytes0; 375ba282f50SJed Brown ptr[0] = (void**)ptr0; 376ba282f50SJed Brown sumbytes = (bytes0 + PETSC_MEMALIGN-1) & ~(PETSC_MEMALIGN-1); 377ba282f50SJed Brown va_start(Argp,ptr0); 378ba282f50SJed Brown for (i=1; i<n; i++) { 379ba282f50SJed Brown bytes[i] = va_arg(Argp,size_t); 380ba282f50SJed Brown ptr[i] = va_arg(Argp,void**); 381ba282f50SJed Brown sumbytes += (bytes[i] + PETSC_MEMALIGN-1) & ~(PETSC_MEMALIGN-1); 382ba282f50SJed Brown } 383ba282f50SJed Brown va_end(Argp); 384ba282f50SJed Brown if (petscmalloccoalesce) { 385ba282f50SJed Brown char *p; 386ba282f50SJed Brown ierr = (*PetscTrMalloc)(sumbytes,lineno,function,filename,(void**)&p);CHKERRQ(ierr); 387ba282f50SJed Brown for (i=0; i<n; i++) { 388ba282f50SJed Brown *ptr[i] = bytes[i] ? p : NULL; 389ba282f50SJed Brown p = (char*)PetscAddrAlign(p + bytes[i]); 390ba282f50SJed Brown } 391ba282f50SJed Brown } else { 392ba282f50SJed Brown for (i=0; i<n; i++) { 393ba282f50SJed Brown ierr = (*PetscTrMalloc)(bytes[i],lineno,function,filename,(void**)ptr[i]);CHKERRQ(ierr); 394ba282f50SJed Brown } 395ba282f50SJed Brown } 396ba282f50SJed Brown if (clear) { 397ba282f50SJed Brown for (i=0; i<n; i++) { 398ba282f50SJed Brown ierr = PetscMemzero(*ptr[i],bytes[i]);CHKERRQ(ierr); 399ba282f50SJed Brown } 400ba282f50SJed Brown } 401ba282f50SJed Brown PetscFunctionReturn(0); 402ba282f50SJed Brown } 403ba282f50SJed Brown 404ba282f50SJed Brown /*@C 405ba282f50SJed Brown PetscFreeA - Free one or more objects, possibly allocated using coalesced malloc 406ba282f50SJed Brown 407ba282f50SJed Brown Not Collective 408ba282f50SJed Brown 409ba282f50SJed Brown Input Parameters: 410ba282f50SJed Brown + n - number of objects to free (at least 1) 411ba282f50SJed Brown . lineno - line number to attribute deallocation (typically __LINE__) 412ba282f50SJed Brown . function - function to attribute deallocation (typically PETSC_FUNCTION_NAME) 413ba282f50SJed Brown . filename - file name to attribute deallocation (typically __FILE__) 414ba282f50SJed Brown - ptr0 ... - first of n pointers to free 415ba282f50SJed Brown 416ba282f50SJed Brown Note: 417ba282f50SJed Brown This function is not normally called directly by users, but rather via the macros PetscFree1(), PetscFree2(), etc. 418ba282f50SJed Brown 419ba282f50SJed Brown Level: developer 420ba282f50SJed Brown 421ba282f50SJed Brown .seealso: PetscMallocAlign(), PetscMallocSet(), PetscMallocA(), PetscFree1(), PetscFree2(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6(), PetscFree7() 422ba282f50SJed Brown @*/ 423ba282f50SJed Brown PetscErrorCode PetscFreeA(int n,int lineno,const char *function,const char *filename,void *ptr0,...) 424ba282f50SJed Brown { 425ba282f50SJed Brown PetscErrorCode ierr; 426ba282f50SJed Brown va_list Argp; 427ba282f50SJed Brown void **ptr[8]; 428ba282f50SJed Brown int i; 429ba282f50SJed Brown 430ba282f50SJed Brown PetscFunctionBegin; 431ba282f50SJed Brown if (n > 8) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Attempt to allocate %d objects but only 8 supported",n); 432ba282f50SJed Brown ptr[0] = (void**)ptr0; 433ba282f50SJed Brown va_start(Argp,ptr0); 434ba282f50SJed Brown for (i=1; i<n; i++) { 435ba282f50SJed Brown ptr[i] = va_arg(Argp,void**); 436ba282f50SJed Brown } 437ba282f50SJed Brown va_end(Argp); 438ba282f50SJed Brown if (petscmalloccoalesce) { 439ba282f50SJed Brown for (i=0; i<n; i++) { /* Find first nonempty allocation */ 440ba282f50SJed Brown if (*ptr[i]) break; 441ba282f50SJed Brown } 442ba282f50SJed Brown while (--n > i) { 443ba282f50SJed Brown *ptr[n] = NULL; 444ba282f50SJed Brown } 445c53cf884SJed Brown ierr = (*PetscTrFree)(*ptr[n],lineno,function,filename);CHKERRQ(ierr); 446ba282f50SJed Brown *ptr[n] = NULL; 447ba282f50SJed Brown } else { 448ba282f50SJed Brown while (--n >= 0) { 449c53cf884SJed Brown ierr = (*PetscTrFree)(*ptr[n],lineno,function,filename);CHKERRQ(ierr); 450ba282f50SJed Brown *ptr[n] = NULL; 451ba282f50SJed Brown } 452ba282f50SJed Brown } 453ba282f50SJed Brown PetscFunctionReturn(0); 454ba282f50SJed Brown } 455