elpa2_kernels_real_bgq.f90 21.4 KB
Newer Older
1
2
!    This file is part of ELPA.
!
Andreas Marek's avatar
Andreas Marek committed
3
!    The ELPA library was originally created by the ELPA consortium,
4
5
!    consisting of the following organizations:
!
6
7
!    - Max Planck Computing and Data Facility (MPCDF), formerly known as
!      Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
8
9
10
!    - Bergische Universität Wuppertal, Lehrstuhl für angewandte
!      Informatik,
!    - Technische Universität München, Lehrstuhl für Informatik mit
Andreas Marek's avatar
Andreas Marek committed
11
12
13
14
15
!      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
16
17
18
19
!    - IBM Deutschland GmbH
!
!
!    More information can be found here:
20
!    http://elpa.mpcdf.mpg.de/
21
22
!
!    ELPA is free software: you can redistribute it and/or modify
Andreas Marek's avatar
Andreas Marek committed
23
24
!    it under the terms of the version 3 of the license of the
!    GNU Lesser General Public License as published by the Free
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
!    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 contains the compute intensive kernels for the Householder transformations.
!
! *** Special IBM BlueGene/Q version with QPX intrinsics in Fortran ***
Andreas Marek's avatar
Andreas Marek committed
48
!
49
50
51
52
53
54
! Copyright of the original code rests with the authors inside the ELPA
! consortium. The copyright of any additional modifications shall rest
! with their original authors, but shall adhere to the licensing terms
! distributed along with the original code in the file "COPYING".
!
! --------------------------------------------------------------------------------------------------
Andreas Marek's avatar
Andreas Marek committed
55
module real_bgq_kernel
56

Andreas Marek's avatar
Andreas Marek committed
57
58
59
60
  private
  public double_hh_trafo_bgq
contains
  subroutine double_hh_trafo_bgq(q, hh, nb, nq, ldq, ldh)
Andreas Marek's avatar
Andreas Marek committed
61
    use precision
Andreas Marek's avatar
Andreas Marek committed
62
63
    implicit none

Andreas Marek's avatar
Andreas Marek committed
64
65
66
    integer(kind=ik), intent(in) :: nb, nq, ldq, ldh
    real(kind=rk), intent(inout) :: q(ldq,*)
    real(kind=rk), intent(in)    :: hh(ldh,*)
Andreas Marek's avatar
Andreas Marek committed
67

Andreas Marek's avatar
Andreas Marek committed
68
69
    real(kind=rk)                :: s
    integer(kind=ik)             :: i
Andreas Marek's avatar
Andreas Marek committed
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110

    ! Safety only:

    if(mod(ldq,4) /= 0) STOP 'double_hh_trafo: ldq not divisible by 4!'

    call alignx(32,q)

    ! Calculate dot product of the two Householder vectors

    s = hh(2,2)*1
    do i=3,nb
       s = s+hh(i,2)*hh(i-1,1)
    enddo

    do i=1,nq-20,24
       call hh_trafo_kernel_24_bgq(q(i   ,1), hh, nb, ldq, ldh, s)
    enddo

    if(nq-i+1 > 16) then
       call hh_trafo_kernel_16_bgq(q(i  ,1), hh, nb, ldq, ldh, s)
       call hh_trafo_kernel_4_bgq(q(i+16,1), hh, nb, ldq, ldh, s)
    else if(nq-i+1 > 12) then
       call hh_trafo_kernel_8_bgq(q(i  ,1), hh, nb, ldq, ldh, s)
       call hh_trafo_kernel_8_bgq(q(i+8,1), hh, nb, ldq, ldh, s)
    else if(nq-i+1 > 8) then
       call hh_trafo_kernel_8_bgq(q(i  ,1), hh, nb, ldq, ldh, s)
       call hh_trafo_kernel_4_bgq(q(i+8,1), hh, nb, ldq, ldh, s)
    else if(nq-i+1 > 4) then
       call hh_trafo_kernel_8_bgq(q(i  ,1), hh, nb, ldq, ldh, s)
    else if(nq-i+1 > 0) then
       call hh_trafo_kernel_4_bgq(q(i  ,1), hh, nb, ldq, ldh, s)
    endif

  end subroutine double_hh_trafo_bgq


  ! --------------------------------------------------------------------------------------------------
  ! The following kernels perform the Householder transformation on Q for 24/16/8/4 rows.
  ! --------------------------------------------------------------------------------------------------

  subroutine hh_trafo_kernel_24_bgq(q, hh, nb, ldq, ldh, s)
Andreas Marek's avatar
Andreas Marek committed
111
    use precision
112
    use elpa_mpi
Andreas Marek's avatar
Andreas Marek committed
113
114
    implicit none

Andreas Marek's avatar
Andreas Marek committed
115
    integer(kind=ik), intent(in) :: nb, ldq, ldh
Andreas Marek's avatar
Andreas Marek committed
116

Andreas Marek's avatar
Andreas Marek committed
117
118
    real(kind=rk), intent(inout) :: q(ldq,*)
    real(kind=rk), intent(in)    :: hh(ldh,*), s
Andreas Marek's avatar
Andreas Marek committed
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317

    VECTOR(REAL(8))::QPX_x1, QPX_x2, QPX_x3, QPX_x4, QPX_x5, QPX_x6
    VECTOR(REAL(8))::QPX_y1, QPX_y2, QPX_y3, QPX_y4, QPX_y5, QPX_y6
    VECTOR(REAL(8))::QPX_q1, QPX_q2, QPX_q3, QPX_q4, QPX_q5, QPX_q6
    VECTOR(REAL(8))::QPX_h1, QPX_h2, QPX_tau1, QPX_tau2, QPX_s
    integer i

    call alignx(32,q)

    !--- multiply Householder vectors with matrix q ---

    QPX_x1 = VEC_LD(0,q(1,2))
    QPX_x2 = VEC_LD(0,q(5,2))
    QPX_x3 = VEC_LD(0,q(9,2))
    QPX_x4 = VEC_LD(0,q(13,2))
    QPX_x5 = VEC_LD(0,q(17,2))
    QPX_x6 = VEC_LD(0,q(21,2))

    QPX_h2 = VEC_SPLATS(hh(2,2))
    QPX_q1 = VEC_LD(0,q(1,1))
    QPX_q2 = VEC_LD(0,q(5,1))
    QPX_q3 = VEC_LD(0,q(9,1))
    QPX_q4 = VEC_LD(0,q(13,1))
    QPX_q5 = VEC_LD(0,q(17,1))
    QPX_q6 = VEC_LD(0,q(21,1))
    QPX_y1 = VEC_MADD(QPX_x1, QPX_h2, QPX_q1)
    QPX_y2 = VEC_MADD(QPX_x2, QPX_h2, QPX_q2)
    QPX_y3 = VEC_MADD(QPX_x3, QPX_h2, QPX_q3)
    QPX_y4 = VEC_MADD(QPX_x4, QPX_h2, QPX_q4)
    QPX_y5 = VEC_MADD(QPX_x5, QPX_h2, QPX_q5)
    QPX_y6 = VEC_MADD(QPX_x6, QPX_h2, QPX_q6)

    do i=3,nb,1

       QPX_q1 = VEC_LD(0,q(1,i))
       QPX_q2 = VEC_LD(0,q(5,i))
       QPX_q3 = VEC_LD(0,q(9,i))
       QPX_q4 = VEC_LD(0,q(13,i))
       QPX_q5 = VEC_LD(0,q(17,i))
       QPX_q6 = VEC_LD(0,q(21,i))
       QPX_h1 = VEC_SPLATS(hh(i-1,1))
       QPX_x1 = VEC_MADD(QPX_q1, QPX_h1, QPX_x1)
       QPX_x2 = VEC_MADD(QPX_q2, QPX_h1, QPX_x2)
       QPX_x3 = VEC_MADD(QPX_q3, QPX_h1, QPX_x3)
       QPX_x4 = VEC_MADD(QPX_q4, QPX_h1, QPX_x4)
       QPX_x5 = VEC_MADD(QPX_q5, QPX_h1, QPX_x5)
       QPX_x6 = VEC_MADD(QPX_q6, QPX_h1, QPX_x6)
       QPX_h2 = VEC_SPLATS(hh(i,2))
       QPX_y1 = VEC_MADD(QPX_q1, QPX_h2, QPX_y1)
       QPX_y2 = VEC_MADD(QPX_q2, QPX_h2, QPX_y2)
       QPX_y3 = VEC_MADD(QPX_q3, QPX_h2, QPX_y3)
       QPX_y4 = VEC_MADD(QPX_q4, QPX_h2, QPX_y4)
       QPX_y5 = VEC_MADD(QPX_q5, QPX_h2, QPX_y5)
       QPX_y6 = VEC_MADD(QPX_q6, QPX_h2, QPX_y6)

    enddo

    QPX_h1 = VEC_SPLATS(hh(nb,1))
    QPX_q1 = VEC_LD(0,q(1,nb+1))
    QPX_q2 = VEC_LD(0,q(5,nb+1))
    QPX_q3 = VEC_LD(0,q(9,nb+1))
    QPX_q4 = VEC_LD(0,q(13,nb+1))
    QPX_q5 = VEC_LD(0,q(17,nb+1))
    QPX_q6 = VEC_LD(0,q(21,nb+1))
    QPX_x1 = VEC_MADD(QPX_q1, QPX_h1, QPX_x1)
    QPX_x2 = VEC_MADD(QPX_q2, QPX_h1, QPX_x2)
    QPX_x3 = VEC_MADD(QPX_q3, QPX_h1, QPX_x3)
    QPX_x4 = VEC_MADD(QPX_q4, QPX_h1, QPX_x4)
    QPX_x5 = VEC_MADD(QPX_q5, QPX_h1, QPX_x5)
    QPX_x6 = VEC_MADD(QPX_q6, QPX_h1, QPX_x6)

    !--- multiply T matrix ---

    QPX_tau1 = VEC_SPLATS(-hh(1,1))
    QPX_x1 = VEC_MUL(QPX_x1, QPX_tau1)
    QPX_x2 = VEC_MUL(QPX_x2, QPX_tau1)
    QPX_x3 = VEC_MUL(QPX_x3, QPX_tau1)
    QPX_x4 = VEC_MUL(QPX_x4, QPX_tau1)
    QPX_x5 = VEC_MUL(QPX_x5, QPX_tau1)
    QPX_x6 = VEC_MUL(QPX_x6, QPX_tau1)
    QPX_tau2 = VEC_SPLATS(-hh(1,2))
    QPX_s = VEC_SPLATS(-hh(1,2)*s)
    QPX_y1 = VEC_MUL(QPX_y1, QPX_tau2)
    QPX_y2 = VEC_MUL(QPX_y2, QPX_tau2)
    QPX_y3 = VEC_MUL(QPX_y3, QPX_tau2)
    QPX_y4 = VEC_MUL(QPX_y4, QPX_tau2)
    QPX_y5 = VEC_MUL(QPX_y5, QPX_tau2)
    QPX_y6 = VEC_MUL(QPX_y6, QPX_tau2)
    QPX_y1 = VEC_MADD(QPX_x1, QPX_s, QPX_y1)
    QPX_y2 = VEC_MADD(QPX_x2, QPX_s, QPX_y2)
    QPX_y3 = VEC_MADD(QPX_x3, QPX_s, QPX_y3)
    QPX_y4 = VEC_MADD(QPX_x4, QPX_s, QPX_y4)
    QPX_y5 = VEC_MADD(QPX_x5, QPX_s, QPX_y5)
    QPX_y6 = VEC_MADD(QPX_x6, QPX_s, QPX_y6)

    !--- rank-2 update of q ---

    QPX_q1 = VEC_LD(0,q(1,1))
    QPX_q2 = VEC_LD(0,q(5,1))
    QPX_q3 = VEC_LD(0,q(9,1))
    QPX_q4 = VEC_LD(0,q(13,1))
    QPX_q5 = VEC_LD(0,q(17,1))
    QPX_q6 = VEC_LD(0,q(21,1))
    QPX_q1 = VEC_ADD(QPX_q1, QPX_y1)
    QPX_q2 = VEC_ADD(QPX_q2, QPX_y2)
    QPX_q3 = VEC_ADD(QPX_q3, QPX_y3)
    QPX_q4 = VEC_ADD(QPX_q4, QPX_y4)
    QPX_q5 = VEC_ADD(QPX_q5, QPX_y5)
    QPX_q6 = VEC_ADD(QPX_q6, QPX_y6)
    call VEC_ST(QPX_q1, 0, q(1,1))
    call VEC_ST(QPX_q2, 0, q(5,1))
    call VEC_ST(QPX_q3, 0, q(9,1))
    call VEC_ST(QPX_q4, 0, q(13,1))
    call VEC_ST(QPX_q5, 0, q(17,1))
    call VEC_ST(QPX_q6, 0, q(21,1))

    QPX_h2 = VEC_SPLATS(hh(2,2))
    QPX_q1 = VEC_LD(0,q(1,2))
    QPX_q2 = VEC_LD(0,q(5,2))
    QPX_q3 = VEC_LD(0,q(9,2))
    QPX_q4 = VEC_LD(0,q(13,2))
    QPX_q5 = VEC_LD(0,q(17,2))
    QPX_q6 = VEC_LD(0,q(21,2))
    QPX_q1 = VEC_MADD(QPX_y1, QPX_h2, QPX_q1)
    QPX_q2 = VEC_MADD(QPX_y2, QPX_h2, QPX_q2)
    QPX_q3 = VEC_MADD(QPX_y3, QPX_h2, QPX_q3)
    QPX_q4 = VEC_MADD(QPX_y4, QPX_h2, QPX_q4)
    QPX_q5 = VEC_MADD(QPX_y5, QPX_h2, QPX_q5)
    QPX_q6 = VEC_MADD(QPX_y6, QPX_h2, QPX_q6)
    QPX_q1 = VEC_ADD(QPX_q1, QPX_x1)
    QPX_q2 = VEC_ADD(QPX_q2, QPX_x2)
    QPX_q3 = VEC_ADD(QPX_q3, QPX_x3)
    QPX_q4 = VEC_ADD(QPX_q4, QPX_x4)
    QPX_q5 = VEC_ADD(QPX_q5, QPX_x5)
    QPX_q6 = VEC_ADD(QPX_q6, QPX_x6)
    call VEC_ST(QPX_q1, 0, q(1,2))
    call VEC_ST(QPX_q2, 0, q(5,2))
    call VEC_ST(QPX_q3, 0, q(9,2))
    call VEC_ST(QPX_q4, 0, q(13,2))
    call VEC_ST(QPX_q5, 0, q(17,2))
    call VEC_ST(QPX_q6, 0, q(21,2))

    do i=3,nb,1

       QPX_q1 = VEC_LD(0,q(1,i))
       QPX_q2 = VEC_LD(0,q(5,i))
       QPX_q3 = VEC_LD(0,q(9,i))
       QPX_q4 = VEC_LD(0,q(13,i))
       QPX_q5 = VEC_LD(0,q(17,i))
       QPX_q6 = VEC_LD(0,q(21,i))
       QPX_h1 = VEC_SPLATS(hh(i-1,1))
       QPX_q1 = VEC_MADD(QPX_x1, QPX_h1, QPX_q1)
       QPX_q2 = VEC_MADD(QPX_x2, QPX_h1, QPX_q2)
       QPX_q3 = VEC_MADD(QPX_x3, QPX_h1, QPX_q3)
       QPX_q4 = VEC_MADD(QPX_x4, QPX_h1, QPX_q4)
       QPX_q5 = VEC_MADD(QPX_x5, QPX_h1, QPX_q5)
       QPX_q6 = VEC_MADD(QPX_x6, QPX_h1, QPX_q6)
       QPX_h2 = VEC_SPLATS(hh(i,2))
       QPX_q1 = VEC_MADD(QPX_y1, QPX_h2, QPX_q1)
       QPX_q2 = VEC_MADD(QPX_y2, QPX_h2, QPX_q2)
       QPX_q3 = VEC_MADD(QPX_y3, QPX_h2, QPX_q3)
       QPX_q4 = VEC_MADD(QPX_y4, QPX_h2, QPX_q4)
       QPX_q5 = VEC_MADD(QPX_y5, QPX_h2, QPX_q5)
       QPX_q6 = VEC_MADD(QPX_y6, QPX_h2, QPX_q6)

       call VEC_ST(QPX_q1, 0, q(1,i))
       call VEC_ST(QPX_q2, 0, q(5,i))
       call VEC_ST(QPX_q3, 0, q(9,i))
       call VEC_ST(QPX_q4, 0, q(13,i))
       call VEC_ST(QPX_q5, 0, q(17,i))
       call VEC_ST(QPX_q6, 0, q(21,i))

    enddo

    QPX_h1 = VEC_SPLATS(hh(nb,1))
    QPX_q1 = VEC_LD(0,q(1,nb+1))
    QPX_q2 = VEC_LD(0,q(5,nb+1))
    QPX_q3 = VEC_LD(0,q(9,nb+1))
    QPX_q4 = VEC_LD(0,q(13,nb+1))
    QPX_q5 = VEC_LD(0,q(17,nb+1))
    QPX_q6 = VEC_LD(0,q(21,nb+1))
    QPX_q1 = VEC_MADD(QPX_x1, QPX_h1, QPX_q1)
    QPX_q2 = VEC_MADD(QPX_x2, QPX_h1, QPX_q2)
    QPX_q3 = VEC_MADD(QPX_x3, QPX_h1, QPX_q3)
    QPX_q4 = VEC_MADD(QPX_x4, QPX_h1, QPX_q4)
    QPX_q5 = VEC_MADD(QPX_x5, QPX_h1, QPX_q5)
    QPX_q6 = VEC_MADD(QPX_x6, QPX_h1, QPX_q6)
    call VEC_ST(QPX_q1, 0, q(1,nb+1))
    call VEC_ST(QPX_q2, 0, q(5,nb+1))
    call VEC_ST(QPX_q3, 0, q(9,nb+1))
    call VEC_ST(QPX_q4, 0, q(13,nb+1))
    call VEC_ST(QPX_q5, 0, q(17,nb+1))
    call VEC_ST(QPX_q6, 0, q(21,nb+1))

  end subroutine hh_trafo_kernel_24_bgq

  ! --------------------------------------------------------------------------------------------------

  subroutine hh_trafo_kernel_16_bgq(q, hh, nb, ldq, ldh, s)
Andreas Marek's avatar
Andreas Marek committed
318
    use precision
319
    use elpa_mpi
Andreas Marek's avatar
Andreas Marek committed
320
321
    implicit none

Andreas Marek's avatar
Andreas Marek committed
322
    integer(kind=ik), intent(in) :: nb, ldq, ldh
Andreas Marek's avatar
Andreas Marek committed
323

Andreas Marek's avatar
Andreas Marek committed
324
325
    real(kind=rk), intent(inout) :: q(ldq,*)
    real(kind=rk), intent(in)    :: hh(ldh,*), s
Andreas Marek's avatar
Andreas Marek committed
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

    VECTOR(REAL(8))::QPX_x1, QPX_x2, QPX_x3, QPX_x4
    VECTOR(REAL(8))::QPX_y1, QPX_y2, QPX_y3, QPX_y4
    VECTOR(REAL(8))::QPX_q1, QPX_q2, QPX_q3, QPX_q4
    VECTOR(REAL(8))::QPX_h1, QPX_h2, QPX_tau1, QPX_tau2, QPX_s
    integer i

    call alignx(32,q)

    !--- multiply Householder vectors with matrix q ---

    QPX_x1 = VEC_LD(0,q(1,2))
    QPX_x2 = VEC_LD(0,q(5,2))
    QPX_x3 = VEC_LD(0,q(9,2))
    QPX_x4 = VEC_LD(0,q(13,2))

    QPX_h2 = VEC_SPLATS(hh(2,2))
    QPX_q1 = VEC_LD(0,q(1,1))
    QPX_q2 = VEC_LD(0,q(5,1))
    QPX_q3 = VEC_LD(0,q(9,1))
    QPX_q4 = VEC_LD(0,q(13,1))
    QPX_y1 = VEC_MADD(QPX_x1, QPX_h2, QPX_q1)
    QPX_y2 = VEC_MADD(QPX_x2, QPX_h2, QPX_q2)
    QPX_y3 = VEC_MADD(QPX_x3, QPX_h2, QPX_q3)
    QPX_y4 = VEC_MADD(QPX_x4, QPX_h2, QPX_q4)

    do i=3,nb,1

       QPX_q1 = VEC_LD(0,q(1,i))
       QPX_q2 = VEC_LD(0,q(5,i))
       QPX_q3 = VEC_LD(0,q(9,i))
       QPX_q4 = VEC_LD(0,q(13,i))
       QPX_h1 = VEC_SPLATS(hh(i-1,1))
       QPX_x1 = VEC_MADD(QPX_q1, QPX_h1, QPX_x1)
       QPX_x2 = VEC_MADD(QPX_q2, QPX_h1, QPX_x2)
       QPX_x3 = VEC_MADD(QPX_q3, QPX_h1, QPX_x3)
       QPX_x4 = VEC_MADD(QPX_q4, QPX_h1, QPX_x4)
       QPX_h2 = VEC_SPLATS(hh(i,2))
       QPX_y1 = VEC_MADD(QPX_q1, QPX_h2, QPX_y1)
       QPX_y2 = VEC_MADD(QPX_q2, QPX_h2, QPX_y2)
       QPX_y3 = VEC_MADD(QPX_q3, QPX_h2, QPX_y3)
       QPX_y4 = VEC_MADD(QPX_q4, QPX_h2, QPX_y4)

    enddo

    QPX_h1 = VEC_SPLATS(hh(nb,1))
    QPX_q1 = VEC_LD(0,q(1,nb+1))
    QPX_q2 = VEC_LD(0,q(5,nb+1))
    QPX_q3 = VEC_LD(0,q(9,nb+1))
    QPX_q4 = VEC_LD(0,q(13,nb+1))
    QPX_x1 = VEC_MADD(QPX_q1, QPX_h1, QPX_x1)
    QPX_x2 = VEC_MADD(QPX_q2, QPX_h1, QPX_x2)
    QPX_x3 = VEC_MADD(QPX_q3, QPX_h1, QPX_x3)
    QPX_x4 = VEC_MADD(QPX_q4, QPX_h1, QPX_x4)

    !--- multiply T matrix ---

    QPX_tau1 = VEC_SPLATS(-hh(1,1))
    QPX_x1 = VEC_MUL(QPX_x1, QPX_tau1)
    QPX_x2 = VEC_MUL(QPX_x2, QPX_tau1)
    QPX_x3 = VEC_MUL(QPX_x3, QPX_tau1)
    QPX_x4 = VEC_MUL(QPX_x4, QPX_tau1)
    QPX_tau2 = VEC_SPLATS(-hh(1,2))
    QPX_s = VEC_SPLATS(-hh(1,2)*s)
    QPX_y1 = VEC_MUL(QPX_y1, QPX_tau2)
    QPX_y2 = VEC_MUL(QPX_y2, QPX_tau2)
    QPX_y3 = VEC_MUL(QPX_y3, QPX_tau2)
    QPX_y4 = VEC_MUL(QPX_y4, QPX_tau2)
    QPX_y1 = VEC_MADD(QPX_x1, QPX_s, QPX_y1)
    QPX_y2 = VEC_MADD(QPX_x2, QPX_s, QPX_y2)
    QPX_y3 = VEC_MADD(QPX_x3, QPX_s, QPX_y3)
    QPX_y4 = VEC_MADD(QPX_x4, QPX_s, QPX_y4)

    !--- rank-2 update of q ---

    QPX_q1 = VEC_LD(0,q(1,1))
    QPX_q2 = VEC_LD(0,q(5,1))
    QPX_q3 = VEC_LD(0,q(9,1))
    QPX_q4 = VEC_LD(0,q(13,1))
    QPX_q1 = VEC_ADD(QPX_q1, QPX_y1)
    QPX_q2 = VEC_ADD(QPX_q2, QPX_y2)
    QPX_q3 = VEC_ADD(QPX_q3, QPX_y3)
    QPX_q4 = VEC_ADD(QPX_q4, QPX_y4)
    call VEC_ST(QPX_q1, 0, q(1,1))
    call VEC_ST(QPX_q2, 0, q(5,1))
    call VEC_ST(QPX_q3, 0, q(9,1))
    call VEC_ST(QPX_q4, 0, q(13,1))

    QPX_h2 = VEC_SPLATS(hh(2,2))
    QPX_q1 = VEC_LD(0,q(1,2))
    QPX_q2 = VEC_LD(0,q(5,2))
    QPX_q3 = VEC_LD(0,q(9,2))
    QPX_q4 = VEC_LD(0,q(13,2))
    QPX_q1 = VEC_MADD(QPX_y1, QPX_h2, QPX_q1)
    QPX_q2 = VEC_MADD(QPX_y2, QPX_h2, QPX_q2)
    QPX_q3 = VEC_MADD(QPX_y3, QPX_h2, QPX_q3)
    QPX_q4 = VEC_MADD(QPX_y4, QPX_h2, QPX_q4)
    QPX_q1 = VEC_ADD(QPX_q1, QPX_x1)
    QPX_q2 = VEC_ADD(QPX_q2, QPX_x2)
    QPX_q3 = VEC_ADD(QPX_q3, QPX_x3)
    QPX_q4 = VEC_ADD(QPX_q4, QPX_x4)
    call VEC_ST(QPX_q1, 0, q(1,2))
    call VEC_ST(QPX_q2, 0, q(5,2))
    call VEC_ST(QPX_q3, 0, q(9,2))
    call VEC_ST(QPX_q4, 0, q(13,2))

    do i=3,nb,1

       QPX_q1 = VEC_LD(0,q(1,i))
       QPX_q2 = VEC_LD(0,q(5,i))
       QPX_q3 = VEC_LD(0,q(9,i))
       QPX_q4 = VEC_LD(0,q(13,i))
       QPX_h1 = VEC_SPLATS(hh(i-1,1))
       QPX_q1 = VEC_MADD(QPX_x1, QPX_h1, QPX_q1)
       QPX_q2 = VEC_MADD(QPX_x2, QPX_h1, QPX_q2)
       QPX_q3 = VEC_MADD(QPX_x3, QPX_h1, QPX_q3)
       QPX_q4 = VEC_MADD(QPX_x4, QPX_h1, QPX_q4)
       QPX_h2 = VEC_SPLATS(hh(i,2))
       QPX_q1 = VEC_MADD(QPX_y1, QPX_h2, QPX_q1)
       QPX_q2 = VEC_MADD(QPX_y2, QPX_h2, QPX_q2)
       QPX_q3 = VEC_MADD(QPX_y3, QPX_h2, QPX_q3)
       QPX_q4 = VEC_MADD(QPX_y4, QPX_h2, QPX_q4)

       call VEC_ST(QPX_q1, 0, q(1,i))
       call VEC_ST(QPX_q2, 0, q(5,i))
       call VEC_ST(QPX_q3, 0, q(9,i))
       call VEC_ST(QPX_q4, 0, q(13,i))

    enddo

    QPX_h1 = VEC_SPLATS(hh(nb,1))
    QPX_q1 = VEC_LD(0,q(1,nb+1))
    QPX_q2 = VEC_LD(0,q(5,nb+1))
    QPX_q3 = VEC_LD(0,q(9,nb+1))
    QPX_q4 = VEC_LD(0,q(13,nb+1))
    QPX_q1 = VEC_MADD(QPX_x1, QPX_h1, QPX_q1)
    QPX_q2 = VEC_MADD(QPX_x2, QPX_h1, QPX_q2)
    QPX_q3 = VEC_MADD(QPX_x3, QPX_h1, QPX_q3)
    QPX_q4 = VEC_MADD(QPX_x4, QPX_h1, QPX_q4)
    call VEC_ST(QPX_q1, 0, q(1,nb+1))
    call VEC_ST(QPX_q2, 0, q(5,nb+1))
    call VEC_ST(QPX_q3, 0, q(9,nb+1))
    call VEC_ST(QPX_q4, 0, q(13,nb+1))

  end subroutine hh_trafo_kernel_16_bgq

  ! --------------------------------------------------------------------------------------------------

  subroutine hh_trafo_kernel_8_bgq(q, hh, nb, ldq, ldh, s)
Andreas Marek's avatar
Andreas Marek committed
475
    use precision
476
    use elpa_mpi
Andreas Marek's avatar
Andreas Marek committed
477
478
    implicit none

Andreas Marek's avatar
Andreas Marek committed
479
    integer(kind=ik), intent(in) :: nb, ldq, ldh
Andreas Marek's avatar
Andreas Marek committed
480

Andreas Marek's avatar
Andreas Marek committed
481
482
483
    real(kind=rk), intent(inout) :: q(ldq,*)
    real(kind=rk), intent(in)    :: hh(ldh,*), s
    integer(kind=ik)             :: i
Andreas Marek's avatar
Andreas Marek committed
484
485
486
    VECTOR(REAL(8))::QPX_x1, QPX_x2, QPX_y1, QPX_y2
    VECTOR(REAL(8))::QPX_q1, QPX_q2
    VECTOR(REAL(8))::QPX_h1, QPX_h2, QPX_tau1, QPX_tau2, QPX_s
Andreas Marek's avatar
Andreas Marek committed
487

Andreas Marek's avatar
Andreas Marek committed
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

    call alignx(32,q)

    !--- multiply Householder vectors with matrix q ---

    QPX_x1 = VEC_LD(0,q(1,2))
    QPX_x2 = VEC_LD(0,q(5,2))

    QPX_h2 = VEC_SPLATS(hh(2,2))
    QPX_q1 = VEC_LD(0,q(1,1))
    QPX_q2 = VEC_LD(0,q(5,1))
    QPX_y1 = VEC_MADD(QPX_x1, QPX_h2, QPX_q1)
    QPX_y2 = VEC_MADD(QPX_x2, QPX_h2, QPX_q2)

    do i=3,nb,1

       QPX_q1 = VEC_LD(0,q(1,i))
       QPX_q2 = VEC_LD(0,q(5,i))
       QPX_h1 = VEC_SPLATS(hh(i-1,1))
       QPX_x1 = VEC_MADD(QPX_q1, QPX_h1, QPX_x1)
       QPX_x2 = VEC_MADD(QPX_q2, QPX_h1, QPX_x2)
       QPX_h2 = VEC_SPLATS(hh(i,2))
       QPX_y1 = VEC_MADD(QPX_q1, QPX_h2, QPX_y1)
       QPX_y2 = VEC_MADD(QPX_q2, QPX_h2, QPX_y2)

    enddo

    QPX_h1 = VEC_SPLATS(hh(nb,1))
    QPX_q1 = VEC_LD(0,q(1,nb+1))
    QPX_q2 = VEC_LD(0,q(5,nb+1))
    QPX_x1 = VEC_MADD(QPX_q1, QPX_h1, QPX_x1)
    QPX_x2 = VEC_MADD(QPX_q2, QPX_h1, QPX_x2)

    !--- multiply T matrix ---

    QPX_tau1 = VEC_SPLATS(-hh(1,1))
    QPX_x1 = VEC_MUL(QPX_x1, QPX_tau1)
    QPX_x2 = VEC_MUL(QPX_x2, QPX_tau1)
    QPX_tau2 = VEC_SPLATS(-hh(1,2))
    QPX_s = VEC_SPLATS(-hh(1,2)*s)
    QPX_y1 = VEC_MUL(QPX_y1, QPX_tau2)
    QPX_y2 = VEC_MUL(QPX_y2, QPX_tau2)
    QPX_y1 = VEC_MADD(QPX_x1, QPX_s, QPX_y1)
    QPX_y2 = VEC_MADD(QPX_x2, QPX_s, QPX_y2)
532

Andreas Marek's avatar
Andreas Marek committed
533
    !--- rank-2 update of q ---
534

Andreas Marek's avatar
Andreas Marek committed
535
536
537
538
539
540
    QPX_q1 = VEC_LD(0,q(1,1))
    QPX_q2 = VEC_LD(0,q(5,1))
    QPX_q1 = VEC_ADD(QPX_q1, QPX_y1)
    QPX_q2 = VEC_ADD(QPX_q2, QPX_y2)
    call VEC_ST(QPX_q1, 0, q(1,1))
    call VEC_ST(QPX_q2, 0, q(5,1))
541

Andreas Marek's avatar
Andreas Marek committed
542
543
544
545
546
547
548
549
550
    QPX_h2 = VEC_SPLATS(hh(2,2))
    QPX_q1 = VEC_LD(0,q(1,2))
    QPX_q2 = VEC_LD(0,q(5,2))
    QPX_q1 = VEC_MADD(QPX_y1, QPX_h2, QPX_q1)
    QPX_q2 = VEC_MADD(QPX_y2, QPX_h2, QPX_q2)
    QPX_q1 = VEC_ADD(QPX_q1, QPX_x1)
    QPX_q2 = VEC_ADD(QPX_q2, QPX_x2)
    call VEC_ST(QPX_q1, 0, q(1,2))
    call VEC_ST(QPX_q2, 0, q(5,2))
551

Andreas Marek's avatar
Andreas Marek committed
552
    do i=3,nb,1
553

Andreas Marek's avatar
Andreas Marek committed
554
555
556
557
558
559
560
561
       QPX_q1 = VEC_LD(0,q(1,i))
       QPX_q2 = VEC_LD(0,q(5,i))
       QPX_h1 = VEC_SPLATS(hh(i-1,1))
       QPX_q1 = VEC_MADD(QPX_x1, QPX_h1, QPX_q1)
       QPX_q2 = VEC_MADD(QPX_x2, QPX_h1, QPX_q2)
       QPX_h2 = VEC_SPLATS(hh(i,2))
       QPX_q1 = VEC_MADD(QPX_y1, QPX_h2, QPX_q1)
       QPX_q2 = VEC_MADD(QPX_y2, QPX_h2, QPX_q2)
562

Andreas Marek's avatar
Andreas Marek committed
563
564
       call VEC_ST(QPX_q1, 0, q(1,i))
       call VEC_ST(QPX_q2, 0, q(5,i))
565

Andreas Marek's avatar
Andreas Marek committed
566
    enddo
567

Andreas Marek's avatar
Andreas Marek committed
568
569
570
571
572
573
574
    QPX_h1 = VEC_SPLATS(hh(nb,1))
    QPX_q1 = VEC_LD(0,q(1,nb+1))
    QPX_q2 = VEC_LD(0,q(5,nb+1))
    QPX_q1 = VEC_MADD(QPX_x1, QPX_h1, QPX_q1)
    QPX_q2 = VEC_MADD(QPX_x2, QPX_h1, QPX_q2)
    call VEC_ST(QPX_q1, 0, q(1,nb+1))
    call VEC_ST(QPX_q2, 0, q(5,nb+1))
575

Andreas Marek's avatar
Andreas Marek committed
576
  end subroutine hh_trafo_kernel_8_bgq
577

Andreas Marek's avatar
Andreas Marek committed
578
  ! --------------------------------------------------------------------------------------------------
579

Andreas Marek's avatar
Andreas Marek committed
580
  subroutine hh_trafo_kernel_4_bgq(q, hh, nb, ldq, ldh, s)
Andreas Marek's avatar
Andreas Marek committed
581
    use precision
582
    use elpa_mpi
Andreas Marek's avatar
Andreas Marek committed
583
    implicit none
584

Andreas Marek's avatar
Andreas Marek committed
585
    integer(kind=ik), intent(in) :: nb, ldq, ldh
Andreas Marek's avatar
Andreas Marek committed
586

Andreas Marek's avatar
Andreas Marek committed
587
588
589
    real(kind=rk), intent(inout) :: q(ldq,*)
    real(kind=rk), intent(in)    :: hh(ldh,*), s
    integer(kind=ik)             :: i
Andreas Marek's avatar
Andreas Marek committed
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
    VECTOR(REAL(8))::QPX_x1, QPX_y1
    VECTOR(REAL(8))::QPX_q1
    VECTOR(REAL(8))::QPX_h1, QPX_h2, QPX_tau1, QPX_tau2, QPX_s

    call alignx(32,q)

    !--- multiply Householder vectors with matrix q ---

    QPX_x1 = VEC_LD(0,q(1,2))

    QPX_h2 = VEC_SPLATS(hh(2,2))
    QPX_q1 = VEC_LD(0,q(1,1))
    QPX_y1 = VEC_MADD(QPX_x1, QPX_h2, QPX_q1)

    do i=3,nb,1

       QPX_q1 = VEC_LD(0,q(1,i))
       QPX_h1 = VEC_SPLATS(hh(i-1,1))
       QPX_x1 = VEC_MADD(QPX_q1, QPX_h1, QPX_x1)
       QPX_h2 = VEC_SPLATS(hh(i,2))
       QPX_y1 = VEC_MADD(QPX_q1, QPX_h2, QPX_y1)

    enddo

    QPX_h1 = VEC_SPLATS(hh(nb,1))
    QPX_q1 = VEC_LD(0,q(1,nb+1))
    QPX_x1 = VEC_MADD(QPX_q1, QPX_h1, QPX_x1)

    !--- multiply T matrix ---

    QPX_tau1 = VEC_SPLATS(-hh(1,1))
    QPX_x1 = VEC_MUL(QPX_x1, QPX_tau1)
    QPX_tau2 = VEC_SPLATS(-hh(1,2))
    QPX_s = VEC_SPLATS(-hh(1,2)*s)
    QPX_y1 = VEC_MUL(QPX_y1, QPX_tau2)
    QPX_y1 = VEC_MADD(QPX_x1, QPX_s, QPX_y1)

    !--- rank-2 update of q ---

    QPX_q1 = VEC_LD(0,q(1,1))
    QPX_q1 = VEC_ADD(QPX_q1, QPX_y1)
    call VEC_ST(QPX_q1, 0, q(1,1))

    QPX_h2 = VEC_SPLATS(hh(2,2))
    QPX_q1 = VEC_LD(0,q(1,2))
    QPX_q1 = VEC_MADD(QPX_y1, QPX_h2, QPX_q1)
    QPX_q1 = VEC_ADD(QPX_q1, QPX_x1)
    call VEC_ST(QPX_q1, 0, q(1,2))

    do i=3,nb,1

       QPX_q1 = VEC_LD(0,q(1,i))
       QPX_h1 = VEC_SPLATS(hh(i-1,1))
       QPX_q1 = VEC_MADD(QPX_x1, QPX_h1, QPX_q1)
       QPX_h2 = VEC_SPLATS(hh(i,2))
       QPX_q1 = VEC_MADD(QPX_y1, QPX_h2, QPX_q1)

       call VEC_ST(QPX_q1, 0, q(1,i))

    enddo
650

Andreas Marek's avatar
Andreas Marek committed
651
652
653
654
    QPX_h1 = VEC_SPLATS(hh(nb,1))
    QPX_q1 = VEC_LD(0,q(1,nb+1))
    QPX_q1 = VEC_MADD(QPX_x1, QPX_h1, QPX_q1)
    call VEC_ST(QPX_q1, 0, q(1,nb+1))
655

Andreas Marek's avatar
Andreas Marek committed
656
657
  end subroutine hh_trafo_kernel_4_bgq
end module real_bgq_kernel
658
! --------------------------------------------------------------------------------------------------