The Elk Code
 
Loading...
Searching...
No Matches
rotzflm.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2008 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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: rotzflm
8! !INTERFACE:
9subroutine rotzflm(rot,lmin,lmax,lmmax,n,ld,zflm1,zflm2)
10! !INPUT/OUTPUT PARAMETERS:
11! rot : rotation matrix (in,real(3,3))
12! lmin : minimum angular momentum (in,integer)
13! lmax : maximum angular momentum (in,integer)
14! lmmax : (lmax+1)^2 or larger (in,integer)
15! n : number of functions to rotate (in,integer)
16! ld : leading dimension (in,integer)
17! zflm1 : coefficients of the complex spherical harmonic expansion for each
18! function (in,complex(ld,n))
19! zflm2 : coefficients of rotated functions (out,complex(ld,n))
20! !DESCRIPTION:
21! Rotates a set of complex functions
22! $$ f_i({\bf r})=\sum_{lm}f_{lm}^iY_{lm}(\hat{\bf r}) $$
23! for all $i$, given the coefficients $f_{lm}^i$ and a rotation matrix $R$.
24! This is done by first the computing the Euler angles $(\alpha,\beta,\gamma)$
25! of $R^{-1}$ (see routine {\tt roteuler}) and then applying the spherical
26! harmonic rotation matrix generated by the routine {\tt ylmrot}.
27!
28! !REVISION HISTORY:
29! Created April 2003 (JKD)
30! Modified, December 2008 (JKD)
31!EOP
32!BOC
33implicit none
34! arguments
35real(8), intent(in) :: rot(3,3)
36integer, intent(in) :: lmin,lmax,lmmax,n,ld
37complex(8), intent(in) :: zflm1(ld,n)
38complex(8), intent(out) :: zflm2(ld,n)
39! local variables
40integer l,lm1,lm2,nm,p
41real(8) det,ang(3),angi(3)
42complex(8), parameter :: zzero=(0.d0,0.d0),zone=(1.d0,0.d0)
43! automatic arrays
44complex(8) d(lmmax,lmmax)
45if ((lmin < 0).or.(lmin > lmax).or.(n < 1)) return
46! find the determinant
47det=rot(1,1)*(rot(2,2)*rot(3,3)-rot(3,2)*rot(2,3)) &
48 +rot(2,1)*(rot(3,2)*rot(1,3)-rot(1,2)*rot(3,3)) &
49 +rot(3,1)*(rot(1,2)*rot(2,3)-rot(2,2)*rot(1,3))
50! calculate the Euler angles of the proper rotation
51if (det > 0.d0) then
52 p=1
53 call roteuler(rot,ang)
54else
55 p=-1
56 call roteuler(-rot(:,:),ang)
57end if
58! inverse rotation: the function is to be rotated, not the spherical harmonics
59angi(1)=-ang(3)
60angi(2)=-ang(2)
61angi(3)=-ang(1)
62! determine the rotation matrix for complex spherical harmonics
63call ylmrot(p,angi,lmax,lmmax,d)
64! apply rotation matrix (d and zflm may have different starting indices)
65lm1=lmin**2+1
66lm2=1
67do l=lmin,lmax
68 nm=2*l+1
69 call zgemm('N','N',nm,n,nm,zone,d(lm1,lm1),lmmax,zflm1(lm2,1),ld,zzero, &
70 zflm2(lm2,1),ld)
71 lm1=lm1+nm
72 lm2=lm2+nm
73end do
74end subroutine
75!EOC
76
subroutine roteuler(rot, ang)
Definition roteuler.f90:10
subroutine rotzflm(rot, lmin, lmax, lmmax, n, ld, zflm1, zflm2)
Definition rotzflm.f90:10
subroutine ylmrot(p, ang, lmax, ld, d)
Definition ylmrot.f90:10