diff options
author | Anthony LaTorre <tlatorre9@gmail.com> | 2011-08-09 16:39:40 -0400 |
---|---|---|
committer | Anthony LaTorre <tlatorre9@gmail.com> | 2011-08-09 16:39:40 -0400 |
commit | 4f3e0b7709bb64ffc24f2a003509d5f480848239 (patch) | |
tree | f64f5efd1238a03e76ae280871f490a32b863ac3 /src | |
parent | 135ad2ca5b5ade4c52ae98f8fc7545fcc88fb449 (diff) | |
download | chroma-4f3e0b7709bb64ffc24f2a003509d5f480848239.tar.gz chroma-4f3e0b7709bb64ffc24f2a003509d5f480848239.tar.bz2 chroma-4f3e0b7709bb64ffc24f2a003509d5f480848239.zip |
switch to indexing child nodes by start and stop indices instead of start and length; this reduces a bit of arithmetic when traversing the bounding volume hierarchy and makes the code in Geometry.build() more concise. add an untested cuda kernel to interleave the bits in three uint64 arrays.
Diffstat (limited to 'src')
-rw-r--r-- | src/alpha.h | 38 | ||||
-rw-r--r-- | src/mesh.h | 39 | ||||
-rw-r--r-- | src/tools.cu | 17 |
3 files changed, 56 insertions, 38 deletions
diff --git a/src/alpha.h b/src/alpha.h index 5e3e803..a9963eb 100644 --- a/src/alpha.h +++ b/src/alpha.h @@ -17,14 +17,14 @@ __device__ void swap(T &a, T &b) struct HitList { - int size; - int indices[ALPHA_DEPTH]; + unsigned int size; + unsigned int indices[ALPHA_DEPTH]; float distances[ALPHA_DEPTH]; }; __device__ void add_to_hit_list(const float &distance, const int &index, HitList &h) { - int i; + unsigned int i; if (h.size >= ALPHA_DEPTH) { if (distance > h.distances[ALPHA_DEPTH-1]) @@ -50,7 +50,7 @@ __device__ void add_to_hit_list(const float &distance, const int &index, HitList } } -__device__ __noinline__ int get_color_alpha(const float3 &origin, const float3& direction) +__device__ int get_color_alpha(const float3 &origin, const float3& direction) { HitList h; h.size = 0; @@ -60,34 +60,34 @@ __device__ __noinline__ int get_color_alpha(const float3 &origin, const float3& if (!intersect_node(origin, direction, g_start_node)) return 0; - int stack[STACK_SIZE]; + unsigned int stack[STACK_SIZE]; - int *head = &stack[0]; - int *node = &stack[1]; - int *tail = &stack[STACK_SIZE-1]; + unsigned int *head = &stack[0]; + unsigned int *node = &stack[1]; + unsigned int *tail = &stack[STACK_SIZE-1]; *node = g_start_node; - int i; + unsigned int i; do { - int first_child = tex1Dfetch(node_map, *node); - int child_count = tex1Dfetch(node_length, *node); + unsigned int first_child = tex1Dfetch(node_map, *node); + unsigned int stop = tex1Dfetch(node_map_end, *node); - while (*node >= g_first_node && child_count == 1) + while (*node >= g_first_node && stop == first_child+1) { *node = first_child; first_child = tex1Dfetch(node_map, *node); - child_count = tex1Dfetch(node_length, *node); + stop = tex1Dfetch(node_map_end, *node); } if (*node >= g_first_node) { - for (i=0; i < child_count; i++) + for (i=first_child; i < stop; i++) { - if (intersect_node(origin, direction, first_child+i)) + if (intersect_node(origin, direction, i)) { - *node = first_child+i; + *node = i; node++; } } @@ -96,9 +96,9 @@ __device__ __noinline__ int get_color_alpha(const float3 &origin, const float3& } else // node is a leaf { - for (i=0; i < child_count; i++) + for (i=first_child; i < stop; i++) { - uint4 triangle_data = g_triangles[first_child+i]; + uint4 triangle_data = g_triangles[i]; float3 v0 = g_vertices[triangle_data.x]; float3 v1 = g_vertices[triangle_data.y]; @@ -106,7 +106,7 @@ __device__ __noinline__ int get_color_alpha(const float3 &origin, const float3& if (intersect_triangle(origin, direction, v0, v1, v2, distance)) { - add_to_hit_list(distance, first_child+i, h); + add_to_hit_list(distance, i, h); } } // triangle loop @@ -18,7 +18,8 @@ texture<float4, 1, cudaReadModeElementType> lower_bounds; /* map to child nodes/triangles and the number of child nodes/triangles */ texture<unsigned int, 1, cudaReadModeElementType> node_map; -texture<unsigned int, 1, cudaReadModeElementType> node_length; +texture<unsigned int, 1, cudaReadModeElementType> node_map_end; +//texture<unsigned int, 1, cudaReadModeElementType> node_length; __device__ float3 make_float3(const float4 &a) { @@ -48,7 +49,7 @@ __device__ bool intersect_node(const float3 &origin, const float3 &direction, co 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__ __noinline__ int intersect_mesh(const float3 &origin, const float3& direction, float &min_distance, int last_hit_triangle = -1) +__device__ int intersect_mesh(const float3 &origin, const float3& direction, float &min_distance, int last_hit_triangle = -1) { int triangle_index = -1; @@ -57,34 +58,34 @@ __device__ __noinline__ int intersect_mesh(const float3 &origin, const float3& d if (!intersect_node(origin, direction, g_start_node)) return -1; - int stack[STACK_SIZE]; + unsigned int stack[STACK_SIZE]; - int *head = &stack[0]; - int *node = &stack[1]; - int *tail = &stack[STACK_SIZE-1]; + unsigned int *head = &stack[0]; + unsigned int *node = &stack[1]; + unsigned int *tail = &stack[STACK_SIZE-1]; *node = g_start_node; - int i; + unsigned int i; do { - int first_child = tex1Dfetch(node_map, *node); - int child_count = tex1Dfetch(node_length, *node); + unsigned int first_child = tex1Dfetch(node_map, *node); + unsigned int stop = tex1Dfetch(node_map_end, *node); - while (*node >= g_first_node && child_count == 1) + while (*node >= g_first_node && stop == first_child+1) { *node = first_child; first_child = tex1Dfetch(node_map, *node); - child_count = tex1Dfetch(node_length, *node); + stop = tex1Dfetch(node_map_end, *node); } if (*node >= g_first_node) { - for (i=0; i < child_count; i++) + for (i=first_child; i < stop; i++) { - if (intersect_node(origin, direction, first_child+i)) + if (intersect_node(origin, direction, i)) { - *node = first_child+i; + *node = i; node++; } } @@ -93,12 +94,12 @@ __device__ __noinline__ int intersect_mesh(const float3 &origin, const float3& d } else // node is a leaf { - for (i=0; i < child_count; i++) + for (i=first_child; i < stop; i++) { - if (last_hit_triangle == first_child+i) + if (last_hit_triangle == i) continue; - uint4 triangle_data = g_triangles[first_child+i]; + uint4 triangle_data = g_triangles[i]; float3 v0 = g_vertices[triangle_data.x]; float3 v1 = g_vertices[triangle_data.y]; @@ -108,14 +109,14 @@ __device__ __noinline__ int intersect_mesh(const float3 &origin, const float3& d { if (triangle_index == -1) { - triangle_index = first_child + i; + triangle_index = i; min_distance = distance; continue; } if (distance < min_distance) { - triangle_index = first_child + i; + triangle_index = i; min_distance = distance; } } diff --git a/src/tools.cu b/src/tools.cu new file mode 100644 index 0000000..3d3fed7 --- /dev/null +++ b/src/tools.cu @@ -0,0 +1,17 @@ +//--*-c-*- + +extern "C" +{ + +__global__ void interleave(int nthreads, unsigned long long *x, unsigned long long *y, unsigned long long *z, int bits, unsigned long long *dest) +{ + int id = blockIdx.x*blockDim.x + threadIdx.x; + + if (id >= nthreads) + return; + + for (int i=0; i < bits; i++) + dest[id] |= (x[id] & 1 << i) << (2*i) | (y[id] & 1 << i) << (2*i+1) | (z[id] & 1 << i) << (2*i+2); +} + +} |