read_real.F90 13.7 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 Naturwissenschaften,
!      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
!    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.
!
!
#include "config-f90.h"
44
45
46
47
48
49
50
51
52
53
54
!>
!> Fortran test programm to demonstrates the use of
!> ELPA 1 real case library.
!> This program can read a matrix from an ascii
!> file and computes then the Eigenvectors.
!> If "HAVE_REDIRECT" was defined at build time
!> the stdout and stderr output of each MPI task
!> can be redirected to files if the environment
!> variable "REDIRECT_ELPA_TEST_OUTPUT" is set
!> to "true".
!>
55
56
57
58
59
60
61
62
63
program read_real

!-------------------------------------------------------------------------------
! Standard eigenvalue problem - REAL version
!
! This program demonstrates the use of the ELPA module
! together with standard scalapack routines
!-------------------------------------------------------------------------------

64
   use precision
65
   use elpa1
66
   use elpa_utilities, only : error_unit
67
68
69
#ifdef WITH_OPENMP
   use test_util
#endif
70
71
72
#ifdef HAVE_REDIRECT
   use redirect
#endif
73
#ifdef HAVE_DETAILED_TIMINGS
74
   use timings
75
#endif
76
77
78
79
#ifdef HAVE_MPI_MODULE
   use mpi
   implicit none
#else
80
81
   implicit none
   include 'mpif.h'
82
#endif
83
84
85
86
87
88

   !-------------------------------------------------------------------------------
   ! Please set system size parameters below!
   ! nblk: Blocking factor in block cyclic distribution
   !-------------------------------------------------------------------------------

89
   integer(kind=ik), parameter :: nblk = 16
90
91
92
93

   !-------------------------------------------------------------------------------
   !  Local Variables

94
   integer(kind=ik)            :: na, nev
95

96
   integer(kind=ik)            :: np_rows, np_cols, na_rows, na_cols
97

98
99
   integer(kind=ik)            :: myid, nprocs, my_prow, my_pcol, mpi_comm_rows, mpi_comm_cols
   integer(kind=ik)            :: i, mpierr, my_blacs_ctxt, sc_desc(9), info, nprow, npcol, lenarg
100

101
   integer(kind=ik), external  :: numroc
102

103
104
   real(kind=rk)               :: err, errmax
   real(kind=rk), allocatable  :: a(:,:), z(:,:), tmp1(:,:), tmp2(:,:), as(:,:), ev(:)
105

106
   character*256               :: filename
107
#ifdef WITH_OPENMP
108
   integer(kind=iK)            :: omp_get_max_threads,  required_mpi_thread_level, provided_mpi_thread_level
109
110
111
112
113
114
115
116
117
118
119
120
121
#endif
   !-------------------------------------------------------------------------------
   !  MPI Initialization

#ifndef WITH_OPENMP
   call mpi_init(mpierr)
#else
   required_mpi_thread_level = MPI_THREAD_MULTIPLE

   call mpi_init_thread(required_mpi_thread_level,     &
                        provided_mpi_thread_level, mpierr)

   if (required_mpi_thread_level .ne. provided_mpi_thread_level) then
122
     write(error_unit,*) "MPI ERROR: MPI_THREAD_MULTIPLE is not provided on this system"
123
124
     write(error_unit,*) "           only ", mpi_thread_level_name(provided_mpi_thread_level), " is available"
     call exit(77)
125
126
127
128
129
130
   endif

#endif
   call mpi_comm_rank(mpi_comm_world,myid,mpierr)
   call mpi_comm_size(mpi_comm_world,nprocs,mpierr)

131
132
133
134
135
136
137
138
#ifdef HAVE_REDIRECT
   if (check_redirect_environment_variable()) then
     if (myid .eq. 0) then
       print *," "
       print *,"Redirection of mpi processes is used"
       print *," "
       if (create_directories() .ne. 1) then
         write(error_unit,*) "Unable to create directory for stdout and stderr!"
Pavel Kus's avatar
Pavel Kus committed
139
         stop 1
140
141
142
143
144
145
146
       endif
     endif
     call MPI_BARRIER(MPI_COMM_WORLD, mpierr)
     call redirect_stdout(myid)
   endif
#endif

147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
#ifdef HAVE_DETAILED_TIMINGS

   ! initialise the timing functionality

#ifdef HAVE_LIBPAPI
   call timer%measure_flops(.true.)
#endif

   call timer%measure_allocated_memory(.true.)
   call timer%measure_virtual_memory(.true.)
   call timer%measure_max_allocated_memory(.true.)

   call timer%set_print_options(&
#ifdef HAVE_LIBPAPI
                print_flop_count=.true., &
                print_flop_rate=.true., &
#endif
                print_allocated_memory = .true. , &
                print_virtual_memory=.true., &
                print_max_allocated_memory=.true.)


  call timer%enable()

  call timer%start("program")
#endif
173
174
175
176
177
178
179
180
181
   !-------------------------------------------------------------------------------
   ! Get the name of the input file containing the matrix and open input file
   ! Please note:
   ! get_command_argument is a FORTRAN 2003 intrinsic which may not be implemented
   ! for every Fortran compiler!!!

   if(myid==0) then
      call get_command_argument(1,filename,lenarg,info)
      if(info/=0) then
182
         write(error_unit,*) 'Usage: test_real matrix_file'
Andreas Marek's avatar
Andreas Marek committed
183
         call mpi_abort(mpi_comm_world,1,mpierr)
184
185
186
      endif
      open(10,file=filename,action='READ',status='OLD',iostat=info)
      if(info/=0) then
187
         write(error_unit,*) 'Error: Unable to open ',trim(filename)
Andreas Marek's avatar
Andreas Marek committed
188
         call mpi_abort(mpi_comm_world,1,mpierr)
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
217
218
219
220
221
222
223
224
225
226
227
228
      endif
   endif
   call mpi_barrier(mpi_comm_world, mpierr) ! Just for safety

   !-------------------------------------------------------------------------------
   ! Selection of number of processor rows/columns
   ! We try to set up the grid square-like, i.e. start the search for possible
   ! divisors of nprocs with a number next to the square root of nprocs
   ! and decrement it until a divisor is found.

   do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
      if(mod(nprocs,np_cols) == 0 ) exit
   enddo
   ! at the end of the above loop, nprocs is always divisible by np_cols

   np_rows = nprocs/np_cols

   if(myid==0) then
      print *
      print '(a)','Standard eigenvalue problem - REAL version'
      print *
      print '(3(a,i0))','Number of processor rows=',np_rows,', cols=',np_cols,', total=',nprocs
      print *
   endif

   !-------------------------------------------------------------------------------
   ! Set up BLACS context and MPI communicators
   !
   ! The BLACS context is only necessary for using Scalapack.
   !
   ! For ELPA, the MPI communicators along rows/cols are sufficient,
   ! and the grid setup may be done in an arbitrary way as long as it is
   ! consistent (i.e. 0<=my_prow<np_rows, 0<=my_pcol<np_cols and every
   ! process has a unique (my_prow,my_pcol) pair).

   my_blacs_ctxt = mpi_comm_world
   call BLACS_Gridinit( my_blacs_ctxt, 'C', np_rows, np_cols )
   call BLACS_Gridinfo( my_blacs_ctxt, nprow, npcol, my_prow, my_pcol )

   ! All ELPA routines need MPI communicators for communicating within
229
   ! rows or columns of processes, these are set in elpa_get_communicators
230

231
   call elpa_get_communicators(mpi_comm_world, my_prow, my_pcol, &
232
233
234
235
236
237
238
239
                               mpi_comm_rows, mpi_comm_cols)

   ! Read matrix size
   if(myid==0) read(10,*) na
   call mpi_bcast(na, 1, mpi_integer, 0, mpi_comm_world, mpierr)

   ! Quick check for plausibility
   if(na<=0 .or. na>10000000) then
240
      if(myid==0) write(error_unit,*) 'Illegal value for matrix size: ',na
241
      call mpi_finalize(mpierr)
Pavel Kus's avatar
Pavel Kus committed
242
      stop 1
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
   endif
   if(myid==0) print *,'Matrix size: ',na

   ! Determine the necessary size of the distributed matrices,
   ! we use the Scalapack tools routine NUMROC for that.

   na_rows = numroc(na, nblk, my_prow, 0, np_rows)
   na_cols = numroc(na, nblk, my_pcol, 0, np_cols)

   ! Set up a scalapack descriptor for the checks below.
   ! For ELPA the following restrictions hold:
   ! - block sizes in both directions must be identical (args 4+5)
   ! - first row and column of the distributed matrix must be on row/col 0/0 (args 6+7)

   call descinit( sc_desc, na, na, nblk, nblk, 0, 0, my_blacs_ctxt, na_rows, info )

   !-------------------------------------------------------------------------------
260
   ! Allocate matrices
261
262
263
#ifdef HAVE_DETAILED_TIMINGS
   call timer%start("set up matrix")
#endif
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
   allocate(a (na_rows,na_cols))
   allocate(z (na_rows,na_cols))
   allocate(as(na_rows,na_cols))

   allocate(ev(na))

   !-------------------------------------------------------------------------------
   ! Read matrix

   call read_matrix(10, na, a, ubound(a,1), nblk, my_prow, my_pcol, np_rows, np_cols)
   if(myid==0) close(10)

   nev = na ! all eigenvaules

   ! Save original matrix A for later accuracy checks

   as = a
281
282
283
#ifdef HAVE_DETAILED_TIMINGS
   call timer%stop("set up matrix")
#endif
284
285
286
287
   !-------------------------------------------------------------------------------
   ! Calculate eigenvalues/eigenvectors

   call mpi_barrier(mpi_comm_world, mpierr) ! for correct timings only
288
289
   call elpa_solve_evp_real_1stage(na, nev, a, na_rows, ev, z, na_rows, nblk, &
                                   mpi_comm_rows, mpi_comm_cols)
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355

   if(myid == 0) print *,'Time tridiag_real :',time_evp_fwd
   if(myid == 0) print *,'Time solve_tridi  :',time_evp_solve
   if(myid == 0) print *,'Time trans_ev_real:',time_evp_back

   if(myid == 0) then
      do i=1,nev
         print '(i6,g25.15)',i,ev(i)
      enddo
   endif

   !-------------------------------------------------------------------------------
   ! Test correctness of result (using plain scalapack routines)

   deallocate(a)
   allocate(tmp1(na_rows,na_cols))

   ! 1. Residual (maximum of || A*Zi - Zi*EVi ||)

   ! tmp1 =  A * Z
   call pdgemm('N','N',na,nev,na,1.d0,as,1,1,sc_desc, &
           z,1,1,sc_desc,0.d0,tmp1,1,1,sc_desc)

   deallocate(as)
   allocate(tmp2(na_rows,na_cols))

   ! tmp2 = Zi*EVi
   tmp2(:,:) = z(:,:)
   do i=1,nev
      call pdscal(na,ev(i),tmp2,1,i,sc_desc,1)
   enddo

   !  tmp1 = A*Zi - Zi*EVi
   tmp1(:,:) =  tmp1(:,:) - tmp2(:,:)

   ! Get maximum norm of columns of tmp1
   errmax = 0
   do i=1,nev
      err = 0
      call pdnrm2(na,err,tmp1,1,i,sc_desc,1)
      errmax = max(errmax, err)
   enddo

   ! Get maximum error norm over all processors
   err = errmax
   call mpi_allreduce(err,errmax,1,MPI_REAL8,MPI_MAX,MPI_COMM_WORLD,mpierr)
   if(myid==0) print *
   if(myid==0) print *,'Error Residual     :',errmax

   ! 2. Eigenvector orthogonality

   ! tmp1 = Z**T * Z
   tmp1 = 0
   call pdgemm('T','N',nev,nev,na,1.d0,z,1,1,sc_desc, &
           z,1,1,sc_desc,0.d0,tmp1,1,1,sc_desc)
   ! Initialize tmp2 to unit matrix
   tmp2 = 0
   call pdlaset('A',nev,nev,0.d0,1.d0,tmp2,1,1,sc_desc)

   ! tmp1 = Z**T * Z - Unit Matrix
   tmp1(:,:) =  tmp1(:,:) - tmp2(:,:)

   ! Get maximum error (max abs value in tmp1)
   err = maxval(abs(tmp1))
   call mpi_allreduce(err,errmax,1,MPI_REAL8,MPI_MAX,MPI_COMM_WORLD,mpierr)
   if(myid==0) print *,'Error Orthogonality:',errmax
356

357
358
359
360
   deallocate(z)
   deallocate(tmp1)
   deallocate(tmp2)
   deallocate(ev)
361
362
363
364
365
366
367
368
#ifdef HAVE_DETAILED_TIMINGS
   call timer%stop("program")
   print *," "
   print *,"Timings program:"
   call timer%print("program")
   print *," "
   print *,"End timings program"
#endif
369
370
371
372
373
374
375
   call blacs_gridexit(my_blacs_ctxt)
   call mpi_finalize(mpierr)

end

!-------------------------------------------------------------------------------
subroutine read_matrix(iunit, na, a, lda, nblk, my_prow, my_pcol, np_rows, np_cols)
376
   use precision
377
   implicit none
378
379
380
381
382
#ifdef HAVE_MPI_MODULE
   use mpi
   implicit none
#else
   implicit none
383
   include 'mpif.h'
384
#endif
385

386
387
   integer(kind=ik), intent(in)  :: iunit, na, lda, nblk, my_prow, my_pcol, np_rows, np_cols
   real(kind=rk), intent(out)    :: a(lda, *)
388

389
390
   integer(kind=ik)              :: i, j, lr, lc, myid, mpierr
   integer(kind=ik), allocatable :: l_row(:), l_col(:)
391

392
   real(kind=rk), allocatable    :: col(:)
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443

   ! allocate and set index arrays

   allocate(l_row(na))
   allocate(l_col(na))

   ! Mapping of global rows/cols to local

   l_row(:) = 0
   l_col(:) = 0

   lr = 0 ! local row counter
   lc = 0 ! local column counter

   do i = 1, na

     if( MOD((i-1)/nblk,np_rows) == my_prow) then
       ! row i is on local processor
       lr = lr+1
       l_row(i) = lr
     endif

     if( MOD((i-1)/nblk,np_cols) == my_pcol) then
       ! column i is on local processor
       lc = lc+1
       l_col(i) = lc
     endif

   enddo

   call mpi_comm_rank(mpi_comm_world,myid,mpierr)
   allocate(col(na))

   do i=1,na
      if(myid==0) read(iunit,*) col(1:i)
      call mpi_bcast(col,i,MPI_REAL8,0,MPI_COMM_WORLD,mpierr)
      if(l_col(i) > 0) then
         do j=1,i
            if(l_row(j)>0) a(l_row(j),l_col(i)) = col(j)
         enddo
      endif
      if(l_row(i) > 0) then
         do j=1,i-1
            if(l_col(j)>0) a(l_row(i),l_col(j)) = col(j)
         enddo
      endif
   enddo

   deallocate(l_row, l_col, col)

end subroutine read_matrix