stream_add_omp_st.c 3.92 KB
Newer Older
1 2 3 4 5 6 7
#include <aml.h>
#include <assert.h>
#include <errno.h>
#include <omp.h>
#include <pthread.h>
#include <stdlib.h>

8 9
#include "utils.h"

10 11 12
#define ITER 10
#define CHUNKING 4

13 14 15 16 17 18 19
size_t numthreads, tilesz, esz;
unsigned long *a, *b, *c;
AML_TILING_1D_DECL(tiling);
AML_AREA_LINUX_DECL(slow);
AML_AREA_LINUX_DECL(fast);
AML_SCRATCH_PAR_DECL(sa);
AML_SCRATCH_PAR_DECL(sb);
20 21 22 23

int kernel(unsigned long *a, unsigned long *b, unsigned long *c, size_t n)
{
	size_t i;
24
	debug("%p = %p + %p [%zi]\n",c,a,b,n);
25 26 27 28 29
	for(i = 0; i < n; i++)
		c[i] = a[i] + b[i];
	return 0;
}

30
void do_work(unsigned long tid)
31
{
32 33

	int offset, i, ai, bi, oldai, oldbi;
34
	unsigned long *ap, *bp, *cp;
35
	void *abaseptr, *bbaseptr;
36
	offset = tid*CHUNKING;
37 38
	ap = aml_tiling_tilestart(&tiling, a, offset);
	bp = aml_tiling_tilestart(&tiling, b, offset);
39
	cp = aml_tiling_tilestart(&tiling, c, offset);
40 41 42 43 44 45 46 47
	abaseptr = aml_scratch_baseptr(&sa);
	bbaseptr = aml_scratch_baseptr(&sb);
	ai = -1; bi = -1;
	for(i = 0; i < CHUNKING-1; i++) {
		struct aml_scratch_request *ar, *br;
		oldai = ai; oldbi = bi;
		aml_scratch_async_pull(&sa, &ar, abaseptr, &ai, a, offset+i+1);
		aml_scratch_async_pull(&sb, &br, bbaseptr, &bi, b, offset+i+1);
48
		kernel(ap, bp, cp, esz);
49 50 51 52
		aml_scratch_wait(&sa, ar);
		aml_scratch_wait(&sb, br);
		ap = aml_tiling_tilestart(&tiling, abaseptr, ai);
		bp = aml_tiling_tilestart(&tiling, bbaseptr, bi);
53
		cp = aml_tiling_tilestart(&tiling, c, offset+i+1);
54 55 56
		aml_scratch_release(&sa, oldai);
		aml_scratch_release(&sb, oldbi);
	}
57
	kernel(ap, bp, cp, esz);
58
}
59

60 61
int main(int argc, char *argv[])
{
62 63 64 65
	AML_BINDING_SINGLE_DECL(binding);
	AML_ARENA_JEMALLOC_DECL(arena);
	AML_DMA_LINUX_SEQ_DECL(dma);
	unsigned long nodemask[AML_NODEMASK_SZ];
66
	struct bitmask *slowb, *fastb;
67
	aml_init(&argc, &argv);
68 69 70 71 72 73
	assert(argc == 4);

	log_init(argv[0]);
	fastb = numa_parse_nodestring_all(argv[1]);
	slowb = numa_parse_nodestring_all(argv[2]);
	unsigned long memsize = 1UL << atoi(argv[3]);
74 75 76 77 78 79 80

	/* use openmp env to figure out how many threads we want
	 * (we actually use 3x as much)
	 */
	#pragma omp parallel
	{
		numthreads = omp_get_num_threads();
81
		tilesz = memsize/(numthreads*CHUNKING);
82
		esz = tilesz/sizeof(unsigned long);
83
	}
84 85 86

	/* initialize all the supporting struct */
	assert(!aml_binding_init(&binding, AML_BINDING_TYPE_SINGLE, 0));
87
	assert(!aml_tiling_init(&tiling, AML_TILING_TYPE_1D, tilesz, memsize));
88 89 90 91 92 93
	assert(!aml_arena_jemalloc_init(&arena, AML_ARENA_JEMALLOC_TYPE_REGULAR));

	assert(!aml_area_linux_init(&slow,
				    AML_AREA_LINUX_MANAGER_TYPE_SINGLE,
				    AML_AREA_LINUX_MBIND_TYPE_REGULAR,
				    AML_AREA_LINUX_MMAP_TYPE_ANONYMOUS,
94
				    &arena, MPOL_BIND, slowb->maskp));
95 96 97 98
	assert(!aml_area_linux_init(&fast,
				    AML_AREA_LINUX_MANAGER_TYPE_SINGLE,
				    AML_AREA_LINUX_MBIND_TYPE_REGULAR,
				    AML_AREA_LINUX_MMAP_TYPE_ANONYMOUS,
99
				    &arena, MPOL_BIND, fastb->maskp));
100 101 102 103 104 105 106
	assert(!aml_dma_linux_seq_init(&dma, numthreads*2));
	assert(!aml_scratch_par_init(&sa, &fast, &slow, &dma, &tiling,
				     2*numthreads, numthreads));
	assert(!aml_scratch_par_init(&sb, &fast, &slow, &dma, &tiling,
				     2*numthreads, numthreads));

	/* allocation */
107 108 109
	a = aml_area_malloc(&slow, memsize);
	b = aml_area_malloc(&slow, memsize);
	c = aml_area_malloc(&fast, memsize);
110
	assert(a != NULL && b != NULL && c != NULL);
111

112
	unsigned long esize = memsize/sizeof(unsigned long);
113
	for(unsigned long i = 0; i < esize; i++) {
114 115 116
		a[i] = i;
		b[i] = esize - i;
		c[i] = 0;
117 118 119
	}

	/* run kernel */
120
	#pragma omp parallel for
121
	for(unsigned long i = 0; i < numthreads; i++) {
122
		do_work(i);
123 124 125 126
	}

	/* validate */
	for(unsigned long i = 0; i < esize; i++) {
127
		assert(c[i] == esize);
128 129
	}

130 131 132
	aml_scratch_par_destroy(&sa);
	aml_scratch_par_destroy(&sb);
	aml_dma_linux_seq_destroy(&dma);
133 134
	aml_area_free(&slow, a);
	aml_area_free(&slow, b);
135
	aml_area_free(&fast, c);
136 137 138 139
	aml_area_linux_destroy(&slow);
	aml_area_linux_destroy(&fast);
	aml_tiling_destroy(&tiling, AML_TILING_TYPE_1D);
	aml_binding_destroy(&binding, AML_BINDING_TYPE_SINGLE);
140 141 142
	aml_finalize();
	return 0;
}