The Elk Code
rotrfmt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2015 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 rotrfmt(rot,nr,nri,rfmt1,rfmt2)
7 use modmain
8 implicit none
9 ! arguments
10 real(8), intent(in) :: rot(3,3)
11 integer, intent(in) :: nr,nri
12 real(8), intent(in) :: rfmt1(*)
13 real(8), intent(out) :: rfmt2(*)
14 ! local variables
15 integer i
16 ! inner part of muffin-tin
17 call rotrflm(rot,lmaxi,nri,lmmaxi,rfmt1,rfmt2)
18 ! outer part of muffin-tin
19 i=lmmaxi*nri+1
20 call rotrflm(rot,lmaxo,nr-nri,lmmaxo,rfmt1(i),rfmt2(i))
21 
22 contains
23 
24 !BOP
25 ! !ROUTINE: rotrflm
26 ! !INTERFACE:
27 subroutine rotrflm(rot,lmax,n,ld,rflm1,rflm2)
28 ! !INPUT/OUTPUT PARAMETERS:
29 ! rot : rotation matrix (in,real(3,3))
30 ! lmax : maximum angular momentum (in,integer)
31 ! n : number of functions to rotate (in,integer)
32 ! ld : leading dimension (in,integer)
33 ! rflm1 : coefficients of the real spherical harmonic expansion for each
34 ! function (in,real(ld,n))
35 ! rflm2 : coefficients of rotated functions (out,complex(ld,n))
36 ! !DESCRIPTION:
37 ! Rotates a set of real functions
38 ! $$ f_i({\bf r})=\sum_{lm}f_{lm}^iR_{lm}(\hat{\bf r}) $$
39 ! for all $i$, given the coefficients $f_{lm}^i$ and a rotation matrix $R$.
40 ! This is done by first the computing the Euler angles $(\alpha,\beta,\gamma)$
41 ! of $R^{-1}$ (see routine {\tt roteuler}) and then applying the spherical
42 ! harmonic rotation matrix generated by the routine {\tt rlmrot}.
43 !
44 ! !REVISION HISTORY:
45 ! Created December 2008 (JKD)
46 !EOP
47 !BOC
48 implicit none
49 ! arguments
50 real(8), intent(in) :: rot(3,3)
51 integer, intent(in) :: lmax,n,ld
52 real(8), intent(in) :: rflm1(ld,*)
53 real(8), intent(out) :: rflm2(ld,*)
54 ! local variables
55 integer l,lm,nm,p
56 real(8) det,ang(3),angi(3)
57 ! automatic arrays
58 real(8) d(ld,ld)
59 if (n == 0) return
60 ! find the determinant
61 det=rot(1,1)*(rot(2,2)*rot(3,3)-rot(3,2)*rot(2,3)) &
62  +rot(2,1)*(rot(3,2)*rot(1,3)-rot(1,2)*rot(3,3)) &
63  +rot(3,1)*(rot(1,2)*rot(2,3)-rot(2,2)*rot(1,3))
64 ! calculate the Euler angles of the proper rotation
65 if (det > 0.d0) then
66  p=1
67  call roteuler(rot,ang)
68 else
69  p=-1
70  call roteuler(-rot(:,:),ang)
71 end if
72 ! inverse rotation: the function is to be rotated, not the spherical harmonics
73 angi(1)=-ang(3)
74 angi(2)=-ang(2)
75 angi(3)=-ang(1)
76 ! determine the rotation matrix for real spherical harmonics
77 call rlmrot(p,angi,lmax,ld,d)
78 ! apply rotation matrix
79 do l=0,lmax
80  nm=2*l+1
81  lm=l**2+1
82  call dgemm('N','N',nm,n,nm,1.d0,d(lm,lm),ld,rflm1(lm,1),ld,0.d0,rflm2(lm,1), &
83  ld)
84 end do
85 end subroutine
86 !EOC
87 
88 !BOP
89 ! !ROUTINE: rlmrot
90 ! !INTERFACE:
91 subroutine rlmrot(p,ang,lmax,ld,d)
92 ! !INPUT/OUTPUT PARAMETERS:
93 ! p : if p=-1 then the rotation matrix is improper (in,integer)
94 ! ang : Euler angles; alpha, beta, gamma (in,real(3))
95 ! lmax : maximum angular momentum (in,integer)
96 ! ld : leading dimension (in,integer)
97 ! d : real spherical harmonic rotation matrix (out,real(ld,*))
98 ! !DESCRIPTION:
99 ! Returns the rotation matrix in the basis of real spherical harmonics given
100 ! the three Euler angles, $(\alpha,\beta,\gamma)$, and the parity, $p$, of the
101 ! rotation. The matrix is determined using the formula of V. V. Nechaev,
102 ! [{\it J. Struct. Chem.} {\bf 35}, 115 (1994)], suitably modified for our
103 ! definition of the real spherical harmonics ($m_1>0$, $m_2>0$):
104 ! \begin{align*}
105 ! &\Delta^l_{00}=d^l_{00}, \\
106 ! &\Delta^l_{m_10}=\sqrt{2}\,(-1)^{m_1}d^l_{0m_1}\cos(m_1\alpha), \\
107 ! &\Delta^l_{0m_2}=\sqrt{2}\,(-1)^{m_2}d^l_{m_20}\cos(m_2\gamma), \\
108 ! &\Delta^l_{-m_10}=-\sqrt{2}\,d^l_{0m_1}\sin(m_1\alpha), \\
109 ! &\Delta^l_{0-m_2}=\sqrt{2}\,d^l_{m_20}\sin(m_2\gamma), \\
110 ! &\Delta^l_{m_1m_2}=(-1)^{m_1}(-1)^{m_2}\{\cos(m_1\alpha)\cos(m_2\gamma)
111 ! [d_A+d_B]-\sin(m_1\alpha)\sin(m_2\gamma)[d_A-d_B]\}, \\
112 ! &\Delta^l_{m_1-m_2}=(-1)^{m_1}\{\sin(m_1\alpha)\cos(m_2\gamma)
113 ! [d_A-d_B]+\cos(m_1\alpha)\sin(m_2\gamma)[d_A+d_B]\}, \\
114 ! &\Delta^l_{-m_1m_2}=-(-1)^{m_2}\{\sin(m_1\alpha)\cos(m_2\gamma)
115 ! [d_A+d_B]+\cos(m_1\alpha)\sin(m_2\gamma)[d_A-d_B]\}, \\
116 ! &\Delta^l_{-m_1-m_2}=\cos(m_1\alpha)\cos(m_2\gamma)
117 ! [d_A-d_B]-\sin(m_1\alpha)\sin(m_2\gamma)[d_A+d_B],
118 ! \end{align*}
119 ! where $d_A\equiv d^l_{-m_1-m_2}$, $d_B\equiv(-1)^{m_1}d^l_{m_1-m_2}$ and
120 ! $d$ is the rotation matrix about the $y$-axis for complex spherical
121 ! harmonics. See the routines {\tt genrlm}, {\tt roteuler} and {\tt ylmroty}.
122 !
123 ! !REVISION HISTORY:
124 ! Created December 2008 (JKD)
125 !EOP
126 !BOC
127 implicit none
128 ! arguments
129 integer, intent(in) :: p
130 real(8), intent(in) :: ang(3)
131 integer, intent(in) :: lmax,ld
132 real(8), intent(out) :: d(ld,*)
133 ! local variables
134 integer l,m1,m2,lm0,lm1,lm2
135 real(8), parameter :: sqtwo=1.4142135623730950488d0
136 real(8) s1,s2,t1,t2,t3,t4,t5,t6,t7,t8
137 ! automatic arrays
138 integer lmi(-lmax:lmax)
139 real(8) ca(lmax),sa(lmax),cg(lmax),sg(lmax),dy(ld,ld)
140 ! generate the complex spherical harmonic rotation matrix about the y-axis
141 call ylmroty(ang(2),lmax,ld,dy)
142 do m1=1,lmax
143  ca(m1)=cos(m1*ang(1)); sa(m1)=sin(m1*ang(1))
144  cg(m1)=cos(m1*ang(3)); sg(m1)=sin(m1*ang(3))
145 end do
146 lm1=0
147 do l=0,lmax
148  do m1=-l,l
149  lm1=lm1+1
150  lmi(m1)=lm1
151  end do
152  lm0=lmi(0)
153  d(lm0,lm0)=dy(lm0,lm0)
154  do m1=1,l
155  if (mod(m1,2) == 0) then
156  s1=1.d0
157  else
158  s1=-1.d0
159  end if
160  t1=sqtwo*dy(lm0,lmi(m1))
161  t2=sqtwo*dy(lmi(m1),lm0)
162  d(lmi(m1),lm0)=s1*t1*ca(m1)
163  d(lm0,lmi(m1))=s1*t2*cg(m1)
164  d(lmi(-m1),lm0)=-t1*sa(m1)
165  d(lm0,lmi(-m1))=t2*sg(m1)
166  do m2=1,l
167  if (mod(m2,2) == 0) then
168  s2=1.d0
169  else
170  s2=-1.d0
171  end if
172  t1=ca(m1)*cg(m2)
173  t2=sa(m1)*sg(m2)
174  t3=sa(m1)*cg(m2)
175  t4=ca(m1)*sg(m2)
176  t5=dy(lmi(-m1),lmi(-m2))
177  t6=s1*dy(lmi(m1),lmi(-m2))
178  t7=t5+t6
179  t8=t5-t6
180  d(lmi(m1),lmi(m2))=s1*s2*(t1*t7-t2*t8)
181  d(lmi(m1),lmi(-m2))=s1*(t3*t8+t4*t7)
182  d(lmi(-m1),lmi(m2))=-s2*(t3*t7+t4*t8)
183  d(lmi(-m1),lmi(-m2))=t1*t8-t2*t7
184  end do
185  end do
186 end do
187 ! apply inversion if required
188 if (p == -1) then
189  do l=1,lmax,2
190  lm1=l**2+1
191  lm2=lm1+2*l
192  d(lm1:lm2,lm1:lm2)=-d(lm1:lm2,lm1:lm2)
193  end do
194 end if
195 end subroutine
196 !EOC
197 
198 end subroutine
199 
integer lmmaxo
Definition: modmain.f90:203
integer lmaxo
Definition: modmain.f90:201
subroutine ylmroty(beta, lmax, ld, dy)
Definition: ylmroty.f90:10
subroutine rotrflm(rot, lmax, n, ld, rflm1, rflm2)
Definition: rotrfmt.f90:28
subroutine rotrfmt(rot, nr, nri, rfmt1, rfmt2)
Definition: rotrfmt.f90:7
subroutine rlmrot(p, ang, lmax, ld, d)
Definition: rotrfmt.f90:92
integer lmmaxi
Definition: modmain.f90:207
subroutine roteuler(rot, ang)
Definition: roteuler.f90:10
integer lmaxi
Definition: modmain.f90:205