The Elk Code
tm3todm.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 F. Bultmark, F. Cricchio, L. Nordstrom and J. K. Dewhurst.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: tm3todm
8 ! !INTERFACE:
9 subroutine tm3todm(l,k,p,r,ld,wkpr,dm)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! l : angular momentum quantum number (in,integer)
12 ! k : k-index of tensor moment (in,integer)
13 ! p : p-index of tensor moment (in,integer)
14 ! r : r-index of tensor moment (in,integer)
15 ! ld : leading dimension (in,integer)
16 ! wkpr : 3-index tensor moment components (in,real(-ld:ld))
17 ! dm : complex Hermitian density matrix (out,complex(ld,2,ld,2))
18 ! !DESCRIPTION:
19 ! The 3-index coupled tensor moment matrices are given by
20 ! $$ \Gamma_t^{kpr}=
21 ! \sqrt{2r+1}\sum_{x=-k}^k\sum_{y=-p}^p
22 ! \begin{pmatrix} k & r & p \\ -x & t & -y \end{pmatrix}
23 ! \Gamma_{xy}^{kp}, $$
24 ! where the irreducible representations are labeled by $k\in\{0,\ldots,2l\}$,
25 ! $p\in\{0,1\}$, $r\in\{|k-p|,\ldots,k+p\}$ and $\Gamma_{xy}^{kp}$ are the
26 ! uncoupled tensor moments (note that the phase $(-1)^{x+y}$ in the original
27 ! formula has been removed because of the Wigner $3j$ condition $x+y=t$). The
28 ! coupled tensor moment matrices are real and orthonormal in the sense
29 ! $$ \tr\big(\Gamma_t^{kpr}\Gamma_{t'}^{k'p'r'}\big)=
30 ! \delta_{kk'}\delta_{pp'}\delta_{rr'}\delta_{tt'}. $$
31 ! It can also be shown that the matrices are complete, thus any general
32 ! complex matrix $D$ of dimension $2(2l+1)$ can be expanded as
33 ! $$ D=\sum_{k=0}^{2l}\sum_{p=0}^1\sum_{r=|k-p|}^{k+p}\sum_{t=-r}^r
34 ! z_t^{kpr}\Gamma_t^{kpr} $$
35 ! where $z_t^{kpr}$ are complex numbers. Likewise, any real matrix can be
36 ! expanded in real tensor moments $w_t^{kpr}$. Using the the symmetry
37 ! properties of the Wigner $3j$-symbols, one can show that the transpose
38 ! $$ \big(\Gamma_t^{kpr}\big)^t=(-1)^{k+p+r+t}\,\Gamma_{-t}^{kpr} $$
39 ! and thus both the symmetric and antisymmetric parts of $\Gamma_t^{kpr}$
40 ! transform under rotation within the same irreducible representation.
41 ! Consequently, any complex Hermitian matrix $D$ can be written as
42 ! $$ D=\sum_{k,p,r,t} w_t^{kpr}\big[(\Gamma_t^{kpr})_{\rm S}
43 ! +i(\Gamma_t^{kpr})_{\rm A}\big], $$
44 ! where the subscripts S and A refer to the symmetric and antisymmetric parts
45 ! of the matrix, respectively. This routine generates the Hermitian density
46 ! matrix $D$ as described above from the real tensor moments $w_t^{kpr}$. For
47 ! a detailed derivation see {\it Phys. Rev. B} {\bf 80}, 035121 (2009),
48 ! {\it J. Phys.: Condens. Matter} {\bf 7}, 9947 (1995) and G. van der Laan in
49 ! {\it Spin-Orbit-Influenced Spectroscopies of Magnetic Solids. Lecture Notes
50 ! in Physics}, {\bf 466} (1996). See also the routines {\tt tm2todm} and
51 ! {\tt tm3rtoz}.
52 !
53 ! !REVISION HISTORY:
54 ! Created 2007 (Francesco Cricchio and Lars Nordstrom)
55 ! Changed normalisation, made the moments real and the matrix Hermitian,
56 ! January 2022 (JKD)
57 !EOP
58 !BOC
59 implicit none
60 ! arguments
61 integer, intent(in) :: l,k,p,r,ld
62 real(8), intent(in) :: wkpr(-ld:ld)
63 complex(8), intent(out) :: dm(ld,2,ld,2)
64 ! local variables
65 integer x,y,t
66 real(8) t0,t1
67 ! automatic arrays
68 real(8) wkp(-ld:ld,-1:1),dmr(ld,2,ld,2)
69 ! external functions
70 real(8), external :: wigner3j
71 if (l < 0) then
72  write(*,*)
73  write(*,'("Error(tm3todm): l < 0 : ",I0)') l
74  write(*,*)
75  stop
76 end if
77 if (k < 0) then
78  write(*,*)
79  write(*,'("Error(tm3todm): k < 0 : ",I0)') k
80  write(*,*)
81  stop
82 end if
83 if (k > 2*l) then
84  write(*,*)
85  write(*,'("Error(tm3todm): k > 2*l :",2(X,I0))') k,2*l
86  write(*,*)
87  stop
88 end if
89 if ((p < 0).or.(p > 1)) then
90  write(*,*)
91  write(*,'("Error(tm3todm): p should be 0 or 1 : ",I0)') p
92  write(*,*)
93  stop
94 end if
95 if (r < abs(k-p)) then
96  write(*,*)
97  write(*,'("Error(tm3todm): r < |k-p| :",2(X,I0))') r,abs(k-p)
98  write(*,*)
99  stop
100 end if
101 if (r > (k+p)) then
102  write(*,*)
103  write(*,'("Error(tm3todm): r > k+p :",2(X,I0))') r,k+p
104  write(*,*)
105  stop
106 end if
107 ! compute 2-index tensor moment from 3-index tensor moment
108 wkp(:,:)=0.d0
109 t0=sqrt(dble(2*r+1))
110 do t=-r,r
111  t1=wkpr(t)
112  if (abs(t1) < 1.d-8) cycle
113  t1=t0*t1
114  do x=-k,k
115  do y=-p,p
116  wkp(x,y)=wkp(x,y)+t1*wigner3j(k,r,p,-x,t,-y)
117  end do
118  end do
119 end do
120 ! compute the real matrix from the 2-index tensor moment
121 call tm2todm(l,k,p,ld,wkp,dmr)
122 ! convert to complex Hermitian matrix
123 call dmrtoz(ld*2,dmr,dm)
124 
125 contains
126 
127 pure subroutine dmrtoz(n,dmr,dm)
128 implicit none
129 ! arguments
130 integer, intent(in) :: n
131 real(8), intent(in) :: dmr(n,n)
132 complex(8), intent(out) :: dm(n,n)
133 ! local variables
134 integer i,j
135 real(8) a,b
136 do j=1,n
137  do i=1,j-1
138 ! symmetric part
139  a=0.5d0*(dmr(i,j)+dmr(j,i))
140 ! antisymmetric part
141  b=0.5d0*(dmr(i,j)-dmr(j,i))
142  dm(i,j)=cmplx(a,b,8)
143  dm(j,i)=cmplx(a,-b,8)
144  end do
145  dm(j,j)=dmr(j,j)
146 end do
147 end subroutine
148 
149 end subroutine
150 !EOC
151 
subroutine tm2todm(l, k, p, ld, wkp, dm)
Definition: tm2todm.f90:10
pure subroutine dmrtoz(n, dmr, dm)
Definition: tm3todm.f90:128
subroutine tm3todm(l, k, p, r, ld, wkpr, dm)
Definition: tm3todm.f90:10