The Elk Code
tm3vdl.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2022 Leon Kerber, J. K. Dewhurst and S. Sharma.
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: tm3vdl
8 ! !INTERFACE:
9 subroutine tm3vdl(l,k,p,r,ld,wkpr,wkpr_v)
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 real tensor moment components (in,real(-ld:ld))
17 ! wkpr_v : complex van der Laan components (out,complex(-ld:ld))
18 ! !DESCRIPTION:
19 ! Converts tensor moments from the real form used internally to the original
20 ! complex convention of G. van der Laan and B. T. Thole, Appendix A in
21 ! {\it J. Phys.: Condens. Matter} {\bf 7}, 9947 (1995). The coupled 3-index
22 ! tensors are given by
23 ! $$ \Gamma_t^{kpr}=\sum_{x=-k}^k\sum_{y=-p}^p
24 ! \Gamma_{xy}^{kp}
25 ! \begin{pmatrix} k & r & p \\ -x & t & -y \end{pmatrix}
26 ! (-1)^{k-x+p-y} \underline{n}_{kpr}^{-1}, $$
27 ! where
28 ! $$ \Big(\Gamma_{xy}^{kp}\Big)_{m_1\sigma_1,m_2\sigma_2}=(-1)^{l-m_2}
29 ! \begin{pmatrix} l & k & l \\ -m_2 & x & m_1 \end{pmatrix}
30 ! (-1)^{s-\sigma_2}
31 ! \begin{pmatrix} s & p & s \\ -\sigma_2 & y & \sigma_1 \end{pmatrix}
32 ! n_{lk}^{-1} n_{sp}^{-1}, $$
33 ! $$ n_{lk}=\frac{(2l)!}{\sqrt{(2l-k)!(2l+1+k)!}} $$
34 ! and
35 ! $$ \underline{n}_{kpr}=
36 ! \begin{pmatrix} k & p & r \\ 0 & 0 & 0 \end{pmatrix} $$
37 ! if $g=k+p+r$ is even, and
38 ! $$ \underline{n}_{kpr}=
39 ! i^g\left[\frac{(g-2k)!(g-2p)!(g-2r)!}{(g+1)!}\right]^{\frac{1}{2}}
40 ! \frac{g!!}{(g-2k)!!(g-2p)!!(g-2r)!!} $$
41 ! if $g$ is odd.
42 !
43 ! !REVISION HISTORY:
44 ! Created October 2022 (JKD and Leon Kerber)
45 !EOP
46 !BOC
47 implicit none
48 ! arguments
49 integer, intent(in) :: l,k,p,r,ld
50 real(8), intent(in) :: wkpr(-ld:ld)
51 complex(8), intent(out) :: wkpr_v(-ld:ld)
52 ! local variables
53 integer g,t
54 real(8) t0,t1
55 complex(8), parameter :: zi=(0.d0,1.d0)
56 complex(8) z1
57 ! external functions
58 real(8), external :: wigner3j,factn,factn2,factr
59 ! complex convention normalisation factors
60 g=k+p+r
61 if (mod(g,2) == 0) then
62  t0=1.d0/wigner3j(k,p,r,0,0,0)
63 else
64  t0=sqrt(factr(g+1,g-2*k)/(factn(g-2*p)*factn(g-2*r)))
65  t0=t0*factn2(g-2*k)*factn2(g-2*p)*factn2(g-2*r)/factn2(g)
66 end if
67 t0=t0*sqrt(factn(2+p)*factn(2*l-k)*factn(2*l+k+1))/factn(2*l)
68 if (mod(k+p,2) /= 0) t0=-t0
69 ! remove Hermitian convention normalisation factors
70 t0=t0/sqrt(dble((2*k+1)*(2*p+1)*(2*r+1)))
71 t0=t0/2.d0
72 ! construct complex tensor moments
73 wkpr_v(-r:r)=0.d0
74 do t=-r,r
75  t1=t0*wkpr(t)
76  if (mod(t,2) /= 0) t1=-t1
77  wkpr_v(t)=wkpr_v(t)+t1*(1.d0,1.d0)
78  if (mod(k+p+r+t,2) /= 0) t1=-t1
79  wkpr_v(-t)=wkpr_v(-t)+t1*(1.d0,-1.d0)
80 end do
81 ! multiply by i⁻ᵍ if g is odd
82 if (mod(g,2) /= 0) then
83  z1=zi**(-g)
84  wkpr_v(-r:r)=wkpr_v(-r:r)*z1
85 end if
86 end subroutine
87 !EOC
88 
subroutine tm3vdl(l, k, p, r, ld, wkpr, wkpr_v)
Definition: tm3vdl.f90:10