TABLE OF CONTENTS
ABINIT/symredcart [ Functions ]
NAME
symredcart
FUNCTION
Convert a symmetry operation from reduced coordinates (integers) to cartesian coordinates (reals) Can operate in real or reciprocal space
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (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/doc/developers/contributors.txt.
INPUTS
symred(3,3)=symmetry matrice in reduced coordinates (integers) (real or reciprocal space) aprim(3,3)=real or reciprocal space dimensional primitive translations (see below) bprim(3,3)=real or reciprocal space dimensional primitive translations (see below)
OUTPUT
symcart(3,3)=symmetry matrice in cartesian coordinates (reals)
NOTES
When aprim=rprimd and bprim=gprimd, the routine operates in real space (on a real space symmetry) When aprim=gprimd and bprim=rprimd, the routine operates in reciprocal space (on a real space symmetry)
PARENTS
m_matlu,m_phonons,symrhg
CHILDREN
SOURCE
37 #if defined HAVE_CONFIG_H 38 #include "config.h" 39 #endif 40 41 #include "abi_common.h" 42 43 44 subroutine symredcart(aprim,bprim,symcart,symred) 45 46 use defs_basis 47 use m_profiling_abi 48 49 !This section has been created automatically by the script Abilint (TD). 50 !Do not modify the following lines by hand. 51 #undef ABI_FUNC 52 #define ABI_FUNC 'symredcart' 53 !End of the abilint section 54 55 implicit none 56 57 !Arguments ------------------------------------ 58 !arrays 59 integer,intent(in) :: symred(3,3) 60 real(dp),intent(in) :: aprim(3,3),bprim(3,3) 61 real(dp),intent(out) :: symcart(3,3) 62 63 !Local variables------------------------------- 64 !scalars 65 integer :: ii,jj,kk 66 real(dp) :: symtmp 67 !arrays 68 real(dp) :: work(3,3) 69 70 ! ************************************************************************* 71 72 work=zero 73 do kk=1,3 74 do jj=1,3 75 symtmp=dble(symred(jj,kk)) 76 do ii=1,3 77 work(ii,jj)=work(ii,jj)+bprim(ii,kk)*symtmp 78 end do 79 end do 80 end do 81 82 symcart=zero 83 do kk=1,3 84 do jj=1,3 85 symtmp=work(jj,kk) 86 do ii=1,3 87 symcart(ii,jj)=symcart(ii,jj)+aprim(ii,kk)*symtmp 88 end do 89 end do 90 end do 91 92 end subroutine symredcart