mod_compute_hh_trafo_real.F90 44.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
!         real(kind=rk8)                :: a(stripe_width,a_dim2,stripe_count)
118
         real(kind=rk8), pointer      :: a(:,:,:)
119
#else
120
         integer(kind=ik), intent(in) :: max_threads, l_nev, thread_width
121
!         real(kind=rk8)                :: a(stripe_width,a_dim2,stripe_count,max_threads)
122
         real(kind=rk8), pointer      :: 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
187

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK2 .or. &
188
             THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX2_BLOCK2 .or. &
189
             THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_SSE_BLOCK2 .or. &
190
191
192
193
194
195
196
             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 */

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

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

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

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

#else /* WITH_OPENMP */

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

#else
225
                 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), &
226
                                            nbw, nl, stripe_width, nbw)
227
#endif
228
229
#endif /* WITH_OPENMP */

230
               enddo
231
232

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
233
             endif
234
235
236
237
238
239
#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)
240
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GENERIC_SIMPLE) then
241
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
242
243
244
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
245
#ifdef WITH_OPENMP
246
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
247
                 call double_hh_trafo_generic_simple_double(a(1,j+off+a_off-1,istripe,my_thread), &
248
                                                       w, nbw, nl, stripe_width, nbw)
249
#else
250
                 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), &
251
252
253
254
255
256
                                                     w, nbw, nl, stripe_width, nbw)

#endif

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

263
#endif
264
265
266

#endif /* WITH_OPENMP */

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

294
295
296
297
298
299
300
301
#if defined(WITH_REAL_SSE_BLOCK2_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
           if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_SSE_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
302
               call double_hh_trafo_real_sse_2hv_double(a(1,j+off+a_off-1,istripe,my_thread), &
303
304
                                                       w, nbw, nl, stripe_width, nbw)
#else
305
               call double_hh_trafo_real_sse_2hv_double(a(1,j+off+a_off-1,istripe), &
306
307
308
309
310
311
312
                                                       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_SSE_BLOCK2_KERNEL */
313

314
#if defined(WITH_REAL_AVX_BLOCK2_KERNEL) || defined(WITH_REAL_AVX2_BLOCK2_KERNEL)
315
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
316

317
318
           if ((THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK2) .or. &
               (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX2_BLOCK2))  then
319
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
320
321
322
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
323
#ifdef WITH_OPENMP
324

325
               call double_hh_trafo_real_avx_avx2_2hv_double(a(1,j+off+a_off-1,istripe,my_thread), &
326
327
                                                       w, nbw, nl, stripe_width, nbw)
#else
328
               call double_hh_trafo_real_avx_avx2_2hv_double(a(1,j+off+a_off-1,istripe), &
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
                                                       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)
380
!              call double_hh_trafo_real_avx_avx2_2hv(a(1,j+off+a_off-1,istripe), w, nbw, nl, stripe_width, nbw)
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
!#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
410
               call quad_hh_trafo_real_avx_avx2_4hv_double(a(1,j+off+a_off-3,istripe,my_thread), w, &
411
412
                                                  nbw, nl, stripe_width, nbw)
#else
413
               call quad_hh_trafo_real_avx_avx2_4hv_double(a(1,j+off+a_off-3,istripe), w, &
414
415
416
417
418
419
420
                                                  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
421
               call double_hh_trafo_real_avx_avx2_2hv_double(a(1,jj+off+a_off-1,istripe,my_thread), &
422
423
                                                    w, nbw, nl, stripe_width, nbw)
#else
424
               call double_hh_trafo_real_avx_avx2_2hv_double(a(1,jj+off+a_off-1,istripe), &
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
                                                    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
454
               call hexa_hh_trafo_real_avx_avx2_6hv_double(a(1,j+off+a_off-5,istripe,my_thread), w, &
455
456
                                                  nbw, nl, stripe_width, nbw)
#else
457
               call hexa_hh_trafo_real_avx_avx2_6hv_double(a(1,j+off+a_off-5,istripe), w, &
458
459
460
461
462
463
464
465
466
                                                  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
467
               call quad_hh_trafo_real_avx_avx2_4hv_double(a(1,jj+off+a_off-3,istripe,my_thread), w, &
468
469
                                                  nbw, nl, stripe_width, nbw)
#else
470
               call quad_hh_trafo_real_avx_avx2_4hv_double(a(1,jj+off+a_off-3,istripe), w, &
471
472
473
474
475
476
477
                                                  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
478
               call double_hh_trafo_real_avx_avx2_2hv_double(a(1,jjj+off+a_off-1,istripe,my_thread), &
479
480
                                                    w, nbw, nl, stripe_width, nbw)
#else
481
               call double_hh_trafo_real_avx_avx2_2hv_double(a(1,jjj+off+a_off-1,istripe), &
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
                                                    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)
574
         real(kind=rk4), pointer      :: a(:,:,:)
575
576
577
#else
         integer(kind=ik), intent(in) :: max_threads, l_nev, thread_width
!         real(kind=rk8)                :: a(stripe_width,a_dim2,stripe_count,max_threads)
578
         real(kind=rk4), pointer      :: a(:,:,:,:)
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
#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
736
737
                 call double_hh_trafo_single(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
                                      stripe_width, nbw)
738
#else
739
740
                 call double_hh_trafo_single(a(1,j+off+a_off-1,istripe), w, nbw, nl, &
                                      stripe_width, nbw)
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
#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
757
                 call double_hh_trafo_real_avx_avx2_2hv_single(a(1,j+off+a_off-1,istripe,my_thread), &
758
759
                                                       w, nbw, nl, stripe_width, nbw)
#else
760
                 call double_hh_trafo_real_avx_avx2_2hv_single(a(1,j+off+a_off-1,istripe), &
761
762
                                                       w, nbw, nl, stripe_width, nbw)
#endif
763
               enddo
764
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
765
             endif
766
767
768
769
770
#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)
771
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_BGP) then
772
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
773
774
775
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
776
#ifdef WITH_OPENMP
777
                 call double_hh_trafo_bgp_single(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
778
779
                                          stripe_width, nbw)
#else
780
                 call double_hh_trafo_bgp_single(a(1,j+off+a_off-1,istripe), w, nbw, nl, &
781
782
                                          stripe_width, nbw)
#endif
783
               enddo
784
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
785
             endif
786
787
788
789
790
791
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_BGP_KERNEL */


#if defined(WITH_REAL_BGQ_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
792
             if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_BGQ) then
793
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
794
795
796
               do j = ncols, 2, -2
                 w(:,1) = bcast_buffer(1:nbw,j+off)
                 w(:,2) = bcast_buffer(1:nbw,j+off-1)
797
#ifdef WITH_OPENMP
798
                 call double_hh_trafo_bgq_single(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
799
800
                                          stripe_width, nbw)
#else
801
                 call double_hh_trafo_bgq_single(a(1,j+off+a_off-1,istripe), w, nbw, nl, &
802
803
                                          stripe_width, nbw)
#endif
804
               enddo
805
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
806
             endif
807
808
809
810
811
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_BGQ_KERNEL */


!#if defined(WITH_AVX_SANDYBRIDGE)
812
!              call double_hh_trafo_real_avx_avx2_2hv(a(1,j+off+a_off-1,istripe), w, nbw, nl, stripe_width, nbw)
813
814
815
!#endif

#ifdef WITH_OPENMP
816
817
             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), &
818
819
820
                                      bcast_buffer(1:nbw,off+1), nbw, nl,     &
                                      stripe_width)
#else
821
822
             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),           &
823
824
825
826
827
828
                                      bcast_buffer(1:nbw,off+1), nbw, nl,     &
                                      stripe_width)
#endif


#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
829
           endif !
830
831
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */

832
#if defined(WITH_REAL_SSE_BLOCK4_KERNEL)
833
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
834
835
836
837
838
839
840
841
842
         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_SSE_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
843
             call quad_hh_trafo_real_sse_4hv_single(a(1,j+off+a_off-3,istripe,my_thread), w, &
844
845
                                                  nbw, nl, stripe_width, nbw)
#else
846
             call quad_hh_trafo_real_sse_4hv_single(a(1,j+off+a_off-3,istripe), w, &
847
848
849
850
851
852
853
                                                  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
854
             call double_hh_trafo_real_sse_2hv_single(a(1,jj+off+a_off-1,istripe,my_thread), &
855
856
                                                    w, nbw, nl, stripe_width, nbw)
#else
857
             call double_hh_trafo_real_sse_2hv_single(a(1,jj+off+a_off-1,istripe), &
858
859
860
861
                                                    w, nbw, nl, stripe_width, nbw)
#endif
           enddo
#ifdef WITH_OPENMP
862
           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), &
863
864
                                          bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#else
865
           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), &
866
867
868
869
870
871
872
                                          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_SSE_BLOCK4_KERNEL */

873
#if defined(WITH_REAL_AVX_BLOCK4_KERNEL) || defined(WITH_REAL_AVX2_BLOCK4_KERNEL)
874
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
875
876
         if ((THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK4) .or. &
             (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX2_BLOCK4)) then
877
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
878
879
880
881
882
883
             ! 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)
884
#ifdef WITH_OPENMP
885

886
             call quad_hh_trafo_real_avx_avx2_4hv_single(a(1,j+off+a_off-3,istripe,my_thread), w, &
887
888
                                                  nbw, nl, stripe_width, nbw)
#else
889
             call quad_hh_trafo_real_avx_avx2_4hv_single(a(1,j+off+a_off-3,istripe), w, &
890
891
                                                  nbw, nl, stripe_width, nbw)
#endif
892
893
894
895
             enddo
             do jj = j, 2, -2
               w(:,1) = bcast_buffer(1:nbw,jj+off)
               w(:,2) = bcast_buffer(1:nbw,jj+off-1)
896
#ifdef WITH_OPENMP
897

898
             call double_hh_trafo_real_avx_avx2_2hv_single(a(1,jj+off+a_off-1,istripe,my_thread), &
899
900
                                                    w, nbw, nl, stripe_width, nbw)
#else
901
             call double_hh_trafo_real_avx_avx2_2hv_single(a(1,jj+off+a_off-1,istripe), &
902
903
                                                    w, nbw, nl, stripe_width, nbw)
#endif
904
             enddo
905
#ifdef WITH_OPENMP
906
             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), &
907
908
                                          bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#else
909
             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), &
910
911
912
                                          bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#endif
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
913
           endif
914
915
916
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_AVX_BLOCK4_KERNEL */

917
918
919
920
921
922
923
924
925
926
927
928
929
#if defined(WITH_REAL_SSE_BLOCK6_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_SSE_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
930
             call hexa_hh_trafo_real_sse_6hv_single(a(1,j+off+a_off-5,istripe,my_thread), w, &
931
932
                                                  nbw, nl, stripe_width, nbw)
#else
933
             call hexa_hh_trafo_real_sse_6hv_single(a(1,j+off+a_off-5,istripe), w, &
934
935
936
937
938
939
940
941
942
                                                  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
943
             call quad_hh_trafo_real_sse_4hv_single(a(1,jj+off+a_off-3,istripe,my_thread), w, &
944
945
                                                  nbw, nl, stripe_width, nbw)
#else
946
             call quad_hh_trafo_real_sse_4hv_single(a(1,jj+off+a_off-3,istripe), w, &
947
948
949
950
951
952
953
                                                  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
954
             call double_hh_trafo_real_sse_2hv_single(a(1,jjj+off+a_off-1,istripe,my_thread), &
955
956
                                                    w, nbw, nl, stripe_width, nbw)
#else
957
             call double_hh_trafo_real_sse_2hv_single(a(1,jjj+off+a_off-1,istripe), &
958
959
960
961
                                                    w, nbw, nl, stripe_width, nbw)
#endif
           enddo
#ifdef WITH_OPENMP
962
           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), &
963
964
                                           bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#else
965
           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), &
966
967
968
969
970
971
                                           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_SSE_BLOCK4_KERNEL */
972

973
#if defined(WITH_REAL_AVX_BLOCK6_KERNEL) || defined(WITH_REAL_AVX2_BLOCK6_KERNEL)
974
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
975

976
977
         if ((THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK6) .or. &
             (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX2_BLOCK6)) then
978
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
979
980
981
982
983
984
985
986
             ! 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)
987
#ifdef WITH_OPENMP
988

989
             call hexa_hh_trafo_real_avx_avx2_6hv_single(a(1,j+off+a_off-5,istripe,my_thread), w, &
990
991
                                                  nbw, nl, stripe_width, nbw)
#else
992
             call hexa_hh_trafo_real_avx_avx2_6hv_single(a(1,j+off+a_off-5,istripe), w, &
993
994
                                                  nbw, nl, stripe_width, nbw)
#endif
995
996
997
998
999
1000
             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)
1001
#ifdef WITH_OPENMP
1002

1003
             call quad_hh_trafo_real_avx_avx2_4hv_single(a(1,jj+off+a_off-3,istripe,my_thread), w, &
1004
1005
                                                  nbw, nl, stripe_width, nbw)
#else
1006
             call quad_hh_trafo_real_avx_avx2_4hv_single(a(1,jj+off+a_off-3,istripe), w, &
1007
1008
                                                  nbw, nl, stripe_width, nbw)
#endif
1009
1010
1011
1012
             enddo
             do jjj = jj, 2, -2
               w(:,1) = bcast_buffer(1:nbw,jjj+off)
               w(:,2) = bcast_buffer(1:nbw,jjj+off-1)
1013
#ifdef WITH_OPENMP
1014

1015
             call double_hh_trafo_real_avx_avx2_2hv_single(a(1,jjj+off+a_off-1,istripe,my_thread), &
1016
1017
                                                    w, nbw, nl, stripe_width, nbw)
#else
1018
             call double_hh_trafo_real_avx_avx2_2hv_single(a(1,jjj+off+a_off-1,istripe), &
1019
1020
                                                    w, nbw, nl, stripe_width, nbw)
#endif
1021
             enddo
1022
#ifdef WITH_OPENMP
1023
             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), &
1024
1025
                                           bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#else
1026
             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), &
1027
1028
1029
                                           bcast_buffer(1:nbw,off+1), nbw, nl, stripe_width)
#endif
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
1030
           endif
1031
1032
1033
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_AVX_BLOCK4_KERNEL */

1034
1035
         endif ! GPU_KERNEL

1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
#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
1046
         call timer%stop("compute_hh_trafo_real_cpu_openmp_single")
1047
#else
1048
         call timer%stop("compute_hh_trafo_real_cpu_single")
1049
1050
1051
1052
#endif
#endif

#ifdef WITH_OPENMP
1053
       end subroutine compute_hh_trafo_real_cpu_openmp_single
1054
#else
1055
       end subroutine compute_hh_trafo_real_cpu_single
1056
#endif
1057
1058
#endif /* WANT_SINGLE_PRECISION_REAL */

1059
1060

end module