TABLE OF CONTENTS


ABINIT/read_plowannier [ Modules ]

[ Top ] [ Modules ]

NAME

 read_plowannier

FUNCTION

  Read Wannier coefficient in the file forlb.ovlp for ucrpa calculation
  this file was typically created in a DFT run with usedmft=1 and nbandkss -1

COPYRIGHT

 Copyright (C) 2006-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 .

INPUTS

 Cryst<cryst_t>= data type gathering info on symmetries and unit cell
    %natom=number of atoms
    %nsym=number of symmetries
    %xred(3,natom)=reduced coordinated of atoms
    %typat(natom)=type of each atom
    %rprimd(3,3)=dimensional primitive translations in real space (bohr)
    %timrev= 2 if time reversal can be used, 1 otherwise
 Kmesh <kmesh_t>= datatype gathering parameters related to the k-point sampling
    %nibz=number of k-points in the IBZ
    %nbz=number of k-points in the BZ
    %bz(3,nbz)=reduced coordinates for k-points in the full Brillouin zone
    %ibz(3,nibz)=reduced coordinates for k-points in the irreducible wedge
    %tab(nbz)=mapping between a kpt in the BZ (array bz) and the irred point in the array ibz
    %tabi(nbz)= -1 if inversion is needed to obtain this particular kpt in the BZ, 1 means identity
    %tabo(nbz)= for each point in the BZ, the index of the symmetry operation S in reciprocal
      space which rotates k_IBZ onto \pm k_BZ (depending on tabi)
    %tabp(nbz)= For each k_BZ, it gives the phase factors associated to non-symmorphic operations, i.e
      e^{-i 2 \pi k_IBZ \cdot R{^-1}t} == e{-i 2\pi k_BZ cdot t} where :
      \transpose R{-1}=S and (S k_IBZ) = \pm k_BZ (depending on ktabi)
    %tabr(nfftot,nbz) For each point r on the real mesh and for each k-point in the BZ, tabr
      gives the index of (R^-1 (r-t)) in the FFT array where R=\transpose S^{-1} and k_BZ=S k_IBZ.
      t is the fractional translation associated to R
  luwindow: T if ucrpa_window is activated, F if not.
  prtvol: integer to give the amount of printing.
 nsppol : number of spin polarization.
 Pawang<pawang_type> angular mesh discretization and related data:

OUTPUT

 bandinf, bandsup : lower and upper bands for define Wannier functions
 coeffW_BZ(nsppol,bandinf:bandsup,nvz,2*lcor+1)
 lcor : angular momentum for correlated orbitals
 itypatcor : correlated species

PARENTS

      cchi0,cchi0q0,prep_calc_ucrpa

CHILDREN

      get_bz_item,wrtout

SOURCE

 59 #if defined HAVE_CONFIG_H
 60 #include "config.h"
 61 #endif
 62 
 63 #include "abi_common.h"
 64 
 65 
 66 subroutine read_plowannier(cryst,bandinf,bandsup,coeffW_BZ,itypatcor,Kmesh,lcor,luwindow,nspinor,nsppol,pawang,prtvol,ucrpa_bands)
 67 
 68  use defs_basis
 69  use defs_datatypes
 70  use defs_abitypes
 71  use m_profiling_abi
 72  use m_errors
 73 
 74  use m_io_tools,      only : open_file
 75  use m_crystal,       only : crystal_t
 76  use m_bz_mesh,       only : kmesh_t, get_BZ_item
 77  use m_pawang,        only : pawang_type
 78 
 79 !This section has been created automatically by the script Abilint (TD).
 80 !Do not modify the following lines by hand.
 81 #undef ABI_FUNC
 82 #define ABI_FUNC 'read_plowannier'
 83  use interfaces_14_hidewrite
 84 !End of the abilint section
 85 
 86  implicit none
 87 !Arguments ------------------------------------
 88 !types and arrays
 89  type(kmesh_t),intent(in) :: Kmesh
 90  type(crystal_t),intent(in) :: Cryst
 91  complex(dpc), allocatable, intent(inout) :: coeffW_BZ(:,:,:,:,:,:)
 92  type(Pawang_type),intent(in) :: Pawang
 93 !scalars
 94  logical, intent(inout) :: luwindow
 95  integer, intent(out) :: bandinf,bandsup,itypatcor,lcor
 96  integer, intent(in) :: nspinor,nsppol,prtvol
 97  integer, intent(in) :: ucrpa_bands(2)
 98 
 99 !Local variables-------------------------------
100  character(len=500) :: message,msg
101  integer :: at_indx,ik_ibz,band1,m1,m2,spin,ik_bz,dummy,isym,itim,iat,indx,ispinor,unt
102  real(dp) :: xx,yy 
103  real(dp) :: kbz(3)
104  complex(dpc), allocatable :: coeffW_IBZ(:,:,:,:,:,:)
105 ! *********************************************************************
106  write(message,*) "Read wannier in iBZ"
107  call wrtout(std_out,message,'COLL')
108 
109  if (open_file('forlb.ovlp',msg,newunit=unt,form='formatted',status='unknown') /= 0) then
110    MSG_ERROR(msg)
111  end if
112  rewind(unt)
113  read(unt,*) message, lcor,itypatcor
114  read(unt,*) message, bandinf,bandsup
115  write(std_out,*) 'read from forlb.ovlp',lcor, bandinf,bandsup
116  write(std_out,*) "for wannier", bandinf,bandsup
117  if(prtvol>0) then
118  endif
119 
120  if(.not.luwindow.and.(ucrpa_bands(1)/=bandinf.or.ucrpa_bands(2)/=bandsup)) then
121    write(msg,'(a,a)')' Bands used for Wannier construction and cRPA construction differ',&
122 & 'It might be physically correct and it is possible with the current implementation'
123    call wrtout(std_out,msg,'COLL')
124    call wrtout(ab_out,msg,'COLL')
125 !   MSG_ERROR(msg)
126  end if
127 
128 !Do not dead the bandinf, bandinf redondance information
129 
130  ABI_ALLOCATE(coeffW_IBZ,(Cryst%nattyp(itypatcor),nsppol,bandinf:bandsup,Kmesh%nibz,nspinor,2*lcor+1))
131  coeffW_IBZ=czero
132  do spin=1,nsppol
133    do ik_ibz=1,Kmesh%nibz
134      !read k
135       read(unt,*)
136       do band1=bandinf,bandsup
137      !read band
138         read(unt,*)
139      !read projection
140         do ispinor=1,nspinor
141           do iat=1,Cryst%nattyp(itypatcor)
142             do m1=1,2*lcor+1
143               read(unt,*) dummy,dummy,dummy,xx,yy
144               coeffW_IBZ(iat,spin,band1,ik_ibz,ispinor,m1)=cmplx(xx,yy)
145             end do
146           end do
147         end do
148       end do
149    end do
150  end do
151  close(unt)
152 
153  ABI_ALLOCATE(coeffW_BZ,(Cryst%nattyp(itypatcor),nsppol,bandinf:bandsup,Kmesh%nbz,nspinor,2*lcor+1))
154  coeffW_BZ=czero
155 
156  if (Kmesh%nbz==Kmesh%nibz) then
157    coeffW_BZ=coeffW_IBZ
158  else if (Cryst%nsym==1) then
159    write(message,*) "Reconstruct in full BZ"
160    call wrtout(std_out,message,'COLL')
161    call wrtout(ab_out,message,'COLL')
162    do ik_bz=1,Kmesh%nbz
163 !     if(prtvol>=10) write(6,*) "ik",ik_bz,Kmesh%tab(ik_bz),Kmesh%tabi(ik_bz),Kmesh%tabo(ik_bz)
164      if(Kmesh%tabi(ik_bz)==1) then
165        coeffW_BZ(:,:,:,ik_bz,:,:)=coeffW_IBZ(:,:,:,Kmesh%tab(ik_bz),:,:)
166      else if(Kmesh%tabi(ik_bz)==-1) then
167        coeffW_BZ(:,:,:,ik_bz,:,:)=conjg(coeffW_IBZ(:,:,:,Kmesh%tab(ik_bz),:,:))
168      endif
169 !     if(prtvol>=10)write(6,*) "coeffW 21",coeffW_BZ(1,1,bandinf,ik_bz,:)
170 !     if(prtvol>=10)write(6,*) "coeffW 21",coeffW_BZ(1,1,bandinf+1,ik_bz,:)
171 !     if(prtvol>=10)write(6,*) "coeffW 25",coeffW_BZ(1,1,bandsup,ik_bz,:)
172    enddo
173  else if (Cryst%nsym>1) then
174    write(message,*) "Reconstruct in full BZ (nsym=0)"
175    call wrtout(std_out,message,'COLL')
176    do ik_bz=1,Kmesh%nbz
177    !write(6,*) ik_bz,Kmesh%nbz
178      call get_BZ_item(Kmesh,ik_bz,kbz,ik_ibz,isym,itim)
179      do indx=1,cryst%nattyp(itypatcor)
180        iat=cryst%atindx1(indx) ! correct index for the full atom list
181        at_indx=cryst%indsym(4,isym,iat) !! see eg sym_matlu and m_crystal
182        do spin=1,nsppol
183          do ispinor=1,nspinor
184            do m1=1,2*lcor+1
185              do m2=1,2*lcor+1
186              coeffW_BZ(indx,spin,:,ik_bz,ispinor,m1)= coeffW_BZ(indx,spin,:,ik_bz,ispinor,m1)&
187 &                     +coeffW_IBZ(at_indx,spin,:,ik_ibz,ispinor,m2)*pawang%zarot(m2,m1,lcor+1,isym)
188              enddo
189 !             write(6,'(20f7.3)') (pawang%zarot(m1,m2,lcor+1,isym),m2=1,2*lcor+1)
190            enddo
191          enddo
192        enddo
193      enddo
194 !     do m1=1,2*lcor+1
195 !      write(6,*) "coeffW_IBZ",coeffW_IBZ(1,8,ik_ibz,m1)
196 !     enddo
197 !     do m1=1,2*lcor+1
198 !      write(6,*) "coeffW_ BZ",coeffW_BZ(1,8,ik_bz,m1)
199 !     enddo
200    enddo
201  end if
202  ABI_DEALLOCATE(coeffW_IBZ)
203 
204 
205 end subroutine read_plowannier