cuUtils_template.Xcu 34.2 KB
Newer Older
Andreas Marek's avatar
Andreas Marek committed
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
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
#if 0
//    This file is part of ELPA.
//
//    The ELPA library was originally created by the ELPA consortium,
//    consisting of the following organizations:
//
//    - Max Planck Computing and Data Facility (MPCDF), formerly known as
//      Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
//    - Bergische Universität Wuppertal, Lehrstuhl für angewandte
//      Informatik,
//    - Technische Universität München, Lehrstuhl für Informatik mit
//      Schwerpunkt Wissenschaftliches Rechnen ,
//    - Fritz-Haber-Institut, Berlin, Abt. Theorie,
//    - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
//      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
//      and
//    - IBM Deutschland GmbH
//
//    This particular source code file contains additions, changes and
//    enhancements authored by Intel Corporation which is not part of
//    the ELPA consortium.
//
//    More information can be found here:
//    http://elpa.mpcdf.mpg.de/
//
//    ELPA is free software: you can redistribute it and/or modify
//    it under the terms of the version 3 of the license of the
//    GNU Lesser General Public License as published by the Free
//    Software Foundation.
//
//    ELPA is distributed in the hope that it will be useful,
//    but WITHOUT ANY WARRANTY; without even the implied warranty of
//    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
//    GNU Lesser General Public License for more details.
//
//    You should have received a copy of the GNU Lesser General Public License
//    along with ELPA.  If not, see <http://www.gnu.org/licenses/>
//
//    ELPA reflects a substantial effort on the part of the original
//    ELPA consortium, and we ask you to respect the spirit of the
//    license that we chose: i.e., please contribute any changes you
//    may have back to the original ELPA library distribution, and keep
//    any derivatives of ELPA under the same license that we chose for
//    the original distribution, the GNU Lesser General Public License.
//
//
// --------------------------------------------------------------------------------------------------
//
// This file was originally written by NVIDIA
// and re-written by A. Marek, MPCDF
#endif


#include "config-f90.h"

#include <cuda_runtime.h>
#include <stdlib.h>
#include <stdio.h>

#if COMPLEXCASE == 1
Andreas Marek's avatar
Andreas Marek committed
61
#include <cuComplex.h>
Andreas Marek's avatar
Andreas Marek committed
62
63
#endif

64
65
66
67
#define BLOCK_CYCLIC_BLOCKSIZE 128
#define GLOBAL_STRIPE_WIDTH 256
#define WARP_SIZE 32

Andreas Marek's avatar
Andreas Marek committed
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
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
// Reset a reduction block
// Limitation: the thread-block size must be a divider of the reduction block's size

#if REALCASE == 1

#ifdef DOUBLE_PRECISION_REAL
__device__ void reset_shared_block_c_double ( double * s_block, int b_size)
#else
__device__ void reset_shared_block_c_single ( float * s_block, int b_size)
#endif

#endif

#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
__device__ void reset_shared_block_c_complex_double ( cuDoubleComplex * s_block, int b_size)
#else
__device__ void reset_shared_block_c_complex_single ( cuFloatComplex * s_block, int b_size)
#endif
#endif

{
    int i, t_idx, s_chunk ;
    t_idx = threadIdx.x;
    s_chunk = b_size / blockDim.x;
    for(i = ((t_idx - 1) * s_chunk + 1) ; i < (t_idx * s_chunk); i++) {
#if  REALCASE == 1
      s_block[i] = 0.0 ;
#endif
#if COMPLEXCASE == 1
      s_block[i].x = 0.0 ;
      s_block[i].y = 0.0 ;
#endif
    }
    __syncthreads();
}

// Reset 2 reduction blocks without an explicit synchronization at the end
// Limitation: : the thread-block size must be a divider of the reduction block's size
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
__device__ void reset_shared_block_pair_c_double( double *s_block_1, double *s_block_2, int b_size)
#else
__device__ void reset_shared_block_pair_c_single( float *s_block_1, float *s_block_2, int b_size)
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
__device__ void reset_shared_block_pair_c_complex_double( cuDoubleComplex *s_block_1, cuDoubleComplex *s_block_2, int b_size)
#else
__device__ void reset_shared_block_pair_c_complex_single( cuFloatComplex *s_block_1, cuFloatComplex *s_block_2, int b_size)
#endif
#endif
{
    int i, t_idx, s_chunk;

    t_idx = threadIdx.x;
    s_chunk = b_size / blockDim.x;
    for(i = ((t_idx - 1) * s_chunk + 1); i < (t_idx * s_chunk); i++)
    {
#if REALCASE == 1
        s_block_1[i] = 0.0 ;
        s_block_2[i] = 0.0 ;
#endif
#if COMPLEXCASE == 1
        s_block_1[i].x = 0.0 ;
        s_block_2[i].x= 0.0 ;
        s_block_1[i].y = 0.0 ;
        s_block_2[i].y= 0.0 ;
#endif
    }
}
// Reset a reduction block
// Limitation: the thread-block size must be a divider of the reduction block's size

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
__device__ void warp_reduce_c_double( double *s_block)
#else
__device__ void warp_reduce_c_single( float *s_block)
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
__device__ void warp_reduce_complex_double( cuDoubleComplex *s_block)
#else
__device__ void warp_reduce_complex_single( cuFloatComplex *s_block)
#endif
#endif
{
    int t_idx ;
    t_idx = threadIdx.x;
    __syncthreads();

#if REALCASE == 1
163
        // attention
Andreas Marek's avatar
Andreas Marek committed
164
165
166
167
168
169
170
171
172
173
174
175
        if (t_idx < 32)
	{
                s_block[t_idx] = s_block[t_idx] + s_block[t_idx + 32] + s_block[t_idx + 64] + s_block[t_idx + 96] ;
        if (t_idx < 8)
                s_block[t_idx] = s_block[t_idx] + s_block[t_idx + 8] + s_block[t_idx + 16] + s_block[t_idx + 24];
        if (t_idx < 4)
                s_block[t_idx] = s_block[t_idx] + s_block[t_idx + 4];
        if (t_idx < 1)
                s_block[t_idx] = s_block[t_idx] + s_block[t_idx + 1] + s_block[t_idx + 2] + s_block[t_idx + 3];
	}
#endif
#if COMPLEXCASE == 1
176
        // attention
Andreas Marek's avatar
Andreas Marek committed
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
205
206
207
208
209
210
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
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
	if (t_idx < 32)
        {
#ifdef DOUBLE_PRECISION_COMPLEX
        s_block[t_idx] = cuCadd(cuCadd(s_block[t_idx],s_block[t_idx + 32]) , cuCadd( s_block[t_idx + 64], s_block[t_idx + 96]) );
        if (t_idx < 8)
        {
        s_block[t_idx] = cuCadd(cuCadd(s_block[t_idx],s_block[t_idx + 8] ) , cuCadd( s_block[t_idx + 16] , s_block[t_idx + 24] ) );

        }
        if (t_idx < 4)
        {
        s_block[t_idx] = cuCadd(s_block[t_idx] , s_block[t_idx + 4]) ;
        }
        if (t_idx < 1)
        {
        s_block[t_idx] = cuCadd(cuCadd(s_block[t_idx],s_block[t_idx + 1] ) , cuCadd( s_block[t_idx +2] , s_block[t_idx + 3] ) );
        }
        }
#else
        s_block[t_idx] = cuCaddf(cuCaddf(s_block[t_idx],s_block[t_idx + 32]) , cuCaddf( s_block[t_idx + 64], s_block[t_idx + 96]) );
        if (t_idx < 8)
        {
        s_block[t_idx] = cuCaddf(cuCaddf(s_block[t_idx],s_block[t_idx + 8] ) , cuCaddf( s_block[t_idx + 16] , s_block[t_idx + 24] ) );

        }
        if (t_idx < 4)
        {
        s_block[t_idx] = cuCaddf(s_block[t_idx] , s_block[t_idx + 4]) ;
        }
        if (t_idx < 1)
        {
        s_block[t_idx] = cuCaddf(cuCaddf(s_block[t_idx],s_block[t_idx + 1] ) , cuCaddf( s_block[t_idx +2] , s_block[t_idx + 3] ) );
        }
        }
#endif
#endif /* COMPLEXCASE == 1 */
}

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
__global__ void my_pack_c_kernel_double(const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, double* src, double* dst)
#else
__global__ void my_pack_c_kernel_single(const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, float* src, float* dst)
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
__global__ void my_pack_c_kernel_complex_double(const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, cuDoubleComplex* src, cuDoubleComplex* dst)
#else
__global__ void my_pack_c_kernel_complex_single(const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, cuFloatComplex* src, cuFloatComplex* dst)
#endif
#endif
{
    int b_id, t_id ;
    int dst_ind ;
    b_id = blockIdx.y;
    t_id = threadIdx.x;

    dst_ind = b_id * stripe_width + t_id;
    if (dst_ind < max_idx)
    {
	// dimension of dst - lnev, nblk
	// dimension of src - stripe_width,a_dim2,stripe_count
#if REALCASE == 1
    	*(dst + dst_ind + (l_nev*blockIdx.x) ) = *(src + t_id + (stripe_width*(n_offset + blockIdx.x)) + ( b_id *stripe_width*a_dim2 ));
#endif
#if COMPLEXCASE == 1
	dst[dst_ind + (l_nev*blockIdx.x)].x = src[t_id + (stripe_width*(n_offset + blockIdx.x)) + ( b_id *stripe_width*a_dim2)].x;
        dst[dst_ind + (l_nev*blockIdx.x)].y = src[t_id + (stripe_width*(n_offset + blockIdx.x)) + ( b_id *stripe_width*a_dim2)].y;
#endif
     }

}

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
__global__ void  my_unpack_c_kernel_double( const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, double* src, double* dst)
#else
__global__ void  my_unpack_c_kernel_single( const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, float* src, float* dst)
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
__global__ void  my_unpack_c_kernel_complex_double( const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, cuDoubleComplex* src, cuDoubleComplex* dst)
#else
__global__ void  my_unpack_c_kernel_complex_single( const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, cuFloatComplex* src, cuFloatComplex* dst)
#endif
#endif
{
    int b_id, t_id ;
    int src_ind;

    b_id = blockIdx.y;
    t_id = threadIdx.x;

    src_ind = b_id * stripe_width + t_id;
    if (src_ind < max_idx)
#if REALCASE == 1
	*(dst + (t_id + ((n_offset + blockIdx.x) * stripe_width) + (b_id * stripe_width * a_dim2 ))) = *(src + src_ind  + (blockIdx.x) *l_nev );
#endif
#if COMPLEXCASE == 1
        dst[ t_id + ((n_offset + blockIdx.x) * stripe_width) + (b_id * stripe_width * a_dim2 )].x = src[ src_ind  + (blockIdx.x) *l_nev ].x;
	dst[ t_id + ((n_offset + blockIdx.x) * stripe_width) + (b_id * stripe_width * a_dim2 )].y = src[ src_ind  + (blockIdx.x) *l_nev ].y;
#endif

}

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
__global__ void extract_hh_tau_c_kernel_double(double* hh, double* hh_tau, const int nbw, const int n, int val)
#else
__global__ void extract_hh_tau_c_kernel_single(float* hh, float* hh_tau, const int nbw, const int n, int val)
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
__global__ void extract_hh_tau_c_kernel_complex_double(cuDoubleComplex* hh, cuDoubleComplex* hh_tau, const int nbw, const int n, int val)
#else
__global__ void extract_hh_tau_c_kernel_complex_single(cuFloatComplex* hh, cuFloatComplex* hh_tau, const int nbw, const int n, int val)
#endif
#endif
{
    int h_idx ;
    h_idx = (blockIdx.x) * blockDim.x + threadIdx.x;

    if (h_idx < n)
    {
	//dimension of hh - (nbw, max_blk_size)
	//dimension of hh_tau - max_blk_size
#if REALCASE == 1
        *(hh_tau + h_idx ) = *(hh +  (h_idx * nbw)) ;
#endif
#if COMPLEXCASE == 1
        hh_tau[h_idx] = hh[h_idx * nbw] ;
#endif
        //  Replace the first element in the HH reflector with 1.0 or 0.0
#if REALCASE == 1
	if( val == 0)
        *(hh + (h_idx * nbw)) = 1.0;
	else
	*(hh + (h_idx * nbw)) = 0.0;
#endif
#if COMPLEXCASE == 1
        if( val == 0)
        {
         hh[(h_idx * nbw)].x = 1.0;
	 hh[h_idx *nbw].y= 0.0;
        }
        else
        {
        hh[(h_idx * nbw)].x = 0.0;
	hh[h_idx*nbw].y =0.0;
        }
#endif
     }
}

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
__global__ void  compute_hh_dotp_c_kernel_double(double* hh, double* v_dot, const int nbw, const int n)
{

339
   __shared__ double hh_s[BLOCK_CYCLIC_BLOCKSIZE] ;
Andreas Marek's avatar
Andreas Marek committed
340
341
342
343
#else
__global__ void  compute_hh_dotp_c_kernel_single(float* hh, float* v_dot, const int nbw, const int n)
{

344
   __shared__ float hh_s[BLOCK_CYCLIC_BLOCKSIZE] ;
Andreas Marek's avatar
Andreas Marek committed
345
346
347
348
349
350
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
__global__ void  compute_hh_dotp_c_kernel_complex_double(cuDoubleComplex* hh, cuDoubleComplex* v_dot, const int nbw, const int n)
{
351
   __shared__ cuDoubleComplex hh_s[BLOCK_CYCLIC_BLOCKSIZE] ;
Andreas Marek's avatar
Andreas Marek committed
352
353
354
355

#else
__global__ void  compute_hh_dotp_c_kernel_complex_single(cuFloatComplex* hh, cuFloatComplex* v_dot, const int nbw, const int n)
{
356
   __shared__ cuFloatComplex hh_s[BLOCK_CYCLIC_BLOCKSIZE] ;
Andreas Marek's avatar
Andreas Marek committed
357
358
359
360
361
362
363
364
365
366
367
#endif
#endif
    int t_idx, v_idx;

    //  The vector index (v_idx) identifies the pair of HH reflectors from which the dot product is computed
    v_idx = blockIdx.x  ;

    //  The thread index indicates the position within the two HH reflectors
    t_idx = threadIdx.x ;

//    //  The contents of the shared memory must be fully reset
368
//     reset_shared_block_c(hh_s, BLOCK_CYCLIC_BLOCKSIZE);
Andreas Marek's avatar
Andreas Marek committed
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521

    //  Initialize the contents of the shared buffer (preparing for reduction)
    if (t_idx  > 0)
    {
#if REALCASE == 1
        *(hh_s + t_idx) = *(hh + t_idx + v_idx * nbw ) *  (*(hh + (t_idx - 1) +  (v_idx +1)* nbw)) ;
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
       hh_s[t_idx] = cuCmul(cuConj(hh[t_idx + v_idx * nbw]),   hh[ (t_idx - 1) +  (v_idx +1)* nbw]) ;
#else
       hh_s[t_idx] = cuCmulf(cuConjf(hh[t_idx + v_idx * nbw]),   hh[ (t_idx - 1) +  (v_idx +1)* nbw]) ;
#endif
#endif
    }
    else
    {
#if REALCASE == 1
        *(hh_s + t_idx) = 0.0 ;
#endif
#if COMPLEXCASE == 1
        hh_s[t_idx].x = 0.0 ;
        hh_s[t_idx].y = 0.0;
#endif
    }
     //  Compute the dot product using a fast reduction

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
     warp_reduce_c_double(hh_s);
#else
     warp_reduce_c_single(hh_s);
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
     warp_reduce_complex_double(hh_s);
#else
     warp_reduce_complex_single(hh_s);
#endif
     __syncthreads();
#endif

      if(t_idx == 0)
      {
#if REALCASE == 1
      *(v_dot + v_idx) = *(hh_s) ;
#endif
#if COMPLEXCASE == 1
      v_dot[v_idx] = hh_s[0] ;
#endif
      }

}

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
extern "C" void launch_my_pack_c_kernel_double(const int row_count, const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, double* a_dev, double* row_group_dev)
#else
extern "C" void launch_my_pack_c_kernel_single(const int row_count, const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, float* a_dev, float* row_group_dev)
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_my_pack_c_kernel_complex_double(const int row_count, const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, cuDoubleComplex* a_dev, cuDoubleComplex* row_group_dev)
#else
extern "C" void launch_my_pack_c_kernel_complex_single(const int row_count, const int n_offset, const int max_idx, const int stripe_width, const int a_dim2, const int stripe_count, const int l_nev, cuFloatComplex* a_dev, cuFloatComplex* row_group_dev)
#endif
#endif
{

	dim3  grid_size;
        grid_size = dim3(row_count, stripe_count, 1);
        cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to mypack kernel: %s, %d\n",cudaGetErrorString(err), err);
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
	my_pack_c_kernel_double<<<grid_size, stripe_width>>>(n_offset, max_idx, stripe_width, a_dim2, stripe_count, l_nev, a_dev, row_group_dev);
#else
	my_pack_c_kernel_single<<<grid_size, stripe_width>>>(n_offset, max_idx, stripe_width, a_dim2, stripe_count, l_nev, a_dev, row_group_dev);
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
        my_pack_c_kernel_complex_double<<<grid_size, stripe_width>>>(n_offset, max_idx, stripe_width, a_dim2, stripe_count, l_nev, a_dev, row_group_dev);
#else
        my_pack_c_kernel_complex_single<<<grid_size, stripe_width>>>(n_offset, max_idx, stripe_width, a_dim2, stripe_count, l_nev, a_dev, row_group_dev);
#endif
#endif
	 err = cudaGetLastError();
        if ( err!= cudaSuccess)
        {
                printf("\n my pack_kernel failed  %s \n",cudaGetErrorString(err) );
        }

}
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
extern "C" void launch_compute_hh_dotp_c_kernel_double(double* bcast_buffer_dev, double* hh_dot_dev,const int nbw,const int n)
#else
extern "C" void launch_compute_hh_dotp_c_kernel_single(float* bcast_buffer_dev, float* hh_dot_dev,const int nbw,const int n)
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_compute_hh_dotp_c_kernel_complex_double(cuDoubleComplex* bcast_buffer_dev, cuDoubleComplex* hh_dot_dev,const int nbw,const int n)
#else
extern "C" void launch_compute_hh_dotp_c_kernel_complex_single(cuFloatComplex* bcast_buffer_dev, cuFloatComplex* hh_dot_dev,const int nbw,const int n)
#endif
#endif
{
	cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to compute_hh kernel: %s, %d\n",cudaGetErrorString(err), err);

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
        compute_hh_dotp_c_kernel_double<<< n-1, nbw >>>(bcast_buffer_dev, hh_dot_dev, nbw, n);
#else
        compute_hh_dotp_c_kernel_single<<< n-1, nbw >>>(bcast_buffer_dev, hh_dot_dev, nbw, n);
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
        compute_hh_dotp_c_kernel_complex_double<<< n-1, nbw >>>(bcast_buffer_dev, hh_dot_dev, nbw, n);
#else
        compute_hh_dotp_c_kernel_complex_single<<< n-1, nbw >>>(bcast_buffer_dev, hh_dot_dev, nbw, n);
#endif
#endif
	err = cudaGetLastError();
        if ( err!= cudaSuccess)
        {
                printf("\n compute _kernel failed  %s \n",cudaGetErrorString(err) );
        }

}
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
extern "C" void launch_extract_hh_tau_c_kernel_double(double* bcast_buffer_dev, double* hh_tau_dev, const int nbw, const int n , const int is_zero)
#else
extern "C" void launch_extract_hh_tau_c_kernel_single(float* bcast_buffer_dev, float* hh_tau_dev, const int nbw, const int n , const int is_zero)
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_extract_hh_tau_c_kernel_complex_double(cuDoubleComplex* bcast_buffer_dev, cuDoubleComplex* hh_tau_dev, const int nbw, const int n , const int is_zero)
#else
extern "C" void launch_extract_hh_tau_c_kernel_complex_single(cuFloatComplex* bcast_buffer_dev, cuFloatComplex* hh_tau_dev, const int nbw, const int n , const int is_zero)
#endif
#endif
{
	int grid_size;
522
	grid_size = 1 + (n - 1) / GLOBAL_STRIPE_WIDTH;
Andreas Marek's avatar
Andreas Marek committed
523
524
525
526
527
	cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to extract kernel: %s, %d\n",cudaGetErrorString(err), err);
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
528
	extract_hh_tau_c_kernel_double<<<grid_size,GLOBAL_STRIPE_WIDTH>>>(bcast_buffer_dev,hh_tau_dev, nbw, n, is_zero);
Andreas Marek's avatar
Andreas Marek committed
529
#else
530
	extract_hh_tau_c_kernel_single<<<grid_size,GLOBAL_STRIPE_WIDTH>>>(bcast_buffer_dev,hh_tau_dev, nbw, n, is_zero);
Andreas Marek's avatar
Andreas Marek committed
531
532
533
534
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
535
        extract_hh_tau_c_kernel_complex_double<<<grid_size,GLOBAL_STRIPE_WIDTH>>>(bcast_buffer_dev,hh_tau_dev, nbw, n, is_zero);
Andreas Marek's avatar
Andreas Marek committed
536
#else
537
        extract_hh_tau_c_kernel_complex_single<<<grid_size,GLOBAL_STRIPE_WIDTH>>>(bcast_buffer_dev,hh_tau_dev, nbw, n, is_zero);
Andreas Marek's avatar
Andreas Marek committed
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
#endif
#endif
	err = cudaGetLastError();
	if ( err!= cudaSuccess)
       	{
		printf("\n  extract _kernel failed  %s \n",cudaGetErrorString(err) );
        }

}

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
extern "C" void launch_my_unpack_c_kernel_double( const int row_count, const int n_offset, const int max_idx, const int stripe_width,const int a_dim2, const int stripe_count, const int l_nev, double* row_group_dev, double* a_dev)
#else
extern "C" void launch_my_unpack_c_kernel_single( const int row_count, const int n_offset, const int max_idx, const int stripe_width,const int a_dim2, const int stripe_count, const int l_nev, float* row_group_dev, float* a_dev)
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_my_unpack_c_kernel_complex_double( const int row_count, const int n_offset, const int max_idx, const int stripe_width,const int a_dim2, const int stripe_count, const int l_nev, cuDoubleComplex* row_group_dev, cuDoubleComplex* a_dev)
#else
extern "C" void launch_my_unpack_c_kernel_complex_single( const int row_count, const int n_offset, const int max_idx, const int stripe_width,const int a_dim2, const int stripe_count, const int l_nev, cuFloatComplex* row_group_dev, cuFloatComplex* a_dev)
#endif
#endif
{
        dim3  grid_size;
	grid_size = dim3(row_count, stripe_count, 1);
        cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to unpack kernel: %s, %d\n",cudaGetErrorString(err), err);
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
        my_unpack_c_kernel_double<<<grid_size, stripe_width>>>(n_offset, max_idx, stripe_width, a_dim2, stripe_count, l_nev, row_group_dev , a_dev);
#else
        my_unpack_c_kernel_single<<<grid_size, stripe_width>>>(n_offset, max_idx, stripe_width, a_dim2, stripe_count, l_nev, row_group_dev , a_dev);
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
        my_unpack_c_kernel_complex_double<<<grid_size, stripe_width>>>(n_offset, max_idx, stripe_width, a_dim2, stripe_count, l_nev, row_group_dev , a_dev);
#else
        my_unpack_c_kernel_complex_single<<<grid_size, stripe_width>>>(n_offset, max_idx, stripe_width, a_dim2, stripe_count, l_nev, row_group_dev , a_dev);
#endif
#endif
        err = cudaGetLastError();
        if ( err != cudaSuccess)
        {
	    printf("\n  my_unpack_c_kernel failed  %s \n",cudaGetErrorString(err) );
        }
}

#if COMPLEXCASE == 1

#ifdef DOUBLE_PRECISION_COMPLEX
__global__ void compute_kernel_reduce_double( cuDoubleComplex* a_dev, int lda , int n ,int nbw ,  cuDoubleComplex *h1_dev )
#else
__global__ void compute_kernel_reduce_single( cuFloatComplex* a_dev, int lda , int n ,int nbw ,  cuFloatComplex *h1_dev )
#endif
{
    int  t_id ;
    int st_ind;

    t_id = threadIdx.x;

    st_ind = (t_id*(t_id+1))/2;
    if(t_id< n)
    {
	for(int i =0;i<=t_id;i++)
        {
         h1_dev[st_ind + i] = a_dev[t_id *lda + i ] ;
	}
    }
    __syncthreads();
}

#ifdef DOUBLE_PRECISION_COMPLEX
__global__ void compute_kernel_reduce_1_double( cuDoubleComplex* a_dev, int lda , int n, cuDoubleComplex *h1_dev )
#else
__global__ void compute_kernel_reduce_1_single( cuFloatComplex* a_dev, int lda , int n, cuFloatComplex *h1_dev )
#endif
{
    int  t_id ;
    int st_ind;

    t_id = threadIdx.x;

    st_ind = (t_id*(t_id+1))/2;
    if(t_id< n)
    {
        for(int i =0;i<=t_id;i++)
         {
	  a_dev[t_id *lda + i ] = h1_dev[st_ind + i];
#ifdef DOUBLE_PRECISION_COMPLEX
	  a_dev[ (i-1)*lda + t_id ] = cuConj(a_dev[ t_id *lda + i-1]) ;
#else
	  a_dev[ (i-1)*lda + t_id ] = cuConjf(a_dev[ t_id *lda + i-1]) ;
#endif
	}
    }
    __syncthreads();
}

#ifdef DOUBLE_PRECISION_COMPLEX
__global__ void  dot_product_c_kernel_double( cuDoubleComplex* hs_dev, cuDoubleComplex* hv_new_dev, cuDoubleComplex tau_new_dev, cuDoubleComplex*  x_dev, cuDoubleComplex *h_dev, cuDoubleComplex *hv_dev, int nr)
#else
__global__ void  dot_product_c_kernel_single( cuFloatComplex* hs_dev, cuFloatComplex* hv_new_dev, cuFloatComplex tau_new_dev, cuFloatComplex*  x_dev, cuFloatComplex *h_dev, cuFloatComplex *hv_dev, int nr)
#endif
{
    int t_id ;

#ifdef DOUBLE_PRECISION_COMPLEX
    __shared__ cuDoubleComplex x_dev_temp[128];
    __shared__ cuDoubleComplex x_val;
#else
    __shared__ cuFloatComplex x_dev_temp[128];
    __shared__ cuFloatComplex x_val;
#endif
    //b_id = blockIdx.y;
    t_id = threadIdx.x;

    if(t_id<nr)
#ifdef DOUBLE_PRECISION_COMPLEX
	 x_dev_temp[t_id] = cuCmul( cuConj(hs_dev[t_id]), hv_new_dev[t_id]) ;
#else
	 x_dev_temp[t_id] = cuCmulf( cuConjf(hs_dev[t_id]), hv_new_dev[t_id]) ;
#endif
    __syncthreads();

    if(t_id==0)
    {
        for(int i=1;i<nr;i++)
#ifdef DOUBLE_PRECISION_COMPLEX
	x_dev_temp[t_id] = cuCadd(x_dev_temp[t_id],x_dev_temp[t_id +i]);
#else
	x_dev_temp[t_id] = cuCaddf(x_dev_temp[t_id],x_dev_temp[t_id +i]);
#endif
    }
    __syncthreads();
     if(t_id ==0)
    {
#ifdef DOUBLE_PRECISION_COMPLEX
      x_val =  cuCmul(x_dev_temp[t_id], tau_new_dev);
#else
      x_val =  cuCmulf(x_dev_temp[t_id], tau_new_dev);
#endif
      x_dev[0] = x_val;
    }
	__syncthreads();
}

#ifdef DOUBLE_PRECISION_COMPLEX
__global__ void  dot_product_c_kernel_1_double(   cuDoubleComplex*  ab_dev, cuDoubleComplex *hs_dev,  cuDoubleComplex*  hv_new_dev, cuDoubleComplex*  x_dev, cuDoubleComplex *h_dev, cuDoubleComplex *hv_dev,  int nb, int nr , int ns )
#else
__global__ void  dot_product_c_kernel_1_single(   cuFloatComplex*  ab_dev, cuFloatComplex *hs_dev,  cuFloatComplex*  hv_new_dev, cuFloatComplex*  x_dev, cuFloatComplex *h_dev, cuFloatComplex *hv_dev,  int nb, int nr , int ns )
#endif
{
    int t_id = threadIdx.x;
        int i;

    if((t_id>0 )&& (t_id < nb))
    {
#ifdef DOUBLE_PRECISION_COMPLEX
	h_dev[t_id] = cuCsub(h_dev[t_id], cuCmul(x_dev[0],hv_dev[t_id]));
#else
	h_dev[t_id] = cuCsubf(h_dev[t_id], cuCmulf(x_dev[0],hv_dev[t_id]));
#endif
        for(i=0;i<nr;i++)
	{
#ifdef DOUBLE_PRECISION_COMPLEX
	 ab_dev[ i+nb-t_id + (t_id+ns-1)*2*nb ] = cuCsub(cuCsub(ab_dev[ i+nb-t_id + (t_id+ns-1)*2*nb],cuCmul(hv_new_dev[i],cuConj(h_dev[t_id])) ),cuCmul(hs_dev[i], cuConj(hv_dev[t_id])));
#else
	 ab_dev[ i+nb-t_id + (t_id+ns-1)*2*nb ] = cuCsubf(cuCsubf(ab_dev[ i+nb-t_id + (t_id+ns-1)*2*nb],cuCmulf(hv_new_dev[i],cuConjf(h_dev[t_id])) ),cuCmulf(hs_dev[i], cuConjf(hv_dev[t_id])));
#endif
 	}
    }
   __syncthreads();
}

#ifdef DOUBLE_PRECISION_COMPLEX
__global__ void  double_hh_transform_kernel_double( cuDoubleComplex*  ab_dev, cuDoubleComplex *hs_dev, cuDoubleComplex *hv_dev,  int nb,  int ns )
#else
__global__ void  double_hh_transform_kernel_single( cuFloatComplex*  ab_dev, cuFloatComplex *hs_dev, cuFloatComplex *hv_dev,  int nb,  int ns )
#endif
{
    int t_id = threadIdx.x;
    if((t_id>0 )&& (t_id < nb))
    {
#ifdef DOUBLE_PRECISION_COMPLEX
         ab_dev[ nb-t_id + (t_id+ns-1)*2*nb ] = cuCsub(ab_dev[ nb-t_id + (t_id+ns-1)*2*nb],cuCmul(hs_dev[0], cuConj(hv_dev[t_id])));
#else
         ab_dev[ nb-t_id + (t_id+ns-1)*2*nb ] = cuCsubf(ab_dev[ nb-t_id + (t_id+ns-1)*2*nb],cuCmulf(hs_dev[0], cuConjf(hv_dev[t_id])));
#endif
    }
   __syncthreads();
}

#ifdef DOUBLE_PRECISION_COMPLEX
__global__ void  double_hh_transform_kernel_2_double( cuDoubleComplex*  ab_dev, cuDoubleComplex *hd_dev, cuDoubleComplex *hv_dev,  int nc,  int ns , int nb )
#else
__global__ void  double_hh_transform_kernel_2_single( cuFloatComplex*  ab_dev, cuFloatComplex *hd_dev, cuFloatComplex *hv_dev,  int nc,  int ns , int nb )
#endif
{
    int t_id = threadIdx.x;
    if(t_id < nc)
    {
#ifdef DOUBLE_PRECISION_COMPLEX
         ab_dev[ t_id + (ns-1)*2*nb ] = cuCsub(cuCsub(ab_dev[ t_id + (ns-1)*2*nb],cuCmul(hd_dev[ t_id], cuConj(hv_dev[0]))) , cuCmul(hv_dev[ t_id], cuConj(hd_dev[0])));
#else
         ab_dev[ t_id + (ns-1)*2*nb ] = cuCsubf(cuCsubf(ab_dev[ t_id + (ns-1)*2*nb],cuCmulf(hd_dev[ t_id], cuConjf(hv_dev[0]))) , cuCmulf(hv_dev[ t_id], cuConjf(hd_dev[0])));

#endif
    }
   __syncthreads();
}

#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_dot_product_kernel_double( cuDoubleComplex* hs_dev, cuDoubleComplex* hv_new_dev, cuDoubleComplex tau_new_dev, cuDoubleComplex*  x_dev, cuDoubleComplex*  h_dev ,cuDoubleComplex*  hv_dev,int nr )
#else
extern "C" void launch_dot_product_kernel_single( cuFloatComplex* hs_dev, cuFloatComplex* hv_new_dev, cuFloatComplex tau_new_dev, cuFloatComplex*  x_dev, cuFloatComplex*  h_dev ,cuFloatComplex*  hv_dev,int nr )
#endif
{
        dim3  grid_size;
        grid_size = dim3(1,1, 1);
        cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to launch_dot_product kernel: %s, %d\n",cudaGetErrorString(err), err);
#ifdef DOUBLE_PRECISION_COMPLEX
        dot_product_c_kernel_double<<<grid_size, nr>>>(hs_dev, hv_new_dev, tau_new_dev, x_dev, h_dev, hv_dev, nr );
#else
        dot_product_c_kernel_single<<<grid_size, nr>>>(hs_dev, hv_new_dev, tau_new_dev, x_dev, h_dev, hv_dev, nr );
#endif
        err = cudaGetLastError();
        if ( err != cudaSuccess)
        {
            printf("\n dot product kernel failed  %s \n",cudaGetErrorString(err) );
        }
}

#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_dot_product_kernel_1_double(  cuDoubleComplex*  ab_dev, cuDoubleComplex *hs_dev,  cuDoubleComplex*  hv_new_dev,cuDoubleComplex*  x_dev, cuDoubleComplex*  h_dev ,cuDoubleComplex*  hv_dev, int nb ,int nr , int ns)
#else
extern "C" void launch_dot_product_kernel_1_single(  cuFloatComplex*  ab_dev, cuFloatComplex *hs_dev,  cuFloatComplex*  hv_new_dev,cuFloatComplex*  x_dev, cuFloatComplex*  h_dev ,cuFloatComplex*  hv_dev, int nb ,int nr , int ns)
#endif
{
        dim3  grid_size;
        grid_size = dim3(1,1, 1);
        cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to launch_dot_product kernel: %s, %d\n",cudaGetErrorString(err), err);
#ifdef DOUBLE_PRECISION_COMPLEX
        dot_product_c_kernel_1_double<<<grid_size, nb>>>( ab_dev, hs_dev, hv_new_dev, x_dev, h_dev, hv_dev, nb, nr, ns );
#else
        dot_product_c_kernel_1_single<<<grid_size, nb>>>( ab_dev, hs_dev, hv_new_dev, x_dev, h_dev, hv_dev, nb, nr, ns );
#endif
        err = cudaGetLastError();
        if ( err != cudaSuccess)
        {
            printf("\n dot product kernel failed  %s \n",cudaGetErrorString(err) );
        }
}

#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_dot_product_kernel_2_double(  cuDoubleComplex*  ab_dev, cuDoubleComplex *hs_dev,  cuDoubleComplex*  hv_dev,cuDoubleComplex*  hd_dev, int nb ,int nr , int ne)
#else
extern "C" void launch_dot_product_kernel_2_single(  cuFloatComplex*  ab_dev, cuFloatComplex *hs_dev,  cuFloatComplex*  hv_dev,cuFloatComplex*  hd_dev, int nb ,int nr , int ne)
#endif
{
        dim3  grid_size;
        grid_size = dim3(1,1, 1);
        cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to launch_dot_product kernel: %s, %d\n",cudaGetErrorString(err), err);
        err = cudaGetLastError();
        if ( err != cudaSuccess)
        {
            printf("\n dot product kernel failed  %s \n",cudaGetErrorString(err) );
        }
}

#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_double_hh_transform_1_double( cuDoubleComplex*  ab_dev, cuDoubleComplex *hs_dev,cuDoubleComplex*  hv_dev, int nb , int ns)
#else
extern "C" void launch_double_hh_transform_1_single( cuFloatComplex*  ab_dev, cuFloatComplex *hs_dev,cuFloatComplex*  hv_dev, int nb , int ns)
#endif
{
        dim3  grid_size;
        grid_size = dim3(1,1, 1);
        cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to launch_double_hh_transform kernel: %s, %d\n",cudaGetErrorString(err), err);
#ifdef DOUBLE_PRECISION_COMPLEX
        double_hh_transform_kernel_double<<<grid_size, nb>>>( ab_dev, hs_dev, hv_dev, nb,  ns );
#else
        double_hh_transform_kernel_single<<<grid_size, nb>>>( ab_dev, hs_dev, hv_dev, nb,  ns );
#endif
        err = cudaGetLastError();
        if ( err != cudaSuccess)
        {
            printf("\n dot product kernel failed  %s \n",cudaGetErrorString(err) );
        }
}

#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_double_hh_transform_2_double( cuDoubleComplex*  ab_dev, cuDoubleComplex *hd_dev,cuDoubleComplex*  hv_dev, int nc , int ns , int nb )
#else
extern "C" void launch_double_hh_transform_2_single( cuFloatComplex*  ab_dev, cuFloatComplex *hd_dev,cuFloatComplex*  hv_dev, int nc , int ns , int nb )
#endif
{
        dim3  grid_size;
        grid_size = dim3(1,1, 1);
        cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to launch_double_hh_transform kernel: %s, %d\n",cudaGetErrorString(err), err);
#ifdef DOUBLE_PRECISION_COMPLEX
        double_hh_transform_kernel_2_double<<<grid_size, nc>>>( ab_dev, hd_dev, hv_dev, nc,  ns, nb );
#else
        double_hh_transform_kernel_2_single<<<grid_size, nc>>>( ab_dev, hd_dev, hv_dev, nc,  ns, nb );
#endif
        err = cudaGetLastError();
        if ( err != cudaSuccess)
        {
            printf("\n dot product kernel failed  %s \n",cudaGetErrorString(err) );
        }
}

#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_compute_kernel_reduce_double( cuDoubleComplex* a_dev, int lda, int n,int nbw, cuDoubleComplex* h_dev)
#else
extern "C" void launch_compute_kernel_reduce_single( cuFloatComplex* a_dev, int lda, int n,int nbw, cuFloatComplex* h_dev)
#endif
{
        dim3  grid_size;
        grid_size = dim3(1,1, 1);
        cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to launch_dot_product kernel: %s, %d\n",cudaGetErrorString(err), err);
#ifdef DOUBLE_PRECISION_COMPLEX
        compute_kernel_reduce_double<<<grid_size,n>>>(a_dev, lda, n, nbw,h_dev);
#else
        compute_kernel_reduce_single<<<grid_size,n>>>(a_dev, lda, n, nbw,h_dev);
#endif
	cudaDeviceSynchronize();
        err = cudaGetLastError();
        if ( err != cudaSuccess)
        {
            printf("\n dot product kernel failed  %s \n",cudaGetErrorString(err) );
        }
}

#ifdef DOUBLE_PRECISION_COMPLEX
extern "C" void launch_compute_kernel_reduce_1_double( cuDoubleComplex* a_dev, int lda, int n , cuDoubleComplex* h_dev)
#else
extern "C" void launch_compute_kernel_reduce_1_single( cuFloatComplex* a_dev, int lda, int n , cuFloatComplex* h_dev)
#endif
{
        dim3  grid_size;
        grid_size = dim3(1,1, 1);
        cudaDeviceSynchronize();
        cudaError_t err = cudaGetLastError();
        if(err != cudaSuccess) printf("error prior to launch_dot_product kernel: %s, %d\n",cudaGetErrorString(err), err);
#ifdef DOUBLE_PRECISION_COMPLEX
        compute_kernel_reduce_1_double<<<grid_size,n>>>(a_dev, lda, n, h_dev);
#else
        compute_kernel_reduce_1_single<<<grid_size,n>>>(a_dev, lda, n, h_dev);
#endif
	cudaDeviceSynchronize();
        err = cudaGetLastError();
        if ( err != cudaSuccess)
        {
            printf("\n dot product kernel failed  %s \n",cudaGetErrorString(err) );
        }

}
#endif /* COMPLEXCASE == 1 */

#ifndef MEMCPY_ALREADY_DEFINED
extern "C" int cuda_MemcpyDeviceToDevice(int val)
{
      val = cudaMemcpyDeviceToDevice;
      return val;
}
#define MEMCPY_ALREADY_DEFINED 1
#endif