aboutsummaryrefslogtreecommitdiff
path: root/boltzgen/kernel/template/pattern/AA.cpp.mako
blob: a61bb41831fe696ce85834f5121d3988efbdaab6 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
<%def name="operator(name, params = None)">
void ${name}_tick(
      ${float_type}* f
    , std::size_t gid
% if params is not None:
% for param_type, param_name in params:
    , ${param_type} ${param_name}
% endfor
% endif
) {
    ${float_type}* preshifted_f = f + ${layout.cell_preshift('gid')};

% for i, c_i in enumerate(descriptor.c):
    const ${float_type} f_curr_${i} = preshifted_f[${layout.pop_offset(i)}];
% endfor

    ${caller.body()}

% for i, c_i in enumerate(descriptor.c):
    preshifted_f[${layout.pop_offset(i)}] = f_next_${descriptor.c.index(-c_i)};
% endfor
}

void ${name}_tock(
      ${float_type}* f
    , std::size_t gid
% if params is not None:
% for param_type, param_name in params:
    , ${param_type} ${param_name}
% endfor
% endif
) {
    ${float_type}* preshifted_f = f + ${layout.cell_preshift('gid')};

% for i, c_i in enumerate(descriptor.c):
    const ${float_type} f_curr_${descriptor.c.index(-c_i)} = preshifted_f[${layout.pop_offset(i) + layout.neighbor_offset(c_i)}];
% endfor

    ${caller.body()}

% for i, c_i in enumerate(descriptor.c):
    preshifted_f[${layout.pop_offset(i) + layout.neighbor_offset(c_i)}] = f_next_${i};
% endfor
}
</%def>

<%def name="functor(name, params = None)">
void ${name}_tick(
      const ${float_type}* f
    , std::size_t gid
% if params is not None:
% for param_type, param_name in params:
    , ${param_type} ${param_name}
% endfor
% endif
) {
    const ${float_type}* preshifted_f = f + ${layout.cell_preshift('gid')};

% for i, c_i in enumerate(descriptor.c):
    const ${float_type} f_curr_${i} = preshifted_f[${layout.pop_offset(descriptor.c.index(-c_i))}];
% endfor

    ${caller.body()}
}

void ${name}_tock(
      const ${float_type}* f
    , std::size_t gid
% if params is not None:
% for param_type, param_name in params:
    , ${param_type} ${param_name}
% endfor
% endif
) {
    const ${float_type}* preshifted_f = f + ${layout.cell_preshift('gid')};

% for i, c_i in enumerate(descriptor.c):
    const ${float_type} f_curr_${i} = preshifted_f[${layout.pop_offset(i) + layout.neighbor_offset(c_i)}];
% endfor

    ${caller.body()}
}
</%def>