1 | //#define REARRANGED_DOMAIN |
---|
2 | #include "hmpp_fun.h" |
---|
3 | |
---|
4 | |
---|
5 | |
---|
6 | #ifdef USING_LOCAL_DIRECTIVES |
---|
7 | #pragma hmpp evaRef codelet, target=CUDA args[*].transfer=atcall |
---|
8 | #endif |
---|
9 | void evaluate_segment_reflective( |
---|
10 | int N1, // Nids |
---|
11 | int N2, |
---|
12 | int N3, |
---|
13 | int N6, |
---|
14 | |
---|
15 | long ids[N1], |
---|
16 | long vol_ids[N2], // domain.boundary_cells |
---|
17 | long edge_ids[N2], // domain.boundary_edges |
---|
18 | double normals[N6], |
---|
19 | |
---|
20 | double stage_edge_values[N3], |
---|
21 | double bed_edge_values[N3], |
---|
22 | double height_edge_values[N3], |
---|
23 | double xmom_edge_values[N3], |
---|
24 | double ymom_edge_values[N3], |
---|
25 | double xvel_edge_values[N3], |
---|
26 | double yvel_edge_values[N3], |
---|
27 | |
---|
28 | double stage_boundary_values[N2], |
---|
29 | double bed_boundary_values[N2], |
---|
30 | double height_boundary_values[N2], |
---|
31 | double xmom_boundary_values[N2], |
---|
32 | double ymom_boundary_values[N2], |
---|
33 | double xvel_boundary_values[N2], |
---|
34 | double yvel_boundary_values[N2] |
---|
35 | ) |
---|
36 | { |
---|
37 | int k, id, id_vol, id_edge; |
---|
38 | |
---|
39 | double n1, n2, q1, q2, r1, r2; |
---|
40 | |
---|
41 | |
---|
42 | #pragma hmppcg gridify(k), & |
---|
43 | #pragma hmppcg & private(id, id_vol, id_edge, n1, n2, q1, q2, r1, r2), & |
---|
44 | #pragma hmppcg & global( edge_ids, normals, & |
---|
45 | #pragma hmppcg & stage_edge_values, bed_edge_values, & |
---|
46 | #pragma hmppcg & height_edge_values, xmom_edge_values, & |
---|
47 | #pragma hmppcg & ymom_edge_values, xvel_edge_values, & |
---|
48 | #pragma hmppcg & yvel_edge_values, stage_boundary_values, & |
---|
49 | #pragma hmppcg & bed_boundary_values, height_boundary_values, & |
---|
50 | #pragma hmppcg & xmom_boundary_values, ymom_boundary_values, & |
---|
51 | #pragma hmppcg & xvel_boundary_values, yvel_boundary_values) |
---|
52 | for(k=0; k<N1; k++) |
---|
53 | { |
---|
54 | //id = ids[k]; |
---|
55 | id = k; |
---|
56 | id_vol = vol_ids[ id ]; |
---|
57 | id_edge= edge_ids[ id ]; |
---|
58 | |
---|
59 | #ifndef REARRANGED_DOMAIN |
---|
60 | n1 = normals[id_vol*6 + id_edge*2]; |
---|
61 | n2 = normals[id_vol*6 + id_edge*2 + 1]; |
---|
62 | q1 = xmom_edge_values[id_vol*3 + id_edge]; |
---|
63 | q2 = ymom_edge_values[id_vol*3 + id_edge]; |
---|
64 | #else |
---|
65 | n1 = normals[id_vol + id_edge*2*N]; |
---|
66 | n2 = normals[id_vol + (id_edge*2+1)*N]; |
---|
67 | q1 = xmom_edge_values[id_vol + id_edge*N]; |
---|
68 | q2 = ymom_edge_values[id_vol + id_edge*N]; |
---|
69 | #endif |
---|
70 | |
---|
71 | r1 = -q1*n1 - q2*n2; |
---|
72 | r2 = -q1*n2 + q2*n1; |
---|
73 | |
---|
74 | #ifndef REARRANGED_DOMAIN |
---|
75 | stage_boundary_values[id] = stage_edge_values[id_vol*3 + id_edge]; |
---|
76 | bed_boundary_values[id] = bed_edge_values[id_vol*3 + id_edge]; |
---|
77 | height_boundary_values[id] = height_edge_values[id_vol*3 + id_edge]; |
---|
78 | |
---|
79 | q1 = xvel_edge_values[id_vol*3 + id_edge]; |
---|
80 | q2 = yvel_edge_values[id_vol*3 + id_edge]; |
---|
81 | #else |
---|
82 | stage_boundary_values[id] = stage_edge_values[id_vol + id_edge*N]; |
---|
83 | bed_boundary_values[id] = bed_edge_values[id_vol + id_edge*N]; |
---|
84 | height_boundary_values[id] = height_edge_values[id_vol + id_edge*N]; |
---|
85 | |
---|
86 | q1 = xvel_edge_values[id_vol + id_edge*N]; |
---|
87 | q2 = yvel_edge_values[id_vol + id_edge*N]; |
---|
88 | #endif |
---|
89 | |
---|
90 | xmom_boundary_values[id] = n1*r1 - n2*r2; |
---|
91 | ymom_boundary_values[id] = n2*r1 + n1*r2; |
---|
92 | |
---|
93 | |
---|
94 | r1 = q1*n1 + q2*n2; |
---|
95 | r2 = q1*n2 - q2*n1; |
---|
96 | |
---|
97 | xvel_boundary_values[id] = n1*r1 - n2*r2; |
---|
98 | yvel_boundary_values[id] = n2*r1 + n1*r2; |
---|
99 | |
---|
100 | vol_ids[k] = ids[k]; |
---|
101 | } |
---|
102 | } |
---|
103 | |
---|
104 | |
---|
105 | /* |
---|
106 | void evaluate_segment_dirichlet_1( |
---|
107 | int N, |
---|
108 | long * ids, |
---|
109 | long * vol_ids, |
---|
110 | long * edge_ids, |
---|
111 | double * boundary_values, |
---|
112 | double *edge_values) |
---|
113 | { |
---|
114 | const int k = |
---|
115 | threadIdx.x+threadIdx.y*blockDim.x+ |
---|
116 | (blockIdx.x+blockIdx.y*gridDim.x)*blockDim.x*blockDim.y; |
---|
117 | |
---|
118 | if (k >= N) |
---|
119 | return; |
---|
120 | |
---|
121 | int id = ids[k], |
---|
122 | id_vol = vol_ids[k], |
---|
123 | id_edge = edge_ids[k]; |
---|
124 | |
---|
125 | boundary_values[id] = edge_values[id_vol*3 + id_edge]; |
---|
126 | } |
---|
127 | |
---|
128 | |
---|
129 | void evaluate_segment_dirichlet_2( |
---|
130 | int N, |
---|
131 | double q_bdry, |
---|
132 | long * ids, |
---|
133 | double * boundary_values) |
---|
134 | { |
---|
135 | const int k = |
---|
136 | threadIdx.x+threadIdx.y*blockDim.x+ |
---|
137 | (blockIdx.x+blockIdx.y*gridDim.x)*blockDim.x*blockDim.y; |
---|
138 | |
---|
139 | int id = ids[k]; |
---|
140 | |
---|
141 | if (k >= N) |
---|
142 | return; |
---|
143 | |
---|
144 | boundary_values[id] = q_bdry; |
---|
145 | } |
---|
146 | |
---|
147 | */ |
---|