mod_compute_hh_trafo_real.F90 25.3 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
!    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 Naturwissenschaftrn,
!      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!      and
!    - IBM Deutschland GmbH
!
!
!    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.
!
! Author: Andreas Marek, MPCDF

44
45
module compute_hh_trafo_real
#include "config-f90.h"
46
  use elpa_mpi
47
48
49
50
51
52
53
54
55
56
57
  implicit none

#ifdef WITH_OPENMP
  public compute_hh_trafo_real_cpu_openmp
#else
  public compute_hh_trafo_real_cpu
#endif

  contains

#ifdef WITH_OPENMP
58
       subroutine compute_hh_trafo_real_cpu_openmp(a, stripe_width, a_dim2, stripe_count, max_threads, l_nev,         &
59
60
                                                   a_off, nbw, max_blk_size, bcast_buffer, kernel_flops, kernel_time, &
                                                   off, ncols, istripe,                                               &
61
                                                   my_thread, thread_width,  THIS_REAL_ELPA_KERNEL)
62
#else
63
       subroutine compute_hh_trafo_real_cpu       (a, stripe_width, a_dim2, stripe_count,                              &
64
65
66
67
68
69
70
71
72
73
74
75
76
                                                   a_off, nbw, max_blk_size, bcast_buffer,  kernel_flops, kernel_time, &
                                                   off, ncols, istripe, last_stripe_width,                             &
                                                   THIS_REAL_ELPA_KERNEL)
#endif


         use precision
         use elpa2_utilities
         use single_hh_trafo_real
#if defined(WITH_REAL_GENERIC_SIMPLE_KERNEL)
         use real_generic_simple_kernel, only : double_hh_trafo_generic_simple
#endif

77
78
79
#if defined(WITH_REAL_GENERIC_KERNEL) && !(defined(DESPERATELY_WANT_ASSUMED_SIZE))
        use real_generic_kernel, only : double_hh_trafo_generic
#endif
80
81
82
83
84
85
86
87
88
89
90

#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
91

92
#if defined(HAVE_AVX) || defined(HAVE_SSE_INTRINSICS) || defined(HAVE_SSE_ASSEMBLY)
93
94
         use kernel_interfaces
#endif
95
96
97
98
99
100
101
102
103
104
105
106
107
         implicit none
         real(kind=rk), intent(inout) :: kernel_time
         integer(kind=lik)            :: kernel_flops
         integer(kind=ik), intent(in) :: nbw, max_blk_size
         real(kind=rk)                :: 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=rk)                :: a(stripe_width,a_dim2,stripe_count)
#else
108
         integer(kind=ik), intent(in) :: max_threads, l_nev, thread_width
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
         real(kind=rk)                :: a(stripe_width,a_dim2,stripe_count,max_threads)
#endif
         integer(kind=ik), intent(in) :: THIS_REAL_ELPA_KERNEL

         ! 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=rk)                :: w(nbw,6), ttt

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

#ifndef WITH_OPENMP
         nl = merge(stripe_width, last_stripe_width, istripe<stripe_count)
#else

         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")
#else
             call timer%stop("compute_hh_trafo_real_cpu")
#endif
#endif
             return
           endif
         endif
#endif

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK2 .or. &
154
             THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX2_BLOCK2 .or. &
155
             THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_SSE_BLOCK2 .or. &
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
             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
174
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
175
               call double_hh_trafo_generic(a(1,j+off+a_off-1,istripe,my_thread), w, &
176
177
                                            nbw, nl, stripe_width, nbw)

178
#else
179
180
               call double_hh_trafo_generic(a(1:stripe_width,j+off+a_off-1:j+off+a_off+nbw-1, &
                                              istripe,my_thread), w(1:nbw,1:6), &
181
                                              nbw, nl, stripe_width, nbw)
182
183
184
185
186
187
188
189
190
#endif

#else /* WITH_OPENMP */

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

#else
191
192
               call double_hh_trafo_generic(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)
193
#endif
194
195
#endif /* WITH_OPENMP */

196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
             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
212
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
213
214
215
               call double_hh_trafo_generic_simple(a(1,j+off+a_off-1,istripe,my_thread), &
                                                     w, nbw, nl, stripe_width, nbw)
#else
216
217
218
219
220
221
222
               call double_hh_trafo_generic_simple(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
223
224
               call double_hh_trafo_generic_simple(a(1,j+off+a_off-1,istripe), &
                                                     w, nbw, nl, stripe_width, nbw)
225
226
227
228
#else
               call double_hh_trafo_generic_simple(a(1:stripe_width,j+off+a_off-1:j+off+a_off-1+nbw,istripe), &
                                                     w, nbw, nl, stripe_width, nbw)

229
#endif
230
231
232

#endif /* WITH_OPENMP */

233
234
235
236
237
238
239
             enddo
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
           endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_GENERIC_SIMPLE_KERNEL */


240
#if defined(WITH_REAL_SSE_ASSEMBLY_KERNEL)
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
#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
               call double_hh_trafo(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
                                      stripe_width, nbw)
#else
               call double_hh_trafo(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 */
258
#endif /* WITH_REAL_SSE_ASSEMBLY_KERNEL */
259

260
261
262
263
#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 */
264
265

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL) || (defined(WITH_ONE_SPECIFIC_REAL_KERNEL) && !defined(WITH_REAL_SSE_BLOCK6_KERNEL) && !defined(WITH_REAL_SSE_BLOCK4_KERNEL))
266
267
268
269
270
271
272
273
274
275
276
             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_2hv(a(1,j+off+a_off-1,istripe,my_thread), &
                                                       w, nbw, nl, stripe_width, nbw)
#else
               call double_hh_trafo_real_sse_2hv(a(1,j+off+a_off-1,istripe), &
                                                       w, nbw, nl, stripe_width, nbw)
#endif
             enddo
277
278
#endif /* defined(WITH_NO_SPECIFIC_REAL_KERNEL) || (defined(WITH_ONE_SPECIFIC_REAL_KERNEL) && !defined(WITH_REAL_SSE_BLOCK6_KERNEL) && !defined(WITH_REAL_SSE_BLOCK4_KERNEL)) */

279
280
281
282
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
           endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_SSE_BLOCK2_KERNEL */
283

284
#if defined(WITH_REAL_AVX_BLOCK2_KERNEL) || defined(WITH_REAL_AVX2_BLOCK2_KERNEL)
285
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
286
287
           if ((THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK2) .or. &
               (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX2_BLOCK2))  then
288
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
289
290

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL) || (defined(WITH_ONE_SPECIFIC_REAL_KERNEL) && !defined(WITH_REAL_AVX_BLOCK6_KERNEL) && !defined(WITH_REAL_AVX_BLOCK4_KERNEL) && !defined(WITH_REAL_AVX2_BLOCK6_KERNEL) && !defined(WITH_REAL_AVX2_BLOCK4_KERNEL))
291
292
293
294
             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
295
               call double_hh_trafo_real_avx_avx2_2hv(a(1,j+off+a_off-1,istripe,my_thread), &
296
297
                                                       w, nbw, nl, stripe_width, nbw)
#else
298
               call double_hh_trafo_real_avx_avx2_2hv(a(1,j+off+a_off-1,istripe), &
299
300
301
                                                       w, nbw, nl, stripe_width, nbw)
#endif
             enddo
302
303
#endif /* defined(WITH_NO_SPECIFIC_REAL_KERNEL) || (defined(WITH_ONE_SPECIFIC_REAL_KERNEL) ... */

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
#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(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
                                          stripe_width, nbw)
#else
               call double_hh_trafo_bgp(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(a(1,j+off+a_off-1,istripe,my_thread), w, nbw, nl, &
                                          stripe_width, nbw)
#else
               call double_hh_trafo_bgq(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)
352
!              call double_hh_trafo_real_avx_avx2_2hv(a(1,j+off+a_off-1,istripe), w, nbw, nl, stripe_width, nbw)
353
354
355
!#endif

#ifdef WITH_OPENMP
356
           if (j==1) call single_hh_trafo_real_cpu_openmp(a(1:stripe_width,1+off+a_off:1+off+a_off+nbw-1,istripe,my_thread), &
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
                                      bcast_buffer(1:nbw,off+1), nbw, nl,     &
                                      stripe_width)
#else
           if (j==1) call single_hh_trafo_real_cpu(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 */



372
373
374
375
#if defined(WITH_REAL_SSE_BLOCK4_KERNEL)
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
         if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_SSE_BLOCK4) then
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
376
377

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL) || (defined(WITH_ONE_SPECIFIC_REAL_KERNEL) && !defined(WITH_REAL_SSE_BLOCK6_KERNEL))
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
           ! 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_4hv(a(1,j+off+a_off-3,istripe,my_thread), w, &
                                                  nbw, nl, stripe_width, nbw)
#else
             call quad_hh_trafo_real_sse_4hv(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_2hv(a(1,jj+off+a_off-1,istripe,my_thread), &
                                                    w, nbw, nl, stripe_width, nbw)
#else
             call double_hh_trafo_real_sse_2hv(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(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(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
410
411
412

#endif /* defined(WITH_NO_SPECIFIC_REAL_KERNEL) || (defined(WITH_ONE_SPECIFIC_REAL_KERNEL) && !defined(WITH_REAL_SSE_BLOCK6_KERNEL)) */

413
414
415
416
417
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
         endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_SSE_BLOCK4_KERNEL */

418
#if defined(WITH_REAL_AVX_BLOCK4_KERNEL) || defined(WITH_REAL_AVX2_BLOCK4_KERNEL)
419
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
420
421
         if ((THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK4) .or. &
             (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX2_BLOCK4)) then
422
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
423
424

#if defined(WITH_NO_SPECIFIC_REAL_KERNEL) || (defined(WITH_ONE_SPECIFIC_REAL_KERNEL) && !defined(WITH_REAL_AVX_BLOCK6_KERNEL) && !defined(WITH_REAL_AVX2_BLOCK6_KERNEL))
425
426
427
428
429
430
431
           ! 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
432
             call quad_hh_trafo_real_avx_avx2_4hv(a(1,j+off+a_off-3,istripe,my_thread), w, &
433
434
                                                  nbw, nl, stripe_width, nbw)
#else
435
             call quad_hh_trafo_real_avx_avx2_4hv(a(1,j+off+a_off-3,istripe), w, &
436
437
438
439
440
441
442
                                                  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
443
             call double_hh_trafo_real_avx_avx2_2hv(a(1,jj+off+a_off-1,istripe,my_thread), &
444
445
                                                    w, nbw, nl, stripe_width, nbw)
#else
446
             call double_hh_trafo_real_avx_avx2_2hv(a(1,jj+off+a_off-1,istripe), &
447
448
449
450
451
452
453
454
455
456
                                                    w, nbw, nl, stripe_width, nbw)
#endif
           enddo
#ifdef WITH_OPENMP
           if (jj==1) call single_hh_trafo_real_cpu_openmp(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(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
457
458
459

#endif /* defined(WITH_NO_SPECIFIC_REAL_KERNEL) || (defined(WITH_ONE_SPECIFIC_REAL_KERNEL) && !defined(WITH_REAL_AVX_BLOCK6_KERNEL) && !defined(WITH_REAL_AVX2_BLOCK6_KERNEL)) */

460
461
462
463
464
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
         endif
#endif /* WITH_NO_SPECIFIC_REAL_KERNEL */
#endif /* WITH_REAL_AVX_BLOCK4_KERNEL */

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
#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
             call hexa_hh_trafo_real_sse_6hv(a(1,j+off+a_off-5,istripe,my_thread), w, &
                                                  nbw, nl, stripe_width, nbw)
#else
             call hexa_hh_trafo_real_sse_6hv(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_4hv(a(1,jj+off+a_off-3,istripe,my_thread), w, &
                                                  nbw, nl, stripe_width, nbw)
#else
             call quad_hh_trafo_real_sse_4hv(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_2hv(a(1,jjj+off+a_off-1,istripe,my_thread), &
                                                    w, nbw, nl, stripe_width, nbw)
#else
             call double_hh_trafo_real_sse_2hv(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(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(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_SSE_BLOCK4_KERNEL */
520

521
#if defined(WITH_REAL_AVX_BLOCK6_KERNEL) || defined(WITH_REAL_AVX2_BLOCK6_KERNEL)
522
#if defined(WITH_NO_SPECIFIC_REAL_KERNEL)
523
524
         if ((THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX_BLOCK6) .or. &
             (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_AVX2_BLOCK6)) then
525
526
527
528
529
530
531
532
533
534
#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
535
             call hexa_hh_trafo_real_avx_avx2_6hv(a(1,j+off+a_off-5,istripe,my_thread), w, &
536
537
                                                  nbw, nl, stripe_width, nbw)
#else
538
             call hexa_hh_trafo_real_avx_avx2_6hv(a(1,j+off+a_off-5,istripe), w, &
539
540
541
542
543
544
545
546
547
                                                  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
548
             call quad_hh_trafo_real_avx_avx2_4hv(a(1,jj+off+a_off-3,istripe,my_thread), w, &
549
550
                                                  nbw, nl, stripe_width, nbw)
#else
551
             call quad_hh_trafo_real_avx_avx2_4hv(a(1,jj+off+a_off-3,istripe), w, &
552
553
554
555
556
557
558
                                                  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
559
             call double_hh_trafo_real_avx_avx2_2hv(a(1,jjj+off+a_off-1,istripe,my_thread), &
560
561
                                                    w, nbw, nl, stripe_width, nbw)
#else
562
             call double_hh_trafo_real_avx_avx2_2hv(a(1,jjj+off+a_off-1,istripe), &
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
                                                    w, nbw, nl, stripe_width, nbw)
#endif
           enddo
#ifdef WITH_OPENMP
           if (jjj==1) call single_hh_trafo_real_cpu_openmp(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(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 */

#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")
#else
         call timer%stop("compute_hh_trafo_real_cpu")
#endif
#endif

#ifdef WITH_OPENMP
       end subroutine compute_hh_trafo_real_cpu_openmp
#else
       end subroutine compute_hh_trafo_real_cpu
#endif

end module