Modularize the Isomalloc stress test
[charm.git] / tests / ampi / isomalloc / test.C
blobf25e3914a3f88381987135b4b8e34c987363d086
1 #ifndef __STDC_FORMAT_MACROS
2 # define __STDC_FORMAT_MACROS
3 #endif
4 #ifndef __STDC_LIMIT_MACROS
5 # define __STDC_LIMIT_MACROS
6 #endif
8 #include <cstdint>
9 #include <cstdlib>
10 #include <cstdio>
11 #include <random>
12 #include <chrono>
13 #include <algorithm>
14 #include <type_traits>
16 // for getrusage
17 #include <sys/time.h>
18 #include <sys/resource.h>
20 #include "mpi.h"
22 #define STRINGIZE(x) #x
23 #define STRINGIZE_VALUE_OF(x) STRINGIZE(x)
24 #define isomalloc_method_str STRINGIZE_VALUE_OF(isomalloc_method)
26 static const uint32_t primes[] =
28   2, 3, 5, 7, 11, 13, 17, 19, 23, 29, 31, 37, 41, 43, 47, 53, 59, 61, 67, 71,
29   73, 79, 83, 89, 97, 101, 103, 107, 109, 113, 127, 131, 137, 139, 149, 151,
30   157, 163, 167, 173, 179, 181, 191, 193, 197, 199, 211, 223, 227, 229, 233,
31   239, 241, 251, 257, 263, 269, 271, 277, 281, 283, 293, 307, 311, 313, 317,
32   331, 337, 347, 349, 353, 359, 367, 373, 379, 383, 389, 397, 401, 409, 419,
33   421, 431, 433, 439, 443, 449, 457, 461, 463, 467, 479, 487, 491, 499, 503,
34   509, 521, 523, 541, 547, 557, 563, 569, 571, 577, 587, 593, 599, 601, 607,
35   613, 617, 619, 631, 641, 643, 647, 653, 659, 661, 673, 677, 683, 691, 701,
36   709, 719, 727, 733, 739, 743, 751, 757, 761, 769, 773, 787, 797, 809, 811,
37   821, 823, 827, 829, 839, 853, 857, 859, 863, 877, 881, 883, 887, 907, 911,
38   919, 929, 937, 941, 947, 953, 967, 971, 977, 983, 991, 997, 1009, 1013, 1019,
39   1021, 1031, 1033, 1039, 1049, 1051, 1061, 1063, 1069, 1087, 1091, 1093, 1097,
40   1103, 1109, 1117, 1123, 1129, 1151, 1153, 1163, 1171, 1181, 1187, 1193, 1201,
41   1213, 1217, 1223, 1229, 1231, 1237, 1249, 1259, 1277, 1279, 1283, 1289, 1291,
42   1297, 1301, 1303, 1307, 1319, 1321, 1327, 1361, 1367, 1373, 1381, 1399, 1409,
43   1423, 1427, 1429, 1433, 1439, 1447, 1451, 1453, 1459, 1471, 1481, 1483, 1487,
44   1489, 1493, 1499, 1511, 1523, 1531, 1543, 1549, 1553, 1559, 1567, 1571, 1579,
45   1583, 1597, 1601, 1607, 1609, 1613, 1619, 1621, 1627, 1637, 1657, 1663, 1667,
46   1669, 1693, 1697, 1699, 1709, 1721, 1723, 1733, 1741, 1747, 1753, 1759, 1777,
47   1783, 1787, 1789, 1801, 1811, 1823, 1831, 1847, 1861, 1867, 1871, 1873, 1877,
48   1879, 1889, 1901, 1907, 1913, 1931, 1933, 1949, 1951, 1973, 1979, 1987, 1993,
49   1997, 1999, 2003, 2011, 2017, 2027, 2029, 2039, 2053, 2063, 2069, 2081, 2083,
50   2087, 2089, 2099, 2111, 2113, 2129, 2131, 2137, 2141, 2143, 2153, 2161, 2179,
51   2203, 2207, 2213, 2221, 2237, 2239, 2243, 2251, 2267, 2269, 2273, 2281, 2287,
52   2293, 2297, 2309, 2311, 2333, 2339, 2341, 2347, 2351, 2357, 2371, 2377, 2381,
53   2383, 2389, 2393, 2399, 2411, 2417, 2423, 2437, 2441, 2447, 2459, 2467, 2473,
54   2477, 2503, 2521, 2531, 2539, 2543, 2549, 2551, 2557, 2579, 2591, 2593, 2609,
55   2617, 2621, 2633, 2647, 2657, 2659, 2663, 2671, 2677, 2683, 2687, 2689, 2693,
56   2699, 2707, 2711, 2713, 2719, 2729, 2731, 2741, 2749, 2753, 2767, 2777, 2789,
57   2791, 2797, 2801, 2803, 2819, 2833, 2837, 2843, 2851, 2857, 2861, 2879, 2887,
58   2897, 2903, 2909, 2917, 2927, 2939, 2953, 2957, 2963, 2969, 2971, 2999, 3001,
59   3011, 3019, 3023, 3037, 3041, 3049, 3061, 3067, 3079, 3083, 3089, 3109, 3119,
60   3121, 3137, 3163, 3167, 3169, 3181, 3187, 3191, 3203, 3209, 3217, 3221, 3229,
61   3251, 3253, 3257, 3259, 3271, 3299, 3301, 3307, 3313, 3319, 3323, 3329, 3331,
62   3343, 3347, 3359, 3361, 3371, 3373, 3389, 3391, 3407, 3413, 3433, 3449, 3457,
63   3461, 3463, 3467, 3469, 3491, 3499, 3511, 3517, 3527, 3529, 3533, 3539, 3541,
64   3547, 3557, 3559, 3571, 3581, 3583, 3593, 3607, 3613, 3617, 3623, 3631, 3637,
65   3643, 3659, 3671, 3673, 3677, 3691, 3697, 3701, 3709, 3719, 3727, 3733, 3739,
66   3761, 3767, 3769, 3779, 3793, 3797, 3803, 3821, 3823, 3833, 3847, 3851, 3853,
67   3863, 3877, 3881, 3889, 3907, 3911, 3917, 3919, 3923, 3929, 3931, 3943, 3947,
68   3967, 3989, 4001, 4003, 4007, 4013, 4019, 4021, 4027, 4049, 4051, 4057, 4073,
69   4079, 4091, 4093, 4099, 4111, 4127, 4129, 4133, 4139, 4153, 4157, 4159, 4177,
70   4201, 4211, 4217, 4219, 4229, 4231, 4241, 4243, 4253, 4259, 4261, 4271, 4273,
71   4283, 4289, 4297, 4327, 4337, 4339, 4349, 4357, 4363, 4373, 4391, 4397, 4409,
72   4421, 4423, 4441, 4447, 4451, 4457, 4463, 4481, 4483, 4493, 4507, 4513, 4517,
73   4519, 4523, 4547, 4549, 4561, 4567, 4583, 4591, 4597, 4603, 4621, 4637, 4639,
74   4643, 4649, 4651, 4657, 4663, 4673, 4679, 4691, 4703, 4721, 4723, 4729, 4733,
75   4751, 4759, 4783, 4787, 4789, 4793, 4799, 4801, 4813, 4817, 4831, 4861, 4871,
76   4877, 4889, 4903, 4909, 4919, 4931, 4933, 4937, 4943, 4951, 4957, 4967, 4969,
77   4973, 4987, 4993, 4999, 5003, 5009, 5011, 5021, 5023, 5039, 5051, 5059, 5077,
78   5081, 5087, 5099, 5101, 5107, 5113, 5119, 5147, 5153, 5167, 5171, 5179, 5189,
79   5197, 5209, 5227, 5231, 5233, 5237, 5261, 5273, 5279, 5281, 5297, 5303, 5309,
80   5323, 5333, 5347, 5351, 5381, 5387, 5393, 5399, 5407, 5413, 5417, 5419, 5431,
81   5437, 5441, 5443, 5449, 5471, 5477, 5479, 5483, 5501, 5503, 5507, 5519, 5521,
82   5527, 5531, 5557, 5563, 5569, 5573, 5581, 5591, 5623, 5639, 5641, 5647, 5651,
83   5653, 5657, 5659, 5669, 5683, 5689, 5693, 5701, 5711, 5717, 5737, 5741, 5743,
84   5749, 5779, 5783, 5791, 5801, 5807, 5813, 5821, 5827, 5839, 5843, 5849, 5851,
85   5857, 5861, 5867, 5869, 5879, 5881, 5897, 5903, 5923, 5927, 5939, 5953, 5981,
86   5987, 6007, 6011, 6029, 6037, 6043, 6047, 6053, 6067, 6073, 6079, 6089, 6091,
87   6101, 6113, 6121, 6131, 6133, 6143, 6151, 6163, 6173, 6197, 6199, 6203, 6211,
88   6217, 6221, 6229, 6247, 6257, 6263, 6269, 6271, 6277, 6287, 6299, 6301, 6311,
89   6317, 6323, 6329, 6337, 6343, 6353, 6359, 6361, 6367, 6373, 6379, 6389, 6397,
90   6421, 6427, 6449, 6451, 6469, 6473, 6481, 6491, 6521, 6529, 6547, 6551, 6553,
91   6563, 6569, 6571, 6577, 6581, 6599, 6607, 6619, 6637, 6653, 6659, 6661, 6673,
92   6679, 6689, 6691, 6701, 6703, 6709, 6719, 6733, 6737, 6761, 6763, 6779, 6781,
93   6791, 6793, 6803, 6823, 6827, 6829, 6833, 6841, 6857, 6863, 6869, 6871, 6883,
94   6899, 6907, 6911, 6917, 6947, 6949, 6959, 6961, 6967, 6971, 6977, 6983, 6991,
95   6997, 7001, 7013, 7019, 7027, 7039, 7043, 7057, 7069, 7079, 7103, 7109, 7121,
96   7127, 7129, 7151, 7159, 7177, 7187, 7193, 7207, 7211, 7213, 7219, 7229, 7237,
97   7243, 7247, 7253, 7283, 7297, 7307, 7309, 7321, 7331, 7333, 7349, 7351, 7369,
98   7393, 7411, 7417, 7433, 7451, 7457, 7459, 7477, 7481, 7487, 7489, 7499, 7507,
99   7517, 7523, 7529, 7537, 7541, 7547, 7549, 7559, 7561, 7573, 7577, 7583, 7589,
100   7591, 7603, 7607, 7621, 7639, 7643, 7649, 7669, 7673, 7681, 7687, 7691, 7699,
101   7703, 7717, 7723, 7727, 7741, 7753, 7757, 7759, 7789, 7793, 7817, 7823, 7829,
102   7841, 7853, 7867, 7873, 7877, 7879, 7883, 7901, 7907, 7919, 7927, 7933, 7937,
103   7949, 7951, 7963, 7993, 8009, 8011, 8017, 8039, 8053, 8059, 8069, 8081, 8087,
104   8089, 8093, 8101, 8111, 8117, 8123, 8147, 8161, 8167, 8171, 8179, 8191, 8209,
105   8219, 8221, 8231, 8233, 8237, 8243, 8263, 8269, 8273, 8287, 8291, 8293, 8297,
106   8311, 8317, 8329, 8353, 8363, 8369, 8377, 8387, 8389, 8419, 8423, 8429, 8431,
107   8443, 8447, 8461, 8467, 8501, 8513, 8521, 8527, 8537, 8539, 8543, 8563, 8573,
108   8581, 8597, 8599, 8609, 8623, 8627, 8629, 8641, 8647, 8663, 8669, 8677, 8681,
109   8689, 8693, 8699, 8707, 8713, 8719, 8731, 8737, 8741, 8747, 8753, 8761, 8779,
110   8783, 8803, 8807, 8819, 8821, 8831, 8837, 8839, 8849, 8861, 8863, 8867, 8887,
111   8893, 8923, 8929, 8933, 8941, 8951, 8963, 8969, 8971, 8999, 9001, 9007, 9011,
112   9013, 9029, 9041, 9043, 9049, 9059, 9067, 9091, 9103, 9109, 9127, 9133, 9137,
113   9151, 9157, 9161, 9173, 9181, 9187, 9199, 9203, 9209, 9221, 9227, 9239, 9241,
114   9257, 9277, 9281, 9283, 9293, 9311, 9319, 9323, 9337, 9341, 9343, 9349, 9371,
115   9377, 9391, 9397, 9403, 9413, 9419, 9421, 9431, 9433, 9437, 9439, 9461, 9463,
116   9467, 9473, 9479, 9491, 9497, 9511, 9521, 9533, 9539, 9547, 9551, 9587, 9601,
117   9613, 9619, 9623, 9629, 9631, 9643, 9649, 9661, 9677, 9679, 9689, 9697, 9719,
118   9721, 9733, 9739, 9743, 9749, 9767, 9769, 9781, 9787, 9791, 9803, 9811, 9817,
119   9829, 9833, 9839, 9851, 9857, 9859, 9871, 9883, 9887, 9901, 9907, 9923, 9929,
120   9931, 9941, 9949, 9967, 9973,
123 static const uint32_t largeprimes[] =
125   143687, 266083, 393073, 523417, 656683, 791447, 927853, 1066909, 1205653,
126   1346567, 1488133, 1630987, 1774937, 1918439, 2063573, 2209499, 2356573,
127   2503597, 2650931, 2799451, 2948161, 3097141, 3246119, 3396661, 3548471,
128   3698881, 3849943, 4002547, 4154791, 4307453, 4459951, 4613237, 4766999,
129   4920613, 5074633, 5229331, 5383601, 5539363, 5695831, 5851717, 6007247,
130   6163567, 6319483, 6477301, 6633877, 6791437, 6947753, 7104341, 7263143,
131   7421383, 7580459, 7738559, 7896949, 8056423, 8215633, 8375273, 8534489,
132   8693653, 8853359, 9014101, 9174793, 9334487, 9495601, 9655759, 9817057,
133   9979111, 10140307,
136 template <typename T, size_t N>
137 static inline constexpr size_t countof(T const (&)[N]) noexcept
139   return N;
142 static inline void request_migration()
144   AMPI_Migrate(AMPI_INFO_LB_SYNC);
147 class malloc_stress
149   // This class avoids heap allocations other than the test itself.
151   const int rank, p;
153   static constexpr size_t maxallocs = 2048;
154   static constexpr size_t maxsegments = 8;
155   static constexpr size_t maxtests = 32;
157   struct alloc_type
158   {
159     unsigned char *ptr;
160     unsigned int size;
161   };
162   alloc_type allocations[maxallocs]{};
163   alloc_type * allocation_queue[maxallocs*2]{};
165   unsigned int segments[maxsegments]{}; // segment boundaries
167   struct test_type
168   {
169     const char *description1, *description2;
170     unsigned long long duration_count;
171   };
172   test_type tests[maxtests]{};
174   unsigned long long overall_duration_count = 0;
176   // ideally these variables could be marked const after a certain point in execution
177   size_t numallocs = 0;
178   size_t numsegments = 1;
179   size_t queuesize;
180 #ifdef MIGRATION_OVERKILL
181   size_t migrationinterval;
182 #endif
184 public:
186   using clock_type = typename std::conditional<
187     std::chrono::high_resolution_clock::is_steady,
188     std::chrono::high_resolution_clock, std::chrono::steady_clock>::type;
189   using time_point = typename std::chrono::time_point<clock_type>;
190   using duration_type = typename std::chrono::duration<unsigned long long, std::nano>;
192   malloc_stress(int rank_, int p_)
193     : rank{rank_}, p{p_}
194   {
195     // this assertion needs to match the way the structures are populated below
196     static_assert(countof(primes) + 3*countof(largeprimes) <= maxallocs, "increase maxallocs");
198     segments[0] = 0;
200     for (size_t i = rank; i < countof(primes); i += p)
201     {
202       allocations[numallocs++].size = primes[i];
203     }
204     segments[numsegments++] = numallocs;
206     for (size_t i = rank; i < countof(largeprimes); i += p)
207     {
208       allocations[numallocs++].size = largeprimes[i];
209     }
210     segments[numsegments++] = numallocs;
212     for (size_t i = rank; i < countof(largeprimes); i += p)
213     {
214       allocations[numallocs++].size = primes[countof(primes)-1-i];
215       allocations[numallocs++].size = largeprimes[i];
216     }
217     segments[numsegments++] = numallocs;
219     queuesize = numallocs * 2;
220 #ifdef MIGRATION_OVERKILL
221     migrationinterval = numallocs / 31; // prime number divisor
222 #endif
223   }
225 private:
227   void alloc_loop(const size_t a_begin, const size_t a_end)
228   {
229     for (size_t a = a_begin; a < a_end; ++a)
230     {
231       alloc_type & allocation = *(allocation_queue[a]);
232       const size_t size = allocation.size;
233       const unsigned char rankbyte = rank & 0xFF;
235       const auto ptr = allocation.ptr = (unsigned char*)malloc(size);
236       if (ptr == nullptr)
237       {
238         printf("Error: malloc() failure in rank %d, test %zu, allocation %zu\n", rank, testnum, a);
239         MPI_Abort(MPI_COMM_WORLD, 1);
240       }
241       ptr[0] = ptr[size-1] = rankbyte;
243 #ifdef MIGRATION_OVERKILL
244       ++eventnum;
245       if (eventnum == migrationinterval)
246       {
247         eventnum = 0;
248         request_migration();
249       }
250 #endif
251     }
252   }
254   void free_loop(const size_t a_begin, const size_t a_end)
255   {
256     for (size_t a = a_begin; a < a_end; ++a)
257     {
258       alloc_type & allocation = *(allocation_queue[a]);
259       const size_t size = allocation.size;
260       const unsigned char rankbyte = rank & 0xFF;
262       const auto ptr = allocation.ptr;
263       if (ptr[0] != rankbyte || ptr[size-1] != rankbyte)
264       {
265         printf("Error: Correctness failure in rank %d, test %zu, allocation %zu\n", rank, testnum, a);
266         MPI_Abort(MPI_COMM_WORLD, 2);
267       }
268       free(ptr);
269       allocation.ptr = nullptr;
271 #ifdef MIGRATION_OVERKILL
272       ++eventnum;
273       if (eventnum == migrationinterval)
274       {
275         eventnum = 0;
276         request_migration();
277       }
278 #endif
279     }
280   }
282   void undetermined_loop(const size_t a_begin, const size_t a_end)
283   {
284     for (size_t a = a_begin; a < a_end; ++a)
285     {
286       alloc_type & allocation = *(allocation_queue[a]);
287       const size_t size = allocation.size;
288       const unsigned char rankbyte = rank & 0xFF;
290       auto ptr = allocation.ptr;
291       if (ptr == nullptr)
292       {
293         ptr = allocation.ptr = (unsigned char*)malloc(size);
294         if (ptr == nullptr)
295         {
296           printf("Error: malloc() failure in rank %d, test %zu, allocation %zu\n", rank, testnum, a);
297           MPI_Abort(MPI_COMM_WORLD, 1);
298         }
299         ptr[0] = ptr[size-1] = rankbyte;
300       }
301       else
302       {
303         if (ptr[0] != rankbyte || ptr[size-1] != rankbyte)
304         {
305           printf("Error: Correctness failure in rank %d, test %zu, allocation %zu\n", rank, testnum, a);
306           MPI_Abort(MPI_COMM_WORLD, 2);
307         }
308         allocation.ptr = nullptr;
309         free(ptr);
310       }
312 #ifdef MIGRATION_OVERKILL
313       ++eventnum;
314       if (eventnum == migrationinterval)
315       {
316         eventnum = 0;
317         request_migration();
318       }
319 #endif
320     }
321   }
323   void reset_queue()
324   {
325     alloc_type ** q = allocation_queue;
326     for (alloc_type * al = allocations, * al_end = al + numallocs; al < al_end; ++al)
327       *(q++) = al;
328   }
330   // tests freeing after every segment
331   using operate_segment = void (malloc_stress::*)(size_t);
333   void test_by_segment(operate_segment between = &malloc_stress::identity_segment, operate_segment after = &malloc_stress::identity_segment)
334   {
335 #ifdef MIGRATION_OVERKILL
336     eventnum = 0;
337 #endif
339     if (testnum == maxtests)
340     {
341       printf("Error: Increase maxtests.\n");
342       MPI_Abort(MPI_COMM_WORLD, 10);
343     }
345     unsigned long long test_duration_count = 0;
347     MPI_Barrier(MPI_COMM_WORLD);
349     for (size_t s = 1; s < numsegments; ++s)
350     {
351       const time_point t_0 = clock_type::now();
352       alloc_loop(segments[s-1], segments[s]);
353       const time_point t_1 = clock_type::now();
355       request_migration();
356       (this->*between)(s);
358       const time_point t_2 = clock_type::now();
359       free_loop(segments[s-1], segments[s]);
360       const time_point t_3 = clock_type::now();
362       request_migration();
363       (this->*after)(s);
365       const duration_type test_duration = (t_3 - t_2) + (t_1 - t_0);
366       test_duration_count += test_duration.count();
367     }
369     unsigned long long duration_count = 0;
370     MPI_Reduce(&test_duration_count, &duration_count, 1, MPI_UNSIGNED_LONG_LONG, MPI_MAX, 0, MPI_COMM_WORLD);
371     tests[testnum].duration_count += duration_count;
373     ++testnum;
374   }
376   // tests freeing at the end
377   using operate_all = void (malloc_stress::*)();
379   void test_all(operate_all between = &malloc_stress::identity_all, operate_all after = &malloc_stress::identity_all)
380   {
381 #ifdef MIGRATION_OVERKILL
382     eventnum = 0;
383 #endif
385     if (testnum == maxtests)
386     {
387       printf("Error: Increase maxtests.\n");
388       MPI_Abort(MPI_COMM_WORLD, 10);
389     }
391     MPI_Barrier(MPI_COMM_WORLD);
393     const time_point t_0 = clock_type::now();
394     alloc_loop(0, numallocs);
395     const time_point t_1 = clock_type::now();
397     request_migration();
398     (this->*between)();
400     const time_point t_2 = clock_type::now();
401     free_loop(0, numallocs);
402     const time_point t_3 = clock_type::now();
404     request_migration();
405     (this->*after)();
407     const duration_type test_duration = (t_3 - t_2) + (t_1 - t_0);
408     unsigned long long test_duration_count = test_duration.count();
410     unsigned long long duration_count = 0;
411     MPI_Reduce(&test_duration_count, &duration_count, 1, MPI_UNSIGNED_LONG_LONG, MPI_MAX, 0, MPI_COMM_WORLD);
412     tests[testnum].duration_count += duration_count;
414     ++testnum;
415   }
417   // tests with alloc and free interspersed
418   void test_mixed()
419   {
420 #ifdef MIGRATION_OVERKILL
421     eventnum = 0;
422 #endif
424     if (testnum == maxtests)
425     {
426       printf("Error: Increase maxtests.\n");
427       MPI_Abort(MPI_COMM_WORLD, 10);
428     }
430     MPI_Barrier(MPI_COMM_WORLD);
432     const time_point t_0 = clock_type::now();
433     undetermined_loop(0, numallocs*2);
434     const time_point t_1 = clock_type::now();
436     request_migration();
438     const duration_type test_duration = t_1 - t_0;
439     unsigned long long test_duration_count = test_duration.count();
441     unsigned long long duration_count = 0;
442     MPI_Reduce(&test_duration_count, &duration_count, 1, MPI_UNSIGNED_LONG_LONG, MPI_MAX, 0, MPI_COMM_WORLD);
443     tests[testnum].duration_count += duration_count;
445     ++testnum;
446   }
448   void identity_segment(size_t) { }
449   void identity_all() { }
451   void shuffle_segment(size_t s)
452   {
453     std::shuffle(allocation_queue + segments[s-1], allocation_queue + segments[s], generator);
454   }
455   void shuffle_segments()
456   {
457     for (size_t s = 1; s < numsegments; ++s)
458       shuffle_segment(s);
459   }
460   void shuffle_all()
461   {
462     std::shuffle(allocation_queue, allocation_queue + numallocs, generator);
463   }
465   void reverse_segment(size_t s)
466   {
467     std::reverse(allocation_queue + segments[s-1], allocation_queue + segments[s]);
468   }
469   void reverse_segments()
470   {
471     for (size_t s = 1; s < numsegments; ++s)
472       reverse_segment(s);
473   }
474   void reverse_all()
475   {
476     std::reverse(allocation_queue, allocation_queue + numallocs);
477   }
479   void shuffle_segments_and_reverse_all()
480   {
481     shuffle_segments();
482     reverse_all();
483   }
485   void tests_static(const char *desc1)
486   {
487     // allocate and free in the original order
489     tests[testnum].description1 = desc1;
490     tests[testnum].description2 = "FIFO free, batched by segment";
491     test_by_segment();
493     tests[testnum].description1 = desc1;
494     tests[testnum].description2 = "FIFO free, one batch";
495     test_all();
498     // allocate in the original order, free in the reverse order
500     tests[testnum].description1 = desc1;
501     tests[testnum].description2 = "LIFO free, batched by segment";
502     test_by_segment(&malloc_stress::reverse_segment, &malloc_stress::reverse_segment);
504     tests[testnum].description1 = desc1;
505     tests[testnum].description2 = "LIFO free within segments, one batch";
506     test_all(&malloc_stress::reverse_segments, &malloc_stress::reverse_segments);
508     tests[testnum].description1 = desc1;
509     tests[testnum].description2 = "LIFO free overall, one batch";
510     test_all(&malloc_stress::reverse_all, &malloc_stress::reverse_all);
511   }
513   void tests_dynamic(const char *desc1)
514   {
515     tests[testnum].description1 = desc1;
516     tests[testnum].description2 = "within segments, batched by segment";
517     test_by_segment(&malloc_stress::shuffle_segment);
519     tests[testnum].description1 = desc1;
520     tests[testnum].description2 = "within segments, FIFO segment order, one batch";
521     test_all(&malloc_stress::shuffle_segments);
523     tests[testnum].description1 = desc1;
524     tests[testnum].description2 = "within segments, LIFO segment order, one batch";
525     test_all(&malloc_stress::shuffle_segments_and_reverse_all); // reverse the segment order
527     tests[testnum].description1 = desc1;
528     tests[testnum].description2 = "overall, one batch";
529     test_all(&malloc_stress::shuffle_all);
530   }
532   void tests_mixed(const char *desc1)
533   {
534     {
535       alloc_type ** q = allocation_queue;
536       for (size_t s = 1; s < numsegments; ++s)
537       {
538         const size_t a_begin = segments[s-1], a_end = segments[s];
540         alloc_type ** const q_begin = q;
542         alloc_type * const al_begin = allocations + a_begin, * const al_end = allocations + a_end;
543         for (alloc_type * al = al_begin; al < al_end; ++al)
544           *(q++) = al;
545         for (alloc_type * al = al_begin; al < al_end; ++al)
546           *(q++) = al;
548         std::shuffle(q_begin, q, generator);
549       }
550     }
551     tests[testnum].description1 = desc1;
552     tests[testnum].description2 = "within segments, batched by segment";
553     test_mixed();
555     {
556       alloc_type ** q = allocation_queue;
557       for (alloc_type * al = allocations, * al_end = al + numallocs; al < al_end; ++al)
558         *(q++) = al;
559       for (alloc_type * al = allocations, * al_end = al + numallocs; al < al_end; ++al)
560         *(q++) = al;
561     }
562     std::shuffle(allocation_queue, allocation_queue + queuesize, generator);
563     tests[testnum].description1 = desc1;
564     tests[testnum].description2 = "overall, one batch";
565     test_mixed();
566   }
568   // these variables should be reset for each run
569   std::mt19937 generator;
570 #ifdef MIGRATION_OVERKILL
571   size_t eventnum;
572 #endif
573   size_t testnum;
575 public:
577   unsigned long long run(uint32_t randomseed)
578   {
579     generator.seed(randomseed);
580 #ifdef MIGRATION_OVERKILL
581     eventnum = 0;
582 #endif
583     testnum = 0;
585     MPI_Barrier(MPI_COMM_WORLD);
586     const time_point t_0 = clock_type::now();
588     reset_queue();
590     tests_static("distinct phases, without shuffling, ");
592     shuffle_segments();
593     tests_static("distinct phases, static shuffle within segments, ");
595     shuffle_all();
596     tests_static("distinct phases, static shuffle overall, ");
598     shuffle_all();
599     tests_dynamic("distinct phases, shuffling between phases ");
601     tests_mixed("one phase, shuffled ");
603     MPI_Barrier(MPI_COMM_WORLD);
604     const time_point t_1 = clock_type::now();
607     const duration_type my_overall_duration = t_1 - t_0;
608     const unsigned long long my_overall_duration_count = my_overall_duration.count();
609     overall_duration_count += my_overall_duration_count;
611     return my_overall_duration_count;
612   }
614   void average(size_t iterations)
615   {
616     for (size_t t = 0; t < testnum; ++t)
617     {
618       test_type & test = tests[t];
619       test.duration_count /= iterations;
620     }
622     overall_duration_count /= iterations;
623   }
625   void print(size_t iterations) const
626   {
627     unsigned long long total_duration_count = 0;
629     printf("detailed results:\n");
630     printf("              test # -  total duration (ns) -   avg. duration (ns) - description\n");
631     for (size_t t = 0; t < testnum; ++t)
632     {
633       const test_type & test = tests[t];
635       total_duration_count += test.duration_count;
637       printf("%20zu - %20llu - %20llu - %s%s\n", t+1,
638              test.duration_count, test.duration_count / iterations,
639              test.description1, test.description2);
640     }
641     printf("                 all - %20llu - %20llu - total\n",
642            total_duration_count, total_duration_count / iterations);
644     printf("overall test duration (ns) including migration, shuffling, barriers, etc: { total: %llu, average: %llu }\n",
645            overall_duration_count, overall_duration_count / iterations);
647     struct rusage usage;
648     getrusage(RUSAGE_SELF, &usage);
650     printf("memory usage high water mark ("
651 #ifdef __APPLE__
652       "bytes"
653 #else
654       "kilobytes"
655 #endif
656       "): %lu\n", usage.ru_maxrss);
657   }
660 static uint32_t generate_randomseed(int rank)
662   uint32_t randomseed;
664   if (rank == 0)
665     randomseed = std::chrono::system_clock::now().time_since_epoch().count();
667   MPI_Bcast(&randomseed, 1, MPI_UINT32_T, 0, MPI_COMM_WORLD);
669   return randomseed;
672 int main(int argc, char** argv)
674   MPI_Init(&argc, &argv);
676   int rank;
677   int p;
678   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
679   MPI_Comm_size(MPI_COMM_WORLD, &p);
681   malloc_stress test{rank, p};
683   bool doIterations = false;
684   size_t iterations = 1;
685   uint32_t randomseed = generate_randomseed(rank);
687   if (argc > 2)
688   {
689     if (strcmp(argv[1], "iterations") == 0)
690     {
691       doIterations = true;
692       iterations = (size_t)atoll(argv[2]);
693       if (iterations == 0)
694         iterations = 1;
695     }
696     else if (strcmp(argv[1], "randomseed") == 0)
697     {
698       randomseed = (uint32_t)atol(argv[2]);
699     }
700   }
702   if (doIterations)
703   {
704     if (rank == 0)
705     {
706       printf("malloc stress test: ./" isomalloc_method_str " (...) +vp%i iterations %zu\n", p, iterations);
707       printf("             round # -        duration (ns) - command\n");
708     }
710     for (size_t i = 0; i < iterations; ++i)
711     {
712       randomseed = generate_randomseed(rank);
714       const unsigned long long duration = test.run(randomseed);
716       if (rank == 0)
717         printf("%20zu - %20llu - ./" isomalloc_method_str " (...) +vp%i randomseed %" PRIu32 "\n",
718                i+1, duration, p, randomseed);
719     }
720   }
721   else
722   {
723     if (rank == 0)
724       printf("malloc stress test: ./" isomalloc_method_str " (...) +vp%i randomseed %" PRIu32 "\n",
725              p, randomseed);
727     test.run(randomseed);
728   }
730   if (rank == 0)
731     test.print(iterations);
733   MPI_Finalize();
735   return 0;