aboutsummaryrefslogtreecommitdiff
path: root/boltzgen/kernel/template/bounce_back_boundary.cl.mako
blob: 8bb4f6c97eaca2e2b30161a1fa00b6d3de2aec0e (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
__kernel void bounce_back_boundary_gid(__global ${float_type}* f_next,
                                       __global ${float_type}* f_prev,
                                       unsigned int gid)
{
    __global ${float_type}* preshifted_f_next = f_next + ${layout.cell_preshift('gid')};
    __global ${float_type}* preshifted_f_prev = f_prev + ${layout.cell_preshift('gid')};

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

<%
    subexpr, assignment = model.bgk(f_eq = model.equilibrium(resolve_moments = True))
%>

% for i, expr in enumerate(subexpr):
    const ${float_type} ${expr[0]} = ${ccode(expr[1])};
% endfor

% for i, expr in enumerate(assignment):
    const ${float_type} ${ccode(expr)}
% endfor

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

% if 'cell_list_dispatch' in extras:
__kernel void bounce_back_boundary_cells(__global ${float_type}* f_next,
                                         __global ${float_type}* f_prev,
                                         __global unsigned int*  cells)
{
    bounce_back_boundary_gid(f_next, f_prev, cells[get_global_id(0)]);
}
% endif