Skip to content
GitLab
Menu
Projects
Groups
Snippets
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
elpa
elpa
Commits
46c851e8
Commit
46c851e8
authored
Aug 30, 2016
by
Andreas Marek
Browse files
Unify real double single precision generic simple kernel
parent
55b6158f
Changes
3
Pipelines
1
Hide whitespace changes
Inline
Side-by-side
Makefile.am
View file @
46c851e8
...
...
@@ -52,6 +52,7 @@ EXTRA_libelpa@SUFFIX@_private_la_DEPENDENCIES = \
src/elpa2_compute_complex_template.X90
\
src/elpa2_kernels/elpa2_kernels_real_template.X90
\
src/elpa2_kernels/elpa2_kernels_complex_template.X90
\
src/elpa2_kernels/elpa2_kernels_real_simple_template.X90
\
src/redist_band.X90
lib_LTLIBRARIES
=
libelpa@SUFFIX@.la
...
...
@@ -773,6 +774,7 @@ EXTRA_DIST = \
src/elpa2_compute_complex_template.X90
\
src/elpa2_kernels/elpa2_kernels_real_template.X90
\
src/elpa2_kernels/elpa2_kernels_complex_template.X90
\
src/elpa2_kernels/elpa2_kernels_real_simple_template.X90
\
src/redist_band.X90
\
src/elpa_qr/elpa_qrkernels.X90
\
src/ev_tridi_band_gpu_c_v2_complex_template.Xcu
\
...
...
src/elpa2_kernels/elpa2_kernels_real_simple.F90
View file @
46c851e8
...
...
@@ -10,7 +10,7 @@
! - 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 Naturwissenschaft
e
n,
! - Max-Plack-Institut für Mathematik in den Naturwissenschaft
r
n,
! Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
! and
! - IBM Deutschland GmbH
...
...
@@ -68,152 +68,21 @@ module real_generic_simple_kernel
#endif
contains
! the intel compiler creates an temp array copy of array q
! This should be prevented if possible without using assumed size arrays
subroutine
double_hh_trafo_generic_simple_double
(
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
=
rk8
),
intent
(
inout
)
::
q
(
ldq
,
*
)
real
(
kind
=
rk8
),
intent
(
in
)
::
hh
(
ldh
,
*
)
#else
real
(
kind
=
rk8
),
intent
(
inout
)
::
q
(
ldq
,
1
:
nb
+1
)
real
(
kind
=
rk8
),
intent
(
in
)
::
hh
(
ldh
,
2
)
#endif
real
(
kind
=
rk8
)
::
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_double"
)
#endif
! 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
#define DOUBLE_PRECISION_REAL 1
#define REAL_DATATYPE rk8
#include "elpa2_kernels_real_simple_template.X90"
#undef DOUBLE_PRECISION_REAL
#undef REAL_DATATYPE
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_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
#undef DOUBLE_PRECISION_REAL
#define REAL_DATATYPE rk4
#include "elpa2_kernels_real_simple_template.X90"
#undef DOUBLE_PRECISION_REAL
#undef REAL_DATATYPE
#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"
)
#endif
! 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 */
end
module
real_generic_simple_kernel
! --------------------------------------------------------------------------------------------------
src/elpa2_kernels/elpa2_kernels_real_simple_template.X90
0 → 100644
View file @
46c851e8
#if 0
! 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 Naturwissenschaften,
! 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.
!
!
! --------------------------------------------------------------------------------------------------
!
! 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.
!
! 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".
!
! --------------------------------------------------------------------------------------------------
#endif
#ifdef DOUBLE_PRECISION_REAL
subroutine double_hh_trafo_generic_simple_double(q, hh, nb, nq, ldq, ldh)
#else
subroutine double_hh_trafo_generic_simple_single(q, hh, nb, nq, ldq, ldh)
#endif
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=REAL_DATATYPE), intent(inout) :: q(ldq,*)
real(kind=REAL_DATATYPE), intent(in) :: hh(ldh,*)
#else
real(kind=REAL_DATATYPE), intent(inout) :: q(ldq,1:nb+1)
real(kind=REAL_DATATYPE), intent(in) :: hh(ldh,2)
#endif
real(kind=REAL_DATATYPE) :: s, h1, h2, tau1, tau2, x(nq), y(nq)
integer(kind=ik) :: i
#ifdef DOUBLE_PRECISION_REAL
#ifdef HAVE_DETAILED_TIMINGS
call timer%start("kernel generic simple: double_hh_trafo_generic_simple_double")
#endif
#else
#ifdef HAVE_DETAILED_TIMINGS
call timer%start("kernel generic simple: double_hh_trafo_generic_simple_single")
#endif
#endif
! 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 DOUBLE_PRECISION_REAL
#ifdef HAVE_DETAILED_TIMINGS
call timer%stop("kernel generic simple: double_hh_trafo_generic_simple_double")
#endif
#else
#ifdef HAVE_DETAILED_TIMINGS
call timer%stop("kernel generic simple: double_hh_trafo_generic_simple_single")
#endif
#endif
#ifdef DOUBLE_PRECISION_REAL
end subroutine double_hh_trafo_generic_simple_double
#else
end subroutine double_hh_trafo_generic_simple_single
#endif
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment