The Elk Code
writetm3.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2008 F. Bultmark, F. Cricchio and L. Nordstrom.
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: writetm3
8 ! !INTERFACE:
9 subroutine writetm3
10 ! !USES:
11 use modmain
12 use moddftu
13 use modtest
14 use modvars
15 ! !DESCRIPTION:
16 ! Decompose the density matrix into 3-index tensor moments and write to
17 ! {\tt TENSMOM.OUT}. See {\it Phys. Rev. B} {\bf 80}, 035121 (2009) and
18 ! {\it J. Phys.: Condens. Matter} {\bf 7} 9947 (1995). See also the routines
19 ! {\tt tm2todm} and {\tt tm3todm}.
20 !
21 ! !REVISION HISTORY:
22 ! Created April 2008 (F. Cricchio and L. Nordstrom)
23 ! Updated, December 2021 (JKD)
24 !EOP
25 !BOC
26 implicit none
27 ! local variables
28 integer is,ia,ias,idu
29 integer l,p,k,r,t
30 real(8) ehx,t0,t1
31 ! automatic arrays
32 real(8) wkpr(-lmmaxdm:lmmaxdm),wkpr_t(-lmmaxdm:lmmaxdm)
33 complex(8) wkpr_v(-lmmaxdm:lmmaxdm)
34 complex(8) gamma(lmmaxdm,2,lmmaxdm,2)
35 ! external functions
36 real(8), external :: trzhmm
37 wkpr(:)=0.d0
38 ! open TENSMOM.OUT file
39 open(50,file='TENSMOM'//trim(filext),form='FORMATTED',action='WRITE')
40 write(50,'("Density matrix decomposition in coupled tensor moments")')
41 write(50,'("Tensor moment type :")')
42 select case(tm3type)
43 case(0)
44  write(50,'(" real, corresponding to Hermitian Γ matrices")')
45 case(1)
46  write(50,'(" complex, see Appendix A of G. van der Laan and B. T. Thole,")')
47  write(50,'(" J. Phys.: Condens. Matter 7, 9947 (1995)")')
48 case(2)
49  write(50,'(" real, corresponding to tesseral spherical harmonics")')
50 end select
51 ! loop over DFT+U entries
52 do idu=1,ndftu
53  is=isldu(1,idu)
54  l=isldu(2,idu)
55 ! scale factor for conventional normalisation
56  t0=sqrt(dble((2*l+1)*2))
57  do ia=1,natoms(is)
58  ias=idxas(ia,is)
59  write(50,*)
60  write(50,'("Species : ",I4," (",A,"), atom : ",I4)') is,trim(spsymb(is)),ia
61  write(50,'(" l = ",I1)') l
62  do k=0,2*l
63  do p=0,1
64  do r=abs(k-p),k+p
65 ! decompose density matrix in 3-index tensor moment components
66  call dmtotm3(l,k,p,r,lmmaxdm,dmatmt(:,:,:,:,ias),wkpr)
67 ! determine the contribution to muffin-tin Hartree + exchange energy
68  call tm3todm(l,k,p,r,lmmaxdm,wkpr,gamma)
69  ehx=0.5d0*trzhmm(lmmaxdm*2,gamma,vmatmt(:,:,:,:,ias))
70 ! write out tensor moment components and vector magnitude
71  write(50,*)
72  write(50,'(" k = ",I1,", p = ",I1,", r = ",I1)') k,p,r
73  select case(tm3type)
74  case(0)
75 ! real convention corresponding to Hermitian Γ matrices
76  do t=-r,r
77  write(50,'(" t = ",I2," : ",F14.8)') t,t0*wkpr(t)
78  end do
79  t1=t0*norm2(wkpr(-r:r))
80  case(1)
81 ! complex convention of G. van der Laan
82  call tm3vdl(l,k,p,r,lmmaxdm,wkpr,wkpr_v)
83  do t=-r,r
84  write(50,'(" t = ",I2," : ",2F14.8)') t,wkpr_v(t)
85  end do
86  t1=sqrt(dot_product(wkpr_v(-r:r),wkpr_v(-r:r)))
87  case(2)
88 ! real convention corresponding to tesseral spherical harmonics
89  call tm3vdl(l,k,p,r,lmmaxdm,wkpr,wkpr_v)
90  call tm3tsrl(r,lmmaxdm,wkpr_v,wkpr_t)
91  do t=-r,r
92  write(50,'(" t = ",I2," : ",F14.8)') t,wkpr_t(t)
93  end do
94  t1=norm2(wkpr_t(-r:r))
95  end select
96  write(50,'(" magnitude : ",F14.8)') t1
97  write(50,'(" Hartree + exchange energy : ",F14.8)') ehx
98 ! write to VARIABLES.OUT if required, but only on the last iteration
99  if (wrtvars.and.tlast) then
100  call writevars('wkpr',n1=is,n2=ia,n3=l,n4=k,n5=p,n6=r,nv=2*r+1, &
101  rva=wkpr(-r:r))
102  end if
103  end do
104  end do
105  end do
106 ! end loop over atoms and species
107  end do
108 end do
109 close(50)
110 ! write last entry of tensor moment components to test file if required
111 if (test) then
112  call writetest(820,'Coupled tensor moments',nv=size(wkpr),tol=2.d-2,rva=wkpr)
113 end if
114 end subroutine
115 !EOC
116 
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
character(256) filext
Definition: modmain.f90:1288
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
pure subroutine tm3tsrl(r, ld, wkpr_v, wkpr_t)
Definition: tm3tsrl.f90:10
complex(8), dimension(:,:,:,:,:), allocatable dmatmt
Definition: moddftu.f90:17
integer, dimension(2, maxdftu) isldu
Definition: moddftu.f90:49
logical tlast
Definition: modmain.f90:1042
integer ndftu
Definition: moddftu.f90:47
subroutine writetm3
Definition: writetm3.f90:10
logical test
Definition: modtest.f90:11
complex(8), dimension(:,:,:,:,:), allocatable vmatmt
Definition: moddftu.f90:21
subroutine dmtotm3(l, k, p, r, ld, dm, wkpr)
Definition: dmtotm3.f90:10
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
logical wrtvars
Definition: modvars.f90:9
subroutine tm3vdl(l, k, p, r, ld, wkpr, wkpr_v)
Definition: tm3vdl.f90:10
subroutine tm3todm(l, k, p, r, ld, wkpr, dm)
Definition: tm3todm.f90:10
character(64), dimension(maxspecies) spsymb
Definition: modmain.f90:78
subroutine writevars(vname, n1, n2, n3, n4, n5, n6, nv, iv, iva, rv, rva, zv, zva, sv, sva)
Definition: modvars.f90:16
integer tm3type
Definition: moddftu.f90:101