TABLE OF CONTENTS


ABINIT/calc_sig_ppm_comp [ Functions ]

[ Top ] [ Functions ]

NAME

 calc_sig_ppm_comp

FUNCTION

 Calculating contributions to self-energy operator using a plasmon-pole model

COPYRIGHT

 Copyright (C) 1999-2018 ABINIT group (FB, GMR, VO, LR, RWG, RShaltaf)
 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

  nomega=number of frequencies to consider
  npwc= number of G vectors in the plasmon pole
  npwc1= 1 if ppmodel==3, =npwc if ppmodel== 4, 1 for all the other cases
  npwc2= 1 if ppmodel==3, =1    if ppmodel== 4, 1 for all the other cases
  npwx=number of G vectors in rhotwgp
  ppmodel=plasmon pole model
  theta_mu_minus_e0i= $\theta(\mu-\epsilon_{k-q,b1,s}), defines if the state is occupied or not
  zcut=small imaginary part to avoid the divergence. (see related input variable)
  omegame0i(nomega)=frequencies where evaluate \Sigma_c ($\omega$ - $\epsilon_i$
  otq(npwc,npwc2)=plasmon pole parameters for this q-point
  botsq(npwc,npwc1)=plasmon pole parameters for this q-point
  eig(npwc,npwc)=the eigvectors of the symmetrized inverse dielectric matrix for this q point
   (first index for G, second index for bands)
  rhotwgp(npwx)=oscillator matrix elements divided by |q+G| i.e
    $\frac{\langle b1 k-q s | e^{-i(q+G)r | b2 k s \rangle}{|q+G|}$

OUTPUT

  sigcme(nomega) (to be described), only relevant if ppm3 or ppm4

  ket(npwc,nomega):

  In case of ppmodel==1,2 it contains

   ket(G,omega) = Sum_G2       conjg(rhotw(G)) * Omega(G,G2) * rhotw(G2)
                          ---------------------------------------------------
                            2 omegatw(G,G2) (omega-E_i + omegatw(G,G2)(2f-1))

NOTES

 Taken from old routine

PARENTS

      calc_sigc_me

CHILDREN

SOURCE

 55 #if defined HAVE_CONFIG_H
 56 #include "config.h"
 57 #endif
 58 
 59 #include "abi_common.h"
 60 
 61 subroutine calc_sig_ppm_comp(npwc,nomega,rhotwgp,botsq,otq,omegame0i_io,zcut,theta_mu_minus_e0i,ket,ppmodel,npwx,npwc1,npwc2)
 62 
 63  use m_profiling_abi
 64  use defs_basis
 65  use m_errors
 66 
 67 !This section has been created automatically by the script Abilint (TD).
 68 !Do not modify the following lines by hand.
 69 #undef ABI_FUNC
 70 #define ABI_FUNC 'calc_sig_ppm_comp'
 71 !End of the abilint section
 72 
 73  implicit none
 74 
 75 !Arguments ------------------------------------
 76 !scalars
 77  integer,intent(in) :: nomega,npwc,npwc1,npwc2,npwx,ppmodel
 78  real(dp),intent(in) :: omegame0i_io,theta_mu_minus_e0i,zcut
 79 !arrays
 80  complex(gwpc),intent(in) :: botsq(npwc,npwc1),rhotwgp(npwx),otq(npwc,npwc2)
 81  complex(gwpc),intent(inout) :: ket(npwc,nomega)
 82 
 83 !Local variables-------------------------------
 84 !scalars
 85  integer :: ig,igp,io
 86  real(dp) :: den,otw,twofm1_zcut
 87  complex(gwpc) :: num,rhotwgdp_igp
 88  logical :: fully_occupied,totally_empty
 89  character(len=500) :: msg
 90 !arrays
 91  complex(gwpc),allocatable :: ket_comp(:)
 92 
 93 !*************************************************************************
 94 
 95  if (ppmodel/=1.and.ppmodel/=2) then
 96    write(msg,'(a,i0,a)')' The completeness trick cannot be used when ppmodel is ',ppmodel,' It should be set to 1 or 2. '
 97    MSG_ERROR(msg)
 98  end if
 99 
100  ABI_ALLOCATE(ket_comp,(npwc))
101  ket_comp(:)=0.d0
102 
103  fully_occupied=(abs(theta_mu_minus_e0i-1.)<0.001)
104  totally_empty=(abs(theta_mu_minus_e0i)<0.001)
105 
106  if(.not.(totally_empty)) then ! not totally empty
107    twofm1_zcut=zcut
108    do igp=1,npwc
109      rhotwgdp_igp=rhotwgp(igp)
110      do ig=1,npwc
111        otw=DBLE(otq(ig,igp)) ! in principle otw -> otw - ieta
112        num = botsq(ig,igp)*rhotwgdp_igp
113 
114        den = omegame0i_io-otw
115        if (den**2>zcut**2) then
116          ket_comp(ig) = ket_comp(ig) - num/(den*otw)*theta_mu_minus_e0i
117        end if
118      end do !ig
119    end do !igp
120  end if ! not totally empty
121 
122  if(.not.(fully_occupied)) then ! not fully occupied
123    twofm1_zcut=-zcut
124 
125    do igp=1,npwc
126      rhotwgdp_igp=rhotwgp(igp)
127      do ig=1,npwc
128        otw=DBLE(otq(ig,igp)) ! in principle otw -> otw - ieta
129        num = botsq(ig,igp)*rhotwgdp_igp
130 
131        den = omegame0i_io-otw
132        if (den**2>zcut**2) then
133          ket_comp(ig) = ket_comp(ig) - num/(den*otw)*(1.-theta_mu_minus_e0i)
134        end if
135      end do !ig
136    end do !igp
137  end if ! not fully occupied
138 
139  do io=1,nomega
140    ket(:,io)=ket(:,io)+0.5*ket_comp(:)
141  end do
142 
143  ABI_DEALLOCATE(ket_comp)
144 
145 end subroutine calc_sig_ppm_comp