The Elk Code
 
Loading...
Searching...
No Matches
ylmrot.f90
Go to the documentation of this file.
1
2! Copyright (C) 2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross
3! This file is distributed under the terms of the GNU Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: ylmrot
8! !INTERFACE:
9subroutine ylmrot(p,ang,lmax,ld,d)
10! !INPUT/OUTPUT PARAMETERS:
11! p : if p=-1 then the rotation matrix is improper (in,integer)
12! ang : Euler angles; alpha, beta, gamma (in,real(3))
13! lmax : maximum angular momentum (in,integer)
14! ld : leading dimension (in,integer)
15! d : complex spherical harmonic rotation matrix (out,complex(ld,*))
16! !DESCRIPTION:
17! Returns the rotation matrix in the basis of complex spherical harmonics
18! given the three Euler angles, $(\alpha,\beta,\gamma)$, and the parity, $p$,
19! of the rotation. The matrix is given by the formula
20! $$ D^l_{m_1m_2}(\alpha,\beta,\gamma)=d^l_{m_1m_2}(\beta)
21! e^{-i(m_1\alpha+m_2\gamma)}, $$
22! where $d$ is the rotation matrix about the $y$-axis. For improper rotations,
23! i.e. those which are a combination of a rotation and inversion, $D$ is
24! modified with $D^l_{m_1m_2}\rightarrow(-1)^l D^l_{m_1m_2}$. See the routines
25! {\tt roteuler} and {\tt ylmroty}.
26!
27! !REVISION HISTORY:
28! Created December 2008 (JKD)
29!EOP
30!BOC
31implicit none
32! arguments
33integer, intent(in) :: p
34real(8), intent(in) :: ang(3)
35integer, intent(in) :: lmax,ld
36complex(8), intent(out) :: d(ld,*)
37! local variables
38integer l,m1,m2,lm1,lm2,n
39real(8) t1
40! automatic arrays
41real(8) dy(ld,ld)
42! generate the rotation matrix about the y-axis
43call ylmroty(ang(2),lmax,ld,dy)
44! apply inversion if required
45if (p == -1) then
46 do l=1,lmax,2
47 lm1=l**2+1
48 lm2=lm1+2*l
49 dy(lm1:lm2,lm1:lm2)=-dy(lm1:lm2,lm1:lm2)
50 end do
51end if
52! rotation by alpha and gamma
53do l=0,lmax
54 n=l*(l+1)+1
55 do m1=-l,l
56 lm1=n+m1
57 do m2=-l,l
58 lm2=n+m2
59 t1=-dble(m1)*ang(1)-dble(m2)*ang(3)
60 d(lm1,lm2)=dy(lm1,lm2)*cmplx(cos(t1),sin(t1),8)
61 end do
62 end do
63end do
64end subroutine
65!EOC
66
subroutine ylmrot(p, ang, lmax, ld, d)
Definition ylmrot.f90:10
subroutine ylmroty(beta, lmax, ld, dy)
Definition ylmroty.f90:10