Blame view

stim/cuda/ivote/update_dir.cuh 7.36 KB
13fe3c84   Laila Saadatifard   update the stimli...
1
2
3
  #ifndef STIM_CUDA_UPDATE_DIR_H
  #define STIM_CUDA_UPDATE_DIR_H
  
13fe3c84   Laila Saadatifard   update the stimli...
4
5
  # include <iostream>
  # include <cuda.h>
96f9b10f   Laila Saadatifard   change the header...
6
  #include <stim/cuda/cudatools.h>
13fe3c84   Laila Saadatifard   update the stimli...
7
8
9
10
11
12
13
  #include <stim/cuda/sharedmem.cuh>
  
  namespace stim{
  	namespace cuda{
  	
  		// this kernel calculates the voting direction for the next iteration based on the angle between the location of this voter and the maximum vote value in its voting area.
  		template<typename T>
84ca9bba   Laila Saadatifard   fix some bugs in ...
14
  		__global__ void cuda_update_dir(T* gpuDir, cudaTextureObject_t in, T* gpuGrad, T* gpuTable, T phi, int rmax,  int x,  int y){
13fe3c84   Laila Saadatifard   update the stimli...
15
16
17
18
19
20
21
  
  			//generate a pointer to shared memory (size will be specified as a kernel parameter)
  			extern __shared__ float s_vote[];
  
  			//calculate the start point for this block
  			int bxi = blockIdx.x * blockDim.x;
  			
13fe3c84   Laila Saadatifard   update the stimli...
22
23
  			// calculate the 2D coordinates for this current thread.
  			int xi = bxi + threadIdx.x;
84ca9bba   Laila Saadatifard   fix some bugs in ...
24
  			int yi = blockIdx.y * blockDim.y + threadIdx.y;
13fe3c84   Laila Saadatifard   update the stimli...
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
  
  			// convert 2D coordinates to 1D
  			int i = yi * x + xi;
  
  			// calculate the voting direction based on the grtadient direction
  			float theta = gpuGrad[2*i];
  
  			//initialize the vote direction to zero
  			gpuDir[i] = 0;
  
  			// define a local variable to maximum value of the vote image in the voting area for this voter
  			float max = 0;
  
  			// define two local variables for the x and y coordinations where the maximum happened
  			int id_x = 0;
  			int id_y = 0;
  
84ca9bba   Laila Saadatifard   fix some bugs in ...
42
43
44
  			//calculate the width of the shared memory block
  			int swidth = 2 * rmax + blockDim.x;
  
13fe3c84   Laila Saadatifard   update the stimli...
45
  			// compute the size of window which will be checked for finding the voting area for this voter
84ca9bba   Laila Saadatifard   fix some bugs in ...
46
47
  			int x_table = 2*rmax +1;
  			int rmax_sq = rmax * rmax;
13fe3c84   Laila Saadatifard   update the stimli...
48
49
  			int tx_rmax = threadIdx.x + rmax;
  			int bxs = bxi - rmax;
bf731970   Laila Saadatifard   fix some bugs in ...
50
  						
84ca9bba   Laila Saadatifard   fix some bugs in ...
51
  			for(int yr = -rmax; yr <= rmax; yr++){
13fe3c84   Laila Saadatifard   update the stimli...
52
53
54
55
56
57
58
59
60
  
  				//copy the portion of the image necessary for this block to shared memory
  				__syncthreads();
  				stim::cuda::sharedMemcpy_tex2D<float>(s_vote, in, bxs, yi + yr , swidth, 1, threadIdx, blockDim);
  				__syncthreads();
  				
  				//if the current thread is outside of the image, it doesn't have to be computed
  				if(xi < x && yi < y){
  
84ca9bba   Laila Saadatifard   fix some bugs in ...
61
  					for(int xr = -rmax; xr <= rmax; xr++){
13fe3c84   Laila Saadatifard   update the stimli...
62
63
64
65
66
  
  						unsigned int ind_t = (rmax - yr) * x_table + rmax - xr;
  
  						// calculate the angle between the voter and the current pixel in x and y directions
  						float atan_angle = gpuTable[ind_t];
bf731970   Laila Saadatifard   fix some bugs in ...
67
  						
13fe3c84   Laila Saadatifard   update the stimli...
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
  						// calculate the voting direction based on the grtadient direction
  						int idx_share_update = xr + tx_rmax ;
  						float share_vote = s_vote[idx_share_update];
  						
  						// check if the current pixel is located in the voting area of this voter.
  						if (((xr * xr + yr *yr)< rmax_sq) && (abs(atan_angle - theta) <phi)){
  
  						// compare the vote value of this pixel with the max value to find the maxima and its index.
  							if  (share_vote>max) {
  
  								max = share_vote;
  								id_x =  xr;
  								id_y =  yr;
  							}
  						}
  					}
  				}
  			}
bf731970   Laila Saadatifard   fix some bugs in ...
86
  							
13fe3c84   Laila Saadatifard   update the stimli...
87
  		unsigned int ind_m = (rmax - id_y) * x_table + (rmax - id_x);
13fe3c84   Laila Saadatifard   update the stimli...
88
89
  		float new_angle = gpuTable[ind_m];
  
84ca9bba   Laila Saadatifard   fix some bugs in ...
90
91
  		if(xi < x && yi < y)
  			gpuDir[i] = new_angle;
13fe3c84   Laila Saadatifard   update the stimli...
92
93
94
95
96
  
  		}
  
  		// this kernel updates the gradient direction by the calculated voting direction.
  		template<typename T>
84ca9bba   Laila Saadatifard   fix some bugs in ...
97
  		__global__ void cuda_update_grad(T* gpuGrad, T* gpuDir, int x, int y){
13fe3c84   Laila Saadatifard   update the stimli...
98
  
13fe3c84   Laila Saadatifard   update the stimli...
99
100
  			// calculate the 2D coordinates for this current thread.
  			int xi = blockIdx.x * blockDim.x + threadIdx.x;
84ca9bba   Laila Saadatifard   fix some bugs in ...
101
102
  			int yi = blockIdx.y * blockDim.y + threadIdx.y;
  		
13fe3c84   Laila Saadatifard   update the stimli...
103
104
  			// convert 2D coordinates to 1D
  			int i = yi * x + xi;
bf731970   Laila Saadatifard   fix some bugs in ...
105
  			
13fe3c84   Laila Saadatifard   update the stimli...
106
107
108
  			//update the gradient image with the vote direction
  			gpuGrad[2*i] = gpuDir[i];
  		}
bf731970   Laila Saadatifard   fix some bugs in ...
109
  		
13fe3c84   Laila Saadatifard   update the stimli...
110
111
112
113
114
115
  		template<typename T>
  		void gpu_update_dir(T* gpuVote, T* gpuGrad, T* gpuTable, T phi, unsigned int rmax, unsigned int x, unsigned int y){
  
  			//get the number of pixels in the image
  			unsigned int pixels = x * y;
  			unsigned int bytes = sizeof(T) * pixels;
bf731970   Laila Saadatifard   fix some bugs in ...
116
  						
13fe3c84   Laila Saadatifard   update the stimli...
117
118
119
  			unsigned int max_threads = stim::maxThreadsPerBlock();
  			dim3 threads(max_threads, 1);
  			dim3 blocks(x/threads.x + (x %threads.x == 0 ? 0:1) , y);
13fe3c84   Laila Saadatifard   update the stimli...
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
  			
  			//define a channel descriptor for a single 32-bit channel
  			cudaChannelFormatDesc channelDesc =
  					   cudaCreateChannelDesc(32, 0, 0, 0,
  											 cudaChannelFormatKindFloat);
  			cudaArray* cuArray;												//declare the cuda array
  			cudaMallocArray(&cuArray, &channelDesc, x, y);			//allocate the cuda array
  
  			// Copy the image data from global memory to the array
  			cudaMemcpyToArray(cuArray, 0, 0, gpuVote, bytes,
  							  cudaMemcpyDeviceToDevice);
  
  			// Specify texture
  			struct cudaResourceDesc resDesc;				//create a resource descriptor
  			memset(&resDesc, 0, sizeof(resDesc));			//set all values to zero
  			resDesc.resType = cudaResourceTypeArray;		//specify the resource descriptor type
  			resDesc.res.array.array = cuArray;				//add a pointer to the cuda array
  
  			// Specify texture object parameters
  			struct cudaTextureDesc texDesc;							//create a texture descriptor
  			memset(&texDesc, 0, sizeof(texDesc));					//set all values in the texture descriptor to zero
  			texDesc.addressMode[0]   = cudaAddressModeWrap;			//use wrapping (around the edges)
  			texDesc.addressMode[1]   = cudaAddressModeWrap;
  			texDesc.filterMode       = cudaFilterModePoint;		//use linear filtering
  			texDesc.readMode         = cudaReadModeElementType;		//reads data based on the element type (32-bit floats)
  			texDesc.normalizedCoords = 0;							//not using normalized coordinates
  
  			// Create texture object
  			cudaTextureObject_t texObj = 0;
  			cudaCreateTextureObject(&texObj, &resDesc, &texDesc, NULL);
  
  			// specify  share memory
  			unsigned int share_bytes = (2*rmax + threads.x)*(1)*4;
  			
  			// allocate space on the GPU for the updated vote direction
  			T* gpuDir;
  			cudaMalloc(&gpuDir, bytes);	
  
  			//call the kernel to calculate the new voting direction
  			cuda_update_dir <<< blocks, threads, share_bytes >>>(gpuDir, texObj, gpuGrad, gpuTable, phi, rmax, x , y);
  
  			//call the kernel to update the gradient direction
  			cuda_update_grad <<< blocks, threads >>>(gpuGrad, gpuDir, x , y);
bf731970   Laila Saadatifard   fix some bugs in ...
163
  			
13fe3c84   Laila Saadatifard   update the stimli...
164
165
166
  			//free allocated memory
  			cudaFree(gpuDir);
  
59781ee3   Pavel Govyadinov   fixed a stask bug...
167
  			cudaDestroyTextureObject(texObj);
5de3a9c2   Pavel Govyadinov   CHECKPOINT: befo...
168
  			cudaFreeArray(cuArray);
59781ee3   Pavel Govyadinov   fixed a stask bug...
169
  
13fe3c84   Laila Saadatifard   update the stimli...
170
  		}
bf731970   Laila Saadatifard   fix some bugs in ...
171
  		
13fe3c84   Laila Saadatifard   update the stimli...
172
173
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
  		template<typename T>
  		void cpu_update_dir(T* cpuVote, T* cpuGrad,T* cpuTable, T phi, unsigned int rmax, unsigned int x, unsigned int y){
  
  			//calculate the number of bytes in the array
  			unsigned int bytes = x * y * sizeof(T);
  
  			//calculate the number of bytes in the atan2 table
  			unsigned int bytes_table = (2*rmax+1) * (2*rmax+1) * sizeof(T);
  
  			//allocate space on the GPU for the Vote Image
  			T* gpuVote;
  			cudaMalloc(&gpuVote, bytes);
  
  			//copy the input vote image to the GPU
  			HANDLE_ERROR(cudaMemcpy(gpuVote, cpuVote, bytes, cudaMemcpyHostToDevice));	
  
  			//allocate space on the GPU for the input Gradient image
  			T* gpuGrad;
  			HANDLE_ERROR(cudaMalloc(&gpuGrad, bytes*2));
  
  			//copy the Gradient data to the GPU
  			HANDLE_ERROR(cudaMemcpy(gpuGrad, cpuGrad, bytes*2, cudaMemcpyHostToDevice));
  
  			//allocate space on the GPU for the atan2 table
  			T* gpuTable;
  			HANDLE_ERROR(cudaMalloc(&gpuTable, bytes_table));
  
  			//copy the atan2 values to the GPU
  			HANDLE_ERROR(cudaMemcpy(gpuTable, cpuTable, bytes_table, cudaMemcpyHostToDevice));
bf731970   Laila Saadatifard   fix some bugs in ...
201
  						
13fe3c84   Laila Saadatifard   update the stimli...
202
203
  			//call the GPU version of the update direction function
  			gpu_update_dir<T>(gpuVote, gpuGrad, gpuTable, phi, rmax, x , y);
bf731970   Laila Saadatifard   fix some bugs in ...
204
  							
13fe3c84   Laila Saadatifard   update the stimli...
205
206
207
208
209
210
211
212
213
214
215
216
  			//copy the new gradient image back to the CPU
  			cudaMemcpy(cpuGrad, gpuGrad, bytes*2, cudaMemcpyDeviceToHost) ;
  
  			//free allocated memory
  			cudaFree(gpuTable);
  			cudaFree(gpuVote);
  			cudaFree(gpuGrad);
  		}
  		
  	}
  }
  
59781ee3   Pavel Govyadinov   fixed a stask bug...
217
  #endif