The Elk Code
genlmirep.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 General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: genlmirep
8 ! !INTERFACE:
9 subroutine genlmirep(elm,ulm)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! elm : eigenvalues of the symmetrised pseudorandom matrix H
14 ! (out,real(lmmaxdb,natmtot))
15 ! ulm : unitary matrix which will transform a density matrix from the (l,m)
16 ! basis to the irreducible basis (in,complex(lmmaxdb,lmmaxdb,natmtot))
17 ! !DESCRIPTION:
18 ! First generates a symmetric, pseudorandom matrix $H$ which is then
19 ! symmetrised by applying all the group symmetries $\{S_i^{\alpha}\}$ at
20 ! atomic site $\alpha$ as
21 ! $$ \tilde{H}=\sum_i S_i^{\alpha} H {S_i^{\alpha}}^{\dag}. $$
22 ! This matrix is then diagonalised. By Schur's second lemma the eigenvalues
23 ! have degeneracies of the same dimension as an irreducible representation
24 ! (IR) of the symmetry group, and can be used to identify the IR. These
25 ! eigenvalues are returned in the matrix {\tt elm}. The conjugate transpose of
26 ! the eigenvector array forms a unitary matrix $U$ which can be applied
27 ! directly to a density matrix $\gamma$ in the $Y_{lm}$ basis as
28 ! $U\gamma U^{\dag}$. This will transform $\gamma$ into into the IR basis. The
29 ! matrix $U$ is returned in the array {\tt ulm}. See the routines
30 ! {\tt bandstr} and {\tt dos}.
31 !
32 ! !REVISION HISTORY:
33 ! Created August 2002 (JKD)
34 !EOP
35 !BOC
36 implicit none
37 ! arguments
38 real(8), intent(out) :: elm(lmmaxdb,natmtot)
39 complex(8), intent(out) :: ulm(lmmaxdb,lmmaxdb,natmtot)
40 ! local variables
41 integer isym,lspl,ias
42 integer i,j,l,lm,n,p,info
43 ! automatic arrays
44 real(8) rwork(3*lmmaxdb)
45 complex(8) dlat(lmmaxdb,lmmaxdb,nsymlat)
46 complex(8) a(lmmaxdb,lmmaxdb),b(lmmaxdb,lmmaxdb)
47 complex(8) h(lmmaxdb,lmmaxdb),work(2*lmmaxdb)
48 ! construct Yₗₘ rotation matrix for each lattice symmetry
49 a(:,:)=0.d0
50 do i=1,lmmaxdb
51  a(i,i)=1.d0
52 end do
53 do isym=1,nsymlat
54  call rotzflm(symlatc(:,:,isym),0,lmaxdb,lmmaxdb,lmmaxdb,lmmaxdb,a, &
55  dlat(:,:,isym))
56 end do
57 ! set up pseudorandom symmetric matrix H
58 h(:,:)=0.d0
59 p=0
60 do l=0,lmaxdb
61  n=2*l
62  lm=l**2+1
63  do i=lm,lm+n
64  do j=i,lm+n
65  p=p+1
66  h(i,j)=p
67  h(j,i)=p
68  end do
69  end do
70 end do
71 ! loop over species and atoms
72 do ias=1,natmtot
73 ! symmetrise H with site symmetries
74  b(:,:)=0.d0
75  do isym=1,nsymsite(ias)
76 ! spatial rotation element in lattice point group
77  lspl=lsplsyms(isym,ias)
78 ! apply lattice symmetry as S H S†
79  call zgemm('N','N',lmmaxdb,lmmaxdb,lmmaxdb,zone,dlat(:,:,lspl),lmmaxdb,h, &
80  lmmaxdb,zzero,a,lmmaxdb)
81  call zgemm('N','C',lmmaxdb,lmmaxdb,lmmaxdb,zone,a,lmmaxdb,dlat(:,:,lspl), &
82  lmmaxdb,zone,b,lmmaxdb)
83  end do
84 ! block diagonalise symmetrised H
85  do l=0,lmaxdb
86  n=2*l+1
87  lm=l**2+1
88  call zheev('V','U',n,b(lm,lm),lmmaxdb,elm(lm,ias),work,2*lmmaxdb,rwork,info)
89  end do
90 ! the unitary matrix U is the conjugate transpose of the eigenvector array since
91 ! this will act on the density matrix as U γ U†
92  do i=1,lmmaxdb
93  do j=1,lmmaxdb
94  ulm(i,j,ias)=conjg(b(j,i))
95  end do
96  end do
97 end do
98 end subroutine
99 !EOC
100 
complex(8), parameter zone
Definition: modmain.f90:1240
subroutine genlmirep(elm, ulm)
Definition: genlmirep.f90:10
integer lmaxdb
Definition: modmain.f90:1081
real(8), dimension(3, 3, 48) symlatc
Definition: modmain.f90:350
complex(8), parameter zzero
Definition: modmain.f90:1240
integer, dimension(:), allocatable nsymsite
Definition: modmain.f90:374
subroutine rotzflm(rot, lmin, lmax, lmmax, n, ld, zflm1, zflm2)
Definition: rotzflm.f90:10
integer, dimension(:,:), allocatable lsplsyms
Definition: modmain.f90:376