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,sm,t0,t1
31 ! automatic arrays
32 real(8) wkpr(-lmmaxdm:lmmaxdm)
33 complex(8) zkpr(-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,'("Components are in the spherical basis")')
42 write(50,'(" (see Phys. Rev. B. 80, 035121 (2009))")')
43 ! loop over DFT+U entries
44 do idu=1,ndftu
45  is=isldu(1,idu)
46  l=isldu(2,idu)
47 ! scale factor for conventional normalisation
48  t0=sqrt(dble((2*l+1)*2))
49  do ia=1,natoms(is)
50  ias=idxas(ia,is)
51  write(50,*)
52  write(50,'("Species : ",I4," (",A,"), atom : ",I4)') is,trim(spsymb(is)),ia
53  write(50,'(" l = ",I1)') l
54  do k=0,2*l
55  do p=0,1
56  do r=abs(k-p),k+p
57 ! decompose density matrix in 3-index tensor moment components
58  call dmtotm3(l,k,p,r,lmmaxdm,dmatmt(:,:,:,:,ias),wkpr)
59 ! determine the contribution to muffin-tin Hartree + exchange energy
60  call tm3todm(l,k,p,r,lmmaxdm,wkpr,gamma)
61  ehx=0.5d0*trzhmm(lmmaxdm*2,gamma,vmatmt(:,:,:,:,ias))
62 ! write out tensor moment components and vector magnitude
63  write(50,*)
64  write(50,'(" k = ",I1,", p = ",I1,", r = ",I1)') k,p,r
65  if (tm3vdl) then
66 ! G. van der Laan complex convention
67  call tm3rtoz(l,k,p,r,lmmaxdm,wkpr,zkpr)
68  do t=-r,r
69  write(50,'(" t = ",I2," : ",2F14.8)') t,zkpr(t)
70  end do
71  else
72 ! new real convention
73  sm=0.d0
74  do t=-r,r
75  t1=t0*wkpr(t)
76  write(50,'(" t = ",I2," : ",F14.8)') t,t1
77  sm=sm+t1**2
78  end do
79  sm=sqrt(sm)
80  write(50,'(" magnitude : ",F14.8)') sm
81  end if
82  write(50,'(" Hartree + exchange energy : ",F14.8)') ehx
83 ! write to VARIABLES.OUT if required, but only on the last iteration
84  if (wrtvars.and.tlast) then
85  call writevars('wkpr',n1=is,n2=ia,n3=l,n4=k,n5=p,n6=r,nv=2*r+1, &
86  rva=wkpr(-r:r))
87  end if
88  end do
89  end do
90  end do
91 ! end loop over atoms and species
92  end do
93 end do
94 close(50)
95 ! write last entry of tensor moment components to test file if required
96 if (test) then
97  call writetest(820,'Coupled tensor moments',nv=size(wkpr),tol=2.d-2,rva=wkpr)
98 end if
99 end subroutine
100 !EOC
101 
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
character(256) filext
Definition: modmain.f90:1301
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
complex(8), dimension(:,:,:,:,:), allocatable dmatmt
Definition: moddftu.f90:16
integer, dimension(2, maxdftu) isldu
Definition: moddftu.f90:40
logical tlast
Definition: modmain.f90:1053
integer ndftu
Definition: moddftu.f90:38
subroutine tm3rtoz(l, k, p, r, ld, wkpr, zkpr)
Definition: tm3rtoz.f90:7
subroutine writetm3
Definition: writetm3.f90:10
logical test
Definition: modtest.f90:11
complex(8), dimension(:,:,:,:,:), allocatable vmatmt
Definition: moddftu.f90:20
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
logical tm3vdl
Definition: moddftu.f90:89
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