/tomo/pyhst

To get this branch, use:
bzr branch http://darksoft.org/webbzr/tomo/pyhst

« back to all changes in this revision

Viewing changes to hst_cuda/hst_cuda_kernels.h

  • Committer: Suren A. Chilingaryan
  • Date: 2017-09-28 10:34:47 UTC
  • Revision ID: csa@suren.me-20170928103447-dggjgnuxmymgew2a
Quality fixes (tex-based)

Show diffs side-by-side

added added

removed removed

Lines of Context:
97
97
    tmp1 = data[pos];
98
98
    tmp2 = data[pos + spitch];
99
99
 
 
100
#ifdef HST_FILTER2
100
101
    pos = 2 * (idy * dpitch + idx);
101
 
 
102
102
    out[pos] = tmp1;
103
103
    out[pos + 1] = tmp2;
 
104
#else
 
105
    pos = (2 * idy * dpitch + idx);
 
106
    out[pos] = tmp1;
 
107
    out[pos + dpitch] = tmp2;
 
108
#endif
104
109
}
105
110
 
106
111
__global__ static void  hst_cuda_pack_kernel_zpad(cufftReal *out, int dpitch, cufftReal *data, int spitch, int num_bins) {
121
126
        tmp2 = 0;
122
127
    }
123
128
 
 
129
#ifdef HST_FILTER2
124
130
    pos = 2 * (idy * dpitch + idx);
125
 
 
126
131
    out[pos] = tmp1;
127
132
    out[pos + 1] = tmp2;
 
133
#else
 
134
    pos = (2 * idy * dpitch + idx);
 
135
    out[pos] = tmp1;
 
136
    out[pos + dpitch] = tmp2;
 
137
#endif
128
138
}
129
139
 
130
140
__global__ static void  hst_cuda_pack_kernel_epad(cufftReal *out, int dpitch, cufftReal *data, int spitch, int num_bins, int mid) {
151
161
        tmp2 = data[pos + spitch];
152
162
    }
153
163
 
 
164
 
 
165
#ifdef HST_FILTER2
154
166
    pos = 2 * (idy * dpitch + idx);
155
 
 
156
167
    out[pos] = tmp1;
157
168
    out[pos + 1] = tmp2;
 
169
#else
 
170
    pos = (2 * idy * dpitch + idx);
 
171
    out[pos] = tmp1;
 
172
    out[pos + dpitch] = tmp2;
 
173
#endif
158
174
}
159
175
 
160
176
 
254
270
 
255
271
*/
256
272
 
257
 
__global__ static void  hst_cuda_unpack_kernel(cufftReal *out, int dpitch, cufftReal *data, int spitch) {
 
273
__global__ static void  hst_cuda_unpack_kernel(cufftReal *out, int dpitch, cufftReal *data, int spitch, float scaling) {
258
274
    float tmp1, tmp2;
259
275
 
260
276
    int  idx = blockIdx.x * BLOCK_SIZE + threadIdx.x;
261
277
    int  idy = blockIdx.y * BLOCK_SIZE + threadIdx.y;
262
278
 
 
279
#ifdef HST_FILTER2
263
280
    int pos = 2 * (idy * spitch + idx);
264
 
 
265
281
    tmp1 = data[pos];
266
282
    tmp2 = data[pos + 1];
 
283
#else
 
284
    int pos = (2 * idy * spitch + idx);
 
285
    tmp1 = data[pos];
 
286
    tmp2 = data[pos + spitch];
 
287
#endif
267
288
 
268
289
    pos = 2 * idy * dpitch + idx ;
269
290
 
270
 
    out[pos] = tmp1;
271
 
    out[pos + dpitch] = tmp2;
 
291
    out[pos] = tmp1 * scaling;
 
292
    out[pos + dpitch] = tmp2 * scaling;
272
293
}
273
294
 
274
295