TABLE OF CONTENTS


ABINIT/gammapositron_fft [ Functions ]

[ Top ] [ Functions ]

NAME

 gammapositron_fft

FUNCTION

 Compute positron electron-positron enhancement factor on a real space FFT grid.

COPYRIGHT

 Copyright (C) 2013-2018 ABINIT group (MT)
 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

  electronpositron <type(electronpositron_type)>=quantities for the electron-positron annihilation
  gprimd(3,3)=dimensional primitive translations in reciprocal space (bohr^-1)
  igamma=type of enhancement factor:
    -1:  gamma=one (test)
     1:  Boronski and Nieminen [1]
     2:  Boronski and Nieminen, RPA limit [1]
     3:  Sterne and Kaiser [2]
     4:  Puska, Seitsonen and Nieminen [3]
  mpi_enreg=information about MPI parallelization
  n3xccc=dimension of the xccc3d array (0 or nfft).
  nfft=number of FFT grid points
  ngfft(18)=contain all needed information about 3D FFT
  rhor_e(nfft)=real space density electronic density (total density)
  rhor_p(nfft)=real space density positronic density
  xccc3d(n3xccc)=3D core electron density for XC core correction

OUTPUT

  gamma(nfft,2)=electron-positron enhancement factor,
                gamma(:,1): using total   electronic density
                gamma(:,2): using valence electronic density

NOTES

  The input densities (rhor_e and rhor_p) should be positive everywhere
  (call mkdenpos routine before entering this one)

PARENTS

      posdoppler,poslifetime

CHILDREN

      gammapositron,xcden

SOURCE

 50 #if defined HAVE_CONFIG_H
 51 #include "config.h"
 52 #endif
 53 
 54 #include "abi_common.h"
 55 
 56 subroutine gammapositron_fft(electronpositron,gamma,gprimd,igamma,mpi_enreg,&
 57 &                            n3xccc,nfft,ngfft,rhor_e,rhor_p,xccc3d)
 58 
 59  use defs_basis
 60  use defs_abitypes
 61  use m_profiling_abi
 62  use m_errors
 63  
 64  use m_electronpositron
 65 
 66 !This section has been created automatically by the script Abilint (TD).
 67 !Do not modify the following lines by hand.
 68 #undef ABI_FUNC
 69 #define ABI_FUNC 'gammapositron_fft'
 70  use interfaces_56_xc, except_this_one => gammapositron_fft
 71 !End of the abilint section
 72 
 73  implicit none
 74 
 75 !Arguments ------------------------------------
 76 !scalars
 77  integer,intent(in) :: igamma,n3xccc,nfft
 78  type(electronpositron_type),pointer :: electronpositron
 79  type(MPI_type),intent(in) :: mpi_enreg
 80 !arrays
 81  integer,intent(in) :: ngfft(18)
 82  real(dp),intent(in) :: gprimd(3,3),rhor_e(nfft),rhor_p(nfft),xccc3d(n3xccc)
 83  real(dp),intent(out) :: gamma(nfft,2)
 84 
 85 !Local variables-------------------------------
 86 !scalars
 87  integer :: cplex,ishift,ngr,ngrad,nspden_ep,usecore
 88 !arrays
 89  real(dp),parameter :: qphon(3)=(/zero,zero,zero/)
 90  real(dp),allocatable :: grhocore2(:),grhoe2(:),rhoc(:,:,:),rhoe(:,:,:)
 91 
 92 ! *************************************************************************
 93 
 94 !Several useful constants
 95  usecore=n3xccc/nfft
 96  cplex=1;ishift=0;ngrad=1;nspden_ep=1
 97  if (electronpositron%ixcpositron==3.or.electronpositron%ixcpositron==31) ngrad=2
 98  ngr=0;if (ngrad==2) ngr=nfft
 99  
100 !Allocate several arrays
101  ABI_ALLOCATE(rhoe,(nfft,nspden_ep,ngrad**2))
102  ABI_ALLOCATE(grhoe2,(ngr))
103  ABI_ALLOCATE(grhocore2,(ngr*usecore))
104 
105 !Store electronic density and its gradients
106  call xcden(cplex,gprimd,ishift,mpi_enreg,nfft,ngfft,ngrad,nspden_ep,&
107 & mpi_enreg%paral_kgb,qphon,rhor_e,rhoe)
108 
109 !Compute squared gradient of the electronic density
110  if (ngrad==2) then
111    grhoe2(:)=rhoe(:,1,2)**2+rhoe(:,1,3)**2+rhoe(:,1,4)**2
112    if (usecore>0) then
113      ABI_ALLOCATE(rhoc,(nfft,1,ngrad**2))
114      call xcden(cplex,gprimd,ishift,mpi_enreg,nfft,ngfft,ngrad,nspden_ep,&
115 &     mpi_enreg%paral_kgb,qphon,xccc3d,rhoc)
116      grhocore2(:)=rhoc(:,1,2)**2+rhoc(:,1,3)**2+rhoc(:,1,4)**2
117      ABI_DEALLOCATE(rhoc)
118    end if
119  end if
120 
121 !Compute enhancement factor on FFT grid
122  call gammapositron(gamma,grhocore2,grhoe2,igamma,ngr,nfft,xccc3d,&
123 & rhoe(:,1,1),rhor_p,usecore)
124 
125 !Release temporary memory
126  ABI_DEALLOCATE(rhoe)
127  ABI_DEALLOCATE(grhoe2)
128  ABI_DEALLOCATE(grhocore2)
129 
130 end subroutine gammapositron_fft