_reg_tools_gpu.cu 19 KB
Newer Older
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
/*
 *  _reg_tools_gpu.cu
 *  
 *
 *  Created by Marc Modat and Pankaj Daga on 24/03/2009.
 *  Copyright (c) 2009, University College London. All rights reserved.
 *  Centre for Medical Image Computing (CMIC)
 *  See the LICENSE.txt file in the nifty_reg root folder
 *
 */

#ifndef _REG_TOOLS_GPU_CU
#define _REG_TOOLS_GPU_CU

#include "_reg_blocksize_gpu.h"
#include "_reg_tools_kernels.cu"


void reg_voxelCentric2NodeCentric_gpu(	nifti_image *targetImage,
					nifti_image *controlPointImage,
					float4 **voxelNMIGradientArray_d,
					float4 **nodeNMIGradientArray_d)
{
	const int nodeNumber = controlPointImage->nx * controlPointImage->ny * controlPointImage->nz;
	const int voxelNumber = targetImage->nx * targetImage->ny * targetImage->nz;
	const int3 targetImageDim = make_int3(targetImage->nx, targetImage->ny, targetImage->nz);
	const int3 gridSize = make_int3(controlPointImage->nx, controlPointImage->ny, controlPointImage->nz);
	const float3 voxelNodeRatio_h = make_float3(
		controlPointImage->dx / targetImage->dx,
		controlPointImage->dy / targetImage->dy,
		controlPointImage->dz / targetImage->dz);

	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_NodeNumber,&nodeNumber,sizeof(int)));
	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_TargetImageDim,&targetImageDim,sizeof(int3)));
	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_ControlPointImageDim,&gridSize,sizeof(int3)));
	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_VoxelNodeRatio,&voxelNodeRatio_h,sizeof(float3)));

	CUDA_SAFE_CALL(cudaBindTexture(0, gradientImageTexture, *voxelNMIGradientArray_d, voxelNumber*sizeof(float4)));

	const unsigned int Grid_reg_voxelCentric2NodeCentric = (unsigned int)ceil((float)nodeNumber/(float)Block_reg_voxelCentric2NodeCentric);
	dim3 B1(Block_reg_voxelCentric2NodeCentric,1,1);
	dim3 G1(Grid_reg_voxelCentric2NodeCentric,1,1);

	reg_voxelCentric2NodeCentric_kernel <<< G1, B1 >>> (*nodeNMIGradientArray_d);
	CUDA_SAFE_CALL(cudaThreadSynchronize());
46 47
#if _VERBOSE
	printf("[VERBOSE] reg_voxelCentric2NodeCentric_gpu kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
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
	       cudaGetErrorString(cudaGetLastError()),G1.x,G1.y,G1.z,B1.x,B1.y,B1.z);
#endif
}

void reg_convertNMIGradientFromVoxelToRealSpace_gpu(	mat44 *sourceMatrix_xyz,
							nifti_image *controlPointImage,
							float4 **nodeNMIGradientArray_d)
{
	const int nodeNumber = controlPointImage->nx * controlPointImage->ny * controlPointImage->nz;
	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_NodeNumber,&nodeNumber,sizeof(int)));

	float4 *matrix_h;CUDA_SAFE_CALL(cudaMallocHost((void **)&matrix_h, 3*sizeof(float4)));
	matrix_h[0] = make_float4(sourceMatrix_xyz->m[0][0], sourceMatrix_xyz->m[0][1], sourceMatrix_xyz->m[0][2], sourceMatrix_xyz->m[0][3]);
	matrix_h[1] = make_float4(sourceMatrix_xyz->m[1][0], sourceMatrix_xyz->m[1][1], sourceMatrix_xyz->m[1][2], sourceMatrix_xyz->m[1][3]);
	matrix_h[2] = make_float4(sourceMatrix_xyz->m[2][0], sourceMatrix_xyz->m[2][1], sourceMatrix_xyz->m[2][2], sourceMatrix_xyz->m[2][3]);
	float4 *matrix_d;
	CUDA_SAFE_CALL(cudaMalloc((void **)&matrix_d, 3*sizeof(float4)));
	CUDA_SAFE_CALL(cudaMemcpy(matrix_d, matrix_h, 3*sizeof(float4), cudaMemcpyHostToDevice));
    CUDA_SAFE_CALL(cudaFreeHost((void *)matrix_h));
	CUDA_SAFE_CALL(cudaBindTexture(0, matrixTexture, matrix_d, 3*sizeof(float4)));
	
	const unsigned int Grid_reg_convertNMIGradientFromVoxelToRealSpace =
		(unsigned int)ceil((float)nodeNumber/(float)Block_reg_convertNMIGradientFromVoxelToRealSpace);
	dim3 B1(Grid_reg_convertNMIGradientFromVoxelToRealSpace,1,1);
	dim3 G1(Block_reg_convertNMIGradientFromVoxelToRealSpace,1,1);

	_reg_convertNMIGradientFromVoxelToRealSpace_kernel <<< G1, B1 >>> (*nodeNMIGradientArray_d);
	CUDA_SAFE_CALL(cudaThreadSynchronize());
76 77
#if _VERBOSE
	printf("[VERBOSE] reg_convertNMIGradientFromVoxelToRealSpace: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98
	       cudaGetErrorString(cudaGetLastError()),G1.x,G1.y,G1.z,B1.x,B1.y,B1.z);
#endif
	CUDA_SAFE_CALL(cudaFree(matrix_d));
}


void reg_initialiseConjugateGradient(	float4 **nodeNMIGradientArray_d,
					float4 **conjugateG_d,
					float4 **conjugateH_d,
					int nodeNumber)
{
	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_NodeNumber,&nodeNumber,sizeof(int)));
	CUDA_SAFE_CALL(cudaBindTexture(0, gradientImageTexture, *nodeNMIGradientArray_d, nodeNumber*sizeof(float4)));

	const unsigned int Grid_reg_initialiseConjugateGradient =
		(unsigned int)ceil((float)nodeNumber/(float)Block_reg_initialiseConjugateGradient);
	dim3 B1(Grid_reg_initialiseConjugateGradient,1,1);
	dim3 G1(Block_reg_initialiseConjugateGradient,1,1);

	reg_initialiseConjugateGradient_kernel <<< G1, B1 >>> (*conjugateG_d);
	CUDA_SAFE_CALL(cudaThreadSynchronize());
99 100
#if _VERBOSE
	printf("[VERBOSE] reg_initialiseConjugateGradient: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124
	       cudaGetErrorString(cudaGetLastError()),G1.x,G1.y,G1.z,B1.x,B1.y,B1.z);
#endif
	CUDA_SAFE_CALL(cudaMemcpy(*conjugateH_d, *conjugateG_d, nodeNumber*sizeof(float4), cudaMemcpyDeviceToDevice));
}

void reg_GetConjugateGradient(	float4 **nodeNMIGradientArray_d,
				float4 **conjugateG_d,
				float4 **conjugateH_d,
				int nodeNumber)
{
	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_NodeNumber,&nodeNumber,sizeof(int)));
	CUDA_SAFE_CALL(cudaBindTexture(0, conjugateGTexture, *conjugateG_d, nodeNumber*sizeof(float4)));
	CUDA_SAFE_CALL(cudaBindTexture(0, conjugateHTexture, *conjugateH_d, nodeNumber*sizeof(float4)));
	CUDA_SAFE_CALL(cudaBindTexture(0, gradientImageTexture, *nodeNMIGradientArray_d, nodeNumber*sizeof(float4)));

	// gam = sum((grad+g)*grad)/sum(HxG);
	const unsigned int Grid_reg_GetConjugateGradient1 = (unsigned int)ceil((float)nodeNumber/(float)Block_reg_GetConjugateGradient1);
	dim3 B1(Block_reg_GetConjugateGradient1,1,1);
	dim3 G1(Grid_reg_GetConjugateGradient1,1,1);

	float2 *sum_d;
	CUDA_SAFE_CALL(cudaMalloc((void **)&sum_d, nodeNumber*sizeof(float2)));
	reg_GetConjugateGradient1_kernel <<< G1, B1 >>> (sum_d);
	CUDA_SAFE_CALL(cudaThreadSynchronize());
125 126
#if _VERBOSE
	printf("[VERBOSE] reg_GetConjugateGradient1 kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146
	       cudaGetErrorString(cudaGetLastError()),G1.x,G1.y,G1.z,B1.x,B1.y,B1.z);
#endif
	float2 *sum_h;CUDA_SAFE_CALL(cudaMallocHost((void **)&sum_h, nodeNumber*sizeof(float2)));
	CUDA_SAFE_CALL(cudaMemcpy(sum_h,sum_d, nodeNumber*sizeof(float2),cudaMemcpyDeviceToHost));
	CUDA_SAFE_CALL(cudaFree(sum_d));
	double dgg = 0.0;
	double gg = 0.0;
	for(int i=0; i<nodeNumber; i++){
		dgg += sum_h[i].x;
		gg += sum_h[i].y;
	}
	float gam = dgg / gg;
	CUDA_SAFE_CALL(cudaFreeHost((void *)sum_h));

	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_ScalingFactor,&gam,sizeof(float)));
	const unsigned int Grid_reg_GetConjugateGradient2 = (unsigned int)ceil((float)nodeNumber/(float)Block_reg_GetConjugateGradient2);
	dim3 B2(Block_reg_GetConjugateGradient2,1,1);
	dim3 G2(Grid_reg_GetConjugateGradient2,1,1);
	reg_GetConjugateGradient2_kernel <<< G2, B2 >>> (*nodeNMIGradientArray_d, *conjugateG_d, *conjugateH_d);
	CUDA_SAFE_CALL(cudaThreadSynchronize());
147 148
#if _VERBOSE
	printf("[VERBOSE] reg_GetConjugateGradient2 kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171
	       cudaGetErrorString(cudaGetLastError()),G1.x,G1.y,G1.z,B1.x,B1.y,B1.z);
#endif


}

float reg_getMaximalLength_gpu(	float4 **nodeNMIGradientArray_d,
				int nodeNumber)
{

	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_NodeNumber,&nodeNumber,sizeof(int)));
	CUDA_SAFE_CALL(cudaBindTexture(0, gradientImageTexture, *nodeNMIGradientArray_d, nodeNumber*sizeof(float4)));

	// each thread extract the maximal value out of 128
	const int threadNumber = (int)ceil((float)nodeNumber/128.0f);
	const unsigned int Grid_reg_getMaximalLength = (unsigned int)ceil((float)threadNumber/(float)Block_reg_getMaximalLength);
	dim3 B1(Block_reg_getMaximalLength,1,1);
	dim3 G1(Grid_reg_getMaximalLength,1,1);

	float *all_d;
	CUDA_SAFE_CALL(cudaMalloc((void **)&all_d, threadNumber*sizeof(float)));
	reg_getMaximalLength_kernel <<< G1, B1 >>> (all_d);
	CUDA_SAFE_CALL(cudaThreadSynchronize());
172 173
#if _VERBOSE
	printf("[VERBOSE] reg_getMaximalLength kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204
	       cudaGetErrorString(cudaGetLastError()),G1.x,G1.y,G1.z,B1.x,B1.y,B1.z);
#endif
	float *all_h;CUDA_SAFE_CALL(cudaMallocHost((void **)&all_h, nodeNumber*sizeof(float)));
	CUDA_SAFE_CALL(cudaMemcpy(all_h, all_d, threadNumber*sizeof(float),cudaMemcpyDeviceToHost));
	CUDA_SAFE_CALL(cudaFree(all_d));
	double maxDistance = 0.0f;
	for(int i=0; i<threadNumber; i++) maxDistance = all_h[i]>maxDistance?all_h[i]:maxDistance;
	CUDA_SAFE_CALL(cudaFreeHost((void *)all_h));

	return maxDistance;
}

void reg_updateControlPointPosition_gpu(nifti_image *controlPointImage,
					float4 **controlPointImageArray_d,
					float4 **bestControlPointPosition_d,
					float4 **nodeNMIGradientArray_d,
					float currentLength)
{
	const int nodeNumber = controlPointImage->nx * controlPointImage->ny * controlPointImage->nz;
	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_NodeNumber,&nodeNumber,sizeof(int)));
	CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_ScalingFactor,&currentLength,sizeof(float)));

	CUDA_SAFE_CALL(cudaBindTexture(0, controlPointTexture, *bestControlPointPosition_d, nodeNumber*sizeof(float4)));
	CUDA_SAFE_CALL(cudaBindTexture(0, gradientImageTexture, *nodeNMIGradientArray_d, nodeNumber*sizeof(float4)));

	const unsigned int Grid_reg_updateControlPointPosition = (unsigned int)ceil((float)nodeNumber/(float)Block_reg_updateControlPointPosition);
	dim3 B1(Block_reg_updateControlPointPosition,1,1);
	dim3 G1(Grid_reg_updateControlPointPosition,1,1);

	reg_updateControlPointPosition_kernel <<< G1, B1 >>> (*controlPointImageArray_d);
	CUDA_SAFE_CALL(cudaThreadSynchronize());
205 206
#if _VERBOSE
	printf("[VERBOSE] reg_updateControlPointPosition kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
207 208 209 210
	       cudaGetErrorString(cudaGetLastError()),G1.x,G1.y,G1.z,B1.x,B1.y,B1.z);
#endif
}

211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322
void reg_gaussianSmoothing_gpu( nifti_image *image,
                                float4 **imageArray_d,
                                float sigma,
                                bool smoothXYZ[8])

{
    const int voxelNumber = image->nx * image->ny * image->nz;
    const int3 imageDim = make_int3(image->nx, image->ny, image->nz);

    CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_ImageDim, &imageDim,sizeof(int3)));
    CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_VoxelNumber, &voxelNumber,sizeof(int3)));

    bool axisToSmooth[8];
    if(smoothXYZ==NULL){
        for(int i=0; i<8; i++) axisToSmooth[i]=true;
    }
    else{
        for(int i=0; i<8; i++) axisToSmooth[i]=smoothXYZ[i];
    }

    for(int n=1; n<4; n++){
        if(axisToSmooth[n]==true){
            float currentSigma;
            if(sigma>0) currentSigma=sigma/image->pixdim[n];
            else currentSigma=fabs(sigma); // voxel based if negative value
            int radius=(int)ceil(currentSigma*3.0f);
            if(radius>0){
                int kernelSize = 1+radius*2;
                float *kernel_h;
                CUDA_SAFE_CALL(cudaMallocHost((void **)&kernel_h, kernelSize*sizeof(float)));
                float kernelSum=0;
                for(int i=-radius; i<=radius; i++){
                    kernel_h[radius+i]=(float)(exp( -(i*i)/(2.0*currentSigma*currentSigma)) / (currentSigma*2.506628274631)); // 2.506... = sqrt(2*pi)
                    kernelSum += kernel_h[radius+i];
                }
                for(int i=0; i<kernelSize; i++)
                    kernel_h[i] /= kernelSum;
                float *kernel_d;
                CUDA_SAFE_CALL(cudaMalloc((void **)&kernel_d, kernelSize*sizeof(float)));
                CUDA_SAFE_CALL(cudaMemcpy(kernel_d, kernel_h, kernelSize*sizeof(float), cudaMemcpyHostToDevice));
                CUDA_SAFE_CALL(cudaFreeHost(kernel_h));

                float4 *smoothedImage;
                CUDA_SAFE_CALL(cudaMalloc((void **)&smoothedImage,voxelNumber*sizeof(float4)));

                CUDA_SAFE_CALL(cudaBindTexture(0, convolutionKernelTexture, kernel_d, kernelSize*sizeof(float)));
                CUDA_SAFE_CALL(cudaBindTexture(0, gradientImageTexture, *imageArray_d, voxelNumber*sizeof(float4)));
                unsigned int Grid_reg_ApplyConvolutionWindow;
                dim3 B,G;
                switch(n){
                    case 1:
                        Grid_reg_ApplyConvolutionWindow =
                            (unsigned int)ceil((float)voxelNumber/(float)Block_reg_ApplyConvolutionWindowAlongX);
                        B=dim3(Block_reg_ApplyConvolutionWindowAlongX,1,1);
                        G=dim3(Grid_reg_ApplyConvolutionWindow,1,1);
                        _reg_ApplyConvolutionWindowAlongX_kernel <<< G, B >>> (smoothedImage, kernelSize);
                        CUDA_SAFE_CALL(cudaThreadSynchronize());
#if _VERBOSE
                        printf("[VERBOSE] reg_ApplyConvolutionWindowAlongX_kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
                            cudaGetErrorString(cudaGetLastError()),G.x,G.y,G.z,B.x,B.y,B.z);
#endif
                        break;
                    case 2:
                        Grid_reg_ApplyConvolutionWindow =
                            (unsigned int)ceil((float)voxelNumber/(float)Block_reg_ApplyConvolutionWindowAlongY);
                        B=dim3(Block_reg_ApplyConvolutionWindowAlongY,1,1);
                        G=dim3(Grid_reg_ApplyConvolutionWindow,1,1);
                        _reg_ApplyConvolutionWindowAlongY_kernel <<< G, B >>> (smoothedImage, kernelSize);
                        CUDA_SAFE_CALL(cudaThreadSynchronize());
#if _VERBOSE
                        printf("[VERBOSE] reg_ApplyConvolutionWindowAlongY_kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
                            cudaGetErrorString(cudaGetLastError()),G.x,G.y,G.z,B.x,B.y,B.z);
#endif
                        break;
                    case 3:
                        Grid_reg_ApplyConvolutionWindow =
                            (unsigned int)ceil((float)voxelNumber/(float)Block_reg_ApplyConvolutionWindowAlongZ);
                        B=dim3(Block_reg_ApplyConvolutionWindowAlongZ,1,1);
                        G=dim3(Grid_reg_ApplyConvolutionWindow,1,1);
                        _reg_ApplyConvolutionWindowAlongZ_kernel <<< G, B >>> (smoothedImage, kernelSize);
                        CUDA_SAFE_CALL(cudaThreadSynchronize());
#if _VERBOSE
                        printf("[VERBOSE] reg_ApplyConvolutionWindowAlongZ_kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
                            cudaGetErrorString(cudaGetLastError()),G.x,G.y,G.z,B.x,B.y,B.z);
#endif
                        break;
                }
                CUDA_SAFE_CALL(cudaFree(kernel_d));
                CUDA_SAFE_CALL(cudaMemcpy(*imageArray_d, smoothedImage, voxelNumber*sizeof(float4), cudaMemcpyDeviceToDevice));
                CUDA_SAFE_CALL(cudaFree(smoothedImage));
            }
        }
    }

}


void reg_smoothImageForCubicSpline_gpu( nifti_image *image,
                                        float4 **imageArray_d,
                                        int *smoothingRadius)
{
    const int voxelNumber = image->nx * image->ny * image->nz;
    const int3 imageDim = make_int3(image->nx, image->ny, image->nz);

    CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_ImageDim, &imageDim,sizeof(int3)));
    CUDA_SAFE_CALL(cudaMemcpyToSymbol(c_VoxelNumber, &voxelNumber,sizeof(int3)));

    for(int n=0; n<3; n++){
        if(smoothingRadius[n]>0){
            int kernelSize = 1+smoothingRadius[n]*2;
            float *kernel_h;
            CUDA_SAFE_CALL(cudaMallocHost((void **)&kernel_h, kernelSize*sizeof(float)));
323
//             float kernelSum=0;
324 325 326 327
            for(int i=-smoothingRadius[n]; i<=smoothingRadius[n]; i++){
                float coeff = fabs(2.0f*(float)(i)/smoothingRadius[n]);
                if(coeff<1.0f)  kernel_h[smoothingRadius[n]+i] = 2.0f/3.0f - coeff*coeff + 0.5f*coeff*coeff*coeff;
                else        kernel_h[smoothingRadius[n]+i] = -(coeff-2.0f)*(coeff-2.0f)*(coeff-2.0f)/6.0f;
328
//                 kernelSum += kernel_h[smoothingRadius[n]+i];
329
            }
330
//             for(int i=0; i<kernelSize; i++) kernel_h[i] /= kernelSum;
331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388
            float *kernel_d;
            CUDA_SAFE_CALL(cudaMalloc((void **)&kernel_d, kernelSize*sizeof(float)));
            CUDA_SAFE_CALL(cudaMemcpy(kernel_d, kernel_h, kernelSize*sizeof(float), cudaMemcpyHostToDevice));
            CUDA_SAFE_CALL(cudaFreeHost(kernel_h));
            CUDA_SAFE_CALL(cudaBindTexture(0, convolutionKernelTexture, kernel_d, kernelSize*sizeof(float)));

            float4 *smoothedImage_d;
            CUDA_SAFE_CALL(cudaMalloc((void **)&smoothedImage_d,voxelNumber*sizeof(float4)));

            CUDA_SAFE_CALL(cudaBindTexture(0, gradientImageTexture, *imageArray_d, voxelNumber*sizeof(float4)));

            unsigned int Grid_reg_ApplyConvolutionWindow;
            dim3 B,G;
            switch(n){
                case 0:
                    Grid_reg_ApplyConvolutionWindow =
                        (unsigned int)ceil((float)voxelNumber/(float)Block_reg_ApplyConvolutionWindowAlongX);
                    B=dim3(Block_reg_ApplyConvolutionWindowAlongX,1,1);
                    G=dim3(Grid_reg_ApplyConvolutionWindow,1,1);
                    _reg_ApplyConvolutionWindowAlongX_kernel <<< G, B >>> (smoothedImage_d, kernelSize);
                    CUDA_SAFE_CALL(cudaThreadSynchronize());
#if _VERBOSE
                    printf("[VERBOSE] reg_ApplyConvolutionWindowAlongX_kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
                        cudaGetErrorString(cudaGetLastError()),G.x,G.y,G.z,B.x,B.y,B.z);
#endif
                    break;
                case 1:
                    Grid_reg_ApplyConvolutionWindow =
                        (unsigned int)ceil((float)voxelNumber/(float)Block_reg_ApplyConvolutionWindowAlongY);
                    B=dim3(Block_reg_ApplyConvolutionWindowAlongY,1,1);
                    G=dim3(Grid_reg_ApplyConvolutionWindow,1,1);
                    _reg_ApplyConvolutionWindowAlongY_kernel <<< G, B >>> (smoothedImage_d, kernelSize);
                    CUDA_SAFE_CALL(cudaThreadSynchronize());
#if _VERBOSE
                    printf("[VERBOSE] reg_ApplyConvolutionWindowAlongY_kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
                        cudaGetErrorString(cudaGetLastError()),G.x,G.y,G.z,B.x,B.y,B.z);
#endif
                    break;
                case 2:
                    Grid_reg_ApplyConvolutionWindow =
                        (unsigned int)ceil((float)voxelNumber/(float)Block_reg_ApplyConvolutionWindowAlongZ);
                    B=dim3(Block_reg_ApplyConvolutionWindowAlongZ,1,1);
                    G=dim3(Grid_reg_ApplyConvolutionWindow,1,1);
                    _reg_ApplyConvolutionWindowAlongZ_kernel <<< G, B >>> (smoothedImage_d, kernelSize);
                    CUDA_SAFE_CALL(cudaThreadSynchronize());
#if _VERBOSE
                    printf("[VERBOSE] reg_ApplyConvolutionWindowAlongZ_kernel: %s - Grid size [%i %i %i] - Block size [%i %i %i]\n",
                        cudaGetErrorString(cudaGetLastError()),G.x,G.y,G.z,B.x,B.y,B.z);
#endif
                    break;
            }
            CUDA_SAFE_CALL(cudaFree(kernel_d));
            CUDA_SAFE_CALL(cudaMemcpy(*imageArray_d, smoothedImage_d, voxelNumber*sizeof(float4), cudaMemcpyDeviceToDevice));
            CUDA_SAFE_CALL(cudaFree(smoothedImage_d));
        }
    }
}

389 390
#endif