The Elk Code
 
Loading...
Searching...
No Matches
rotdmat.f90
Go to the documentation of this file.
1
2! Copyright (C) 2014 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
6subroutine rotdmat(rspl,rspn,lmax,nspinor,ld,dmat1,dmat2)
7implicit none
8! arguments
9real(8), intent(in) :: rspl(3,3),rspn(3,3)
10integer, intent(in) :: lmax,nspinor,ld
11complex(8), intent(in) :: dmat1(ld,nspinor,ld,nspinor)
12complex(8), intent(inout) :: dmat2(ld,nspinor,ld,nspinor)
13! local variables
14integer lds,ispn,jspn,p
15integer lmmax,l,lm1,lm2,nm
16real(8), parameter :: eps=1.d-8
17real(8) det,ang(3),angi(3)
18real(8) v(3),th
19complex(8), parameter :: zzero=(0.d0,0.d0),zone=(1.d0,0.d0)
20complex(8) su2(2,2),a(2,2),b(2,2)
21! allocatable arrays
22complex(8), allocatable :: dm(:,:,:,:),c(:,:),d(:,:)
23! external functions
24real(8), external :: r3mdet
25lmmax=(lmax+1)**2
26allocate(dm(ld,nspinor,ld,nspinor))
27allocate(c(lmmax,lmmax),d(lmmax,lmmax))
28! find the determinant of the spatial rotation matrix
29det=r3mdet(rspl)
30! calculate the Euler angles of the proper rotation
31if (det > 0.d0) then
32 p=1
33 call roteuler(rspl,ang)
34else
35 p=-1
36 call roteuler(-rspl(:,:),ang)
37end if
38! inverse rotation: the matrix is to be rotated, not the spherical harmonics
39angi(1)=-ang(3)
40angi(2)=-ang(2)
41angi(3)=-ang(1)
42! determine the rotation matrix for complex spherical harmonics
43call ylmrot(p,angi,lmax,lmmax,d)
44! apply (l,m) rotation matrix as U*D*conjg(U')
45lds=ld*nspinor
46do ispn=1,nspinor
47 do jspn=1,nspinor
48 lm1=1
49 do l=0,lmax
50 nm=2*l+1
51 call zgemm('N','N',nm,lmmax,nm,zone,d(lm1,lm1),lmmax, &
52 dmat1(lm1,ispn,1,jspn),lds,zzero,c(lm1,1),lmmax)
53 lm1=lm1+nm
54 end do
55 lm1=1
56 do l=0,lmax
57 nm=2*l+1
58 call zgemm('N','C',lmmax,nm,nm,zone,c(1,lm1),lmmax,d(lm1,lm1),lmmax, &
59 zzero,dm(1,ispn,lm1,jspn),lds)
60 lm1=lm1+nm
61 end do
62 end do
63end do
64! spin rotation if required
65if (nspinor == 2) then
66! convert spin rotation matrix to axis-angle form
67 call rotaxang(eps,rspn,det,v,th)
68! find the SU(2) representation of the rotation matrix
69 call axangsu2(v,th,su2)
70! apply SU(2) symmetry matrix as U*D*U† and add to dmat2
71 do lm1=1,lmmax
72 do lm2=1,lmmax
73 a(:,:)=dm(lm1,:,lm2,:)
74 call z2mm(su2,a,b)
75 call z2mmct(b,su2,a)
76 dmat2(lm1,:,lm2,:)=dmat2(lm1,:,lm2,:)+a(:,:)
77 end do
78 end do
79else
80 dmat2(1:lmmax,1,1:lmmax,1)=dmat2(1:lmmax,1,1:lmmax,1)+dm(1:lmmax,1,1:lmmax,1)
81end if
82deallocate(dm,c,d)
83end subroutine
84
pure subroutine axangsu2(v, th, su2)
Definition axangsu2.f90:10
pure real(8) function r3mdet(a)
Definition r3mdet.f90:10
subroutine rotaxang(eps, rot, det, v, th)
Definition rotaxang.f90:10
subroutine rotdmat(rspl, rspn, lmax, nspinor, ld, dmat1, dmat2)
Definition rotdmat.f90:7
subroutine roteuler(rot, ang)
Definition roteuler.f90:10
subroutine ylmrot(p, ang, lmax, ld, d)
Definition ylmrot.f90:10
pure subroutine z2mm(a, b, c)
Definition z2mm.f90:10
pure subroutine z2mmct(a, b, c)
Definition z2mmct.f90:10