The Elk Code
writedftu.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 F. Bultmark, F. Cricchio and L. Nordstrom.
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 subroutine writedftu
7 use modmain
8 use moddftu
9 implicit none
10 ! local variables
11 integer ispn,jspn,idu
12 integer is,ia,ias
13 integer k,l,ll,m1,m2,lm1,lm2
14 if (dftu == 0) return
15 ! machine readable density matrix file
16 open(50,file='DMATMT'//trim(filext),form='FORMATTED',action='WRITE')
17 do idu=1,ndftu
18  is=isldu(1,idu)
19  l=isldu(2,idu)
20  ll=l*(l+1)+1
21  do ia=1,natoms(is)
22  ias=idxas(ia,is)
23  write(50,*)
24  write(50,*)
25  write(50,'(3I4," : species, atom, l")') is,ia,l
26  do ispn=1,nspinor
27  do jspn=1,nspinor
28  write(50,*)
29  write(50,'(2I4," : ispn, jspn; m1, m2, dmatmt below")') ispn,jspn
30  do m1=-l,l
31  lm1=ll+m1
32  do m2=-l,l
33  lm2=ll+m2
34  write(50,'(2I6," ",2G18.10)') m1,m2,dmatmt(lm1,ispn,lm2,jspn,ias)
35  end do
36  end do
37  end do
38  end do
39  end do
40 end do
41 close(50)
42 ! machine readable potential matrix file
43 open(50,file='VMATMT'//trim(filext),form='FORMATTED',action='WRITE')
44 do idu=1,ndftu
45  is=isldu(1,idu)
46  l=isldu(2,idu)
47  ll=l*(l+1)+1
48  do ia=1,natoms(is)
49  ias=idxas(ia,is)
50  write(50,*)
51  write(50,*)
52  write(50,'(3I4," : species, atom, l")') is,ia,l
53  do ispn=1,nspinor
54  do jspn=1,nspinor
55  write(50,*)
56  write(50,'(2I4," : ispn, jspn; m1, m2, vmatmt below")') ispn,jspn
57  do m1=-l,l
58  lm1=ll+m1
59  do m2=-l,l
60  lm2=ll+m2
61  write(50,'(2I6," ",2G18.10)') m1,m2,vmatmt(lm1,ispn,lm2,jspn,ias)
62  end do
63  end do
64  end do
65  end do
66  end do
67 end do
68 close(50)
69 ! Slater parameters
70 open(50,file='FDU'//trim(filext),form='FORMATTED',action='WRITE')
71 do idu=1,ndftu
72  is=isldu(1,idu)
73  l=isldu(2,idu)
74  write(50,*)
75  write(50,'(2I4," : species, l")') is,l
76  do k=0,2*l,2
77  write(50,'(G18.10," : F^(",I1,")")') fdu(k,idu),k
78  end do
79  write(50,'(G18.10," : U")') ujdu(1,idu)
80  write(50,'(G18.10," : J")') ujdu(2,idu)
81  if (inpdftu >= 4) write(50,'(G18.10," : screening length λ")') lamdu(idu)
82 end do
83 close(50)
84 end subroutine
85 
character(256) filext
Definition: modmain.f90:1301
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
real(8), dimension(0:2 *lmaxdm, maxdftu) fdu
Definition: moddftu.f90:48
complex(8), dimension(:,:,:,:,:), allocatable dmatmt
Definition: moddftu.f90:16
integer, dimension(2, maxdftu) isldu
Definition: moddftu.f90:40
integer ndftu
Definition: moddftu.f90:38
real(8), dimension(maxdftu) lamdu
Definition: moddftu.f90:52
integer inpdftu
Definition: moddftu.f90:34
integer nspinor
Definition: modmain.f90:267
complex(8), dimension(:,:,:,:,:), allocatable vmatmt
Definition: moddftu.f90:20
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer dftu
Definition: moddftu.f90:32
real(8), dimension(2, maxdftu) ujdu
Definition: moddftu.f90:42
subroutine writedftu
Definition: writedftu.f90:7