The Elk Code
symrvfmt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine symrvfmt(tspin,tnc,nrmt_,nrmti_,npmt_,ld,rvfmt)
7 use modmain
8 implicit none
9 ! arguments
10 logical, intent(in) :: tspin,tnc
11 integer, intent(in) :: nrmt_(nspecies),nrmti_(nspecies),npmt_(nspecies),ld
12 real(8), intent(inout) :: rvfmt(ld,natmtot,*)
13 ! local variables
14 integer nd,isym,lspl,lspn
15 integer is,ia,ja,ias,jas
16 integer nr,nri,np,i
17 real(8) sc(3,3),t0,t1
18 real(8) x1,x2,x3,y1,y2,y3
19 ! allocatable arrays
20 real(8), allocatable :: rvfmt1(:,:,:),rvfmt2(:,:)
21 ! dimension of the vector field
22 if (tnc) then
23  nd=3
24 else
25  nd=1
26 end if
27 allocate(rvfmt1(ld,natmmax,nd),rvfmt2(ld,nd))
28 t0=1.d0/dble(nsymcrys)
29 do is=1,nspecies
30  nr=nrmt_(is)
31  nri=nrmti_(is)
32  np=npmt_(is)
33 ! make copy of vector field for all atoms of current species
34  do i=1,nd
35  do ia=1,natoms(is)
36  ias=idxas(ia,is)
37  rvfmt1(1:np,ia,i)=rvfmt(1:np,ias,i)
38  end do
39  end do
40  do ia=1,natoms(is)
41 ! only symmetrise first equivalent atom (rotate into others)
42  if (.not.tfeqat(ia,is)) cycle
43  ias=idxas(ia,is)
44  rvfmt(1:np,ias,1:nd)=0.d0
45 ! begin loop over crystal symmetries
46  do isym=1,nsymcrys
47 ! equivalent atom
48  ja=ieqatom(ia,is,isym)
49 ! parallel transport of vector field
50  lspl=lsplsymc(isym)
51  do i=1,nd
52  call rotrfmt(symlatc(:,:,lspl),nr,nri,rvfmt1(:,ja,i),rvfmt2(:,i))
53  end do
54  if (tspin) then
55 ! global spin proper rotation matrix in Cartesian coordinates
56  lspn=lspnsymc(isym)
57  sc(1:3,1:3)=symlatd(lspn)*symlatc(1:3,1:3,lspn)
58  else
59 ! set spin rotation equal to spatial rotation
60  lspn=lspl
61  sc(1:3,1:3)=symlatc(1:3,1:3,lspl)
62  end if
63 ! global spin rotation of vector field
64  if (tnc) then
65 ! non-collinear case
66  do i=1,np
67  x1=rvfmt2(i,1); x2=rvfmt2(i,2); x3=rvfmt2(i,3)
68  y1=sc(1,1)*x1+sc(1,2)*x2+sc(1,3)*x3
69  y2=sc(2,1)*x1+sc(2,2)*x2+sc(2,3)*x3
70  y3=sc(3,1)*x1+sc(3,2)*x2+sc(3,3)*x3
71  rvfmt(i,ias,1)=rvfmt(i,ias,1)+y1
72  rvfmt(i,ias,2)=rvfmt(i,ias,2)+y2
73  rvfmt(i,ias,3)=rvfmt(i,ias,3)+y3
74  end do
75  else
76 ! collinear case
77  t1=sc(3,3)
78  rvfmt(1:np,ias,1)=rvfmt(1:np,ias,1)+t1*rvfmt2(1:np,1)
79  end if
80 ! end loop over crystal symmetries
81  end do
82 ! normalise
83  do i=1,nd
84  rvfmt(1:np,ias,i)=t0*rvfmt(1:np,ias,i)
85  end do
86 ! rotate into equivalent atoms
87  do ja=1,natoms(is)
88  if (eqatoms(ia,ja,is).and.(ia /= ja)) then
89  isym=findloc(ieqatom(ia,is,1:nsymcrys),ja,1)
90  jas=idxas(ja,is)
91 ! parallel transport of vector field (using operation inverse)
92  lspl=isymlat(lsplsymc(isym))
93  do i=1,nd
94  call rotrfmt(symlatc(:,:,lspl),nr,nri,rvfmt(:,ias,i),rvfmt(:,jas,i))
95  end do
96  if (tspin) then
97 ! inverse of global proper rotation matrix in Cartesian coordinates
98  lspn=isymlat(lspnsymc(isym))
99  sc(1:3,1:3)=symlatd(lspn)*symlatc(1:3,1:3,lspn)
100  else
101 ! set spin rotation equal to spatial rotation
102  lspn=lspl
103  sc(1:3,1:3)=symlatc(1:3,1:3,lspl)
104  end if
105 ! global spin rotation of vector field
106  if (tnc) then
107 ! non-collinear case
108  do i=1,np
109  x1=rvfmt(i,jas,1); x2=rvfmt(i,jas,2); x3=rvfmt(i,jas,3)
110  y1=sc(1,1)*x1+sc(1,2)*x2+sc(1,3)*x3
111  y2=sc(2,1)*x1+sc(2,2)*x2+sc(2,3)*x3
112  y3=sc(3,1)*x1+sc(3,2)*x2+sc(3,3)*x3
113  rvfmt(i,jas,1)=y1; rvfmt(i,jas,2)=y2; rvfmt(i,jas,3)=y3
114  end do
115  else
116 ! collinear case
117  t1=sc(3,3)
118  rvfmt(1:np,jas,1)=t1*rvfmt(1:np,jas,1)
119  end if
120  end if
121  end do
122 ! end loop over atoms and species
123  end do
124 end do
125 deallocate(rvfmt1,rvfmt2)
126 end subroutine
127 
integer, dimension(maxsymcrys) lspnsymc
Definition: modmain.f90:366
integer natmmax
Definition: modmain.f90:38
logical, dimension(:,:), allocatable tfeqat
Definition: modmain.f90:372
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
subroutine symrvfmt(tspin, tnc, nrmt_, nrmti_, npmt_, ld, rvfmt)
Definition: symrvfmt.f90:7
integer nsymcrys
Definition: modmain.f90:358
integer, dimension(48) isymlat
Definition: modmain.f90:348
logical, dimension(:,:,:), allocatable eqatoms
Definition: modmain.f90:370
integer, dimension(:,:,:), allocatable ieqatom
Definition: modmain.f90:368
integer, dimension(maxsymcrys) lsplsymc
Definition: modmain.f90:364
real(8), dimension(3, 3, 48) symlatc
Definition: modmain.f90:350
integer, dimension(48) symlatd
Definition: modmain.f90:346
subroutine rotrfmt(rot, nr, nri, rfmt1, rfmt2)
Definition: rotrfmt.f90:7
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36