TABLE OF CONTENTS
m_abi_linalg/abi_chpgv [ Functions ]
[ Top ] [ m_abi_linalg ] [ Functions ]
NAME
abi_chpgv
FUNCTION
INPUTS
SOURCE
139 subroutine abi_chpgv(itype,jobz,uplo,n,a,b,w,z,ldz) 140 141 !Arguments ------------------------------------ 142 integer,intent(in) :: itype 143 character(len=1), intent(in) :: jobz 144 character(len=1), intent(in) :: uplo 145 integer, intent(in) :: n,ldz 146 complex(spc), intent(inout) :: a(:,:) 147 complex(spc), intent(inout) :: b(:,:) 148 complex(spc), intent(out) :: z(:,:) 149 real(sp), intent(out) :: w(:) 150 151 !Local variables------------------------------- 152 integer :: info 153 real(sp),pointer :: rwork(:) 154 complex(spc),pointer :: work(:) 155 156 ! ********************************************************************* 157 158 ABI_CHECK(lapack_packed_storage,"BUG(1) in abi_chpgv (storage)!") 159 ABI_CHECK(lapack_single_precision,"BUG(2) in abi_chpgv (precision)!") 160 ABI_CHECK(n<=eigen_c_maxsize,"BUG(3) in abi_chpgv (maxsize)!") 161 162 work => eigen_c_work ; rwork => eigen_c_rwork 163 164 !===== LAPACK 165 if (eigen_c_lwork==0) then 166 ABI_MALLOC(work,(2*n-1)) 167 end if 168 if (eigen_c_lrwork==0) then 169 ABI_MALLOC(rwork,(3*n-2)) 170 end if 171 call chpgv(itype,jobz,uplo,n,a,b,w,z,ldz,work,rwork,info) 172 if (eigen_c_lwork==0) then 173 ABI_FREE(work) 174 end if 175 if (eigen_c_lrwork==0) then 176 ABI_FREE(rwork) 177 end if 178 179 ABI_CHECK(info==0,"abi_chpgv returned info!=0!") 180 181 end subroutine abi_chpgv
m_abi_linalg/abi_dhpgv [ Functions ]
[ Top ] [ m_abi_linalg ] [ Functions ]
NAME
abi_dhpgv
FUNCTION
INPUTS
SOURCE
33 subroutine abi_dhpgv(itype,jobz,uplo,n,a,b,w,z,ldz,istwf_k,use_slk,use_gpu_elpa) 34 35 use m_fstrings, only : sjoin, itoa 36 37 !Arguments ------------------------------------ 38 integer :: itype 39 character(len=1), intent(in) :: jobz 40 character(len=1), intent(in) :: uplo 41 integer, intent(in) :: n,ldz 42 real(dp), intent(inout) :: a(:) 43 real(dp), intent(inout) :: b(:) 44 real(dp), intent(out) :: z(:,:) 45 real(dp), intent(out) :: w(:) 46 integer, optional, intent(in) :: istwf_k 47 integer, optional, intent(in) :: use_slk,use_gpu_elpa 48 49 !Local variables------------------------------- 50 integer :: info,use_slk_,use_gpu_elpa_,istwf_k_ 51 #ifdef HAVE_LINALG_SCALAPACK 52 type(matrix_scalapack) :: sca_a,sca_b,sca_ev 53 integer :: ierr 54 #endif 55 56 ! ********************************************************************* 57 58 ABI_CHECK(lapack_packed_storage,"BUG(1) in abi_dhpgv (storage)!") 59 ABI_CHECK(lapack_double_precision,"BUG(2) in abi_dhpgv (precision)!") 60 ABI_CHECK(n<=eigen_d_maxsize,"BUG(3) in abi_dhpgv (maxsize)!") 61 62 info = 0 !to avoid unwanted warning when info is not set by scalapack 63 64 use_slk_ = 0; if (present(use_slk)) use_slk_ = use_slk 65 istwf_k_ = 1; if (present(istwf_k)) istwf_k_ = istwf_k 66 use_gpu_elpa_=0 67 #ifdef HAVE_LINALG_ELPA 68 if (present(use_gpu_elpa)) use_gpu_elpa_=use_gpu_elpa 69 #endif 70 71 !===== SCALAPACK 72 if (ABI_LINALG_SCALAPACK_ISON.and.use_slk_==1.and.n>slk_minsize) then 73 #if defined HAVE_LINALG_SCALAPACK 74 z = zero 75 call sca_a%init(n,n,slk_processor,istwf_k_) 76 call sca_b%init(n,n,slk_processor,istwf_k_) 77 call sca_ev%init(n,n,slk_processor,istwf_k_) 78 #ifdef HAVE_LINALG_ELPA 79 call matrix_from_global_sym(sca_a,a,istwf_k_) 80 call matrix_from_global_sym(sca_b,b,istwf_k_) 81 #else 82 call matrix_from_global(sca_a,a,istwf_k_) 83 call matrix_from_global(sca_b,b,istwf_k_) 84 #endif 85 call compute_generalized_eigen_problem(slk_processor,sca_a,sca_b,& 86 & sca_ev,w,slk_communicator,istwf_k_,use_gpu_elpa=use_gpu_elpa_) 87 call matrix_to_global(sca_a,a,istwf_k_) 88 call matrix_to_global(sca_b,b,istwf_k_) 89 call matrix_to_reference(sca_ev,z,istwf_k_) 90 call xmpi_sum(z,slk_communicator,ierr) 91 call sca_a%free() 92 call sca_ev%free() 93 #endif 94 95 !===== LAPACK 96 else 97 if (istwf_k_/=2) then 98 call zhpgv(itype,jobz,uplo,n,a,b,w,z,ldz,eigen_z_work,eigen_z_rwork,info) 99 else 100 call dspgv(itype,jobz,uplo,n,a,b,w,z,ldz,eigen_d_work,info) 101 endif 102 end if 103 104 if (info < 0) then 105 ABI_COMMENT(sjoin("argument #", itoa(-info), "had an illegal value")) 106 end if 107 108 if (info > 0) then 109 ABI_COMMENT("DSPEV failed to converge") 110 if (info <= n) then 111 ABI_COMMENT(sjoin("DSPEV failed to converge;", itoa(info), " off-diagonal elements of")) 112 ABI_COMMENT(" an intermediate tridiagonal form did not converge to zero.") 113 else 114 ABI_COMMENT("The factorization of B could not be completed and no eigenvalues or eigenvectors were computed.") 115 endif 116 end if 117 118 ABI_CHECK(info==0,"abi_dhpgv returned info!=0!") 119 120 #ifndef HAVE_LINALG_ELPA 121 ABI_UNUSED(use_gpu_elpa) 122 #endif 123 124 end subroutine abi_dhpgv
m_abi_linalg/abi_xhpgv [ Functions ]
[ Top ] [ m_abi_linalg ] [ Functions ]
NAME
abi_xhpgv
FUNCTION
abi_xhpgv 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), stored in packed format 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_zhpgv [ Functions ]
[ Top ] [ m_abi_linalg ] [ Functions ]
NAME
abi_zhpgv
FUNCTION
INPUTS
SOURCE
196 subroutine abi_zhpgv(itype,jobz,uplo,n,a,b,w,z,ldz) 197 198 !Arguments ------------------------------------ 199 integer,intent(in) :: itype 200 integer, intent(in) :: n,ldz 201 character(len=1), intent(in) :: jobz 202 character(len=1), intent(in) :: uplo 203 complex(dpc), intent(inout) :: a(:,:) 204 complex(dpc), intent(inout) :: b(:,:) 205 complex(dpc), intent(out) :: z(:,:) 206 real(dp), intent(out) :: w(:) 207 208 !Local variables------------------------------- 209 integer :: info 210 real(dp),pointer :: rwork(:) 211 complex(dpc),pointer :: work(:) 212 213 ! ********************************************************************* 214 215 ABI_CHECK(lapack_packed_storage,"BUG(1) in abi_zhpgv (storage)!") 216 ABI_CHECK(lapack_double_precision,"BUG(2) in abi_zhpgv (precision)!") 217 ABI_CHECK(n<=eigen_z_maxsize,"BUG(3) in abi_zhpgv (maxsize)!") 218 219 work => eigen_z_work ; rwork => eigen_z_rwork 220 221 !===== LAPACK 222 if (eigen_z_lwork==0) then 223 ABI_MALLOC(work,(2*n-1)) 224 end if 225 if (eigen_z_lrwork==0) then 226 ABI_MALLOC(rwork,(3*n-2)) 227 end if 228 call zhpgv(itype,jobz,uplo,n,a,b,w,z,ldz,work,rwork,info) 229 if (eigen_z_lwork==0) then 230 ABI_FREE(work) 231 end if 232 if (eigen_z_lrwork==0) then 233 ABI_FREE(rwork) 234 end if 235 236 ABI_CHECK(info==0,"abi_zhpgv returned info!=0!") 237 238 end subroutine abi_zhpgv