Phonon part pf last commit was missing!

git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@13229 c92efa57-630b-4861-b058-cf58834340f0
This commit is contained in:
giannozz 2017-01-08 14:38:26 +00:00
parent f86660892c
commit 1695899967
3 changed files with 33 additions and 35 deletions

View File

@ -136,7 +136,7 @@ SUBROUTINE elph_tetra_lambda()
CALL symdyn_munu_new (el_ph_sum, u, xq, s, invs, rtau, irt, at, bg, nsymq, nat, irotmq, minus_q)
!
dosef(1:2) = 0.0_dp
CALL opt_tetra_dos_t (et_col, nspin, nbnd, nkstot, ntetra, tetra, ef, dosef)
CALL opt_tetra_dos_t (et_col, nspin, nbnd, nkstot, ef, dosef)
!
dosef(1:2) = 0.5_dp * dosef(1:2)
!
@ -227,7 +227,7 @@ SUBROUTINE elph_tetra_delta1(nbnd_fs,iq,tfst,tlst,et_col,wght)
USE wvfct, ONLY: nbnd
USE klist, ONLY: nkstot
USE lsda_mod, ONLY : nspin
USE ktetra, ONLY : tetra, ntetra, wlsm
USE ktetra, ONLY : tetra, ntetra, nntetra, wlsm
!
INTEGER,INTENT(IN) :: nbnd_fs, iq, tfst, tlst
REAL(dp),INTENT(IN) :: et_col(nbnd, nkstot)
@ -259,7 +259,7 @@ SUBROUTINE elph_tetra_delta1(nbnd_fs,iq,tfst,tlst,et_col,wght)
!
ei0(1:4, 1:nbnd_fs) = 0.0_dp
ej0(1:4, 1:nbnd_fs) = 0.0_dp
DO ii = 1, 20
DO ii = 1, nntetra
!
ik = tetra(ii, nt) + nk
DO ibnd = 1, nbnd_fs
@ -357,7 +357,7 @@ SUBROUTINE elph_tetra_delta1(nbnd_fs,iq,tfst,tlst,et_col,wght)
!
END DO
!
DO ii = 1, 20
DO ii = 1, nntetra
!
ik = tetra(ii, nt) + nk
DO jj = 1, 4
@ -371,8 +371,10 @@ SUBROUTINE elph_tetra_delta1(nbnd_fs,iq,tfst,tlst,et_col,wght)
!
END DO ! ns
!
wght(1:nbnd_fs,1:nbnd_fs,1:nkstot) = wght(1:nbnd_fs,1:nbnd_fs,1:nkstot) / REAL(ntetra, dp)
IF(nspin == 1) wght(1:nbnd_fs,1:nbnd_fs,1:nkstot) = wght(1:nbnd_fs,1:nbnd_fs,1:nkstot) * 2.0_dp
wght(1:nbnd_fs,1:nbnd_fs,1:nkstot) = wght(1:nbnd_fs,1:nbnd_fs,1:nkstot) &
/ REAL(ntetra, dp)
IF(nspin == 1) wght(1:nbnd_fs,1:nbnd_fs,1:nkstot) = 2.0_dp * &
wght(1:nbnd_fs,1:nbnd_fs,1:nkstot)
!
CALL mp_sum(wght, intra_image_comm)
!
@ -542,7 +544,7 @@ SUBROUTINE elph_tetra_gamma()
END DO
!
dosef(1:2) = 0.0_dp
CALL opt_tetra_dos_t (et_col, nspin, nbnd, nkstot, ntetra, tetra, ef, dosef)
CALL opt_tetra_dos_t (et_col, nspin, nbnd, nkstot, ef, dosef)
!
dosef(1:2) = 0.5_dp * dosef(1:2)
!
@ -635,7 +637,7 @@ SUBROUTINE elph_tetra_step1(nbnd_fs,iq,tfst,tlst,et_col,wght)
USE wvfct, ONLY: nbnd
USE klist, ONLY: nkstot
USE lsda_mod, ONLY : nspin
USE ktetra, ONLY : ntetra, tetra, wlsm
USE ktetra, ONLY : ntetra, tetra, nntetra, wlsm
!
INTEGER,INTENT(IN) :: nbnd_fs, iq, tfst, tlst
REAL(dp),INTENT(IN) :: et_col(nbnd, nkstot)
@ -667,7 +669,7 @@ SUBROUTINE elph_tetra_step1(nbnd_fs,iq,tfst,tlst,et_col,wght)
!
ei0(1:4, 1:nbnd_fs) = 0.0_dp
ej0(1:4, 1:nbnd_fs) = 0.0_dp
DO ii = 1, 20
DO ii = 1, nntetra
!
ik = tetra(ii, nt) + nk
DO ibnd = 1, nbnd_fs
@ -874,7 +876,7 @@ SUBROUTINE elph_tetra_step1(nbnd_fs,iq,tfst,tlst,et_col,wght)
!
END DO ! ibnd
!
DO ii = 1, 20
DO ii = 1, nntetra
!
ik = tetra(ii, nt) + nk
!

View File

@ -829,6 +829,7 @@ matdyn.o : ../../Modules/parameters.o
matdyn.o : ../../Modules/parser.o
matdyn.o : ../../PW/src/pwcom.o
matdyn.o : ../../PW/src/symm_base.o
matdyn.o : ../../PW/src/tetra.o
matdyn.o : io_dyn_mat.o
mix_pot.o : ../../Modules/io_files.o
mix_pot.o : ../../Modules/kind.o

View File

@ -143,6 +143,7 @@ PROGRAM matdyn
USE rap_point_group, ONLY : code_group
USE bz_form, ONLY : transform_label_coord
USE parser, ONLY : read_line
USE ktetra, ONLY : tetra_dos_t
USE ifconstants, ONLY : frc, atm, zeu, tau_blk, ityp_blk, m_loc
!
@ -156,14 +157,14 @@ PROGRAM matdyn
INTEGER:: nax, nax_blk
INTEGER, PARAMETER:: ntypx=10, nrwsx=200
REAL(DP), PARAMETER :: eps=1.0d-6
INTEGER :: nr1, nr2, nr3, nsc, nk1, nk2, nk3, ntetra, ibrav
INTEGER :: nr1, nr2, nr3, nsc, nk1, nk2, nk3, ibrav
CHARACTER(LEN=256) :: flfrc, flfrq, flvec, fltau, fldos, filename, fldyn, fleig
CHARACTER(LEN=10) :: asr
LOGICAL :: dos, has_zstar, q_in_cryst_coord, eigen_similarity
COMPLEX(DP), ALLOCATABLE :: dyn(:,:,:,:), dyn_blk(:,:,:,:), frc_ifc(:,:,:,:)
COMPLEX(DP), ALLOCATABLE :: z(:,:)
REAL(DP), ALLOCATABLE:: tau(:,:), q(:,:), w2(:,:), freq(:,:), wq(:)
INTEGER, ALLOCATABLE:: tetra(:,:), ityp(:), itau_blk(:)
INTEGER, ALLOCATABLE:: ityp(:), itau_blk(:)
REAL(DP) :: omega,alat, &! cell parameters and volume
at_blk(3,3), bg_blk(3,3), &! original cell
omega_blk, &! original cell volume
@ -182,7 +183,7 @@ PROGRAM matdyn
!
LOGICAL :: readtau, la2F, xmlifc, lo_to_split, na_ifc, fd, nosym
!
REAL(DP) :: qhat(3), qh, DeltaE, Emin=0._dp, Emax, E, DOSofE(1), qq
REAL(DP) :: qhat(3), qh, DeltaE, Emin=0._dp, Emax, E, DOSofE(2), qq
REAL(DP) :: delta, pathL
REAL(DP), ALLOCATABLE :: xqaux(:,:)
INTEGER, ALLOCATABLE :: nqb(:)
@ -391,11 +392,10 @@ PROGRAM matdyn
IF (dos) THEN
IF (nk1 < 1 .OR. nk2 < 1 .OR. nk3 < 1) &
CALL errore ('matdyn','specify correct q-point grid!',1)
ntetra = 6 * nk1 * nk2 * nk3
nqx = nk1*nk2*nk3
ALLOCATE ( tetra(4,ntetra), q(3,nqx) )
ALLOCATE ( q(3,nqx) )
CALL gen_qpoints (ibrav, at, bg, nat, tau, ityp, nk1, nk2, nk3, &
ntetra, nqx, nq, q, tetra, nosym)
nqx, nq, q, nosym)
ELSE
!
! read q-point list
@ -403,7 +403,6 @@ PROGRAM matdyn
IF (ionode) READ (5,*) nq
CALL mp_bcast(nq, ionode_id, world_comm)
ALLOCATE ( q(3,nq) )
ALLOCATE( tetra(1,1) )
IF (.NOT.q_in_band_form) THEN
DO n = 1,nq
IF (ionode) READ (5,*) (q(i,n),i=1,3)
@ -711,7 +710,7 @@ PROGRAM matdyn
IF (ionode) OPEN (unit=2,file=fldos,status='unknown',form='formatted')
DO n= 1, ndos
E = Emin + (n - 1) * DeltaE
CALL dos_t(freq, 1, 3*nat, nq, ntetra, tetra, E, DOSofE)
CALL tetra_dos_t(freq, 1, 3*nat, nq, E, DOSofE)
!
! The factor 0.5 corrects for the factor 2 in dos_t,
! that accounts for the spin in the electron DOS.
@ -745,7 +744,7 @@ PROGRAM matdyn
call a2Fdos (nat, nq, nr1, nr2, nr3, ibrav, at, bg, tau, alat, &
nsc, nat_blk, at_blk, bg_blk, itau_blk, omega_blk, &
rws, nrws, dos, Emin, DeltaE, ndos, &
ntetra, tetra, asr, q, freq,fd)
asr, q, freq,fd)
!
IF (.NOT.dos) THEN
DO isig=1,10
@ -1959,7 +1958,7 @@ END SUBROUTINE write_tau
!
!-----------------------------------------------------------------------
SUBROUTINE gen_qpoints (ibrav, at_, bg_, nat, tau, ityp, nk1, nk2, nk3, &
ntetra, nqx, nq, q, tetra, nosym)
nqx, nq, q, nosym)
!-----------------------------------------------------------------------
!
USE kinds, ONLY : DP
@ -1970,11 +1969,11 @@ SUBROUTINE gen_qpoints (ibrav, at_, bg_, nat, tau, ityp, nk1, nk2, nk3, &
!
IMPLICIT NONE
! input
INTEGER :: ibrav, nat, nk1, nk2, nk3, ntetra, ityp(*)
INTEGER :: ibrav, nat, nk1, nk2, nk3, ityp(*)
REAL(DP) :: at_(3,3), bg_(3,3), tau(3,nat)
LOGICAL :: nosym
! output
INTEGER :: nqx, nq, tetra(4,ntetra)
INTEGER :: nqx, nq
REAL(DP) :: q(3,nqx)
! local
REAL(DP) :: xqq(3), wk(nqx), mdum(3,nat)
@ -2000,11 +1999,8 @@ SUBROUTINE gen_qpoints (ibrav, at_, bg_, nat, tau, ityp, nk1, nk2, nk3, &
CALL irreducible_BZ (nrot, s, nsym, time_reversal, magnetic_sym, &
at, bg, nqx, nq, q, wk, t_rev)
!
IF (ntetra /= 6 * nk1 * nk2 * nk3) &
CALL errore ('gen_qpoints','inconsistent ntetra',1)
!
CALL tetra_init (nsym, s, time_reversal, t_rev, at, bg, nqx, 0, 0, 0, &
nk1, nk2, nk3, nq, q, ntetra, tetra)
nk1, nk2, nk3, nq, q)
!
RETURN
END SUBROUTINE gen_qpoints
@ -2013,7 +2009,7 @@ END SUBROUTINE gen_qpoints
SUBROUTINE a2Fdos &
(nat, nq, nr1, nr2, nr3, ibrav, at, bg, tau, alat, &
nsc, nat_blk, at_blk, bg_blk, itau_blk, omega_blk, rws, nrws, &
dos, Emin, DeltaE, ndos, ntetra, tetra, asr, q, freq,fd )
dos, Emin, DeltaE, ndos, asr, q, freq,fd )
!-----------------------------------------------------------------------
!
USE kinds, ONLY : DP
@ -2023,11 +2019,11 @@ SUBROUTINE a2Fdos &
USE mp_images, ONLY : intra_image_comm
USE ifconstants, ONLY : zeu, tau_blk
USE constants, ONLY : pi, RY_TO_THZ
USE ktetra, ONLY : tetra_init
!
IMPLICIT NONE
!
INTEGER, INTENT(in) :: nat, nq, nr1, nr2, nr3, ibrav, ndos, ntetra, &
tetra(4, ntetra)
INTEGER, INTENT(in) :: nat, nq, nr1, nr2, nr3, ibrav, ndos
LOGICAL, INTENT(in) :: dos,fd
CHARACTER(LEN=*), INTENT(IN) :: asr
REAL(DP), INTENT(in) :: freq(3*nat,nq), q(3,nq), at(3,3), bg(3,3), &
@ -2170,8 +2166,7 @@ SUBROUTINE a2Fdos &
dos_tot = 0.0d0
do j=1,nmodes
!
dos_a2F(j) = dos_gam(nmodes, nq, j, ntetra, tetra, &
gamma, freq, E)
dos_a2F(j) = dos_gam(nmodes, nq, j, gamma, freq, E)
dos_a2F(j) = dos_a2F(j) / dos_ee(isig) / 2.d0 / pi
dos_tot = dos_tot + dos_a2F(j)
!
@ -2297,7 +2292,7 @@ subroutine setgam (q, gam, nat, at,bg,tau,itau_blk,nsc,alat, &
end subroutine setgam
!
!--------------------------------------------------------------------
function dos_gam (nbndx, nq, jbnd, ntetra, tetra, gamma, et, ef)
function dos_gam (nbndx, nq, jbnd, gamma, et, ef)
!--------------------------------------------------------------------
! calculates weights with the tetrahedron method (Bloechl version)
! this subroutine is based on tweights.f90 belonging to PW
@ -2306,11 +2301,11 @@ function dos_gam (nbndx, nq, jbnd, ntetra, tetra, gamma, et, ef)
! and "et" means the frequency(mode,q-point)
!
USE kinds, ONLY: DP
use parameters
! USE ifconstants, ONLY : gamma
USE parameters
USE ktetra, ONLY : ntetra, tetra
implicit none
!
integer :: nq, nbndx, ntetra, tetra(4,ntetra), jbnd
integer :: nq, nbndx, jbnd
real(DP) :: et(nbndx,nq), gamma(nbndx,nq), func
real(DP) :: ef