The Elk Code
initeph.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2020 Chung-Yu Wang, J. K. Dewhurst, S. Sharma and
3 ! E. K. U. Gross. This file is distributed under the terms of the GNU General
4 ! Public License. See the file COPYING for license details.
5 
6 subroutine initeph
7 use modmain
8 use modphonon
9 use modbog
10 implicit none
11 ! local variables
12 integer iq,ik,n,i
13 ! automatic arrays
14 complex(8) u(nstsv,nstsv),v(nstsv,nstsv)
15 complex(8) w(nbph,nbph),x(nbph,nbph),y(nbph)
16 complex(8) ev(nbph,nbph)
17 ! allocatable arrays
18 complex(8), allocatable :: ephmat(:,:,:)
19 
20 ! combined target array for fermionic and bosonic density matrices
21 if (allocated(duvwx)) deallocate(duvwx)
22 n=2*nstsv*nstsv*nkpt+2*nbph*nbph*nqpt
23 allocate(duvwx(n))
24 
25 !------------------------------!
26 ! electronic variables !
27 !------------------------------!
28 if (allocated(evaluv)) deallocate(evaluv)
29 allocate(evaluv(nstsv,nkpt))
30 if (allocated(vnorm)) deallocate(vnorm)
31 allocate(vnorm(nstsv,nkpt))
32 ! associate the electronic density matrices with target
33 dvv(1:nstsv,1:nstsv,1:nkpt) => duvwx(1:)
34 i=size(dvv)+1
35 duv(1:nstsv,1:nstsv,1:nkpt) => duvwx(i:)
36 i=i+size(duv)
37 if (task == 270) then
38 ! initialise the density matrices to random numbers
39  dvv(:,:,:)=0.d0
40  duv(:,:,:)=0.d0
41  do ik=1,nkpt
42  call rndevsv(1.d0,dvv(:,:,ik))
43  call rndevsv(1.d0,duv(:,:,ik))
44  end do
45 else
46  do ik=1,nkpt
47 ! get the eigenvalues from file
48  call getevaluv(ik,evaluv(:,ik))
49 ! get the eigenvectors from file
50  call getevecuv(ik,vkl(:,ik),u,v)
51 ! calculate the density matrices
52  call dmatuv(nstsv,efermi,evalsv(:,ik),u,v,dvv(:,:,ik),duv(:,:,ik), &
53  vnorm(:,ik))
54  end do
55 end if
56 
57 !----------------------------!
58 ! phononic variables !
59 !----------------------------!
60 ! initialise the phonon dynamical matrices
61 call initph
62 ! find the eigenvalues of the dynamical matrices and store in global array
63 do iq=1,nqpt
64  call dynev(dynq(:,:,iq),wphq(:,iq),ev)
65 end do
66 if (allocated(evalwx)) deallocate(evalwx)
67 allocate(evalwx(nbph,nqpt))
68 if (allocated(xnorm)) deallocate(xnorm)
69 allocate(xnorm(nbph,nqpt))
70 ! associate the phononic density matrices with target
71 dxx(1:nbph,1:nbph,1:nqpt) => duvwx(i:)
72 i=i+size(dxx)
73 dwx(1:nbph,1:nbph,1:nqpt) => duvwx(i:)
74 if (task == 270) then
75 ! zero the density matrices
76  dxx(:,:,:)=0.d0
77  dwx(:,:,:)=0.d0
78 else
79  do iq=1,nqpt
80 ! get the eigenvalues from file
81  call getevalwx(iq,evalwx(:,iq))
82 ! get the eigenvectors from file
83  call getevecwxy(iq,w,x,y)
84 ! calculate the density matrices
85  call dmatwx(nbph,w,x,dxx(:,:,iq),dwx(:,:,iq),xnorm(:,iq))
86  end do
87 end if
88 
89 !-----------------------------------!
90 ! electron-phonon variables !
91 !-----------------------------------!
92 if (any(task == [270,271])) then
93 ! allocate the electron-phonon matrix elements array
94  if (allocated(ephmkq)) deallocate(ephmkq)
95  allocate(ephmkq(nstsv,nstsv,nbph,nkptnr,nqpt))
96 ! read the matrix elements from file and store in global array
97  allocate(ephmat(nstsv,nstsv,nbph))
98  do iq=1,nqpt
99  do ik=1,nkptnr
100  call getephmat(iq,ik,ephmat)
101  ephmkq(:,:,:,ik,iq)=ephmat(:,:,:)
102  end do
103 ! zero the electron-phonon coupling for phonon small phonon frequencies
104  do i=1,nbph
105  if (wphq(i,iq) < wphcut) ephmkq(:,:,i,:,iq)=0.d0
106  end do
107  end do
108  deallocate(ephmat)
109 end if
110 
111 end subroutine
112 
real(8) efermi
Definition: modmain.f90:907
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
integer task
Definition: modmain.f90:1299
real(8), dimension(:,:), allocatable evalwx
Definition: modbog.f90:37
integer nqpt
Definition: modmain.f90:525
integer nkpt
Definition: modmain.f90:461
real(8), dimension(:,:), allocatable evaluv
Definition: modbog.f90:15
subroutine getevalwx(iq, evalwxp)
Definition: getevalwx.f90:7
subroutine rndevsv(rndm, evecsv)
Definition: rndevsv.f90:7
subroutine getevecuv(ikp, vpl, u, v)
Definition: getevecuv.f90:7
integer nkptnr
Definition: modmain.f90:463
subroutine dmatwx(n, w, x, dxx, dwx, xn)
Definition: dmatwx.f90:7
complex(8), dimension(:,:,:), pointer, contiguous dwx
Definition: modbog.f90:43
subroutine getevecwxy(iq, w, x, y)
Definition: getevecwxy.f90:7
real(8), dimension(:,:), allocatable xnorm
Definition: modbog.f90:41
subroutine getevaluv(ik, evaluvp)
Definition: getevaluv.f90:7
subroutine dynev(dq, wq, ev)
Definition: dynev.f90:7
complex(8), dimension(:,:,:), allocatable dynq
Definition: modphonon.f90:27
Definition: modbog.f90:6
subroutine dmatuv(n, ef, e, u, v, dvv, duv, vn)
Definition: dmatuv.f90:7
real(8) wphcut
Definition: modphonon.f90:148
subroutine initph
Definition: initph.f90:7
complex(8), dimension(:), allocatable, target duvwx
Definition: modbog.f90:9
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
complex(8), dimension(:,:,:), pointer, contiguous duv
Definition: modbog.f90:25
complex(8), dimension(:,:,:), pointer, contiguous dxx
Definition: modbog.f90:43
real(8), dimension(:,:), allocatable vnorm
Definition: modbog.f90:17
subroutine getephmat(iqp, ikp, ephmat)
Definition: getephmat.f90:7
subroutine initeph
Definition: initeph.f90:7
complex(8), dimension(:,:,:), pointer, contiguous dvv
Definition: modbog.f90:25
real(8), dimension(:,:), allocatable wphq
Definition: modphonon.f90:31
complex(4), dimension(:,:,:,:,:), allocatable ephmkq
Definition: modphonon.f90:157