Merge branch 'develop01' into 'develop'

Debug of USPP in EPW

See merge request QEF/q-e!533
This commit is contained in:
giannozz 2019-07-17 11:42:16 +00:00
commit 39ca8859df
13 changed files with 416 additions and 372 deletions

View File

@ -73,7 +73,7 @@
COMPLEX(kind=DP) :: sum_nc(npol)
!! auxiliary sum variable non-collinear case
!
IF ( .NOT. okvan) RETURN
IF (.NOT. okvan) RETURN
!
CALL start_clock('adddvscf2')
!

View File

@ -94,7 +94,7 @@
!! e^{-i q * \tau} * conjg(e^{-i q * \tau})
COMPLEX(kind=DP) :: fact1
!! -i * omega
COMPLEX(kind=DP), EXTERNAL :: zdotc
COMPLEX(kind=DP), EXTERNAL :: ZDOTC
!! the scalar product function
COMPLEX(kind=DP), ALLOCATABLE :: aux1(:), aux2(:), &
aux3(:), aux5(:), sk(:)
@ -105,7 +105,7 @@
COMPLEX(kind=DP), POINTER :: qgmq(:)
!! the augmentation function at q+G
!
IF ( .NOT. okvan) RETURN
IF (.NOT. okvan) RETURN
!
CALL start_clock('dvanqq2')
!
@ -113,6 +113,7 @@
int2(:,:,:,:,:) = czero
int4(:,:,:,:,:) = czero
int5(:,:,:,:,:) = czero
!
ALLOCATE ( sk(ngm) )
ALLOCATE ( aux1(ngm) )
ALLOCATE ( aux2(ngm) )
@ -124,17 +125,32 @@
ALLOCATE ( qgm(ngm))
ALLOCATE ( ylmk0(ngm, lmaxq * lmaxq) )
ALLOCATE ( ylmkq(ngm, lmaxq * lmaxq) )
sk(:) = czero
aux1(:) = czero
aux2(:) = czero
aux3(:) = czero
aux5(:) = czero
qmodg(:) = zero
qmod(:) = zero
qgmq(:) = czero
qgm(:) = czero
ylmk0(:,:) = zero
ylmkq(:,:) = zero
!
! compute spherical harmonics
!
CALL ylmr2( lmaxq * lmaxq, ngm, g, gg, ylmk0 )
!
DO ig = 1, ngm
qmodg(ig) = sqrt( gg(ig) )
ENDDO
!
ALLOCATE ( qpg(3, ngm) )
qpg(:,:) = zero
!
CALL setqmod( ngm, xq, g, qmod, qpg )
CALL ylmr2(lmaxq * lmaxq, ngm, qpg, qmod, ylmkq)
!
DEALLOCATE (qpg)
DO ig = 1, ngm
qmod(ig) = sqrt( qmod(ig) )
@ -143,8 +159,10 @@
! we start by computing the FT of the effective potential
!
ALLOCATE (veff(dfftp%nnr,nspin_mag))
veff(:,:) = czero
!
DO is = 1, nspin_mag
IF (nspin_mag.ne.4 .or. is==1) THEN
IF (nspin_mag /= 4 .OR. is == 1) THEN
DO ir = 1, dfftp%nnr
veff(ir,is) = CMPLX(vltot(ir) + v%of_r(ir,is), zero, kind=DP)
ENDDO
@ -200,7 +218,7 @@
aux5(ig) = sk(ig) * ( g(ipol,ig) + xq(ipol) )
ENDDO
int2(ih,jh,ipol,na,nb) = fact * fact1 * &
zdotc(ngm, aux1, 1, aux5, 1)
ZDOTC(ngm, aux1, 1, aux5, 1)
!
DO jpol = 1, 3
IF (jpol >= ipol) THEN
@ -210,7 +228,7 @@
ENDDO
int5(ijh,ipol,jpol,na,nb) = &
conjg(fact) * tpiba2 * omega * &
zdotc(ngm, aux3, 1, aux1, 1)
ZDOTC(ngm, aux3, 1, aux1, 1)
ELSE
int5(ijh,ipol,jpol,na,nb) = &
int5(ijh,jpol,ipol,na,nb)
@ -232,14 +250,14 @@
aux2(ig) = veff(dfftp%nl(ig),is) * g(ipol,ig)
ENDDO
int1(ih,jh,ipol,nb,is) = - fact1 * &
zdotc(ngm, aux1, 1, aux2, 1)
ZDOTC(ngm, aux1, 1, aux2, 1)
DO jpol = 1, 3
IF (jpol >= ipol) THEN
DO ig = 1, ngm
aux3(ig) = aux2(ig) * g(jpol,ig)
ENDDO
int4(ijh,ipol,jpol,nb,is) = - tpiba2 * &
omega * zdotc(ngm, aux3, 1, aux1, 1)
omega * ZDOTC(ngm, aux3, 1, aux1, 1)
ELSE
int4(ijh,ipol,jpol,nb,is) = &
int4(ijh,jpol,ipol,nb,is)
@ -304,14 +322,14 @@
ENDIF
!
!DBRM
!write(*,'(a,e20.12)') 'int1 = ', &
!SUM((REAL(REAL(int1(:,:,:,:,:))))**2)+SUM((REAL(AIMAG(int1(:,:,:,:,:))))**2)
!write(*,'(a,e20.12)') 'int2 = ', &
!SUM((REAL(REAL(int2(:,:,:,:,:))))**2)+SUM((REAL(AIMAG(int2(:,:,:,:,:))))**2)
!write(*,'(a,e20.12)') 'int4 = ', &
!SUM((REAL(REAL(int4(:,:,:,:,:))))**2)+SUM((REAL(AIMAG(int4(:,:,:,:,:))))**2)
!write(*,'(a,e20.12)') 'int5 = ', &
!SUM((REAL(REAL(int5(:,:,:,:,:))))**2)+SUM((REAL(AIMAG(int5(:,:,:,:,:))))**2)
! write(*,'(a,e20.12)') 'int1 = ', &
! SUM((REAL(REAL(int1(:,:,:,:,:))))**2)+SUM((REAL(AIMAG(int1(:,:,:,:,:))))**2)
! write(*,'(a,e20.12)') 'int2 = ', &
! SUM((REAL(REAL(int2(:,:,:,:,:))))**2)+SUM((REAL(AIMAG(int2(:,:,:,:,:))))**2)
! write(*,'(a,e20.12)') 'int4 = ', &
! SUM((REAL(REAL(int4(:,:,:,:,:))))**2)+SUM((REAL(AIMAG(int4(:,:,:,:,:))))**2)
! write(*,'(a,e20.12)') 'int5 = ', &
! SUM((REAL(REAL(int5(:,:,:,:,:))))**2)+SUM((REAL(AIMAG(int5(:,:,:,:,:))))**2)
!END
!
DEALLOCATE (sk)

View File

@ -175,13 +175,14 @@
!
CALL invfft('Rho', drhoc, dfftp)
!
IF ( .NOT. lsda) THEN
aux(:) = czero
IF (.NOT. lsda) THEN
DO ir = 1, dfftp%nnr
aux(ir) = drhoc(ir) * dmuxc(ir,1,1)
ENDDO
ELSE
is = isk_loc(ik)
DO ir=1, dfftp%nnr
DO ir = 1, dfftp%nnr
aux(ir) = drhoc(ir) * 0.5d0 * (dmuxc(ir, is, 1) + dmuxc(ir, is, 2))
ENDDO
ENDIF

View File

@ -10,7 +10,7 @@
! adapted from PH/dvqpsi_us_only (QE)
!
!----------------------------------------------------------------------
subroutine dvqpsi_us_only3 (ik, uact, xxkq, igkq, npwq)
SUBROUTINE dvqpsi_us_only3 (ik, uact, xxkq, igkq, npwq)
!----------------------------------------------------------------------
!!
!! This routine calculates dV_bare/dtau * psi for one perturbation
@ -35,7 +35,7 @@
USE eqv, ONLY : dvpsi
USE elph2, ONLY : lower_band, upper_band
USE noncollin_module, ONLY : noncolin, npol
USE constants_epw, ONLY : czero, cone, eps12
USE constants_epw, ONLY : czero, zero, cone, eps12
USE klist_epw, ONLY : isk_loc
!
IMPLICIT NONE
@ -100,24 +100,24 @@
ALLOCATE ( ps1_nc(nkb, npol, lower_band:upper_band) )
ALLOCATE ( ps2_nc(nkb, npol, lower_band:upper_band, 3) )
ALLOCATE ( deff_nc(nhm, nhm, nat, nspin) )
ps1_nc(:,:,:) = czero
ps2_nc(:,:,:,:) = czero
deff_nc(:,:,:,:) = czero
ELSE
ALLOCATE ( ps1(nkb, lower_band:upper_band) )
ALLOCATE ( ps2(nkb, lower_band:upper_band, 3) )
ALLOCATE ( deff(nhm, nhm, nat) )
ps1(:,:) = czero
ps2(:,:,:) = czero
deff(:,:,:) = zero
ENDIF
ALLOCATE ( aux(npwx) )
aux(:) = czero
!
IF (lsda) current_spin = isk_loc(ik)
!
! we first compute the coefficients of the vectors
!
IF (noncolin) THEN
ps1_nc(:,:,:) = czero
ps2_nc(:,:,:,:) = czero
ELSE
ps1(:,:) = czero
ps2(:,:,:) = czero
ENDIF
!
DO ibnd = lower_band, upper_band
IF (noncolin) THEN
CALL compute_deff_nc( deff_nc, et(ibnd,ik) )
@ -142,7 +142,7 @@
DO js = 1, npol
ijs = ijs + 1
ps1_nc(ikb,is,ibnd) = ps1_nc(ikb,is,ibnd) + &
deff_nc(ih,jh,na,ijs) * &
deff_nc(ih,jh,na,ijs) * &
alphap(ipol,ik)%nc(jkb,js,ibnd) * uact(mu+ipol)
ps2_nc(ikb,is,ibnd,ipol) = ps2_nc(ikb,is,ibnd,ipol) + &
deff_nc(ih,jh,na,ijs) * becp1(ik)%nc(jkb,js,ibnd) * &
@ -157,52 +157,52 @@
deff(ih,jh,na) * becp1(ik)%k(jkb,ibnd) * &
(0.d0,-1.d0) * uact(mu+ipol) * tpiba
ENDIF
! IF (okvan) THEN
! IF (noncolin) THEN
! ijs = 0
! DO is = 1, npol
! DO js = 1, npol
! ijs = ijs + 1
! ps1_nc(ikb,is,ibnd) = ps1_nc(ikb,is,ibnd) + &
! int1_nc(ih,jh,ipol,na,ijs) * &
! becp1(ik)%nc(jkb,js,ibnd) * uact(mu+ipol)
! ENDDO
! ENDDO
! ELSE
! ps1(ikb,ibnd) = ps1(ikb, ibnd) + &
! int1(ih,jh,ipol,na,current_spin) * &
! becp1(ik)%k(jkb,ibnd) * uact(mu+ipol)
! ENDIF
! ENDIF ! okvan
IF (okvan) THEN
IF (noncolin) THEN
ijs = 0
DO is = 1, npol
DO js = 1, npol
ijs = ijs + 1
ps1_nc(ikb,is,ibnd) = ps1_nc(ikb,is,ibnd) + &
int1_nc(ih,jh,ipol,na,ijs) * &
becp1(ik)%nc(jkb,js,ibnd) * uact(mu+ipol)
ENDDO
ENDDO
ELSE
ps1(ikb,ibnd) = ps1(ikb, ibnd) + &
int1(ih,jh,ipol,na,current_spin) * &
becp1(ik)%k(jkb,ibnd) * uact(mu+ipol)
ENDIF
ENDIF ! okvan
ENDIF ! uact>0
! IF (okvan) THEN
! DO nb = 1, nat
! nu = 3 * (nb - 1)
! IF (noncolin) THEN
! IF (lspinorb) THEN
! ijs = 0
! DO is = 1, npol
! DO js = 1, npol
! ijs = ijs + 1
! ps1_nc(ikb,is,ibnd) = ps1_nc(ikb,is,ibnd) + &
! int2_so(ih,jh,ipol,nb,na,ijs) * &
! becp1(ik)%nc(jkb,js,ibnd) * uact(nu+ipol)
! ENDDO
! ENDDO
! ELSE
! DO is = 1, npol
! ps1_nc(ikb,is,ibnd) = ps1_nc(ikb,is,ibnd) + &
! int2(ih,jh,ipol,nb,na) * &
! becp1(ik)%nc(jkb,is,ibnd) * uact(nu+ipol)
! ENDDO
! ENDIF
! ELSE
! ps1(ikb,ibnd) = ps1(ikb,ibnd) + &
! int2(ih,jh,ipol,nb,na) * &
! becp1(ik)%k(jkb,ibnd) * uact(nu+ipol)
! ENDIF
! ENDDO
! ENDIF ! okvan
IF (okvan) THEN
DO nb = 1, nat
nu = 3 * (nb - 1)
IF (noncolin) THEN
IF (lspinorb) THEN
ijs = 0
DO is = 1, npol
DO js = 1, npol
ijs = ijs + 1
ps1_nc(ikb,is,ibnd) = ps1_nc(ikb,is,ibnd) + &
int2_so(ih,jh,ipol,nb,na,ijs) * &
becp1(ik)%nc(jkb,js,ibnd) * uact(nu+ipol)
ENDDO
ENDDO
ELSE
DO is = 1, npol
ps1_nc(ikb,is,ibnd) = ps1_nc(ikb,is,ibnd) + &
int2(ih,jh,ipol,nb,na) * &
becp1(ik)%nc(jkb,is,ibnd) * uact(nu+ipol)
ENDDO
ENDIF
ELSE
ps1(ikb,ibnd) = ps1(ikb,ibnd) + &
int2(ih,jh,ipol,nb,na) * &
becp1(ik)%k(jkb,ibnd) * uact(nu+ipol)
ENDIF
ENDDO
ENDIF ! okvan
ENDDO ! ipol
ENDDO ! jh
ENDDO ! ih
@ -216,27 +216,27 @@
!
IF (nkb > 0) THEN
IF (noncolin) THEN
CALL zgemm('n', 'n', npwq, (upper_band-lower_band+1)*npol, nkb, &
cone, vkb, npwx, ps1_nc, nkb, cone, dvpsi, npwx)
CALL ZGEMM( 'n', 'n', npwq, (upper_band-lower_band+1)*npol, nkb, &
cone, vkb, npwx, ps1_nc, nkb, cone, dvpsi, npwx )
ELSE
CALL zgemm('n', 'n', npwq, (upper_band-lower_band+1), nkb, &
cone, vkb, npwx, ps1, nkb, cone, dvpsi, npwx)
CALL ZGEMM( 'n', 'n', npwq, (upper_band-lower_band+1), nkb, &
cone, vkb, npwx, ps1, nkb, cone, dvpsi, npwx )
ENDIF
ENDIF
!
! This term is proportional to (k+q+G)_\alpha*beta(k+q+G)
!
DO ikb=1, nkb
DO ipol=1, 3
ok = .false.
DO ikb = 1, nkb
DO ipol = 1, 3
ok = .FALSE.
IF (noncolin) THEN
DO ibnd = lower_band, upper_band
ok = ok .OR. (ABS(ps2_nc(ikb,1,ibnd,ipol) ) > eps12 ) .OR. &
(ABS(ps2_nc(ikb,2,ibnd,ipol) ) > eps12 )
ok = ok .OR. ( ABS(ps2_nc(ikb,1,ibnd,ipol) ) > eps12 ) .OR. &
( ABS(ps2_nc(ikb,2,ibnd,ipol) ) > eps12 )
ENDDO
ELSE
DO ibnd=lower_band, upper_band
ok = ok .OR. (ABS(ps2(ikb,ibnd,ipol) ) > eps12)
DO ibnd = lower_band, upper_band
ok = ok .OR. ( ABS(ps2(ikb,ibnd,ipol) ) > eps12)
ENDDO
ENDIF
IF (ok) THEN
@ -247,10 +247,10 @@
ENDDO
DO ibnd = lower_band, upper_band
IF (noncolin) THEN
CALL zaxpy( npwq, ps2_nc(ikb,1,ibnd,ipol), aux, 1, dvpsi(1,ibnd), 1 )
CALL zaxpy( npwq, ps2_nc(ikb,2,ibnd,ipol), aux, 1, dvpsi(1+npwx,ibnd), 1 )
CALL ZAXPY( npwq, ps2_nc(ikb,1,ibnd,ipol), aux, 1, dvpsi(1,ibnd), 1 )
CALL ZAXPY( npwq, ps2_nc(ikb,2,ibnd,ipol), aux, 1, dvpsi(1+npwx,ibnd), 1 )
ELSE
CALL zaxpy( npwq, ps2(ikb,ibnd,ipol), aux, 1, dvpsi(1,ibnd), 1 )
CALL ZAXPY( npwq, ps2(ikb,ibnd,ipol), aux, 1, dvpsi(1,ibnd), 1 )
ENDIF
ENDDO
ENDIF
@ -259,12 +259,12 @@
!
DEALLOCATE (aux)
IF (noncolin) THEN
DEALLOCATE (ps2_nc)
DEALLOCATE (ps1_nc)
DEALLOCATE (ps2_nc)
DEALLOCATE (deff_nc)
ELSE
DEALLOCATE (ps2)
DEALLOCATE (ps1)
DEALLOCATE (ps2)
DEALLOCATE (deff)
ENDIF
!

View File

@ -56,11 +56,12 @@
USE constants_epw, ONLY : ryd2ev, zero, czero
USE fft_base, ONLY : dfftp
USE control_ph, ONLY : u_from_file
USE noncollin_module, ONLY : m_loc, npol
USE noncollin_module, ONLY : m_loc, npol, noncolin
USE iotk_module, ONLY : iotk_open_read, iotk_scan_dat, iotk_free_unit, &
iotk_close_read
USE division, ONLY : fkbounds
USE uspp, ONLY : okvan
USE spin_orb, ONLY : lspinorb
USE lrus, ONLY : becp1
USE becmod, ONLY : becp, deallocate_bec_type
USE phus, ONLY : int1, int1_nc, int2, int2_so, &
@ -85,7 +86,7 @@
INTEGER :: maxvalue
!! Temporary integer for max value
INTEGER :: nqxq_tmp
!! Maximum G+q length ?
!! Maximum G+q length
INTEGER :: ibnd
!! Band index
INTEGER :: ik
@ -230,9 +231,13 @@
nqxq_tmp = INT(((SQRT(gcutm) + qnorm_tmp) / dq + 4) * cell_factor)
IF (nqxq_tmp > maxvalue) maxvalue = nqxq_tmp
ENDDO
!
IF (maxvalue > nqxq) THEN
IF (ALLOCATED(qrad)) DEALLOCATE (qrad)
ALLOCATE (qrad(maxvalue, nbetam * (nbetam + 1) / 2, lmaxq, nsp))
qrad(:,:,:,:) = zero
! RM - need to call init_us_1 to re-calculate qrad
CALL init_us_1
ENDIF
!
! do not perform the check if restart
@ -321,11 +326,11 @@
!
dynq(:, :, :) = czero
epmatq(:, :, :, :, :) = czero
epsi(:, :) = zero
zstar(:, :, :) = zero
bmat(:, :, :, :) = czero
cu(:, :, :) = czero
cuq(:, :, :) = czero
epsi(:, :) = zero
zstar(:, :, :) = zero
!
! read interatomic force constat matrix from q2r
IF (lifc) CALL read_ifc
@ -473,7 +478,7 @@
DO iq = 1, nq
! SP: First the vlocq needs to be initialized properly with the first
! q in the star
xq = xq0
xq = xq0
CALL epw_init(.false.)
!
! retrieve the q in the star
@ -705,18 +710,22 @@
DEALLOCATE (int2)
DEALLOCATE (int4)
DEALLOCATE (int5)
DEALLOCATE (int1_nc)
DEALLOCATE (int4_nc)
DEALLOCATE (int2_so)
DEALLOCATE (int5_so)
IF (noncolin) THEN
DEALLOCATE (int1_nc)
DEALLOCATE (int4_nc)
IF (lspinorb) THEN
DEALLOCATE (int2_so)
DEALLOCATE (int5_so)
ENDIF
ENDIF
ENDIF
DO ik=1, nks
DO ipol=1, 3
DO ik = 1, nks
DO ipol = 1, 3
CALL deallocate_bec_type( alphap(ipol,ik) )
ENDDO
ENDDO
DEALLOCATE (alphap)
DO ik=1, size(becp1)
DO ik = 1, size(becp1)
CALL deallocate_bec_type( becp1(ik) )
ENDDO
DEALLOCATE (becp1)

View File

@ -76,7 +76,7 @@
USE division, ONLY : fkbounds
USE mp, ONLY : mp_barrier, mp_bcast, mp_sum
USE io_global, ONLY : ionode_id
USE mp_global, ONLY : inter_pool_comm, intra_pool_comm, root_pool
USE mp_global, ONLY : inter_pool_comm
USE mp_world, ONLY : mpime, world_comm
#if defined(__MPI)
USE parallel_include, ONLY : MPI_MODE_RDONLY, MPI_INFO_NULL, MPI_OFFSET_KIND, &
@ -1885,7 +1885,7 @@
#endif
USE io_global, ONLY : ionode_id
USE mp, ONLY : mp_barrier, mp_bcast
USE mp_global, ONLY : intra_pool_comm, inter_pool_comm, root_pool, world_comm
USE mp_global, ONLY : inter_pool_comm, world_comm
USE mp_world, ONLY : mpime
!
implicit none

View File

@ -77,7 +77,7 @@
USE division, ONLY : fkbounds
USE mp, ONLY : mp_barrier, mp_bcast, mp_sum
USE io_global, ONLY : ionode_id
USE mp_global, ONLY : inter_pool_comm, intra_pool_comm, root_pool
USE mp_global, ONLY : inter_pool_comm
USE mp_world, ONLY : mpime, world_comm
#if defined(__MPI)
USE parallel_include, ONLY : MPI_MODE_RDONLY, MPI_INFO_NULL, MPI_OFFSET_KIND, &

View File

@ -20,9 +20,9 @@
!
USE kinds, ONLY : DP
USE ions_base, ONLY : nat, ntyp => nsp, tau
USE becmod, ONLY : calbec
USE becmod, ONLY : calbec, becp, allocate_bec_type
USE lrus, ONLY : becp1
USE uspp, ONLY : vkb
USE uspp, ONLY : vkb, nlcc_any, okvan, nkb
USE pwcom, ONLY : npwx, nbnd, nks
USE klist_epw, ONLY : xk_loc, isk_loc
USE constants, ONLY : tpi
@ -33,21 +33,17 @@
USE atom, ONLY : msh, rgrid
USE wavefunctions, ONLY : evc
USE noncollin_module, ONLY : noncolin, npol, nspin_mag
USE uspp_param, ONLY : upf
USE uspp_param, ONLY : upf, nhm
USE m_gth, ONLY : setlocq_gth
USE units_lr, ONLY : lrwfc, iuwfc
USE phcom, ONLY : vlocq
USE qpoint, ONLY : xq, eigqts
USE nlcc_ph, ONLY : drc
USE uspp, ONLY : nlcc_any
USE elph2, ONLY : igk_k_all, ngk_all
USE mp, ONLY : mp_barrier
USE mp_global, ONLY : inter_pool_comm, my_pool_id
USE spin_orb, ONLY : lspinorb
USE uspp_param, ONLY : nhm
USE uspp, ONLY : okvan, nkb
USE lsda_mod, ONLY : nspin, lsda, current_spin
USE becmod, ONLY : becp, allocate_bec_type
USE phus, ONLY : int1, int1_nc, int2, int2_so, &
int4, int4_nc, int5, int5_so, &
alphap
@ -110,7 +106,7 @@
CALL allocate_bec_type(nkb, nbnd, becp)
ENDIF
!
DO na=1, nat
DO na = 1, nat
!
! xq here is the first q of the star
arg = (xq(1) * tau(1, na) + &
@ -144,7 +140,7 @@
ALLOCATE (aux1(npwx*npol, nbnd))
!ALLOCATE (evc(npwx*npol, nbnd))
!
DO ik=1, nks
DO ik = 1, nks
!
!
IF (lsda) current_spin = isk_loc(ik)
@ -180,16 +176,15 @@
ENDDO
ENDIF
ENDDO
CALL calbec( ngk(ik), vkb, aux1, alphap(ipol,ik) )
CALL calbec( ngk(ik), vkb, aux1, alphap(ipol,ik) )
ENDDO
!
!
ENDDO
!
DEALLOCATE (aux1)
!
IF( .NOT. ALLOCATED(igk_k_all)) ALLOCATE (igk_k_all(npwx,nkstot))
IF( .NOT. ALLOCATED(ngk_all)) ALLOCATE (ngk_all(nkstot))
IF( .NOT. ALLOCATED(igk_k_all) ) ALLOCATE (igk_k_all(npwx,nkstot))
IF( .NOT. ALLOCATED(ngk_all) ) ALLOCATE (ngk_all(nkstot))
!
#if defined(__MPI)
!
@ -204,7 +199,7 @@
!
#endif
!
IF ( .NOT. first_run) CALL dvanqq2()
IF ( .NOT. first_run ) CALL dvanqq2()
!
CALL stop_clock( 'epw_init' )
!

View File

@ -73,7 +73,6 @@
USE io_global, ONLY : meta_ionode, meta_ionode_id, ionode_id
USE io_epw, ONLY : iunkf, iunqf
USE noncollin_module, ONLY : npol
USE wavefunctions, ONLY : evc
USE wvfct, ONLY : npwx
#if defined(__NAG)
USE F90_UNIX_ENV, ONLY : iargc, getarg

View File

@ -33,7 +33,7 @@
USE mp, ONLY : mp_sum
USE lrus, ONLY : int3
USE qpoint, ONLY : eigqts
USE constants_epw, ONLY : czero
USE constants_epw, ONLY : czero, zero
!
IMPLICIT NONE
!
@ -72,15 +72,15 @@
REAL(kind=DP), ALLOCATABLE :: ylmk0(:,:)
!! the spherical harmonics at q+G
!
COMPLEX(kind=DP), EXTERNAL :: zdotc
COMPLEX(kind=DP), EXTERNAL :: ZDOTC
!! the scalar product function
COMPLEX(kind=DP), ALLOCATABLE :: aux1(:), aux2(:,:)
COMPLEX(kind=DP), ALLOCATABLE :: qgm(:)
!! the augmentation function at G
!! the augmentation function at q+G
COMPLEX(kind=DP), ALLOCATABLE :: veff(:)
!! effective potential
!
IF ( .NOT. okvan) RETURN
IF (.NOT. okvan) RETURN
!
CALL start_clock('newdq2')
!
@ -93,11 +93,19 @@
ALLOCATE ( qgm(ngm) )
ALLOCATE ( qmod(ngm) )
ALLOCATE ( qg(3,ngm) )
aux1(:) = czero
aux2(:,:) = czero
veff(:) = czero
ylmk0(:,:) = zero
qgm(:) = czero
qmod(:) = zero
qg(:,:) = zero
!
! first compute the spherical harmonics
!
CALL setqmod( ngm, xq0, g, qmod, qg )
CALL ylmr2( lmaxq * lmaxq, ngm, qg, qmod, ylmk0 )
!
DO ig = 1, ngm
qmod(ig) = sqrt( qmod(ig) )
ENDDO
@ -127,6 +135,7 @@
!
DO ih = 1, nh(nt)
DO jh = ih, nh(nt)
!
CALL qvan2( ngm, ih, jh, nt, qmod, qgm, ylmk0 )
!
DO na = 1, nat
@ -139,7 +148,7 @@
ENDDO
DO is = 1, nspin_mag
int3(ih,jh,na,is,ipert) = omega * &
zdotc(ngm,aux1,1,aux2(1,is),1)
ZDOTC(ngm,aux1,1,aux2(1,is),1)
ENDDO
ENDIF
ENDDO

View File

@ -1567,8 +1567,8 @@
!
any_uspp = ANY( upf(:)%tvanp )
!
IF (any_uspp .and. noncolin) CALL errore('pw2wan90epw',&
'noncolin calculation not implimented with USP',1)
IF ( any_uspp ) CALL errore('pw2wan90epw',&
'dipole matrix calculation not implimented with USP - set vme=.true.',1)
!
ALLOCATE (dmec(3, nbnd, nbnd, nks))
!

View File

@ -27,7 +27,7 @@
Comput. Phys. Commun. 209, 116 (2016)
Program EPW v.5.0.0 starts on 18Jan2019 at 19: 6:16
Program EPW v.5.1.0 starts on 17Jul2019 at 12:13:57
This program is part of the open-source Quantum ESPRESSO suite
for quantum simulation of materials; please cite
@ -45,7 +45,8 @@
./sic.save/
IMPORTANT: XC functional enforced from input :
Exchange-correlation = PBE ( 1 4 3 4 0 0)
Exchange-correlation= PBE
( 1 4 3 4 0 0 0)
Any further DFT definition will be discarded
Please, verify this is what you really want
@ -66,7 +67,8 @@
number of atomic types = 2
kinetic-energy cut-off = 30.0000 Ry
charge density cut-off = 120.0000 Ry
Exchange-correlation = PBE ( 1 4 3 4 0 0)
Exchange-correlation= PBE
( 1 4 3 4 0 0 0)
celldm(1)= 8.23700 celldm(2)= 0.00000 celldm(3)= 0.00000
@ -147,9 +149,9 @@
l(4) = 1
Q(r) pseudized with 0 coefficients
EPW : 0.26s CPU 0.26s WALL
EPW : 0.27s CPU 0.28s WALL
EPW : 0.31s CPU 0.32s WALL
EPW : 0.32s CPU 0.33s WALL
No wavefunction gauge setting applied
-------------------------------------------------------------------
@ -251,16 +253,13 @@
( -0.16235 -0.16235 0.16235) : 0.87049
-------------------------------------------------------------------
WANNIER : 2.04s CPU 2.05s WALL ( 1 calls)
WANNIER : 1.99s CPU 1.99s WALL ( 1 calls)
-------------------------------------------------------------------
Dipole matrix elements calculated
Calculating kgmap
Progress kgmap: ########################################
kmaps : 0.22s CPU 0.22s WALL ( 1 calls)
kmaps : 0.20s CPU 0.21s WALL ( 1 calls)
Symmetries of Bravais lattice: 48
Symmetries of crystal: 24
@ -372,23 +371,23 @@
The .epb files have been correctly written
Computes the analytic long-range interaction for polar materials [lpolar]
Construct the Wigner-Seitz cell using Wannier centers and atomic positions
Number of WS vectors for electrons 79
Number of WS vectors for phonons 63
Number of WS vectors for electron-phonon 53
Maximum number of cores for efficient parallelization 106
Velocity matrix elements calculated
Writing Hamiltonian, Dynamical matrix and EP vertex in Wann rep to file
Reading Hamiltonian, Dynamical matrix and EP vertex in Wann rep from file
Finished reading Wann rep data from file
===================================================================
Memory usage: VmHWM = 31Mb
VmPeak = 301Mb
Memory usage: VmHWM = 30Mb
VmPeak = 302Mb
===================================================================
Using q-mesh file: pathq.dat
@ -419,60 +418,60 @@
ik = 1 coord.: 0.0000000 0.0000000 0.0000000
ibnd jbnd imode enk[eV] enk+q[eV] omega(q)[meV] |g|[meV]
------------------------------------------------------------------------------
2 2 1 9.5362 9.5362 0.0000 0.3356727864E+02
2 2 2 9.5362 9.5362 0.0000 0.3356727864E+02
2 2 3 9.5362 9.5362 0.0000 0.3356727864E+02
2 2 4 9.5362 9.5362 96.9696 0.1152036546E+03
2 2 5 9.5362 9.5362 96.9696 0.1152036546E+03
2 2 6 9.5362 9.5362 96.9696 0.1152036546E+03
2 3 1 9.5362 9.5362 0.0000 0.3356727864E+02
2 3 2 9.5362 9.5362 0.0000 0.3356727864E+02
2 3 3 9.5362 9.5362 0.0000 0.3356727864E+02
2 3 4 9.5362 9.5362 96.9696 0.1152036546E+03
2 3 5 9.5362 9.5362 96.9696 0.1152036546E+03
2 3 6 9.5362 9.5362 96.9696 0.1152036546E+03
2 4 1 9.5362 9.5362 0.0000 0.3356727864E+02
2 4 2 9.5362 9.5362 0.0000 0.3356727864E+02
2 4 3 9.5362 9.5362 0.0000 0.3356727864E+02
2 4 4 9.5362 9.5362 96.9696 0.1152036546E+03
2 4 5 9.5362 9.5362 96.9696 0.1152036546E+03
2 4 6 9.5362 9.5362 96.9696 0.1152036546E+03
3 2 1 9.5362 9.5362 0.0000 0.3356727864E+02
3 2 2 9.5362 9.5362 0.0000 0.3356727864E+02
3 2 3 9.5362 9.5362 0.0000 0.3356727864E+02
3 2 4 9.5362 9.5362 96.9696 0.1152036546E+03
3 2 5 9.5362 9.5362 96.9696 0.1152036546E+03
3 2 6 9.5362 9.5362 96.9696 0.1152036546E+03
3 3 1 9.5362 9.5362 0.0000 0.3356727864E+02
3 3 2 9.5362 9.5362 0.0000 0.3356727864E+02
3 3 3 9.5362 9.5362 0.0000 0.3356727864E+02
3 3 4 9.5362 9.5362 96.9696 0.1152036546E+03
3 3 5 9.5362 9.5362 96.9696 0.1152036546E+03
3 3 6 9.5362 9.5362 96.9696 0.1152036546E+03
3 4 1 9.5362 9.5362 0.0000 0.3356727864E+02
3 4 2 9.5362 9.5362 0.0000 0.3356727864E+02
3 4 3 9.5362 9.5362 0.0000 0.3356727864E+02
3 4 4 9.5362 9.5362 96.9696 0.1152036546E+03
3 4 5 9.5362 9.5362 96.9696 0.1152036546E+03
3 4 6 9.5362 9.5362 96.9696 0.1152036546E+03
4 2 1 9.5362 9.5362 0.0000 0.3356727864E+02
4 2 2 9.5362 9.5362 0.0000 0.3356727864E+02
4 2 3 9.5362 9.5362 0.0000 0.3356727864E+02
4 2 4 9.5362 9.5362 96.9696 0.1152036546E+03
4 2 5 9.5362 9.5362 96.9696 0.1152036546E+03
4 2 6 9.5362 9.5362 96.9696 0.1152036546E+03
4 3 1 9.5362 9.5362 0.0000 0.3356727864E+02
4 3 2 9.5362 9.5362 0.0000 0.3356727864E+02
4 3 3 9.5362 9.5362 0.0000 0.3356727864E+02
4 3 4 9.5362 9.5362 96.9696 0.1152036546E+03
4 3 5 9.5362 9.5362 96.9696 0.1152036546E+03
4 3 6 9.5362 9.5362 96.9696 0.1152036546E+03
4 4 1 9.5362 9.5362 0.0000 0.3356727864E+02
4 4 2 9.5362 9.5362 0.0000 0.3356727864E+02
4 4 3 9.5362 9.5362 0.0000 0.3356727864E+02
4 4 4 9.5362 9.5362 96.9696 0.1152036546E+03
4 4 5 9.5362 9.5362 96.9696 0.1152036546E+03
4 4 6 9.5362 9.5362 96.9696 0.1152036546E+03
2 2 1 9.5362 9.5362 0.0000 0.3519101343E+02
2 2 2 9.5362 9.5362 0.0000 0.3519101343E+02
2 2 3 9.5362 9.5362 0.0000 0.3519101343E+02
2 2 4 9.5362 9.5362 96.9696 0.1152035784E+03
2 2 5 9.5362 9.5362 96.9696 0.1152035784E+03
2 2 6 9.5362 9.5362 96.9696 0.1152035784E+03
2 3 1 9.5362 9.5362 0.0000 0.3519101343E+02
2 3 2 9.5362 9.5362 0.0000 0.3519101343E+02
2 3 3 9.5362 9.5362 0.0000 0.3519101343E+02
2 3 4 9.5362 9.5362 96.9696 0.1152035784E+03
2 3 5 9.5362 9.5362 96.9696 0.1152035784E+03
2 3 6 9.5362 9.5362 96.9696 0.1152035784E+03
2 4 1 9.5362 9.5362 0.0000 0.3519101343E+02
2 4 2 9.5362 9.5362 0.0000 0.3519101343E+02
2 4 3 9.5362 9.5362 0.0000 0.3519101343E+02
2 4 4 9.5362 9.5362 96.9696 0.1152035784E+03
2 4 5 9.5362 9.5362 96.9696 0.1152035784E+03
2 4 6 9.5362 9.5362 96.9696 0.1152035784E+03
3 2 1 9.5362 9.5362 0.0000 0.3519101343E+02
3 2 2 9.5362 9.5362 0.0000 0.3519101343E+02
3 2 3 9.5362 9.5362 0.0000 0.3519101343E+02
3 2 4 9.5362 9.5362 96.9696 0.1152035784E+03
3 2 5 9.5362 9.5362 96.9696 0.1152035784E+03
3 2 6 9.5362 9.5362 96.9696 0.1152035784E+03
3 3 1 9.5362 9.5362 0.0000 0.3519101343E+02
3 3 2 9.5362 9.5362 0.0000 0.3519101343E+02
3 3 3 9.5362 9.5362 0.0000 0.3519101343E+02
3 3 4 9.5362 9.5362 96.9696 0.1152035784E+03
3 3 5 9.5362 9.5362 96.9696 0.1152035784E+03
3 3 6 9.5362 9.5362 96.9696 0.1152035784E+03
3 4 1 9.5362 9.5362 0.0000 0.3519101343E+02
3 4 2 9.5362 9.5362 0.0000 0.3519101343E+02
3 4 3 9.5362 9.5362 0.0000 0.3519101343E+02
3 4 4 9.5362 9.5362 96.9696 0.1152035784E+03
3 4 5 9.5362 9.5362 96.9696 0.1152035784E+03
3 4 6 9.5362 9.5362 96.9696 0.1152035784E+03
4 2 1 9.5362 9.5362 0.0000 0.3519101343E+02
4 2 2 9.5362 9.5362 0.0000 0.3519101343E+02
4 2 3 9.5362 9.5362 0.0000 0.3519101343E+02
4 2 4 9.5362 9.5362 96.9696 0.1152035784E+03
4 2 5 9.5362 9.5362 96.9696 0.1152035784E+03
4 2 6 9.5362 9.5362 96.9696 0.1152035784E+03
4 3 1 9.5362 9.5362 0.0000 0.3519101343E+02
4 3 2 9.5362 9.5362 0.0000 0.3519101343E+02
4 3 3 9.5362 9.5362 0.0000 0.3519101343E+02
4 3 4 9.5362 9.5362 96.9696 0.1152035784E+03
4 3 5 9.5362 9.5362 96.9696 0.1152035784E+03
4 3 6 9.5362 9.5362 96.9696 0.1152035784E+03
4 4 1 9.5362 9.5362 0.0000 0.3519101343E+02
4 4 2 9.5362 9.5362 0.0000 0.3519101343E+02
4 4 3 9.5362 9.5362 0.0000 0.3519101343E+02
4 4 4 9.5362 9.5362 96.9696 0.1152035784E+03
4 4 5 9.5362 9.5362 96.9696 0.1152035784E+03
4 4 6 9.5362 9.5362 96.9696 0.1152035784E+03
------------------------------------------------------------------------------
Electron-phonon vertex |g| (meV)
@ -481,60 +480,60 @@
ik = 1 coord.: 0.0000000 0.0000000 0.0000000
ibnd jbnd imode enk[eV] enk+q[eV] omega(q)[meV] |g|[meV]
------------------------------------------------------------------------------
2 2 1 9.5362 3.3767 29.0464 0.3098790956E+02
2 2 2 9.5362 3.3767 29.0464 0.3098790956E+02
2 2 3 9.5362 3.3767 62.4896 0.9280967135E+02
2 2 4 9.5362 3.3767 94.2352 0.7542390365E+02
2 2 5 9.5362 3.3767 94.2352 0.7542390365E+02
2 2 6 9.5362 3.3767 109.2915 0.1251759224E+03
2 3 1 9.5362 8.6927 29.0464 0.2601275765E+02
2 3 2 9.5362 8.6927 29.0464 0.2601275765E+02
2 3 3 9.5362 8.6927 62.4896 0.1258616594E+03
2 3 4 9.5362 8.6927 94.2352 0.1244187173E+03
2 3 5 9.5362 8.6927 94.2352 0.1244187173E+03
2 3 6 9.5362 8.6927 109.2915 0.1464359089E+03
2 4 1 9.5362 8.6927 29.0464 0.2601275765E+02
2 4 2 9.5362 8.6927 29.0464 0.2601275765E+02
2 4 3 9.5362 8.6927 62.4896 0.1258616594E+03
2 4 4 9.5362 8.6927 94.2352 0.1244187173E+03
2 4 5 9.5362 8.6927 94.2352 0.1244187173E+03
2 4 6 9.5362 8.6927 109.2915 0.1464359089E+03
3 2 1 9.5362 3.3767 29.0464 0.3098790956E+02
3 2 2 9.5362 3.3767 29.0464 0.3098790956E+02
3 2 3 9.5362 3.3767 62.4896 0.9280967135E+02
3 2 4 9.5362 3.3767 94.2352 0.7542390365E+02
3 2 5 9.5362 3.3767 94.2352 0.7542390365E+02
3 2 6 9.5362 3.3767 109.2915 0.1251759224E+03
3 3 1 9.5362 8.6927 29.0464 0.2601275765E+02
3 3 2 9.5362 8.6927 29.0464 0.2601275765E+02
3 3 3 9.5362 8.6927 62.4896 0.1258616594E+03
3 3 4 9.5362 8.6927 94.2352 0.1244187173E+03
3 3 5 9.5362 8.6927 94.2352 0.1244187173E+03
3 3 6 9.5362 8.6927 109.2915 0.1464359089E+03
3 4 1 9.5362 8.6927 29.0464 0.2601275765E+02
3 4 2 9.5362 8.6927 29.0464 0.2601275765E+02
3 4 3 9.5362 8.6927 62.4896 0.1258616594E+03
3 4 4 9.5362 8.6927 94.2352 0.1244187173E+03
3 4 5 9.5362 8.6927 94.2352 0.1244187173E+03
3 4 6 9.5362 8.6927 109.2915 0.1464359089E+03
4 2 1 9.5362 3.3767 29.0464 0.3098790956E+02
4 2 2 9.5362 3.3767 29.0464 0.3098790956E+02
4 2 3 9.5362 3.3767 62.4896 0.9280967135E+02
4 2 4 9.5362 3.3767 94.2352 0.7542390365E+02
4 2 5 9.5362 3.3767 94.2352 0.7542390365E+02
4 2 6 9.5362 3.3767 109.2915 0.1251759224E+03
4 3 1 9.5362 8.6927 29.0464 0.2601275765E+02
4 3 2 9.5362 8.6927 29.0464 0.2601275765E+02
4 3 3 9.5362 8.6927 62.4896 0.1258616594E+03
4 3 4 9.5362 8.6927 94.2352 0.1244187173E+03
4 3 5 9.5362 8.6927 94.2352 0.1244187173E+03
4 3 6 9.5362 8.6927 109.2915 0.1464359089E+03
4 4 1 9.5362 8.6927 29.0464 0.2601275765E+02
4 4 2 9.5362 8.6927 29.0464 0.2601275765E+02
4 4 3 9.5362 8.6927 62.4896 0.1258616594E+03
4 4 4 9.5362 8.6927 94.2352 0.1244187173E+03
4 4 5 9.5362 8.6927 94.2352 0.1244187173E+03
4 4 6 9.5362 8.6927 109.2915 0.1464359089E+03
2 2 1 9.5362 3.3767 29.0461 0.2853406535E+02
2 2 2 9.5362 3.3767 29.0461 0.2853406535E+02
2 2 3 9.5362 3.3767 62.4893 0.1121901696E+03
2 2 4 9.5362 3.3767 94.2351 0.1032444000E+03
2 2 5 9.5362 3.3767 94.2351 0.1032444000E+03
2 2 6 9.5362 3.3767 109.2914 0.1372711279E+03
2 3 1 9.5362 8.6927 29.0461 0.2853406535E+02
2 3 2 9.5362 8.6927 29.0461 0.2853406535E+02
2 3 3 9.5362 8.6927 62.4893 0.1121901696E+03
2 3 4 9.5362 8.6927 94.2351 0.1032444000E+03
2 3 5 9.5362 8.6927 94.2351 0.1032444000E+03
2 3 6 9.5362 8.6927 109.2914 0.1372711279E+03
2 4 1 9.5362 8.6927 29.0461 0.2853406535E+02
2 4 2 9.5362 8.6927 29.0461 0.2853406535E+02
2 4 3 9.5362 8.6927 62.4893 0.1121901696E+03
2 4 4 9.5362 8.6927 94.2351 0.1032444000E+03
2 4 5 9.5362 8.6927 94.2351 0.1032444000E+03
2 4 6 9.5362 8.6927 109.2914 0.1372711279E+03
3 2 1 9.5362 3.3767 29.0461 0.2738064309E+02
3 2 2 9.5362 3.3767 29.0461 0.2738064309E+02
3 2 3 9.5362 3.3767 62.4893 0.1177059697E+03
3 2 4 9.5362 3.3767 94.2351 0.1139944117E+03
3 2 5 9.5362 3.3767 94.2351 0.1139944117E+03
3 2 6 9.5362 3.3767 109.2914 0.1409130214E+03
3 3 1 9.5362 8.6927 29.0461 0.2738064309E+02
3 3 2 9.5362 8.6927 29.0461 0.2738064309E+02
3 3 3 9.5362 8.6927 62.4893 0.1177059697E+03
3 3 4 9.5362 8.6927 94.2351 0.1139944117E+03
3 3 5 9.5362 8.6927 94.2351 0.1139944117E+03
3 3 6 9.5362 8.6927 109.2914 0.1409130214E+03
3 4 1 9.5362 8.6927 29.0461 0.2738064309E+02
3 4 2 9.5362 8.6927 29.0461 0.2738064309E+02
3 4 3 9.5362 8.6927 62.4893 0.1177059697E+03
3 4 4 9.5362 8.6927 94.2351 0.1139944117E+03
3 4 5 9.5362 8.6927 94.2351 0.1139944117E+03
3 4 6 9.5362 8.6927 109.2914 0.1409130214E+03
4 2 1 9.5362 3.3767 29.0461 0.2738064309E+02
4 2 2 9.5362 3.3767 29.0461 0.2738064309E+02
4 2 3 9.5362 3.3767 62.4893 0.1177059697E+03
4 2 4 9.5362 3.3767 94.2351 0.1139944117E+03
4 2 5 9.5362 3.3767 94.2351 0.1139944117E+03
4 2 6 9.5362 3.3767 109.2914 0.1409130214E+03
4 3 1 9.5362 8.6927 29.0461 0.2738064309E+02
4 3 2 9.5362 8.6927 29.0461 0.2738064309E+02
4 3 3 9.5362 8.6927 62.4893 0.1177059697E+03
4 3 4 9.5362 8.6927 94.2351 0.1139944117E+03
4 3 5 9.5362 8.6927 94.2351 0.1139944117E+03
4 3 6 9.5362 8.6927 109.2914 0.1409130214E+03
4 4 1 9.5362 8.6927 29.0461 0.2738064309E+02
4 4 2 9.5362 8.6927 29.0461 0.2738064309E+02
4 4 3 9.5362 8.6927 62.4893 0.1177059697E+03
4 4 4 9.5362 8.6927 94.2351 0.1139944117E+03
4 4 5 9.5362 8.6927 94.2351 0.1139944117E+03
4 4 6 9.5362 8.6927 109.2914 0.1409130214E+03
------------------------------------------------------------------------------
Electron-phonon vertex |g| (meV)
@ -543,60 +542,60 @@
ik = 1 coord.: 0.0000000 0.0000000 0.0000000
ibnd jbnd imode enk[eV] enk+q[eV] omega(q)[meV] |g|[meV]
------------------------------------------------------------------------------
2 2 1 9.5362 4.5130 40.9084 0.4194847423E+02
2 2 2 9.5362 4.5130 40.9084 0.4194847423E+02
2 2 3 9.5362 4.5130 62.8687 0.3170165291E+02
2 2 4 9.5362 4.5130 93.7701 0.1197897914E+03
2 2 5 9.5362 4.5130 93.7701 0.1197897914E+03
2 2 6 9.5362 4.5130 109.2933 0.5188756740E+02
2 3 1 9.5362 6.9263 40.9084 0.1141115167E+02
2 3 2 9.5362 6.9263 40.9084 0.1141115167E+02
2 3 3 9.5362 6.9263 62.8687 0.1314917363E+03
2 3 4 9.5362 6.9263 93.7701 0.7353510442E+02
2 3 5 9.5362 6.9263 93.7701 0.7353510442E+02
2 3 6 9.5362 6.9263 109.2933 0.1643855028E+03
2 4 1 9.5362 6.9263 40.9084 0.1141115167E+02
2 4 2 9.5362 6.9263 40.9084 0.1141115167E+02
2 4 3 9.5362 6.9263 62.8687 0.1314917363E+03
2 4 4 9.5362 6.9263 93.7701 0.7353510442E+02
2 4 5 9.5362 6.9263 93.7701 0.7353510442E+02
2 4 6 9.5362 6.9263 109.2933 0.1643855028E+03
3 2 1 9.5362 4.5130 40.9084 0.4194847423E+02
3 2 2 9.5362 4.5130 40.9084 0.4194847423E+02
3 2 3 9.5362 4.5130 62.8687 0.3170165291E+02
3 2 4 9.5362 4.5130 93.7701 0.1197897914E+03
3 2 5 9.5362 4.5130 93.7701 0.1197897914E+03
3 2 6 9.5362 4.5130 109.2933 0.5188756740E+02
3 3 1 9.5362 6.9263 40.9084 0.1141115167E+02
3 3 2 9.5362 6.9263 40.9084 0.1141115167E+02
3 3 3 9.5362 6.9263 62.8687 0.1314917363E+03
3 3 4 9.5362 6.9263 93.7701 0.7353510442E+02
3 3 5 9.5362 6.9263 93.7701 0.7353510442E+02
3 3 6 9.5362 6.9263 109.2933 0.1643855028E+03
3 4 1 9.5362 6.9263 40.9084 0.1141115167E+02
3 4 2 9.5362 6.9263 40.9084 0.1141115167E+02
3 4 3 9.5362 6.9263 62.8687 0.1314917363E+03
3 4 4 9.5362 6.9263 93.7701 0.7353510442E+02
3 4 5 9.5362 6.9263 93.7701 0.7353510442E+02
3 4 6 9.5362 6.9263 109.2933 0.1643855028E+03
4 2 1 9.5362 4.5130 40.9084 0.4194847423E+02
4 2 2 9.5362 4.5130 40.9084 0.4194847423E+02
4 2 3 9.5362 4.5130 62.8687 0.3170165291E+02
4 2 4 9.5362 4.5130 93.7701 0.1197897914E+03
4 2 5 9.5362 4.5130 93.7701 0.1197897914E+03
4 2 6 9.5362 4.5130 109.2933 0.5188756740E+02
4 3 1 9.5362 6.9263 40.9084 0.1141115167E+02
4 3 2 9.5362 6.9263 40.9084 0.1141115167E+02
4 3 3 9.5362 6.9263 62.8687 0.1314917363E+03
4 3 4 9.5362 6.9263 93.7701 0.7353510442E+02
4 3 5 9.5362 6.9263 93.7701 0.7353510442E+02
4 3 6 9.5362 6.9263 109.2933 0.1643855028E+03
4 4 1 9.5362 6.9263 40.9084 0.1141115167E+02
4 4 2 9.5362 6.9263 40.9084 0.1141115167E+02
4 4 3 9.5362 6.9263 62.8687 0.1314917363E+03
4 4 4 9.5362 6.9263 93.7701 0.7353510442E+02
4 4 5 9.5362 6.9263 93.7701 0.7353510442E+02
4 4 6 9.5362 6.9263 109.2933 0.1643855028E+03
2 2 1 9.5362 4.5130 40.9082 0.2477307077E+02
2 2 2 9.5362 4.5130 40.9082 0.2477307077E+02
2 2 3 9.5362 4.5130 62.8657 0.1014416910E+03
2 2 4 9.5362 4.5130 93.7699 0.9348130446E+02
2 2 5 9.5362 4.5130 93.7699 0.9348130446E+02
2 2 6 9.5362 4.5130 109.2888 0.1287308185E+03
2 3 1 9.5362 6.9263 40.9082 0.2477307077E+02
2 3 2 9.5362 6.9263 40.9082 0.2477307077E+02
2 3 3 9.5362 6.9263 62.8657 0.1014416910E+03
2 3 4 9.5362 6.9263 93.7699 0.9348130446E+02
2 3 5 9.5362 6.9263 93.7699 0.9348130446E+02
2 3 6 9.5362 6.9263 109.2888 0.1287308185E+03
2 4 1 9.5362 6.9263 40.9082 0.2477307077E+02
2 4 2 9.5362 6.9263 40.9082 0.2477307077E+02
2 4 3 9.5362 6.9263 62.8657 0.1014416910E+03
2 4 4 9.5362 6.9263 93.7699 0.9348130446E+02
2 4 5 9.5362 6.9263 93.7699 0.9348130446E+02
2 4 6 9.5362 6.9263 109.2888 0.1287308185E+03
3 2 1 9.5362 4.5130 40.9082 0.2651779195E+02
3 2 2 9.5362 4.5130 40.9082 0.2651779195E+02
3 2 3 9.5362 4.5130 62.8657 0.1124393853E+03
3 2 4 9.5362 4.5130 93.7699 0.9062483731E+02
3 2 5 9.5362 4.5130 93.7699 0.9062483731E+02
3 2 6 9.5362 4.5130 109.2888 0.1417295666E+03
3 3 1 9.5362 6.9263 40.9082 0.2651779195E+02
3 3 2 9.5362 6.9263 40.9082 0.2651779195E+02
3 3 3 9.5362 6.9263 62.8657 0.1124393853E+03
3 3 4 9.5362 6.9263 93.7699 0.9062483731E+02
3 3 5 9.5362 6.9263 93.7699 0.9062483731E+02
3 3 6 9.5362 6.9263 109.2888 0.1417295666E+03
3 4 1 9.5362 6.9263 40.9082 0.2651779195E+02
3 4 2 9.5362 6.9263 40.9082 0.2651779195E+02
3 4 3 9.5362 6.9263 62.8657 0.1124393853E+03
3 4 4 9.5362 6.9263 93.7699 0.9062483731E+02
3 4 5 9.5362 6.9263 93.7699 0.9062483731E+02
3 4 6 9.5362 6.9263 109.2888 0.1417295666E+03
4 2 1 9.5362 4.5130 40.9082 0.2651779195E+02
4 2 2 9.5362 4.5130 40.9082 0.2651779195E+02
4 2 3 9.5362 4.5130 62.8657 0.1124393853E+03
4 2 4 9.5362 4.5130 93.7699 0.9062483731E+02
4 2 5 9.5362 4.5130 93.7699 0.9062483731E+02
4 2 6 9.5362 4.5130 109.2888 0.1417295666E+03
4 3 1 9.5362 6.9263 40.9082 0.2651779195E+02
4 3 2 9.5362 6.9263 40.9082 0.2651779195E+02
4 3 3 9.5362 6.9263 62.8657 0.1124393853E+03
4 3 4 9.5362 6.9263 93.7699 0.9062483731E+02
4 3 5 9.5362 6.9263 93.7699 0.9062483731E+02
4 3 6 9.5362 6.9263 109.2888 0.1417295666E+03
4 4 1 9.5362 6.9263 40.9082 0.2651779195E+02
4 4 2 9.5362 6.9263 40.9082 0.2651779195E+02
4 4 3 9.5362 6.9263 62.8657 0.1124393853E+03
4 4 4 9.5362 6.9263 93.7699 0.9062483731E+02
4 4 5 9.5362 6.9263 93.7699 0.9062483731E+02
4 4 6 9.5362 6.9263 109.2888 0.1417295666E+03
------------------------------------------------------------------------------
Electron-phonon vertex |g| (meV)
@ -605,85 +604,86 @@
ik = 1 coord.: 0.0000000 0.0000000 0.0000000
ibnd jbnd imode enk[eV] enk+q[eV] omega(q)[meV] |g|[meV]
------------------------------------------------------------------------------
2 2 1 9.5362 2.5291 42.3398 0.5933004589E+02
2 2 2 9.5362 2.5291 58.0658 0.5280358644E+02
2 2 3 9.5362 2.5291 65.8200 0.1093263246E+03
2 2 4 9.5362 2.5291 92.0117 0.6149939971E+01
2 2 5 9.5362 2.5291 92.2193 0.4953068343E+02
2 2 6 9.5362 2.5291 105.8096 0.1336299247E+03
2 3 1 9.5362 3.9023 42.3398 0.2307923449E+02
2 3 2 9.5362 3.9023 58.0658 0.4875566957E+02
2 3 3 9.5362 3.9023 65.8200 0.1115183127E+03
2 3 4 9.5362 3.9023 92.0117 0.9111482842E+02
2 3 5 9.5362 3.9023 92.2193 0.5178673920E+02
2 3 6 9.5362 3.9023 105.8096 0.3556518912E+02
2 4 1 9.5362 7.3625 42.3398 0.5010294468E+02
2 4 2 9.5362 7.3625 58.0658 0.2949120021E+02
2 4 3 9.5362 7.3625 65.8200 0.1480592816E+03
2 4 4 9.5362 7.3625 92.0117 0.1078754791E+03
2 4 5 9.5362 7.3625 92.2193 0.1869269993E+02
2 4 6 9.5362 7.3625 105.8096 0.1278261567E+03
3 2 1 9.5362 2.5291 42.3398 0.5933004589E+02
3 2 2 9.5362 2.5291 58.0658 0.5280358644E+02
3 2 3 9.5362 2.5291 65.8200 0.1093263246E+03
3 2 4 9.5362 2.5291 92.0117 0.6149939971E+01
3 2 5 9.5362 2.5291 92.2193 0.4953068343E+02
3 2 6 9.5362 2.5291 105.8096 0.1336299247E+03
3 3 1 9.5362 3.9023 42.3398 0.2307923449E+02
3 3 2 9.5362 3.9023 58.0658 0.4875566957E+02
3 3 3 9.5362 3.9023 65.8200 0.1115183127E+03
3 3 4 9.5362 3.9023 92.0117 0.9111482842E+02
3 3 5 9.5362 3.9023 92.2193 0.5178673920E+02
3 3 6 9.5362 3.9023 105.8096 0.3556518912E+02
3 4 1 9.5362 7.3625 42.3398 0.5010294468E+02
3 4 2 9.5362 7.3625 58.0658 0.2949120021E+02
3 4 3 9.5362 7.3625 65.8200 0.1480592816E+03
3 4 4 9.5362 7.3625 92.0117 0.1078754791E+03
3 4 5 9.5362 7.3625 92.2193 0.1869269993E+02
3 4 6 9.5362 7.3625 105.8096 0.1278261567E+03
4 2 1 9.5362 2.5291 42.3398 0.5933004589E+02
4 2 2 9.5362 2.5291 58.0658 0.5280358644E+02
4 2 3 9.5362 2.5291 65.8200 0.1093263246E+03
4 2 4 9.5362 2.5291 92.0117 0.6149939971E+01
4 2 5 9.5362 2.5291 92.2193 0.4953068343E+02
4 2 6 9.5362 2.5291 105.8096 0.1336299247E+03
4 3 1 9.5362 3.9023 42.3398 0.2307923449E+02
4 3 2 9.5362 3.9023 58.0658 0.4875566957E+02
4 3 3 9.5362 3.9023 65.8200 0.1115183127E+03
4 3 4 9.5362 3.9023 92.0117 0.9111482842E+02
4 3 5 9.5362 3.9023 92.2193 0.5178673920E+02
4 3 6 9.5362 3.9023 105.8096 0.3556518912E+02
4 4 1 9.5362 7.3625 42.3398 0.5010294468E+02
4 4 2 9.5362 7.3625 58.0658 0.2949120021E+02
4 4 3 9.5362 7.3625 65.8200 0.1480592816E+03
4 4 4 9.5362 7.3625 92.0117 0.1078754791E+03
4 4 5 9.5362 7.3625 92.2193 0.1869269993E+02
4 4 6 9.5362 7.3625 105.8096 0.1278261567E+03
2 2 1 9.5362 2.5291 42.3397 0.4664118925E+02
2 2 2 9.5362 2.5291 58.0655 0.4993970589E+02
2 2 3 9.5362 2.5291 65.8200 0.1206035829E+03
2 2 4 9.5362 2.5291 92.0116 0.9092749278E+02
2 2 5 9.5362 2.5291 92.2178 0.4644773796E+02
2 2 6 9.5362 2.5291 105.8068 0.1046925861E+03
2 3 1 9.5362 3.9023 42.3397 0.4664118925E+02
2 3 2 9.5362 3.9023 58.0655 0.4993970589E+02
2 3 3 9.5362 3.9023 65.8200 0.1206035829E+03
2 3 4 9.5362 3.9023 92.0116 0.9092749278E+02
2 3 5 9.5362 3.9023 92.2178 0.4644773796E+02
2 3 6 9.5362 3.9023 105.8068 0.1046925861E+03
2 4 1 9.5362 7.3625 42.3397 0.4664118925E+02
2 4 2 9.5362 7.3625 58.0655 0.4993970589E+02
2 4 3 9.5362 7.3625 65.8200 0.1206035829E+03
2 4 4 9.5362 7.3625 92.0116 0.9092749278E+02
2 4 5 9.5362 7.3625 92.2178 0.4644773796E+02
2 4 6 9.5362 7.3625 105.8068 0.1046925861E+03
3 2 1 9.5362 2.5291 42.3397 0.4180598853E+02
3 2 2 9.5362 2.5291 58.0655 0.5300248882E+02
3 2 3 9.5362 2.5291 65.8200 0.1115028654E+03
3 2 4 9.5362 2.5291 92.0116 0.8015134537E+02
3 2 5 9.5362 2.5291 92.2178 0.4693935033E+02
3 2 6 9.5362 2.5291 105.8068 0.9008050165E+02
3 3 1 9.5362 3.9023 42.3397 0.4180598853E+02
3 3 2 9.5362 3.9023 58.0655 0.5300248882E+02
3 3 3 9.5362 3.9023 65.8200 0.1115028654E+03
3 3 4 9.5362 3.9023 92.0116 0.8015134537E+02
3 3 5 9.5362 3.9023 92.2178 0.4693935033E+02
3 3 6 9.5362 3.9023 105.8068 0.9008050165E+02
3 4 1 9.5362 7.3625 42.3397 0.4180598853E+02
3 4 2 9.5362 7.3625 58.0655 0.5300248882E+02
3 4 3 9.5362 7.3625 65.8200 0.1115028654E+03
3 4 4 9.5362 7.3625 92.0116 0.8015134537E+02
3 4 5 9.5362 7.3625 92.2178 0.4693935033E+02
3 4 6 9.5362 7.3625 105.8068 0.9008050165E+02
4 2 1 9.5362 2.5291 42.3397 0.5137983042E+02
4 2 2 9.5362 2.5291 58.0655 0.2705001339E+02
4 2 3 9.5362 2.5291 65.8200 0.1390387254E+03
4 2 4 9.5362 2.5291 92.0116 0.7269688252E+02
4 2 5 9.5362 2.5291 92.2178 0.3352457384E+02
4 2 6 9.5362 2.5291 105.8068 0.1280165544E+03
4 3 1 9.5362 3.9023 42.3397 0.5137983042E+02
4 3 2 9.5362 3.9023 58.0655 0.2705001339E+02
4 3 3 9.5362 3.9023 65.8200 0.1390387254E+03
4 3 4 9.5362 3.9023 92.0116 0.7269688252E+02
4 3 5 9.5362 3.9023 92.2178 0.3352457384E+02
4 3 6 9.5362 3.9023 105.8068 0.1280165544E+03
4 4 1 9.5362 7.3625 42.3397 0.5137983042E+02
4 4 2 9.5362 7.3625 58.0655 0.2705001339E+02
4 4 3 9.5362 7.3625 65.8200 0.1390387254E+03
4 4 4 9.5362 7.3625 92.0116 0.7269688252E+02
4 4 5 9.5362 7.3625 92.2178 0.3352457384E+02
4 4 6 9.5362 7.3625 105.8068 0.1280165544E+03
------------------------------------------------------------------------------
===================================================================
Memory usage: VmHWM = 31Mb
VmPeak = 305Mb
VmPeak = 306Mb
===================================================================
Unfolding on the coarse grid
dvanqq2 : 0.26s CPU 0.26s WALL ( 28 calls)
elphon_wrap : 14.45s CPU 14.67s WALL ( 1 calls)
dvanqq2 : 0.24s CPU 0.25s WALL ( 27 calls)
elphon_wrap : 13.97s CPU 14.14s WALL ( 1 calls)
INITIALIZATION:
init_vloc : 0.00s CPU 0.00s WALL ( 1 calls)
init_us_1 : 0.29s CPU 0.29s WALL ( 4 calls)
init_us_1 : 0.28s CPU 0.28s WALL ( 4 calls)
newd : 0.00s CPU 0.00s WALL ( 1 calls)
dvanqq2 : 0.24s CPU 0.25s WALL ( 27 calls)
Electron-Phonon interpolation
ephwann : 0.72s CPU 0.75s WALL ( 1 calls)
ephwann : 0.71s CPU 0.73s WALL ( 1 calls)
ep-interp : 0.02s CPU 0.02s WALL ( 4 calls)
Ham: step 1 : 0.00s CPU 0.00s WALL ( 1 calls)
Ham: step 2 : 0.00s CPU 0.00s WALL ( 1 calls)
Ham: step 2 : 0.00s CPU 0.01s WALL ( 1 calls)
ep: step 1 : 0.00s CPU 0.00s WALL ( 162 calls)
ep: step 2 : 0.05s CPU 0.06s WALL ( 162 calls)
DynW2B : 0.00s CPU 0.00s WALL ( 4 calls)
@ -692,7 +692,7 @@
Total program execution
EPW : 17.52s CPU 17.80s WALL
EPW : 16.99s CPU 17.19s WALL
Please consider citing:

View File

@ -1,8 +1,6 @@
--
&inputepw
prefix = 'sic'
amass(1) = 28.0855
amass(2) = 12.0107
outdir = './'
dvscf_dir = './save'
@ -17,10 +15,25 @@
wannierize = .true.
nbndsub = 4
nbndskip = 0
num_iter = 300
proj(1) = 'Si:sp3'
dis_win_max = 12
dis_froz_max= 9.2
proj(1) = 'f= 0.000, 0.000, 0.000, : sp3'
wdata(1) = 'bands_plot = .true.'
wdata(2) = 'begin kpoint_path'
wdata(3) = 'G 0.000 0.000 0.000 X 0.500 0.000 0.500'
wdata(4) = 'X 0.500 0.000 0.500 W 0.500 0.250 0.750'
wdata(5) = 'W 0.500 0.250 0.750 L 0.500 0.500 0.500'
wdata(6) = 'L 0.500 0.500 0.500 K 0.375 0.375 0.750'
wdata(7) = 'K 0.375 0.375 0.750 G 0.000 0.000 0.000'
wdata(8) = 'G 0.000 0.000 0.000 L 0.500 0.500 0.500'
wdata(9) = 'end kpoint_path'
wdata(10) = 'bands_plot_format = gnuplot'
wdata(11) = 'guiding_centres = .true.'
wdata(12) = 'dis_num_iter = 500'
wdata(13) = 'num_print_cycles = 50'
use_ws = .true.
vme = .true.
elecselfen = .false.
phonselfen = .false.