1!
2! Copyright (C) 2001-2009 Quantum ESPRESSO group
3! This file is distributed under the terms of the
4! GNU General Public License. See the file `License'
5! in the root directory of the present distribution,
6! or http://www.gnu.org/copyleft/gpl.txt .
7!
8!---------------------------------------------------------------------------------
9SUBROUTINE symmetrize_at( nsym, s, invs, ft, irt, nat, tau, at, bg, alat, omega )
10  !-------------------------------------------------------------------------------
11  !! Forces atomic coordinates to have the symmetry of a given point group.
12  !
13  USE io_global,  ONLY : stdout
14  USE cellmd,     ONLY : at_old, lmovecell
15  USE kinds
16  !
17  IMPLICIT NONE
18  !
19  INTEGER, INTENT(IN) :: nsym
20  !! total number of crystal symmetries
21  INTEGER, INTENT(IN) :: s(3,3,48)
22  !! symmetry matrices, in crystal axis
23  INTEGER, INTENT(IN) :: invs(48)
24  !! index of inverse operation: \(S^{-1}_i=S(\text{invs}(i))\)
25  INTEGER, INTENT(IN) :: nat
26  !! number of atoms
27  INTEGER, INTENT(IN) :: irt(48,nat)
28  !! symmetric atom for each atom and sym.op.
29  REAL(DP), INTENT(IN) :: ft(3,48)
30  !! fractional translations, in crystal axis
31  REAL(DP), INTENT(INOUT) :: tau(3,nat)
32  !! atomic positions
33  REAL(DP), INTENT(INOUT) :: at(3,3)
34  !! lattice vectors of the simulation cell
35  REAL(DP), INTENT(INOUT) :: bg(3,3)
36  !! reciprocal lattice vectors
37  REAL(DP), INTENT(INOUT) :: alat
38  !! lattice parameter
39  REAL(DP), INTENT(INOUT) :: omega
40  !! volume of the simulation cell
41  !
42  ! ... local variables
43  !
44  integer :: na, icar, ipol, jpol, kpol, lpol, irot
45  real(DP) , allocatable :: xau (:,:)
46  ! atomic coordinates in crystal axis
47  real(DP) :: work, bg_old(3,3), sat(3,3), wrk(3,3), ba(3,3)
48  !
49  allocate(xau(3,nat))
50  !
51  !     Compute the coordinates of each atom in the basis of
52  !     the direct lattice vectors
53  !
54  xau = tau
55  tau = 0.d0
56
57  call cryst_to_cart( nat, xau, bg, -1 )
58
59  do irot = 1, nsym
60     do na = 1, nat
61        do kpol = 1, 3
62           work =  s (1, kpol, irot) * xau (1, na) + &
63                   s (2, kpol, irot) * xau (2, na) + &
64                   s (3, kpol, irot) * xau (3, na) - &
65                   ft(kpol,irot)
66           tau (kpol, irt(irot,na)) = tau (kpol, irt(irot,na)) + work &
67                                    - nint(work-xau(kpol,irt(irot,na)))
68        enddo
69     enddo
70  enddo
71  tau (:,:) = tau(:,:)/nsym
72  !
73  !  If the cell is moving then the lattice vectors has to be
74  !  symmetrized as well
75  !
76  if (lmovecell) then
77     CALL recips( at_old(1,1), at_old(1,2), at_old(1,3), &
78                  bg_old(1,1), bg_old(1,2), bg_old(1,3) )
79     do ipol=1,3
80        do jpol=1,3
81           ba(ipol,jpol) = bg_old(1,ipol) * at(1,jpol) + &
82                           bg_old(2,ipol) * at(2,jpol) + &
83                           bg_old(3,ipol) * at(3,jpol)
84        end do
85     end do
86     at = 0.d0
87     !
88     !  at(i) = 1/nsym sum_S  at_old(m) S(l,m) <bg_old(l)|at(k)> invS(i,k)
89     !
90     do irot=1,nsym
91        do icar = 1, 3
92           do lpol = 1, 3
93              sat(icar,lpol) = at_old(icar,1) * s(lpol,1,irot) &
94                             + at_old(icar,2) * s(lpol,2,irot) &
95                             + at_old(icar,3) * s(lpol,3,irot)
96           end do
97        end do
98        do icar = 1, 3
99           do kpol =1, 3
100              wrk(icar,kpol) = sat(icar,1) * ba(1,kpol) &
101                             + sat(icar,2) * ba(2,kpol) &
102                             + sat(icar,3) * ba(3,kpol)
103           end do
104        end do
105        do icar = 1, 3
106           do ipol =1, 3
107              at(icar,ipol) = at(icar,ipol)  &
108                            + wrk(icar,1) * s(ipol,1,invs(irot)) &
109                            + wrk(icar,2) * s(ipol,2,invs(irot)) &
110                            + wrk(icar,3) * s(ipol,3,invs(irot))
111           end do
112        end do
113     end do
114     at(:,:) = at(:,:) / nsym
115     CALL volume( alat, at(1,1), at(1,2), at(1,3), omega )
116     CALL recips( at(1,1), at(1,2), at(1,3), bg(1,1), bg(1,2), bg(1,3) )
117  end if
118  !
119  !   deallocate work space
120  !
121  deallocate (xau)
122  !
123  call cryst_to_cart(nat, tau, at, 1)
124
125  write (stdout,*) " SYMMETRIZED ATOMIC COORDINATES "
126
127  call output_tau(lmovecell, .FALSE.)
128  !
129  return
130  !
131END SUBROUTINE symmetrize_at
132
133