1 2!------------------------------------------------------------------------------- 3 4subroutine solve_evp_complex(na, nev, a, lda, ev, q, ldq, nblk, mpi_comm_rows, mpi_comm_cols) 5 6!------------------------------------------------------------------------------- 7! solve_evp_complex: Solves the complex eigenvalue problem 8! 9! Parameters 10! 11! na Order of matrix a 12! 13! nev Number of eigenvalues needed 14! The smallest nev eigenvalues/eigenvectors are calculated. 15! 16! a(lda,*) Distributed matrix for which eigenvalues are to be computed. 17! Distribution is like in Scalapack. 18! The full matrix must be set (not only one half like in scalapack). 19! Destroyed on exit (upper and lower half). 20! 21! lda Leading dimension of a 22! 23! ev(na) On output: eigenvalues of a, every processor gets the complete set 24! 25! q(ldq,*) On output: Eigenvectors of a 26! Distribution is like in Scalapack. 27! Must be always dimensioned to the full size (corresponding to (na,na)) 28! even if only a part of the eigenvalues is needed. 29! 30! ldq Leading dimension of q 31! 32! nblk blocksize of cyclic distribution, must be the same in both directions! 33! 34! mpi_comm_rows 35! mpi_comm_cols 36! MPI-Communicators for rows/columns 37! 38!------------------------------------------------------------------------------- 39 40 use ELPA1 41 42 implicit none 43 44 include 'mpif.h' 45 46 integer, intent(in) :: na, nev, lda, ldq, nblk, mpi_comm_rows, mpi_comm_cols 47 complex*16 :: a(lda,*), q(ldq,*) 48 real*8 :: ev(na) 49 50 integer my_prow, my_pcol, np_rows, np_cols, mpierr 51 integer l_rows, l_cols, l_cols_nev 52 real*8, allocatable :: q_real(:,:), e(:) 53 complex*16, allocatable :: tau(:) 54 real*8 ttt0, ttt1 55 56 call mpi_comm_rank(mpi_comm_rows,my_prow,mpierr) 57 call mpi_comm_size(mpi_comm_rows,np_rows,mpierr) 58 call mpi_comm_rank(mpi_comm_cols,my_pcol,mpierr) 59 call mpi_comm_size(mpi_comm_cols,np_cols,mpierr) 60 61 l_rows = local_index(na, my_prow, np_rows, nblk, -1) ! Local rows of a and q 62 l_cols = local_index(na, my_pcol, np_cols, nblk, -1) ! Local columns of q 63 64 l_cols_nev = local_index(nev, my_pcol, np_cols, nblk, -1) ! Local columns corresponding to nev 65 66 allocate(e(na), tau(na)) 67 allocate(q_real(l_rows,l_cols)) 68 69 ttt0 = MPI_Wtime() 70 call tridiag_complex(na, a, lda, nblk, mpi_comm_rows, mpi_comm_cols, ev, e, tau) 71 ttt1 = MPI_Wtime() 72 if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) print *,'Time tridiag_complex :',ttt1-ttt0 73 time_evp_fwd = ttt1-ttt0 74 75 ttt0 = MPI_Wtime() 76 call solve_tridi(na, nev, ev, e, q_real, l_rows, nblk, mpi_comm_rows, mpi_comm_cols) 77 ttt1 = MPI_Wtime() 78 if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) print *,'Time solve_tridi :',ttt1-ttt0 79 time_evp_solve = ttt1-ttt0 80 81 ttt0 = MPI_Wtime() 82 q(1:l_rows,1:l_cols_nev) = q_real(1:l_rows,1:l_cols_nev) 83 84 call trans_ev_complex(na, nev, a, lda, tau, q, ldq, nblk, mpi_comm_rows, mpi_comm_cols) 85 ttt1 = MPI_Wtime() 86 if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) print *,'Time trans_ev_complex:',ttt1-ttt0 87 time_evp_back = ttt1-ttt0 88 89 deallocate(q_real) 90 deallocate(e, tau) 91 92end subroutine solve_evp_complex 93 94!------------------------------------------------------------------------------- 95 96