TABLE OF CONTENTS


ABINIT/update_cprj [ Functions ]

[ Top ] [ Functions ]

NAME

 update_cprj

FUNCTION

  Update the matrix elements of the PAW projectors in case of self-consistent GW.

COPYRIGHT

  Copyright (C) 2008-2018 ABINIT group (MG)
  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

  dimlmn(natom)=number of (l,m,n) components for each atom (only for PAW)
  nkibz=number of k-points
  nsppol=number of spin
  nbnds=number of bands in the present GW calculation
  m_lda_to_qp(nbnds,nbnds,nkibz,nsppol)= expansion of the QP amplitudes in terms of KS wavefunctions
  natom=number of atomd in unit cell

OUTPUT

  Cprj_ibz(natom,nspinor*nkibz*nbnds*nsppol) <type(pawcprj_type)>=projected wave functions 
   <Proj_i|Cnk> with all NL projectors. On exit, it contains the projections onto the 
   QP amplitudes.

TODO

 To be moved to cprj_utils, although here we use complex variables.

PARENTS

      mlwfovlp_qp

CHILDREN

      pawcprj_alloc,pawcprj_free

SOURCE

 40 #if defined HAVE_CONFIG_H
 41 #include "config.h"
 42 #endif
 43 
 44 #include "abi_common.h"
 45 
 46 subroutine update_cprj(natom,nkibz,nbnds,nsppol,nspinor,m_lda_to_qp,dimlmn,Cprj_ibz)
 47 
 48  use defs_basis
 49  use m_profiling_abi
 50  use m_errors
 51  use m_pawcprj, only : pawcprj_type, pawcprj_alloc, pawcprj_free
 52 
 53 !This section has been created automatically by the script Abilint (TD).
 54 !Do not modify the following lines by hand.
 55 #undef ABI_FUNC
 56 #define ABI_FUNC 'update_cprj'
 57 !End of the abilint section
 58 
 59  implicit none
 60 
 61 !Arguments ------------------------------------
 62 !scalars
 63  integer,intent(in) :: natom,nbnds,nkibz,nsppol,nspinor
 64 !arrays
 65  integer,intent(in) :: dimlmn(natom)
 66  complex(dpc),intent(in) :: m_lda_to_qp(nbnds,nbnds,nkibz,nsppol)
 67  type(pawcprj_type),intent(inout) :: Cprj_ibz(natom,nspinor*nbnds*nkibz*nsppol)
 68 
 69 !Local variables-------------------------------
 70 !scalars
 71  integer :: iat,ib,ik,is,shift,indx_kibz,ilmn,nlmn,ispinor,ibsp,spad,ibdx
 72 !arrays
 73  real(dp),allocatable :: re_p(:),im_p(:),vect(:,:),umat(:,:,:)
 74  type(pawcprj_type),allocatable :: Cprj_ks(:,:)
 75 
 76 !************************************************************************
 77 
 78  DBG_ENTER("COLL")
 79 
 80  ABI_DATATYPE_ALLOCATE(Cprj_ks,(natom,nspinor*nbnds))
 81  call pawcprj_alloc(Cprj_ks,0,dimlmn)
 82 
 83  ABI_ALLOCATE(re_p,(nbnds))
 84  ABI_ALLOCATE(im_p,(nbnds))
 85  ABI_ALLOCATE(vect,(2,nbnds))
 86  ABI_ALLOCATE(umat,(2,nbnds,nbnds))
 87  !
 88  ! $ \Psi^{QP}_{r,b} = \sum_n \Psi^{KS}_{r,n} M_{n,b} $ 
 89  !
 90  ! therefore the updated PAW projections are given by:
 91  !
 92  ! $ \<\tprj_j|\Psi^{QP}_a\> = sum_b M_{b,a} <\tprj_j|\Psi^{KS}_b\> $.
 93  !
 94  do is=1,nsppol
 95    do ik=1,nkibz
 96 
 97     shift=nspinor*nbnds*nkibz*(is-1)
 98     indx_kibz=nspinor*nbnds*(ik-1)+shift
 99     ibsp=0
100     do ib=1,nbnds
101       do ispinor=1,nspinor
102         ibsp=ibsp+1
103         do iat=1,natom
104           Cprj_ks(iat,ibsp)%cp(:,:)=Cprj_ibz(iat,indx_kibz+ibsp)%cp(:,:)
105         end do
106       end do
107     end do
108     
109     umat(1,:,:)=TRANSPOSE( REAL (m_lda_to_qp(:,:,ik,is)) )
110     umat(2,:,:)=TRANSPOSE( AIMAG(m_lda_to_qp(:,:,ik,is)) )
111 
112     do iat=1,natom
113       nlmn=dimlmn(iat)
114       do ilmn=1,nlmn
115 
116         do ispinor=1,nspinor
117            ! * Retrieve projections for this spinor component, at fixed atom and ilmn.
118            spad=(ispinor-1)
119            ibdx=0
120            do ib=1,nbnds*nspinor,nspinor 
121             ibdx=ibdx+1
122             vect(1,ibdx)=Cprj_ks(iat,ib+spad)%cp(1,ilmn)
123             vect(2,ibdx)=Cprj_ks(iat,ib+spad)%cp(2,ilmn)
124            end do
125 
126            re_p(:)= &
127 &            MATMUL(umat(1,:,:),vect(1,:)) &
128 &           -MATMUL(umat(2,:,:),vect(2,:))
129 
130            im_p(:)= &
131 &            MATMUL(umat(1,:,:),vect(2,:)) &
132 &           +MATMUL(umat(2,:,:),vect(1,:))
133            
134            ! === Save values ===
135            ibdx=0
136            do ib=1,nbnds*nspinor,nspinor 
137             ibdx=ibdx+1
138             Cprj_ibz(iat,indx_kibz+spad+ib)%cp(1,ilmn)=re_p(ibdx)
139             Cprj_ibz(iat,indx_kibz+spad+ib)%cp(2,ilmn)=im_p(ibdx)
140            end do
141         end do !ispinor
142 
143       end do !ilmn
144     end do !iat
145 
146    end do !ik
147  end do !is
148 
149  ABI_DEALLOCATE(re_p)
150  ABI_DEALLOCATE(im_p)
151  ABI_DEALLOCATE(vect)
152  ABI_DEALLOCATE(umat)
153 
154  call pawcprj_free(Cprj_ks) 
155  ABI_DATATYPE_DEALLOCATE(Cprj_ks)
156 
157  DBG_EXIT("COLL")
158 
159 end subroutine update_cprj