TABLE OF CONTENTS
ABINIT/wvl_nl_gradient [ Functions ]
NAME
wvl_nl_gradient
FUNCTION
Compute the non local part of the wavefunction gradient.
COPYRIGHT
Copyright (C) 2005-2018 ABINIT group (DC) 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
OUTPUT
SIDE EFFECTS
NOTES
PARENTS
forstr,vtorho
CHILDREN
nonlocal_forces,wrtout,xmpi_sum
SOURCE
31 #if defined HAVE_CONFIG_H 32 #include "config.h" 33 #endif 34 35 #include "abi_common.h" 36 37 subroutine wvl_nl_gradient(grnl, mpi_enreg, natom, rprimd, wvl, xcart) 38 39 use defs_basis 40 use defs_abitypes 41 use defs_wvltypes 42 use m_profiling_abi 43 use m_errors 44 use m_xmpi 45 46 #if defined HAVE_BIGDFT 47 use BigDFT_API, only: nonlocal_forces 48 #endif 49 50 !This section has been created automatically by the script Abilint (TD). 51 !Do not modify the following lines by hand. 52 #undef ABI_FUNC 53 #define ABI_FUNC 'wvl_nl_gradient' 54 use interfaces_14_hidewrite 55 !End of the abilint section 56 57 implicit none 58 59 !Arguments ------------------------------------ 60 !scalars 61 integer, intent(in) :: natom 62 type(MPI_type),intent(in) :: mpi_enreg 63 type(wvl_data),intent(inout) :: wvl 64 !arrays 65 real(dp),intent(in) :: xcart(3,natom),rprimd(3,3) 66 real(dp),intent(inout) :: grnl(3,natom) 67 68 !Local variables------------------------------- 69 #if defined HAVE_BIGDFT 70 !scalars 71 integer :: ia,ierr,igeo,me,nproc,spaceComm 72 character(len=500) :: message 73 !arrays 74 real(dp),allocatable :: gxyz(:,:) 75 real(dp)::strtens(6,4) 76 #endif 77 78 ! ************************************************************************* 79 80 #if defined HAVE_BIGDFT 81 82 !Compute forces 83 write(message, '(a,a)' ) ' wvl_nl_gradient(): compute non-local part to gradient.' 84 call wrtout(std_out,message,'COLL') 85 86 !Nullify output arrays. 87 grnl(:, :) = zero 88 strtens(:,:)=zero 89 90 ABI_ALLOCATE(gxyz,(3, natom)) 91 gxyz(:,:) = zero 92 93 !Add the nonlocal part of the forces to grtn (BigDFT routine) 94 spaceComm=mpi_enreg%comm_wvl 95 me=xmpi_comm_rank(spaceComm) 96 nproc=xmpi_comm_size(spaceComm) 97 call nonlocal_forces(wvl%descr%Glr, & 98 & wvl%descr%h(1), wvl%descr%h(2), wvl%descr%h(3), wvl%descr%atoms, & 99 & xcart, wvl%wfs%ks%orbs, wvl%projectors%nlpsp, wvl%wfs%ks%Lzd%Glr%wfd, & 100 & wvl%wfs%ks%psi, gxyz, .true.,strtens(1,2), & 101 & proj_G=wvl%projectors%G,paw=wvl%descr%paw) 102 103 if (nproc > 1) then 104 call xmpi_sum(gxyz, spaceComm, ierr) 105 end if 106 107 !Forces should be in reduced coordinates. 108 do ia = 1, natom, 1 109 do igeo = 1, 3, 1 110 grnl(igeo, ia) = - rprimd(1, igeo) * gxyz(1, ia) - & 111 & rprimd(2, igeo) * gxyz(2, ia) - & 112 & rprimd(3, igeo) * gxyz(3, ia) 113 end do 114 end do 115 ABI_DEALLOCATE(gxyz) 116 117 #else 118 BIGDFT_NOTENABLED_ERROR() 119 if (.false.) write(std_out,*) natom,mpi_enreg%nproc,wvl%wfs%ks,xcart(1,1),rprimd(1,1),grnl(1,1) 120 #endif 121 122 end subroutine wvl_nl_gradient