TABLE OF CONTENTS


ABINIT/m_read_plowannier [ Modules ]

[ Top ] [ Modules ]

NAME

  m_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-2024 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

OUTPUT

SOURCE

21 #if defined HAVE_CONFIG_H
22 #include "config.h"
23 #endif
24 
25 #include "abi_common.h"
26 
27 MODULE m_read_plowannier
28 
29  use defs_basis
30  use m_abicore
31  use m_errors
32 
33  use m_io_tools,      only : open_file
34  use m_crystal,       only : crystal_t
35  use m_bz_mesh,       only : kmesh_t
36  use m_pawang,        only : pawang_type
37 
38  implicit none
39 
40  private
41 
42  public :: read_plowannier

m_read_plowannier/read_plowannier [ Modules ]

[ Top ] [ m_read_plowannier ] [ 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-2024 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

SOURCE

 98 subroutine read_plowannier(cryst,bandinf,bandsup,coeffW_BZ,itypatcor,Kmesh,lcor,luwindow,nspinor,nsppol,pawang,prtvol,ucrpa_bands)
 99 
100 !Arguments ------------------------------------
101 !types and arrays
102  type(kmesh_t),intent(in) :: Kmesh
103  type(crystal_t),intent(in) :: Cryst
104  complex(dpc), allocatable, intent(inout) :: coeffW_BZ(:,:,:,:,:,:)
105  type(Pawang_type),intent(in) :: Pawang
106 !scalars
107  logical, intent(inout) :: luwindow
108  integer, intent(out) :: bandinf,bandsup,itypatcor,lcor
109  integer, intent(in) :: nspinor,nsppol,prtvol
110  integer, intent(in) :: ucrpa_bands(2)
111 
112 !Local variables-------------------------------
113  character(len=500) :: message,msg
114  integer :: at_indx,ik_ibz,band1,m1,m2,spin,ik_bz,dummy,isym,itim,iat,indx,ispinor,unt
115  real(dp) :: xx,yy
116  real(dp) :: kbz(3)
117  complex(dpc), allocatable :: coeffW_IBZ(:,:,:,:,:,:)
118 ! *********************************************************************
119  write(message,*) "Read wannier in iBZ"
120  call wrtout(std_out,message,'COLL')
121 
122  if (open_file('forlb.ovlp',msg,newunit=unt,form='formatted',status='unknown') /= 0) then
123    ABI_ERROR(msg)
124  end if
125  rewind(unt)
126  read(unt,*) message
127  read(unt,*) message, lcor,itypatcor
128  read(unt,*) message, bandinf,bandsup
129  write(std_out,*) 'read from forlb.ovlp',lcor, bandinf,bandsup
130  write(std_out,*) "for wannier", bandinf,bandsup
131  if(prtvol>0) then
132  endif
133 
134  if(.not.luwindow.and.(ucrpa_bands(1)/=bandinf.or.ucrpa_bands(2)/=bandsup)) then
135    write(msg,'(a,a)')' Bands used for Wannier construction and cRPA construction differ',&
136 & 'It might be physically correct and it is possible with the current implementation'
137    call wrtout(std_out,msg,'COLL')
138    call wrtout(ab_out,msg,'COLL')
139 !   ABI_ERROR(msg)
140  end if
141 
142 !Do not dead the bandinf, bandinf redondance information
143 
144  ABI_MALLOC(coeffW_IBZ,(Cryst%nattyp(itypatcor),nsppol,bandinf:bandsup,Kmesh%nibz,nspinor,2*lcor+1))
145  coeffW_IBZ=czero
146  do spin=1,nsppol
147    do ik_ibz=1,Kmesh%nibz
148      !read k
149       read(unt,*)
150       do band1=bandinf,bandsup
151      !read band
152         read(unt,*)
153      !read projection
154         do ispinor=1,nspinor
155           do iat=1,Cryst%nattyp(itypatcor)
156             do m1=1,2*lcor+1
157               read(unt,*) dummy,dummy,dummy,dummy,xx,yy
158               coeffW_IBZ(iat,spin,band1,ik_ibz,ispinor,m1)=cmplx(xx,yy)
159             end do
160           end do
161         end do
162       end do
163    end do
164  end do
165  close(unt)
166 
167  ABI_MALLOC(coeffW_BZ,(Cryst%nattyp(itypatcor),nsppol,bandinf:bandsup,Kmesh%nbz,nspinor,2*lcor+1))
168  coeffW_BZ=czero
169 
170  if (Kmesh%nbz==Kmesh%nibz) then
171    coeffW_BZ=coeffW_IBZ
172  else if (Cryst%nsym==1) then
173    write(message,*) "Reconstruct in full BZ"
174    call wrtout(std_out,message,'COLL')
175    call wrtout(ab_out,message,'COLL')
176    do ik_bz=1,Kmesh%nbz
177 !     if(prtvol>=10) write(6,*) "ik",ik_bz,Kmesh%tab(ik_bz),Kmesh%tabi(ik_bz),Kmesh%tabo(ik_bz)
178      if(Kmesh%tabi(ik_bz)==1) then
179        coeffW_BZ(:,:,:,ik_bz,:,:)=coeffW_IBZ(:,:,:,Kmesh%tab(ik_bz),:,:)
180      else if(Kmesh%tabi(ik_bz)==-1) then
181        coeffW_BZ(:,:,:,ik_bz,:,:)=conjg(coeffW_IBZ(:,:,:,Kmesh%tab(ik_bz),:,:))
182      endif
183 !     if(prtvol>=10)write(6,*) "coeffW 21",coeffW_BZ(1,1,bandinf,ik_bz,:)
184 !     if(prtvol>=10)write(6,*) "coeffW 21",coeffW_BZ(1,1,bandinf+1,ik_bz,:)
185 !     if(prtvol>=10)write(6,*) "coeffW 25",coeffW_BZ(1,1,bandsup,ik_bz,:)
186    enddo
187  else if (Cryst%nsym>1) then
188    write(message,*) "Reconstruct in full BZ (nsym=0)"
189    call wrtout(std_out,message,'COLL')
190    do ik_bz=1,Kmesh%nbz
191    !write(6,*) ik_bz,Kmesh%nbz
192      call kmesh%get_BZ_item(ik_bz,kbz,ik_ibz,isym,itim)
193      do indx=1,cryst%nattyp(itypatcor)
194        iat=cryst%atindx1(indx) ! correct index for the full atom list
195        at_indx=cryst%indsym(4,isym,iat) !! see eg sym_matlu and m_crystal
196        do spin=1,nsppol
197          do ispinor=1,nspinor
198            do m1=1,2*lcor+1
199              do m2=1,2*lcor+1
200              coeffW_BZ(indx,spin,:,ik_bz,ispinor,m1)= coeffW_BZ(indx,spin,:,ik_bz,ispinor,m1)&
201 &                     +coeffW_IBZ(at_indx,spin,:,ik_ibz,ispinor,m2)*pawang%zarot(m2,m1,lcor+1,isym)
202              enddo
203 !             write(6,'(20f7.3)') (pawang%zarot(m1,m2,lcor+1,isym),m2=1,2*lcor+1)
204            enddo
205          enddo
206        enddo
207      enddo
208 !     do m1=1,2*lcor+1
209 !      write(6,*) "coeffW_IBZ",coeffW_IBZ(1,8,ik_ibz,m1)
210 !     enddo
211 !     do m1=1,2*lcor+1
212 !      write(6,*) "coeffW_ BZ",coeffW_BZ(1,8,ik_bz,m1)
213 !     enddo
214    enddo
215  end if
216  ABI_FREE(coeffW_IBZ)
217 
218 
219 end subroutine read_plowannier