elpa_invert_trm.X90 9.51 KB
Newer Older
Andreas Marek's avatar
Andreas Marek committed
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
!    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".

53
#include "../sanity.X90"
Andreas Marek's avatar
Andreas Marek committed
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112

       use precision
       use elpa1_compute
       use elpa_utilities
       use elpa_mpi
#ifdef HAVE_DETAILED_TIMINGS
       use timings
#else
       use timings_dummy
#endif
       implicit none

       integer(kind=ik)             :: na, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
#if  REALCASE ==  1
#ifdef USE_ASSUMED_SIZE
       real(kind=REAL_DATATYPE)     :: a(lda,*)
#else
       real(kind=REAL_DATATYPE)     :: a(lda,matrixCols)
#endif
#endif

#if COMPLEXCASE == 1
#ifdef USE_ASSUMED_SIZE
       complex(kind=COMPLEX_DATATYPE)   :: a(lda,*)
#else
       complex(kind=COMPLEX_DATATYPE)   :: a(lda,matrixCols)
#endif
#endif

       integer(kind=ik)             :: my_prow, my_pcol, np_rows, np_cols, mpierr
       integer(kind=ik)             :: l_cols, l_rows, l_col1, l_row1, l_colx, l_rowx
       integer(kind=ik)             :: n, nc, i, info, ns, nb
#if REALCASE ==  1
       real(kind=REAL_DATATYPE), allocatable   :: tmp1(:), tmp2(:,:), tmat1(:,:), tmat2(:,:)
#endif
#if COMPLEXCASE == 1
       complex(kind=COMPLEX_DATATYPE), allocatable    :: tmp1(:), tmp2(:,:), tmat1(:,:), tmat2(:,:)
#endif
       logical, intent(in)          :: wantDebug
       logical                      :: success
       integer(kind=ik)             :: istat
       character(200)               :: errorMessage

       call timer%start("mpi_communication")
       call mpi_comm_rank(mpi_comm_rows,my_prow,mpierr)
       call mpi_comm_size(mpi_comm_rows,np_rows,mpierr)
       call mpi_comm_rank(mpi_comm_cols,my_pcol,mpierr)
       call mpi_comm_size(mpi_comm_cols,np_cols,mpierr)
       call timer%stop("mpi_communication")
       success = .true.

       l_rows = local_index(na, my_prow, np_rows, nblk, -1) ! Local rows of a
       l_cols = local_index(na, my_pcol, np_cols, nblk, -1) ! Local cols of a

       allocate(tmp1(nblk*nblk), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"elpa_invert_trm_&
	 &MATH_DATATYPE&
	 &: error when allocating tmp1 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
113
         stop 1
Andreas Marek's avatar
Andreas Marek committed
114
115
116
117
118
119
120
       endif

       allocate(tmp2(nblk,nblk), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"elpa_invert_trm_&
	 &MATH_DATATYPE&
	 &: error when allocating tmp2 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
121
         stop 1
Andreas Marek's avatar
Andreas Marek committed
122
123
124
125
126
127
128
129
130
131
       endif

       tmp1 = 0
       tmp2 = 0

       allocate(tmat1(l_rows,nblk), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"elpa_invert_trm_&
	 &MATH_DATATYPE&
	 &: error when allocating tmat1 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
132
         stop 1
Andreas Marek's avatar
Andreas Marek committed
133
134
135
136
137
138
139
       endif

       allocate(tmat2(nblk,l_cols), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"elpa_invert_trm_&
	 &MATH_DATATYPE&
	 &: error when allocating tmat2 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
140
         stop 1
Andreas Marek's avatar
Andreas Marek committed
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
       endif

       tmat1 = 0
       tmat2 = 0


       ns = ((na-1)/nblk)*nblk + 1

       do n = ns,1,-nblk

         l_row1 = local_index(n, my_prow, np_rows, nblk, +1)
         l_col1 = local_index(n, my_pcol, np_cols, nblk, +1)

         nb = nblk
         if (na-n+1 < nblk) nb = na-n+1

         l_rowx = local_index(n+nb, my_prow, np_rows, nblk, +1)
         l_colx = local_index(n+nb, my_pcol, np_cols, nblk, +1)

         if (my_prow==prow(n, nblk, np_rows)) then

           if (my_pcol==pcol(n, nblk, np_cols)) then
             call timer%start("blas")
164
165
166

             call PRECISION_TRTRI('U', 'N', nb, a(l_row1,l_col1), lda, info)

Andreas Marek's avatar
Andreas Marek committed
167
168
169
170
171
             call timer%stop("blas")

             if (info/=0) then
               if (wantDebug) write(error_unit,*) "elpa_invert_trm_&
	       &MATH_DATATYPE&
172
173

#if REALCASE == 1
Andreas Marek's avatar
Andreas Marek committed
174
	       &: Error in DTRTRI"
175
176
177
178
179
#endif
#if COMPLEXCASE == 1
	       &: Error in ZTRTRI"
#endif

Andreas Marek's avatar
Andreas Marek committed
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
               success = .false.
               return
             endif

             nc = 0
             do i=1,nb
               tmp1(nc+1:nc+i) = a(l_row1:l_row1+i-1,l_col1+i-1)
               nc = nc+i
             enddo
           endif
#ifdef WITH_MPI
           call timer%start("mpi_communication")
           call MPI_Bcast(tmp1, nb*(nb+1)/2,        &
#if REALCASE == 1
                          MPI_REAL_PRECISION,       &
#endif
#if COMPLEXCASE == 1
                          MPI_COMPLEX_PRECISION,       &
#endif
                          pcol(n, nblk, np_cols), mpi_comm_cols, mpierr)
           call timer%stop("mpi_communication")
#endif /* WITH_MPI */
           nc = 0
           do i=1,nb
             tmp2(1:i,i) = tmp1(nc+1:nc+i)
             nc = nc+i
           enddo

           call timer%start("blas")
           if (l_cols-l_colx+1>0) &
210
	      call PRECISION_TRMM ('L', 'U', 'N', 'N', nb, l_cols-l_colx+1,  &
Andreas Marek's avatar
Andreas Marek committed
211
#if REALCASE == 1
212
                                   CONST_1_0,   &
Andreas Marek's avatar
Andreas Marek committed
213
214
#endif
#if COMPLEXCASE == 1
215
                                   CONST_COMPLEX_PAIR_1_0, &
Andreas Marek's avatar
Andreas Marek committed
216
#endif
217
                                   tmp2, ubound(tmp2,dim=1), a(l_row1,l_colx), lda)
Andreas Marek's avatar
Andreas Marek committed
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
           call timer%stop("blas")
           if (l_colx<=l_cols)   tmat2(1:nb,l_colx:l_cols) = a(l_row1:l_row1+nb-1,l_colx:l_cols)
           if (my_pcol==pcol(n, nblk, np_cols)) tmat2(1:nb,l_col1:l_col1+nb-1) = tmp2(1:nb,1:nb) ! tmp2 has the lower left triangle 0

         endif

         if (l_row1>1) then
           if (my_pcol==pcol(n, nblk, np_cols)) then
             tmat1(1:l_row1-1,1:nb) = a(1:l_row1-1,l_col1:l_col1+nb-1)
             a(1:l_row1-1,l_col1:l_col1+nb-1) = 0
           endif

           do i=1,nb
#ifdef WITH_MPI
             call timer%start("mpi_communication")
             call MPI_Bcast(tmat1(1,i), l_row1-1,       &
#if REALCASE == 1
                             MPI_REAL_PRECISION,       &
#endif
#if COMPLEXCASE == 1
                             MPI_COMPLEX_PRECISION,       &
#endif
                             pcol(n, nblk, np_cols), mpi_comm_cols, mpierr)

             call timer%stop("mpi_communication")
#endif /* WITH_MPI */
           enddo
         endif
#ifdef WITH_MPI
         call timer%start("mpi_communication")
         if (l_cols-l_col1+1>0) &
	    call MPI_Bcast(tmat2(1,l_col1), (l_cols-l_col1+1)*nblk,      &
#if REALCASE == 1
                             MPI_REAL_PRECISION,       &
#endif
#if COMPLEXCASE == 1
                             MPI_COMPLEX_PRECISION,       &
#endif
                             prow(n, nblk, np_rows), mpi_comm_rows, mpierr)

          call timer%stop("mpi_communication")
#endif /* WITH_MPI */

         call timer%start("blas")
         if (l_row1>1 .and. l_cols-l_col1+1>0) &
263
264
265
266
267
268
269
270
           call PRECISION_GEMM('N', 'N', l_row1-1, l_cols-l_col1+1, nb,    &
#if REALCASE == 1
                               -CONST_1_0,   &
#endif
#if COMPLEXCASE == 1
                               -CONST_COMPLEX_PAIR_1_0,        &
#endif
                                tmat1, ubound(tmat1,dim=1), tmat2(1,l_col1), ubound(tmat2,dim=1), &
Andreas Marek's avatar
Andreas Marek committed
271
#if REALCASE == 1
272
                                CONST_1_0,   &
Andreas Marek's avatar
Andreas Marek committed
273
274
#endif
#if COMPLEXCASE == 1
275
                                CONST_COMPLEX_PAIR_1_0,        &
Andreas Marek's avatar
Andreas Marek committed
276
#endif
277
278
                                 a(1,l_col1), lda)

Andreas Marek's avatar
Andreas Marek committed
279
280
281
282
283
284
285
286
287
         call timer%stop("blas")

       enddo

       deallocate(tmp1, tmp2, tmat1, tmat2, stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"elpa_invert_trm_&
	 &MATH_DATATYPE&
	 &: error when deallocating tmp1 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
288
         stop 1
Andreas Marek's avatar
Andreas Marek committed
289
290
291
292
293
294
       endif

#undef REALCASE
#undef COMPLEXCASE
#undef DOUBLE_PRECISION
#undef SINGLE_PRECISION