The Elk Code
 
Loading...
Searching...
No Matches
dmtotm3.f90
Go to the documentation of this file.
1
2! Copyright (C) 2008 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: dmtotm3
8! !INTERFACE:
9subroutine dmtotm3(l,k,p,r,ld,dm,wkpr)
10! !INPUT/OUTPUT PARAMETERS:
11! l : angular momentum (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! dm : density matrix (in,complex(ld,2,ld,2))
17! wkpr : 3-index spherical tensor moments (out,real(-ld:ld))
18! !DESCRIPTION:
19! Determines the 3-index spherical tensor moments of a density matrix $D$ with
20! $$ w_t^{kpr}=\tr\big(\Gamma_t^{kpr}D\big). $$
21! This exploits the orthonormality of the $\Gamma_t^{kpr}$ matrices. See the
22! routines {\tt tm2todm} and {\tt tm3todm} for more details.
23!
24! !REVISION HISTORY:
25! Created April 2008 (F. Cricchio and L. Nordstrom)
26! Modified, January 2014 (JKD)
27! Changed to real tensor moments, December 2021 (JKD)
28!EOP
29!BOC
30implicit none
31integer, intent(in) :: l,k,p,r,ld
32complex(8), intent(in) :: dm(ld,2,ld,2)
33real(8), intent(out) :: wkpr(-ld:ld)
34! local variables
35integer n,t
36! automatic arrays
37real(8) w(-ld:ld)
38complex(8) gamma(ld,2,ld,2)
39! external functions
40real(8), external :: trzhmm
41if (l < 0) then
42 write(*,*)
43 write(*,'("Error(dmtotm3): l < 0 : ",I8)') l
44 write(*,*)
45 stop
46end if
47if (k < 0) then
48 write(*,*)
49 write(*,'("Error(dmtotm3): k < 0 : ",I8)') k
50 write(*,*)
51 stop
52end if
53if (k > 2*l) then
54 write(*,*)
55 write(*,'("Error(dmtotm3): k > 2*l : ",2I8)') k,2*l
56 write(*,*)
57 stop
58end if
59if ((p < 0).or.(p > 1)) then
60 write(*,*)
61 write(*,'("Error(dmtotm3): p should be 0 or 1 : ",I8)') p
62 write(*,*)
63 stop
64end if
65if (r < abs(k-p)) then
66 write(*,*)
67 write(*,'("Error(dmtotm3): r < |k-p| : ",2I8)') r,abs(k-p)
68 write(*,*)
69 stop
70end if
71if (r > (k+p)) then
72 write(*,*)
73 write(*,'("Error(dmtotm3): r > k+p : ",2I8)') r,k+p
74 write(*,*)
75 stop
76end if
77n=ld*2
78wkpr(:)=0.d0
79do t=-r,r
80 w(:)=0.d0
81 w(t)=1.d0
82 call tm3todm(l,k,p,r,ld,w,gamma)
83 wkpr(t)=trzhmm(n,gamma,dm)
84end do
85end subroutine
86!EOC
87
subroutine dmtotm3(l, k, p, r, ld, dm, wkpr)
Definition dmtotm3.f90:10
subroutine tm3todm(l, k, p, r, ld, wkpr, dm)
Definition tm3todm.f90:10
pure real(8) function trzhmm(n, a, b)
Definition trzhmm.f90:10