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
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
|
//-*-c-*-
#include <math_constants.h>
#include "linalg.h"
#include "matrix.h"
#include "rotate.h"
#include "intersect.h"
#define STACK_SIZE 100
/* flattened triangle mesh */
texture<float4, 1, cudaReadModeElementType> mesh;
/* lower/upper bounds for the bounding box associated with each node/leaf */
texture<float4, 1, cudaReadModeElementType> upper_bounds;
texture<float4, 1, cudaReadModeElementType> lower_bounds;
/* map to child nodes/triangles and the number of child nodes/triangles */
texture<uint, 1, cudaReadModeElementType> child_map_arr;
texture<uint, 1, cudaReadModeElementType> child_len_arr;
__device__ float3 make_float3(const float4 &a)
{
return make_float3(a.x, a.y, a.z);
}
/* Test the intersection between a ray starting at `origin` traveling in the
direction `direction` and the bounding box around node `i`. If the ray
intersects the bounding box return true, else return false. */
__device__ bool intersect_node(const float3 &origin, const float3 &direction, const int &i)
{
float3 lower_bound = make_float3(tex1Dfetch(lower_bounds, i));
float3 upper_bound = make_float3(tex1Dfetch(upper_bounds, i));
return intersect_box(origin, direction, lower_bound, upper_bound);
}
/* Find the intersection between a ray starting at `origin` traveling in the
direction `direction` and the global mesh texture. If the ray intersects
the texture return the index of the triangle which the ray intersected,
else return -1. */
__device__ int intersect_mesh(const float3 &origin, const float3& direction, const int first_leaf)
{
int triangle_idx = -1;
float distance;
float min_distance;
if (!intersect_node(origin, direction, 0))
return -1;
int stack[STACK_SIZE];
int *head = &stack[0];
int *node = &stack[1];
int *tail = &stack[STACK_SIZE-1];
*node = 0;
int i;
do
{
int child_map = tex1Dfetch(child_map_arr, *node);
int child_len = tex1Dfetch(child_len_arr, *node);
if (*node < first_leaf)
{
for (i=0; i < child_len; i++)
if (intersect_node(origin, direction, child_map+i))
*node++ = child_map+i;
if (node == head)
break;
node--;
}
else // node is a leaf
{
for (i=0; i < child_len; i++)
{
int mesh_idx = 3*(child_map + i);
float3 v0 = make_float3(tex1Dfetch(mesh, mesh_idx));
float3 v1 = make_float3(tex1Dfetch(mesh, mesh_idx+1));
float3 v2 = make_float3(tex1Dfetch(mesh, mesh_idx+2));
if (intersect_triangle(origin, direction, v0, v1, v2, distance))
{
if (triangle_idx == -1)
{
triangle_idx = child_map + i;
min_distance = distance;
continue;
}
if (distance < min_distance)
{
triangle_idx = child_map + i;
min_distance = distance;
}
}
} // triangle loop
node--;
} // node is a leaf
} // while loop
while (node != head);
return triangle_idx;
}
extern "C"
{
/* Translate `points` by the vector `v` */
__global__ void translate(int max_idx, float3 *points, float3 v)
{
int idx = blockIdx.x*blockDim.x + threadIdx.x;
if (idx > max_idx)
return;
*(points+idx) += v;
}
/* Rotate `points` through an angle `phi` counter-clockwise about the
axis `axis` (when looking towards +infinity). */
__global__ void rotate(int max_idx, float3 *points, float phi, float3 axis)
{
int idx = blockIdx.x*blockDim.x + threadIdx.x;
if (idx > max_idx)
return;
*(points+idx) = rotate(*(points+idx), phi, axis);
}
/* Trace the rays starting at `origins` traveling in the direction `directions`
to their intersection with the global mesh. If the ray intersects the mesh
set the pixel associated with the ray to a 32 bit color whose brightness is
determined by the cosine of the angle between the ray and the normal of the
triangle it intersected, else set the pixel to 0. */
__global__ void ray_trace(int max_idx, float3 *origins, float3 *directions, int first_leaf, int *pixels)
{
int idx = blockIdx.x*blockDim.x + threadIdx.x;
if (idx > max_idx)
return;
float3 origin = *(origins+idx);
float3 direction = *(directions+idx);
direction /= norm(direction);
int intersection_idx = intersect_mesh(origin, direction, first_leaf);
if (intersection_idx == -1)
{
*(pixels+idx) = 0;
}
else
{
int mesh_idx = 3*intersection_idx;
float3 v0 = make_float3(tex1Dfetch(mesh, mesh_idx));
float3 v1 = make_float3(tex1Dfetch(mesh, mesh_idx+1));
float3 v2 = make_float3(tex1Dfetch(mesh, mesh_idx+2));
*(pixels+idx) = get_color(direction, v0, v1, v2);
}
} // ray_trace
/* Propagate the photons starting at `origins` traveling in the direction
`directions` to their intersection with the global mesh. If the ray
intersects the mesh set the hit_solid array value associated with the
photon to the triangle index of the triangle the photon intersected, else
set the hit_solid array value to -1. */
__global__ void propagate(int max_idx, float3 *origins, float3 *directions, int first_leaf, int *hit_solids)
{
int idx = blockIdx.x*blockDim.x + threadIdx.x;
if (idx > max_idx)
return;
float3 origin = *(origins+idx);
float3 direction = *(directions+idx);
direction /= norm(direction);
*(hit_solids+idx) = intersect_mesh(origin, direction, first_leaf);
} // propagate
} // extern "c"
|