sub_gwdos.f90

Go to the documentation of this file.
00001 subroutine calc_dos_GW(nkb1,nkb2,nkb3,NTK,NWF,nsgm,idlt,dmna,dmnr,shift_value,b1,b2,b3,sgmw,SK0,EMK,AW) 
00002  !
00003  use m_tetrahedron
00004  !
00005  implicit none 
00006  integer::nkb1,nkb2,nkb3,NTK,NWF,nsgm 
00007  real(8)::idlt,dmna,dmnr,shift_value 
00008  real(8)::b1(3),b2(3),b3(3),SK0(3,NTK),sgmw(nsgm)  
00009  complex(8)::EMK(NWF,NTK,nsgm)           
00010  integer::jb,ie,ikb1,ikb2,ikb3,ik   
00011  real(8)::delta,ReZ,ImZ,SUM_REAL 
00012  complex(8)::en  
00013  real(8)::AW(nsgm) 
00014  !
00015  !ttrhdrn
00016  !
00017  integer,allocatable::imt1(:)!imt1(4*nkb1*nkb2*nkb3*6)  
00018  integer,allocatable::index_kpt(:,:,:)!index_kpt(nkb1,nkb2,nkb3)    
00019  complex(8),allocatable::fk_1D(:)!fk_1D(NTK)
00020  complex(8),allocatable::gk_1D(:)!gk_1D(NTK)
00021  complex(8),allocatable::fk_3D(:,:,:)!fk_3D(nkb1,nkb2,nkb3) 
00022  complex(8),allocatable::gk_3D(:,:,:)!gk_3D(nkb1,nkb2,nkb3)
00023  complex(8),allocatable::xo(:,:,:)!xo(nkb1,nkb2,nkb3) 
00024  complex(8),allocatable::xow(:,:,:,:)!xow(nsgm,nkb1,nkb2,nkb3)
00025  !
00026  real(8),parameter::pi=dacos(-1.0d0)
00027  real(8),parameter::au=27.21151d0
00028  !
00029  allocate(imt1(4*nkb1*nkb2*nkb3*6)); imt1=0   
00030  allocate(index_kpt(nkb1,nkb2,nkb3)); index_kpt=0  
00031  call make_index_kpt(NTK,nkb1,nkb2,nkb3,SK0(1,1),index_kpt(1,1,1)) 
00032  call ttrhdrn_mkidx(nkb1,nkb2,nkb3,imt1(1),b1(1),b2(1),b3(1))
00033  !
00034  allocate(xow(nsgm,nkb1,nkb2,nkb3)); xow=0.0d0 
00035 !$OMP PARALLEL PRIVATE(ie,jb,delta,fk_1D,gk_1D,ik,en,ReZ,ImZ,fk_3D,gk_3D,ikb3,ikb2,ikb1,xo) 
00036  allocate(fk_1D(NTK)); fk_1D=0.0d0 
00037  allocate(gk_1D(NTK)); gk_1D=0.0d0 
00038  allocate(fk_3D(nkb1,nkb2,nkb3)); fk_3D=0.0d0 
00039  allocate(gk_3D(nkb1,nkb2,nkb3)); gk_3D=0.0d0 
00040  allocate(xo(nkb1,nkb2,nkb3)); xo=0.0d0  
00041 !$OMP DO 
00042  do ie=1,nsgm 
00043   do jb=1,NWF 
00044    !
00045    delta=sgmw(ie)+shift_value 
00046    !
00047    fk_1D=0.0d0
00048    gk_1D=0.0d0 
00049    do ik=1,NTK 
00050     en=EMK(jb,ik,ie)
00051     ReZ=delta-dble(en) 
00052     ImZ=-(idlt+dabs(dimag(en)))!retarded 
00053     gk_1D(ik)=cmplx(ReZ,ImZ)
00054     fk_1D(ik)=1.0d0 
00055    enddo!ik 
00056    !
00057    fk_3D=0.0d0
00058    gk_3D=0.0d0 
00059    do ikb3=1,nkb3
00060     do ikb2=1,nkb2
00061      do ikb1=1,nkb1
00062       ik=index_kpt(ikb1,ikb2,ikb3)
00063       fk_3D(ikb1,ikb2,ikb3)=fk_1D(ik)
00064       gk_3D(ikb1,ikb2,ikb3)=gk_1D(ik)
00065      enddo!ikb1 
00066     enddo!ikb2 
00067    enddo!ikb3 
00068    !
00069    xo=0.0d0 
00070    call ttrhdrn_simple(dmna,dmnr,nkb1,nkb2,nkb3,imt1(1),fk_3D(1,1,1),gk_3D(1,1,1),xo(1,1,1))
00071    xow(ie,:,:,:)=xow(ie,:,:,:)+xo(:,:,:) 
00072    !
00073   enddo!jb 
00074  enddo!ie
00075 !$OMP END DO
00076  deallocate(fk_1D,gk_1D,fk_3D,gk_3D,xo) 
00077 !$OMP END PARALLEL
00078  !
00079  AW=0.0d0 
00080  do ie=1,nsgm
00081   SUM_REAL=0.0d0 
00082   do ikb3=1,nkb3
00083    do ikb2=1,nkb2
00084     do ikb1=1,nkb1
00085      SUM_REAL=SUM_REAL+dabs(dimag(xow(ie,ikb1,ikb2,ikb3))) 
00086     enddo!ikb1 
00087    enddo!ikb2 
00088   enddo!ikb3 
00089  AW(ie)=2.0d0*SUM_REAL/dble(NTK)/pi 
00090  enddo!ie  
00091  !
00092  deallocate(imt1,index_kpt,xow) 
00093  !
00094 RETURN  
00095 END 

Generated on 17 Nov 2020 for respack by  doxygen 1.6.1