mirror of https://gitlab.com/QEF/q-e.git
79 lines
2.5 KiB
Fortran
79 lines
2.5 KiB
Fortran
!
|
|
! Copyright (C) 2001-2018 Quantum ESPRESSO group
|
|
! This file is distributed under the terms of the
|
|
! GNU General Public License. See the file `License'
|
|
! in the root directory of the present distribution,
|
|
! or http://www.gnu.org/copyleft/gpl.txt .
|
|
!
|
|
!-----------------------------------------------------------------------
|
|
subroutine sgam_lr (at, bg, nsym, s, irt, tau, rtau, nat)
|
|
!-----------------------------------------------------------------------
|
|
!
|
|
! This routine computes the vector rtau which contains for each
|
|
! atom and each rotation the vector S\tau_a - \tau_b, where
|
|
! b is the rotated a atom, given by the array irt. These rtau are
|
|
! non zero only if fractional translations are present.
|
|
!
|
|
USE kinds, ONLY : DP
|
|
implicit none
|
|
!
|
|
! first the dummy variables
|
|
!
|
|
integer, intent(in) :: nsym, s (3, 3, 48), nat, irt (48, nat)
|
|
! nsym: number of symmetries of the point group
|
|
! s: matrices of symmetry operations
|
|
! nat : number of atoms in the unit cell
|
|
! irt(n,m) = transformed of atom m for symmetry n
|
|
real(DP), intent(in) :: at (3, 3), bg (3, 3), tau (3, nat)
|
|
! at: direct lattice vectors
|
|
! bg: reciprocal lattice vectors
|
|
! tau: coordinates of the atoms
|
|
real(DP), intent(out):: rtau (3, 48, nat)
|
|
! rtau: the direct translations
|
|
!
|
|
! here the local variables
|
|
!
|
|
integer :: na, nb, isym, ipol
|
|
! counters on: atoms, symmetry operations, polarization
|
|
real(DP) , allocatable :: xau (:,:)
|
|
real(DP) :: ft (3)
|
|
!
|
|
allocate (xau(3,nat))
|
|
!
|
|
! compute the atomic coordinates in crystal axis, xau
|
|
!
|
|
do na = 1, nat
|
|
do ipol = 1, 3
|
|
xau (ipol, na) = bg (1, ipol) * tau (1, na) + &
|
|
bg (2, ipol) * tau (2, na) + &
|
|
bg (3, ipol) * tau (3, na)
|
|
enddo
|
|
enddo
|
|
!
|
|
! for each symmetry operation, compute the atomic coordinates
|
|
! of the rotated atom, ft, and calculate rtau = Stau'-tau
|
|
!
|
|
rtau(:,:,:) = 0.0_dp
|
|
do isym = 1, nsym
|
|
do na = 1, nat
|
|
nb = irt (isym, na)
|
|
do ipol = 1, 3
|
|
ft (ipol) = s (1, ipol, isym) * xau (1, na) + &
|
|
s (2, ipol, isym) * xau (2, na) + &
|
|
s (3, ipol, isym) * xau (3, na) - xau (ipol, nb)
|
|
enddo
|
|
do ipol = 1, 3
|
|
rtau (ipol, isym, na) = at (ipol, 1) * ft (1) + &
|
|
at (ipol, 2) * ft (2) + &
|
|
at (ipol, 3) * ft (3)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
!
|
|
! deallocate workspace
|
|
!
|
|
deallocate(xau)
|
|
return
|
|
end subroutine sgam_lr
|
|
!
|