xref: /petsc/src/sys/memory/mal.c (revision de1d6c17da65edfefa08fc4c71a1bba4ab15b007)
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*/
6e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H)
7e5c89e4eSSatish Balay #include <malloc.h>
8e5c89e4eSSatish Balay #endif
9*de1d6c17SHong Zhang #if defined(PETSC_HAVE_MEMKIND)
10*de1d6c17SHong Zhang #include <memkind.h>
11*de1d6c17SHong Zhang typedef enum {PETSC_MK_DEFAULT=0,PETSC_MK_HBW_PREFERRED=1} PetscMemkindType;
12*de1d6c17SHong Zhang PetscMemkindType currentmktype = PETSC_MK_HBW_PREFERRED;
13*de1d6c17SHong Zhang PetscMemkindType previousmktype = PETSC_MK_HBW_PREFERRED;
14*de1d6c17SHong Zhang #endif
15e5c89e4eSSatish Balay /*
16e5c89e4eSSatish Balay         We want to make sure that all mallocs of double or complex numbers are complex aligned.
17e5c89e4eSSatish Balay     1) on systems with memalign() we call that routine to get an aligned memory location
18e5c89e4eSSatish Balay     2) on systems without memalign() we
19e5c89e4eSSatish Balay        - allocate one sizeof(PetscScalar) extra space
20e5c89e4eSSatish Balay        - we shift the pointer up slightly if needed to get PetscScalar aligned
210700a824SBarry Smith        - if shifted we store at ptr[-1] the amount of shift (plus a classid)
22e5c89e4eSSatish Balay */
230700a824SBarry Smith #define SHIFT_CLASSID 456123
24e5c89e4eSSatish Balay 
25efca3c55SSatish Balay PetscErrorCode  PetscMallocAlign(size_t mem,int line,const char func[],const char file[],void **result)
26e5c89e4eSSatish Balay {
27f0ba7cfcSLisandro Dalcin   if (!mem) { *result = NULL; return 0; }
28e5c89e4eSSatish Balay #if defined(PETSC_HAVE_DOUBLE_ALIGN_MALLOC) && (PETSC_MEMALIGN == 8)
29*de1d6c17SHong Zhang #  if defined(PETSC_HAVE_MEMKIND)
30*de1d6c17SHong Zhang   {
31*de1d6c17SHong Zhang     int ierr;
32*de1d6c17SHong Zhang     if (!currentmktype) ierr = memkind_posix_memalign(MEMKIND_DEFAULT,result,PETSC_MEMALIGN,mem);
33*de1d6c17SHong Zhang     else ierr = memkind_posix_memalign(MEMKIND_HBW_PREFERRED,result,PETSC_MEMALIGN,mem);
34*de1d6c17SHong Zhang     if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MEM,"Memory requested with memkind %.0f",(PetscLogDouble)mem);
35*de1d6c17SHong Zhang   }
36*de1d6c17SHong Zhang #  else
37e5c89e4eSSatish Balay   *result = malloc(mem);
38*de1d6c17SHong Zhang #  endif
39e5c89e4eSSatish Balay #elif defined(PETSC_HAVE_MEMALIGN)
40*de1d6c17SHong Zhang #  if defined(PETSC_HAVE_MEMKIND)
41*de1d6c17SHong Zhang   {
42*de1d6c17SHong Zhang     int ierr;
43*de1d6c17SHong Zhang     if (!currentmktype) ierr = memkind_posix_memalign(MEMKIND_DEFAULT,result,PETSC_MEMALIGN,mem);
44*de1d6c17SHong Zhang     else ierr = memkind_posix_memalign(MEMKIND_HBW_PREFERRED,result,PETSC_MEMALIGN,mem);
45*de1d6c17SHong Zhang     if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MEM,"Memory requested with memkind %.0f",(PetscLogDouble)mem);
46*de1d6c17SHong Zhang   }
47*de1d6c17SHong Zhang #  else
48e5c89e4eSSatish Balay   *result = memalign(PETSC_MEMALIGN,mem);
49*de1d6c17SHong Zhang #  endif
50e5c89e4eSSatish Balay #else
51e5c89e4eSSatish Balay   {
52e5c89e4eSSatish Balay     /*
53e5c89e4eSSatish Balay       malloc space for two extra chunks and shift ptr 1 + enough to get it PetscScalar aligned
54e5c89e4eSSatish Balay     */
55*de1d6c17SHong Zhang #  if defined(PETSC_HAVE_MEMKIND)
56*de1d6c17SHong Zhang     int *ptr,ierr;
57*de1d6c17SHong Zhang     if (!currentmktype) ierr = memkind_posix_memalign(MEMKIND_DEFAULT,&ptr,PETSC_MEMALIGN,mem+2*PETSC_MEMALIGN);
58*de1d6c17SHong Zhang     else ierr = memkind_posix_memalign(MEMKIND_HBW_PREFERRED,&ptr,PETSC_MEMALIGN,mem+2*PETSC_MEMALIGN);
59*de1d6c17SHong Zhang     if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MEM,"Memory requested with memkind %.0f",(PetscLogDouble)(mem+2*PETSC_MEMALIGN));
60*de1d6c17SHong Zhang #  else
61f0ba7cfcSLisandro Dalcin     int *ptr = (int*)malloc(mem + 2*PETSC_MEMALIGN);
62*de1d6c17SHong Zhang #  endif
63e5c89e4eSSatish Balay     if (ptr) {
64f0ba7cfcSLisandro Dalcin       int shift    = (int)(((PETSC_UINTPTR_T) ptr) % PETSC_MEMALIGN);
65e5c89e4eSSatish Balay       shift        = (2*PETSC_MEMALIGN - shift)/sizeof(int);
660700a824SBarry Smith       ptr[shift-1] = shift + SHIFT_CLASSID;
67e5c89e4eSSatish Balay       ptr         += shift;
68e5c89e4eSSatish Balay       *result      = (void*)ptr;
69e5c89e4eSSatish Balay     } else {
70f0ba7cfcSLisandro Dalcin       *result      = NULL;
71e5c89e4eSSatish Balay     }
72e5c89e4eSSatish Balay   }
73e5c89e4eSSatish Balay #endif
74f0ba7cfcSLisandro Dalcin   if (!*result) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_MEM,PETSC_ERROR_INITIAL,"Memory requested %.0f",(PetscLogDouble)mem);
75e5c89e4eSSatish Balay   return 0;
76e5c89e4eSSatish Balay }
77e5c89e4eSSatish Balay 
78efca3c55SSatish Balay PetscErrorCode  PetscFreeAlign(void *ptr,int line,const char func[],const char file[])
79e5c89e4eSSatish Balay {
80f0ba7cfcSLisandro Dalcin   if (!ptr) return 0;
81e5c89e4eSSatish Balay #if (!(defined(PETSC_HAVE_DOUBLE_ALIGN_MALLOC) && (PETSC_MEMALIGN == 8)) && !defined(PETSC_HAVE_MEMALIGN))
82f0ba7cfcSLisandro Dalcin   {
83e5c89e4eSSatish Balay     /*
84e5c89e4eSSatish Balay       Previous int tells us how many ints the pointer has been shifted from
85e5c89e4eSSatish Balay       the original address provided by the system malloc().
86e5c89e4eSSatish Balay     */
87f0ba7cfcSLisandro Dalcin     int shift = *(((int*)ptr)-1) - SHIFT_CLASSID;
88efca3c55SSatish 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");
89efca3c55SSatish Balay     if (shift < 0) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"Likely memory corruption in heap");
90e5c89e4eSSatish Balay     ptr = (void*)(((int*)ptr) - shift);
91e5c89e4eSSatish Balay   }
92f0ba7cfcSLisandro Dalcin #endif
93e5c89e4eSSatish Balay 
94e5c89e4eSSatish Balay #if defined(PETSC_HAVE_FREE_RETURN_INT)
95e5c89e4eSSatish Balay   int err = free(ptr);
96efca3c55SSatish Balay   if (err) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"System free returned error %d\n",err);
97*de1d6c17SHong Zhang #elif defined(PETSC_HAVE_MEMKIND)
98*de1d6c17SHong Zhang   memkind_free(0,ptr); /* specify the kind to 0 so that memkind will look up for the right type */
99e5c89e4eSSatish Balay #else
100e5c89e4eSSatish Balay   free(ptr);
101e5c89e4eSSatish Balay #endif
102e5c89e4eSSatish Balay   return 0;
103e5c89e4eSSatish Balay }
104e5c89e4eSSatish Balay 
1053221ece2SMatthew G. Knepley PetscErrorCode PetscReallocAlign(size_t mem, int line, const char func[], const char file[], void **result)
1063221ece2SMatthew G. Knepley {
107c22f1541SToby Isaac   PetscErrorCode ierr;
108c22f1541SToby Isaac 
109c22f1541SToby Isaac   if (!mem) {
110c22f1541SToby Isaac     ierr = PetscFreeAlign(*result, line, func, file);
111c22f1541SToby Isaac     if (ierr) return ierr;
112c22f1541SToby Isaac     *result = NULL;
113c22f1541SToby Isaac     return 0;
114c22f1541SToby Isaac   }
1153221ece2SMatthew G. Knepley #if (!(defined(PETSC_HAVE_DOUBLE_ALIGN_MALLOC) && (PETSC_MEMALIGN == 8)) && !defined(PETSC_HAVE_MEMALIGN))
1163221ece2SMatthew G. Knepley   {
1173221ece2SMatthew G. Knepley     /*
1183221ece2SMatthew G. Knepley       Previous int tells us how many ints the pointer has been shifted from
1193221ece2SMatthew G. Knepley       the original address provided by the system malloc().
1203221ece2SMatthew G. Knepley     */
1213221ece2SMatthew G. Knepley     int shift = *(((int*)*result)-1) - SHIFT_CLASSID;
1223221ece2SMatthew 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");
1233221ece2SMatthew G. Knepley     if (shift < 0) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"Likely memory corruption in heap");
1243221ece2SMatthew G. Knepley     *result = (void*)(((int*)*result) - shift);
1253221ece2SMatthew G. Knepley   }
1263221ece2SMatthew G. Knepley #endif
1273221ece2SMatthew G. Knepley 
128c22f1541SToby Isaac #if (defined(PETSC_HAVE_DOUBLE_ALIGN_MALLOC) && (PETSC_MEMALIGN == 8)) || defined(PETSC_HAVE_MEMALIGN)
129*de1d6c17SHong Zhang #  if defined(PETSC_HAVE_MEMKIND)
130*de1d6c17SHong Zhang   if (!currentmktype) *result = memkind_realloc(MEMKIND_DEFAULT,*result,mem);
131*de1d6c17SHong Zhang   else *result = memkind_realloc(MEMKIND_HBW_PREFERRED,*result,mem);
132*de1d6c17SHong Zhang #  else
1333221ece2SMatthew G. Knepley   *result = realloc(*result, mem);
134*de1d6c17SHong Zhang #  endif
1353221ece2SMatthew G. Knepley #else
1363221ece2SMatthew G. Knepley   {
1373221ece2SMatthew G. Knepley     /*
1383221ece2SMatthew G. Knepley       malloc space for two extra chunks and shift ptr 1 + enough to get it PetscScalar aligned
1393221ece2SMatthew G. Knepley     */
140*de1d6c17SHong Zhang #  if defined(PETSC_HAVE_MEMKIND)
141*de1d6c17SHong Zhang     int *ptr;
142*de1d6c17SHong Zhang     if (!currentmktype) ptr = (int*)memkind_realloc(MEMKIND_DEFAULT,*result,mem+2*PETSC_MEMALIGN);
143*de1d6c17SHong Zhang     else ptr = (int*)memkind_realloc(MEMKIND_HBW_PREFERRED,*result,mem+2*PETSC_MEMALIGN);
144*de1d6c17SHong Zhang #  else
1453221ece2SMatthew G. Knepley     int *ptr = (int *) realloc(*result, mem + 2*PETSC_MEMALIGN);
146*de1d6c17SHong Zhang #  endif
1473221ece2SMatthew G. Knepley     if (ptr) {
1483221ece2SMatthew G. Knepley       int shift    = (int)(((PETSC_UINTPTR_T) ptr) % PETSC_MEMALIGN);
1493221ece2SMatthew G. Knepley       shift        = (2*PETSC_MEMALIGN - shift)/sizeof(int);
1503221ece2SMatthew G. Knepley       ptr[shift-1] = shift + SHIFT_CLASSID;
1513221ece2SMatthew G. Knepley       ptr         += shift;
1523221ece2SMatthew G. Knepley       *result      = (void*)ptr;
1533221ece2SMatthew G. Knepley     } else {
1543221ece2SMatthew G. Knepley       *result      = NULL;
1553221ece2SMatthew G. Knepley     }
1563221ece2SMatthew G. Knepley   }
1573221ece2SMatthew G. Knepley #endif
1583221ece2SMatthew G. Knepley   if (!*result) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_MEM,PETSC_ERROR_INITIAL,"Memory requested %.0f",(PetscLogDouble)mem);
159c22f1541SToby Isaac #if defined(PETSC_HAVE_MEMALIGN)
160c22f1541SToby Isaac   /* There are no standard guarantees that realloc() maintains the alignment of memalign(), so I think we have to
161c22f1541SToby Isaac    * realloc and, if the alignment is wrong, malloc/copy/free. */
162c22f1541SToby Isaac   if (((size_t) (*result)) % PETSC_MEMALIGN) {
163c22f1541SToby Isaac     void *newResult;
164c22f1541SToby Isaac 
165c22f1541SToby Isaac     newResult = memalign(PETSC_MEMALIGN,mem);
166c22f1541SToby Isaac     if (!newResult) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_MEM,PETSC_ERROR_INITIAL,"Memory requested %.0f",(PetscLogDouble)mem);
167c22f1541SToby Isaac     ierr = PetscMemcpy(newResult,*result,mem);
168c22f1541SToby Isaac     if (ierr) return ierr;
169c22f1541SToby Isaac #  if defined(PETSC_HAVE_FREE_RETURN_INT)
170c22f1541SToby Isaac     {
171c22f1541SToby Isaac       int err = free(*result);
172c22f1541SToby Isaac       if (err) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"System free returned error %d\n",err);
173c22f1541SToby Isaac     }
174c22f1541SToby Isaac #  else
175*de1d6c17SHong Zhang #    if defined(PETSC_HAVE_MEMKIND)
176*de1d6c17SHong Zhang     memkind_free(0,*result);
177*de1d6c17SHong Zhang #    else
178c22f1541SToby Isaac     free(*result);
179c22f1541SToby Isaac #    endif
180*de1d6c17SHong Zhang #  endif
181c22f1541SToby Isaac     *result = newResult;
182c22f1541SToby Isaac   }
183c22f1541SToby Isaac #endif
1843221ece2SMatthew G. Knepley   return 0;
1853221ece2SMatthew G. Knepley }
1863221ece2SMatthew G. Knepley 
187efca3c55SSatish Balay PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**) = PetscMallocAlign;
188efca3c55SSatish Balay PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[])           = PetscFreeAlign;
1893221ece2SMatthew G. Knepley PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**) = PetscReallocAlign;
190e5c89e4eSSatish Balay 
191ace3abfcSBarry Smith PetscBool petscsetmallocvisited = PETSC_FALSE;
192e5c89e4eSSatish Balay 
193e5c89e4eSSatish Balay /*@C
1941d1a0024SBarry Smith    PetscMallocSet - Sets the routines used to do mallocs and frees.
195e5c89e4eSSatish Balay    This routine MUST be called before PetscInitialize() and may be
196e5c89e4eSSatish Balay    called only once.
197e5c89e4eSSatish Balay 
198e5c89e4eSSatish Balay    Not Collective
199e5c89e4eSSatish Balay 
200e5c89e4eSSatish Balay    Input Parameters:
201e5c89e4eSSatish Balay +  malloc - the malloc routine
202e5c89e4eSSatish Balay -  free - the free routine
203e5c89e4eSSatish Balay 
204e5c89e4eSSatish Balay    Level: developer
205e5c89e4eSSatish Balay 
206e5c89e4eSSatish Balay    Concepts: malloc
207e5c89e4eSSatish Balay    Concepts: memory^allocation
208e5c89e4eSSatish Balay 
209e5c89e4eSSatish Balay @*/
210efca3c55SSatish Balay PetscErrorCode  PetscMallocSet(PetscErrorCode (*imalloc)(size_t,int,const char[],const char[],void**),
211efca3c55SSatish Balay                                               PetscErrorCode (*ifree)(void*,int,const char[],const char[]))
212e5c89e4eSSatish Balay {
213e5c89e4eSSatish Balay   PetscFunctionBegin;
214e32f2f54SBarry Smith   if (petscsetmallocvisited && (imalloc != PetscTrMalloc || ifree != PetscTrFree)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"cannot call multiple times");
215e5c89e4eSSatish Balay   PetscTrMalloc         = imalloc;
216e5c89e4eSSatish Balay   PetscTrFree           = ifree;
217e5c89e4eSSatish Balay   petscsetmallocvisited = PETSC_TRUE;
218e5c89e4eSSatish Balay   PetscFunctionReturn(0);
219e5c89e4eSSatish Balay }
220e5c89e4eSSatish Balay 
221e5c89e4eSSatish Balay /*@C
2221d1a0024SBarry Smith    PetscMallocClear - Resets the routines used to do mallocs and frees to the
223e5c89e4eSSatish Balay         defaults.
224e5c89e4eSSatish Balay 
225e5c89e4eSSatish Balay    Not Collective
226e5c89e4eSSatish Balay 
227e5c89e4eSSatish Balay    Level: developer
228e5c89e4eSSatish Balay 
229e5c89e4eSSatish Balay    Notes:
230e5c89e4eSSatish Balay     In general one should never run a PETSc program with different malloc() and
231e5c89e4eSSatish Balay     free() settings for different parts; this is because one NEVER wants to
232e5c89e4eSSatish Balay     free() an address that was malloced by a different memory management system
233e5c89e4eSSatish Balay 
234e5c89e4eSSatish Balay @*/
2357087cfbeSBarry Smith PetscErrorCode  PetscMallocClear(void)
236e5c89e4eSSatish Balay {
237e5c89e4eSSatish Balay   PetscFunctionBegin;
238e5c89e4eSSatish Balay   PetscTrMalloc         = PetscMallocAlign;
239e5c89e4eSSatish Balay   PetscTrFree           = PetscFreeAlign;
240e5c89e4eSSatish Balay   petscsetmallocvisited = PETSC_FALSE;
241e5c89e4eSSatish Balay   PetscFunctionReturn(0);
242e5c89e4eSSatish Balay }
243b44d5720SBarry Smith 
244b44d5720SBarry Smith PetscErrorCode PetscMemoryTrace(const char label[])
245b44d5720SBarry Smith {
246b44d5720SBarry Smith   PetscErrorCode        ierr;
247b44d5720SBarry Smith   PetscLogDouble        mem,mal;
248b44d5720SBarry Smith   static PetscLogDouble oldmem = 0,oldmal = 0;
249b44d5720SBarry Smith 
250b44d5720SBarry Smith   PetscFunctionBegin;
251b44d5720SBarry Smith   ierr = PetscMemoryGetCurrentUsage(&mem);CHKERRQ(ierr);
252b44d5720SBarry Smith   ierr = PetscMallocGetCurrentUsage(&mal);CHKERRQ(ierr);
253b44d5720SBarry Smith 
254b44d5720SBarry 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);
255b44d5720SBarry Smith   oldmem = mem;
256b44d5720SBarry Smith   oldmal = mal;
257b44d5720SBarry Smith   PetscFunctionReturn(0);
258b44d5720SBarry Smith }
25913850c04SHong Zhang 
26050a41461SHong Zhang static PetscErrorCode (*PetscTrMallocOld)(size_t,int,const char[],const char[],void**) = PetscMallocAlign;
26150a41461SHong Zhang static PetscErrorCode (*PetscTrFreeOld)(void*,int,const char[],const char[])           = PetscFreeAlign;
262*de1d6c17SHong Zhang 
263*de1d6c17SHong Zhang /*@C
264*de1d6c17SHong Zhang    PetscMallocSetDRAM - Set PetscMalloc to use DRAM.
265*de1d6c17SHong Zhang      If memkind is available, change the memkind type. Otherwise, switch the
266*de1d6c17SHong Zhang      current malloc and free routines to the PetscMallocAlign and
267*de1d6c17SHong Zhang      PetscFreeAlign (PETSc default).
268*de1d6c17SHong Zhang 
269*de1d6c17SHong Zhang    Not Collective
270*de1d6c17SHong Zhang 
271*de1d6c17SHong Zhang    Level: developer
272*de1d6c17SHong Zhang 
273*de1d6c17SHong Zhang    Notes:
274*de1d6c17SHong Zhang      This provides a way to do the allocation on DRAM temporarily. One
275*de1d6c17SHong Zhang      can switch back to the previous choice by calling PetscMallocReset().
276*de1d6c17SHong Zhang 
277*de1d6c17SHong Zhang .seealso: PetscMallocReset()
278*de1d6c17SHong Zhang @*/
27913850c04SHong Zhang PetscErrorCode PetscMallocSetDRAM(void)
28013850c04SHong Zhang {
28113850c04SHong Zhang   PetscFunctionBegin;
282*de1d6c17SHong Zhang   if (PetscTrMalloc == PetscMallocAlign) {
283*de1d6c17SHong Zhang #if defined(PETSC_HAVE_MEMKIND)
284*de1d6c17SHong Zhang     previousmktype = currentmktype;
285*de1d6c17SHong Zhang     currentmktype  = PETSC_MK_DEFAULT;
286*de1d6c17SHong Zhang #endif
287*de1d6c17SHong Zhang   } else {
28813850c04SHong Zhang     /* Save the previous choice */
28913850c04SHong Zhang     PetscTrMallocOld = PetscTrMalloc;
29013850c04SHong Zhang     PetscTrFreeOld   = PetscTrFree;
29113850c04SHong Zhang     PetscTrMalloc    = PetscMallocAlign;
29213850c04SHong Zhang     PetscTrFree      = PetscFreeAlign;
293*de1d6c17SHong Zhang   }
29413850c04SHong Zhang   PetscFunctionReturn(0);
29513850c04SHong Zhang }
29613850c04SHong Zhang 
297*de1d6c17SHong Zhang /*@C
298*de1d6c17SHong Zhang    PetscMallocResetDRAM - Reset the changes made by PetscMallocSetDRAM
299*de1d6c17SHong Zhang 
300*de1d6c17SHong Zhang    Not Collective
301*de1d6c17SHong Zhang 
302*de1d6c17SHong Zhang    Level: developer
303*de1d6c17SHong Zhang 
304*de1d6c17SHong Zhang .seealso: PetscMallocSetDRAM()
305*de1d6c17SHong Zhang @*/
30613850c04SHong Zhang PetscErrorCode PetscMallocResetDRAM(void)
30713850c04SHong Zhang {
30813850c04SHong Zhang   PetscFunctionBegin;
309*de1d6c17SHong Zhang   if (PetscTrMalloc == PetscMallocAlign) {
310*de1d6c17SHong Zhang #if defined(PETSC_HAVE_MEMKIND)
311*de1d6c17SHong Zhang     currentmktype = previousmktype;
312*de1d6c17SHong Zhang #endif
313*de1d6c17SHong Zhang   } else {
31413850c04SHong Zhang     /* Reset to the previous choice */
31513850c04SHong Zhang     PetscTrMalloc = PetscTrMallocOld;
31613850c04SHong Zhang     PetscTrFree   = PetscTrFreeOld;
317*de1d6c17SHong Zhang   }
31813850c04SHong Zhang   PetscFunctionReturn(0);
31913850c04SHong Zhang }
320