TABLE OF CONTENTS
ABINIT/sygrad [ Functions ]
NAME
sygrad
FUNCTION
Symmetrize derivatives of energy with respect to coordinates. Unsymmetrized gradients are input as dedt; symmetrized grads are then placed in fred. If nsym=1 simply copy dedt into fred (only symmetry is identity).
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR) 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
natom=number of atoms in cell dedt(3,natom)=unsymmetrized gradients wrt dimensionless tn (hartree) nsym=number of symmetry operators in group symrec(3,3,nsym)=symmetries of group in terms of operations on reciprocal space primitive translations--see comments below indsym(4,nsym,natom)=label given by subroutine symatm, indicating atom label which gets rotated into given atom by given symmetry (first three elements are related primitive translation-- see symatm where this is computed)
OUTPUT
fred(3,3,natom)=symmetrized gradients wrt reduced coordinates (hartree)
NOTES
symmetrization of gradients with respect to reduced coordinates tn is conducted according to the expression $[d(e)/d(t(n,a))]_{symmetrized} = (1/Nsym)*Sum(S)*symrec(n,m,S)* [d(e)/d(t(m,b))]_{unsymmetrized}$ where $t(m,b)= (symrel^{-1})(m,n)*(t(n,a)-tnons(n))$ and tnons is a possible nonsymmorphic translation. The label "b" here refers to the atom which gets rotated into "a" under symmetry "S". symrel is the symmetry matrix in real space, which is the inverse transpose of symrec. symrec is the symmetry matrix in reciprocal space. $sym_{cartesian} = R * symrel * R^{-1} = G * symrec * G^{-1}$ where the columns of R and G are the dimensional primitive translations in real and reciprocal space respectively. Note the use of "symrec" in the symmetrization expression above.
PARENTS
forces
CHILDREN
SOURCE
56 #if defined HAVE_CONFIG_H 57 #include "config.h" 58 #endif 59 60 #include "abi_common.h" 61 62 63 subroutine sygrad(fred,natom,dedt,nsym,symrec,indsym) 64 65 use m_profiling_abi 66 67 use defs_basis 68 69 !This section has been created automatically by the script Abilint (TD). 70 !Do not modify the following lines by hand. 71 #undef ABI_FUNC 72 #define ABI_FUNC 'sygrad' 73 !End of the abilint section 74 75 implicit none 76 77 !Arguments ------------------------------------ 78 !scalars 79 integer,intent(in) :: natom,nsym 80 !arrays 81 integer,intent(in) :: indsym(4,nsym,natom),symrec(3,3,nsym) 82 real(dp),intent(in) :: dedt(3,natom) 83 real(dp),intent(out) :: fred(3,natom) 84 85 !Local variables------------------------------- 86 !scalars 87 integer :: ia,ind,isym,mu 88 real(dp),parameter :: tol=1.0d-30 89 real(dp) :: summ 90 91 ! ************************************************************************* 92 ! 93 if (nsym==1) then 94 ! only symmetry is identity so simply copy 95 fred(:,:)=dedt(:,:) 96 else 97 ! actually conduct symmetrization 98 do ia=1,natom 99 do mu=1,3 100 summ=0._dp 101 do isym=1,nsym 102 ind=indsym(4,isym,ia) 103 summ=summ+dble(symrec(mu,1,isym))*dedt(1,ind)+& 104 & dble(symrec(mu,2,isym))*dedt(2,ind)+& 105 & dble(symrec(mu,3,isym))*dedt(3,ind) 106 end do 107 fred(mu,ia)=summ/dble(nsym) 108 if(abs(fred(mu,ia))<tol)fred(mu,ia)=0.0_dp 109 end do 110 end do 111 end if 112 113 end subroutine sygrad