TABLE OF CONTENTS


m_abi_linalg/abi_chegv [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_chegv

FUNCTION

INPUTS

PARENTS

SOURCE

239 subroutine abi_chegv(itype,jobz,uplo,n,a,lda,b,ldb,w,work,lwork,rwork,info)
240 
241 !This section has been created automatically by the script Abilint (TD).
242 !Do not modify the following lines by hand.
243 #undef ABI_FUNC
244 #define ABI_FUNC 'abi_chegv'
245 !End of the abilint section
246 
247  implicit none
248 
249  !Arguments ------------------------------------
250  integer,intent(in) :: itype
251  character(len=1), intent(in) :: jobz
252  character(len=1), intent(in) :: uplo
253  integer, intent(in) :: n
254  integer, intent(in) :: lda
255  integer, intent(in) :: ldb
256  integer, intent(in) :: lwork
257  complex(spc),target, intent(inout) :: a(lda,*)
258  complex(spc),target, intent(inout) :: b(ldb,*)
259  complex(spc),target, intent(inout) :: work(*)
260  real(sp), intent(inout) :: rwork(lwork)
261  real(sp), intent(out) :: w(n)
262  integer, intent(out) :: info
263 
264 #ifdef HAVE_LINALG_PLASMA
265  integer :: jobz_plasma_a
266  type(c_ptr) :: plasma_work
267 #endif
268     ! *********************************************************************
269 
270 #ifdef HAVE_LINALG_PLASMA
271  ! FDahm & LNGuyen  (November 2012) :
272  !  In Plasma v 2.4.6, eigen routines support only
273  !the eigenvalues computation (jobz=N) and not the
274  !full eigenvectors bases determination (jobz=V)
275  if (LSAME(jobz,'N')) then
276     jobz_plasma_a = jobz_plasma(jobz)
277 
278     call PLASMA_Alloc_Workspace_chegv(n,n,plasma_work,info)
279 
280     MSG_ERROR("This code is broken")
281     !call PLASMA_chegv(itype,jobz_plasma_a,uplo_plasma(uplo),n,a,lda,b,ldb,w,plasma_work,rwork,lwork,info)
282 
283     call PLASMA_Dealloc_handle(plasma_work,info)
284  else
285 #endif
286    call chegv(itype,jobz,uplo,n,a,lda,b,ldb,w,work,lwork,rwork,info)
287 #ifdef HAVE_LINALG_PLASMA
288  end if
289 #endif
290 
291 end subroutine abi_chegv

m_abi_linalg/abi_chegv_alloc [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_chegv_alloc

FUNCTION

INPUTS

PARENTS

SOURCE

307  subroutine abi_chegv_alloc(itype,jobz,uplo,n,a,b,w)
308 
309 !This section has been created automatically by the script Abilint (TD).
310 !Do not modify the following lines by hand.
311 #undef ABI_FUNC
312 #define ABI_FUNC 'abi_chegv_alloc'
313 
314 !End of the abilint section
315 
316  implicit none
317 
318 !Arguments ------------------------------------
319  integer :: itype
320  character(len=1), intent(in) :: jobz
321  character(len=1), intent(in) :: uplo
322  integer, intent(in) :: n
323  complex(spc), DEV_CONTARRD intent(inout) :: a(:,:) ! here I need lda to be consistent
324  complex(spc), DEV_CONTARRD intent(inout) :: b(:,:)
325  real(sp), intent(out) :: w(n)
326  integer :: info
327 
328 ! *********************************************************************
329 
330  call abi_chegv(itype,jobz,uplo,n,a,n,b,n,w,eigen_c_work,eigen_c_lwork,eigen_c_rwork,info)
331 
332 end subroutine abi_chegv_alloc

m_abi_linalg/abi_dhegv [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_dhegv

FUNCTION

INPUTS

PARENTS

SOURCE

 38   subroutine abi_dhegv(itype,jobz,uplo,n,a,lda,b,ldb,w,work,lwork,rwork,info, &
 39 &       x_cplx,istwf_k,timopt,tim_xeigen,use_slk,use_gpu)
 40 
 41 !This section has been created automatically by the script Abilint (TD).
 42 !Do not modify the following lines by hand.
 43 #undef ABI_FUNC
 44 #define ABI_FUNC 'abi_dhegv'
 45 !End of the abilint section
 46 
 47  implicit none
 48 !Arguments ------------------------------------
 49  integer :: itype
 50  character(len=1), intent(in) :: jobz
 51  character(len=1), intent(in) :: uplo
 52  integer, intent(in) :: n
 53  integer, intent(in) :: lda
 54  integer, intent(in) :: ldb
 55  integer, intent(in) :: lwork
 56  real(dp), target, intent(inout) :: a(lda,*) ! FIXME should be cplex * lda
 57  real(dp), target, intent(inout) :: b(ldb,*)
 58  real(dp), target, intent(inout) :: work(*)
 59  real(dp),target,optional, intent(inout) :: rwork(lwork)
 60  real(dp),target,intent(out) :: w(n)
 61  integer, intent(out) :: info
 62 
 63 !Optional Arguments ------------------------------------
 64  integer, optional, intent(in) :: x_cplx
 65  integer, optional, intent(in) :: istwf_k
 66  integer, optional, intent(in) :: use_slk
 67  integer, optional, intent(in) :: timopt,tim_xeigen,use_gpu
 68 
 69 !Local variables ------------------------------------
 70  character(len=500) :: msg
 71  integer :: istwf_k_ , usegpu_, use_slk_, cplx_
 72  real(dp) :: tsec(2)
 73 
 74 #ifdef HAVE_LINALG_MAGMA
 75  integer :: liwork,lwork_, lrwork_
 76  integer, dimension(:), allocatable :: iwork
 77 #endif
 78 #ifdef HAVE_LINALG_PLASMA
 79  integer :: jobz_plasma_a
 80  type(c_ptr) :: plasma_work
 81 #endif
 82 
 83 ! *********************************************************************
 84 
 85  if (present(tim_xeigen).and.present(timopt)) then
 86     if(abs(timopt)==3) then
 87       call timab(tim_xeigen,1,tsec)
 88     end if
 89  end if
 90 
 91  usegpu_=0; if (present(use_gpu)) usegpu_=use_gpu
 92  use_slk_=0; if (present(use_slk)) use_slk_=use_slk
 93  istwf_k_=1; if (present(istwf_k)) istwf_k_=istwf_k
 94  cplx_=1 ; if(PRESENT(x_cplx)) cplx_ = x_cplx
 95 
 96  if( n > eigen_d_maxsize ) then
 97     write(msg,'(a,2i3)')' Eigen size higher than max size set!!',n,eigen_d_maxsize
 98     MSG_ERROR(msg)
 99  endif
100 
101 #ifdef HAVE_LINALG_MAGMA
102  if (usegpu_==1) then
103     lrwork_= 1 + 5*n + 2*(n**2)
104     liwork= 3 + 5*n
105     ABI_ALLOCATE(iwork,(liwork))
106     if (cplx_ == 2 .and. present(rwork)) then
107        lwork_=33*n + (n**2)
108        call magmaf_zhegvd(itype,jobz,uplo,n,a,lda,b,ldb,w, &
109 &            work(1:2*lwork_),lwork_,rwork(1:lrwork_),lrwork_,iwork(1:liwork),liwork,info)
110     else
111        lwork_= 1+n*(6 + 2*n)
112        call magmaf_dsygvd(itype,jobz,uplo,n,a,lda,b,ldb,w,&
113 &            work(1:lwork_),lwork_,iwork(1:liwork),liwork,info)
114     endif
115     ABI_DEALLOCATE(iwork)
116  else
117 #endif
118 
119 #ifdef HAVE_LINALG_SCALAPACK
120  if(use_slk_ == 1) then
121    ABI_CHECK(PRESENT(x_cplx),"x_cplx must be present")
122    call compute_eigen2(abi_communicator,abi_processor,cplx_,n,n,a,b,w,istwf_k_)
123    info = 0 ! This is to avoid unwanted warning but it's not clean
124  else
125 #endif
126 
127 #ifdef HAVE_LINALG_PLASMA
128  ! FDahm & LNGuyen  (November 2012) :
129  !  In Plasma v 2.4.6, eigen routines support only
130  !the eigenvalues computation (jobz=N) and not the
131  !full eigenvectors bases determination (jobz=V)
132  if (LSAME(jobz,'N')) then
133     jobz_plasma_a = jobz_plasma(jobz)
134 
135     MSG_ERROR("This code is broken")
136     if (cplx_ == 2 .and. present(rwork)) then
137        call PLASMA_Alloc_Workspace_zhegv(n,n,plasma_work,info)
138 
139 !       info = PLASMA_zhegv_c(itype,jobz_plasma_a,uplo_plasma(uplo),&
140 !&        n,c_loc(a),lda,c_loc(b),ldb,c_loc(w),plasma_work,c_loc(rwork),lwork)
141 
142        call PLASMA_Dealloc_handle(plasma_work,info)
143     else
144        call PLASMA_Alloc_Workspace_dsygv(n,n,plasma_work,info)
145 
146 !       info = PLASMA_dsygv_c(itype,jobz_plasma_a,uplo_plasma(uplo),&
147 !&        n,c_loc(a),lda,c_loc(b),ldb,c_loc(w),plasma_work,c_loc(lwork),lwork)
148 
149        call PLASMA_Dealloc_handle(plasma_work,info)
150     endif
151  else
152 #endif
153  if ( cplx_ == 2 .and. present(rwork)) then
154     call zhegv(itype,jobz,uplo,n,a,lda,b,ldb,w,work,lwork/2,rwork,info)
155  else
156     call dsygv(itype,jobz,uplo,n,a,lda,b,ldb,w,work,lwork,info)
157  endif
158 #ifdef HAVE_LINALG_PLASMA
159  end if
160 #endif
161 #ifdef HAVE_LINALG_SCALAPACK
162  end if
163 #endif
164 #ifdef HAVE_LINALG_MAGMA
165  end if
166 #endif
167 
168  if(info/=0) then
169     write(msg,'(a,i3)')' Problem in abi_xhegv, info= ',info
170     MSG_ERROR(msg)
171  endif
172 
173  if (present(tim_xeigen).and.present(timopt)) then
174     if(abs(timopt)==3) then
175       call timab(tim_xeigen,2,tsec)
176     end if
177  end if
178 
179 end subroutine abi_dhegv

m_abi_linalg/abi_dhegv_alloc [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_dhegv_alloc

FUNCTION

INPUTS

PARENTS

SOURCE

194  subroutine abi_dhegv_alloc(itype,jobz,uplo,n,a,b,w, &
195 &       x_cplx,istwf_k,timopt,tim_xeigen,use_slk,use_gpu)
196 
197 !This section has been created automatically by the script Abilint (TD).
198 !Do not modify the following lines by hand.
199 #undef ABI_FUNC
200 #define ABI_FUNC 'abi_dhegv_alloc'
201 
202 !End of the abilint section
203 
204  implicit none
205 
206 !Arguments ------------------------------------
207  integer,intent(in) :: itype
208  character(len=1), intent(in) :: jobz
209  character(len=1), intent(in) :: uplo
210  integer, intent(in) :: n
211  real(dp), DEV_CONTARRD intent(inout) :: a(:,:)  ! here I need lda to be consistent
212  real(dp), DEV_CONTARRD intent(inout) :: b(:,:)
213  real(dp), intent(out) :: w(n)
214  integer :: info
215 
216 !Optional Arguments ------------------------------------
217  integer, optional, intent(in) :: x_cplx
218  integer, optional, intent(in) :: istwf_k
219  integer, optional, intent(in) :: timopt,tim_xeigen,use_slk, use_gpu
220 
221  call  abi_dhegv(itype,jobz,uplo,n,a,n,b,n,w,eigen_d_work,eigen_d_lwork,eigen_z_rwork,info, &
222 &         x_cplx,istwf_k,timopt,tim_xeigen,use_slk,use_gpu)
223 
224 end subroutine abi_dhegv_alloc

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-2018 ABINIT group (LNguyen,FDahm (CS))
  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

PARENTS

SOURCE

348 subroutine abi_zhegv(itype,jobz,uplo,n,a,lda,b,ldb,w,work,lwork,rwork,info)
349 
350 !This section has been created automatically by the script Abilint (TD).
351 !Do not modify the following lines by hand.
352 #undef ABI_FUNC
353 #define ABI_FUNC 'abi_zhegv'
354 !End of the abilint section
355 
356  implicit none
357 
358 !Arguments ------------------------------------
359  integer,intent(in) :: itype
360  character(len=1), intent(in) :: jobz
361  character(len=1), intent(in) :: uplo
362  integer, intent(in) :: n
363  integer, intent(in) :: lda
364  integer, intent(in) :: ldb
365  integer, intent(in) :: lwork
366  complex(dpc), target,intent(inout) :: a(lda,*)
367  complex(dpc), target,intent(inout) :: b(ldb,*)
368  complex(dpc), target,intent(inout) :: work(*)
369  real(dp), intent(inout) :: rwork(lwork)
370  real(dp), intent(out) :: w(n)
371  integer, intent(out) :: info
372 
373 #ifdef HAVE_LINALG_PLASMA
374  !Optional Arguments ------------------------------------
375  integer :: jobz_plasma_a
376  type(c_ptr) :: plasma_work
377 #endif
378 
379  ! *********************************************************************
380 
381 #ifdef HAVE_LINALG_PLASMA
382  ! FDahm & LNGuyen  (November 2012) :
383  !  In Plasma v 2.4.6, eigen routines support only
384  !the eigenvalues computation (jobz=N) and not the
385  !full eigenvectors bases determination (jobz=V)
386  if (LSAME(jobz,'N')) then
387     jobz_plasma_a = jobz_plasma(jobz)
388 
389     call PLASMA_Alloc_Workspace_zhegv(n,n,plasma_work,info)
390 
391     MSG_ERROR("This code is broken")
392     !call PLASMA_zhegv(itype,jobz_plasma_a,uplo_plasma(uplo),n,a,lda,b,ldb,w,plasma_work,rwork,lwork,info)
393 
394     call PLASMA_Dealloc_handle(plasma_work,info)
395  else
396 #endif
397    call zhegv(itype,jobz,uplo,n,a,lda,b,ldb,w,work,lwork,rwork,info)
398 #ifdef HAVE_LINALG_PLASMA
399  end if
400 #endif
401 
402 end subroutine abi_zhegv

m_abi_linalg/abi_zhegv_alloc [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_zhegv_alloc

FUNCTION

INPUTS

PARENTS

SOURCE

417 subroutine abi_zhegv_alloc(itype,jobz,uplo,n,a,b,w)
418 
419 !This section has been created automatically by the script Abilint (TD).
420 !Do not modify the following lines by hand.
421 #undef ABI_FUNC
422 #define ABI_FUNC 'abi_zhegv_alloc'
423 
424 !End of the abilint section
425 
426  implicit none
427 
428  !Arguments ------------------------------------
429  integer,intent(in) :: itype
430  character(len=1), intent(in) :: jobz
431  character(len=1), intent(in) :: uplo
432  integer, intent(in) :: n
433  complex(dpc), DEV_CONTARRD intent(inout) :: a(:,:) ! here I need lda to be consistent
434  complex(dpc), DEV_CONTARRD intent(inout) :: b(:,:)
435  real(dp), intent(out) :: w(n)
436  integer :: info
437 
438 ! *********************************************************************
439  call abi_zhegv(itype,jobz,uplo,n,a,n,b,n,w,eigen_z_work,eigen_z_lwork,eigen_z_rwork,info)
440 
441 end subroutine abi_zhegv_alloc