The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine initeph
7use modmain
8use modphonon
9use modbog
10implicit none
11! local variables
12integer iq,ik,n,i
13! automatic arrays
14complex(8) u(nstsv,nstsv),v(nstsv,nstsv)
15complex(8) w(nbph,nbph),x(nbph,nbph),y(nbph)
16complex(8) ev(nbph,nbph)
17! allocatable arrays
18complex(8), allocatable :: ephmat(:,:,:)
19
20! combined target array for fermionic and bosonic density matrices
21if (allocated(duvwx)) deallocate(duvwx)
22n=2*nstsv*nstsv*nkpt+2*nbph*nbph*nqpt
23allocate(duvwx(n))
24
25!------------------------------!
26! electronic variables !
27!------------------------------!
28if (allocated(evaluv)) deallocate(evaluv)
29allocate(evaluv(nstsv,nkpt))
30if (allocated(vnorm)) deallocate(vnorm)
31allocate(vnorm(nstsv,nkpt))
32! associate the electronic density matrices with target
33dvv(1:nstsv,1:nstsv,1:nkpt) => duvwx(1:)
34i=nstsv*nstsv*nkpt+1
35duv(1:nstsv,1:nstsv,1:nkpt) => duvwx(i:)
36i=i+nstsv*nstsv*nkpt
37if (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
45else
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
55end if
56
57!----------------------------!
58! phononic variables !
59!----------------------------!
60! initialise the phonon dynamical matrices
61call initph
62! find the eigenvalues of the dynamical matrices and store in global array
63do iq=1,nqpt
64 call dynev(dynq(:,:,iq),wphq(:,iq),ev)
65end do
66if (allocated(evalwx)) deallocate(evalwx)
67allocate(evalwx(nbph,nqpt))
68if (allocated(xnorm)) deallocate(xnorm)
69allocate(xnorm(nbph,nqpt))
70! associate the phononic density matrices with target
71dxx(1:nbph,1:nbph,1:nqpt) => duvwx(i:)
72i=i+nbph*nbph*nqpt
73dwx(1:nbph,1:nbph,1:nqpt) => duvwx(i:)
74if (task == 270) then
75! zero the density matrices
76 dxx(:,:,:)=0.d0
77 dwx(:,:,:)=0.d0
78else
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
87end if
88
89!-----------------------------------!
90! electron-phonon variables !
91!-----------------------------------!
92if (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)
109end if
110
111end subroutine
112
subroutine dmatuv(n, ef, e, u, v, dvv, duv, vn)
Definition dmatuv.f90:7
subroutine dmatwx(n, w, x, dxx, dwx, xn)
Definition dmatwx.f90:7
subroutine dynev(dq, wq, ev)
Definition dynev.f90:7
subroutine getephmat(iqp, ikp, ephmat)
Definition getephmat.f90:7
subroutine getevaluv(ik, evaluvp)
Definition getevaluv.f90:7
subroutine getevalwx(iq, evalwxp)
Definition getevalwx.f90:7
subroutine getevecuv(ikp, vpl, u, v)
Definition getevecuv.f90:7
subroutine getevecwxy(iq, w, x, y)
Definition getevecwxy.f90:7
subroutine initeph
Definition initeph.f90:7
subroutine initph
Definition initph.f90:7
complex(8), dimension(:,:,:), pointer, contiguous duv
Definition modbog.f90:25
complex(8), dimension(:), allocatable, target duvwx
Definition modbog.f90:9
complex(8), dimension(:,:,:), pointer, contiguous dvv
Definition modbog.f90:25
real(8), dimension(:,:), allocatable evalwx
Definition modbog.f90:37
complex(8), dimension(:,:,:), pointer, contiguous dxx
Definition modbog.f90:43
complex(8), dimension(:,:,:), pointer, contiguous dwx
Definition modbog.f90:43
real(8), dimension(:,:), allocatable xnorm
Definition modbog.f90:41
real(8), dimension(:,:), allocatable evaluv
Definition modbog.f90:15
real(8), dimension(:,:), allocatable vnorm
Definition modbog.f90:17
real(8) efermi
Definition modmain.f90:904
integer nkptnr
Definition modmain.f90:463
integer nqpt
Definition modmain.f90:525
integer nkpt
Definition modmain.f90:461
integer task
Definition modmain.f90:1298
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
real(8) wphcut
complex(4), dimension(:,:,:,:,:), allocatable ephmkq
complex(8), dimension(:,:,:), allocatable dynq
Definition modphonon.f90:27
real(8), dimension(:,:), allocatable wphq
Definition modphonon.f90:31
subroutine rndevsv(rndm, evecsv)
Definition rndevsv.f90:7