diff options
Diffstat (limited to 'src')
-rw-r--r-- | src/kernel.cu | 28 |
1 files changed, 25 insertions, 3 deletions
diff --git a/src/kernel.cu b/src/kernel.cu index a26460a..cc567b8 100644 --- a/src/kernel.cu +++ b/src/kernel.cu @@ -292,7 +292,24 @@ __global__ void ray_trace_alpha(int nthreads, float3 *positions, float3 *directi } // ray_trace -__global__ void propagate(int first_photon, int nthreads, int *photon_offsets, curandState *rng_states, +__global__ void swap(int *values, int nswap, int *offset_a, int *offset_b) +{ + int id = blockIdx.x*blockDim.x + threadIdx.x; + + if (id < nswap) { + int a = offset_a[id]; + int b = offset_b[id]; + + int tmp = values[a]; + values[a] = values[b]; + values[b] = tmp; + } +} + +__global__ void propagate(int first_photon, int nthreads, + unsigned int *input_queue, + unsigned int *output_queue, + curandState *rng_states, float3 *positions, float3 *directions, float *wavelengths, float3 *polarizations, float *times, unsigned int *histories, int *last_hit_triangles, int max_steps) { @@ -303,7 +320,7 @@ __global__ void propagate(int first_photon, int nthreads, int *photon_offsets, c curandState rng = rng_states[id]; - int photon_id = photon_offsets[first_photon + id]; + int photon_id = input_queue[first_photon + id]; Photon p; p.position = positions[photon_id]; @@ -364,7 +381,12 @@ __global__ void propagate(int first_photon, int nthreads, int *photon_offsets, c times[photon_id] = p.time; histories[photon_id] = p.history; last_hit_triangles[photon_id] = p.last_hit_triangle; - + + // Not done, put photon in output queue + if ( (p.history & (NO_HIT | BULK_ABSORB | SURFACE_DETECT | SURFACE_ABSORB)) == 0) { + int out_idx = atomicAdd(output_queue, 1); + output_queue[out_idx] = photon_id; + } } // propagate } // extern "c" |