! !----------------------------------------------------------------------- subroutine ephbloch2wane ( nbnd, nbndsub, nks, nksq, nkbl, lgamma, xk, & cu, epmatk, nrr, irvec, wslen, epmatw) !----------------------------------------------------------------------- ! ! From the electron-phonon matrix elements in Bloch representation (coarse mesh), ! find the corresponding matrix elements in Wannier representation ! ! input : ! ! output : ! ! Feliciano Giustino, UCB ! !----------------------------------------------------------------------- ! #include "f_defs.h" ! USE kinds, only : DP use pwcom, only : at, bg, celldm use control_flags, ONLY : iverbosity #ifdef __PARA use para USE io_global, ONLY : ionode_id USE mp_global, ONLY : mpime #endif implicit none ! ! input variables ! integer :: nbnd, nbndsub, nks, nksq, nkbl, nrr, irvec (3, nrr) ! number of bands ! number of bands in the optimal subspace ! number of kpoints ! number of kpoint blocks, in the pool ! number of kpoint blocks, total ! number of WS points and coordinates logical :: lgamma ! true if we are working with q=0 real(kind=DP) :: xk (3, nks), wslen (nrr) ! kpoint coordinates (cartesian in units of 2piba) ! WS vectors length (alat units) complex(kind=DP) :: cu (nbnd, nbndsub, nks), epmatk ( nbnd, nbnd, nks) ! rotation matrix from wannier code ! e-p matrix in bloch representation, coarse mesh ! ! output variables ! complex(kind=DP) :: epmatw ( nbndsub, nbndsub, nrr) ! e-p matrix in wannier basis ! ! work variables ! complex(kind=DP) :: epmats (nbndsub, nbndsub, nksq), eptmp(nbndsub, nbnd) ! e-p matrix in smooth Bloch basis, coarse mesh ! e-p matrix, temporary ! integer :: ik, ikk, ikq, ibnd, jbnd, ir, pbnd, qbnd real(kind=DP), parameter :: bohr2ang = 0.5291772108, twopi = 6.28318530717959 real(kind=DP) :: rdotk, tmp complex(kind=DP), parameter :: ci = (0.d0,1.d0), czero = (0.d0, 0.d0), & cone = (1.d0, 0.d0) complex(kind=DP) :: cfac ! ! !---------------------------------------------------------- ! STEP 1: rotation to optimally smooth Bloch states !---------------------------------------------------------- ! ! g~ = U_k+q^\dagger g U_k ! ! g is epmatk (ibnd, jbnd, ik) ! g~ is epmats (ibnd, jbnd, ik) ! call start_clock ( 'ep: step 1' ) ! do ik = 1, nksq ! if (lgamma) then ikk = ik ikq = ik else ikk = 2 * ik - 1 ikq = ikk + 1 endif ! ! the two zgemm calls perform the following ops: ! epmats = [ cu(ikq)^\dagger * epmatk ] * cu(ikk) ! [here we have a size-reduction from nbnd*nbnd to nbndsub*nbndsub] ! call zgemm ('c', 'n', nbndsub, nbnd, nbnd, cone, cu(:,:,ikq), & nbnd, epmatk(:,:,ik), nbnd, czero, eptmp, nbndsub) call zgemm ('n', 'n', nbndsub, nbndsub, nbnd, cone, eptmp, & nbndsub, cu(:,:,ikk), nbnd, czero, epmats(:,:,ik), nbndsub) ! enddo ! call stop_clock ( 'ep: step 1' ) ! !---------------------------------------------------------------------- ! STEP 2: Fourier transform to obtain matrix elements in wannier basis !---------------------------------------------------------------------- ! ! g (R) = (1/nkc) sum_k e^{-ikR} g~(k) ! ! epmatw (nbndsub,nbndsub,ir) is g(R) ! call start_clock ( 'ep: step 2' ) ! epmatw (:, :, :) = czero ! ! bring xk in crystal coordinates ! call cryst_to_cart (nks, xk, at, -1) ! do ir = 1, nrr ! do ik = 1, nksq ! if (lgamma) then ikk = ik else ikk = 2 * ik - 1 endif ! rdotk = twopi * dot_product( xk ( :, ikk), float(irvec( :, ir) )) cfac = exp( -ci*rdotk ) / float(nkbl) epmatw ( :, :, ir) = epmatw ( :, :, ir) + cfac * epmats ( :, :, ik) ! enddo ! enddo ! #ifdef __PARA call poolreduce (2 * nrr * nbndsub * nbndsub, epmatw ) #endif ! ! bring xk back into cart coord ! call cryst_to_cart (nks, xk, bg, 1) ! !@ if (iverbosity.eq.1) then ! ! check spatial decay of matrix elements in Wannier basis ! the unit in r-space is angstrom, and I am plotting ! the matrix for the first mode only ! #ifdef __PARA if (mpime.eq.ionode_id) then #endif !@ write(301, '(/3x,a/)') 'Spatial decay of e-p matrix elements in Wannier basis' do ir = 1, nrr ! tmp = maxval ( abs(epmatw(:,:,ir)) ) write(301, *) wslen(ir) * celldm (1) * bohr2ang, tmp ! enddo #ifdef __PARA endif #endif !@ endif ! call stop_clock ( 'ep: step 2' ) ! return end subroutine ephbloch2wane !