mem_mgr.txx 7.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267
  1. /**
  2. * \file mem_mgr.txx
  3. * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
  4. * \date 9-21-2014
  5. * \brief This file contains the definition of a simple memory manager which
  6. * uses a pre-allocated buffer of size defined in call to the constructor.
  7. */
  8. #include <omp.h>
  9. #include <algorithm>
  10. #include <cstring>
  11. #include <cassert>
  12. namespace pvfmm{
  13. namespace mem{
  14. template <class T>
  15. uintptr_t TypeTraits<T>::ID(){
  16. return (uintptr_t)&ID;
  17. }
  18. template <class T>
  19. bool TypeTraits<T>::IsPOD(){
  20. return false;
  21. }
  22. MemoryManager::MemHead* MemoryManager::GetMemHead(void* p){
  23. static uintptr_t alignment=MEM_ALIGN-1;
  24. static uintptr_t header_size=(uintptr_t)(sizeof(MemoryManager::MemHead)+alignment) & ~(uintptr_t)alignment;
  25. return (MemHead*)(((char*)p)-header_size);
  26. }
  27. void* MemoryManager::malloc(const size_t& n_elem, const size_t& type_size, const uintptr_t& type_id) const{
  28. if(!n_elem) return NULL;
  29. static uintptr_t alignment=MEM_ALIGN-1;
  30. static uintptr_t header_size=(uintptr_t)(sizeof(MemHead)+alignment) & ~(uintptr_t)alignment;
  31. size_t size=n_elem*type_size+header_size;
  32. size=(uintptr_t)(size+alignment) & ~(uintptr_t)alignment;
  33. char* base=NULL;
  34. omp_set_lock(&omp_lock);
  35. std::multimap<size_t, size_t>::iterator it=free_map.lower_bound(size);
  36. size_t n_indx=(it!=free_map.end()?it->second:0);
  37. if(n_indx){ // Allocate from buff
  38. size_t n_free_indx=(it->first>size?new_node():0);
  39. MemNode& n=node_buff[n_indx-1];
  40. assert(n.size==it->first);
  41. assert(n.it==it);
  42. assert(n.free);
  43. if(n_free_indx){ // Create a node for the remaining free part.
  44. MemNode& n_free=node_buff[n_free_indx-1];
  45. n_free=n;
  46. n_free.size-=size;
  47. n_free.mem_ptr=(char*)n_free.mem_ptr+size;
  48. { // Insert n_free to the link list
  49. n_free.prev=n_indx;
  50. if(n_free.next){
  51. size_t n_next_indx=n_free.next;
  52. MemNode& n_next=node_buff[n_next_indx-1];
  53. n_next.prev=n_free_indx;
  54. }
  55. n.next=n_free_indx;
  56. }
  57. assert(n_free.free); // Insert n_free to free map
  58. n_free.it=free_map.insert(std::make_pair(n_free.size,n_free_indx));
  59. n.size=size; // Update n
  60. }
  61. n.free=false;
  62. free_map.erase(it);
  63. base = n.mem_ptr;
  64. }
  65. omp_unset_lock(&omp_lock);
  66. if(!base){ // Use system malloc
  67. size+=2+alignment;
  68. char* p = (char*)::malloc(size);
  69. base = (char*)((uintptr_t)(p+2+alignment) & ~(uintptr_t)alignment);
  70. ((uint16_t*)base)[-1] = (uint16_t)(base-p);
  71. }
  72. { // Check out-of-bounds write
  73. #ifndef NDEBUG
  74. if(n_indx){
  75. #pragma omp parallel for
  76. for(size_t i=0;i<size;i++) assert(base[i]==init_mem_val);
  77. }
  78. #endif
  79. }
  80. MemHead* mem_head=(MemHead*)base;
  81. { // Set mem_head
  82. mem_head->n_indx=n_indx;
  83. mem_head->n_elem=n_elem;
  84. mem_head->type_id=type_id;
  85. }
  86. { // Set header check_sum
  87. #ifndef NDEBUG
  88. size_t check_sum=0;
  89. mem_head->check_sum=0;
  90. for(size_t i=0;i<header_size;i++){
  91. check_sum+=base[i];
  92. }
  93. check_sum=check_sum & ((1UL << sizeof(mem_head->check_sum))-1);
  94. mem_head->check_sum=check_sum;
  95. #endif
  96. }
  97. return (void*)(base+header_size);
  98. }
  99. void MemoryManager::free(void* p, const size_t& type_size, const uintptr_t& type_id) const{
  100. if(!p) return;
  101. static uintptr_t alignment=MEM_ALIGN-1;
  102. static uintptr_t header_size=(uintptr_t)(sizeof(MemHead)+alignment) & ~(uintptr_t)alignment;
  103. char* base=(char*)((char*)p-header_size);
  104. MemHead* mem_head=(MemHead*)base;
  105. assert(mem_head->type_id==type_id);
  106. if(base<&buff[0] || base>=&buff[buff_size]){ // Use system free
  107. char* p=(char*)((uintptr_t)base-((uint16_t*)base)[-1]);
  108. return ::free(p);
  109. }
  110. size_t n_indx=mem_head->n_indx;
  111. assert(n_indx>0 && n_indx<=node_buff.size());
  112. { // Verify header check_sum; set array to init_mem_val
  113. #ifndef NDEBUG
  114. { // Verify header check_sum
  115. size_t check_sum=0;
  116. for(size_t i=0;i<header_size;i++){
  117. check_sum+=base[i];
  118. }
  119. check_sum-=mem_head->check_sum;
  120. check_sum=check_sum & ((1UL << sizeof(mem_head->check_sum))-1);
  121. assert(check_sum==mem_head->check_sum);
  122. }
  123. size_t size=mem_head->n_elem*type_size;
  124. #pragma omp parallel for
  125. for(size_t i=0;i<size;i++) ((char*)p)[i]=init_mem_val;
  126. for(size_t i=0;i<sizeof(MemHead);i++) base[i]=init_mem_val;
  127. #endif
  128. }
  129. omp_set_lock(&omp_lock);
  130. MemNode& n=node_buff[n_indx-1];
  131. assert(!n.free && n.size>0 && n.mem_ptr==base);
  132. if(n.prev!=0 && node_buff[n.prev-1].free){
  133. size_t n_prev_indx=n.prev;
  134. MemNode& n_prev=node_buff[n_prev_indx-1];
  135. n.size+=n_prev.size;
  136. n.mem_ptr=n_prev.mem_ptr;
  137. n.prev=n_prev.prev;
  138. free_map.erase(n_prev.it);
  139. delete_node(n_prev_indx);
  140. if(n.prev){
  141. size_t n_prev_indx=n.prev;
  142. MemNode& n_prev=node_buff[n_prev_indx-1];
  143. n_prev.next=n_indx;
  144. }
  145. }
  146. if(n.next!=0 && node_buff[n.next-1].free){
  147. size_t n_next_indx=n.next;
  148. MemNode& n_next=node_buff[n_next_indx-1];
  149. n.size+=n_next.size;
  150. n.next=n_next.next;
  151. free_map.erase(n_next.it);
  152. delete_node(n_next_indx);
  153. if(n.next){
  154. size_t n_next_indx=n.next;
  155. MemNode& n_next=node_buff[n_next_indx-1];
  156. n_next.prev=n_indx;
  157. }
  158. }
  159. n.free=true; // Insert n to free_map
  160. n.it=free_map.insert(std::make_pair(n.size,n_indx));
  161. omp_unset_lock(&omp_lock);
  162. }
  163. size_t MemoryManager::new_node() const{
  164. if(node_stack.empty()){
  165. node_buff.resize(node_buff.size()+1);
  166. node_stack.push(node_buff.size());
  167. }
  168. size_t indx=node_stack.top();
  169. node_stack.pop();
  170. assert(indx);
  171. return indx;
  172. }
  173. void MemoryManager::delete_node(size_t indx) const{
  174. assert(indx);
  175. assert(indx<=node_buff.size());
  176. MemNode& n=node_buff[indx-1];
  177. n.free=false;
  178. n.size=0;
  179. n.prev=0;
  180. n.next=0;
  181. n.mem_ptr=NULL;
  182. node_stack.push(indx);
  183. }
  184. template <class T>
  185. T* aligned_new(size_t n_elem, const MemoryManager* mem_mgr){
  186. if(!n_elem) return NULL;
  187. static MemoryManager def_mem_mgr(0);
  188. if(!mem_mgr) mem_mgr=&def_mem_mgr;
  189. T* A=(T*)mem_mgr->malloc(n_elem, sizeof(T), TypeTraits<T>::ID());
  190. if(!TypeTraits<T>::IsPOD()){ // Call constructors
  191. //printf("%s\n", __PRETTY_FUNCTION__);
  192. #pragma omp parallel for
  193. for(size_t i=0;i<n_elem;i++){
  194. T* Ai=new(A+i) T();
  195. assert(Ai==(A+i));
  196. }
  197. }else{
  198. for(size_t i=0;i<n_elem*sizeof(T);i++){
  199. ((char*)A)[i]=0;
  200. }
  201. }
  202. assert(A);
  203. return A;
  204. }
  205. template <class T>
  206. void aligned_delete(T* A, const MemoryManager* mem_mgr){
  207. if (!A) return;
  208. if(!TypeTraits<T>::IsPOD()){ // Call destructors
  209. //printf("%s\n", __PRETTY_FUNCTION__);
  210. MemoryManager::MemHead* mem_head=MemoryManager::GetMemHead(A);
  211. size_t n_elem=mem_head->n_elem;
  212. for(size_t i=0;i<n_elem;i++){
  213. (A+i)->~T();
  214. }
  215. }else{
  216. MemoryManager::MemHead* mem_head=MemoryManager::GetMemHead(A);
  217. size_t n_elem=mem_head->n_elem;
  218. for(size_t i=0;i<n_elem*sizeof(T);i++){
  219. ((char*)A)[i]=0;
  220. }
  221. }
  222. static MemoryManager def_mem_mgr(0);
  223. if(!mem_mgr) mem_mgr=&def_mem_mgr;
  224. mem_mgr->free(A, sizeof(T), TypeTraits<T>::ID());
  225. }
  226. void* memcopy( void * destination, const void * source, size_t num){
  227. if(destination==source || num==0) return destination;
  228. return memcpy ( destination, source, num );
  229. }
  230. }//end namespace
  231. }//end namespace