TABLE OF CONTENTS


ABINIT/chi0_bbp_mask [ Functions ]

[ Top ] [ Functions ]

NAME

  chi0_bbp_mask

FUNCTION

INPUTS

OUTPUT

PARENTS

      cchi0,cchi0q0

CHILDREN

SOURCE

373 #if defined HAVE_CONFIG_H
374 #include "config.h"
375 #endif
376 
377 #include "abi_common.h"
378 
379 subroutine chi0_bbp_mask(Ep,use_tr,QP_BSt,mband,ikmq_ibz,ik_ibz,spin,spin_fact,bbp_mask)
380 
381  use defs_basis
382  use defs_datatypes
383  use m_profiling_abi
384  use m_errors
385 
386  use m_gwdefs,        only : GW_TOL_DOCC, em1params_t
387 
388 !This section has been created automatically by the script Abilint (TD).
389 !Do not modify the following lines by hand.
390 #undef ABI_FUNC
391 #define ABI_FUNC 'chi0_bbp_mask'
392 !End of the abilint section
393 
394  implicit none
395 
396 !Arguments ------------------------------------
397 !scalars
398  integer,intent(in) :: spin,ik_ibz,ikmq_ibz,mband
399  real(dp),intent(in) :: spin_fact
400  logical,intent(in) :: use_tr
401  type(em1params_t),intent(in) :: Ep
402  type(ebands_t),target,intent(in) :: QP_BSt
403 !arrays
404  logical,intent(out) :: bbp_mask(mband,mband)
405 
406 !Local variables-------------------------------
407 !scalars
408  integer :: ib1,ib2
409  real(dp) :: deltaeGW_b1kmq_b2k,deltaf_b1kmq_b2k,e_b1_kmq,f_b1_kmq 
410  character(len=500) :: msg
411 !arrays
412  real(dp),pointer :: qp_energy(:,:,:),qp_occ(:,:,:)
413 
414 !************************************************************************
415 
416  qp_energy => QP_BSt%eig(:,:,:)
417  qp_occ    => QP_BSt%occ(:,:,:)
418 
419  bbp_mask=.FALSE.
420  !use_tr = (Ep%awtr==1)
421 
422  SELECT CASE (Ep%gwcomp) 
423 
424  CASE (0)
425  
426  do ib1=1,Ep%nbnds ! Loop over "conduction" states.
427    e_b1_kmq = qp_energy(ib1,ikmq_ibz,spin)
428    f_b1_kmq =    qp_occ(ib1,ikmq_ibz,spin)
429 
430    do ib2=1,Ep%nbnds ! Loop over "valence" states.
431      deltaf_b1kmq_b2k   = spin_fact*(f_b1_kmq-qp_occ(ib2,ik_ibz,spin))
432      deltaeGW_b1kmq_b2k = e_b1_kmq-qp_energy(ib2,ik_ibz,spin)
433 
434      SELECT CASE (Ep%spmeth)
435      CASE (0) ! Standard Adler-Wiser expression.
436 
437        if (ABS(deltaf_b1kmq_b2k) >= GW_TOL_DOCC) then 
438          bbp_mask(ib1,ib2)=.TRUE.
439          if (use_tr .and. ib1<ib2) bbp_mask(ib1,ib2)=.FALSE. ! GAIN a factor ~2 thanks to time-reversal.
440        end if
441 
442      CASE (1,2) ! Spectral method, WARNING time-reversal here is always assumed!
443        if (ABS(deltaf_b1kmq_b2k) >= GW_TOL_DOCC) then 
444          bbp_mask(ib1,ib2)=.TRUE.
445          if (deltaeGW_b1kmq_b2k<zero) bbp_mask(ib1,ib2)=.FALSE.  ! Only positive frequencies are needed for the Hilbert transform.
446          !$if (use_tr .and. ib1<ib2) bbp_mask(ib1,ib2)=.FALSE. ! GAIN a factor ~2 thanks to time-reversal.
447        end if
448 
449      CASE  DEFAULT
450        write(msg,'(a,i0)')" Wrong value for spmeth: ",Ep%spmeth
451        MSG_ERROR(msg)
452      END SELECT
453 !       write(6,*) "bbp_mask(ib1,ib2)",bbp_mask(ib1,ib2)
454 
455    end do !ib2
456  end do !ib1
457 
458  CASE (1) 
459    ! Extrapolar technique
460    ABI_CHECK(Ep%spmeth==0,"Hilbert transform and extrapolar method are not compatible")
461 
462    do ib1=1,Ep%nbnds ! Loop over "conduction" states.
463      e_b1_kmq=qp_energy(ib1,ikmq_ibz,spin)
464      f_b1_kmq=   qp_occ(ib1,ikmq_ibz,spin)
465 
466      do ib2=1,Ep%nbnds ! Loop over "valence" states.
467 
468        deltaf_b1kmq_b2k  =spin_fact*(f_b1_kmq-qp_occ(ib2,ik_ibz,spin))
469        deltaeGW_b1kmq_b2k=e_b1_kmq-qp_energy(ib2,ik_ibz,spin)
470 
471        ! * When the completeness correction is used,
472        !   we need to also consider transitions with vanishing deltaf
473        !
474        ! Rangel: This is to compute chi in metals correctly with the extrapolar method.
475        bbp_mask(ib1,ib2)=.TRUE.
476        !if (qp_occ(ib2,ik_ibz,is) < GW_TOL_DOCC) CYCLE
477        if (qp_occ(ib2,ik_ibz,spin) < GW_TOL_DOCC .and. (ABS(deltaf_b1kmq_b2k) < GW_TOL_DOCC .or. ib1<ib2)) then
478          bbp_mask(ib1,ib2)=.FALSE.
479        end if
480 
481      end do 
482    end do
483 
484   CASE DEFAULT
485     write(msg,'(a,i0)')" Wrong value of gwcomp: ",Ep%gwcomp 
486     MSG_ERROR(msg)
487   END SELECT
488 
489 end subroutine chi0_bbp_mask

ABINIT/make_transitions [ Functions ]

[ Top ] [ Functions ]

NAME

 make_transitions

FUNCTION

  Calculate transition energies entering the espression for the irreducible polarizability.

COPYRIGHT

  Copyright (C) 2007-2017 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

  nsspol=1 for spin unpolarized, 2 for spin polarized calculations
  nbnds=total number of bands
  kmesh<kmesh_t>=datatype gathering info on the k-mesh:
   | %nbz=number of k-points in the full BZ
   | %nibz=number of k-points in the IBZ
   | %tab(nkbz)=table giving for each k-point in the BZ, the corresponding irreducible point in the IBZ array
   | %bz(3,nkbz)=reduced coordinated of k-points
  TOL_DELTA_OCC=tolerance on the difference of the occupation numbers
  gw_energy(nbnds,kmesh%nkibz,nsppol)=quasi-particle energies energies 
  occ(nbnds,kmesh%nkibz,nsppol)=occupation numbers
  chi0alg=integer defining the method used to calculate chi0
   0 ==> calculate chi0 using the Adler-Wiser expression
   1 ==> use spectral method 
  timrev=if 2, time-reversal symmetry is considered; 1 otherwise

OUTPUT

 my_max_rest,my_min_rest=Maximum and minimum resonant (posite) transition energy.
 max_rest,min_rest=Maximun and minimum resonant (posite) transition energy treated by this node.

PARENTS

      cchi0,cchi0q0

CHILDREN

SOURCE

 42 #if defined HAVE_CONFIG_H
 43 #include "config.h"
 44 #endif
 45 
 46 #include "abi_common.h"
 47 
 48 subroutine make_transitions(Wfd,chi0alg,nbnds,nbvw,nsppol,symchi,timrev,TOL_DELTA_OCC,&
 49 & max_rest,min_rest,my_max_rest,my_min_rest,Kmesh,Ltg_q,gw_energy,occ,qpoint,bbp_ks_distrb)
 50 
 51  use defs_basis
 52  use defs_datatypes
 53  use defs_abitypes
 54  use m_profiling_abi
 55  use m_errors
 56  use m_hdr
 57 
 58  use m_bz_mesh, only : kmesh_t, has_BZ_item, littlegroup_t
 59  use m_wfd,     only : wfd_t
 60 
 61 !This section has been created automatically by the script Abilint (TD).
 62 !Do not modify the following lines by hand.
 63 #undef ABI_FUNC
 64 #define ABI_FUNC 'make_transitions'
 65  use interfaces_14_hidewrite
 66 !End of the abilint section
 67 
 68  implicit none
 69 
 70 !Arguments ------------------------------------
 71 !scalars
 72  integer,intent(in) :: chi0alg,nbnds,nbvw,nsppol,symchi,timrev
 73  real(dp),intent(in) :: TOL_DELTA_OCC
 74  real(dp),intent(out) :: max_rest,min_rest
 75  real(dp),intent(out) :: my_max_rest,my_min_rest
 76  type(kmesh_t),intent(in) :: Kmesh
 77  type(littlegroup_t),intent(in) :: Ltg_q
 78  type(wfd_t),intent(in) :: Wfd
 79 !arrays
 80  real(dp),intent(in) :: gw_energy(nbnds,Kmesh%nibz,nsppol)
 81  real(dp),intent(in) :: occ(nbnds,Kmesh%nibz,nsppol),qpoint(3)
 82  integer,intent(in) :: bbp_ks_distrb(Wfd%mband,Wfd%mband,Kmesh%nbz,Wfd%nsppol)
 83 
 84 !Local variables-------------------------------
 85 !scalars
 86  integer :: ib1,ib2,ii,ik_bz,ik_ibz,ikmq_bz,ikmq_ibz,is,nt,ntrans,my_ntrans,iloop
 87  real(dp) :: delta_ene,delta_occ,spin_fact
 88  character(len=500) :: msg
 89 !arrays
 90  integer :: G0(3)
 91  real(dp) :: kmq(3)
 92 
 93 !************************************************************************
 94 
 95  DBG_ENTER("COLL")
 96 
 97  if (chi0alg<0 .or. chi0alg>=2) then 
 98    write(msg,'(a,i3,a)')' chi0alg = ',chi0alg,' not allowed '
 99    MSG_BUG(msg)
100  end if 
101  if (timrev/=1 .and. timrev/=2) then 
102    write(msg,'(a,i3,a)')' timrev = ',timrev,' not allowed'
103    MSG_BUG(msg)
104  end if 
105 
106  ABI_UNUSED(nbvw)
107  !
108  ! In the first loop calculate total number of transitions for this q-point
109  ! as well min and max transition without taking into account distribution of bands. 
110  ! In the second iteration calculate min and Max transition for this processor.
111  !
112  spin_fact=half; if (nsppol==2) spin_fact=one
113  my_max_rest=smallest_real; my_min_rest=greatest_real
114     max_rest=smallest_real;    min_rest=greatest_real
115 
116  do iloop=1,2
117    nt=0
118    do ik_bz=1,Kmesh%nbz
119      ik_ibz=Kmesh%tab(ik_bz)
120      kmq(:)=Kmesh%bz(:,ik_bz)-qpoint(:)
121 
122      if (symchi==1) then  
123        if (Ltg_q%ibzq(ik_bz)/=1) cycle ! This point does not belong to the IBZ defined by the little group
124      end if 
125      !
126      ! Find kp=k-q-G0 and also G0 where kp is in the first BZ
127      if (.not.has_BZ_item(Kmesh,kmq,ikmq_bz,g0)) then ! Stop as the weight 1.0/nkbz is wrong.
128        write(msg,'(4a,2(2a,3f12.6),2a)')ch10,&
129 &        ' make_transitions : ERROR - ',ch10,&
130 &        ' kp  = k-q-G0 not found in the BZ mesh',ch10,&
131 &        ' k   = ',(Kmesh%bz(ii,ik_bz),ii=1,3),ch10,&
132 &        ' k-q = ',(kmq(ii),ii=1,3),ch10,&
133 &        ' weight in cchi0/cchi0q is wrong ' 
134        MSG_ERROR(msg)
135      end if 
136 
137      ikmq_ibz=Kmesh%tab(ikmq_bz)
138      do is=1,nsppol
139        do ib1=1,nbnds           
140          do ib2=1,nbnds
141 
142            if (iloop==2) then
143              if (bbp_ks_distrb(ib1,ib2,ik_bz,is)/=Wfd%my_rank) cycle
144            end if
145 
146            if (timrev==2 .and. ib1<ib2) cycle ! Thanks to time-reversal we gain a factor ~2.
147 
148            delta_occ=spin_fact*(occ(ib1,ikmq_ibz,is)-occ(ib2,ik_ibz,is))
149            delta_ene=gw_energy(ib1,ikmq_ibz,is)-gw_energy(ib2,ik_ibz,is)
150 
151            if (chi0alg==0)  then ! Adler-Wiser expression. Skip only if factor due to occupation number is smaller than TOL_DELTA_OCC
152              if (abs(delta_occ) < abs(TOL_DELTA_OCC)) cycle
153            else if (chi0alg==1) then
154              ! Spectral method with time-reversal, only resonant transitions 
155              ! This has to changed to include spectral method without time-reversal
156              if (delta_ene < -abs(TOL_DELTA_OCC) .or. abs(delta_occ) < abs(TOL_DELTA_OCC)) cycle
157            end if
158            !
159            ! We have a new transition
160            nt=nt+1
161 
162            if (iloop==1) then 
163              max_rest=MAX(max_rest,zero,delta_ene)
164              if (delta_ene>=-tol6) min_rest=MIN(min_rest,delta_ene)
165            end if
166            if (iloop==2) then 
167              my_max_rest=MAX(my_max_rest,zero,delta_ene)
168              if (delta_ene>=-tol6) my_min_rest=MIN(my_min_rest,delta_ene)
169            end if
170 
171          end do 
172        end do 
173      end do
174    end do
175    if (iloop==1) ntrans=nt
176    if (iloop==2) my_ntrans=nt
177  end do !iloop
178 
179  write(msg,'(2a,i9,2a,f8.3,3a,f8.3,a)')ch10,&
180 &  ' Total number of transitions = ',ntrans,ch10,&
181 &  ' min resonant     = ',min_rest*Ha_eV,' [eV] ',ch10,&
182 &  ' Max resonant     = ',max_rest*Ha_eV,' [eV] '
183  call wrtout(std_out,msg,'COLL')
184 
185  if (Wfd%nproc/=1) then 
186    write(msg,'(2a,i9,2a,f8.3,3a,f8.3,a)')ch10,&
187 &    ' Total number of transitions for this processor= ',my_ntrans,ch10,&
188 &    ' min resonant     = ',my_min_rest*Ha_eV,' [eV] ',ch10,&
189 &    ' Max resonant     = ',my_max_rest*Ha_eV,' [eV] '
190    call wrtout(std_out,msg,'PERS')
191  end if
192 
193  DBG_EXIT("COLL")
194 
195 end subroutine make_transitions 

ABINIT/sigma_distribute_bks [ Functions ]

[ Top ] [ Functions ]

NAME

  sigma_distribute_bks

FUNCTION

  Distribute the loop over (b,k,s) used to calculate the self-energy matrix elements 
  taking into account the MPI distribution of the wavefunctions and the use of 
  symmetries to reduce the BZ sum to an appropriate irreducible wedge.

INPUTS

 nsppol
 can_symmetrize(nsppol)=.TRUE if symmetries can be used to reduce the number of k-points to be summed.
 Kmesh<kmesh_t>
 Qmesh<kmesh_t>
 Ltg_kgw<littlegroup_t>
 Wfd(wfd_t),intent(inout) :: 
 mg0(3)
 kptgw(3)
 [bks_mask(Wfd%mband,Kmesh%nbz,nsppol)]
 [got(Wfd%nproc)]=The number of tasks already assigned to the nodes.
 [global]=If true, an MPI global communication is performed such that each node will have the same table. Useful
   if for implementing algorithms in which each node needs to know the global distribution of the tasks, not only
   the task it has to complete. Defaults to .FALSE.

OUTPUT

  my_nbks
  proc_distrb(Wfd%mband,Kmesh%nbz,nsppol)

SIDE EFFECTS

  Wfd%bks_tab

PARENTS

      calc_sigc_me,calc_sigx_me,cohsex_me

CHILDREN

SOURCE

238 #if defined HAVE_CONFIG_H
239 #include "config.h"
240 #endif
241 
242 #include "abi_common.h"
243 
244 subroutine sigma_distribute_bks(Wfd,Kmesh,Ltg_kgw,Qmesh,nsppol,can_symmetrize,kptgw,mg0,my_nbks,proc_distrb,got,bks_mask,global)
245 
246  use defs_basis
247  use defs_datatypes
248  use defs_abitypes
249  use m_xmpi
250  use m_errors
251  use m_profiling_abi
252 
253  use m_bz_mesh,  only : kmesh_t, littlegroup_t, findqg0
254  use m_wfd,      only : wfd_t, wfd_distribute_bands, wfd_update_bkstab
255 
256 !This section has been created automatically by the script Abilint (TD).
257 !Do not modify the following lines by hand.
258 #undef ABI_FUNC
259 #define ABI_FUNC 'sigma_distribute_bks'
260 !End of the abilint section
261 
262  implicit none
263 
264 !Arguments ------------------------------------
265 !scalars
266  integer,intent(in) :: nsppol
267  integer,intent(out) :: my_nbks
268  logical,optional,intent(in) :: global
269  type(kmesh_t),intent(in) :: Kmesh,Qmesh
270  type(littlegroup_t),intent(in) :: Ltg_kgw
271  type(wfd_t),intent(inout) :: Wfd
272 !arrays
273  integer,intent(in) :: mg0(3)
274  integer,optional,intent(inout) :: got(Wfd%nproc)
275  integer,intent(out) :: proc_distrb(Wfd%mband,Kmesh%nbz,nsppol)
276  real(dp),intent(in) :: kptgw(3)
277  logical,intent(in) :: can_symmetrize(Wfd%nsppol)
278  logical,optional,intent(in) :: bks_mask(Wfd%mband,Kmesh%nbz,nsppol)
279 
280 !Local variables-------------------------------
281 !scalars
282  integer :: ierr,ik_bz,ik_ibz,spin,iq_bz,my_nband
283  !character(len=500) :: msg
284 !arrays
285  integer :: g0(3)
286  real(dp) :: kgwmk(3)
287  integer :: get_more(Wfd%nproc),my_band_list(Wfd%mband)
288  !integer :: test(Wfd%mband,Kmesh%nbz,nsppol)
289  logical :: bmask(Wfd%mband)
290 
291 !************************************************************************
292 
293  call wfd_update_bkstab(Wfd)
294 
295  get_more=0; if (PRESENT(got)) get_more=got
296 
297  ! Different distribution of tasks depending whether symmetries can be used or not.
298  proc_distrb= xmpi_undefined_rank
299 
300  do spin=1,Wfd%nsppol
301 
302    if (can_symmetrize(spin)) then 
303      do ik_bz=1,Kmesh%nbz
304        ik_ibz = Kmesh%tab(ik_bz)
305        kgwmk= kptgw-Kmesh%bz(:,ik_bz) ! kptgw must be inside the BZ
306        call findqg0(iq_bz,g0,kgwmk,Qmesh%nbz,Qmesh%bz,mG0) ! Identify q_bz and G0 where q_bz+G0=k_gw-k_bz
307        if (Ltg_kgw%ibzq(iq_bz)==1) then
308          bmask=.FALSE.; bmask(1:Wfd%nband(ik_ibz,spin))=.TRUE.
309          if (PRESENT(bks_mask)) bmask = bks_mask(:,ik_bz,spin)
310          call wfd_distribute_bands(Wfd,ik_ibz,spin,my_nband,my_band_list,got=get_more,bmask=bmask)
311          if (my_nband>0) proc_distrb(my_band_list(1:my_nband),ik_bz,spin)=Wfd%my_rank
312        end if
313      end do
314 
315    else ! No symmetries for this spin. Divide the full BZ among procs.
316      do ik_bz=1,Kmesh%nbz
317        ik_ibz = Kmesh%tab(ik_bz)
318        bmask=.FALSE.; bmask(1:Wfd%nband(ik_ibz,spin))=.TRUE.
319        if (PRESENT(bks_mask)) bmask = bks_mask(:,ik_bz,spin)
320        call wfd_distribute_bands(Wfd,ik_ibz,spin,my_nband,my_band_list,got=get_more,bmask=bmask)
321        if (my_nband>0) proc_distrb(my_band_list(1:my_nband),ik_bz,spin)=Wfd%my_rank
322      end do
323    end if
324  end do ! spin
325 
326  if (PRESENT(global)) then 
327    if (global) then ! Each node will have the same table so that it will know how the tasks are distributed.
328      proc_distrb = proc_distrb + 1
329      where (proc_distrb == xmpi_undefined_rank+1)
330        proc_distrb = 0
331      end where
332      call xmpi_sum(proc_distrb,Wfd%comm,ierr)
333      where (proc_distrb == 0)
334        proc_distrb = xmpi_undefined_rank
335      elsewhere
336        proc_distrb = proc_distrb -1
337      end where
338      !where (proc_distrb /= xmpi_undefined_rank) 
339      !  ltest = (ANY(proc_distrb == (/(ii,ii=0,Wfd%nproc-1)/)))
340      !end where
341      !if (.not.ltest) then
342      !  write(std_out,*)proc_distrb
343      !  MSG_BUG("Bug in the generation of proc_distrb table")
344      !end if
345    end if
346  end if
347 
348  my_nbks = COUNT(proc_distrb==Wfd%my_rank)
349  if (PRESENT(got)) got=get_more
350 
351 end subroutine sigma_distribute_bks