mem_mgr.txx 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572
  1. #include <omp.h>
  2. #include <algorithm>
  3. #include <cstring>
  4. #include <cassert>
  5. #include <pvfmm/profile.hpp>
  6. namespace pvfmm {
  7. #ifdef PVFMM_MEMDEBUG
  8. template <class ValueType> inline ConstIterator<ValueType>::ConstIterator(const ValueType* base_, difference_type len_, bool dynamic_alloc) {
  9. this->base = (char*)base_;
  10. this->len = len_ * (Long)sizeof(ValueType);
  11. this->offset = 0;
  12. PVFMM_ASSERT_MSG((uintptr_t)(this->base + this->offset) % alignof(ValueType) == 0, "invalid alignment during pointer type conversion.");
  13. if (dynamic_alloc) {
  14. MemoryManager::MemHead& mh = *&MemoryManager::GetMemHead((char*)this->base);
  15. MemoryManager::CheckMemHead(mh);
  16. alloc_ctr = mh.alloc_ctr;
  17. mem_head = &mh;
  18. } else
  19. mem_head = NULL;
  20. }
  21. template <class ValueType> inline void ConstIterator<ValueType>::IteratorAssertChecks(Long j) const {
  22. const auto& base = this->base;
  23. const auto& offset = this->offset + j * (Long)sizeof(ValueType);
  24. const auto& len = this->len;
  25. const auto& mem_head = this->mem_head;
  26. const auto& alloc_ctr = this->alloc_ctr;
  27. if (*this == NULL) PVFMM_WARN("dereferencing a NULL pointer is undefined.");
  28. PVFMM_ASSERT_MSG(offset >= 0 && offset + (Long)sizeof(ValueType) <= len, "access to pointer [B" << (offset < 0 ? "" : "+") << offset << ",B" << (offset + (Long)sizeof(ValueType) < 0 ? "" : "+") << offset + (Long)sizeof(ValueType) << ") is outside of the range [B,B+" << len << ").");
  29. if (mem_head) {
  30. MemoryManager::MemHead& mh = *(MemoryManager::MemHead*)(mem_head);
  31. PVFMM_ASSERT_MSG(mh.alloc_ctr == alloc_ctr, "invalid memory address or corrupted memory.");
  32. }
  33. }
  34. template <class ValueType> inline typename ConstIterator<ValueType>::reference ConstIterator<ValueType>::operator*() const {
  35. this->IteratorAssertChecks();
  36. return *(ValueType*)(base + offset);
  37. }
  38. template <class ValueType> inline const typename ConstIterator<ValueType>::value_type* ConstIterator<ValueType>::operator->() const {
  39. this->IteratorAssertChecks();
  40. return (ValueType*)(base + offset);
  41. }
  42. template <class ValueType> inline typename ConstIterator<ValueType>::reference ConstIterator<ValueType>::operator[](difference_type j) const {
  43. this->IteratorAssertChecks(j);
  44. return *(ValueType*)(base + offset + j * (Long)sizeof(ValueType));
  45. }
  46. template <class ValueType> inline typename Iterator<ValueType>::reference Iterator<ValueType>::operator*() const {
  47. this->IteratorAssertChecks();
  48. return *(ValueType*)(this->base + this->offset);
  49. }
  50. template <class ValueType> inline typename Iterator<ValueType>::value_type* Iterator<ValueType>::operator->() const {
  51. this->IteratorAssertChecks();
  52. return (ValueType*)(this->base + this->offset);
  53. }
  54. template <class ValueType> inline typename Iterator<ValueType>::reference Iterator<ValueType>::operator[](difference_type j) const {
  55. this->IteratorAssertChecks(j);
  56. return *(ValueType*)(this->base + this->offset + j * (Long)sizeof(ValueType));
  57. }
  58. template <class ValueType, Long DIM> inline StaticArray<ValueType, DIM>::StaticArray() {
  59. arr = aligned_new<ValueType>(DIM);
  60. // arr = PVFMM_PTR2ITR(ValueType, arr_, DIM);
  61. Iterator<ValueType>::operator=(arr);
  62. }
  63. template <class ValueType, Long DIM> inline StaticArray<ValueType, DIM>::~StaticArray() { aligned_delete<ValueType>(arr); }
  64. template <class ValueType, Long DIM> inline StaticArray<ValueType, DIM>::StaticArray(const StaticArray& I) {
  65. arr = aligned_new<ValueType>(DIM);
  66. // arr = PVFMM_PTR2ITR(ValueType, arr_, DIM);
  67. Iterator<ValueType>::operator=(arr);
  68. for (Long i = 0; i < DIM; i++) (*this)[i] = I[i];
  69. }
  70. template <class ValueType, Long DIM> inline StaticArray<ValueType, DIM>& StaticArray<ValueType, DIM>::operator=(const StaticArray& I) {
  71. for (Long i = 0; i < DIM; i++) (*this)[i] = I[i];
  72. return *this;
  73. }
  74. #endif
  75. template <class T> inline uintptr_t TypeTraits<T>::ID() { return (uintptr_t) & ID; }
  76. template <class T> inline bool TypeTraits<T>::IsPOD() { return false; }
  77. #define PVFMMDefinePOD(type) \
  78. template <> bool inline TypeTraits<type>::IsPOD() { \
  79. return true; \
  80. };
  81. PVFMMDefinePOD(char);
  82. PVFMMDefinePOD(float);
  83. PVFMMDefinePOD(double);
  84. PVFMMDefinePOD(int);
  85. PVFMMDefinePOD(long long);
  86. PVFMMDefinePOD(unsigned long);
  87. PVFMMDefinePOD(char*);
  88. PVFMMDefinePOD(float*);
  89. PVFMMDefinePOD(double*);
  90. #undef PVFMMDefinePOD
  91. inline MemoryManager::MemoryManager(Long N) {
  92. buff_size = N;
  93. { // Allocate buff
  94. assert(PVFMM_MEM_ALIGN <= 0x8000);
  95. Long alignment = PVFMM_MEM_ALIGN - 1;
  96. char* base_ptr = (char*)::malloc(N + 2 + alignment);
  97. PVFMM_ASSERT_MSG(base_ptr, "memory allocation failed.");
  98. buff = (char*)((uintptr_t)(base_ptr + 2 + alignment) & ~(uintptr_t)alignment);
  99. ((uint16_t*)buff)[-1] = (uint16_t)(buff - base_ptr);
  100. }
  101. { // Initialize to init_mem_val
  102. #ifdef PVFMM_MEMDEBUG
  103. #pragma omp parallel for
  104. for (Long i = 0; i < buff_size; i++) {
  105. buff[i] = init_mem_val;
  106. }
  107. #endif
  108. }
  109. n_dummy_indx = new_node();
  110. Long n_indx = new_node();
  111. MemNode& n_dummy = node_buff[n_dummy_indx - 1];
  112. MemNode& n = node_buff[n_indx - 1];
  113. n_dummy.size = 0;
  114. n_dummy.free = false;
  115. n_dummy.prev = 0;
  116. n_dummy.next = n_indx;
  117. n_dummy.mem_ptr = &buff[0];
  118. assert(n_indx);
  119. n.size = N;
  120. n.free = true;
  121. n.prev = n_dummy_indx;
  122. n.next = 0;
  123. n.mem_ptr = &buff[0];
  124. n.it = free_map.insert(std::make_pair(N, n_indx));
  125. omp_init_lock(&omp_lock);
  126. }
  127. inline MemoryManager::~MemoryManager() {
  128. Check();
  129. MemNode* n_dummy = &node_buff[n_dummy_indx - 1];
  130. MemNode* n = &node_buff[n_dummy->next - 1];
  131. if (!n->free || n->size != buff_size || node_stack.size() != node_buff.size() - 2 || !system_malloc.empty()) {
  132. PVFMM_WARN("memory leak detected.");
  133. }
  134. omp_destroy_lock(&omp_lock);
  135. { // free buff
  136. assert(buff);
  137. ::free(buff - ((uint16_t*)buff)[-1]);
  138. }
  139. }
  140. inline MemoryManager::MemHead& MemoryManager::GetMemHead(char* I) {
  141. PVFMM_ASSERT_MSG(I != NULL, "NULL pointer exception.");
  142. static uintptr_t alignment = PVFMM_MEM_ALIGN - 1;
  143. static uintptr_t header_size = (uintptr_t)(sizeof(MemHead) + alignment) & ~(uintptr_t)alignment;
  144. return *(MemHead*)(((char*)I) - header_size);
  145. }
  146. inline void MemoryManager::CheckMemHead(const MemHead& mem_head) { // Verify header check_sum
  147. #ifdef PVFMM_MEMDEBUG
  148. Long check_sum = 0;
  149. const unsigned char* base_ = (const unsigned char*)&mem_head;
  150. for (Integer i = 0; i < sizeof(MemHead); i++) {
  151. check_sum += base_[i];
  152. }
  153. check_sum -= mem_head.check_sum;
  154. check_sum = check_sum & ((1UL << (8 * sizeof(mem_head.check_sum))) - 1);
  155. PVFMM_ASSERT_MSG(check_sum == mem_head.check_sum, "invalid memory address or corrupted memory.");
  156. #endif
  157. }
  158. inline Iterator<char> MemoryManager::malloc(const Long n_elem, const Long type_size, const uintptr_t type_id) const {
  159. if (!n_elem) return NULL;
  160. static uintptr_t alignment = PVFMM_MEM_ALIGN - 1;
  161. static uintptr_t header_size = (uintptr_t)(sizeof(MemHead) + alignment) & ~(uintptr_t)alignment;
  162. Long size = n_elem * type_size + header_size;
  163. size = (uintptr_t)(size + alignment) & ~(uintptr_t)alignment;
  164. char* base = NULL;
  165. omp_set_lock(&omp_lock);
  166. static Long alloc_ctr = 0;
  167. alloc_ctr++;
  168. Long head_alloc_ctr = alloc_ctr;
  169. std::multimap<Long, Long>::iterator it = free_map.lower_bound(size);
  170. Long n_indx = (it != free_map.end() ? it->second : 0);
  171. if (n_indx) { // Allocate from buff
  172. Long n_free_indx = (it->first > size ? new_node() : 0);
  173. MemNode& n = node_buff[n_indx - 1];
  174. assert(n.size == it->first);
  175. assert(n.it == it);
  176. assert(n.free);
  177. if (n_free_indx) { // Create a node for the remaining free part.
  178. MemNode& n_free = node_buff[n_free_indx - 1];
  179. n_free = n;
  180. n_free.size -= size;
  181. n_free.mem_ptr = (char*)n_free.mem_ptr + size;
  182. { // Insert n_free to the link list
  183. n_free.prev = n_indx;
  184. if (n_free.next) {
  185. Long n_next_indx = n_free.next;
  186. MemNode& n_next = node_buff[n_next_indx - 1];
  187. n_next.prev = n_free_indx;
  188. }
  189. n.next = n_free_indx;
  190. }
  191. assert(n_free.free); // Insert n_free to free map
  192. n_free.it = free_map.insert(std::make_pair(n_free.size, n_free_indx));
  193. n.size = size; // Update n
  194. }
  195. n.free = false;
  196. free_map.erase(it);
  197. base = n.mem_ptr;
  198. }
  199. omp_unset_lock(&omp_lock);
  200. if (!base) { // Use system malloc
  201. Long end_padding = 8; // to check for out-of-bound writes
  202. char* p = (char*)::malloc(size + 2 + alignment + end_padding);
  203. PVFMM_ASSERT_MSG(p, "memory allocation failed.");
  204. #ifdef PVFMM_MEMDEBUG
  205. { // system_malloc.insert(p)
  206. omp_set_lock(&omp_lock);
  207. system_malloc.insert(p);
  208. omp_unset_lock(&omp_lock);
  209. }
  210. { // set p[*] to init_mem_val
  211. #pragma omp parallel for
  212. for (Long i = 0; i < size + 2 + alignment + end_padding; i++) p[i] = init_mem_val;
  213. }
  214. #endif
  215. { // base <-- align(p)
  216. base = (char*)((uintptr_t)(p + 2 + alignment) & ~(uintptr_t)alignment);
  217. ((uint16_t*)base)[-1] = (uint16_t)(base - p);
  218. }
  219. }
  220. { // Check out-of-bounds write
  221. #ifdef PVFMM_MEMDEBUG
  222. if (n_indx) {
  223. #pragma omp parallel for
  224. for (Long i = 0; i < size; i++) PVFMM_ASSERT_MSG(base[i] == init_mem_val, "memory corruption detected.");
  225. }
  226. #endif
  227. }
  228. MemHead& mem_head = *(MemHead*)base;
  229. { // Set mem_head
  230. #ifdef PVFMM_MEMDEBUG
  231. for (Integer i = 0; i < sizeof(MemHead); i++) base[i] = init_mem_val;
  232. #endif
  233. mem_head.n_indx = n_indx;
  234. mem_head.n_elem = n_elem;
  235. mem_head.type_size = type_size;
  236. mem_head.alloc_ctr = head_alloc_ctr;
  237. mem_head.type_id = type_id;
  238. }
  239. { // Set header check_sum
  240. #ifdef PVFMM_MEMDEBUG
  241. Long check_sum = 0;
  242. unsigned char* base_ = (unsigned char*)base;
  243. mem_head.check_sum = 0;
  244. for (Integer i = 0; i < sizeof(MemHead); i++) check_sum += base_[i];
  245. check_sum = check_sum & ((1UL << (8 * sizeof(mem_head.check_sum))) - 1);
  246. mem_head.check_sum = check_sum;
  247. #endif
  248. }
  249. Profile::Add_MEM(n_elem * type_size);
  250. #ifdef PVFMM_MEMDEBUG
  251. return Iterator<char>(base + header_size, n_elem * type_size, true);
  252. #else
  253. return base + header_size;
  254. #endif
  255. }
  256. inline void MemoryManager::free(Iterator<char> p) const {
  257. if (p == NULL) return;
  258. static uintptr_t alignment = PVFMM_MEM_ALIGN - 1;
  259. static uintptr_t header_size = (uintptr_t)(sizeof(MemHead) + alignment) & ~(uintptr_t)alignment;
  260. MemHead& mem_head = GetMemHead(&p[0]);
  261. Long n_indx = mem_head.n_indx;
  262. Long n_elem = mem_head.n_elem;
  263. Long type_size = mem_head.type_size;
  264. char* base = (char*)&mem_head;
  265. { // Verify header check_sum; set array to init_mem_val
  266. #ifdef PVFMM_MEMDEBUG
  267. CheckMemHead(mem_head);
  268. Long size = mem_head.n_elem * mem_head.type_size;
  269. #pragma omp parallel for
  270. for (Long i = 0; i < size; i++) p[i] = init_mem_val;
  271. for (Integer i = 0; i < sizeof(MemHead); i++) base[i] = init_mem_val;
  272. #endif
  273. }
  274. if (n_indx == 0) { // Use system free
  275. assert(base < &buff[0] || base >= &buff[buff_size]);
  276. char* p_;
  277. { // p_ <-- unalign(base)
  278. p_ = (char*)((uintptr_t)base - ((uint16_t*)base)[-1]);
  279. }
  280. #ifdef PVFMM_MEMDEBUG
  281. { // Check out-of-bounds write
  282. base[-1] = init_mem_val;
  283. base[-2] = init_mem_val;
  284. Long size = n_elem * type_size + header_size;
  285. size = (uintptr_t)(size + alignment) & ~(uintptr_t)alignment;
  286. Long end_padding = 8; // to check for out-of-bound writes
  287. #pragma omp parallel for
  288. for (Long i = 0; i < size + 2 + alignment + end_padding; i++) {
  289. PVFMM_ASSERT_MSG(p_[i] == init_mem_val, "memory corruption detected.");
  290. }
  291. }
  292. { // system_malloc.erase(p_)
  293. omp_set_lock(&omp_lock);
  294. PVFMM_ASSERT_MSG(system_malloc.erase(p_) == 1, "double free or corruption.");
  295. omp_unset_lock(&omp_lock);
  296. }
  297. #endif
  298. ::free(p_);
  299. } else {
  300. assert(n_indx <= node_buff.size());
  301. omp_set_lock(&omp_lock);
  302. MemNode& n = node_buff[n_indx - 1];
  303. assert(!n.free && n.size > 0 && n.mem_ptr == base);
  304. if (n.prev != 0 && node_buff[n.prev - 1].free) {
  305. Long n_prev_indx = n.prev;
  306. MemNode& n_prev = node_buff[n_prev_indx - 1];
  307. n.size += n_prev.size;
  308. n.mem_ptr = n_prev.mem_ptr;
  309. n.prev = n_prev.prev;
  310. free_map.erase(n_prev.it);
  311. delete_node(n_prev_indx);
  312. if (n.prev) {
  313. node_buff[n.prev - 1].next = n_indx;
  314. }
  315. }
  316. if (n.next != 0 && node_buff[n.next - 1].free) {
  317. Long n_next_indx = n.next;
  318. MemNode& n_next = node_buff[n_next_indx - 1];
  319. n.size += n_next.size;
  320. n.next = n_next.next;
  321. free_map.erase(n_next.it);
  322. delete_node(n_next_indx);
  323. if (n.next) {
  324. node_buff[n.next - 1].prev = n_indx;
  325. }
  326. }
  327. n.free = true; // Insert n to free_map
  328. n.it = free_map.insert(std::make_pair(n.size, n_indx));
  329. omp_unset_lock(&omp_lock);
  330. }
  331. Profile::Add_MEM(-n_elem * type_size);
  332. }
  333. inline void MemoryManager::print() const {
  334. if (!buff_size) return;
  335. omp_set_lock(&omp_lock);
  336. Long size = 0;
  337. Long largest_size = 0;
  338. MemNode* n = &node_buff[n_dummy_indx - 1];
  339. std::cout << "\n|";
  340. while (n->next) {
  341. n = &node_buff[n->next - 1];
  342. if (n->free) {
  343. std::cout << ' ';
  344. largest_size = std::max(largest_size, n->size);
  345. } else {
  346. std::cout << '#';
  347. size += n->size;
  348. }
  349. }
  350. std::cout << "| allocated=" << round(size * 1000.0 / buff_size) / 10 << "%";
  351. std::cout << " largest_free=" << round(largest_size * 1000.0 / buff_size) / 10 << "%\n";
  352. omp_unset_lock(&omp_lock);
  353. }
  354. inline void MemoryManager::test() {
  355. Long M = 2000000000;
  356. { // With memory manager
  357. Long N = M * sizeof(double) * 1.1;
  358. double tt;
  359. Iterator<double> tmp;
  360. std::cout << "With memory manager: ";
  361. MemoryManager memgr(N);
  362. for (Integer j = 0; j < 3; j++) {
  363. tmp = (Iterator<double>)memgr.malloc(M * sizeof(double));
  364. PVFMM_ASSERT(tmp != NULL);
  365. tt = omp_get_wtime();
  366. #pragma omp parallel for
  367. for (Long i = 0; i < M; i += 64) tmp[i] = (double)i;
  368. tt = omp_get_wtime() - tt;
  369. std::cout << tt << ' ';
  370. memgr.free((Iterator<char>)tmp);
  371. }
  372. std::cout << '\n';
  373. }
  374. { // Without memory manager
  375. double tt;
  376. double* tmp;
  377. std::cout << "Without memory manager: ";
  378. for (Integer j = 0; j < 3; j++) {
  379. tmp = (double*)::malloc(M * sizeof(double));
  380. PVFMM_ASSERT(tmp != NULL);
  381. tt = omp_get_wtime();
  382. #pragma omp parallel for
  383. for (Long i = 0; i < M; i += 64) tmp[i] = (double)i;
  384. tt = omp_get_wtime() - tt;
  385. std::cout << tt << ' ';
  386. ::free(tmp);
  387. }
  388. std::cout << '\n';
  389. }
  390. }
  391. inline void MemoryManager::Check() const {
  392. #ifdef PVFMM_MEMDEBUG
  393. // print();
  394. omp_set_lock(&omp_lock);
  395. MemNode* curr_node = &node_buff[n_dummy_indx - 1];
  396. while (curr_node->next) {
  397. curr_node = &node_buff[curr_node->next - 1];
  398. if (curr_node->free) {
  399. char* base = curr_node->mem_ptr;
  400. #pragma omp parallel for
  401. for (Long i = 0; i < curr_node->size; i++) {
  402. PVFMM_ASSERT_MSG(base[i] == init_mem_val, "memory corruption detected.");
  403. }
  404. }
  405. }
  406. omp_unset_lock(&omp_lock);
  407. #endif
  408. }
  409. inline Long MemoryManager::new_node() const {
  410. if (node_stack.empty()) {
  411. node_buff.resize(node_buff.size() + 1);
  412. node_stack.push(node_buff.size());
  413. }
  414. Long indx = node_stack.top();
  415. node_stack.pop();
  416. assert(indx);
  417. return indx;
  418. }
  419. inline void MemoryManager::delete_node(Long indx) const {
  420. assert(indx);
  421. assert(indx <= node_buff.size());
  422. MemNode& n = node_buff[indx - 1];
  423. n.free = false;
  424. n.size = 0;
  425. n.prev = 0;
  426. n.next = 0;
  427. n.mem_ptr = NULL;
  428. node_stack.push(indx);
  429. }
  430. template <class ValueType> inline Iterator<ValueType> aligned_new(Long n_elem, const MemoryManager* mem_mgr) {
  431. if (!n_elem) return NULL;
  432. static MemoryManager def_mem_mgr(0);
  433. if (!mem_mgr) mem_mgr = &def_mem_mgr;
  434. Iterator<ValueType> A = (Iterator<ValueType>)mem_mgr->malloc(n_elem, sizeof(ValueType));
  435. PVFMM_ASSERT_MSG(A != NULL, "memory allocation failed.");
  436. if (!TypeTraits<ValueType>::IsPOD()) { // Call constructors
  437. // printf("%s\n", __PRETTY_FUNCTION__);
  438. #pragma omp parallel for schedule(static)
  439. for (Long i = 0; i < n_elem; i++) {
  440. ValueType* Ai = new (&A[i]) ValueType();
  441. assert(Ai == (&A[i]));
  442. }
  443. } else {
  444. #ifdef PVFMM_MEMDEBUG
  445. static Long random_init_val = 1;
  446. Iterator<char> A_ = (Iterator<char>)A;
  447. #pragma omp parallel for schedule(static)
  448. for (Long i = 0; i < n_elem * sizeof(ValueType); i++) {
  449. A_[i] = random_init_val + i;
  450. }
  451. random_init_val += n_elem * sizeof(ValueType);
  452. #endif
  453. }
  454. return A;
  455. }
  456. template <class ValueType> inline void aligned_delete(Iterator<ValueType> A, const MemoryManager* mem_mgr) {
  457. if (A == NULL) return;
  458. if (!TypeTraits<ValueType>::IsPOD()) { // Call destructors
  459. // printf("%s\n", __PRETTY_FUNCTION__);
  460. MemoryManager::MemHead& mem_head = MemoryManager::GetMemHead((char*)&A[0]);
  461. #ifdef PVFMM_MEMDEBUG
  462. MemoryManager::CheckMemHead(mem_head);
  463. // PVFMM_ASSERT_MSG(mem_head.n_elem==1 || mem_head.type_id==TypeTraits<ValueType>::ID(), "pointer to aligned_delete has different type than what was used in aligned_new.");
  464. #endif
  465. Long n_elem = mem_head.n_elem;
  466. for (Long i = 0; i < n_elem; i++) {
  467. A[i].~ValueType();
  468. }
  469. } else {
  470. #ifdef PVFMM_MEMDEBUG
  471. MemoryManager::MemHead& mem_head = MemoryManager::GetMemHead((char*)&A[0]);
  472. MemoryManager::CheckMemHead(mem_head);
  473. // PVFMM_ASSERT_MSG(mem_head.type_id==TypeTraits<ValueType>::ID(), "pointer to aligned_delete has different type than what was used in aligned_new.");
  474. Long size = mem_head.n_elem * mem_head.type_size;
  475. Iterator<char> A_ = (Iterator<char>)A;
  476. #pragma omp parallel for
  477. for (Long i = 0; i < size; i++) {
  478. A_[i] = 0;
  479. }
  480. #endif
  481. }
  482. static MemoryManager def_mem_mgr(0);
  483. if (!mem_mgr) mem_mgr = &def_mem_mgr;
  484. mem_mgr->free((Iterator<char>)A);
  485. }
  486. template <class ValueType> inline Iterator<ValueType> memcopy(Iterator<ValueType> destination, ConstIterator<ValueType> source, Long num) {
  487. if (destination != source && num) {
  488. #ifdef PVFMM_MEMDEBUG
  489. destination[num - 1];
  490. source[num - 1];
  491. #endif
  492. if (TypeTraits<ValueType>::IsPOD()) {
  493. memcpy(&destination[0], &source[0], num * sizeof(ValueType));
  494. } else {
  495. for (Long i = 0; i < num; i++) destination[i] = source[i];
  496. }
  497. }
  498. return destination;
  499. }
  500. template <class ValueType> inline Iterator<ValueType> memset(Iterator<ValueType> ptr, int value, Long num) {
  501. if (num) {
  502. #ifdef PVFMM_MEMDEBUG
  503. ptr[0];
  504. ptr[num - 1];
  505. #endif
  506. ::memset(&ptr[0], value, num * sizeof(ValueType));
  507. }
  508. return ptr;
  509. }
  510. } // end namespace