TABLE OF CONTENTS


ABINIT/calc_ucrpa [ Functions ]

[ Top ] [ Functions ]

NAME

 calc_ucrpa

FUNCTION

 Calculate the effective interaction in the correlated orbitals

COPYRIGHT

 Copyright (C) 1999-2024 ABINIT group (TApplencourt,BA)
 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

 npwe : number of plane wave for the dielectric constant
 npw : number of plane wave
 nomega  : number of frequencis
 bandinf,bandsup : kohn sham band
 optimisation : string for the optimisation
 Wfd:: MPI communicator
 mesh <kmesh_t>
    %nbz=Number of points in the BZ
    %nibz=Number of points in IBZ
    %kibz(,nibz)=k-point coordinates, irreducible Brillouin zone
    %kbz(3,nbz)=k-point coordinates, full Brillouin zone
    %ktab(nbz)= table giving for each k-point in the BZ (kBZ), the corresponding
    %ktabi(nbz)= for each k-point in the BZ defines whether inversion has to be considered
    %ktabp(nbz)= phase factor associated to tnons
 M1_q_m(bandinf:bandsup,bandinf:bandsup,npw,Qmesh%nibz): Oscillator strengh in Wannier basis
 rhot1_q_m(bandinf:bandsup,bandinf:bandsup,npw,Qmesh%nibz): Oscillator strengh
                                          multiplied by coulomb potential in Wannier basis

OUTPUT

NOTES

SOURCE

 83  subroutine calc_ucrpa(itypatcor,cryst,Kmesh,lpawu,M1_q_m,Qmesh,npwe,&
 84 & npw,nsym,nomega,omegamin,omegamax,bandinf,bandsup,optimisation,ucvol,Wfd,fname,plowan_compute,rhot1,wanbz)
 85 
 86  use defs_basis
 87  use m_abicore
 88  use m_xmpi
 89  use m_errors
 90 
 91  use m_io_tools,      only : open_file
 92  use m_wfd,           only : wfd_t
 93  use m_io_screening,  only : read_screening, em1_ncname
 94  use m_bz_mesh,       only : kmesh_t
 95  use m_crystal,       only : crystal_t
 96  use m_plowannier,    only : operwan_realspace_type,plowannier_type
 97  implicit none
 98 !   _____                   _
 99 !  |_   _|                 | |
100 !    | |  _ __  _ __  _   _| |_
101 !    | | | '_ \| '_ \| | | | __|
102 !   _| |_| | | | |_) | |_| | |_
103 !  |_____|_| |_| .__/ \__,_|\__|
104 !              | |
105 !              |_|
106 
107 !Arguments ------------------------------------
108  integer, intent(in)   :: itypatcor,lpawu,npw,npwe,nsym
109  integer, intent(in)   :: nomega
110  integer, intent(in)   :: bandinf
111  integer, intent(in)   :: bandsup
112  integer, intent(in)   :: plowan_compute
113  character(len=fnlen), intent(in) :: fname
114  character(len=*), intent(in) :: optimisation
115  real(dp), intent(in) :: ucvol,omegamin,omegamax
116 
117  class(wfd_t),intent(inout) :: Wfd
118  type(kmesh_t),intent(in) :: Kmesh,Qmesh
119  type(crystal_t),intent(in) :: Cryst
120  type(operwan_realspace_type),intent(in) :: rhot1(npw,Qmesh%nibz)
121  type(plowannier_type),intent(in) :: wanbz
122  complex(dpc), intent(in) :: M1_q_m(cryst%nattyp(itypatcor),Wfd%nspinor,Wfd%nspinor,2*lpawu+1,2*lpawu+1,npw,Qmesh%nibz)
123 
124 !Local variables ------------------------------
125 !scalars
126  real(dp) :: x
127  real(dp) :: t1,t2
128  real(dp):: tol
129  complex(dpc) :: uu,jj
130 
131  complex :: nC,ualter
132 
133  integer :: iatom1,iatom2,iatom3,iatom4,pos1,pos2,pos3,pos4,il1,il2,il3,il4,ispin,one_orbital
134  integer :: im_paral,iqalloc,ib1,ib2,m1,m2,m3,m4,iqibz,mbband1,mbband2,mbband3,mbband4,spin1,spin2
135  integer :: ierr,ik_bz,ik_ibz,iq_ibz,i,iG1,iG2,iG,iiG,iomega,iomega1,ispinor1,ispinor2,ispinor3,ispinor4
136  integer :: lpawu_read,nkibz,nbband,nkbz,nprocs,nqalloc,nqibz,ms1,ms2,ms3,ms4,mbband,nspinor
137  integer :: isym_kgw,iik,unt, cp_paral
138  complex(dpc) ::ph_mkt,cplx1,cplx2
139 
140  logical :: wannier=.TRUE.
141  logical :: verbose=.FALSE.
142  logical :: bug=.FALSE.
143  logical :: lscr_one
144 
145  character(len=500) :: message
146 
147 !arrays
148  complex(dpc), allocatable :: V_m(:,:,:,:)
149  complex(dpc), allocatable :: U_m(:,:,:,:)
150  complex(dpc),allocatable :: uspin(:,:),jspin(:,:)
151 ! complex(dpc), allocatable :: coeffW_BZ(:,:,:),coeffW_IBZ(:,:,:)
152  complex(dpc), allocatable :: rhot_q_m1m3(:,:,:,:,:,:),rhot_q_m2m4(:,:,:,:,:,:)
153  complex(dpc), allocatable :: rhot_q_m1m3_npwe(:,:,:,:,:,:),rhot_q_m2m4_npwe(:,:,:,:,:,:)
154  complex(dpc),allocatable :: trrho(:,:),sumrhorhoeps(:)
155  complex(gwpc), allocatable :: scr(:,:,:,:)
156 
157  real(dp),allocatable :: k_coord(:,:)!,k_coordIBZ(:,:)
158  real(dp),allocatable :: q_coord(:,:)
159  real(dp),allocatable:: normG(:)
160  complex(dpc),allocatable:: uomega(:),jomega(:)
161  real(dp),allocatable:: omega(:)
162  complex(dpc),allocatable:: eiqr(:)
163 
164  integer,allocatable:: ikmq_bz_t(:,:)
165 
166  logical,allocatable :: bijection(:)
167 !************************************************************************
168 
169  write(message,*) ch10, '==== Calculation of the screened interaction ===='
170  call wrtout(std_out,message,'COLL')
171  call wrtout(ab_out,message,'COLL')
172  write(message,*) ""
173  call wrtout(std_out,message,'COLL')
174  call wrtout(ab_out,message,'COLL')
175  nkbz = Kmesh%nbz
176  nqibz= Qmesh%nibz
177  nspinor=Wfd%nspinor
178 
179  nbband=1+bandsup-bandinf
180 !  _  __            ____
181 ! | |/ /   ___     / __ \
182 ! | ' /   ( _ )   | |  | |
183 ! |  <    / _ \/\ | |  | |
184 ! | . \  | (_>  < | |__| |
185 ! |_|\_\  \___/\/  \___\_\
186 
187  write(message,*) "Read K and Q mesh"
188  call wrtout(std_out,message,'COLL')
189  call wrtout(ab_out,message,'COLL')
190 
191  ABI_MALLOC(k_coord,(nkbz,3))
192  ABI_MALLOC(q_coord,(nqibz,4))
193  ABI_MALLOC(eiqr,(nqibz))
194  eiqr=czero
195 !==Read k and q==!
196 !open(unit=2012,file='ikbz_COORD',form='formatted',status='unknown')
197 !read(2012,*) (ik_bz,k_coord(ik_bz,:),i=1,nkbz)
198 !close(2012)
199 
200  do ik_bz=1,nkbz
201    call kmesh%get_BZ_item(ik_bz,k_coord(ik_bz,:),ik_ibz,isym_kgw,iik,ph_mkt)
202  end do
203 
204 ! open(unit=2012,file='iqbz_COORD',form='formatted',status='unknown')
205 ! read(2012,*)
206  do i=1,nqibz
207 !   read(2012,*) iq_ibz,q_coord(iq_ibz,:)
208    q_coord(i,1)=Qmesh%ibz(1,i)
209    q_coord(i,2)=Qmesh%ibz(2,i)
210    q_coord(i,3)=Qmesh%ibz(3,i)
211    q_coord(i,4)=Qmesh%wt(i)
212 !   if (iq_ibz > nqibz) then
213 !     write(message,*) iq_ibz,nqibz," Error on line",i,"Are you in iBZ ?"
214 !     call wrtout(std_out,message,'COLL')
215 !   end if
216  end do
217 ! close(2012)
218 
219 !==Bijection and array for k-q==!
220  ABI_MALLOC(bijection,(nkbz))
221  ABI_MALLOC(ikmq_bz_t,(nkbz,nqibz))
222  bijection(:)=.FALSE.
223  if (nsym==1) then
224  do ik_bz=1,nkbz
225    do iq_ibz=1,nqibz
226      ikmq_bz_t(ik_bz,iq_ibz)=findkmq(ik_bz,k_coord,q_coord(iq_ibz,:),nkbz)
227      if (ikmq_bz_t(ik_bz,iq_ibz)>nkbz.and.nsym==1) then
228        BUG=.TRUE.
229        write(message,*) "No K-Q for K/Q =",ik_bz,iq_ibz
230        ABI_ERROR(message)
231      end if
232      bijection(ikmq_bz_t(ik_bz,iq_ibz))=.TRUE.
233    end do
234 
235    if (count(bijection).NE.nqibz.and.nsym==1) then
236    BUG=.TRUE.
237    write(message,*) 'No bijection ',ik_bz
238    ABI_ERROR(message)
239    end if
240 
241    bijection(:)=.FALSE.
242  end do
243 
244  if (.NOT.BUG.and.nsym==1) then
245    write(message,*)  "Bijection Ok."
246    call wrtout(std_out,message,'COLL')
247  end if
248 endif
249 !                                           _____
250 !                                          / ____|
251 !   _ __   ___  _ __ _ __ ___   ___       | |  __
252 !  | '_ \ / _ \| '__| '_ ` _ \ / _ \      | | |_ |
253 !  | | | | (_) | |  | | | | | |  __/      | |__| |
254 !  |_| |_|\___/|_|  |_| |_| |_|\___|       \_____|
255 
256  ABI_MALLOC(normG,(npw))
257  if (verbose) then
258    write(message,*) 'Read the potential and G norm'
259    call wrtout(std_out,message,'COLL')
260    if (open_file('normeG',message,newunit=unt,form='formatted',status='unknown') /= 0) then
261      ABI_ERROR(message)
262    end if
263    read(unt,*) (iiG,x,normG(iiG),iG=1,npw)
264    close(unt)
265    !!False norme for G=0 idd G is the inverse of the potential inverse du potentiel (q=0)
266    normG(1)=0
267  end if
268 
269 !========================================================================
270 !------------------------------------------------------------------------
271 
272 !  FIRST PART OF THE ROUTINE: USE M_G^(nn')(q,k) to do formal checks.
273 !                             USE rhot_q_n to do compute bare interaction
274 
275 !------------------------------------------------------------------------
276 !========================================================================
277 
278 !                               _
279 !                              ( )
280 !   _ __ ___        _ __  _ __ |/
281 !  | '_ ` _ \      | '_ \| '_ \
282 !  | | | | | |     | | | | | | |
283 !  |_| |_| |_|     |_| |_|_| |_|
284 
285 !==========================================================
286 !==========================================================
287 ! tol=1E-1
288 ! tolerance for the normalization of wfc: should be around 0.01.
289  tol = 1 ! very large for test.
290  write(message,*) 'Check the norm of M'
291  call wrtout(std_out,message,'COLL')
292  write(message,*) 'Tolerance :',tol
293  call wrtout(std_out,message,'COLL')
294 !  __      __
295 !  \ \    / /
296 !   \ \  / /      _ __
297 !    \ \/ /      | '_ \
298 !     \  /       | | | |
299 !      \/        |_| |_|
300 !==========================================================================
301 !==Compute V_{n,n'}: bare interaction in the KS basis <- rhot_q_n -> V_n
302 !==========================================================================
303 !==========================================================
304 !==Compute V_{n,n'}
305 !==========================================================
306  if(verbose) then
307    write(message,*) ""
308    call wrtout(std_out,message,'COLL')
309    call wrtout(ab_out,message,'COLL')
310    write(message,*)  "==Calcul of the bare kohn-sham interaction V n=="
311    call wrtout(std_out,message,'COLL')
312  endif
313  tol=1E+1
314 
315  if (.NOT.wannier) RETURN
316 
317 !========================================================================
318 !------------------------------------------------------------------------
319 
320 !  SECOND PART OF THE ROUTINE: Read Wannier projections
321 
322 !------------------------------------------------------------------------
323 !========================================================================
324 !
325 ! \ \        / /       (_)
326 !  \ \  /\  / /_ _ _ __  _ __  _  ___ _ __
327 !   \ \/  \/ / _` | '_ \| '_ \| |/ _ \ '__|
328 !    \  /\  / (_| | | | | | | | |  __/ |
329 !     \/  \/ \__,_|_| |_|_| |_|_|\___|_|
330 
331 !==========================================================
332 !==========================================================
333 !== Read Wannier coefficient in forlb.ovlp
334 !==========================================================
335  nkibz=Kmesh%nibz
336 
337 !Read "l"
338  if (plowan_compute<10) then
339    write(message,*) ""
340    call wrtout(std_out,message,'COLL')
341    call wrtout(ab_out,message,'COLL')
342    write(message,*) "Read wannier in iBZ"
343    call wrtout(ab_out,message,'COLL')
344    call wrtout(std_out,message,'COLL')
345    if (open_file('forlb.ovlp',message,newunit=unt,form='formatted',status='unknown') /= 0) then
346      ABI_ERROR(message)
347    end if
348    rewind(unt)
349  read(unt,*) message
350    read(unt,*) message, lpawu_read
351    read(unt,*) message, ib1, ib2
352    close(unt)
353    mbband=2*lpawu_read+1
354  else
355    ib1=wanbz%bandi_wan
356    ib2=wanbz%bandf_wan
357    mbband=2*wanbz%latom_wan(1)%lcalc(1)+1
358    write(message,*)"Read l and bands from wanbz",ib1,ib2,mbband
359    call wrtout(ab_out,message,'COLL')
360    call wrtout(std_out,message,'COLL')
361  endif
362  if(ib1/=bandinf.and.ib2/=bandsup) then
363    write(message,*) "Error with bands",ib1,bandinf,ib2,bandsup
364    ABI_ERROR(message)
365  endif
366 !!Read the bandinf, bandinf redondance information
367 
368 
369 !*******************************************************
370 !if (3==4) then
371 !!  USELESS START

ABINIT/m_calc_ucrpa [ Modules ]

[ Top ] [ Modules ]

NAME

  m_calc_ucrpa

FUNCTION

 Calculate the effective interaction in the correlated orbitals

COPYRIGHT

 Copyright (C) 2006-2024 ABINIT group (BAmadon,ROuterovitch)
 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

20 #if defined HAVE_CONFIG_H
21 #include "config.h"
22 #endif
23 
24 
25 #include "abi_common.h"
26 
27 MODULE m_calc_ucrpa
28 
29 #ifndef HAVE_CRPA_OPTIM
30 #ifdef FC_INTEL
31 #warning "optimization of m_calc_ucrpa is deactivated on intel fortran"
32 !DEC$ NOOPTIMIZE
33 #endif
34 #endif
35 
36   use defs_basis
37  implicit none
38 
39  private
40 
41  public :: calc_ucrpa