elpa2_kernels_real_simple.F90 4.38 KB
Newer Older
1 2
!    This file is part of ELPA.
!
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
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
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 48 49 50
!    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.
!
! This is the small and simple version (no hand unrolling of loops etc.) but for some
! compilers this performs better than a sophisticated version with transformed and unrolled loops.
!
! It should be compiled with the highest possible optimization level.
51
!
52 53 54 55 56 57
! 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".
!
! --------------------------------------------------------------------------------------------------
58 59 60

#include "config-f90.h"

Andreas Marek's avatar
Andreas Marek committed
61
module real_generic_simple_kernel
62

Andreas Marek's avatar
Andreas Marek committed
63 64 65 66
  private
  public double_hh_trafo_generic_simple
contains
  subroutine double_hh_trafo_generic_simple(q, hh, nb, nq, ldq, ldh)
Andreas Marek's avatar
Andreas Marek committed
67
    use precision
68 69 70
#ifdef HAVE_DETAILED_TIMINGS
    use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
71
    implicit none
72

Andreas Marek's avatar
Andreas Marek committed
73 74 75
    integer(kind=ik), intent(in) :: nb, nq, ldq, ldh
    real(kind=rk), intent(inout) :: q(ldq,*)
    real(kind=rk), intent(in)    :: hh(ldh,*)
76

Andreas Marek's avatar
Andreas Marek committed
77 78
    real(kind=rk)                :: s, h1, h2, tau1, tau2, x(nq), y(nq)
    integer(kind=ik)             :: i
79

80 81 82
#ifdef HAVE_DETAILED_TIMINGS
    call timer%start("kernel generic simple: double_hh_trafo_generic_simple")
#endif
Andreas Marek's avatar
Andreas Marek committed
83
    ! Calculate dot product of the two Householder vectors
84

Andreas Marek's avatar
Andreas Marek committed
85 86 87 88
    s = hh(2,2)*1
    do i=3,nb
       s = s+hh(i,2)*hh(i-1,1)
    enddo
89

Andreas Marek's avatar
Andreas Marek committed
90
    ! Do the Householder transformations
91

Andreas Marek's avatar
Andreas Marek committed
92
    x(1:nq) = q(1:nq,2)
93

Andreas Marek's avatar
Andreas Marek committed
94
    y(1:nq) = q(1:nq,1) + q(1:nq,2)*hh(2,2)
95

Andreas Marek's avatar
Andreas Marek committed
96 97 98 99 100 101
    do i=3,nb
       h1 = hh(i-1,1)
       h2 = hh(i,2)
       x(1:nq) = x(1:nq) + q(1:nq,i)*h1
       y(1:nq) = y(1:nq) + q(1:nq,i)*h2
    enddo
102

Andreas Marek's avatar
Andreas Marek committed
103
    x(1:nq) = x(1:nq) + q(1:nq,nb+1)*hh(nb,1)
104

Andreas Marek's avatar
Andreas Marek committed
105 106
    tau1 = hh(1,1)
    tau2 = hh(1,2)
107

Andreas Marek's avatar
Andreas Marek committed
108 109 110 111 112
    h1 = -tau1
    x(1:nq) = x(1:nq)*h1
    h1 = -tau2
    h2 = -tau2*s
    y(1:nq) = y(1:nq)*h1 + x(1:nq)*h2
113

Andreas Marek's avatar
Andreas Marek committed
114 115
    q(1:nq,1) = q(1:nq,1) + y(1:nq)
    q(1:nq,2) = q(1:nq,2) + x(1:nq) + y(1:nq)*hh(2,2)
116

Andreas Marek's avatar
Andreas Marek committed
117 118 119 120 121
    do i=3,nb
       h1 = hh(i-1,1)
       h2 = hh(i,2)
       q(1:nq,i) = q(1:nq,i) + x(1:nq)*h1 + y(1:nq)*h2
    enddo
122

Andreas Marek's avatar
Andreas Marek committed
123
    q(1:nq,nb+1) = q(1:nq,nb+1) + x(1:nq)*hh(nb,1)
124

125 126 127 128
#ifdef HAVE_DETAILED_TIMINGS
    call timer%stop("kernel generic simple: double_hh_trafo_generic_simple")
#endif

Andreas Marek's avatar
Andreas Marek committed
129 130
  end subroutine double_hh_trafo_generic_simple
end module real_generic_simple_kernel
131
! --------------------------------------------------------------------------------------------------