elpa2_kernels_real_simple.F90 6.6 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
  private
64 65 66 67 68 69 70
  public double_hh_trafo_generic_simple_double

#ifdef WANT_SINGLE_PRECISION_REAL
 public double_hh_trafo_generic_simple_single
#endif

  contains
71 72
  ! the intel compiler creates an temp array copy of array q
  ! This should be prevented if possible without using assumed size arrays
73
  subroutine double_hh_trafo_generic_simple_double(q, hh, nb, nq, ldq, ldh)
Andreas Marek's avatar
Andreas Marek committed
74
    use precision
75 76 77
#ifdef HAVE_DETAILED_TIMINGS
    use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
78
    implicit none
79

Andreas Marek's avatar
Andreas Marek committed
80
    integer(kind=ik), intent(in) :: nb, nq, ldq, ldh
81
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
82 83
    real(kind=rk8), intent(inout) :: q(ldq,*)
    real(kind=rk8), intent(in)    :: hh(ldh,*)
84
#else
85 86
    real(kind=rk8), intent(inout) :: q(ldq,1:nb+1)
    real(kind=rk8), intent(in)    :: hh(ldh,2)
87
#endif
88

89
    real(kind=rk8)                :: s, h1, h2, tau1, tau2, x(nq), y(nq)
Andreas Marek's avatar
Andreas Marek committed
90
    integer(kind=ik)             :: i
91

92
#ifdef HAVE_DETAILED_TIMINGS
93
    call timer%start("kernel generic simple: double_hh_trafo_generic_simple_double")
94
#endif
Andreas Marek's avatar
Andreas Marek committed
95
    ! Calculate dot product of the two Householder vectors
96

Andreas Marek's avatar
Andreas Marek committed
97 98 99 100
    s = hh(2,2)*1
    do i=3,nb
       s = s+hh(i,2)*hh(i-1,1)
    enddo
101

Andreas Marek's avatar
Andreas Marek committed
102
    ! Do the Householder transformations
103

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

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

Andreas Marek's avatar
Andreas Marek committed
108 109 110 111 112 113
    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
114

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

Andreas Marek's avatar
Andreas Marek committed
117 118
    tau1 = hh(1,1)
    tau2 = hh(1,2)
119

Andreas Marek's avatar
Andreas Marek committed
120 121 122 123 124
    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
125

Andreas Marek's avatar
Andreas Marek committed
126 127
    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)
128

Andreas Marek's avatar
Andreas Marek committed
129 130 131 132 133
    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
134

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

137
#ifdef HAVE_DETAILED_TIMINGS
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
    call timer%stop("kernel generic simple: double_hh_trafo_generic_simple_double")
#endif

  end subroutine double_hh_trafo_generic_simple_double

#ifdef WANT_SINGLE_PRECISION_REAL
  ! single precision implementation at the moment duplicated !!!

  subroutine double_hh_trafo_generic_simple_single(q, hh, nb, nq, ldq, ldh)
    use precision
#ifdef HAVE_DETAILED_TIMINGS
    use timings
#endif
    implicit none

    integer(kind=ik), intent(in) :: nb, nq, ldq, ldh
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
    real(kind=rk4), intent(inout) :: q(ldq,*)
    real(kind=rk4), intent(in)    :: hh(ldh,*)
#else
    real(kind=rk4), intent(inout) :: q(ldq,1:nb+1)
    real(kind=rk4), intent(in)    :: hh(ldh,2)
#endif

    real(kind=rk4)                :: s, h1, h2, tau1, tau2, x(nq), y(nq)
    integer(kind=ik)             :: i

#ifdef HAVE_DETAILED_TIMINGS
    call timer%start("kernel generic simple: double_hh_trafo_generic_simple_single")
167
#endif
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
    ! 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 the Householder transformations

    x(1:nq) = q(1:nq,2)

    y(1:nq) = q(1:nq,1) + q(1:nq,2)*hh(2,2)

    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

    x(1:nq) = x(1:nq) + q(1:nq,nb+1)*hh(nb,1)

    tau1 = hh(1,1)
    tau2 = hh(1,2)

    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

    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)

    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

    q(1:nq,nb+1) = q(1:nq,nb+1) + x(1:nq)*hh(nb,1)

#ifdef HAVE_DETAILED_TIMINGS
    call timer%stop("kernel generic simple: double_hh_trafo_generic_simple_single")
#endif

  end subroutine double_hh_trafo_generic_simple_single

#endif /* WANT_SINGLE_PRECISION_REAL */
217

Andreas Marek's avatar
Andreas Marek committed
218
end module real_generic_simple_kernel
219
! --------------------------------------------------------------------------------------------------