summaryrefslogtreecommitdiffstats
path: root/cuda/3d
diff options
context:
space:
mode:
authorWillem Jan Palenstijn <Willem.Jan.Palenstijn@cwi.nl>2021-11-22 14:01:19 +0100
committerWillem Jan Palenstijn <wjp@usecode.org>2021-11-26 12:09:09 +0100
commit87ca44ca75a966cb6c1be88b201f9132ee176003 (patch)
treef4aebf7e0c20601cc7f61b3f1e9a4212105baee7 /cuda/3d
parent38d53a0900a41d013879168236278b8ba9cff99b (diff)
downloadastra-87ca44ca75a966cb6c1be88b201f9132ee176003.tar.gz
astra-87ca44ca75a966cb6c1be88b201f9132ee176003.tar.bz2
astra-87ca44ca75a966cb6c1be88b201f9132ee176003.tar.xz
astra-87ca44ca75a966cb6c1be88b201f9132ee176003.zip
Replace texref by texobj in par3d_fp
Diffstat (limited to 'cuda/3d')
-rw-r--r--cuda/3d/par3d_fp.cu79
1 files changed, 44 insertions, 35 deletions
diff --git a/cuda/3d/par3d_fp.cu b/cuda/3d/par3d_fp.cu
index 5daddc1..cae75f1 100644
--- a/cuda/3d/par3d_fp.cu
+++ b/cuda/3d/par3d_fp.cu
@@ -35,10 +35,6 @@ along with the ASTRA Toolbox. If not, see <http://www.gnu.org/licenses/>.
#include <cuda.h>
-typedef texture<float, 3, cudaReadModeElementType> texture3D;
-
-static texture3D gT_par3DVolumeTexture;
-
namespace astraCUDA3d {
static const unsigned int g_anglesPerBlock = 4;
@@ -63,21 +59,26 @@ __constant__ float gC_DetVY[g_MaxAngles];
__constant__ float gC_DetVZ[g_MaxAngles];
-static bool bindVolumeDataTexture(const cudaArray* array)
+static bool bindVolumeDataTexture(cudaArray* array, cudaTextureObject_t& texObj)
{
- cudaChannelFormatDesc channelDesc = cudaCreateChannelDesc<float>();
-
- gT_par3DVolumeTexture.addressMode[0] = cudaAddressModeBorder;
- gT_par3DVolumeTexture.addressMode[1] = cudaAddressModeBorder;
- gT_par3DVolumeTexture.addressMode[2] = cudaAddressModeBorder;
- gT_par3DVolumeTexture.filterMode = cudaFilterModeLinear;
- gT_par3DVolumeTexture.normalized = false;
-
- cudaBindTextureToArray(gT_par3DVolumeTexture, array, channelDesc);
-
- // TODO: error value?
-
- return true;
+ cudaChannelFormatDesc channelDesc =
+ cudaCreateChannelDesc(32, 0, 0, 0, cudaChannelFormatKindFloat);
+
+ cudaResourceDesc resDesc;
+ memset(&resDesc, 0, sizeof(resDesc));
+ resDesc.resType = cudaResourceTypeArray;
+ resDesc.res.array.array = array;
+
+ cudaTextureDesc texDesc;
+ memset(&texDesc, 0, sizeof(texDesc));
+ texDesc.addressMode[0] = cudaAddressModeBorder;
+ texDesc.addressMode[1] = cudaAddressModeBorder;
+ texDesc.addressMode[2] = cudaAddressModeBorder;
+ texDesc.filterMode = cudaFilterModeLinear;
+ texDesc.readMode = cudaReadModeElementType;
+ texDesc.normalizedCoords = 0;
+
+ return checkCuda(cudaCreateTextureObject(&texObj, &resDesc, &texDesc, NULL), "par3d_fp texture");
}
@@ -89,7 +90,7 @@ struct DIR_X {
__device__ float c0(float x, float y, float z) const { return x; }
__device__ float c1(float x, float y, float z) const { return y; }
__device__ float c2(float x, float y, float z) const { return z; }
- __device__ float tex(float f0, float f1, float f2) const { return tex3D(gT_par3DVolumeTexture, f0, f1, f2); }
+ __device__ float tex(cudaTextureObject_t tex, float f0, float f1, float f2) const { return tex3D<float>(tex, f0, f1, f2); }
__device__ float x(float f0, float f1, float f2) const { return f0; }
__device__ float y(float f0, float f1, float f2) const { return f1; }
__device__ float z(float f0, float f1, float f2) const { return f2; }
@@ -103,7 +104,7 @@ struct DIR_Y {
__device__ float c0(float x, float y, float z) const { return y; }
__device__ float c1(float x, float y, float z) const { return x; }
__device__ float c2(float x, float y, float z) const { return z; }
- __device__ float tex(float f0, float f1, float f2) const { return tex3D(gT_par3DVolumeTexture, f1, f0, f2); }
+ __device__ float tex(cudaTextureObject_t tex, float f0, float f1, float f2) const { return tex3D<float>(tex, f1, f0, f2); }
__device__ float x(float f0, float f1, float f2) const { return f1; }
__device__ float y(float f0, float f1, float f2) const { return f0; }
__device__ float z(float f0, float f1, float f2) const { return f2; }
@@ -117,7 +118,7 @@ struct DIR_Z {
__device__ float c0(float x, float y, float z) const { return z; }
__device__ float c1(float x, float y, float z) const { return x; }
__device__ float c2(float x, float y, float z) const { return y; }
- __device__ float tex(float f0, float f1, float f2) const { return tex3D(gT_par3DVolumeTexture, f1, f2, f0); }
+ __device__ float tex(cudaTextureObject_t tex, float f0, float f1, float f2) const { return tex3D<float>(tex, f1, f2, f0); }
__device__ float x(float f0, float f1, float f2) const { return f1; }
__device__ float y(float f0, float f1, float f2) const { return f2; }
__device__ float z(float f0, float f1, float f2) const { return f0; }
@@ -145,6 +146,7 @@ struct SCALE_NONCUBE {
template<class COORD, class SCALE>
__global__ void par3D_FP_t(float* D_projData, unsigned int projPitch,
+ cudaTextureObject_t tex,
unsigned int startSlice,
unsigned int startAngle, unsigned int endAngle,
const SDimensions3D dims,
@@ -210,7 +212,7 @@ __global__ void par3D_FP_t(float* D_projData, unsigned int projPitch,
for (int s = startSlice; s < endSlice; ++s)
{
- fVal += c.tex(f0, f1, f2);
+ fVal += c.tex(tex, f0, f1, f2);
f0 += 1.0f;
f1 += a1;
f2 += a2;
@@ -225,6 +227,7 @@ __global__ void par3D_FP_t(float* D_projData, unsigned int projPitch,
// Supersampling version
template<class COORD>
__global__ void par3D_FP_SS_t(float* D_projData, unsigned int projPitch,
+ cudaTextureObject_t tex,
unsigned int startSlice,
unsigned int startAngle, unsigned int endAngle,
const SDimensions3D dims, int iRaysPerDetDim,
@@ -301,7 +304,7 @@ __global__ void par3D_FP_SS_t(float* D_projData, unsigned int projPitch,
for (int s = startSlice; s < endSlice; ++s)
{
- fVal += c.tex(f0, f1, f2);
+ fVal += c.tex(tex, f0, f1, f2);
f0 += 1.0f;
f1 += a1;
f2 += a2;
@@ -415,7 +418,9 @@ __global__ void par3D_FP_SumSqW_t(float* D_projData, unsigned int projPitch,
bool Par3DFP_Array_internal(cudaPitchedPtr D_projData,
- const SDimensions3D& dims, unsigned int angleCount, const SPar3DProjection* angles,
+ cudaTextureObject_t D_texObj,
+ const SDimensions3D& dims,
+ unsigned int angleCount, const SPar3DProjection* angles,
const SProjectorParams3D& params)
{
// transfer angles to constant memory
@@ -520,29 +525,29 @@ bool Par3DFP_Array_internal(cudaPitchedPtr D_projData,
for (unsigned int i = 0; i < dims.iVolX; i += g_blockSlices)
if (params.iRaysPerDetDim == 1)
if (cube)
- par3D_FP_t<DIR_X><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), i, blockStart, blockEnd, dims, scube);
+ par3D_FP_t<DIR_X><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), D_texObj, i, blockStart, blockEnd, dims, scube);
else
- par3D_FP_t<DIR_X><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), i, blockStart, blockEnd, dims, snoncubeX);
+ par3D_FP_t<DIR_X><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), D_texObj,i, blockStart, blockEnd, dims, snoncubeX);
else
- par3D_FP_SS_t<DIR_X><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), i, blockStart, blockEnd, dims, params.iRaysPerDetDim, snoncubeX);
+ par3D_FP_SS_t<DIR_X><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), D_texObj,i, blockStart, blockEnd, dims, params.iRaysPerDetDim, snoncubeX);
} else if (blockDirection == 1) {
for (unsigned int i = 0; i < dims.iVolY; i += g_blockSlices)
if (params.iRaysPerDetDim == 1)
if (cube)
- par3D_FP_t<DIR_Y><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), i, blockStart, blockEnd, dims, scube);
+ par3D_FP_t<DIR_Y><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), D_texObj,i, blockStart, blockEnd, dims, scube);
else
- par3D_FP_t<DIR_Y><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), i, blockStart, blockEnd, dims, snoncubeY);
+ par3D_FP_t<DIR_Y><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), D_texObj,i, blockStart, blockEnd, dims, snoncubeY);
else
- par3D_FP_SS_t<DIR_Y><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), i, blockStart, blockEnd, dims, params.iRaysPerDetDim, snoncubeY);
+ par3D_FP_SS_t<DIR_Y><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), D_texObj,i, blockStart, blockEnd, dims, params.iRaysPerDetDim, snoncubeY);
} else if (blockDirection == 2) {
for (unsigned int i = 0; i < dims.iVolZ; i += g_blockSlices)
if (params.iRaysPerDetDim == 1)
if (cube)
- par3D_FP_t<DIR_Z><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), i, blockStart, blockEnd, dims, scube);
+ par3D_FP_t<DIR_Z><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), D_texObj,i, blockStart, blockEnd, dims, scube);
else
- par3D_FP_t<DIR_Z><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), i, blockStart, blockEnd, dims, snoncubeZ);
+ par3D_FP_t<DIR_Z><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), D_texObj,i, blockStart, blockEnd, dims, snoncubeZ);
else
- par3D_FP_SS_t<DIR_Z><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), i, blockStart, blockEnd, dims, params.iRaysPerDetDim, snoncubeZ);
+ par3D_FP_SS_t<DIR_Z><<<dimGrid, dimBlock, 0, stream>>>((float*)D_projData.ptr, D_projData.pitch/sizeof(float), D_texObj,i, blockStart, blockEnd, dims, params.iRaysPerDetDim, snoncubeZ);
}
}
@@ -569,10 +574,12 @@ bool Par3DFP(cudaPitchedPtr D_volumeData,
const SDimensions3D& dims, const SPar3DProjection* angles,
const SProjectorParams3D& params)
{
+ cudaTextureObject_t D_texObj;
+
// transfer volume to array
cudaArray* cuArray = allocateVolumeArray(dims);
transferVolumeToArray(D_volumeData, cuArray, dims);
- bindVolumeDataTexture(cuArray);
+ bindVolumeDataTexture(cuArray, D_texObj);
bool ret;
@@ -584,7 +591,7 @@ bool Par3DFP(cudaPitchedPtr D_volumeData,
cudaPitchedPtr D_subprojData = D_projData;
D_subprojData.ptr = (char*)D_projData.ptr + iAngle * D_projData.pitch;
- ret = Par3DFP_Array_internal(D_subprojData,
+ ret = Par3DFP_Array_internal(D_subprojData, D_texObj,
dims, iEndAngle - iAngle, angles + iAngle,
params);
if (!ret)
@@ -593,6 +600,8 @@ bool Par3DFP(cudaPitchedPtr D_volumeData,
cudaFreeArray(cuArray);
+ cudaDestroyTextureObject(D_texObj);
+
return ret;
}