elpa1_template.X90 9.1 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
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
51
52
53
54
#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
!
!    This particular source code file contains additions, changes and
!    enhancements authored by Intel Corporation which is not part of
!    the ELPA consortium.
!
!    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.
!
!
! ELPA1 -- Faster replacements for ScaLAPACK symmetric eigenvalue routines
!
! 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

55
#include "../general/sanity.X90"
56
57
58
59
60

function elpa_solve_evp_&
         &MATH_DATATYPE&
	 &_1stage_&
	 &PRECISION&
61
	 &_impl (obj, a, ev, q) result(success)
62
63
64
65
   use precision
   use cuda_functions
   use mod_check_for_gpu
   use iso_c_binding
66
   use elpa_abstract_impl
67
68
69
70
   use elpa_mpi
   use elpa1_compute
   implicit none

71
   class(elpa_abstract_impl_t), intent(inout) :: obj
72
   real(kind=REAL_DATATYPE), intent(out)           :: ev(obj%na)
73
74
#if REALCASE == 1
#ifdef USE_ASSUMED_SIZE
75
76
   real(kind=C_DATATYPE_KIND), intent(inout)       :: a(obj%local_nrows,*)
   real(kind=C_DATATYPE_KIND), intent(out)         :: q(obj%local_nrows,*)
77
#else
Andreas Marek's avatar
Andreas Marek committed
78
79
   real(kind=C_DATATYPE_KIND), intent(inout)       :: a(obj%local_nrows,obj%local_ncols)
   real(kind=C_DATATYPE_KIND), intent(out)         :: q(obj%local_nrows,obj%local_ncols)
80
81
82
83
84
85
#endif
   real(kind=C_DATATYPE_KIND), allocatable         :: tau(:)
#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef USE_ASSUMED_SIZE
86
87
   complex(kind=C_DATATYPE_KIND), intent(inout)    :: a(obj%local_nrows,*)
   complex(kind=C_DATATYPE_KIND), intent(out)      :: q(obj%local_nrows,*)
88
#else
Andreas Marek's avatar
Andreas Marek committed
89
90
   complex(kind=C_DATATYPE_KIND), intent(inout)    :: a(obj%local_nrows,obj%local_ncols)
   complex(kind=C_DATATYPE_KIND), intent(out)      :: q(obj%local_nrows,obj%local_ncols)
91
92
93
94
95
96
97
#endif

   real(kind=REAL_DATATYPE), allocatable           :: q_real(:,:)
   complex(kind=C_DATATYPE_KIND), allocatable      :: tau(:)
   integer(kind=c_int)                             :: l_cols, l_rows, l_cols_nev, np_rows, np_cols
#endif /* COMPLEXCASE */

98
   logical                                         :: useGPU
99
100
101
102
103
104
105
106
107
108
   logical                                         :: success

   logical                                         :: do_useGPU
   integer(kind=ik)                                :: numberOfGPUDevices

   integer(kind=c_int)                             :: my_pe, n_pes, my_prow, my_pcol, mpierr
   real(kind=C_DATATYPE_KIND), allocatable         :: e(:)
   logical                                         :: wantDebug
   integer(kind=c_int)                             :: istat
   character(200)                                  :: errorMessage
109
110
   integer(kind=ik)                                :: na, nev, lda, ldq, nblk, matrixCols, &
                                                      mpi_comm_rows, mpi_comm_cols, mpi_comm_all
111

112
   call obj%timer%start("elpa_solve_evp_&
113
114
115
   &MATH_DATATYPE&
   &_1stage_&
   &PRECISION&
116
   &")
117

118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
   na         = obj%na
   nev        = obj%nev
   lda        = obj%local_nrows
   ldq        = obj%local_nrows
   nblk       = obj%nblk
   matrixCols = obj%local_ncols

   mpi_comm_rows = obj%get("mpi_comm_rows")
   mpi_comm_cols = obj%get("mpi_comm_cols")
   mpi_comm_all  = obj%get("mpi_comm_parent")

   if (obj%get("gpu") .eq. 1) then
     useGPU =.true.
   else
     useGPU = .false.
   endif

135
   call obj%timer%start("mpi_communication")
136
137
138
139
140
141
142
143
144
145
146
147

   call mpi_comm_rank(mpi_comm_all,my_pe,mpierr)
   call mpi_comm_size(mpi_comm_all,n_pes,mpierr)

   call mpi_comm_rank(mpi_comm_rows,my_prow,mpierr)
   call mpi_comm_rank(mpi_comm_cols,my_pcol,mpierr)

#if COMPLEXCASE == 1
   call mpi_comm_size(mpi_comm_rows,np_rows,mpierr)
   call mpi_comm_size(mpi_comm_cols,np_cols,mpierr)
#endif

148
   call obj%timer%stop("mpi_communication")
149
150
   success = .true.

151
152
   wantDebug = obj%get("debug") == 1
   do_useGPU = .false.
153

154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171

   if (useGPU) then
     if (check_for_gpu(my_pe,numberOfGPUDevices, wantDebug=wantDebug)) then

       do_useGPU = .true.
       ! set the neccessary parameters
       cudaMemcpyHostToDevice   = cuda_memcpyHostToDevice()
       cudaMemcpyDeviceToHost   = cuda_memcpyDeviceToHost()
       cudaMemcpyDeviceToDevice = cuda_memcpyDeviceToDevice()
       cudaHostRegisterPortable = cuda_hostRegisterPortable()
       cudaHostRegisterMapped   = cuda_hostRegisterMapped()
     else
       print *,"GPUs are requested but not detected! Aborting..."
       success = .false.
       return
     endif
   else
     ! check whether set by environment variable
172
     do_useGPU = obj%get("gpu") == 1
173
174

     if (do_useGPU) then
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
       if (check_for_gpu(my_pe,numberOfGPUDevices, wantDebug=wantDebug)) then

         ! set the neccessary parameters
         cudaMemcpyHostToDevice   = cuda_memcpyHostToDevice()
         cudaMemcpyDeviceToHost   = cuda_memcpyDeviceToHost()
         cudaMemcpyDeviceToDevice = cuda_memcpyDeviceToDevice()
         cudaHostRegisterPortable = cuda_hostRegisterPortable()
         cudaHostRegisterMapped   = cuda_hostRegisterMapped()
       else
         print *,"GPUs are requested but not detected! Aborting..."
         success = .false.
         return
       endif
     endif
   endif
#if COMPLEXCASE == 1
   l_rows = local_index(na, my_prow, np_rows, nblk, -1) ! Local rows of a and q
   l_cols = local_index(na, my_pcol, np_cols, nblk, -1) ! Local columns of q

   l_cols_nev = local_index(nev, my_pcol, np_cols, nblk, -1) ! Local columns corresponding to nev

   allocate(q_real(l_rows,l_cols), stat=istat, errmsg=errorMessage)
   if (istat .ne. 0) then
     print *,"solve_evp_&
     &MATH_DATATYPE&
     &_1stage_&
     &PRECISION&
     &" // ": error when allocating q_real "//errorMessage
     stop 1
   endif
#endif
   allocate(e(na), tau(na), stat=istat, errmsg=errorMessage)
   if (istat .ne. 0) then
     print *,"solve_evp_&
     &MATH_DATATYPE&
     &_1stage_&
     &PRECISION&
     &" // ": error when allocating e, tau "//errorMessage
     stop 1
   endif

216
   call obj%timer%start("forward")
217
218
219
220
   call tridiag_&
   &MATH_DATATYPE&
   &_&
   &PRECISION&
221
   & (obj, na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau, do_useGPU)
222
   call obj%timer%stop("forward")
223

224
   call obj%timer%start("solve")
225
226
   call solve_tridi_&
   &PRECISION&
227
   & (obj, na, nev, ev, e,  &
228
229
230
231
232
233
234
#if REALCASE == 1
      q, ldq,          &
#endif
#if COMPLEXCASE == 1
      q_real, l_rows,  &
#endif
      nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, wantDebug, success)
235
   call obj%timer%stop("solve")
236
237
238
239
240
   if (.not.(success)) return

#if COMPLEXCASE == 1
   q(1:l_rows,1:l_cols_nev) = q_real(1:l_rows,1:l_cols_nev)
#endif
241
242

   call obj%timer%start("back")
243
244
245
246
   call trans_ev_&
   &MATH_DATATYPE&
   &_&
   &PRECISION&
247
   & (obj, na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, do_useGPU)
248
   call obj%timer%stop("back")
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
#if COMPLEXCASE == 1
   deallocate(q_real, stat=istat, errmsg=errorMessage)
   if (istat .ne. 0) then
     print *,"solve_evp_&
     &MATH_DATATYPE&
     &_1stage_&
     &PRECISION&
     &" // ": error when deallocating q_real "//errorMessage
     stop 1
   endif
#endif
   deallocate(e, tau, stat=istat, errmsg=errorMessage)
   if (istat .ne. 0) then
     print *,"solve_evp_&
     &MATH_DATATYPE&
     &_1stage_&
     &PRECISION&
     &" // ": error when deallocating e, tau "//errorMessage
     stop 1
   endif

270
   call obj%timer%stop("elpa_solve_evp_&
271
272
273
   &MATH_DATATYPE&
   &_1stage_&
   &PRECISION&
274
   &")
275
276
277
end function