The Elk Code
findsym.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007-2008 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: findsym
8 ! !INTERFACE:
9 subroutine findsym(apl1,apl2,nsym,lspl,lspn,iea)
10 ! !USES:
11 use modmain
12 use moddftu
13 use modtddft
14 ! !INPUT/OUTPUT PARAMETERS:
15 ! apl1 : first set of atomic positions in lattice coordinates
16 ! (in,real(3,maxatoms,maxspecies))
17 ! apl2 : second set of atomic positions in lattice coordinates
18 ! (in,real(3,maxatoms,maxspecies))
19 ! nsym : number of symmetries (out,integer)
20 ! lspl : spatial rotation element in lattice point group for each symmetry
21 ! (out,integer(48))
22 ! lspn : spin rotation element in lattice point group for each symmetry
23 ! (out,integer(48))
24 ! iea : equivalent atom index for each symmetry
25 ! (out,integer(iea(natmmax,nspecies,48))
26 ! !DESCRIPTION:
27 ! Finds the symmetries which rotate one set of atomic positions into another.
28 ! Both sets of positions differ only by a translation vector and have the same
29 ! muffin-tin magnetic fields (stored in the global array {\tt bfcmt}). Any
30 ! symmetry element consists of a spatial rotation of the atomic position
31 ! vectors followed by a global magnetic rotation: $\{\alpha_S|\alpha_R\}$. In
32 ! the case of spin-orbit coupling $\alpha_S=\alpha_R$. The symmetries are
33 ! returned as indices of elements in the Bravais lattice point group. An
34 ! index to equivalent atoms is stored in the array {\tt iea}.
35 !
36 ! !REVISION HISTORY:
37 ! Created April 2007 (JKD)
38 ! Fixed use of proper rotations for spin, February 2008 (L. Nordstrom)
39 !EOP
40 !BOC
41 implicit none
42 ! arguments
43 real(8), intent(in) :: apl1(3,maxatoms,maxspecies)
44 real(8), intent(in) :: apl2(3,maxatoms,maxspecies)
45 integer, intent(out) :: nsym,lspl(48),lspn(48)
46 integer, intent(out) :: iea(natmmax,nspecies,48)
47 ! local variables
48 integer isym,jsym,jsym0,jsym1
49 integer is,ia,ias,ja,jas,its,n
50 real(8) sl(3,3),scp(3,3)
51 real(8) c(3,3),d(3,3),v(3),t1
52 ! automatic arrays
53 integer jea(natmmax,nspecies)
54 real(8) apl3(3,natmmax)
55 complex(8) dmat(lmmaxdm,nspinor,lmmaxdm,nspinor)
56 ! external functions
57 real(8), external :: dznrm2
58 nsym=0
59 ! loop over lattice symmetries (spatial rotations)
60 do isym=1,nsymlat
61 ! make real copy of lattice rotation symmetry
62  sl(1:3,1:3)=dble(symlat(1:3,1:3,isym))
63 ! loop over species
64  do is=1,nspecies
65 ! map apl1 coordinates to [0,1) and store in apl3
66  do ia=1,natoms(is)
67  apl3(1:3,ia)=apl1(1:3,ia,is)
68  call r3frac(epslat,apl3(:,ia))
69  end do
70  do ja=1,natoms(is)
71 ! apply lattice symmetry to atomic positions
72  v(1:3)=sl(1:3,1)*apl2(1,ja,is) &
73  +sl(1:3,2)*apl2(2,ja,is) &
74  +sl(1:3,3)*apl2(3,ja,is)
75 ! map coordinates to [0,1)
76  call r3frac(epslat,v)
77 ! check if atomic positions are invariant
78  do ia=1,natoms(is)
79  t1=abs(apl3(1,ia)-v(1))+abs(apl3(2,ia)-v(2))+abs(apl3(3,ia)-v(3))
80  if (t1 < epslat) then
81 ! equivalent atom index
82  jea(ia,is)=ja
83  goto 10
84  end if
85  end do
86 ! not invariant so try new spatial rotation
87  goto 40
88 10 continue
89  end do
90  end do
91 ! all atomic positions invariant at this point
92  jsym=1
93 ! spin polarised case
94  if (spinpol) then
95 ! check invariance of magnetic fields under global spin rotation
96  if (spinorb) then
97 ! with spin-orbit coupling spin rotation equals spatial rotation
98  jsym0=isym
99  jsym1=isym
100  else
101 ! without spin-orbit coupling spin rotation independent of spatial rotation
102  jsym0=1
103  jsym1=nsymlat
104  end if
105  do jsym=jsym0,jsym1
106 ! proper part of symmetry matrix
107  scp(1:3,1:3)=dble(symlatd(jsym))*symlatc(1:3,1:3,jsym)
108 ! rotate global field and check invariance using proper part of symmetry matrix
109  v(1:3)=scp(1:3,1)*bfieldc0(1) &
110  +scp(1:3,2)*bfieldc0(2) &
111  +scp(1:3,3)*bfieldc0(3)
112  t1=abs(bfieldc0(1)-v(1))+abs(bfieldc0(2)-v(2))+abs(bfieldc0(3)-v(3))
113 ! if not invariant try a different global spin rotation
114  if (t1 > epslat) goto 20
115 ! rotate muffin-tin magnetic fields and check invariance
116  do is=1,nspecies
117  do ia=1,natoms(is)
118 ! equivalent atom
119  ja=jea(ia,is)
120  v(1:3)=scp(1:3,1)*bfcmt0(1,ja,is) &
121  +scp(1:3,2)*bfcmt0(2,ja,is) &
122  +scp(1:3,3)*bfcmt0(3,ja,is)
123  t1=abs(bfcmt0(1,ia,is)-v(1)) &
124  +abs(bfcmt0(2,ia,is)-v(2)) &
125  +abs(bfcmt0(3,ia,is)-v(3))
126 ! if not invariant try a different global spin rotation
127  if (t1 > epslat) goto 20
128  end do
129  end do
130 ! all fields invariant
131  goto 30
132 20 continue
133 ! end loop over global spin rotations
134  end do
135 ! magnetic fields not invariant so try different spatial rotation
136  goto 40
137  end if
138 30 continue
139 ! check invariance of density matrices for fixed tensor moment calculations
140  if (ftmtype /= 0) then
141  n=(lmmaxdm*nspinor)**2
142  do is=1,nspecies
143  do ia=1,natoms(is)
144  ias=idxas(ia,is)
145 ! equivalent atom
146  ja=jea(ia,is)
147  jas=idxas(ja,is)
148 ! rotate the fixed tensor moment density matrix
149  dmat(:,:,:,:)=0.d0
150  call rotdmat(symlatc(:,:,isym),symlatc(:,:,jsym),lmaxdm,nspinor, &
151  lmmaxdm,dmftm(:,:,:,:,jas),dmat)
152 ! check invariance
153  dmat(:,:,:,:)=dmat(:,:,:,:)-dmftm(:,:,:,:,ias)
154  t1=dznrm2(n,dmat,1)/dble(n)
155  if (t1 > epsdmat) goto 40
156  end do
157  end do
158  end if
159 ! check invariance of static spin-dependent vector potential
160  if (tafsp) then
161  call r3mm(symlatc(:,:,isym),afspc,c)
162  call r3mmt(c,scp,d)
163  t1=sum(abs(afspc(1:3,1:3)-d(1:3,1:3)))
164  if (t1 > epslat) goto 40
165  end if
166 ! check invariance of time- and spin-dependent vector potential
167  if (tafspt) then
168  do its=1,ntimes
169  call r3mm(symlatc(:,:,isym),afspt(:,:,its),c)
170  call r3mmt(c,scp,d)
171  t1=sum(abs(afspt(1:3,1:3,its)-d(1:3,1:3)))
172  if (t1 > epslat) goto 40
173  end do
174  end if
175 ! everything invariant so add symmetry to set
176  nsym=nsym+1
177  lspl(nsym)=isym
178  lspn(nsym)=jsym
179  do is=1,nspecies
180  do ia=1,natoms(is)
181  iea(ia,is,nsym)=jea(ia,is)
182  end do
183  end do
184 40 continue
185 ! end loop over spatial rotations
186 end do
187 end subroutine
188 !EOC
189 
pure subroutine r3mmt(a, b, c)
Definition: r3mmt.f90:10
real(8), dimension(3, maxatoms, maxspecies) bfcmt0
Definition: modmain.f90:275
real(8), dimension(3, 3) afspc
Definition: modmain.f90:331
logical spinpol
Definition: modmain.f90:228
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
complex(8), dimension(:,:,:,:,:), allocatable dmftm
Definition: moddftu.f90:80
integer ntimes
Definition: modtddft.f90:42
integer, dimension(3, 3, 48) symlat
Definition: modmain.f90:344
integer ftmtype
Definition: moddftu.f90:70
integer nsymlat
Definition: modmain.f90:342
real(8), dimension(3, 3, 48) symlatc
Definition: modmain.f90:350
pure subroutine r3frac(eps, v)
Definition: r3frac.f90:10
integer, dimension(48) symlatd
Definition: modmain.f90:346
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
subroutine rotdmat(rspl, rspn, lmax, nspinor, ld, dmat1, dmat2)
Definition: rotdmat.f90:7
real(8) epslat
Definition: modmain.f90:24
subroutine findsym(apl1, apl2, nsym, lspl, lspn, iea)
Definition: findsym.f90:10
logical tafspt
Definition: modtddft.f90:60
logical spinorb
Definition: modmain.f90:230
real(8) epsdmat
Definition: moddftu.f90:18
real(8), dimension(3) bfieldc0
Definition: modmain.f90:271
logical tafsp
Definition: modmain.f90:329
integer, parameter lmaxdm
Definition: moddftu.f90:13
pure subroutine r3mm(a, b, c)
Definition: r3mm.f90:10
real(8), dimension(:,:,:), allocatable afspt
Definition: modtddft.f90:62