mem_mgr.txx 7.5 KB

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