TABLE OF CONTENTS
ABINIT/moddiel_csrb [ Functions ]
NAME
moddiel_csrb
FUNCTION
- Compute a model real-space dielectric functions from the input density and dielar. Based on formula (1) of PRB 47, nb 15, p. 9892 (by C, S, R and B) - formulas for kf and qtf are from J. phys.: condens. matter 2 (1990) 7597-7611 - The plasmon frequency/pulsation is from Theory of the inhomogeneous electron gas (ed. Lundqvist and march) P327
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DCA, XG, 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/Infos/contributors .
INPUTS
OUTPUT
rdiemac : real space, diagonal, dielectric permittivity
SIDE EFFECTS
WARNINGS
This is experimental code : input, ouptput, results and any other feature may vary greatly.
NOTES
PARENTS
CHILDREN
fourdp,laplacian
SOURCE
40 #if defined HAVE_CONFIG_H 41 #include "config.h" 42 #endif 43 44 #include "abi_common.h" 45 46 47 subroutine moddiel_csrb(dielar,dtset,gprimd,mpi_enreg,rdiemac,rhor_in) 48 49 use defs_basis 50 use defs_abitypes 51 use m_profiling_abi 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 'moddiel_csrb' 57 use interfaces_53_ffts 58 use interfaces_56_recipspace 59 !End of the abilint section 60 61 implicit none 62 63 !Arguments ------------------------------- 64 !scalars 65 type(MPI_type),intent(in) :: mpi_enreg 66 type(dataset_type),intent(in) :: dtset 67 !arrays 68 real(dp),intent(in) :: dielar(7),gprimd(3,3),rhor_in(dtset%nfft,dtset%nspden) 69 real(dp),intent(out) :: rdiemac(dtset%nfft,dtset%nspden) 70 71 !Local variables ------------------------- 72 !real(dp) :: invqtf2(2) 73 !real(dp) :: kf(2),wp2(2) 74 !scalars 75 integer :: ifft 76 real(dp) :: alpha 77 complex(dpc) :: invqtf2,kf,wp2 78 ! character(len=fnlen) :: filapp 79 !arrays 80 real(dp) :: g2cart(dtset%nfft),qdiemac(2,dtset%nfft,dtset%nspden) 81 real(dp) :: rhog(2,dtset%nfft,dtset%nspden),rhor(dtset%nfft,dtset%nspden) 82 real(dp) :: xred2(size(dtset%xred_orig,1),size(dtset%xred_orig,2)) 83 complex(dpc) :: cqdiemac(dtset%nfft,dtset%nspden) 84 complex(dpc) :: crhog(dtset%nfft,dtset%nspden) 85 86 ! ********************************************************************* 87 88 xred2=dtset%xred_orig(:,:,1) 89 alpha=0.0d0 90 !presently works only with nspden=1 91 rhor=rhor_in 92 call laplacian(gprimd,mpi_enreg,dtset%nfft,dtset%nspden,dtset%ngfft,dtset%paral_kgb,g2cart_out=g2cart) 93 94 95 96 call fourdp(1, rhog(:,:,1), rhor(:,1),-1,mpi_enreg,dtset%nfft,dtset%ngfft,dtset%paral_kgb,0) ! 97 98 crhog(:,:)=cmplx(rhog(1,:,:),rhog(2,:,:)) 99 100 cqdiemac=zero 101 do ifft=1,dtset%nfft 102 ! wp2(:)=((4.0d0*pi*rhog(:,ifft,1))) 103 ! kf(:)=((3.0d0*pi*pi*rhog(:,ifft,1))**2)**(1.0d0/6.0d0) 104 ! invqtf2(:)=(4.0d0*kf(:)/pi)**(-1) 105 ! qdiemac(:,ifft,1)=one & 106 ! & + ((dielar(3)-1)**(-1) & 107 ! & + alpha*g2cart(ifft)*invqtf2(:) & 108 ! & + g2cart(ifft)*g2cart(ifft)*4.0d0*pi*pi/(4.0d0*wp2(:)))**(-1) 109 if(crhog(ifft,1)/=(zero,zero)) then 110 wp2=((4.0d0*pi*crhog(ifft,1))) 111 kf=((3.0d0*pi*pi*crhog(ifft,1)))**(1.0d0/3.0d0) 112 invqtf2=(4.0d0*kf/pi)**(-1) 113 cqdiemac(ifft,1)=one & 114 & + ((dielar(3)-1)**(-1) & 115 & + alpha*g2cart(ifft)*invqtf2 & 116 & + g2cart(ifft)*g2cart(ifft)*4.0d0*pi*pi/(4.0d0*wp2))**(-1) 117 write(212,*) cqdiemac(ifft,1),wp2,kf,invqtf2,crhog(ifft,1) 118 end if 119 end do 120 qdiemac(1,:,1)=real(cqdiemac(:,1),dp) 121 qdiemac(2,:,1)=real(aimag(cqdiemac(:,1)),dp) 122 call fourdp(1, qdiemac(:,:,1), rdiemac(:,1),1,mpi_enreg,dtset%nfft,dtset%ngfft,dtset%paral_kgb,0) 123 124 end subroutine moddiel_csrb