TABLE OF CONTENTS
ABINIT/fock_ACE_getghc [ Functions ]
NAME
fock_ACE_getghc
FUNCTION
Compute the matrix elements <G|Vx|psi> of the Fock operator in the ACE context.
COPYRIGHT
Copyright (C) 2013-2018 ABINIT group (CMartins,FJ,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 .
INPUTS
cwavef(2,npw*nspinor*ndat)= planewave coefficients of wavefunctions on which Fock operator is applied. gs_ham <type(gs_hamiltonian_type)>=all data for the Hamiltonian to be applied mpi_enreg= information about MPI parallelization
SIDE EFFECTS
ghc(2,npw*ndat)= matrix elements <G|H|C> or <G|H-lambda.S|C> (if sij_opt>=0 or =-1 in getghc) contains the fock exchange term for cwavef at the end.
NOTES
The current version assumes that : * nspinor = 1 * no "my_nspinor" * no restriction to the value of istwfk_bz (but must be tested in all case) * all the data for the occupied states (cgocc_bz) are the same as those for the current states (cg)
PARENTS
getghc
CHILDREN
dotprod_g
SOURCE
39 #if defined HAVE_CONFIG_H 40 #include "config.h" 41 #endif 42 43 #include "abi_common.h" 44 45 subroutine fock_ACE_getghc(cwavef,ghc,gs_ham,mpi_enreg) 46 47 use defs_basis 48 use defs_abitypes 49 use m_errors 50 use m_xmpi 51 use m_cgtools, only :dotprod_g 52 use m_fock 53 use m_hamiltonian, only : gs_hamiltonian_type 54 use defs_datatypes, only: pseudopotential_type 55 56 !This section has been created automatically by the script Abilint (TD). 57 !Do not modify the following lines by hand. 58 #undef ABI_FUNC 59 #define ABI_FUNC 'fock_ACE_getghc' 60 !End of the abilint section 61 62 implicit none 63 64 !Arguments ------------------------------------ 65 ! Scalars 66 type(MPI_type),intent(in) :: mpi_enreg 67 type(gs_hamiltonian_type),target,intent(inout) :: gs_ham 68 ! Arrays 69 real(dp),intent(inout) :: cwavef(:,:)!,ghc(2,gs_ham%npw_k) 70 real(dp),intent(inout) :: ghc(:,:) 71 72 !Local variables------------------------------- 73 ! Scalars 74 integer :: iband,ikpt,ipw,my_nspinor,nband_k,npw 75 real(dp) :: doti,dotr,eigen 76 type(fock_common_type),pointer :: fockcommon 77 ! Arrays 78 real(dp), allocatable :: ghc1(:,:),xi(:,:) 79 80 81 ! ************************************************************************* 82 83 ABI_CHECK(associated(gs_ham%fockcommon),"fock must be associated!") 84 fockcommon => gs_ham%fockcommon 85 86 ABI_CHECK(gs_ham%nspinor==1,"only allowed for nspinor=1!") 87 ABI_CHECK(gs_ham%npw_k==gs_ham%npw_kp,"only allowed for npw_k=npw_kp (ground state)!") 88 89 ikpt=fockcommon%ikpt 90 npw=gs_ham%npw_k 91 nband_k=fockcommon%nband(ikpt) 92 my_nspinor=max(1,gs_ham%nspinor/mpi_enreg%nproc_spinor) 93 !*Initialization of the array ghc1 94 !*ghc1 will contain the exact exchange contribution to the Hamiltonian 95 ABI_ALLOCATE(ghc1,(2,npw*my_nspinor)) 96 ghc1=zero 97 ABI_ALLOCATE(xi,(2,npw*my_nspinor)) 98 99 do iband=1, nband_k 100 xi(1,:)=gs_ham%fockACE_k%xi(1,:,iband) 101 xi(2,:)=gs_ham%fockACE_k%xi(2,:,iband) 102 103 call dotprod_g(dotr,doti,gs_ham%istwf_k,npw*my_nspinor,2,xi,cwavef,mpi_enreg%me_g0,mpi_enreg%comm_fft) 104 105 ghc1(1,:)=ghc1(1,:)-(dotr*gs_ham%fockACE_k%xi(1,:,iband)-doti*gs_ham%fockACE_k%xi(2,:,iband)) 106 ghc1(2,:)=ghc1(2,:)-(dotr*gs_ham%fockACE_k%xi(2,:,iband)+doti*gs_ham%fockACE_k%xi(1,:,iband)) 107 end do 108 ABI_DEALLOCATE(xi) 109 110 !* If the calculation is parallelized, perform an MPI_allreduce to sum all the contributions in the array ghc 111 ! ghc(:,:)=ghc(:,:)/mpi_enreg%nproc_kpt + ghc1(:,:) 112 ghc(:,:)=ghc(:,:) + ghc1(:,:) 113 114 ! call xmpi_sum(ghc,mpi_enreg%comm_kpt,ier) 115 116 ! ============================================ 117 ! === Calculate the contribution to energy === 118 ! ============================================ 119 !* Only the contribution when cwavef=cgocc_bz are calculated, in order to cancel exactly the self-interaction 120 !* at each convergence step. (consistent definition with the defintion of hartree energy) 121 if (fockcommon%ieigen/=0) then 122 eigen=zero 123 !* Dot product of cwavef and ghc 124 !* inspired from the routine 53_spacepar/meanvalue_g but without the reference to parallelism and filtering 125 if(gs_ham%istwf_k==2) then 126 eigen=half*cwavef(1,1)*ghc1(1,1) 127 else 128 eigen=cwavef(1,1)*ghc1(1,1)+cwavef(2,1)*ghc1(2,1) 129 end if 130 do ipw=2,npw 131 eigen=eigen+cwavef(1,ipw)*ghc1(1,ipw)+cwavef(2,ipw)*ghc1(2,ipw) 132 end do 133 if(gs_ham%istwf_k>=2) eigen=two*eigen 134 ! call xmpi_sum(eigen,mpi_enreg%comm_kpt,ier) 135 fockcommon%eigen_ikpt(fockcommon%ieigen)= eigen 136 fockcommon%ieigen = 0 137 end if 138 139 ! =============================== 140 ! === Deallocate local arrays === 141 ! =============================== 142 143 ABI_DEALLOCATE(ghc1) 144 145 end subroutine fock_ACE_getghc