m_rd_dat_eps.f90

Go to the documentation of this file.
00001 module m_rd_dat_eps 
00002 use m_rd_dat_wfn
00003 implicit none
00004 public::rd_dat_chi_cutoff
00005 public::rd_dat_ttrhdrn 
00006 public::rd_dat_wgrid 
00007 public::rd_dat_sq 
00008 public::rd_dat_eps 
00009 !chi_cutoff(300)  
00010 real(8),public::Ecut_for_eps
00011 !ttrhdrn(302)  
00012 real(8),public::idlt!Green's function delt (au)!ttrhdrn
00013 real(8),public::dmna!dmna (au)!ttrhdrn
00014 real(8),public::dmnr!dmnr (au)!ttrhdrn
00015 !wgrid(135)  
00016 integer,public::Num_freq_grid
00017 integer,public::ne 
00018 complex(8),public,allocatable::em(:)!em(ne)
00019 complex(8),allocatable::pole_of_chi(:)!pole_of_chi(ne)  
00020 complex(8),allocatable::mat_b(:,:)!mat_b(ne,ne) 
00021 !sq(301)  
00022 integer,public::Nq_irr!TOTAL NUMBER OF irreducible q POINTS 
00023 integer,public::NTQ!TOTAL NUMBER OF Q POINTS IN MK MESHES 
00024 integer,public::NTGQ!TOTAL NUMBER OF G VECTORS FOR EPSILON
00025 integer,public::Lq1,Lq2,Lq3  
00026 real(8),public,allocatable::SQI(:,:)!SQI(3,Nq_irr) 
00027 real(8),public,allocatable::SQ(:,:)!SQ(3,NTQ)  
00028 integer,public,allocatable::numirrq(:)!numirrq(NTQ)  
00029 integer,public,allocatable::numrotq(:)!numrotq(NTQ)  
00030 integer,public,allocatable::trsq(:)!trsq(NTQ)  
00031 integer,public,allocatable::RWq(:,:)!RWq(3,NTQ)
00032 integer,public,allocatable::LG0(:,:,:)!LG0(3,NTG,NTQ)    
00033 integer,public,allocatable::NGQ_eps(:)!NGQ_eps(NTQ)
00034 integer,public,allocatable::NGQ_psi(:)!NGQ_psi(NTQ)  
00035 integer,public,allocatable::packingq(:,:,:,:)!packingq(-Lq1:Lq1,-Lq2:Lq2,-Lq3:Lq3,Nq_irr) 
00036 !epsqw(600-) 
00037 !complex(8),public,allocatable::epstmp(:,:,:)!epstmp(NTGQ,NTGQ,ne)
00038 !complex(8),public,allocatable::epstmpgm(:,:,:,:)!epstmpgm(NTGQ,NTGQ,ne,3)
00039 complex(4),public,allocatable::epsirr(:,:,:,:)!epsirr(NTGQ,NTGQ,ne,Nq_irr) 
00040 contains
00041 !--
00042 subroutine rd_dat_chi_cutoff
00043 implicit none 
00044 OPEN(300,FILE='./dir-eps/dat.chi_cutoff')
00045 read(300,*) Ecut_for_eps
00046 end subroutine
00047 !--
00048 !20180423 
00049 !--
00050 subroutine rd_dat_ttrhdrn
00051 implicit none 
00052 OPEN(302,FILE='./dir-eps/dat.ttrhdrn')
00053 read(302,*) idlt!ttrhdrn Green's function delt (au)
00054 read(302,*) dmna!(au) 
00055 read(302,*) dmnr!(au) 
00056 end subroutine
00057 !--
00058 subroutine rd_dat_wgrid 
00059 implicit none 
00060 integer::ie,je,ke 
00061 real(8)::x,y
00062 complex(8),allocatable::emp(:)!emp(ne+1)
00063 complex(8),allocatable::mat_c(:,:)!mat_c(ne,ne) 
00064 complex(8)::sum_cmpx 
00065 !
00066 OPEN(135,FILE='./dir-eps/dat.wgrid') 
00067 rewind(135) 
00068 read(135,*) Num_freq_grid
00069 ne=Num_freq_grid
00070 allocate(em(ne)); em(:)=0.0d0 
00071 do ie=1,ne 
00072  read(135,'(2f15.10)') em(ie)
00073 enddo 
00074 !
00075 allocate(emp(ne+1));emp(:)=0.0d0 
00076 do ie=1,ne
00077  emp(ie)=em(ie)
00078 enddo
00079 emp(ne+1)=em(ne)+1.3d0*(em(ne)-em(ne-1)) 
00080 !
00081 !pole of chi
00082 !
00083 allocate(pole_of_chi(ne)); pole_of_chi(:)=0.0d0 
00084 do ie=1,ne 
00085  x=dble((emp(ie+1)+emp(ie))/2.0d0) 
00086  y=-dble(1.5d0*(emp(ie+1)-emp(ie))) 
00087  pole_of_chi(ie)=cmplx(x,y) 
00088 enddo 
00089 !
00090 !mat_b
00091 !
00092 allocate(mat_b(ne,ne)); mat_b=0.0d0 
00093 do ie=1,ne
00094  do je=1,ne 
00095   mat_b(ie,je)=1.0d0/(em(ie)-pole_of_chi(je))-1.0d0/(em(ie)+pole_of_chi(je)) 
00096  enddo 
00097 enddo 
00098 !
00099 !check: mat_b * mat_c = 1
00100 !
00101 allocate(mat_c(ne,ne)); mat_c=0.0d0 
00102 mat_c(:,:)=mat_b(:,:)
00103 !
00104 call invmat_complex(ne,mat_b) 
00105 !
00106 do ie=1,ne
00107  do je=1,ne
00108   sum_cmpx=0.0d0
00109   do ke=1,ne
00110    sum_cmpx=sum_cmpx+mat_b(ie,ke)*mat_c(ke,je)
00111   enddo 
00112   !write(6,'(2I5,2F15.10)') ie,je,sum_cmpx 
00113  enddo 
00114 enddo 
00115 end subroutine
00116 !--
00117 subroutine rd_dat_sq 
00118 implicit none 
00119 integer::iq,i,iik,ierr,chdir,jk,iop,iqir,NG_for_eps,NG_for_psi,ig,i1,j1,k1  
00120 character(99)::dirname 
00121 logical::file_e 
00122 real(8)::q1,q2,q3,ktmp(3)  
00123 integer::RWtmp(3) 
00124 integer,allocatable::KGtmp(:,:)!KGtmp(3,NTG)
00125 !
00126 Nq_irr=Nk_irr 
00127 allocate(SQI(3,Nq_irr));SQI(:,:)=0.0D0 
00128 ierr=CHDIR("./dir-eps") 
00129 do iq=1,Nq_irr 
00130  write(dirname,"('q',i3.3)")iq 
00131  ierr=CHDIR(dirname) 
00132  inquire(file='dat.log.400',exist=file_e) 
00133  if(file_e)then 
00134   write(6,*)'dat.log.400 exists in ',trim(dirname) 
00135  else
00136   write(6,*)'error: no dat.log.400 in ',trim(dirname) 
00137  endif 
00138  OPEN(301,FILE='dat.sq') 
00139  rewind(301)
00140  inquire(file='dat.sq',exist=file_e) 
00141  if(file_e)then 
00142   read(301,*)(SQI(i,iq),i=1,3) 
00143  else
00144   write(6,*)'error: no dat.sq in',trim(dirname) 
00145  endif 
00146  close(301) 
00147  ierr=CHDIR("..") 
00148 enddo!iq 
00149 ierr=CHDIR("..") 
00150 !call system('pwd') 
00151 !--
00152 !gen(SQ,numirrq,numrotq,trsq,RWq) 
00153 NTQ=NTK
00154 allocate(SQ(3,NTQ));SQ(:,:)=0.0d0 
00155 allocate(numirrq(NTQ));numirrq(:)=0 
00156 allocate(numrotq(NTQ));numrotq(:)=0 
00157 allocate(trsq(NTQ));trsq(:)=0
00158 allocate(RWq(3,NTQ));RWq(:,:)=0      
00159 !
00160 do iq=1,Nq_irr 
00161  SQ(:,iq)=SQI(:,iq) 
00162  numirrq(iq)=iq
00163  numrotq(iq)=1
00164  trsq(iq)=1
00165  RWq(1:3,iq)=0
00166 enddo 
00167 jk=Nq_irr 
00168 do iq=1,Nq_irr
00169  do iop=1,Nsymq
00170   !sym
00171   ktmp(:)=0.0d0;RWtmp(:)=0  
00172   ktmp(1)=dble(rg(1,1,iop))*SQI(1,iq)+dble(rg(1,2,iop))*SQI(2,iq)+dble(rg(1,3,iop))*SQI(3,iq)
00173   ktmp(2)=dble(rg(2,1,iop))*SQI(1,iq)+dble(rg(2,2,iop))*SQI(2,iq)+dble(rg(2,3,iop))*SQI(3,iq)
00174   ktmp(3)=dble(rg(3,1,iop))*SQI(1,iq)+dble(rg(3,2,iop))*SQI(2,iq)+dble(rg(3,3,iop))*SQI(3,iq)
00175   !
00176   call kcheck(ktmp(1),RWtmp(1))!rewind check 
00177   !
00178   do iik=1,jk 
00179    if(abs(SQ(1,iik)-ktmp(1))<1.0d-4.and.abs(SQ(2,iik)-ktmp(2))<1.0d-4.and.abs(SQ(3,iik)-ktmp(3))<1.0d-4) goto 1100
00180   enddo!iik
00181   jk=jk+1
00182   SQ(:,jk)=ktmp(:)
00183   numirrq(jk)=iq
00184   numrotq(jk)=iop
00185   trsq(jk)=1
00186   RWq(:,jk)=RWtmp(:)
00187 !time-reversal
00188 1100 ktmp(:)=0.0d0;RWtmp(:)=0  
00189   ktmp(1)=dble(rg(1,1,iop))*SQI(1,iq)+dble(rg(1,2,iop))*SQI(2,iq)+dble(rg(1,3,iop))*SQI(3,iq)
00190   ktmp(2)=dble(rg(2,1,iop))*SQI(1,iq)+dble(rg(2,2,iop))*SQI(2,iq)+dble(rg(2,3,iop))*SQI(3,iq) 
00191   ktmp(3)=dble(rg(3,1,iop))*SQI(1,iq)+dble(rg(3,2,iop))*SQI(2,iq)+dble(rg(3,3,iop))*SQI(3,iq) 
00192   !
00193   call kcheck_trs(ktmp(1),RWtmp(1))!rewind check 20170321 
00194   !
00195   do iik=1,jk
00196    if(abs(SQ(1,iik)-(-ktmp(1)))<1.0d-4.and.abs(SQ(2,iik)-(-ktmp(2)))<1.0d-4.and.abs(SQ(3,iik)-(-ktmp(3)))<1.0d-4) goto 2100
00197   enddo!iik
00198   jk=jk+1
00199   SQ(:,jk)=-ktmp(:)
00200   numirrq(jk)=iq
00201   numrotq(jk)=iop
00202   trsq(jk)=-1
00203   RWq(:,jk)=RWtmp(:)
00204 2100 enddo!iop 
00205 enddo!iq  
00206 !--
00207 if(NTQ/=jk)then 
00208  write(6,*)'ERROR;STOP;NTQ should be jk'   
00209  write(6,*)'NTQ=',NTQ,'jk=',jk
00210  stop 
00211 endif 
00212 !--
00213 !gen(NGQ_eps,NGQ_psi,LG0)
00214 allocate(KGtmp(3,NTG));KGtmp(:,:)=0 
00215 allocate(LG0(3,NTG,NTQ));LG0(:,:,:)=0
00216 allocate(NGQ_eps(NTQ));NGQ_eps(:)=0
00217 allocate(NGQ_psi(NTQ));NGQ_psi(:)=0
00218 !
00219 do iq=1,Nq_irr 
00220  q1=SQI(1,iq); q2=SQI(2,iq); q3=SQI(3,iq)
00221  call make_LG0(NTG,b1(1),b2(1),b3(1),Ecut_for_eps,Ecut_for_psi,q1,q2,q3,LG0(1,1,iq),NG_for_eps,NG_for_psi) 
00222  NGQ_eps(iq)=NG_for_eps
00223  NGQ_psi(iq)=NG_for_psi  
00224  write(6,'(i8,3f10.5,a8,i8,a8,i10)') iq,q1,q2,q3,'NGeps',NG_for_eps,'NGpsi',NG_for_psi  
00225 enddo 
00226 do iq=Nq_irr+1,NTQ 
00227  if(trsq(iq)==1)then 
00228   iqir=numirrq(iq)
00229   iop=numrotq(iq) 
00230   !
00231   !q1=dble(rg(1,1,iop))*SQI(1,iqir)+dble(rg(1,2,iop))*SQI(2,iqir)+dble(rg(1,3,iop))*SQI(3,iqir)+dble(RWq(1,iq)) 
00232   !q2=dble(rg(2,1,iop))*SQI(1,iqir)+dble(rg(2,2,iop))*SQI(2,iqir)+dble(rg(2,3,iop))*SQI(3,iqir)+dble(RWq(2,iq))  
00233   !q3=dble(rg(3,1,iop))*SQI(1,iqir)+dble(rg(3,2,iop))*SQI(2,iqir)+dble(rg(3,3,iop))*SQI(3,iqir)+dble(RWq(3,iq))  
00234   !
00235   q1=rg(1,1,iop)*SQI(1,iqir)+rg(1,2,iop)*SQI(2,iqir)+rg(1,3,iop)*SQI(3,iqir)+dble(RWq(1,iq)) 
00236   q2=rg(2,1,iop)*SQI(1,iqir)+rg(2,2,iop)*SQI(2,iqir)+rg(2,3,iop)*SQI(3,iqir)+dble(RWq(2,iq))  
00237   q3=rg(3,1,iop)*SQI(1,iqir)+rg(3,2,iop)*SQI(2,iqir)+rg(3,3,iop)*SQI(3,iqir)+dble(RWq(3,iq))  
00238   !
00239   call make_LG0(NTG,b1(1),b2(1),b3(1),Ecut_for_eps,Ecut_for_psi,q1,q2,q3,LG0(1,1,iq),NG_for_eps,NG_for_psi) 
00240   NGQ_eps(iq)=NG_for_eps
00241   NGQ_psi(iq)=NG_for_psi  
00242   write(6,'(i8,3f10.5,a8,i8,a8,i10)') iq,q1,q2,q3,'NGeps',NG_for_eps,'NGpsi',NG_for_psi  
00243  elseif(trsq(iq)==-1)then  
00244   iqir=numirrq(iq)
00245   iop=numrotq(iq) 
00246   !
00247   !q1=dble(rg(1,1,iop))*SQI(1,iqir)+dble(rg(1,2,iop))*SQI(2,iqir)+dble(rg(1,3,iop))*SQI(3,iqir)+dble(RWq(1,iq)) 
00248   !q2=dble(rg(2,1,iop))*SQI(1,iqir)+dble(rg(2,2,iop))*SQI(2,iqir)+dble(rg(2,3,iop))*SQI(3,iqir)+dble(RWq(2,iq))  
00249   !q3=dble(rg(3,1,iop))*SQI(1,iqir)+dble(rg(3,2,iop))*SQI(2,iqir)+dble(rg(3,3,iop))*SQI(3,iqir)+dble(RWq(3,iq))  
00250   !
00251   q1=rg(1,1,iop)*SQI(1,iqir)+rg(1,2,iop)*SQI(2,iqir)+rg(1,3,iop)*SQI(3,iqir)+dble(RWq(1,iq)) 
00252   q2=rg(2,1,iop)*SQI(1,iqir)+rg(2,2,iop)*SQI(2,iqir)+rg(2,3,iop)*SQI(3,iqir)+dble(RWq(2,iq))  
00253   q3=rg(3,1,iop)*SQI(1,iqir)+rg(3,2,iop)*SQI(2,iqir)+rg(3,3,iop)*SQI(3,iqir)+dble(RWq(3,iq))  
00254   !
00255   KGtmp(:,:)=0 
00256   call make_LG0(NTG,b1(1),b2(1),b3(1),Ecut_for_eps,Ecut_for_psi,q1,q2,q3,KGtmp(1,1),NG_for_eps,NG_for_psi) 
00257   NGQ_eps(iq)=NG_for_eps 
00258   NGQ_psi(iq)=NG_for_psi  
00259   write(6,'(i8,3f10.5,a8,i8,a8,i10)') iq,q1,q2,q3,'NGeps',NG_for_eps,'NGpsi',NG_for_psi  
00260   LG0(:,:,iq)=-KGtmp(:,:) 
00261  endif 
00262 enddo 
00263 !--
00264 Lq1=maxval(abs(LG0(1,:,:)))+1
00265 Lq2=maxval(abs(LG0(2,:,:)))+1
00266 Lq3=maxval(abs(LG0(3,:,:)))+1
00267 allocate(packingq(-Lq1:Lq1,-Lq2:Lq2,-Lq3:Lq3,Nq_irr));packingq(:,:,:,:)=0 
00268 do iq=1,Nq_irr 
00269  do ig=1,NGQ_eps(iq) 
00270   i1=LG0(1,ig,iq)
00271   j1=LG0(2,ig,iq)
00272   k1=LG0(3,ig,iq) 
00273   packingq(i1,j1,k1,iq)=ig 
00274   !write(6,'(a,5i10)')'iq,ig,i1,j1,k1',iq,ig,i1,j1,k1
00275  enddo 
00276 enddo 
00277 !-- 
00278 NTGQ=maxval(NGQ_eps(:))
00279 write(6,*)'Lq1=',Lq1 
00280 write(6,*)'Lq2=',Lq2 
00281 write(6,*)'Lq3=',Lq3 
00282 write(6,*)'NTGQ=',NTGQ  
00283 !--
00284 do iq=1,NTQ 
00285  write(6,'(i5,3f15.10,6i5)')iq,(SQ(i1,iq),i1=1,3),numirrq(iq),numrotq(iq),trsq(iq),(RWq(j1,iq),j1=1,3) 
00286 enddo  
00287 !--
00288 end subroutine
00289 !--
00290 subroutine rd_dat_eps 
00291 implicit none 
00292 integer::ierr,chdir,iq,iqgm,ix,file_num,NG_for_eps,ig,jg,ie
00293 character(99)::filename,dirname 
00294 complex(8),allocatable::epstmp(:,:,:)!epstmp(NTGQ,NTGQ,ne)
00295 complex(8),allocatable::epstmpgm(:,:,:,:)!epstmpgm(NTGQ,NTGQ,ne,3)
00296 !--
00297 !OPEN(600-,R,FILE='dat.epsqw',FORM='unformatted') 
00298 allocate(epstmp(NTGQ,NTGQ,ne));epstmp(:,:,:)=0.0d0!real8
00299 allocate(epstmpgm(NTGQ,NTGQ,ne,3));epstmpgm(:,:,:,:)=0.0d0!real8
00300 allocate(epsirr(NTGQ,NTGQ,ne,Nq_irr));epsirr(:,:,:,:)=0.0d0!real4 
00301 !
00302 ierr=CHDIR("./dir-eps") 
00303 call system('pwd') 
00304 do iq=1,Nq_irr 
00305  if(abs(SKI(1,iq))<1.0d-5.and.abs(SKI(2,iq))<1.0d-5.and.abs(SKI(3,iq))<1.0d-5)then 
00306   write(dirname,"('q',i3.3)")iq
00307   iqgm=iq 
00308   ierr=CHDIR(dirname) 
00309   do ix=1,3 
00310    file_num=600+(ix-1)  
00311    write(filename,'("dat.epsqw.",i3.3)')file_num 
00312    open(file_num,file=filename,form='unformatted') 
00313    write(6,'(a10,a15,a5,a10)')'read: ',trim(filename),'in ',trim(dirname)  
00314    rewind(file_num) 
00315    NG_for_eps=NGQ_eps(iq)
00316    read(file_num)(((epstmpgm(ig,jg,ie,ix),ig=1,NG_for_eps),jg=1,NG_for_eps),ie=1,ne)
00317   enddo!ix 
00318   !epsirr: real4
00319   !epstmp: real8 
00320   epsirr(:,:,:,iq)=(epstmpgm(:,:,:,1)+epstmpgm(:,:,:,2)+epstmpgm(:,:,:,3))/3.0d0
00321  else
00322   write(dirname,"('q',i3.3)")iq
00323   ierr=CHDIR(dirname) 
00324   file_num=600 
00325   write(filename,'("dat.epsqw.",i3.3)')file_num 
00326   open(file_num,file=filename,form='unformatted') 
00327   write(6,'(a10,a15,a5,a10)')'read: ',trim(filename),'in ',trim(dirname)  
00328   rewind(file_num) 
00329   NG_for_eps=NGQ_eps(iq)
00330   read(file_num)(((epstmp(ig,jg,ie),ig=1,NG_for_eps),jg=1,NG_for_eps),ie=1,ne)
00331   !epsirr: real4 
00332   !epstmp: real8 
00333   epsirr(:,:,:,iq)=epstmp(:,:,:)
00334  endif!q=0 or not 
00335  ierr=CHDIR("..") 
00336 enddo!iq   
00337 ierr=CHDIR("..") 
00338 call system('pwd') 
00339 !
00340 !do iq=1,Nq_irr 
00341 ! NG_for_eps=NGQ_eps(iq)
00342 ! write(6,*) NG_for_eps 
00343 ! do ig=1,NG_for_eps 
00344 !  write(6,*) epsirr(ig,ig,100,iq) 
00345 ! enddo  
00346 !enddo  
00347 !
00348 !do iq=1,Nq_irr 
00349 ! NG_for_eps=NGQ_eps(iq) 
00350 ! do ie=1,nen 
00351 !  do igL=1,NG_for_eps 
00352 !   epsirr(igL,igL,ie,iq)=epsirr(igL,igL,ie,iq)-1.0d0    
00353 !  enddo 
00354 ! enddo 
00355 !enddo 
00356 !write(6,*) 'FINISH READING EPSILON'
00357 !
00358 end subroutine
00359 !--
00360 end module 

Generated on 17 Nov 2020 for respack by  doxygen 1.6.1