TABLE OF CONTENTS
ABINIT/fred2fdeloc [ Functions ]
NAME
fred2fdeloc
FUNCTION
calculate delocalized forces from reduced coordinate ones
COPYRIGHT
Copyright (C) 2003-2018 ABINIT group (MVer) 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
btinv(3*(natom-1),3*natom)= inverse transpose of B matrix (see delocint) natom = number of atoms gprimd(3,3)=dimensional translations in reciprocal space (bohr-1)
OUTPUT
deloc_force(3*(natom-1))=delocalized forces from reduced coordinate ones fred(3,natom)=delocalized forces in reduced coordinates
SIDE EFFECTS
NOTES
PARENTS
pred_delocint,xfh_recover_deloc
CHILDREN
dgemm,dgemv,wrtout
SOURCE
37 #if defined HAVE_CONFIG_H 38 #include "config.h" 39 #endif 40 41 #include "abi_common.h" 42 43 44 subroutine fred2fdeloc(btinv,deloc_force,fred,natom,gprimd) 45 46 use defs_basis 47 use m_profiling_abi 48 use m_linalg_interfaces 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 'fred2fdeloc' 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 !arrays 63 real(dp),intent(in) :: btinv(3*(natom-1),3*natom),gprimd(3,3),fred(3,natom) 64 real(dp),intent(out) :: deloc_force(3*(natom-1)) 65 66 !Local variables------------------------------- 67 integer :: ii 68 !arrays 69 real(dp) :: fcart(3,natom) 70 character(len=500) :: message 71 72 ! ****************************************************************** 73 74 !make cartesian forces 75 76 call dgemm('N','N',3,natom,3,one,& 77 & gprimd,3,fred,3,zero,fcart,3) 78 79 !turn cartesian to delocalized forces 80 call dgemv('N',3*(natom-1),3*natom,one,& 81 & btinv,3*(natom-1),fcart,1,zero,deloc_force,1) 82 83 write (message,'(a)') 'fred2fdeloc : deloc_force = ' 84 call wrtout(std_out,message,'COLL') 85 86 do ii = 1, 3*(natom-1) 87 write (message,'(I6,E16.6)') ii, deloc_force(ii) 88 call wrtout(std_out,message,'COLL') 89 end do 90 91 end subroutine fred2fdeloc