TABLE OF CONTENTS


ABINIT/local_ks_green [ Functions ]

[ Top ] [ Functions ]

NAME

 local_ks_green

FUNCTION

 Compute the sum over k-point of ks green function.
 do the fourier transformation and print it

COPYRIGHT

 Copyright (C) 1999-2018 ABINIT group (BAmadon)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  cryst_struc
  istep    =  step of iteration for LDA.
  lda_occup
  mpi_enreg=informations about MPI parallelization
  paw_dmft =  data for self-consistent LDA+DMFT calculations.

OUTPUT

  paw_dmft =  data for self-consistent LDA+DMFT calculations.

NOTES

PARENTS

      dmft_solve

CHILDREN

      fourier_fct,wrtout

SOURCE

 37 #if defined HAVE_CONFIG_H
 38 #include "config.h"
 39 #endif
 40 
 41 
 42 #include "abi_common.h"
 43 
 44 subroutine local_ks_green(green,paw_dmft,prtopt)
 45 
 46  use defs_basis
 47  use m_errors
 48  use m_profiling_abi
 49 
 50  use m_crystal, only : crystal_t
 51  use m_green, only : green_type,fourier_fct
 52  use m_paw_dmft, only : paw_dmft_type
 53 
 54 !This section has been created automatically by the script Abilint (TD).
 55 !Do not modify the following lines by hand.
 56 #undef ABI_FUNC
 57 #define ABI_FUNC 'local_ks_green'
 58  use interfaces_14_hidewrite
 59 !End of the abilint section
 60 
 61  implicit none
 62 
 63 !Arguments ------------------------------------
 64 !scalars
 65  type(green_type), intent(in) :: green
 66  type(paw_dmft_type), intent(in)  :: paw_dmft
 67  integer, intent(in) :: prtopt
 68 
 69 !Local variables ------------------------------
 70  character(len=500) :: message
 71  integer :: iband,ifreq,ikpt,isppol,itau,lsub,ltau,mbandc,nkpt,nsppol
 72  character(len=1) :: tag_is
 73  character(len=fnlen) :: tmpfil
 74  integer,allocatable :: unitgreenlocks_arr(:)
 75  real(dp) :: beta
 76  real(dp), allocatable :: tau(:)
 77  complex(dpc), allocatable :: loc_ks(:,:,:)
 78  complex(dpc), allocatable :: loc_ks_tau(:,:,:),fw(:),ft(:)
 79 !scalars
 80 !************************************************************************
 81  mbandc=paw_dmft%mbandc
 82  nkpt=paw_dmft%nkpt
 83  nsppol=paw_dmft%nsppol
 84  ltau=128
 85  ABI_ALLOCATE(tau,(ltau))
 86  do itau=1,ltau
 87    tau(itau)=float(itau-1)/float(ltau)/paw_dmft%temp
 88  end do
 89  beta=one/paw_dmft%temp
 90 
 91 !Only imaginary frequencies here
 92  if(green%w_type=="real") then
 93    message = ' compute_energy not implemented for real frequency'
 94    MSG_BUG(message)
 95  end if
 96 
 97 !=========================================
 98 !Compute local band ks green function
 99 ! should be computed in compute_green: it would be less costly in memory.
100 !=========================================
101  ABI_ALLOCATE(loc_ks,(nsppol,mbandc,paw_dmft%dmft_nwlo))
102  if(green%oper(1)%has_operks==1) then
103    loc_ks(:,:,:)=czero
104    do isppol=1,nsppol
105      do iband=1,mbandc
106        do ifreq=1,paw_dmft%dmft_nwlo
107          do ikpt=1,nkpt
108            loc_ks(isppol,iband,ifreq)=loc_ks(isppol,iband,ifreq)+  &
109 &           green%oper(ifreq)%ks(isppol,ikpt,iband,iband)*paw_dmft%wtk(ikpt)
110          end do
111        end do
112      end do
113    end do
114  else
115    message = ' green fct is not computed in ks space'
116    MSG_BUG(message)
117  end if
118 
119 !=========================================
120 !Compute fourier transformation 
121 !=========================================
122 
123  ABI_ALLOCATE(loc_ks_tau,(nsppol,mbandc,ltau))
124  ABI_ALLOCATE(fw,(paw_dmft%dmft_nwlo))
125  ABI_ALLOCATE(ft,(ltau))
126  loc_ks_tau(:,:,:)=czero
127  do isppol=1,nsppol
128    do iband=1,mbandc
129      do ifreq=1,paw_dmft%dmft_nwlo
130        fw(ifreq)=loc_ks(isppol,iband,ifreq)
131      end do
132      call fourier_fct(fw,ft,.true.,ltau,-1,paw_dmft) ! inverse fourier
133      do itau=1,ltau
134        loc_ks_tau(isppol,iband,itau)=ft(itau)
135      end do
136    end do
137  end do
138  ABI_DEALLOCATE(fw)
139  ABI_DEALLOCATE(ft)
140  do isppol=1,nsppol
141    do iband=1,mbandc
142      do itau=1,ltau
143        loc_ks_tau(isppol,iband,itau)=(loc_ks_tau(isppol,iband,itau)+conjg(loc_ks_tau(isppol,iband,itau)))/two
144      end do
145    end do
146  end do
147 
148 !=========================================
149 !Print out ksloc green function
150 !=========================================
151  if(abs(prtopt)==1) then
152    ABI_ALLOCATE(unitgreenlocks_arr,(nsppol))
153    do isppol=1,nsppol
154      write(tag_is,'(i1)')isppol
155      tmpfil = trim(paw_dmft%filapp)//'Gtau_locks_isppol'//tag_is
156      write(message,'(3a)') ch10," == Print green function on file ",tmpfil
157      call wrtout(std_out,message,'COLL')
158      unitgreenlocks_arr(isppol)=500+isppol-1
159      open (unit=unitgreenlocks_arr(isppol),file=trim(tmpfil),status='unknown',form='formatted')
160      rewind(unitgreenlocks_arr(isppol))
161      write(message,'(a,a,a,i4)') 'opened file : ', trim(tmpfil), ' unit', unitgreenlocks_arr(isppol)
162      write(message,'(a,a)') ch10,"# New record : First 40 bands"
163      call wrtout(unitgreenlocks_arr(isppol),message,'COLL')
164      do lsub=1,mbandc/40+1
165        do itau=1, ltau
166          write(message,'(2x,50(e10.3,2x))') tau(itau), &
167 &         (real(loc_ks_tau(isppol,iband,itau)),iband=40*(lsub-1)+1,min(40*lsub,mbandc))
168          call wrtout(unitgreenlocks_arr(isppol),message,'COLL')
169        end do
170        write(message,'(2x,50(e10.3,2x))') beta, &
171 &       ((-one-real(loc_ks_tau(isppol,iband,1))),iband=40*(lsub-1)+1,min(40*lsub,mbandc))
172        call wrtout(unitgreenlocks_arr(isppol),message,'COLL')
173        if(40*lsub<mbandc) then
174          write(message,'(a,a,i5,a,i5)')    &
175 &         ch10,"# Same record, Following bands : From ",    &
176 &         40*(lsub),"  to ",min(40*(lsub+1),mbandc)
177          call wrtout(unitgreenlocks_arr(isppol),message,'COLL')
178        end if
179      end do
180 !    call flush(unitgreenlocks_arr(isppol))
181    end do
182    ABI_DEALLOCATE(unitgreenlocks_arr)
183  end if
184 
185 !Deallocations
186  ABI_DEALLOCATE(loc_ks)
187  ABI_DEALLOCATE(loc_ks_tau)
188  ABI_DEALLOCATE(tau)
189 
190 end subroutine local_ks_green