elpa2_kernels_real_simple.F90 4.52 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
    integer(kind=ik), intent(in) :: nb, nq, ldq, ldh
74
#ifdef DESPERATELY_WANT_ASSUMED_SIZE
Andreas Marek's avatar
Andreas Marek committed
75 76
    real(kind=rk), intent(inout) :: q(ldq,*)
    real(kind=rk), intent(in)    :: hh(ldh,*)
77 78 79 80
#else
    real(kind=rk), intent(inout) :: q(ldq,1:nb+1)
    real(kind=rk), intent(in)    :: hh(ldh,2)
#endif
81

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

85 86 87
#ifdef HAVE_DETAILED_TIMINGS
    call timer%start("kernel generic simple: double_hh_trafo_generic_simple")
#endif
Andreas Marek's avatar
Andreas Marek committed
88
    ! Calculate dot product of the two Householder vectors
89

Andreas Marek's avatar
Andreas Marek committed
90 91 92 93
    s = hh(2,2)*1
    do i=3,nb
       s = s+hh(i,2)*hh(i-1,1)
    enddo
94

Andreas Marek's avatar
Andreas Marek committed
95
    ! Do the Householder transformations
96

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

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

Andreas Marek's avatar
Andreas Marek committed
101 102 103 104 105 106
    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
107

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

Andreas Marek's avatar
Andreas Marek committed
110 111
    tau1 = hh(1,1)
    tau2 = hh(1,2)
112

Andreas Marek's avatar
Andreas Marek committed
113 114 115 116 117
    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
118

Andreas Marek's avatar
Andreas Marek committed
119 120
    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)
121

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

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

130 131 132 133
#ifdef HAVE_DETAILED_TIMINGS
    call timer%stop("kernel generic simple: double_hh_trafo_generic_simple")
#endif

Andreas Marek's avatar
Andreas Marek committed
134 135
  end subroutine double_hh_trafo_generic_simple
end module real_generic_simple_kernel
136
! --------------------------------------------------------------------------------------------------