/normxcorr/trunk

To get this branch, use:
bzr branch http://suren.me/webbzr/normxcorr/trunk

« back to all changes in this revision

Viewing changes to cuda/normxcorr_hw.cu

  • Committer: Suren A. Chilingaryan
  • Date: 2009-11-26 18:37:24 UTC
  • Revision ID: csa@dside.dyndns.org-20091126183724-dck2wcyt2ookjtab
Add error checking to initialization process

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
#if defined(_WIN32) || defined(_WIN64)
 
2
# include <windows.h>
 
3
    typedef UINT8 uint8_t;
 
4
    typedef UINT16 uint16_t;
 
5
    typedef UINT32 uint32_t;
 
6
    typedef INT8 int8_t;
 
7
    typedef INT16 int16_t;
 
8
    typedef INT32 int32_t;
 
9
#else
 
10
# include <stdint.h>
 
11
#endif
 
12
 
1
13
#include <stdio.h>
2
 
#include <stdint.h>
 
14
#include <stdlib.h>
3
15
 
4
16
#include <cuda.h>
5
17
#include <cuda_runtime.h>
24
36
    ACTION_COMPUTE_FRAGMENT = 11,
25
37
} TAction;
26
38
 
 
39
typedef enum {
 
40
    ERROR_CUFFT = 1,
 
41
    ERROR_CUDA_MALLOC = 2,
 
42
    ERROR_MALLOC = 3
 
43
} TError;
27
44
 
28
45
struct STProcessingState {
29
46
    cufftComplex *cuda_base_buffer;     // Stored FFT's of the template image
84
101
    }
85
102
}
86
103
 
87
 
static void fftInit(TProcessingState *ps) {
 
104
#include <unistd.h>
 
105
static int fftInit(TProcessingState *ps) {
 
106
    cufftResult cufft_err;
 
107
    cudaError cuda_err;
 
108
    
88
109
    fftFree(ps);
89
110
 
90
111
//    cublasInit();
91
112
//    cufftPlan2d(&ps->cufft_plan, ps->fft_size, ps->fft_size, CUFFT_C2C);
92
113
 
93
 
    cufftPlan2d(&ps->cufft_r2c_plan, ps->fft_size, ps->fft_size, CUFFT_R2C);
94
 
    cufftPlan2d(&ps->cufft_c2r_plan, ps->fft_size, ps->fft_size, CUFFT_C2R);
 
114
    cufft_err = cufftPlan2d(&ps->cufft_r2c_plan, ps->fft_size, ps->fft_size, CUFFT_R2C);
 
115
    if (cufft_err) {
 
116
        reportError("Problem initializing c2r plan, cufft code: %i", cufft_err);
 
117
        return ERROR_CUFFT;
 
118
    }   
 
119
    
 
120
    cufft_err = cufftPlan2d(&ps->cufft_c2r_plan, ps->fft_size, ps->fft_size, CUFFT_C2R);
 
121
    if (cufft_err) {
 
122
        reportError("Problem initializing r2c plan, cufft code: %i", cufft_err);
 
123
        cufftDestroy(ps->cufft_r2c_plan);
 
124
        return ERROR_CUFFT;
 
125
    }
95
126
 
96
127
    ps->fft_initialized = true;
97
128
 
98
 
    cudaMalloc((void**)&ps->cuda_base_buffer, ps->ncp * ps->fft_alloc_size * sizeof(cufftComplex));
99
 
 
100
 
    cudaMalloc((void**)&ps->cuda_data_buffer, ps->fft_alloc_size * sizeof(cufftComplex));
101
 
    cudaMalloc((void**)&ps->cuda_input_buffer, ps->ncp * ps->fft_alloc_size * sizeof(uint8_t));
102
 
    cudaMalloc((void**)&ps->cuda_final_buffer, ps->ncp * ps->fft_alloc_size * sizeof(float));
 
129
    cuda_err = cudaMalloc((void**)&ps->cuda_base_buffer, ps->ncp * ps->fft_alloc_size * sizeof(cufftComplex));
 
130
    if (cuda_err) {
 
131
        reportError("Device memory allocation of %u*%u*cufftComplex bytes for cuda_base_buffer is failed", ps->ncp, ps->fft_alloc_size);
 
132
        fftFree(ps);
 
133
        return ERROR_CUDA_MALLOC;
 
134
    }
 
135
 
 
136
    cuda_err = cudaMalloc((void**)&ps->cuda_data_buffer, ps->fft_alloc_size * sizeof(cufftComplex));
 
137
    if (cuda_err) {
 
138
        reportError("Device memory allocation of %u*cufftComplex bytes for cuda_data_buffer is failed", ps->fft_alloc_size);
 
139
        fftFree(ps);
 
140
        return ERROR_CUDA_MALLOC;
 
141
    }
 
142
 
 
143
    cuda_err = cudaMalloc((void**)&ps->cuda_input_buffer, ps->ncp * ps->fft_alloc_size * sizeof(uint8_t));
 
144
    if (cuda_err) {
 
145
        reportError("Device memory allocation of %u*%u*uint8 bytes for cuda_input_buffer is failed", ps->ncp, ps->fft_alloc_size);
 
146
        fftFree(ps);
 
147
        return ERROR_CUDA_MALLOC;
 
148
    }
 
149
 
 
150
    cuda_err = cudaMalloc((void**)&ps->cuda_final_buffer, ps->ncp * ps->fft_alloc_size * sizeof(float));
 
151
    if (cuda_err) {
 
152
        reportError("Device memory allocation of %u*%u*float bytes for cuda_final_buffer is failed", ps->ncp, ps->fft_alloc_size);
 
153
        fftFree(ps);
 
154
        return ERROR_CUDA_MALLOC;
 
155
    }
103
156
    cudaMemset((void*)ps->cuda_final_buffer, 0, ps->ncp * ps->fft_alloc_size * sizeof(float));
104
 
    cudaMalloc((void**)&ps->cuda_result_buffer, ps->fft_alloc_size * sizeof(cufftReal));
105
 
    cudaMalloc((void**)&ps->cuda_temp_buffer, ps->fft_alloc_size * sizeof(cufftReal));
 
157
 
 
158
    cuda_err = cudaMalloc((void**)&ps->cuda_result_buffer, ps->fft_alloc_size * sizeof(cufftReal));
 
159
    if (cuda_err) {
 
160
        reportError("Device memory allocation of %u*cufftReal bytes for cuda_result_buffer is failed", ps->fft_alloc_size);
 
161
        fftFree(ps);
 
162
        return ERROR_CUDA_MALLOC;
 
163
    }
 
164
 
 
165
    cuda_err = cudaMalloc((void**)&ps->cuda_temp_buffer, ps->fft_alloc_size * sizeof(cufftReal));
 
166
    if (cuda_err) {
 
167
        reportError("Device memory allocation of %u*cufftReal bytes for cuda_temp_buffer is failed", ps->fft_alloc_size);
 
168
        fftFree(ps);
 
169
        return ERROR_CUDA_MALLOC;
 
170
    }
106
171
    cudaMemset((void*)ps->cuda_temp_buffer, 0, ps->fft_alloc_size * sizeof(cufftReal));
107
172
 
108
 
    cudaMalloc((void**)&ps->cuda_lsum_buffer, ps->ncp * ps->fft_alloc_size * sizeof(float));
109
 
    cudaMalloc((void**)&ps->cuda_denom_buffer, ps->ncp * ps->fft_alloc_size * sizeof(float));
110
 
    cudaMalloc((void**)&ps->cuda_nonzero_buffer, ps->ncp * ps->fft_alloc_size * sizeof(uint16_t));
 
173
    cuda_err = cudaMalloc((void**)&ps->cuda_lsum_buffer, ps->ncp * ps->fft_alloc_size * sizeof(float));
 
174
    if (cuda_err) {
 
175
        reportError("Device memory allocation of %u*%u*float bytes for cuda_lsum_buffer is failed", ps->ncp, ps->fft_alloc_size);
 
176
        fftFree(ps);
 
177
        return ERROR_CUDA_MALLOC;
 
178
    }
 
179
 
 
180
    cuda_err = cudaMalloc((void**)&ps->cuda_denom_buffer, ps->ncp * ps->fft_alloc_size * sizeof(float));
 
181
    if (cuda_err) {
 
182
        reportError("Device memory allocation of %u*%u*float bytes for cuda_denom_buffer is failed", ps->ncp, ps->fft_alloc_size);
 
183
        fftFree(ps);
 
184
        return ERROR_CUDA_MALLOC;
 
185
    }
 
186
 
 
187
    cuda_err = cudaMalloc((void**)&ps->cuda_nonzero_buffer, ps->ncp * ps->fft_alloc_size * sizeof(uint16_t));
 
188
    if (cuda_err) {
 
189
        reportError("Device memory allocation of %u*%u*uint16 bytes for cuda_nonzero_buffer is failed", ps->ncp, ps->fft_alloc_size);
 
190
        fftFree(ps);
 
191
        return ERROR_CUDA_MALLOC;
 
192
    }
111
193
    cudaMemset((void*)ps->cuda_nonzero_buffer, 0, ps->ncp * ps->fft_alloc_size * sizeof(uint16_t));
112
194
    
113
195
    ps->cuda_nonzero_items = (uint16_t*)malloc(ps->ncp * sizeof(uint16_t));
 
196
    if (!ps->cuda_nonzero_items) {
 
197
        reportError("Host memory allocation of %u*uint16 bytes for cuda_nonzero_items is failed", ps->ncp);
 
198
        fftFree(ps);
 
199
        return ERROR_MALLOC;
 
200
    }
114
201
    
115
202
    ps->grid_size = (int*)malloc(ps->ncp*sizeof(int));
 
203
    if (!ps->grid_size) {
 
204
        reportError("Host memory allocation of %u*int bytes for grid_size is failed", ps->ncp);
 
205
        fftFree(ps);
 
206
        return ERROR_MALLOC;
 
207
    }
 
208
    
 
209
    return 0;
116
210
}
117
211
 
118
212
static void fftPrepare(TProcessingState *ps) {
168
262
    dataPtr = (uint8_t*)mxGetData(data);
169
263
    cudaMemcpy(cudaInputPtr, dataPtr, N2*sizeof(uint8_t), cudaMemcpyHostToDevice);
170
264
    vecPack<<<input_grid_dim, input_block_dim>>>(cudaRealPtr, size, cudaInputPtr, N);
 
265
 
171
266
    cufftExecR2C(ps->cufft_r2c_plan, cudaRealPtr, cudaPtr);
172
267
 
173
268
// Loading various stuff
213
308
    dataPtr = (uint8_t*)mxGetData(data);
214
309
    cudaMemcpy(cudaInputPtr, dataPtr, N2*sizeof(uint8_t), cudaMemcpyHostToDevice);
215
310
    vecPack<<<input_grid_dim, input_block_dim>>>(cudaRealPtr, size, cudaInputPtr, N);
 
311
 
216
312
    cufftExecR2C(ps->cufft_r2c_plan, cudaRealPtr, cudaDataPtr);
217
313
 
218
314
    vecMul<<<grid_dim,block_dim>>>(cudaDataPtr, cudaPtr, size/2+1);
271
367
 
272
368
 
273
369
void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[]) {
 
370
    int err;
 
371
    int64_t *errPtr;
 
372
    
274
373
    int deviceCount;
275
374
    cudaDeviceProp deviceProp;
276
375
    
386
485
 
387
486
    ps = pstate;
388
487
    
389
 
    action = (TAction)mxGetScalar((mxArray*)prhs[1]);
 
488
    action = (TAction)int(mxGetScalar((mxArray*)prhs[1]));
390
489
 
391
490
//    reportMessage("Executing normxcorr_hw action: %u", action);
392
491
 
492
591
            ps->fft_alloc_size = ps->fft_size2;
493
592
        }
494
593
 
495
 
        fftInit(ps);
 
594
        err = fftInit(ps);
 
595
 
 
596
        if (nlhs == 1) {
 
597
            idMatrix = mxCreateNumericMatrix(1, 1, mxINT64_CLASS, mxREAL);
 
598
            if (idMatrix) {
 
599
                errPtr = (int64_t*)mxGetData(idMatrix);
 
600
                errPtr[0] = err;
 
601
                plhs[0] = idMatrix;
 
602
            } else {
 
603
                reportError("Initialization of result matrix is failed");
 
604
                return;
 
605
            }
 
606
        }
496
607
     break;
497
608
     case ACTION_PREPARE:
498
609
        fftPrepare(ps);