TABLE OF CONTENTS


m_abi_linalg/abi_cheev [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_cheev

FUNCTION

INPUTS

SOURCE

146 subroutine abi_cheev(jobz,uplo,n,a,lda,w)
147 
148 !Arguments ------------------------------------
149  character(len=1), intent(in) :: jobz
150  character(len=1), intent(in) :: uplo
151  integer, intent(in) :: n,lda
152  complex(spc), intent(inout) :: a(lda,*)
153  real(sp), intent(out) :: w(n)
154 
155 !Local variables-------------------------------
156  integer :: info,lwork
157  real(sp),pointer :: rwork(:)
158  complex(spc),pointer :: work(:)
159 
160 ! *********************************************************************
161 
162  ABI_CHECK(lapack_full_storage,"BUG(1) in abi_cheev (storage)!")
163  ABI_CHECK(lapack_single_precision,"BUG(2) in abi_cheev (precision)!")
164  ABI_CHECK(n<=eigen_c_maxsize,"BUG(3) in abi_cheev (maxsize)!")
165 
166  work => eigen_c_work ; rwork => eigen_c_rwork
167  lwork=eigen_c_lwork
168 
169 !===== PLASMA
170  !FDahm & LNGuyen  (November 2012) :
171  !  In Plasma v 2.4.6, eigen routines support only
172  !  the eigenvalues computation (jobz=N) and not the
173  ! full eigenvectors bases determination (jobz=V)
174  if (ABI_LINALG_PLASMA_ISON.and.LSAME(jobz,'N')) then
175 #if defined HAVE_LINALG_PLASMA
176    if (eigen_c_lwork==0) then
177      ABI_MALLOC(work,(n**2))
178    end if
179    call PLASMA_Alloc_Workspace_cheev(n,n,plasma_work,info)
180    info = PLASMA_cheev_c(jobz_plasma(jobz),uplo_plasma(uplo),n,c_loc(a),lda,c_loc(w),&
181 &                        plasma_work,c_loc(work),n)
182    call PLASMA_Dealloc_handle(plasma_work,info)
183    if (eigen_c_lwork==0) then
184      ABI_FREE(work)
185    end if
186 #endif
187 
188 !===== LAPACK
189  else
190    if (eigen_c_lwork==0) then
191      lwork=2*n-1
192      ABI_MALLOC(work,(lwork))
193    end if
194    if (eigen_c_lrwork==0) then
195      ABI_MALLOC(rwork,(3*n-2))
196    end if
197    call cheev(jobz,uplo,n,a,lda,w,work,lwork,rwork,info)
198    if (eigen_c_lwork==0) then
199      ABI_FREE(work)
200    end if
201    if (eigen_c_lrwork==0) then
202      ABI_FREE(rwork)
203    end if
204  end if
205 
206  ABI_CHECK(info==0,"abi_cheev returned info!=!0")
207 
208 end subroutine abi_cheev

m_abi_linalg/abi_dheev [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_dheev

FUNCTION

INPUTS

SOURCE

 31   subroutine abi_dheev(jobz,uplo,n,a,lda,w,&
 32 &            x_cplx,istwf_k,timopt,tim_xeigen,&
 33 &            use_gpu_elpa,use_gpu_magma,use_slk)
 34 
 35 !Arguments ------------------------------------
 36  character(len=1), intent(in) :: jobz
 37  character(len=1), intent(in) :: uplo
 38  integer, intent(in) :: n,lda
 39  real(dp), intent(inout) :: a(n,*)  ! FIXME should be cplex * lda
 40  real(dp), intent(out) :: w(n)
 41  integer, optional, intent(in) :: istwf_k
 42  integer, optional, intent(in) :: x_cplx
 43  integer, optional, intent(in) :: timopt,tim_xeigen
 44  integer, optional, intent(in) :: use_gpu_elpa,use_gpu_magma,use_slk
 45 
 46 !Local variables-------------------------------
 47  integer :: cplx_,istwf_k_,use_gpu_elpa_,use_gpu_magma_,use_slk_
 48  integer :: info
 49  real(dp) :: tsec(2)
 50 
 51 ! *********************************************************************
 52 
 53  ABI_CHECK(lapack_full_storage,"BUG(1) in abi_dheev (storage)!")
 54  ABI_CHECK(lapack_double_precision,"BUG(2) in abi_dheev (precision)!")
 55  ABI_CHECK(n<=eigen_d_maxsize,"BUG(3) in abi_dheev (maxsize)!")
 56 
 57  if (present(tim_xeigen).and.present(timopt)) then
 58    if(abs(timopt)==3) call timab(tim_xeigen,1,tsec)
 59  end if
 60 
 61  cplx_=1 ; if(present(x_cplx)) cplx_ = x_cplx
 62  use_gpu_magma_=0;if (present(use_gpu_magma)) use_gpu_magma_=use_gpu_magma
 63  istwf_k_=1;if (present(istwf_k)) istwf_k_=istwf_k
 64  use_slk_ = 0; if(present(use_slk)) use_slk_ = use_slk
 65  use_gpu_elpa_=0
 66 #ifdef HAVE_LINALG_ELPA
 67  if (present(use_gpu_elpa)) use_gpu_elpa_=use_gpu_elpa
 68 #endif
 69 
 70 !===== MAGMA
 71  if (ABI_LINALG_MAGMA_ISON.and.use_gpu_magma_==1) then
 72 #if defined HAVE_LINALG_MAGMA
 73    ABI_CHECK((lapack_divide_conquer),"BUG(4) in abi_dheev (d&c)!")
 74    if (cplx_ == 2) then
 75      call magmaf_zheevd(jobz,uplo,n,a,lda,w,eigen_z_work,eigen_z_lwork, &
 76 &           eigen_z_rwork,eigen_z_lrwork,eigen_iwork,eigen_liwork,info)
 77    else
 78      call magmaf_dsyevd(jobz,uplo,n,a,lda,w,eigen_d_work,eigen_d_lwork,&
 79 &                       eigen_iwork,eigen_liwork,info)
 80    endif
 81 #endif
 82 
 83 !===== SCALAPACK
 84  else if (ABI_LINALG_SCALAPACK_ISON.and.use_slk_==1.and.n>slk_minsize) then
 85 #if defined HAVE_LINALG_SCALAPACK
 86    ABI_CHECK(present(x_cplx),"BUG(5) in abi_dheev (x_cplx)!")
 87    call compute_eigen1(slk_communicator,slk_processor,cplx_,n,n,a,w,istwf_k_,&
 88 &                      use_gpu_elpa=use_gpu_elpa_)
 89    info = 0 ! This is to avoid unwanted warning but it's not clean
 90 #endif
 91 
 92 !===== PLASMA
 93  !FDahm & LNGuyen  (November 2012) :
 94  !  In Plasma v 2.4.6, eigen routines support only
 95  !  the eigenvalues computation (jobz=N) and not the
 96  !  full eigenvectors bases determination (jobz=V)
 97  else if (ABI_LINALG_PLASMA_ISON.and.LSAME(jobz,'N')) then
 98 #if defined HAVE_LINALG_PLASMA
 99    jobz_plasma_a = jobz_plasma(jobz)
100    if (cplx_ == 2) then
101      call PLASMA_Alloc_Workspace_zheev(n,n,plasma_work,info)
102      info = PLASMA_zheev_c(jobz_plasma(jobz),uplo_plasma(uplo),n,c_loc(a),lda,c_loc(w),&
103 &                          plasma_work,c_loc(eigen_z_work),n)
104    else
105      call PLASMA_Alloc_Workspace_dsyev(n,n,plasma_work,info)
106      info = PLASMA_dsyev_c(jobz_plasma(jobz),uplo_plasma(uplo),n,c_loc(a),lda,c_loc(w),&
107 &                          plasma_work,c_loc(eigen_d_work),n)
108    endif
109    call PLASMA_Dealloc_handle(plasma_work,info)
110 #endif
111 
112 !===== LAPACK
113  else
114    if (cplx_ == 2) then
115      call zheev(jobz,uplo,n,a,lda,w,eigen_z_work,eigen_z_lwork,eigen_z_rwork,info)
116    else
117      call dsyev(jobz,uplo,n,a,lda,w,eigen_d_work,eigen_d_lwork,info)
118    endif
119  end if
120 
121  if (present(tim_xeigen).and.present(timopt)) then
122    if(abs(timopt)==3) call timab(tim_xeigen,2,tsec)
123  end if
124 
125  ABI_CHECK(info==0,"abi_dheev returned info!=0!")
126 
127 #ifndef HAVE_LINALG_ELPA
128  ABI_UNUSED(use_gpu_elpa)
129 #endif
130 
131 end subroutine abi_dheev

m_abi_linalg/abi_xheev [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

  abi_xheev

FUNCTION

  abi_xheev is the generic function that compute
  all eigenvalues and, optionally, eigenvectors of a
  symmetric or hermitian matrix A.

COPYRIGHT

  Copyright (C) 2001-2024 ABINIT group (LNguyen,FDahm,MT)
  This file is distributed under the terms of the
  GNU General Public License, see ~ABINIT/Infos/copyright
  or http://www.gnu.org/copyleft/gpl.txt .

SOURCE


m_abi_linalg/abi_zheev [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_zheev

FUNCTION

INPUTS

SOURCE

223 subroutine abi_zheev(jobz,uplo,n,a,lda,w)
224 
225 !Arguments ------------------------------------
226  character(len=1), intent(in) :: jobz
227  character(len=1), intent(in) :: uplo
228  integer, intent(in) :: n,lda
229  complex(dpc), intent(inout) :: a(lda,*)
230  real(dp), intent(out) :: w(n)
231 
232 !Local variables-------------------------------
233  integer :: info,lwork
234  real(dp),pointer :: rwork(:)
235  complex(dpc),pointer :: work(:)
236 
237 ! *********************************************************************
238 
239  ABI_CHECK(lapack_full_storage,"BUG(1) in abi_zheev (storage)!")
240  ABI_CHECK(lapack_double_precision,"BUG(2) in abi_zheev (precision)!")
241  ABI_CHECK(n<=eigen_z_maxsize,"BUG(3) in abi_zheev (maxsize)!")
242 
243  work => eigen_z_work ; rwork => eigen_z_rwork
244  lwork=eigen_z_lwork
245 
246 !===== PLASMA
247  !FDahm & LNGuyen  (November 2012) :
248  !  In Plasma v 2.4.6, eigen routines support only
249  ! the eigenvalues computation (jobz=N) and not the
250  ! full eigenvectors bases determination (jobz=V)
251  if (ABI_LINALG_PLASMA_ISON.and.LSAME(jobz,'N')) then
252 #if defined HAVE_LINALG_PLASMA
253    if (eigen_z_lwork==0) then
254      ABI_MALLOC(work,(n**2))
255    end if
256    call PLASMA_Alloc_Workspace_zheev(n,n,plasma_work,info)
257    info = PLASMA_zheev_c(jobz_plasma(jobz),uplo_plasma(uplo),&
258 &                        plasma_work,c_loc(work),n)
259    call PLASMA_Dealloc_handle(plasma_work,info)
260    if (eigen_z_lwork==0) then
261      ABI_FREE(work)
262    end if
263 #endif
264 
265 !===== LAPACK
266  else
267    if (eigen_z_lwork==0) then
268      lwork=2*n-1
269      ABI_MALLOC(work,(lwork))
270    end if
271    if (eigen_z_lrwork==0) then
272      ABI_MALLOC(rwork,(3*n-2))
273    end if
274   call zheev(jobz,uplo,n,a,lda,w,work,lwork,rwork,info)
275    if (eigen_z_lwork==0) then
276      ABI_FREE(work)
277    end if
278    if (eigen_z_lrwork==0) then
279      ABI_FREE(rwork)
280    end if
281  end if
282 
283  ABI_CHECK(info==0,"abi_zheev returned info !=0!")
284 
285 end subroutine abi_zheev