File: mode_group.f90

package info (click to toggle)
espresso 6.7-4
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 311,068 kB
  • sloc: f90: 447,429; ansic: 52,566; sh: 40,631; xml: 37,561; tcl: 20,077; lisp: 5,923; makefile: 4,503; python: 4,379; perl: 1,219; cpp: 761; fortran: 618; java: 568; awk: 128
file content (117 lines) | stat: -rw-r--r-- 3,781 bytes parent folder | download | duplicates (6)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
!
! Copyright (C) 2001-2008 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 mode_group &
    (modenum, xq, at, bg, nat, nrot, s, irt, minus_q, rtau, sym)
  !-----------------------------------------------------------------------
  !
  ! This routine selects, among the symmetry matrices of the point group
  ! of a crystal, the symmetry operations which leave a given mode unchanged
  ! For the moment it assumes that the mode modenum displaces the atom
  ! modenum/3 in the direction mod(modenum,3)+1
  !
  USE kinds, ONLY : DP
  USE constants, ONLY : tpi
  implicit none

  integer, intent(in) :: nat, s (3, 3, 48), irt (48, nat), nrot, modenum
  ! nat : the number of atoms of the system
  ! s   : the symmetry matrices
  ! irt : the rotated atom
  ! nrot: number of symmetry operations
  ! modenum: which displacement pattern

  real(DP), intent(in) :: xq (3), rtau (3, 48, nat), bg (3, 3), at (3, 3)
  ! xq  : the q point
  ! rtau: the translations of each atom
  ! bg  : the reciprocal lattice vectors
  ! at  : the direct lattice vectors
  logical, intent(in) :: minus_q
  ! if true Sq=>-q+G  symmetry is used
  logical, intent(inout) :: sym (48)
  ! on  input: .true. if symm. op. has to be tested
  ! on output: .true. if symm. op. does not change mode modenum
  !
  integer :: isym, nas, ipols, na, sna, ipol, jpol
  ! counters
  real(DP) :: arg
  ! auxiliary
  complex(DP), allocatable :: u (:,:)
  ! the original pattern
  complex(DP)              :: fase, sum
  ! the phase of the mode
  ! check for orthogonality
  complex(DP), allocatable :: work_u (:,:), work_ru (:,:)
  ! the working pattern
  ! the rotated working pattern

  allocate(u(3, nat), work_u(3, nat), work_ru (3, nat))

  if (modenum > 3*nat .or. modenum < 1) call errore ('mode_group', &
       'wrong modenum', 1)
  nas = (modenum - 1) / 3 + 1
  ipols = mod (modenum - 1, 3) + 1
  u (:,:) = (0.d0, 0.d0)
  u (ipols, nas) = (1.d0, 0.d0)
  do na = 1, nat
     call trnvecc (u (1, na), at, bg, - 1)
  enddo
  do isym = 1, nrot
     if (sym (isym) ) then
        do na = 1, nat
           do ipol = 1, 3
              work_u (ipol, na) = u (ipol, na)
           enddo
        enddo
        work_ru (:,:) = (0.d0, 0.d0)
        do na = 1, nat
           sna = irt (isym, na)
           arg = 0.d0
           do ipol = 1, 3
              arg = arg + xq (ipol) * rtau (ipol, isym, na)
           enddo
           arg = arg * tpi
           if (isym.eq.nrot.and.minus_q) then
              fase = CMPLX(cos (arg), sin (arg) ,kind=DP)
           else
              fase = CMPLX(cos (arg), - sin (arg) ,kind=DP)
           endif
           do ipol = 1, 3
              do jpol = 1, 3
                 work_ru (ipol, sna) = work_ru (ipol, sna) + s (jpol, ipol, &
                      isym) * work_u (jpol, na) * fase
              enddo
           enddo
        enddo
        !
        !    Transform back the rotated pattern
        !
        do na = 1, nat
           call trnvecc (work_ru (1, na), at, bg, 1)
           call trnvecc (work_u (1, na), at, bg, 1)
        enddo
        !
        !   only if the pattern remain the same up to a phase we keep
        !   the symmetry
        !
        sum = (0.d0, 0.d0)
        do na = 1, nat
           do ipol = 1, 3
              sum = sum + CONJG(work_u (ipol, na) ) * work_ru (ipol, na)
           enddo
        enddo
        sum = abs (sum)
        if (abs (sum - 1.d0) .gt.1.d-7) sym (isym) = .false.
     endif

  enddo
  deallocate ( work_ru, work_u, u)
  return

end subroutine mode_group