mod_compute_hh_trafo_real.F90 38.1 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
!    This file is part of ELPA.
!
!    The ELPA library was originally created by the ELPA consortium,
!    consisting of the following organizations:
!
!    - 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 Naturwissenschaftrn,
!      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!      and
!    - IBM Deutschland GmbH
!
!
!    More information can be found here:
!    http://elpa.rzg.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 written by A. Marek, MPCDF



45
46
module compute_hh_trafo_real
#include "config-f90.h"
47
  use elpa_mpi
48
49
50
  implicit none

#ifdef WITH_OPENMP
51
52
53
54
55
56
57
58
  public compute_hh_trafo_real_cpu_openmp_double
#else
  public compute_hh_trafo_real_cpu_double
#endif

#ifdef WANT_SINGLE_PRECISION_REAL
#ifdef WITH_OPENMP
  public compute_hh_trafo_real_cpu_openmp_single
59
#else
60
61
62
63
  public compute_hh_trafo_real_cpu_single
#endif


64
65
66
67
68
#endif

  contains

#ifdef WITH_OPENMP
69
       subroutine compute_hh_trafo_real_cpu_openmp_double(a, a_dev, stripe_width, a_dim2, stripe_count, max_threads, l_nev, &
70
                                                   a_off, nbw, max_blk_size, bcast_buffer, bcast_buffer_dev, hh_dot_dev,    &
71
                                                   hh_tau_dev, kernel_flops, kernel_time, off, ncols, istripe,              &
72
                                                   my_thread, thread_width, THIS_REAL_ELPA_KERNEL)
73
#else
74
       subroutine compute_hh_trafo_real_cpu_double       (a, a_dev, stripe_width, a_dim2, stripe_count,                       &
75
76
77
                                                   a_off, nbw, max_blk_size, bcast_buffer,  bcast_buffer_dev, hh_dot_dev,     &
                                                   hh_tau_dev, kernel_flops, kernel_time, off, ncols, istripe,                &
                                                   last_stripe_width, THIS_REAL_ELPA_KERNEL)
78
79
80
81
#endif


         use precision
82
         use iso_c_binding
83
84
         use elpa2_utilities
         use single_hh_trafo_real
85
86
87
         use cuda_c_kernel
         use cuda_functions

88
#if defined(WITH_REAL_GENERIC_SIMPLE_KERNEL)
89
         use real_generic_simple_kernel !, only : double_hh_trafo_generic_simple
90
91
#endif

92
#if defined(WITH_REAL_GENERIC_KERNEL) && !(defined(DESPERATELY_WANT_ASSUMED_SIZE))
93
        use real_generic_kernel !, only : double_hh_trafo_generic
94
#endif
95
96

#if defined(WITH_REAL_BGP_KERNEL)
97
         use real_bgp_kernel !, only : double_hh_trafo_bgp
98
99
100
#endif

#if defined(WITH_REAL_BGQ_KERNEL)
101
         use real_bgq_kernel !, only : double_hh_trafo_bgq
102
103
104
105
106
#endif
#ifdef HAVE_DETAILED_TIMINGS
         use timings
#endif
         implicit none
107
         real(kind=c_double), intent(inout) :: kernel_time  ! MPI_WTIME always needs double
108
109
         integer(kind=lik)            :: kernel_flops
         integer(kind=ik), intent(in) :: nbw, max_blk_size
110
         real(kind=rk8)                :: bcast_buffer(nbw,max_blk_size)
111
112
113
114
115
116
         integer(kind=ik), intent(in) :: a_off

         integer(kind=ik), intent(in) :: stripe_width,a_dim2,stripe_count

#ifndef WITH_OPENMP
         integer(kind=ik), intent(in) :: last_stripe_width
117
118
!         real(kind=rk8)                :: a(stripe_width,a_dim2,stripe_count)
         real(kind=rk8), allocatable   :: a(:,:,:)
119
#else
120
         integer(kind=ik), intent(in) :: max_threads, l_nev, thread_width
121
122
!         real(kind=rk8)                :: a(stripe_width,a_dim2,stripe_count,max_threads)
         real(kind=rk8), allocatable   :: a(:,:,:,:)
123
124
125
#endif
         integer(kind=ik), intent(in) :: THIS_REAL_ELPA_KERNEL

126
127
128
129
130
         integer(kind=c_intptr_t)     :: a_dev
         integer(kind=c_intptr_t)     :: bcast_buffer_dev
         integer(kind=c_size_t)       :: dev_offset
         integer(kind=c_intptr_t)     :: hh_dot_dev
         integer(kind=c_intptr_t)     :: hh_tau_dev
131
132
133
134
135
136
         ! Private variables in OMP regions (my_thread) should better be in the argument list!
         integer(kind=ik)             :: off, ncols, istripe
#ifdef WITH_OPENMP
         integer(kind=ik)             :: my_thread, noff
#endif
         integer(kind=ik)             :: j, nl, jj, jjj
137
         real(kind=rk8)                :: w(nbw,6)
138
         real(kind=c_double)          :: ttt ! MPI_WTIME always needs double
139

140
141
142
143
144
         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GPU) then
           ! ncols - indicates the number of HH reflectors to apply; at least 1 must be available
           if (ncols < 1) return
         endif

145
146
#ifdef HAVE_DETAILED_TIMINGS
#ifdef WITH_OPENMP
147
         call timer%start("compute_hh_trafo_real_cpu_openmp_double")
148
#else
149
         call timer%start("compute_hh_trafo_real_cpu_double")
150
151
152
153
154
155
156
#endif
#endif
         ttt = mpi_wtime()

#ifndef WITH_OPENMP
         nl = merge(stripe_width, last_stripe_width, istripe<stripe_count)
#else
157
158
159
160
         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GPU) then
           print *,"compute_hh_trafo_real GPU OPENMP: not yet implemented"
           stop
         endif
161
162
163
164
165
166
167
168
169

         if (istripe<stripe_count) then
           nl = stripe_width
         else
           noff = (my_thread-1)*thread_width + (istripe-1)*stripe_width
           nl = min(my_thread*thread_width-noff, l_nev-noff)
           if (nl<=0) then
#ifdef HAVE_DETAILED_TIMINGS
#ifdef WITH_OPENMP
170
             call timer%stop("compute_hh_trafo_real_cpu_openmp_double")
171
#else
172
             call timer%stop("compute_hh_trafo_real_cpu_double")
173
174
175
176
177
#endif
#endif
             return
           endif
         endif
178
179
180
#endif /* not WITH_OPENMP */

         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GPU) then
181
182
           dev_offset = (0 + (a_off * stripe_width) + ( (istripe - 1) * stripe_width *a_dim2 )) *size_of_double_real_datatype
           call launch_compute_hh_trafo_c_kernel_real_double(a_dev + dev_offset, bcast_buffer_dev, hh_dot_dev, &
183
184
                                                      hh_tau_dev, nl, nbw, stripe_width, off, ncols)
         else ! not CUDA kernel
185
186

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
187
188
189
190
191
192
           if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK2 .or. &
               THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GENERIC    .or. &
               THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GENERIC_SIMPLE .or. &
               THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_SSE .or.        &
               THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_BGP .or.        &
               THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_BGQ) then
193
194
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */

195
             !FORTRAN CODE / X86 INRINISIC CODE / BG ASSEMBLER USING 2 HOUSEHOLDER VECTORS
196
197
#if defined(WITH_REAL_GENERIC_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
198
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GENERIC) then
199
200
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */

201
202
203
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
204
205

#ifdef WITH_OPENMP
206
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
207
                 call double_hh_trafo_generic_double(a(1,j+off+a_off-1,istripe,my_thread), w, &
208
209
                                            nbw, nl, stripe_width, nbw)

210
#else
211
                 call double_hh_trafo_generic_double(a(1:stripe_width,j+off+a_off-1:j+off+a_off+nbw-1, &
212
                                              istripe,my_thread), w(1:nbw,1:6), &
213
                                              nbw, nl, stripe_width, nbw)
214
215
216
217
218
#endif

#else /* WITH_OPENMP */

#ifdef DESPERATELY_WANT_ASSUMED_SIZE
219
                 call double_hh_trafo_generic_double(a(1,j+off+a_off-1,istripe),w, &
220
221
222
                                            nbw, nl, stripe_width, nbw)

#else
223
                 call double_hh_trafo_generic_double(a(1:stripe_width,j+off+a_off-1:j+off+a_off+nbw-1,istripe),w(1:nbw,1:6), &
224
                                            nbw, nl, stripe_width, nbw)
225
#endif
226
227
#endif /* WITH_OPENMP */

228
               enddo
229
230

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
231
             endif
232
233
234
235
236
237
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_GENERIC_KERNEL */


#if defined(WITH_REAL_GENERIC_SIMPLE_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
238
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GENERIC_SIMPLE) then
239
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
240
241
242
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
243
#ifdef WITH_OPENMP
244
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
245
                 call double_hh_trafo_generic_simple_double(a(1,j+off+a_off-1,istripe,my_thread), &
246
                                                       w, nbw, nl, stripe_width, nbw)
247
#else
248
                 call double_hh_trafo_generic_simple_double(a(1:stripe_width,j+off+a_off-1:j+off+a_off-1+nbw,istripe,my_thread), &
249
250
251
252
253
254
                                                     w, nbw, nl, stripe_width, nbw)

#endif

#else /* WITH_OPENMP */
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
255
                 call double_hh_trafo_generic_simple_double(a(1,j+off+a_off-1,istripe), &
256
                                                       w, nbw, nl, stripe_width, nbw)
257
#else
258
                 call double_hh_trafo_generic_simple_double(a(1:stripe_width,j+off+a_off-1:j+off+a_off-1+nbw,istripe), &
259
260
                                                     w, nbw, nl, stripe_width, nbw)

261
#endif
262
263
264

#endif /* WITH_OPENMP */

265
               enddo
266
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
267
             endif
268
269
270
271
272
273
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_GENERIC_SIMPLE_KERNEL */


#if defined(WITH_REAL_SSE_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
274
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_SSE) then
275
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
276
277
278
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
279
#ifdef WITH_OPENMP
280
                 call double_hh_trafo_double(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
281
282
                                      stripe_width, nbw)
#else
283
                 call double_hh_trafo_double(a(1,j+off+a_off-1,istripe), w, nbw, nl, &
284
285
                                      stripe_width, nbw)
#endif
286
               enddo
287
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
288
             endif
289
290
291
292
293
294
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_SSE_KERNEL */


#if defined(WITH_REAL_AVX_BLOCK2_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
295
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK2) then
296
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
297
298
299
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
300
#ifdef WITH_OPENMP
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
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
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
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
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
                 call double_hh_trafo_real_sse_avx_2hv_double(a(1,j+off+a_off-1,istripe,my_thread), &
                                                       w, nbw, nl, stripe_width, nbw)
#else
                 call double_hh_trafo_real_sse_avx_2hv_double(a(1,j+off+a_off-1,istripe), &
                                                       w, nbw, nl, stripe_width, nbw)
#endif
               enddo
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_AVX_BLOCK2_KERNEL */

#if defined(WITH_REAL_BGP_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_BGP) then
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
#ifdef WITH_OPENMP
                 call double_hh_trafo_bgp_double(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
                                          stripe_width, nbw)
#else
                 call double_hh_trafo_bgp_double(a(1,j+off+a_off-1,istripe), w, nbw, nl, &
                                          stripe_width, nbw)
#endif
               enddo
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_BGP_KERNEL */


#if defined(WITH_REAL_BGQ_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_BGQ) then
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
#ifdef WITH_OPENMP
                 call double_hh_trafo_bgq_double(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
                                          stripe_width, nbw)
#else
                 call double_hh_trafo_bgq_double(a(1,j+off+a_off-1,istripe), w, nbw, nl, &
                                          stripe_width, nbw)
#endif
               enddo
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_BGQ_KERNEL */


!#if defined(WITH_AVX_SANDYBRIDGE)
!              call double_hh_trafo_real_sse_avx_2hv(a(1,j+off+a_off-1,istripe), w, nbw, nl, stripe_width, nbw)
!#endif

#ifdef WITH_OPENMP
             if (j==1) call single_hh_trafo_real_cpu_openmp_double(a(1:stripe_width, &
                                      1+off+a_off:1+off+a_off+nbw-1,istripe,my_thread), &
                                      bcast_buffer(1:nbw,off+1), nbw, nl,     &
                                      stripe_width)
#else
             if (j==1) call single_hh_trafo_real_cpu_double(a(1:stripe_width,1+off+a_off:1+off+a_off+nbw-1,istripe),           &
                                      bcast_buffer(1:nbw,off+1), nbw, nl,     &
                                      stripe_width)
#endif


#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
           endif !
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */

#if defined(WITH_REAL_AVX_BLOCK4_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
           if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK4) then
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
             ! X86 INTRINSIC CODE, USING 4 HOUSEHOLDER VECTORS
             do j = ncols, 4, -4
               w(:,1) = bcast_buffer(1:nbw,j+off)
               w(:,2) = bcast_buffer(1:nbw,j+off-1)
               w(:,3) = bcast_buffer(1:nbw,j+off-2)
               w(:,4) = bcast_buffer(1:nbw,j+off-3)
#ifdef WITH_OPENMP
               call quad_hh_trafo_real_sse_avx_4hv_double(a(1,j+off+a_off-3,istripe,my_thread), w, &
                                                  nbw, nl, stripe_width, nbw)
#else
               call quad_hh_trafo_real_sse_avx_4hv_double(a(1,j+off+a_off-3,istripe), w, &
                                                  nbw, nl, stripe_width, nbw)
#endif
             enddo
             do jj = j, 2, -2
               w(:,1) = bcast_buffer(1:nbw,jj+off)
               w(:,2) = bcast_buffer(1:nbw,jj+off-1)
#ifdef WITH_OPENMP
               call double_hh_trafo_real_sse_avx_2hv_double(a(1,jj+off+a_off-1,istripe,my_thread), &
                                                    w, nbw, nl, stripe_width, nbw)
#else
               call double_hh_trafo_real_sse_avx_2hv_double(a(1,jj+off+a_off-1,istripe), &
                                                    w, nbw, nl, stripe_width, nbw)
#endif
             enddo
#ifdef WITH_OPENMP
             if (jj==1) call single_hh_trafo_real_cpu_openmp_double(a(1:stripe_width,1+off+a_off:1+off+a_off+nbw-1,istripe,my_thread), &
                                          bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#else
             if (jj==1) call single_hh_trafo_real_cpu_double(a(1:stripe_width,1+off+a_off:1+off+a_off+nbw-1,istripe), &
                                          bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#endif
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
           endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_AVX_BLOCK4_KERNEL */


#if defined(WITH_REAL_AVX_BLOCK6_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
           if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK6) then
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
             ! X86 INTRINSIC CODE, USING 6 HOUSEHOLDER VECTORS
             do j = ncols, 6, -6
               w(:,1) = bcast_buffer(1:nbw,j+off)
               w(:,2) = bcast_buffer(1:nbw,j+off-1)
               w(:,3) = bcast_buffer(1:nbw,j+off-2)
               w(:,4) = bcast_buffer(1:nbw,j+off-3)
               w(:,5) = bcast_buffer(1:nbw,j+off-4)
               w(:,6) = bcast_buffer(1:nbw,j+off-5)
#ifdef WITH_OPENMP
               call hexa_hh_trafo_real_sse_avx_6hv_double(a(1,j+off+a_off-5,istripe,my_thread), w, &
                                                  nbw, nl, stripe_width, nbw)
#else
               call hexa_hh_trafo_real_sse_avx_6hv_double(a(1,j+off+a_off-5,istripe), w, &
                                                  nbw, nl, stripe_width, nbw)
#endif
             enddo
             do jj = j, 4, -4
               w(:,1) = bcast_buffer(1:nbw,jj+off)
               w(:,2) = bcast_buffer(1:nbw,jj+off-1)
               w(:,3) = bcast_buffer(1:nbw,jj+off-2)
               w(:,4) = bcast_buffer(1:nbw,jj+off-3)
#ifdef WITH_OPENMP
               call quad_hh_trafo_real_sse_avx_4hv_double(a(1,jj+off+a_off-3,istripe,my_thread), w, &
                                                  nbw, nl, stripe_width, nbw)
#else
               call quad_hh_trafo_real_sse_avx_4hv_double(a(1,jj+off+a_off-3,istripe), w, &
                                                  nbw, nl, stripe_width, nbw)
#endif
             enddo
             do jjj = jj, 2, -2
               w(:,1) = bcast_buffer(1:nbw,jjj+off)
               w(:,2) = bcast_buffer(1:nbw,jjj+off-1)
#ifdef WITH_OPENMP
               call double_hh_trafo_real_sse_avx_2hv_double(a(1,jjj+off+a_off-1,istripe,my_thread), &
                                                    w, nbw, nl, stripe_width, nbw)
#else
               call double_hh_trafo_real_sse_avx_2hv_double(a(1,jjj+off+a_off-1,istripe), &
                                                    w, nbw, nl, stripe_width, nbw)
#endif
             enddo
#ifdef WITH_OPENMP
             if (jjj==1) call single_hh_trafo_real_cpu_openmp_double(a(1:stripe_width,1+off+a_off:1+off+a_off+nbw-1,istripe,my_thread), &
                                           bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#else
             if (jjj==1) call single_hh_trafo_real_cpu_double(a(1:stripe_width,1+off+a_off:1+off+a_off+nbw-1,istripe), &
                                           bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#endif
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
           endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_AVX_BLOCK4_KERNEL */

         endif ! GPU_KERNEL

#ifdef WITH_OPENMP
         if (my_thread==1) then
#endif
           kernel_flops = kernel_flops + 4*int(nl,8)*int(ncols,8)*int(nbw,8)
           kernel_time = kernel_time + mpi_wtime()-ttt
#ifdef WITH_OPENMP
         endif
#endif
#ifdef HAVE_DETAILED_TIMINGS
#ifdef WITH_OPENMP
         call timer%stop("compute_hh_trafo_real_cpu_openmp_double")
#else
         call timer%stop("compute_hh_trafo_real_cpu_double")
#endif
#endif

#ifdef WITH_OPENMP
       end subroutine compute_hh_trafo_real_cpu_openmp_double
#else
       end subroutine compute_hh_trafo_real_cpu_double
#endif

#ifdef WANT_SINGLE_PRECISION_REAL
! SINGLE PRECISION IMPLEMENTATION AT THE MOMENT DUPLICATED !!!

#ifdef WITH_OPENMP
       subroutine compute_hh_trafo_real_cpu_openmp_single(a, a_dev, stripe_width, a_dim2, stripe_count, max_threads, l_nev, &
                                                   a_off, nbw, max_blk_size, bcast_buffer, bcast_buffer_dev, hh_dot_dev,    &
                                                   hh_tau_dev, kernel_flops, kernel_time, off, ncols, istripe,              &
                                                   my_thread, thread_width, THIS_REAL_ELPA_KERNEL)
#else
       subroutine compute_hh_trafo_real_cpu_single       (a, a_dev, stripe_width, a_dim2, stripe_count,                      &
                                                   a_off, nbw, max_blk_size, bcast_buffer,  bcast_buffer_dev, hh_dot_dev,     &
                                                   hh_tau_dev, kernel_flops, kernel_time, off, ncols, istripe,                &
                                                   last_stripe_width, THIS_REAL_ELPA_KERNEL)
#endif


         use precision
         use iso_c_binding
         use elpa2_utilities
         use single_hh_trafo_real
         use cuda_c_kernel
         use cuda_functions

#if defined(WITH_REAL_GENERIC_SIMPLE_KERNEL)
         use real_generic_simple_kernel !, only : double_hh_trafo_generic_simple
#endif

#if defined(WITH_REAL_GENERIC_KERNEL) && !(defined(DESPERATELY_WANT_ASSUMED_SIZE))
        use real_generic_kernel !, only : double_hh_trafo_generic
#endif

#if defined(WITH_REAL_BGP_KERNEL)
         use real_bgp_kernel !, only : double_hh_trafo_bgp
#endif

#if defined(WITH_REAL_BGQ_KERNEL)
         use real_bgq_kernel !, only : double_hh_trafo_bgq
#endif
#ifdef HAVE_DETAILED_TIMINGS
         use timings
#endif
         implicit none
         real(kind=c_double), intent(inout) :: kernel_time  ! MPI_WTIME always needs double
         integer(kind=lik)            :: kernel_flops
         integer(kind=ik), intent(in) :: nbw, max_blk_size
         real(kind=rk4)                :: bcast_buffer(nbw,max_blk_size)
         integer(kind=ik), intent(in) :: a_off

         integer(kind=ik), intent(in) :: stripe_width,a_dim2,stripe_count

#ifndef WITH_OPENMP
         integer(kind=ik), intent(in) :: last_stripe_width
!         real(kind=rk4)                :: a(stripe_width,a_dim2,stripe_count)
         real(kind=rk4), allocatable   :: a(:,:,:)
#else
         integer(kind=ik), intent(in) :: max_threads, l_nev, thread_width
!         real(kind=rk8)                :: a(stripe_width,a_dim2,stripe_count,max_threads)
         real(kind=rk4), allocatable   :: a(:,:,:,:)
#endif
         integer(kind=ik), intent(in) :: THIS_REAL_ELPA_KERNEL

         integer(kind=c_intptr_t)     :: a_dev
         integer(kind=c_intptr_t)     :: bcast_buffer_dev
         integer(kind=c_size_t)       :: dev_offset
         integer(kind=c_intptr_t)     :: hh_dot_dev
         integer(kind=c_intptr_t)     :: hh_tau_dev
         ! Private variables in OMP regions (my_thread) should better be in the argument list!
         integer(kind=ik)             :: off, ncols, istripe
#ifdef WITH_OPENMP
         integer(kind=ik)             :: my_thread, noff
#endif
         integer(kind=ik)             :: j, nl, jj, jjj
         real(kind=rk4)                :: w(nbw,6)
         real(kind=c_double)          :: ttt ! MPI_WTIME always needs double

         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GPU) then
           ! ncols - indicates the number of HH reflectors to apply; at least 1 must be available
           if (ncols < 1) return
         endif

#ifdef HAVE_DETAILED_TIMINGS
#ifdef WITH_OPENMP
         call timer%start("compute_hh_trafo_real_cpu_openmp_single")
#else
         call timer%start("compute_hh_trafo_real_cpu_single")
#endif
#endif
         ttt = mpi_wtime()

#ifndef WITH_OPENMP
         nl = merge(stripe_width, last_stripe_width, istripe<stripe_count)
#else
         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GPU) then
           print *,"compute_hh_trafo_real GPU OPENMP: not yet implemented"
           stop
         endif

         if (istripe<stripe_count) then
           nl = stripe_width
         else
           noff = (my_thread-1)*thread_width + (istripe-1)*stripe_width
           nl = min(my_thread*thread_width-noff, l_nev-noff)
           if (nl<=0) then
#ifdef HAVE_DETAILED_TIMINGS
#ifdef WITH_OPENMP
             call timer%stop("compute_hh_trafo_real_cpu_openmp_single")
#else
             call timer%stop("compute_hh_trafo_real_cpu_single")
#endif
#endif
             return
           endif
         endif
#endif /* not WITH_OPENMP */

         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GPU) then
           dev_offset = (0 + (a_off * stripe_width) + ( (istripe - 1) * stripe_width *a_dim2 )) *size_of_single_real_datatype
           call launch_compute_hh_trafo_c_kernel_real_single(a_dev + dev_offset, bcast_buffer_dev, hh_dot_dev, &
                                                      hh_tau_dev, nl, nbw, stripe_width, off, ncols)
         else ! not CUDA kernel

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
           if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK2 .or. &
               THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GENERIC    .or. &
               THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GENERIC_SIMPLE .or. &
               THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_SSE .or.        &
               THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_BGP .or.        &
               THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_BGQ) then
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */

             !FORTRAN CODE / X86 INRINISIC CODE / BG ASSEMBLER USING 2 HOUSEHOLDER VECTORS
#if defined(WITH_REAL_GENERIC_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GENERIC) then
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */

               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)

#ifdef WITH_OPENMP
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
                 call double_hh_trafo_generic_single(a(1,j+off+a_off-1,istripe,my_thread), w, &
                                            nbw, nl, stripe_width, nbw)

#else
                 call double_hh_trafo_generic_single(a(1:stripe_width,j+off+a_off-1:j+off+a_off+nbw-1, &
                                              istripe,my_thread), w(1:nbw,1:6), &
                                              nbw, nl, stripe_width, nbw)
#endif

#else /* WITH_OPENMP */

#ifdef DESPERATELY_WANT_ASSUMED_SIZE
                 call double_hh_trafo_generic_single(a(1,j+off+a_off-1,istripe),w, &
                                            nbw, nl, stripe_width, nbw)

#else
                 call double_hh_trafo_generic_single(a(1:stripe_width,j+off+a_off-1:j+off+a_off+nbw-1,istripe),w(1:nbw,1:6), &
                                            nbw, nl, stripe_width, nbw)
#endif
#endif /* WITH_OPENMP */

               enddo

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_GENERIC_KERNEL */


#if defined(WITH_REAL_GENERIC_SIMPLE_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GENERIC_SIMPLE) then
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
#ifdef WITH_OPENMP
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
                 call double_hh_trafo_generic_simple_single(a(1,j+off+a_off-1,istripe,my_thread), &
                                                       w, nbw, nl, stripe_width, nbw)
#else
                 call double_hh_trafo_generic_simple_single(a(1:stripe_width,j+off+a_off-1:j+off+a_off-1+nbw,istripe,my_thread), &
                                                     w, nbw, nl, stripe_width, nbw)

#endif

#else /* WITH_OPENMP */
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
                 call double_hh_trafo_generic_simple_single(a(1,j+off+a_off-1,istripe), &
                                                       w, nbw, nl, stripe_width, nbw)
#else
                 call double_hh_trafo_generic_simple_single(a(1:stripe_width,j+off+a_off-1:j+off+a_off-1+nbw,istripe), &
                                                     w, nbw, nl, stripe_width, nbw)

#endif

#endif /* WITH_OPENMP */

               enddo
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_GENERIC_SIMPLE_KERNEL */


#if defined(WITH_REAL_SSE_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_SSE) then
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
#ifdef WITH_OPENMP
712
713
                 call double_hh_trafo_single(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
                                      stripe_width, nbw)
714
#else
715
716
                 call double_hh_trafo_single(a(1,j+off+a_off-1,istripe), w, nbw, nl, &
                                      stripe_width, nbw)
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
#endif
               enddo
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_SSE_KERNEL */


#if defined(WITH_REAL_AVX_BLOCK2_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK2) then
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
#ifdef WITH_OPENMP
                 call double_hh_trafo_real_sse_avx_2hv_single(a(1,j+off+a_off-1,istripe,my_thread), &
734
735
                                                       w, nbw, nl, stripe_width, nbw)
#else
736
                 call double_hh_trafo_real_sse_avx_2hv_single(a(1,j+off+a_off-1,istripe), &
737
738
                                                       w, nbw, nl, stripe_width, nbw)
#endif
739
               enddo
740
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
741
             endif
742
743
744
745
746
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_AVX_BLOCK2_KERNEL */

#if defined(WITH_REAL_BGP_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
747
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_BGP) then
748
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
749
750
751
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
752
#ifdef WITH_OPENMP
753
                 call double_hh_trafo_bgp_single(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
754
755
                                          stripe_width, nbw)
#else
756
                 call double_hh_trafo_bgp_single(a(1,j+off+a_off-1,istripe), w, nbw, nl, &
757
758
                                          stripe_width, nbw)
#endif
759
               enddo
760
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
761
             endif
762
763
764
765
766
767
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_BGP_KERNEL */


#if defined(WITH_REAL_BGQ_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
768
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_BGQ) then
769
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
770
771
772
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
773
#ifdef WITH_OPENMP
774
                 call double_hh_trafo_bgq_single(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
775
776
                                          stripe_width, nbw)
#else
777
                 call double_hh_trafo_bgq_single(a(1,j+off+a_off-1,istripe), w, nbw, nl, &
778
779
                                          stripe_width, nbw)
#endif
780
               enddo
781
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
782
             endif
783
784
785
786
787
788
789
790
791
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_BGQ_KERNEL */


!#if defined(WITH_AVX_SANDYBRIDGE)
!              call double_hh_trafo_real_sse_avx_2hv(a(1,j+off+a_off-1,istripe), w, nbw, nl, stripe_width, nbw)
!#endif

#ifdef WITH_OPENMP
792
793
             if (j==1) call single_hh_trafo_real_cpu_openmp_single(a(1:stripe_width, &
                                   1+off+a_off:1+off+a_off+nbw-1,istripe,my_thread), &
794
795
796
                                      bcast_buffer(1:nbw,off+1), nbw, nl,     &
                                      stripe_width)
#else
797
798
             if (j==1) call single_hh_trafo_real_cpu_single(a(1:stripe_width, &
                                    1+off+a_off:1+off+a_off+nbw-1,istripe),           &
799
800
801
802
803
804
                                      bcast_buffer(1:nbw,off+1), nbw, nl,     &
                                      stripe_width)
#endif


#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
805
           endif !
806
807
808
809
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */

#if defined(WITH_REAL_AVX_BLOCK4_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
810
           if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK4) then
811
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
812
813
814
815
816
817
             ! X86 INTRINSIC CODE, USING 4 HOUSEHOLDER VECTORS
             do j = ncols, 4, -4
               w(:,1) = bcast_buffer(1:nbw,j+off)
               w(:,2) = bcast_buffer(1:nbw,j+off-1)
               w(:,3) = bcast_buffer(1:nbw,j+off-2)
               w(:,4) = bcast_buffer(1:nbw,j+off-3)
818
#ifdef WITH_OPENMP
819
               call quad_hh_trafo_real_sse_avx_4hv_single(a(1,j+off+a_off-3,istripe,my_thread), w, &
820
821
                                                  nbw, nl, stripe_width, nbw)
#else
822
               call quad_hh_trafo_real_sse_avx_4hv_single(a(1,j+off+a_off-3,istripe), w, &
823
824
                                                  nbw, nl, stripe_width, nbw)
#endif
825
826
827
828
             enddo
             do jj = j, 2, -2
               w(:,1) = bcast_buffer(1:nbw,jj+off)
               w(:,2) = bcast_buffer(1:nbw,jj+off-1)
829
#ifdef WITH_OPENMP
830
               call double_hh_trafo_real_sse_avx_2hv_single(a(1,jj+off+a_off-1,istripe,my_thread), &
831
832
                                                    w, nbw, nl, stripe_width, nbw)
#else
833
               call double_hh_trafo_real_sse_avx_2hv_single(a(1,jj+off+a_off-1,istripe), &
834
835
                                                    w, nbw, nl, stripe_width, nbw)
#endif
836
             enddo
837
#ifdef WITH_OPENMP
838
             if (jj==1) call single_hh_trafo_real_cpu_openmp_single(a(1:stripe_width,1+off+a_off:1+off+a_off+nbw-1,istripe,my_thread), &
839
840
                                          bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#else
841
             if (jj==1) call single_hh_trafo_real_cpu_single(a(1:stripe_width,1+off+a_off:1+off+a_off+nbw-1,istripe), &
842
843
844
                                          bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#endif
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
845
           endif
846
847
848
849
850
851
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_AVX_BLOCK4_KERNEL */


#if defined(WITH_REAL_AVX_BLOCK6_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
852
           if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK6) then
853
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
854
855
856
857
858
859
860
861
             ! X86 INTRINSIC CODE, USING 6 HOUSEHOLDER VECTORS
             do j = ncols, 6, -6
               w(:,1) = bcast_buffer(1:nbw,j+off)
               w(:,2) = bcast_buffer(1:nbw,j+off-1)
               w(:,3) = bcast_buffer(1:nbw,j+off-2)
               w(:,4) = bcast_buffer(1:nbw,j+off-3)
               w(:,5) = bcast_buffer(1:nbw,j+off-4)
               w(:,6) = bcast_buffer(1:nbw,j+off-5)
862
#ifdef WITH_OPENMP
863
               call hexa_hh_trafo_real_sse_avx_6hv_single(a(1,j+off+a_off-5,istripe,my_thread), w, &
864
865
                                                  nbw, nl, stripe_width, nbw)
#else
866
               call hexa_hh_trafo_real_sse_avx_6hv_single(a(1,j+off+a_off-5,istripe), w, &
867
868
                                                  nbw, nl, stripe_width, nbw)
#endif
869
870
871
872
873
874
             enddo
             do jj = j, 4, -4
               w(:,1) = bcast_buffer(1:nbw,jj+off)
               w(:,2) = bcast_buffer(1:nbw,jj+off-1)
               w(:,3) = bcast_buffer(1:nbw,jj+off-2)
               w(:,4) = bcast_buffer(1:nbw,jj+off-3)
875
#ifdef WITH_OPENMP
876
               call quad_hh_trafo_real_sse_avx_4hv_single(a(1,jj+off+a_off-3,istripe,my_thread), w, &
877
878
                                                  nbw, nl, stripe_width, nbw)
#else
879
               call quad_hh_trafo_real_sse_avx_4hv_single(a(1,jj+off+a_off-3,istripe), w, &
880
881
                                                  nbw, nl, stripe_width, nbw)
#endif
882
883
884
885
             enddo
             do jjj = jj, 2, -2
               w(:,1) = bcast_buffer(1:nbw,jjj+off)
               w(:,2) = bcast_buffer(1:nbw,jjj+off-1)
886
#ifdef WITH_OPENMP
887
               call double_hh_trafo_real_sse_avx_2hv_single(a(1,jjj+off+a_off-1,istripe,my_thread), &
888
889
                                                    w, nbw, nl, stripe_width, nbw)
#else
890
               call double_hh_trafo_real_sse_avx_2hv_single(a(1,jjj+off+a_off-1,istripe), &
891
892
                                                    w, nbw, nl, stripe_width, nbw)
#endif
893
             enddo
894
#ifdef WITH_OPENMP
895
             if (jjj==1) call single_hh_trafo_real_cpu_openmp_single(a(1:stripe_width,1+off+a_off:1+off+a_off+nbw-1,istripe,my_thread), &
896
897
                                           bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#else
898
             if (jjj==1) call single_hh_trafo_real_cpu_single(a(1:stripe_width,1+off+a_off:1+off+a_off+nbw-1,istripe), &
899
900
901
                                           bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#endif
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
902
           endif
903
904
905
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_AVX_BLOCK4_KERNEL */

906
907
         endif ! GPU_KERNEL

908
909
910
911
912
913
914
915
916
917
#ifdef WITH_OPENMP
         if (my_thread==1) then
#endif
           kernel_flops = kernel_flops + 4*int(nl,8)*int(ncols,8)*int(nbw,8)
           kernel_time = kernel_time + mpi_wtime()-ttt
#ifdef WITH_OPENMP
         endif
#endif
#ifdef HAVE_DETAILED_TIMINGS
#ifdef WITH_OPENMP
918
         call timer%stop("compute_hh_trafo_real_cpu_openmp_single")
919
#else
920
         call timer%stop("compute_hh_trafo_real_cpu_single")
921
922
923
924
#endif
#endif

#ifdef WITH_OPENMP
925
       end subroutine compute_hh_trafo_real_cpu_openmp_single
926
#else
927
       end subroutine compute_hh_trafo_real_cpu_single
928
#endif
929
930
#endif /* WANT_SINGLE_PRECISION_REAL */

931
932

end module