TABLE OF CONTENTS


m_abi_linalg/abi_chegv [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_chegv

FUNCTION

INPUTS

SOURCE

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

m_abi_linalg/abi_dhegv [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_dhegv

FUNCTION

INPUTS

SOURCE

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

m_abi_linalg/abi_xhegv [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

  abi_xhegv

FUNCTION

  abi_xhegv is the generic function that compute
  all eigenvalues and, optionally, eigenvectors of a
  generalized symmetric-definite eigenproblem, of the form
  A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.
  Here A and B are assumed to be symmetric (or hermitian) and
  B is also positive definite.

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_zhegv [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_zhegv

FUNCTION

INPUTS

SOURCE

234 subroutine abi_zhegv(itype,jobz,uplo,n,a,lda,b,ldb,w)
235 
236  implicit none
237 
238 !Arguments ------------------------------------
239  integer,intent(in) :: itype
240  character(len=1), intent(in) :: jobz
241  character(len=1), intent(in) :: uplo
242  integer, intent(in) :: n,lda,ldb
243  complex(dpc), intent(inout) :: a(lda,*)
244  complex(dpc), intent(inout) :: b(ldb,*)
245  real(dp), intent(out) :: w(n)
246 
247 !Local variables-------------------------------
248  integer :: info,lwork
249  real(dp),pointer :: rwork(:)
250  complex(dpc),pointer :: work(:)
251 
252 ! *********************************************************************
253 
254  ABI_CHECK(lapack_full_storage,"BUG(1) in abi_zhegv (storage)!")
255  ABI_CHECK(lapack_double_precision,"BUG(2) in abi_zhegv (precision)!")
256  ABI_CHECK(n<=eigen_z_maxsize,"BUG(3) in abi_zhegv (maxsize)!")
257 
258  work => eigen_z_work ; rwork => eigen_z_rwork
259  lwork=eigen_z_lwork
260 
261 !===== PLASMA
262  !FDahm & LNGuyen  (November 2012) :
263  !  In Plasma v 2.4.6, eigen routines support only
264  !  the eigenvalues computation (jobz=N) and not the
265  ! full eigenvectors bases determination (jobz=V)
266  if (ABI_LINALG_PLASMA_ISON.and.LSAME(jobz,'N')) then
267 #if defined HAVE_LINALG_PLASMA
268    if (eigen_z_lwork==0) then
269      ABI_MALLOC(work,(n**2))
270    end if
271    call PLASMA_Alloc_Workspace_zhegv(n,n,plasma_work,info)
272    call PLASMA_zhegv(itype,jobz_plasma(jobz),uplo_plasma(uplo),n,a,lda,b,ldb,w,&
273 &                    plasma_work,c_loc(work),n)
274    call PLASMA_Dealloc_handle(plasma_work,info)
275    if (eigen_z_lwork==0) then
276      ABI_FREE(work)
277    end if
278 #endif
279 
280 !===== LAPACK
281  else
282    if (eigen_z_lwork==0) then
283      lwork=2*n-1
284      ABI_MALLOC(work,(lwork))
285    end if
286    if (eigen_z_lrwork==0) then
287      ABI_MALLOC(rwork,(3*n-2))
288    end if
289    call zhegv(itype,jobz,uplo,n,a,lda,b,ldb,w,work,lwork,rwork,info)
290    if (eigen_z_lwork==0) then
291      ABI_FREE(work)
292    end if
293    if (eigen_z_lrwork==0) then
294      ABI_FREE(rwork)
295    end if
296  end if
297 
298  ABI_CHECK(info==0,"abi_zhegv returned info!=0!")
299 
300 end subroutine abi_zhegv